G+Smo  24.08.0
Geometry + Simulation Modules
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
gsBiharmonicAssembler.hpp
1 
15 #pragma once
16 
17 #include <gsElasticity/gsBiharmonicAssembler.h>
18 
19 #include <gsPde/gsPoissonPde.h>
20 #include <gsElasticity/gsVisitorBiharmonic.h>
21 
22 namespace gismo
23 {
24 template<class T>
26  const gsMultiBasis<T> & basis,
27  const gsBoundaryConditions<T> & bconditions,
28  const gsFunction<T> & body_force)
29 {
30  gsPiecewiseFunction<T> rightHandSides;
31  rightHandSides.addPiece(body_force);
32  typename gsPde<T>::Ptr pde( new gsPoissonPde<T>(patches,bconditions,rightHandSides) );
33  m_bases.push_back(basis);
34  m_bases.push_back(basis);
35  Base::initialize(pde, m_bases, defaultOptions());
36 }
37 
38 template <class T>
40 {
41  gsOptionList opt = Base::defaultOptions();
42  opt.addReal("LocalStiff","Stiffening degree for the Jacobian-based local stiffening",0.);
43  return opt;
44 }
45 
46 template <class T>
48 {
49  std::vector<gsDofMapper> m_dofMappers(2);
50  for (unsigned d = 0; d < 2; d++)
51  m_bases[d].getMapper((dirichlet::strategy)m_options.getInt("DirichletStrategy"),
52  iFace::glue,m_pde_ptr->bc(),m_dofMappers[d],d,true);
53 
54  m_system = gsSparseSystem<T>(m_dofMappers, gsVector<index_t>::Ones(2));
55  reserve();
56  Base::computeDirichletDofs(0);
57  Base::computeDirichletDofs(1);
58 }
59 
60 template <class T>
62 {
63  // Pick up values from options
64  const T bdA = m_options.getReal("bdA");
65  const index_t bdB = m_options.getInt("bdB");
66  const T bdO = m_options.getReal("bdO");
67 
68  index_t dim = m_bases[0][0].dim();
69  index_t deg = 0;
70  for (index_t d = 0; d < dim; ++d )
71  if (m_bases[0][0].degree(d) > deg)
72  deg = m_bases[0][0].degree(d);
73 
74  index_t numElPerColumn = pow((bdA*deg+bdB),dim)*2;
75  m_system.reserve(numElPerColumn*(1+bdO),m_pde_ptr->numRhs());
76 }
77 
78 template<class T>
79 void gsBiharmonicAssembler<T>::assemble(bool saveEliminationMatrix)
80 {
81  m_system.matrix().setZero();
82  reserve();
83  m_system.rhs().setZero(Base::numDofs(),m_pde_ptr->numRhs());
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  gsVisitorBiharmonic<T> visitor(*m_pde_ptr, saveEliminationMatrix ? &eliminationMatrix : nullptr);
93  Base::template push<gsVisitorBiharmonic<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 //--------------------- SOLUTION CONSTRUCTION ----------------------------------//
105 
106 template <class T>
108  const std::vector<gsMatrix<T> > & fixedDoFs,
109  gsMultiPatch<T> & solution) const
110 {
111  Base::constructSolution(solVector,fixedDoFs,solution,gsVector<index_t>::Zero(1));
112 }
113 
114 template <class T>
116  const std::vector<gsMatrix<T> > & fixedDoFs,
117  gsMultiPatch<T> & solutionAux) const
118 {
119  Base::constructSolution(solVector,fixedDoFs,solutionAux,gsVector<index_t>::Ones(1));
120 }
121 
122 template <class T>
124  const std::vector<gsMatrix<T> > & fixedDoFs,
125  gsMultiPatch<T> & solutionMain, gsMultiPatch<T> & solutionAux) const
126 {
127  constructSolution(solVector,fixedDoFs,solutionMain);
128  constructSolutionAux(solVector,fixedDoFs,solutionAux);
129 }
130 
131 }// namespace gismo ends
void addPiece(const gsFunction< T > &func)
Add a piece.
Definition: gsPiecewiseFunction.h:113
Describes a Poisson PDE.
void refresh()
Creates the mappers and setups the sparse system. to be implemented in derived classes, see scalarProblemGalerkinRefresh() for a possible implementation.
Definition: gsBiharmonicAssembler.hpp:23
static gsOptionList defaultOptions()
Returns the list of default options for assembly.
Definition: gsBiharmonicAssembler.hpp:39
#define index_t
Definition: gsConfig.h:32
A function from a n-dimensional domain to an m-dimensional image.
Definition: gsFunction.h:59
virtual void constructSolution(const gsMatrix< T > &solVector, const std::vector< gsMatrix< T > > &fixedDoFs, gsMultiPatch< T > &solution) const
construct the solution of the equation
Definition: gsBiharmonicAssembler.hpp:107
Visitor for the biharmonic equation.
Definition: gsBiharmonicAssembler.h:25
A vector with arbitrary coefficient type and fixed or dynamic size.
Definition: gsVector.h:35
Holds a set of patch-wise bases and their topology information.
Definition: gsMultiBasis.h:36
gsBiharmonicAssembler(gsMultiPatch< T > const &patches, gsMultiBasis< T > const &bases, gsBoundaryConditions< T > const &bconditions, gsBoundaryConditions< T > const &bconditions2, const gsFunction< T > &rhs, dirichlet::strategy dirStrategy, iFace::strategy intStrategy=iFace::glue)
Constructor of the assembler object.
Definition: gsBiharmonicAssembler.h:56
void assemble()
Main assemble routine, to be implemented in derived classes.
Definition: gsBiharmonicAssembler.hpp:31
virtual void constructSolutionAux(const gsMatrix< T > &solVector, const std::vector< gsMatrix< T > > &fixedDoFs, gsMultiPatch< T > &solutionAux) const
construct the Laplacian of the solution
Definition: gsBiharmonicAssembler.hpp:115
A Poisson PDE.
Definition: gsPoissonPde.h:34
Container class for a set of geometry patches and their topology, that is, the interface connections ...
Definition: gsMultiPatch.h:33
Class containing a set of boundary conditions.
Definition: gsBoundaryConditions.h:341
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
A function depending on an index i, typically referring to a patch/sub-domain. On each patch a differ...
Definition: gsPiecewiseFunction.h:28
virtual void reserve()
a custom reserve function to allocate memory for the sparse matrix
Definition: gsBiharmonicAssembler.hpp:61