MAST
transient_assembly.h
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 #ifndef __mast__transient_assembly__
21 #define __mast__transient_assembly__
22 
23 // MAST includes
24 #include "base/assembly_base.h"
25 
26 // libMesh includes
27 #include "libmesh/nonlinear_implicit_system.h"
28 
29 
30 
31 namespace MAST {
32 
33  // Forward declerations
34  class TransientSolverBase;
35  class TransientAssemblyElemOperations;
36 
37 
39  public MAST::AssemblyBase {
40  public:
41 
48 
49  public:
52  virtual void post_assembly(const libMesh::NumericVector<Real>& X,
53  libMesh::NumericVector<Real>* R,
54  libMesh::SparseMatrix<Real>* J,
55  libMesh::NonlinearImplicitSystem& S) = 0;
56  };
57 
62 
63 
68  virtual ~TransientAssembly();
69 
70 
77  void
79 
80 
85  virtual void
86  residual_and_jacobian (const libMesh::NumericVector<Real>& X,
87  libMesh::NumericVector<Real>* R,
88  libMesh::SparseMatrix<Real>* J,
89  libMesh::NonlinearImplicitSystem& S);
90 
99  virtual void
100  linearized_jacobian_solution_product(const libMesh::NumericVector<Real>& X,
101  const libMesh::NumericVector<Real>& dX,
102  libMesh::NumericVector<Real>& JdX,
103  libMesh::NonlinearImplicitSystem& S);
104 
115  virtual bool
117  libMesh::NumericVector<Real>& sensitivity_rhs);
118 
119 
120 
121  protected:
122 
128 
129  };
130 
131 
132 }
133 
134 #endif // __mast__transient_assembly__
TransientAssembly()
constructor associates this assembly object with the system
MAST::TransientAssembly::PostAssemblyOperation * _post_assembly
this object, if non-NULL is user-provided to perform actions after assembly and before returning to t...
PostAssemblyOperation()
virtual ~TransientAssembly()
destructor resets the association of this assembly object with the system
virtual void linearized_jacobian_solution_product(const libMesh::NumericVector< Real > &X, const libMesh::NumericVector< Real > &dX, libMesh::NumericVector< Real > &JdX, libMesh::NonlinearImplicitSystem &S)
calculates the product of the Jacobian and a perturbation in solution vector .
virtual ~PostAssemblyOperation()
virtual void post_assembly(const libMesh::NumericVector< Real > &X, libMesh::NumericVector< Real > *R, libMesh::SparseMatrix< Real > *J, libMesh::NonlinearImplicitSystem &S)=0
virtual bool sensitivity_assemble(const MAST::FunctionBase &f, libMesh::NumericVector< Real > &sensitivity_rhs)
Assembly function.
void set_post_assembly_operation(MAST::TransientAssembly::PostAssemblyOperation &post)
sets the PostAssemblyOperation object for use after assembly.
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 ...
user-provided object to perform actions after assembly and before returning to the solver...