20 #ifndef __mast_libmesh_residual_and_jacobian_h__ 21 #define __mast_libmesh_residual_and_jacobian_h__ 31 #include <libmesh/nonlinear_implicit_system.h> 32 #include <libmesh/dof_map.h> 37 namespace libMeshWrapper {
39 template <
typename ScalarType,
46 "Scalar type of assembly and element operations must be same");
59 template <
typename VecType,
typename MatType,
typename ContextType>
65 Assert0( R || J,
"Atleast one assembled quantity should be specified.");
73 sol_accessor(*c.sys, X);
75 using elem_vector_t =
typename ElemOpsType::vector_t;
76 using elem_matrix_t =
typename ElemOpsType::matrix_t;
82 libMesh::MeshBase::const_element_iterator
83 el = c.mesh->active_local_elements_begin(),
84 end_el = c.mesh->active_local_elements_end();
86 for ( ; el != end_el; ++el) {
91 sol_accessor.
init(*c.elem);
93 res_e.setZero(sol_accessor.
n_dofs());
94 if (J) jac_e.setZero(sol_accessor.
n_dofs(), sol_accessor.
n_dofs());
97 _e_ops->compute(c, sol_accessor, res_e, J?&jac_e:
nullptr);
103 <ScalarType, VecType, MatType, elem_vector_t, elem_matrix_t>
104 (*R, *J, c.sys->get_dof_map(), sol_accessor.
dof_indices(), res_e, jac_e);
107 <ScalarType, VecType, elem_vector_t>
108 (*R, c.sys->get_dof_map(), sol_accessor.
dof_indices(), res_e);
111 <ScalarType, MatType, elem_matrix_t>
112 (*J, c.sys->get_dof_map(), sol_accessor.
dof_indices(), jac_e);
131 #endif // __mast_libmesh_residual_and_jacobian_h__ void assemble(ContextType &c, const VecType &X, VecType *R, MatType *J)
void init(const libMesh::Elem &e)
void set_elem_ops(ElemOpsType &e_ops)
void finalize(ValType &m)
const std::vector< libMesh::dof_id_type > & dof_indices() const
std::enable_if< Dim< 3, ScalarType >::typesource_load_multiplier(const SourceLoadFieldType *f, const SectionAreaType *s, ContextType &c) { Assert0(f, "Invalid pointer");Assert0(s, "Invalid pointer");return f-> value(c) *s -> value(c)
std::enable_if< std::is_same< ScalarType, real_t >::value, void >::type constrain_and_add_matrix_and_vector(VecType &v, MatType &m, const libMesh::DofMap &dof_map, std::vector< libMesh::dof_id_type > &dof_indices, SubVecType &v_sub, SubMatType &m_sub)
virtual ~ResidualAndJacobian()
std::enable_if< std::is_same< ScalarType, real_t >::value, void >::type constrain_and_add_vector(VecType &v, const libMesh::DofMap &dof_map, std::vector< libMesh::dof_id_type > &dof_indices, SubVecType &v_sub)
void set_finalize_jac(bool f)
#define Assert0(cond, msg)
std::enable_if< std::is_same< ScalarType, real_t >::value, void >::type constrain_and_add_matrix(MatType &m, const libMesh::DofMap &dof_map, std::vector< libMesh::dof_id_type > &dof_indices, SubMatType &m_sub)