G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsElasticityAssembler.hpp
Go to the documentation of this file.
1
16#pragma once
17
19
20#include <gsUtils/gsPointGrid.h>
24
25// Element visitors
31
32namespace gismo
33{
34
35template<class T>
37 const gsMultiBasis<T> & basis,
38 const gsBoundaryConditions<T> & bconditions,
39 const gsFunction<T> & body_force)
40{
41 // Originally concieved as a meaningful class, now gsPde is just a container for
42 // the domain, boundary conditions and the right-hand side;
43 // TUDO: change/remove gsPde from gsAssembler logic
44 gsPiecewiseFunction<T> rightHandSides;
45 rightHandSides.addPiece(body_force);
46 typename gsPde<T>::Ptr pde( new gsBasePde<T>(patches,bconditions,rightHandSides) );
47 // gsAssembler<>::initialize requires a vector of bases, one for each unknown;
48 // different bases are used to compute Dirichlet DoFs;
49 // but always the first basis is used for the assembly;
50 // TODO: change gsAssembler logic
51 m_dim = body_force.targetDim();
52 for (short_t d = 0; d < m_dim; ++d)
53 m_bases.push_back(basis);
54
55 Base::initialize(pde, m_bases, defaultOptions());
56}
57
58template<class T>
60 gsMultiBasis<T> const & basisDisp,
61 gsMultiBasis<T> const & basisPres,
62 gsBoundaryConditions<T> const & bconditions,
63 const gsFunction<T> & body_force)
64{
65 // same as above
66 gsPiecewiseFunction<T> rightHandSides;
67 rightHandSides.addPiece(body_force);
68 typename gsPde<T>::Ptr pde( new gsBasePde<T>(patches,bconditions,rightHandSides) );
69 // same as above
70 m_dim = body_force.targetDim();
71 for (short_t d = 0; d < m_dim; ++d)
72 m_bases.push_back(basisDisp);
73 m_bases.push_back(basisPres);
75 Base::initialize(pde, m_bases, defaultOptions());
76 m_options.setInt("MaterialLaw",material_law::mixed_hooke);
77}
78
79template <class T>
81{
82 gsOptionList opt = Base::defaultOptions();
83 opt.addReal("YoungsModulus","Youngs modulus of the material",200e9);
84 opt.addReal("PoissonsRatio","Poisson's ratio of the material",0.33);
85 opt.addReal("ForceScaling","Force scaling parameter",1.);
86 opt.addInt("MaterialLaw","Material law: 0 for St. Venant-Kirchhof, 1 for Neo-Hooke",material_law::hooke);
87 opt.addReal("LocalStiff","Stiffening degree for the Jacobian-based local stiffening",0.);
88 opt.addSwitch("Check","Check bijectivity of the displacement field before matrix assebmly",false);
89 return opt;
90}
91
92template <class T>
94{
95 // Pick up values from options
96 const T bdA = m_options.getReal("bdA");
97 const index_t bdB = m_options.getInt("bdB");
98 const T bdO = m_options.getReal("bdO");
99
100 index_t deg = 0;
101 for (index_t d = 0; d < m_bases[0][0].dim(); ++d )
102 if (m_bases[0][0].degree(d) > deg)
103 deg = m_bases[0][0].degree(d);
104
105 if (m_bases.size() == unsigned(m_dim)) // displacement formulation
106 {
107 // m_dim velocity*velocity blocks
108 index_t numElPerColumn = pow((bdA*deg+bdB),m_dim)*m_dim;
109 m_system.reserve(numElPerColumn*(1+bdO),1);
110 }
111 else // mixed formulation (displacement + pressure)
112 {
113 // m_dim velocity*velocity blocks + 1 pressure*velocity block (additioanal factor 2 for subgrid element)
114 index_t numElPerColumn = pow((bdA*deg+bdB),m_dim)*m_dim + pow((2*bdA*deg+bdB),m_dim);
115 m_system.reserve(numElPerColumn*(1+bdO),1);
116 }
117}
118
119template <class T>
121{
122 GISMO_ENSURE(m_dim == m_pde_ptr->domain().parDim(), "The RHS dimension and the domain dimension don't match!");
123 GISMO_ENSURE(m_dim == 2 || m_dim == 3, "Only two- and three-dimenstion domains are supported!");
124
125 std::vector<gsDofMapper> m_dofMappers(m_bases.size());
126 for (unsigned d = 0; d < m_bases.size(); d++)
127 m_bases[d].getMapper((dirichlet::strategy)m_options.getInt("DirichletStrategy"),
128 iFace::glue,m_pde_ptr->bc(),m_dofMappers[d],d,true);
129
130 m_system = gsSparseSystem<T>(m_dofMappers, gsVector<index_t>::Ones(m_bases.size()));
131 reserve();
132
133 for (unsigned d = 0; d < m_bases.size(); ++d)
134 Base::computeDirichletDofs(d);
135}
136
137//--------------------- SYSTEM ASSEMBLY ----------------------------------//
138
139template<class T>
140void gsElasticityAssembler<T>::assemble(bool saveEliminationMatrix)
141{
142 m_system.matrix().setZero();
143 reserve();
144 m_system.rhs().setZero();
145
146 // Compute volumetric integrals and write to the global linear system
147 if (m_bases.size() == unsigned(m_dim)) // displacement formulation
148 {
149 GISMO_ENSURE(m_options.getInt("MaterialLaw") == material_law::hooke,
150 "Material law not specified OR not supported!");
151 if (saveEliminationMatrix)
152 {
153 eliminationMatrix.resize(Base::numDofs(),Base::numFixedDofs());
154 eliminationMatrix.setZero();
155 eliminationMatrix.reservePerColumn(m_system.numColNz(m_bases[0],m_options));
156 }
157
158 gsVisitorLinearElasticity<T> visitor(*m_pde_ptr, saveEliminationMatrix ? &eliminationMatrix : nullptr);
159 Base::template push<gsVisitorLinearElasticity<T> >(visitor);
160
161 if (saveEliminationMatrix)
162 {
163 Base::rhsWithZeroDDofs = m_system.rhs();
164 eliminationMatrix.makeCompressed();
165 }
166
167 }
168 else // mixed formulation (displacement + pressure)
169 {
170 GISMO_ENSURE(m_options.getInt("MaterialLaw") == material_law::mixed_hooke,
171 "Material law not specified OR not supported!");
172 gsVisitorMixedLinearElasticity<T> visitor(*m_pde_ptr);
173 Base::template push<gsVisitorMixedLinearElasticity<T> >(visitor);
174 }
175
176 // Compute surface integrals and write to the global rhs vector
177 Base::template push<gsVisitorElasticityNeumann<T> >(m_pde_ptr->bc().neumannSides());
178
179 m_system.matrix().makeCompressed();
180}
181
182template <class T>
184 const std::vector<gsMatrix<T> > & fixedDoFs)
185{
186 gsMultiPatch<T> displacement;
187 constructSolution(solutionVector,fixedDoFs,displacement);
188 if (m_options.getSwitch("Check"))
189 if (checkDisplacement(m_pde_ptr->patches(),displacement) != -1)
190 return false;
191
192 if (m_bases.size() == unsigned(m_dim)) // displacement formulation
193 assemble(displacement);
194 else // mixed formulation (displacement + pressure)
195 {
196 gsMultiPatch<T> pressure;
197 constructPressure(solutionVector,fixedDoFs,pressure);
198 assemble(displacement,pressure);
199 }
200 return true;
201}
202
203template<class T>
205{
206 GISMO_ENSURE(m_options.getInt("MaterialLaw") == material_law::saint_venant_kirchhoff ||
207 m_options.getInt("MaterialLaw") == material_law::neo_hooke_ln ||
208 m_options.getInt("MaterialLaw") == material_law::neo_hooke_quad,
209 "Material law not specified OR not supported!");
210 m_system.matrix().setZero();
211 reserve();
212 m_system.rhs().setZero();
213
214 // Compute volumetric integrals and write to the global linear system
215 gsVisitorNonLinearElasticity<T> visitor(*m_pde_ptr,displacement);
216 Base::template push<gsVisitorNonLinearElasticity<T> >(visitor);
217 // Compute surface integrals and write to the global rhs vector
218 // change to reuse rhs from linear system
219 Base::template push<gsVisitorElasticityNeumann<T> >(m_pde_ptr->bc().neumannSides());
220
221 m_system.matrix().makeCompressed();
222}
223
224template<class T>
226 const gsMultiPatch<T> & pressure)
227{
228 GISMO_ENSURE(m_options.getInt("MaterialLaw") == material_law::mixed_neo_hooke_ln,
229 "Material law not specified OR not supported!");
230 m_options.setInt("MaterialLaw",material_law::mixed_neo_hooke_ln);
231 m_system.matrix().setZero();
232 reserve();
233 m_system.rhs().setZero();
234
235 // Compute volumetric integrals and write to the global linear systemz
236 gsVisitorMixedNonLinearElasticity<T> visitor(*m_pde_ptr,displacement,pressure);
237 Base::template push<gsVisitorMixedNonLinearElasticity<T> >(visitor);
238 // Compute surface integrals and write to the global rhs vector
239 // change to reuse rhs from linear system
240 Base::template push<gsVisitorElasticityNeumann<T> >(m_pde_ptr->bc().neumannSides());
241
242 m_system.matrix().makeCompressed();
243}
244
245//--------------------- SOLUTION CONSTRUCTION ----------------------------------//
246
247template <class T>
249 const std::vector<gsMatrix<T> > & fixedDoFs,
250 gsMultiPatch<T> & result) const
251{
252 gsVector<index_t> unknowns(m_dim);
253 for (short_t d = 0; d < m_dim; ++d)
254 unknowns.at(d) = d;
255 Base::constructSolution(solVector,fixedDoFs,result,unknowns);
256}
257
258template <class T>
260 const std::vector<gsMatrix<T> > & fixedDoFs,
261 gsMultiPatch<T> & displacement, gsMultiPatch<T> & pressure) const
262{
263 GISMO_ENSURE(m_bases.size() == unsigned(m_dim) + 1, "Not a mixed formulation: can't construct pressure.");
264 // construct displacement
265 constructSolution(solVector,fixedDoFs,displacement);
266 // construct pressure
267 constructPressure(solVector,fixedDoFs,pressure);
268}
269
270template <class T>
272 const std::vector<gsMatrix<T> > & fixedDoFs,
273 gsMultiPatch<T>& pressure) const
274{
275 GISMO_ENSURE(m_bases.size() == unsigned(m_dim) + 1, "Not a mixed formulation: can't construct pressure.");
276 gsVector<index_t> unknowns(1);
277 unknowns.at(0) = m_dim;
278 Base::constructSolution(solVector,fixedDoFs,pressure,unknowns);
279}
280
281//--------------------- SPECIALS ----------------------------------//
282
283template <class T>
285 gsPiecewiseFunction<T> & result,
287{
289 GISMO_ENSURE(m_dim == 2, "Invalid stress components for a 2D problem");
292 GISMO_ENSURE(m_dim == 3, "Invalid stress type for a 3D problem");
293 GISMO_ENSURE(m_options.getInt("MaterialLaw") == material_law::hooke ||
294 m_options.getInt("MaterialLaw") == material_law::neo_hooke_ln ||
295 m_options.getInt("MaterialLaw") == material_law::saint_venant_kirchhoff ||
296 m_options.getInt("MaterialLaw") == material_law::neo_hooke_quad,
297 "Pressure field not provided! Can't compute stresses with the chosen material law.");
298 result.clear();
299
300 for (size_t p = 0; p < m_pde_ptr->domain().nPatches(); ++p )
301 result.addPiecePointer(new gsCauchyStressFunction<T>(p,comp,m_options,
302 &(m_pde_ptr->domain()),&displacement));
303}
304
305template <class T>
307 const gsMultiPatch<T> & pressure,
308 gsPiecewiseFunction<T> & result,
310{
312 GISMO_ENSURE(m_dim == 2, "Invalid stress components for a 2D problem");
315 GISMO_ENSURE(m_dim == 3, "Invalid stress type for a 3D problem");
316 GISMO_ENSURE(m_options.getInt("MaterialLaw") == material_law::mixed_neo_hooke_ln ||
317 m_options.getInt("MaterialLaw") == material_law::mixed_hooke,
318 "Pressure field is not necessary to compute stresses with the chosen material law.");
319 result.clear();
320
321 for (size_t p = 0; p < m_pde_ptr->domain().nPatches(); ++p )
322 result.addPiecePointer(new gsCauchyStressFunction<T>(p,comp,m_options,
323 &(m_pde_ptr->domain()), &displacement, &pressure));
324}
325
326}// namespace gismo ends
Class containing a set of boundary conditions.
Definition gsBoundaryConditions.h:342
Compute Cauchy stresses for a previously computed/defined displacement field. Can be pushed into gsPi...
Definition gsElasticityFunctions.h:30
virtual void reserve()
a custom reserve function to allocate memory for the sparse matrix
Definition gsElasticityAssembler.hpp:93
virtual void constructSolution(const gsMatrix< T > &solVector, const std::vector< gsMatrix< T > > &fixedDoFs, gsMultiPatch< T > &displacement) const
Construct displacement from computed solution vector and fixed degrees of freedom.
Definition gsElasticityAssembler.hpp:248
gsElasticityAssembler(const gsMultiPatch< T > &patches, const gsMultiBasis< T > &basis, const gsBoundaryConditions< T > &bconditions, const gsFunction< T > &body_force)
Constructor for displacement formulation.
Definition gsElasticityAssembler.hpp:36
virtual void assemble()
Main assemble routine, to be implemented in derived classes.
Definition gsElasticityAssembler.h:62
virtual void refresh()
Refresh routine to set dof-mappers.
Definition gsElasticityAssembler.hpp:120
virtual void constructPressure(const gsMatrix< T > &solVector, const std::vector< gsMatrix< T > > &fixedDoFs, gsMultiPatch< T > &pressure) const
@ brief Construct pressure from computed solution vector
Definition gsElasticityAssembler.hpp:271
virtual void constructCauchyStresses(const gsMultiPatch< T > &displacement, gsPiecewiseFunction< T > &result, stress_components::components component=stress_components::von_mises) const
Construct Cauchy stresses for evaluation or visualization.
Definition gsElasticityAssembler.hpp:284
static gsOptionList defaultOptions()
Returns the list of default options for assembly.
Definition gsElasticityAssembler.hpp:80
A function from a n-dimensional domain to an m-dimensional image.
Definition gsFunction.h:60
virtual short_t targetDim() const
Dimension of the target space.
Definition gsFunctionSet.h:595
A matrix with arbitrary coefficient type and fixed or dynamic size.
Definition gsMatrix.h:41
Holds a set of patch-wise bases and their topology information.
Definition gsMultiBasis.h:37
Container class for a set of geometry patches and their topology, that is, the interface connections ...
Definition gsMultiPatch.h:100
Class which holds a list of parameters/options, and provides easy access to them.
Definition gsOptionList.h:33
void addInt(const std::string &label, const std::string &desc, const index_t &value)
Adds a option named label, with description desc and value value.
Definition gsOptionList.cpp:201
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
void addSwitch(const std::string &label, const std::string &desc, const bool &value)
Adds a option named label, with description desc and value value.
Definition gsOptionList.cpp:235
A function depending on an index i, typically referring to a patch/sub-domain. On each patch a differ...
Definition gsPiecewiseFunction.h:29
void addPiece(const gsFunction< T > &func)
Add a piece.
Definition gsPiecewiseFunction.h:113
void clear()
Clear (delete) all functions.
Definition gsPiecewiseFunction.h:148
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
T at(index_t i) const
Returns the i-th element of the vector.
Definition gsVector.h:177
IMHO, a useless class, but it is necessary to use the gsAssembler class. Contains proper information ...
Provides several simple utility and naming classes.
#define short_t
Definition gsConfig.h:35
#define index_t
Definition gsConfig.h:32
#define GISMO_ENSURE(cond, message)
Definition gsDebug.h:102
Provides linear and nonlinear elasticity systems for 2D plain strain and 3D continua.
Provides isogeometric meshing and modelling routines.
Provides functions to generate structured point data.
Visitor class for the surface load integration.
Visitor class for volumetric integration of the linear elasticity system.
Visitor class for volumetric integration of the mixed linear elasticity system.
Visitor class for volumetric integration of the mixed nonlinear elasticity system.
Element visitor for nonlinear elasticity for 2D plain strain and 3D continua.
The G+Smo namespace, containing all definitions for the library.
index_t checkDisplacement(gsMultiPatch< T > const &domain, gsMultiPatch< T > const &displacement)
Checks whether the deformed configuration is bijective, i.e. det(Jac(geo+disp)) > 0; returns -1 if ye...
Definition gsGeoUtils.hpp:245
@ mixed_neo_hooke_ln
S = lambda/2*(J^2-1)C^-1 + mu*(I-C^-1)
Definition gsBaseUtils.h:135
@ neo_hooke_quad
S = lambda*ln(J)*C^-1 + mu*(I-C^-1)
Definition gsBaseUtils.h:134
@ neo_hooke_ln
S = 2*mu*E + lambda*tr(E)*I.
Definition gsBaseUtils.h:133
@ saint_venant_kirchhoff
sigma = 2*mu*eps + lambda*tr(eps)*I
Definition gsBaseUtils.h:132
@ hooke
sigma = 2*mu*eps + p*I
Definition gsBaseUtils.h:131
components
Definition gsBaseUtils.h:112
@ normal_3D_vector
return all components of the 2D stress tensor as a 2x2 matrix
Definition gsBaseUtils.h:117
@ shear_3D_vector
Definition gsBaseUtils.h:119
@ all_3D_matrix
Definition gsBaseUtils.h:121
@ all_2D_vector
return von Mises stress as a scala
Definition gsBaseUtils.h:114
@ all_2D_matrix
Definition gsBaseUtils.h:116