MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
continuation_solver_base.h
Go to the documentation of this file.
1 /*
2  * MAST: Multidisciplinary-design Adaptation and Sensitivity Toolkit
3  * Copyright (C) 2013-2020 Manav Bhatia and MAST authors
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 #ifndef __mast__continuation_solver_base_h__
21 #define __mast__continuation_solver_base_h__
22 
23 // MAST includes
24 #include "base/mast_data_types.h"
25 
26 // libMesh includes
27 #include "libmesh/numeric_vector.h"
28 
29 
30 namespace MAST {
31 
32  // Forward declerations
33  class AssemblyElemOperations;
34  class AssemblyBase;
35  class Parameter;
36 
54 
55  public:
56 
58 
59 
60  virtual ~ContinuationSolverBase();
61 
65  void
67  MAST::AssemblyBase& assembly,
68  MAST::Parameter& p);
69 
74 
79  virtual void initialize(Real dp) =0;
80 
84  virtual void solve();
85 
90  unsigned int
92 
96  Real
98 
102  Real
104 
108  Real
110 
114  Real
116 
120  Real
122 
126  Real
128 
135  unsigned int
137 
142 
143  protected:
144 
145  virtual void
146  _solve_NR_iterate(libMesh::NumericVector<Real> &X,
147  MAST::Parameter &p) = 0;
148 
159  void
160  _solve(const libMesh::NumericVector<Real> &X,
161  const MAST::Parameter &p,
162  libMesh::NumericVector<Real> &f,
163  bool update_f,
164  libMesh::NumericVector<Real> &dfdp,
165  bool update_dfdp,
166  const libMesh::NumericVector<Real> &dgdX,
167  const Real dgdp,
168  const Real g,
169  libMesh::NumericVector<Real> &dX,
170  Real &dp);
171 
172 
183  void
184  _solve_schur_factorization(const libMesh::NumericVector<Real> &X,
185  const MAST::Parameter &p,
186  libMesh::SparseMatrix<Real> &jac,
187  bool update_jac,
188  libMesh::NumericVector<Real> &f,
189  bool update_f,
190  libMesh::NumericVector<Real> &dfdp,
191  bool update_dfdp,
192  libMesh::NumericVector<Real> &dXdp,
193  bool update_dXdp,
194  const libMesh::NumericVector<Real> &dgdX,
195  const Real dgdp,
196  const Real g,
197  libMesh::NumericVector<Real> &dX,
198  Real &dp);
199 
204  Real _res_norm(const libMesh::NumericVector<Real> &X,
205  const MAST::Parameter &p);
206 
207  virtual Real
208  _g(const libMesh::NumericVector<Real> &X,
209  const MAST::Parameter &p) = 0;
210 
214  virtual void _save_iteration_data() = 0;
215 
219  virtual void _reset_iterations() = 0;
220 
222 
226 
227  Real
229  _X_scale,
230  _p_scale;
231 
232  std::unique_ptr<libMesh::NumericVector<Real>>
234  };
235 }
236 
237 
238 #endif // __mast__continuation_solver_base_h__
Real rel_tol
Relative tolerance for the solver.
Real max_step
maximum step size allowed with adaptivity
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
libMesh::Real Real
virtual void initialize(Real dp)=0
initializes the data structure based on initial load step dp.
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...
the equation set is: the N-R updates are calculated such that This equation is solved using Schur-f...
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_assembly_and_load_parameter(MAST::AssemblyElemOperations &elem_ops, MAST::AssemblyBase &assembly, MAST::Parameter &p)
sets the assembly object for this solver
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
unsigned int step_desired_iters
desired N-R iterations per load-step.
bool schur_factorization
flag to use Schur-factorizaiton (default) or monolithic solver