Go to the documentation of this file.
20 #ifndef LIBMESH_IMPLICIT_SYSTEM_H
21 #define LIBMESH_IMPLICIT_SYSTEM_H
24 #include "libmesh/explicit_system.h"
25 #include "libmesh/auto_ptr.h"
34 template <
typename T>
class LinearSolver;
35 template <
typename T>
class SparseMatrix;
66 const std::string &
name,
67 const unsigned int number);
93 virtual void clear ()
override;
99 virtual void reinit ()
override;
118 virtual std::string
system_type ()
const override {
return "Implicit"; }
141 virtual std::pair<unsigned int, Real>
165 { libmesh_not_implemented(); }
186 { libmesh_not_implemented(); }
195 virtual std::pair<unsigned int, Real>
206 virtual std::pair<unsigned int, Real>
219 virtual std::pair<unsigned int, Real>
234 virtual std::pair<unsigned int, Real>
335 template <
template <
typename>
class>
347 bool have_matrix (
const std::string & mat_name)
const;
386 virtual unsigned int n_matrices ()
const override;
426 std::map<std::string, SparseMatrix<Number> *>
_matrices;
448 return cast_int<unsigned int>(
_matrices.size());
451 template <
template <
typename>
class MatrixType>
458 libmesh_error_msg(
"ERROR: Too late. Cannot add matrices to the system after initialization"
459 <<
"\n any more. You should have done this earlier.");
467 _matrices.insert (std::make_pair (mat_name, buf));
474 #endif // LIBMESH_IMPLICIT_SYSTEM_H
virtual void assemble_residual_derivatives(const ParameterVector ¶meters) override
Residual parameter derivative function.
const SparseMatrix< Number > * request_matrix(const std::string &mat_name) const
Manages consistently variables, degrees of freedom, coefficient vectors, and matrices for implicit sy...
virtual unsigned int n_matrices() const override
virtual std::pair< unsigned int, Real > get_linear_solve_parameters() const
virtual std::pair< unsigned int, Real > adjoint_solve(const QoISet &qoi_indices=QoISet()) override
Assembles & solves the linear system (dR/du)^T*z = dq/du, for those quantities of interest q specifie...
virtual void reinit() override
Reinitializes the member data fields associated with the system, so that, e.g., assemble() may be use...
The libMesh namespace provides an interface to certain functionality in the library.
SparseMatrix< Number > & add_matrix(const std::string &mat_name)
Adds the additional matrix mat_name to this system.
virtual void assemble() override
Prepares matrix and rhs for system assembly, then calls user assembly function.
virtual std::string system_type() const override
std::map< std::string, SparseMatrix< Number > * > _matrices
Some systems need an arbitrary number of matrices.
const Parallel::Communicator & comm() const
virtual void qoi_parameter_hessian_vector_product(const QoISet &qoi_indices, const ParameterVector ¶meters, const ParameterVector &vector, SensitivityData &product) override
For each of the system's quantities of interest q in qoi[qoi_indices], and for a vector of parameters...
Data structure for holding completed parameter sensitivity calculations.
unsigned int number() const
ImplicitSystem sys_type
The type of system.
virtual LinearSolver< Number > * get_linear_solver() const
ExplicitSystem Parent
The type of the parent.
std::map< std::string, SparseMatrix< Number > * >::const_iterator const_matrices_iterator
virtual std::pair< unsigned int, Real > weighted_sensitivity_solve(const ParameterVector ¶meters, const ParameterVector &weights) override
Assembles & solves the linear system(s) (dR/du)*u_w = sum(w_p*-dR/dp), for those parameters p contain...
virtual ~ImplicitSystem()
Destructor.
Data structure for specifying which Quantities of Interest should be calculated in an adjoint or a pa...
virtual void clear() override
Clear all the data structures associated with the system.
bool _can_add_matrices
true when additional matrices may still be added, false otherwise.
virtual void assembly(bool, bool, bool=false, bool=false)
Assembles a residual in rhs and/or a jacobian in matrix, as requested.
virtual void adjoint_qoi_parameter_sensitivity(const QoISet &qoi_indices, const ParameterVector ¶meters, SensitivityData &sensitivities) override
Solves for the derivative of each of the system's quantities of interest q in qoi[qoi_indices] with r...
virtual void init_matrices()
Initializes the matrices associated with this system.
SparseMatrix< Number > * matrix
The system matrix.
Data structure for specifying which Parameters should be independent variables in a parameter sensiti...
virtual std::pair< unsigned int, Real > sensitivity_solve(const ParameterVector ¶meters) override
Assembles & solves the linear system(s) (dR/du)*u_p = -dR/dp, for those parameters contained within p...
bool have_matrix(const std::string &mat_name) const
std::map< std::string, SparseMatrix< Number > * >::iterator matrices_iterator
Matrix iterator typedefs.
virtual void release_linear_solver(LinearSolver< Number > *) const
Releases a pointer to a linear solver acquired by this->get_linear_solver()
This is the EquationSystems class.
virtual void solve() override
For explicit systems, just assemble and solve the system A*x=b.
Manages consistently variables, degrees of freedom, and coefficient vectors for explicit systems.
ImplicitSystem(EquationSystems &es, const std::string &name, const unsigned int number)
Constructor.
const std::string & name() const
void remove_matrix(const std::string &mat_name)
Removes the additional matrix mat_name from this system.
virtual void forward_qoi_parameter_sensitivity(const QoISet &qoi_indices, const ParameterVector ¶meters, SensitivityData &sensitivities) override
Solves for the derivative of each of the system's quantities of interest q in qoi[qoi_indices] with r...
virtual std::pair< unsigned int, Real > weighted_sensitivity_adjoint_solve(const ParameterVector ¶meters, const ParameterVector &weights, const QoISet &qoi_indices=QoISet()) override
Assembles & solves the linear system(s) (dR/du)^T*z_w = sum(w_p*(d^2q/dudp - d^2R/dudp*z)),...
bool zero_out_matrix_and_rhs
By default, the system will zero out the matrix and the right hand side.
virtual void init_data() override
Initializes the member data fields associated with the system, so that, e.g., assemble() may be used.
const SparseMatrix< Number > & get_matrix(const std::string &mat_name) const
void add_system_matrix()
Add the system matrix to the _matrices data structure.
virtual void qoi_parameter_hessian(const QoISet &qoi_indices, const ParameterVector ¶meters, SensitivityData &hessian) override
For each of the system's quantities of interest q in qoi[qoi_indices], and for a vector of parameters...
virtual void disable_cache() override
Avoids use of any cached data that might affect any solve result.