21 #include <gsElasticity/gsVisitorMass.h>
27 gsMassAssembler<T>::gsMassAssembler(
const gsMultiPatch<T> & patches,
28 const gsMultiBasis<T> & basis,
29 const gsBoundaryConditions<T> & bconditions,
30 const gsFunction<T> & body_force)
36 gsPiecewiseFunction<T> rightHandSides;
37 rightHandSides.addPiece(body_force);
38 typename gsPde<T>::Ptr pde(
new gsBasePde<T>(patches,bconditions,rightHandSides) );
43 m_dim = body_force.targetDim();
44 for (
short_t d = 0; d < m_dim; ++d)
45 m_bases.push_back(basis);
47 Base::initialize(pde, m_bases, defaultOptions());
54 opt.
addReal(
"Density",
"Density of the material",1.);
61 GISMO_ENSURE(m_dim == m_pde_ptr->domain().parDim(),
"The RHS dimension and the domain dimension don't match!");
62 GISMO_ENSURE(m_dim == 2 || m_dim == 3,
"Only two- and three-dimenstion domains are supported!");
64 std::vector<gsDofMapper> m_dofMappers(m_bases.size());
65 for (
unsigned d = 0; d < m_bases.size(); d++)
66 m_bases[d].getMapper((dirichlet::strategy)m_options.getInt(
"DirichletStrategy"),
67 iFace::glue,m_pde_ptr->bc(),m_dofMappers[d],d,
true);
71 m_options.setReal(
"bdO",m_bases.size()*(1+m_options.getReal(
"bdO"))-1);
72 m_system.reserve(m_bases[0], m_options, 1);
74 for (
unsigned d = 0; d < m_bases.size(); ++d)
75 Base::computeDirichletDofs(d);
82 m_system.matrix().setZero();
83 m_system.reserve(m_bases[0], m_options, 1);
84 m_system.rhs().setZero(Base::numDofs(),1);
86 if (saveEliminationMatrix)
88 eliminationMatrix.resize(Base::numDofs(),Base::numFixedDofs());
89 eliminationMatrix.setZero();
90 eliminationMatrix.reservePerColumn(m_system.numColNz(m_bases[0],m_options));
93 gsVisitorMass<T> visitor(saveEliminationMatrix ? &eliminationMatrix :
nullptr);
94 Base::template push<gsVisitorMass<T> >(visitor);
96 m_system.matrix().makeCompressed();
98 if (saveEliminationMatrix)
100 Base::rhsWithZeroDDofs = m_system.rhs();
101 eliminationMatrix.makeCompressed();
static gsOptionList defaultOptions()
Returns the list of default options for assembly.
Definition: gsMassAssembler.hpp:51
#define short_t
Definition: gsConfig.h:35
Provides mass matrix for elasticity systems in 2D plain strain and 3D continua.
#define GISMO_ENSURE(cond, message)
Definition: gsDebug.h:102
The visitor computes element mass integrals.
Definition: gsVisitorMass.h:29
A vector with arbitrary coefficient type and fixed or dynamic size.
Definition: gsVector.h:35
virtual void refresh()
Refresh routine to set dof-mappers.
Definition: gsMassAssembler.hpp:59
IMHO, a useless class, but it is necessary to use the gsAssembler class. Contains proper information ...
void addReal(const std::string &label, const std::string &desc, const Real &value)
Adds a option named label, with description desc and value value.
Definition: gsOptionList.cpp:211
Class which holds a list of parameters/options, and provides easy access to them. ...
Definition: gsOptionList.h:32
virtual void assemble()
Main assemble routine, to be implemented in derived classes.
Definition: gsAssembler.hpp:51