MAST
continuation_solver_base.cpp
Go to the documentation of this file.
1 /*
2  * MAST: Multidisciplinary-design Adaptation and Sensitivity Toolkit
3  * Copyright (C) 2013-2019 Manav Bhatia
4  *
5  * This library is free software; you can redistribute it and/or
6  * modify it under the terms of the GNU Lesser General Public
7  * License as published by the Free Software Foundation; either
8  * version 2.1 of the License, or (at your option) any later version.
9  *
10  * This library is distributed in the hope that it will be useful,
11  * but WITHOUT ANY WARRANTY; without even the implied warranty of
12  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13  * Lesser General Public License for more details.
14  *
15  * You should have received a copy of the GNU Lesser General Public
16  * License along with this library; if not, write to the Free Software
17  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
18  */
19 
20 // MAST includes
22 #include "base/nonlinear_system.h"
23 #include "base/assembly_base.h"
25 #include "base/parameter.h"
26 
27 // libMesh includes
28 #include "libmesh/linear_solver.h"
29 #include "libmesh/dof_map.h"
30 #include "libmesh/petsc_matrix.h"
31 #include "libmesh/petsc_vector.h"
32 
33 // PETSc includes
34 #include <petscmat.h>
35 
36 
38 max_it (10),
39 abs_tol (1.e-8),
40 rel_tol (1.e-8),
41 arc_length (0.),
42 min_step (2.),
43 max_step (20.),
44 step_size_change_exponent (0.5),
45 step_desired_iters (5),
46 schur_factorization (true),
47 _initialized (false),
48 _elem_ops (nullptr),
49 _assembly (nullptr),
50 _p (nullptr),
51 _p0 (0.),
52 _X_scale (0.),
53 _p_scale (0.) {
54 
55 }
56 
57 
58 
60 
61 }
62 
63 
64 
65 void
68  MAST::AssemblyBase& assembly,
69  MAST::Parameter& p) {
70 
71  libmesh_assert(!_elem_ops);
72 
73  _elem_ops = &elem_ops;
74  _assembly = &assembly;
75  _p = &p;
76 
77 
78 }
79 
80 
81 void
83 
84  _elem_ops = nullptr;
85  _assembly = nullptr;
86  _p = nullptr;
87 
88 }
89 
90 
91 
92 void
94 
95  libmesh_assert(_initialized);
96 
97  unsigned int
98  iter = 0;
99 
100  libMesh::NumericVector<Real>
101  &X = *_assembly->system().solution;
102 
103  _p0 = (*_p)();
104  _X0.reset(X.clone().release());
105 
106  Real
107  norm0 = _res_norm(X, *_p),
108  norm = norm0;
109 
110  bool
111  cont = true;
112 
113  // save data for possible reuse if the iterations are restarted.
115 
116  while (cont) {
117 
118  libMesh::out
119  << std::setw(10) << "iter: "
120  << std::setw(5) << iter
121  << std::setw(10) << "res-l2: "
122  << std::setw(15) << norm
123  << std::setw(20) << "relative res-l2: "
124  << std::setw(15) << norm/norm0 << std::endl;
125 
126  _solve_NR_iterate(X, *_p);
127  norm = _res_norm(X, *_p);
128  iter++;
129 
130  if (norm < abs_tol) cont = false;
131  if (norm/norm0 < rel_tol) cont = false;
132  if (iter >= max_it) {
133  if (arc_length > min_step) {
134 
135  // reduce step-size if possible, otherwise terminate
136  libMesh::out
137  << " Retrying with smaller step-size" << std::endl;
138 
139  // reanalyze with a smaller step-size
140  iter = 0;
141  arc_length = std::max(min_step, .5*arc_length);
142  X.zero();
143  X.add(1., *_X0);
144  X.close();
145  *_p = _p0;
147  cont = true;
148  }
149  else
150  cont = false;
151  }
152  }
153 
154  libMesh::out
155  << std::setw(10) << "iter: "
156  << std::setw(5) << iter
157  << std::setw(10) << "res-l2: "
158  << std::setw(15) << norm
159  << std::setw(20) << "relative res-l2: "
160  << std::setw(15) << norm/norm0
161  << std::setw(20) << "Terminated" << std::endl;
162 
163  if (iter) {
164  Real
165  factor = std::pow((1.*step_desired_iters)/(1.*iter+1.), step_size_change_exponent);
166  if (factor > 1.) {
167  arc_length = std::min(max_step, factor*arc_length);
168  libMesh::out
169  << std::setw(30) << "increased step-size: "
170  << std::setw(15) << arc_length << std::endl;
171  }
172  else if (factor < 1.) {
173  arc_length = std::max(min_step, factor*arc_length);
174  libMesh::out
175  << std::setw(30) << "reduced step-size: "
176  << std::setw(15) << arc_length << std::endl;
177  }
178  }
179 
180 }
181 
182 
183 void
184 MAST::ContinuationSolverBase::_solve(const libMesh::NumericVector<Real> &X,
185  const MAST::Parameter &p,
186  libMesh::NumericVector<Real> &f,
187  bool update_f,
188  libMesh::NumericVector<Real> &dfdp,
189  bool update_dfdp,
190  const libMesh::NumericVector<Real> &dgdX,
191  const Real dgdp,
192  const Real g,
193  libMesh::NumericVector<Real> &dX,
194  Real &dp) {
195 
196  libmesh_assert(_elem_ops);
197 
199  &system = _assembly->system();
200 
201  libmesh_assert(system.operation() == MAST::NonlinearSystem::NONE);
202 
203  //
204  // { f(X, p) } = { 0 }
205  // { g(X, p) } = { 0 }
206  //
207  // [df/dX df/dp] { dX } = { -f}
208  // [dg/dX dg/dp] { dp } = { -g}
209  //
210 
211  // create matrix, vector and solver objects for the enlarged system
212  // first, the sparsity pattern.
213  // We will append one row and column to the sparsity pattern of the
214  // system. The highest rank cpu will own this new unknown so that
215  // the dof ids do not change
216 
218  // create the new matrix and vector quantities for the
219  // combined system
221  libMesh::DofMap& dof_map = system.get_dof_map();
222  const libMesh::Parallel::Communicator
223  &comm = system.comm();
224 
225  PetscInt
226  rank = comm.rank(),
227  my_m = dof_map.n_dofs(),
228  my_n = my_m,
229  total_m = my_m + 1, // size of the system with the constraint equation
230  total_n = total_m,
231  n_l = dof_map.n_dofs_on_processor(system.processor_id()),
232  m_l = n_l,
233  total_n_l = n_l,
234  total_m_l = n_l;
235 
236  std::vector<libMesh::dof_id_type>
237  n_nz = dof_map.get_n_nz(), // number of non-zeros per row in the diagonal block on this rank
238  n_oz = dof_map.get_n_oz(); // number of non-zeros per row in the off-diagonal blocks on other ranks
239 
240 
241  if (rank == comm.size()-1) {
242 
243  // the last rank will own the new dof, so add one more non-zero to
244  // the rows on this rank
245  total_n_l++;
246  total_m_l++;
247 
248  // add one more non-zero to the block matrix
249  for (unsigned int i=0; i<n_nz.size(); i++)
250  n_nz[i]++;
251  // the final row on this rank will be a full row, that is all
252  // entries of this row will be non-zero
253  n_nz.push_back(total_n_l);
254  // likewise, all entries in the off-diagonal block will be nonzero
255  n_oz.push_back(total_n - total_n_l);
256  }
257  else {
258 
259  // add an extra non-zero in the off-diagonal block
260  for (unsigned int i=0; i<n_oz.size(); i++)
261  n_oz[i]++;
262  }
263 
264  // create the matrix
265  PetscErrorCode ierr;
266  Mat mat;
267 
268  ierr = MatCreate(comm.get(), &mat); CHKERRABORT(comm.get(), ierr);
269  ierr = MatSetSizes(mat,
270  total_m_l, total_n_l,
271  total_m, total_n); CHKERRABORT(comm.get(), ierr);
272 
273  if (libMesh::on_command_line("--solver_system_names")) {
274 
275  std::string nm = _assembly->system().name() + "_continuation_";
276  MatSetOptionsPrefix(mat, nm.c_str());
277  }
278  ierr = MatSetFromOptions(mat); CHKERRABORT(comm.get(), ierr);
279 
280  ierr = MatSeqAIJSetPreallocation(mat,
281  my_m,
282  (PetscInt*)&n_nz[0]); CHKERRABORT(comm.get(), ierr);
283  ierr = MatMPIAIJSetPreallocation(mat,
284  0,
285  (PetscInt*)&n_nz[0],
286  0,
287  (PetscInt*)&n_oz[0]); CHKERRABORT(comm.get(), ierr);
288  ierr = MatSeqBAIJSetPreallocation (mat, 1,
289  0, (PetscInt*)&n_nz[0]); CHKERRABORT(comm.get(), ierr);
290  ierr = MatMPIBAIJSetPreallocation (mat, 1,
291  0, (PetscInt*)&n_nz[0],
292  0, (PetscInt*)&n_oz[0]); CHKERRABORT(comm.get(), ierr);
293 
294  // now add the block entries
295  ierr = MatSetOption(mat,
296  MAT_NEW_NONZERO_ALLOCATION_ERR,
297  PETSC_TRUE); CHKERRABORT(comm.get(), ierr);
298 
299  // now create the vectors
300  Vec res_vec, sol_vec;
301 
302  ierr = MatCreateVecs(mat, &res_vec, PETSC_NULL); CHKERRABORT(comm.get(), ierr);
303  ierr = MatCreateVecs(mat, &sol_vec, PETSC_NULL); CHKERRABORT(comm.get(), ierr);
304 
305 
306  std::unique_ptr<libMesh::SparseMatrix<Real> >
307  jac_mat(new libMesh::PetscMatrix<Real>(mat, comm));
308 
309  std::unique_ptr<libMesh::NumericVector<Real> >
310  res(new libMesh::PetscVector<Real>(res_vec, comm)),
311  sol(new libMesh::PetscVector<Real>(sol_vec, comm));
312 
314  // now, assemble the information into the new matrix/vector
316  // compute the dfdp term if asked for, since that goes into the matrix.
319  if (update_dfdp)
321 
322  // now compute the jacobian
324  _assembly->close_matrix = false;
326  update_f?res.get():nullptr,
327  jac_mat.get(),
328  system);
329  _assembly->close_matrix = true;
330 
333 
334  // first set the data for the last column
335  if (!update_f) {
336  for (unsigned int i=dof_map.first_dof(rank); i<dof_map.end_dof(rank); i++)
337  res->set(i, f.el(i));
338  }
339 
340  // diagonal entry on the last rank
341  if (rank == comm.size()-1)
342  res->set(total_m-1, g);
343 
344  // finish assembling the residual vector if it was not updated by
345  // the residual and jacobian routine
346  res->close();
347  res->scale(-1.);
348  res->close();
349 
350 
351  // first set the data for the last column
352  for (unsigned int i=dof_map.first_dof(rank);
353  i<dof_map.end_dof(rank); i++) {
354 
355  jac_mat->set( i, total_m-1, dfdp.el(i));
356  jac_mat->set(total_m-1, i, dgdX.el(i));
357  }
358 
359  // diagonal entry on the last rank
360  if (rank == comm.size()-1)
361  jac_mat->set(total_m-1, total_m-1, dgdp);
362 
363  jac_mat->close();
364 
365 
367  // now, copy the information to the new matrix/vector
369  KSP ksp;
370  PC pc;
371 
372  // setup the KSP
373  ierr = KSPCreate(comm.get(), &ksp); CHKERRABORT(comm.get(), ierr);
374 
375  if (libMesh::on_command_line("--solver_system_names")) {
376 
377  std::string nm = _assembly->system().name() + "_continuation_";
378  KSPSetOptionsPrefix(ksp, nm.c_str());
379  }
380 
381  ierr = KSPSetOperators(ksp, mat, mat); CHKERRABORT(comm.get(), ierr);
382  ierr = KSPSetFromOptions(ksp); CHKERRABORT(comm.get(), ierr);
383 
384  // setup the PC
385  ierr = KSPGetPC(ksp, &pc); CHKERRABORT(comm.get(), ierr);
386  ierr = PCSetFromOptions(pc); CHKERRABORT(comm.get(), ierr);
387 
388  // now solve
389  ierr = KSPSolve(ksp, res_vec, sol_vec);
390 
391  // copy the solution back to the system
392  dX.zero();
393  for (unsigned int i=dof_map.first_dof(rank);
394  i<dof_map.end_dof(rank); i++)
395  dX.set(i, sol->el(i));
396  dX.close();
397 
398 #ifdef LIBMESH_ENABLE_CONSTRAINTS
399  system.get_dof_map().enforce_constraints_exactly (system,
400  &dX,
401  /* homogeneous = */ true);
402 #endif
403 
404  std::vector<Real> val(1, 0.);
405  std::vector<libMesh::numeric_index_type> dofs(1, total_m-1);
406  sol->localize(val, dofs);
407  dp = val[0];
408 
409  // destroy the objects
410  ierr = KSPDestroy(&ksp);
411  ierr = MatDestroy(&mat);
412  ierr = VecDestroy(&res_vec);
413  ierr = VecDestroy(&sol_vec);
414 }
415 
416 
417 
418 
419 void
421 _solve_schur_factorization(const libMesh::NumericVector<Real> &X,
422  const MAST::Parameter &p,
423  libMesh::SparseMatrix<Real> &jac,
424  bool update_jac,
425  libMesh::NumericVector<Real> &f,
426  bool update_f,
427  libMesh::NumericVector<Real> &dfdp,
428  bool update_dfdp,
429  libMesh::NumericVector<Real> &dXdp,
430  bool update_dXdp,
431  const libMesh::NumericVector<Real> &dgdX,
432  const Real dgdp,
433  const Real g,
434  libMesh::NumericVector<Real> &dX,
435  Real &dp) {
436 
437  libmesh_assert(_elem_ops);
438 
440  &system = _assembly->system();
441 
442  libmesh_assert(system.operation() == MAST::NonlinearSystem::NONE);
443 
444  //
445  // { f(X, p) } = { 0 }
446  // { g(X, p) } = { 0 }
447  //
448  // [df/dX df/dp] { dX } = { -f}
449  // [dg/dX dg/dp] { dp } = { -g}
450  //
451  // Schur factorization:
452  // df/dX dX = -f - df/dp dp
453  //
454  // Substitute in second equation
455  // dg/dp dp - dg/dX inv(df/dX) ( f + df/dp dp) = -g
456  // => [dg/dp - dg/dX inv(df/dX) df/dp ] dp = -g + dg/dX inv(df/dX) f
457  // => [dg/dp + dg/dX dX/dp ] dp = -g + dg/dX r1
458  //
459  // 1. solve r1 = inv(df/dX) f
460  // 2. solve dXdp = -inv(df/dX) df/dp
461  // 3. solve a = dg/dp + dg/dX dXdp
462  // 4. solve a dp = -g + dg/dX r1
463  // 5. solve df/dX dX = -f - df/dp dp
464  //
465 
466 
467  std::pair<unsigned int, Real>
468  solver_params = system.get_linear_solve_parameters(),
469  rval;
470 
471 
472  std::unique_ptr<libMesh::NumericVector<Real>>
473  r1(X.zero_clone().release());
474 
475  libMesh::SparseMatrix<Real>
476  *pc = system.request_matrix("Preconditioner");
477 
478  Real
479  a = 0.;
480 
482  // STEP 1: r1 = inv(df/dX) f
486 
487  if (update_f || update_jac)
489  update_f? &f:nullptr,
490  update_jac? &jac:nullptr,
491  system);
492  rval = system.linear_solver->solve (jac, pc,
493  *r1,
494  f,
495  solver_params.second,
496  solver_params.first);
497 
498 #ifdef LIBMESH_ENABLE_CONSTRAINTS
499  system.get_dof_map().enforce_constraints_exactly (system,
500  r1.get(),
501  /* homogeneous = */ true);
502 #endif
503 
504 
506  // STEP 2: dXdp = - inv(df/dX) df/dp
509  if (update_dfdp)
511 
512  if (update_dfdp || update_dXdp) {
513 
514  rval = system.linear_solver->solve (jac, pc,
515  dXdp,
516  dfdp,
517  solver_params.second,
518  solver_params.first);
519 
520  dXdp.scale(-1.);
521  dXdp.close();
522 
523  // The linear solver may not have fit our constraints exactly
524 #ifdef LIBMESH_ENABLE_CONSTRAINTS
525  system.get_dof_map().enforce_constraints_exactly (system,
526  &dXdp,
527  /* homogeneous = */ true);
528 #endif
529  }
530 
532  // STEP 3: a = dg/dp + dg/dX dXdp
534  a = dgdp + dgdX.dot(dXdp);
535 
537  // STEP 4: a dp = -g + dg/dX r1
539  dp = 1./a * (- g + dgdX.dot(*r1));
540 
542  // STEP 5: df/dX dX = -f - df/dp dp
544  r1->zero();
545  r1->add(1., f);
546  r1->add(dp, dfdp);
547  r1->scale(-1.);
548  r1->close();
549  rval = system.linear_solver->solve (jac, pc,
550  dX,
551  *r1,
552  solver_params.second,
553  solver_params.first);
554 
555 
556  // The linear solver may not have fit our constraints exactly
557 #ifdef LIBMESH_ENABLE_CONSTRAINTS
558  system.get_dof_map().enforce_constraints_exactly (system,
559  &dX,
560  /* homogeneous = */ true);
561 #endif
562 
565 }
566 
567 
568 Real
570 _res_norm(const libMesh::NumericVector<Real> &X,
571  const MAST::Parameter &p) {
572 
573  libmesh_assert(_elem_ops);
574 
576  &system = _assembly->system();
577 
578  libmesh_assert(system.operation() == MAST::NonlinearSystem::NONE);
579 
580  Real
581  val = 0.;
582 
583  std::unique_ptr<libMesh::NumericVector<Real>>
584  f(X.zero_clone().release());
585 
588 
589  _assembly->residual_and_jacobian(X, f.get(), nullptr, system);
590 
593 
594  val = std::pow(std::pow(_g(X, p), 2) + std::pow(f->l2_norm(), 2), 0.5);
595 
596  return val;
597 }
598 
Real rel_tol
Relative tolerance for the solver.
Real max_step
maximum step size allowed with adaptivity
This class implements a system for solution of nonlinear systems.
virtual void clear_elem_operation_object()
clears the association of this object with the assembly element operation object. ...
virtual bool sensitivity_assemble(const MAST::FunctionBase &f, libMesh::NumericVector< Real > &sensitivity_rhs)
Assembly function.
virtual void solve()
solves for the next load step
This is a scalar function whose value can be changed and one that can be used as a design variable in...
Definition: parameter.h:35
bool close_matrix
flag to control the closing fo the Jacobian after assembly
virtual void set_elem_operation_object(MAST::AssemblyElemOperations &elem_ops)
attaches a element operation to this object, and associated this with the element operation object...
libMesh::Real Real
void _solve_schur_factorization(const libMesh::NumericVector< Real > &X, const MAST::Parameter &p, libMesh::SparseMatrix< Real > &jac, bool update_jac, libMesh::NumericVector< Real > &f, bool update_f, libMesh::NumericVector< Real > &dfdp, bool update_dfdp, libMesh::NumericVector< Real > &dXdp, bool update_dXdp, const libMesh::NumericVector< Real > &dgdX, const Real dgdp, const Real g, libMesh::NumericVector< Real > &dX, Real &dp)
solves for the linear system of equation using Schur factorization.
virtual void _solve_NR_iterate(libMesh::NumericVector< Real > &X, MAST::Parameter &p)=0
Real min_step
minimum step size allowed with adaptivity
virtual void _save_iteration_data()=0
method saves any data for possible resuse if the solution step is restarted
Real _res_norm(const libMesh::NumericVector< Real > &X, const MAST::Parameter &p)
void _solve(const libMesh::NumericVector< Real > &X, const MAST::Parameter &p, libMesh::NumericVector< Real > &f, bool update_f, libMesh::NumericVector< Real > &dfdp, bool update_dfdp, const libMesh::NumericVector< Real > &dgdX, const Real dgdp, const Real g, libMesh::NumericVector< Real > &dX, Real &dp)
solves for the linear system of equation as a monolithic system dX and dp are returned from the solu...
virtual std::pair< unsigned int, Real > get_linear_solve_parameters()
calls NonlinearImplicitSystem::set_solver_parameters() before accessing the values.
virtual Real _g(const libMesh::NumericVector< Real > &X, const MAST::Parameter &p)=0
Real arc_length
arc length that the solver is required to satisfy for the update.
Real abs_tol
Absolute tolerance for the solver.
std::unique_ptr< libMesh::NumericVector< Real > > _X0
void set_operation(MAST::NonlinearSystem::Operation op)
sets the current operation of the system
void set_assembly_and_load_parameter(MAST::AssemblyElemOperations &elem_ops, MAST::AssemblyBase &assembly, MAST::Parameter &p)
sets the assembly object for this solver
std::unique_ptr< libMesh::LinearSolver< Real > > linear_solver
The LinearSolver for solution of the linear equations.
virtual void residual_and_jacobian(const libMesh::NumericVector< Real > &X, libMesh::NumericVector< Real > *R, libMesh::SparseMatrix< Real > *J, libMesh::NonlinearImplicitSystem &S)
function that assembles the matrices and vectors quantities for nonlinear solution ...
MAST::NonlinearSystem::Operation operation()
virtual void _reset_iterations()=0
method resets any data if a soltion step is restarted
unsigned int max_it
Maximum number of Newton-Raphson iterations for the solver.
Real step_size_change_exponent
exponent used in step size update.
MAST::AssemblyElemOperations * _elem_ops
void clear_assembly_and_load_parameters()
clears the assembly object from this solver
const MAST::NonlinearSystem & system() const
unsigned int step_desired_iters
desired N-R iterations per load-step.