MAST
Multidisciplinary-design Adaptation and Sensitivity Toolkit (MAST)
assembly_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__system_assembly_base_h__
21 #define __mast__system_assembly_base_h__
22 
23 // C++ includes
24 #include <map>
25 #include <memory>
26 
27 
28 // MAST includes
29 #include "base/mast_data_types.h"
30 
31 
32 // libMesh includes
33 #include "libmesh/system.h"
34 #include "libmesh/nonlinear_implicit_system.h"
35 #include "libmesh/sparse_matrix.h"
36 
37 
38 namespace MAST {
39 
40  // Forward declerations
41  class PhysicsDisciplineBase;
42  class SystemInitialization;
43  class ElementBase;
44  class MeshFieldFunction;
45  class NonlinearSystem;
46  class FEBase;
47  class AssemblyElemOperations;
48  class OutputAssemblyElemOperations;
49  class FunctionBase;
50 
51  class AssemblyBase:
52  public libMesh::NonlinearImplicitSystem::ComputeResidualandJacobian {
53  public:
54 
59  AssemblyBase();
60 
64  virtual ~AssemblyBase();
65 
71  std::set<unsigned int> diagonal_elem_subdomain_id;
72 
73 
74  class SolverMonitor {
75  public:
77  virtual ~SolverMonitor(){}
78  virtual void init(MAST::AssemblyBase& assembly) = 0;
79  virtual void clear() = 0;
80  };
81 
94  public:
95  ElemParameterDependence(bool o_flag): override_flag(o_flag) {}
97  virtual bool if_elem_depends_on_parameter(const libMesh::Elem& e,
98  const MAST::FunctionBase& p) const = 0;
105  };
106 
111 
117 
123 
124 
129 
130 
134  virtual void
137 
144  void
146 
147 
148  void
150 
151 
155  virtual void
157 
158 
163  virtual void
165 
166 
171  virtual void
173 
174 
179  const MAST::NonlinearSystem& system() const;
180 
186 
191 
192 
198 
204 
208  void clear_solver_monitor();
209 
215 
216 
221 
226  virtual void
227  residual_and_jacobian (const libMesh::NumericVector<Real>& X,
228  libMesh::NumericVector<Real>* R,
229  libMesh::SparseMatrix<Real>* J,
230  libMesh::NonlinearImplicitSystem& S) {
231  libmesh_assert(false); // implement in the derived class
232  }
233 
245  virtual bool
246  sensitivity_assemble (const libMesh::NumericVector<Real>& X,
247  bool if_localize_sol,
248  const MAST::FunctionBase& f,
249  libMesh::NumericVector<Real>& sensitivity_rhs,
250  bool close_vector = true) {
251  libmesh_assert(false); // implemented in the derived class
252  }
253 
257  virtual void
258  calculate_output(const libMesh::NumericVector<Real>& X,
259  bool if_localize_sol,
261 
262 
266  virtual void
267  calculate_output_derivative(const libMesh::NumericVector<Real>& X,
268  bool if_localize_sol,
270  libMesh::NumericVector<Real>& dq_dX);
271 
272 
282  virtual void
283  calculate_output_direct_sensitivity(const libMesh::NumericVector<Real>& X,
284  bool if_localize_sol,
285  const libMesh::NumericVector<Real>* dXdp,
286  bool if_localize_sol_sens,
287  const MAST::FunctionBase& p,
289 
290 
296  virtual Real
297  calculate_output_adjoint_sensitivity(const libMesh::NumericVector<Real>& X,
298  bool if_localize_sol,
299  const libMesh::NumericVector<Real>& adj_sol,
300  const MAST::FunctionBase& p,
303  const bool include_partial_sens = true);
304 
305 
313  virtual void
315  (const libMesh::NumericVector<Real>& X,
316  bool if_localize_sol,
317  const libMesh::NumericVector<Real>& adj_sol,
318  const std::vector<const MAST::FunctionBase*>& p_vec,
321  std::vector<Real>& sens);
322 
323 
329  std::unique_ptr<libMesh::NumericVector<Real> >
330  build_localized_vector(const libMesh::System& sys,
331  const libMesh::NumericVector<Real>& global) const;
332 
333 
334  protected:
335 
340 
341 
346 
347 
352 
357 
363 
370  };
371 
372 }
373 
374 #endif // __mast__system_assembly_base_h__
375 
MAST::AssemblyElemOperations * _elem_ops
provides assembly elem operations for use by this class
const MAST::NonlinearSystem & system() const
AssemblyBase()
constructor takes a reference to the discipline that provides the boundary conditions, volume loads, properties, etc.
void clear_solver_monitor()
clears the monitor object
This class implements a system for solution of nonlinear systems.
virtual void calculate_output_derivative(const libMesh::NumericVector< Real > &X, bool if_localize_sol, MAST::OutputAssemblyElemOperations &output, libMesh::NumericVector< Real > &dq_dX)
calculates
bool override_flag
if true, assume zero solution sensitivity when elem does not dependent on parameter.
virtual Real calculate_output_adjoint_sensitivity(const libMesh::NumericVector< Real > &X, bool if_localize_sol, const libMesh::NumericVector< Real > &adj_sol, const MAST::FunctionBase &p, MAST::AssemblyElemOperations &elem_ops, MAST::OutputAssemblyElemOperations &output, const bool include_partial_sens=true)
Evaluates the total sensitivity of output wrt p using the adjoint solution provided in adj_sol for a ...
virtual void clear_elem_operation_object()
clears the association of this object with the assembly element operation object. ...
This provides a wrapper FieldFunction compatible class that interpolates the solution using libMesh&#39;s...
void attach_solution_function(MAST::MeshFieldFunction &f)
tells the assembly object that this function is will need to be initialized before each residual eval...
void attach_elem_parameter_dependence_object(MAST::AssemblyBase::ElemParameterDependence &dep)
This object, if provided by user, will be used to reduce unnecessary computations in sensitivity anal...
void set_solver_monitor(MAST::AssemblyBase::SolverMonitor &monitor)
attaches the solver monitor, which is a user provided routine that is called each time ...
This provides the base class for definitin of element level contribution of output quantity in an ana...
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
std::set< unsigned int > diagonal_elem_subdomain_id
subdomain ids for which residuakl and Jacobian contributions will not be computed.
Definition: assembly_base.h:71
MAST::SystemInitialization * _system
System for which this assembly is performed.
virtual void init(MAST::AssemblyBase &assembly)=0
virtual void calculate_output(const libMesh::NumericVector< Real > &X, bool if_localize_sol, MAST::OutputAssemblyElemOperations &output)
calculates the value of quantity .
void detach_solution_function()
removes the attachment of the solution function
MAST::PhysicsDisciplineBase * _discipline
PhysicsDisciplineBase object for which this class is assembling.
virtual void clear_discipline_and_system()
clears association with a system to this discipline
std::unique_ptr< libMesh::NumericVector< Real > > build_localized_vector(const libMesh::System &sys, const libMesh::NumericVector< Real > &global) const
localizes the parallel vector so that the local copy stores all values necessary for calculation of t...
virtual void calculate_output_direct_sensitivity(const libMesh::NumericVector< Real > &X, bool if_localize_sol, const libMesh::NumericVector< Real > *dXdp, bool if_localize_sol_sens, const MAST::FunctionBase &p, MAST::OutputAssemblyElemOperations &output)
evaluates the sensitivity of the outputs in the attached discipline with respect to the parametrs in ...
void clear_elem_parameter_dependence_object()
MAST::AssemblyBase::SolverMonitor * get_solver_monitor()
Inherited objects from this class can be provided by the user provide assessment of whether or not an...
Definition: assembly_base.h:93
virtual void calculate_output_adjoint_sensitivity_multiple_parameters_no_direct(const libMesh::NumericVector< Real > &X, bool if_localize_sol, const libMesh::NumericVector< Real > &adj_sol, const std::vector< const MAST::FunctionBase *> &p_vec, MAST::AssemblyElemOperations &elem_ops, MAST::OutputAssemblyElemOperations &output, std::vector< Real > &sens)
Evaluates the dot product between adj_sol and sensitivity of residual about X for multiple parameter_...
MAST::SystemInitialization & system_init()
MAST::AssemblyBase::ElemParameterDependence * _param_dependence
If provided by user, this object is used by sensitiivty analysis to check for whether or the current ...
virtual void set_discipline_and_system(MAST::PhysicsDisciplineBase &discipline, MAST::SystemInitialization &system)
attaches a system to this discipline
const MAST::PhysicsDisciplineBase & discipline() const
virtual ~AssemblyBase()
virtual destructor
MAST::AssemblyElemOperations & get_elem_ops()
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::AssemblyBase::SolverMonitor * _solver_monitor
User provided solver monitor is attached to the linear nonlinear solvers, if provided.
MAST::MeshFieldFunction * _sol_function
system solution that will be initialized before each solution
virtual bool sensitivity_assemble(const libMesh::NumericVector< Real > &X, bool if_localize_sol, const MAST::FunctionBase &f, libMesh::NumericVector< Real > &sensitivity_rhs, bool close_vector=true)
Assembly function.