G+Smo  24.08.0
Geometry + Simulation Modules
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
gsMassAssembler.hpp
1 
16 #pragma once
17 
19 #include <gsElasticity/gsBasePde.h>
20 
21 #include <gsElasticity/gsVisitorMass.h>
22 
23 namespace gismo
24 {
25 
26 template<class T>
27 gsMassAssembler<T>::gsMassAssembler(const gsMultiPatch<T> & patches,
28  const gsMultiBasis<T> & basis,
29  const gsBoundaryConditions<T> & bconditions,
30  const gsFunction<T> & body_force)
31 {
32  // Originally concieved as a meaningful class, now gsPde is just a container for
33  // the domain, boundary conditions and the right-hand side;
34  // any derived class can surve this purpuse, for example gsPoissonPde;
35  // TUDO: change/remove gsPde from gsAssembler logic
36  gsPiecewiseFunction<T> rightHandSides;
37  rightHandSides.addPiece(body_force);
38  typename gsPde<T>::Ptr pde( new gsBasePde<T>(patches,bconditions,rightHandSides) );
39  // gsAssembler<>::initialize requires a vector of bases, one for each unknown;
40  // different bases are used to compute Dirichlet DoFs;
41  // but always the first basis is used for the assembly;
42  // TODO: change gsAssembler logic
43  m_dim = body_force.targetDim();
44  for (short_t d = 0; d < m_dim; ++d)
45  m_bases.push_back(basis);
46 
47  Base::initialize(pde, m_bases, defaultOptions());
48 }
49 
50 template <class T>
52 {
53  gsOptionList opt = Base::defaultOptions();
54  opt.addReal("Density","Density of the material",1.);
55  return opt;
56 }
57 
58 template <class T>
60 {
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!");
63 
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);
68 
69  m_system = gsSparseSystem<T>(m_dofMappers,gsVector<index_t>::Ones(m_bases.size()));
70 
71  m_options.setReal("bdO",m_bases.size()*(1+m_options.getReal("bdO"))-1);
72  m_system.reserve(m_bases[0], m_options, 1);
73 
74  for (unsigned d = 0; d < m_bases.size(); ++d)
75  Base::computeDirichletDofs(d);
76 }
77 
78 template<class T>
79 void gsMassAssembler<T>::assemble(bool saveEliminationMatrix)
80 {
81  // allocate space for the linear system
82  m_system.matrix().setZero();
83  m_system.reserve(m_bases[0], m_options, 1);
84  m_system.rhs().setZero(Base::numDofs(),1);
85 
86  if (saveEliminationMatrix)
87  {
88  eliminationMatrix.resize(Base::numDofs(),Base::numFixedDofs());
89  eliminationMatrix.setZero();
90  eliminationMatrix.reservePerColumn(m_system.numColNz(m_bases[0],m_options));
91  }
92 
93  gsVisitorMass<T> visitor(saveEliminationMatrix ? &eliminationMatrix : nullptr);
94  Base::template push<gsVisitorMass<T> >(visitor);
95 
96  m_system.matrix().makeCompressed();
97 
98  if (saveEliminationMatrix)
99  {
100  Base::rhsWithZeroDDofs = m_system.rhs();
101  eliminationMatrix.makeCompressed();
102  }
103 }
104 
105 
106 }// namespace gismo ends
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