G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsMassAssembler.hpp
1
16#pragma once
17
21
22namespace gismo
23{
24
25template<class T>
26gsMassAssembler<T>::gsMassAssembler(const gsMultiPatch<T> & patches,
27 const gsMultiBasis<T> & basis,
28 const gsBoundaryConditions<T> & bconditions,
29 const gsFunction<T> & body_force)
30{
31 // Originally concieved as a meaningful class, now gsPde is just a container for
32 // the domain, boundary conditions and the right-hand side;
33 // any derived class can surve this purpuse, for example gsPoissonPde;
34 // TUDO: change/remove gsPde from gsAssembler logic
35 gsPiecewiseFunction<T> rightHandSides;
36 rightHandSides.addPiece(body_force);
37 typename gsPde<T>::Ptr pde( new gsBasePde<T>(patches,bconditions,rightHandSides) );
38 // gsAssembler<>::initialize requires a vector of bases, one for each unknown;
39 // different bases are used to compute Dirichlet DoFs;
40 // but always the first basis is used for the assembly;
41 // TODO: change gsAssembler logic
42 m_dim = body_force.targetDim();
43 for (short_t d = 0; d < m_dim; ++d)
44 m_bases.push_back(basis);
45
46 Base::initialize(pde, m_bases, defaultOptions());
47}
48
49template <class T>
51{
52 gsOptionList opt = Base::defaultOptions();
53 opt.addReal("Density","Density of the material",1.);
54 return opt;
55}
56
57template <class T>
59{
60 GISMO_ENSURE(m_dim == m_pde_ptr->domain().parDim(), "The RHS dimension and the domain dimension don't match!");
61 GISMO_ENSURE(m_dim == 2 || m_dim == 3, "Only two- and three-dimenstion domains are supported!");
62
63 std::vector<gsDofMapper> m_dofMappers(m_bases.size());
64 for (unsigned d = 0; d < m_bases.size(); d++)
65 m_bases[d].getMapper((dirichlet::strategy)m_options.getInt("DirichletStrategy"),
66 iFace::glue,m_pde_ptr->bc(),m_dofMappers[d],d,true);
67
68 m_system = gsSparseSystem<T>(m_dofMappers,gsVector<index_t>::Ones(m_bases.size()));
69
70 m_options.setReal("bdO",m_bases.size()*(1+m_options.getReal("bdO"))-1);
71 m_system.reserve(m_bases[0], m_options, 1);
72
73 for (unsigned d = 0; d < m_bases.size(); ++d)
74 Base::computeDirichletDofs(d);
75}
76
77template<class T>
78void gsMassAssembler<T>::assemble(bool saveEliminationMatrix)
79{
80 // allocate space for the linear system
81 m_system.matrix().setZero();
82 m_system.reserve(m_bases[0], m_options, 1);
83 m_system.rhs().setZero(Base::numDofs(),1);
84
85 if (saveEliminationMatrix)
86 {
87 eliminationMatrix.resize(Base::numDofs(),Base::numFixedDofs());
88 eliminationMatrix.setZero();
89 eliminationMatrix.reservePerColumn(m_system.numColNz(m_bases[0],m_options));
90 }
91
92 gsVisitorMassElasticity<T> visitor(saveEliminationMatrix ? &eliminationMatrix : nullptr);
93 Base::template push<gsVisitorMassElasticity<T> >(visitor);
94
95 m_system.matrix().makeCompressed();
96
97 if (saveEliminationMatrix)
98 {
99 Base::rhsWithZeroDDofs = m_system.rhs();
100 eliminationMatrix.makeCompressed();
101 }
102}
103
104
105}// namespace gismo ends
virtual void assemble()
Main assemble routine, to be implemented in derived classes.
Definition gsMassAssembler.h:50
virtual void refresh()
Refresh routine to set dof-mappers.
Definition gsMassAssembler.hpp:58
static gsOptionList defaultOptions()
Returns the list of default options for assembly.
Definition gsMassAssembler.hpp:50
Class which holds a list of parameters/options, and provides easy access to them.
Definition gsOptionList.h:33
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
A sparse linear system indexed by sets of degrees of freedom.
Definition gsSparseSystem.h:30
A vector with arbitrary coefficient type and fixed or dynamic size.
Definition gsVector.h:37
IMHO, a useless class, but it is necessary to use the gsAssembler class. Contains proper information ...
#define short_t
Definition gsConfig.h:35
#define GISMO_ENSURE(cond, message)
Definition gsDebug.h:102
Provides mass matrix for elasticity systems in 2D plain strain and 3D continua.
Visitor class for the mass matrix assembly for elasticity problems.
The G+Smo namespace, containing all definitions for the library.