G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsThinShellAssemblerDWR.h
Go to the documentation of this file.
1
16#pragma once
17
21
22namespace gismo
23{
24
25enum class GoalFunction : short_t
26{
27 Displacement = 0,
28 Stretch = 1,
29 PStrain = 2,
30 PStress = 3,
31 MembraneStrain = 10,
32 MembraneStress = 11,
33 MembraneForce = 12,
34 FlexuralStrain = 20,
35 FlexuralStress = 21,
36 FlexuralMoment = 22,
37 Modal = 100,
38 Buckling = 200,
39};
40
41#ifdef GISMO_WITH_PYBIND11
42
46 void pybind11_enum_GoalFunction(pybind11::module &m);
47
48#endif // GISMO_WITH_PYBIND11
49
50
51template<class T> class gsThinShellAssemblerDWRBase;
52
53template <short_t d, class T, bool bending>
54class gsThinShellAssemblerDWR : public gsThinShellAssembler<d,T,bending>,
55 public gsThinShellAssemblerDWRBase<T>
56
57{
58public:
59 typedef gsThinShellAssembler<d,T,bending> Base;
60 typedef typename gsThinShellAssemblerDWRBase<T>::bContainer bContainer;
61
62 virtual ~gsThinShellAssemblerDWR()
63 {
64 delete m_assemblerL;
65 delete m_assemblerH;
66 };
67
68 // empty constructor
69 gsThinShellAssemblerDWR() {};
70
71 gsThinShellAssemblerDWR(
72 const gsMultiPatch<T> & patches,
73 const gsMultiBasis<T> & basisL,
74 const gsMultiBasis<T> & basisH,
75 const gsBoundaryConditions<T> & bconditions,
76 const gsFunction<T> & surface_force,
77 typename gsMaterialMatrixBase<T>::uPtr & materialmatrix
78 );
79
80 gsThinShellAssemblerDWR(
81 const gsMultiPatch<T> & patches,
82 const gsMultiBasis<T> & basisL,
83 const gsMultiBasis<T> & basisH,
84 const gsBoundaryConditions<T> & bconditions,
85 const gsFunction<T> & surface_force,
86 gsMaterialMatrixBase<T> * materialmatrix
87 );
88
90 gsOptionList & optionsL() {return m_assemblerL->options();}
92 gsOptionList & optionsH() {return m_assemblerH->options();}
93
95 gsExprAssembler<T> assemblerL() {return m_assemblerL->assembler(); }
97 gsExprAssembler<T> assemblerH() {return m_assemblerH->assembler(); }
98
100 void setOptions(gsOptionList & options) { m_assemblerL->setOptions(options); m_assemblerH->setOptions(options); }
101
103 void setPointLoads(const gsPointLoads<T> & pLoads) { m_pLoads = pLoads; m_assemblerL->setPointLoads(pLoads); m_assemblerH->setPointLoads(pLoads); }
104
106 void setFoundation(const gsFunction<T> & foundation) { m_assemblerL->setFoundation(foundation); m_assemblerH->setFoundation(foundation); }
107
109 void setPressure(const gsFunction<T> & pressure) { m_assemblerL->setPressure(pressure); m_assemblerH->setPressure(pressure); }
110
112 void updateBCs(const gsBoundaryConditions<T> & bconditions) { m_assemblerL->updateBCs(bconditions); m_assemblerH->updateBCs(bconditions); }
113
115 void setBasisL(const gsMultiBasis<T> & basis) { m_assemblerL->setBasis(basis); }
116 void setBasisH(const gsMultiBasis<T> & basis) { m_assemblerH->setBasis(basis); }
117
119 void setUndeformed(const gsMultiPatch<T> & patches) { m_patches = patches; m_assemblerL->setUndeformed(patches); m_assemblerH->setUndeformed(patches); }
120
122 void homogenizeDirichlet() { m_assemblerL->homogenizeDirichlet(); m_assemblerH->homogenizeDirichlet(); }
123
125 index_t numDofsL() const { return m_assemblerL->numDofs(); }
126 index_t numDofsH() const { return m_assemblerH->numDofs(); }
127
128 void setSpaceBasisL(const gsMultiPatch<T> & spaceBasis) { m_assemblerL->setSpaceBasis(spaceBasis); }
129 void setSpaceBasisH(const gsMultiPatch<T> & spaceBasis) { m_assemblerH->setSpaceBasis(spaceBasis); }
130
131
133 //--------------------- SYSTEM ASSEMBLY ----------------------------------//
135 ThinShellAssemblerStatus assembleL();
136 ThinShellAssemblerStatus assembleH();
137
138 ThinShellAssemblerStatus assembleMassL(bool lumped = false)
139 { return _assembleMass(m_assemblerL, m_massL, lumped); }
140 ThinShellAssemblerStatus assembleMassH(bool lumped = false)
141 { return _assembleMass(m_assemblerH, m_massH, lumped); }
142 ThinShellAssemblerStatus assembleMatrixL()
143 { return _assembleMatrix(m_assemblerL, m_matrixL); }
144 ThinShellAssemblerStatus assembleMatrixH()
145 { return _assembleMatrix(m_assemblerH, m_matrixH); }
146 ThinShellAssemblerStatus assembleMatrixL(const gsMultiPatch<T> & deformed)
147 { return _assembleMatrix(m_assemblerL,deformed, m_matrixL); }
148 ThinShellAssemblerStatus assembleMatrixH(const gsMultiPatch<T> & deformed)
149 { return _assembleMatrix(m_assemblerH,deformed, m_matrixH); }
150
151 ThinShellAssemblerStatus assemblePrimalL()
152 { return _assemblePrimal(m_assemblerL, m_pL); }
153 ThinShellAssemblerStatus assemblePrimalH()
154 { return _assemblePrimal(m_assemblerH, m_pH); }
155 ThinShellAssemblerStatus assemblePrimalL(const gsMultiPatch<T> & deformed)
156 { return _assemblePrimal(m_assemblerL,deformed, m_pL); }
157 ThinShellAssemblerStatus assemblePrimalH(const gsMultiPatch<T> & deformed)
158 { return _assemblePrimal(m_assemblerH,deformed, m_pH); }
159
160 ThinShellAssemblerStatus assembleDualL(const gsMultiPatch<T> & primal)
161 { return _assembleDual(m_assemblerL,primal, m_dL); }
162 ThinShellAssemblerStatus assembleDualH(const gsMultiPatch<T> & primal)
163 { return _assembleDual(m_assemblerH,primal, m_dH); }
164
165 ThinShellAssemblerStatus assembleDualL(const bContainer & bnds, const gsMultiPatch<T> & primal)
166 { return _assembleDual(bnds,m_assemblerL,primal, m_dL); }
167 ThinShellAssemblerStatus assembleDualH(const bContainer & bnds, const gsMultiPatch<T> & primal)
168 { return _assembleDual(bnds,m_assemblerH,primal, m_dH); }
169
170 ThinShellAssemblerStatus assembleDualL(const gsMatrix<T> & points, const gsMultiPatch<T> & primal)
171 { return _assembleDual(points,m_assemblerL,primal, m_dL); }
172 ThinShellAssemblerStatus assembleDualH(const gsMatrix<T> & points, const gsMultiPatch<T> & primal)
173 { return _assembleDual(points,m_assemblerH,primal, m_dH); }
174
175 ThinShellAssemblerStatus assembleDualL(const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed)
176 { return _assembleDual(m_assemblerL,primal,deformed, m_dL); }
177 ThinShellAssemblerStatus assembleDualH(const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed)
178 { return _assembleDual(m_assemblerH,primal,deformed, m_dH); }
179
180 ThinShellAssemblerStatus assembleDualL(const bContainer & bnds, const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed)
181 { return _assembleDual(bnds,m_assemblerL,primal,deformed, m_dL); }
182 ThinShellAssemblerStatus assembleDualH(const bContainer & bnds, const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed)
183 { return _assembleDual(bnds,m_assemblerH,primal,deformed, m_dH); }
184
185 ThinShellAssemblerStatus assembleDualL(const gsMatrix<T> & points, const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed)
186 { return _assembleDual(points,m_assemblerL,primal,deformed, m_dL); }
187 ThinShellAssemblerStatus assembleDualH(const gsMatrix<T> & points, const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed)
188 { return _assembleDual(points,m_assemblerH,primal,deformed, m_dH); }
189
190 const gsSparseMatrix<T> & matrixL() const { return m_matrixL; }
191 const gsSparseMatrix<T> & matrixH() const { return m_matrixH; }
192
193 const gsSparseMatrix<T> & massL() const { return m_massL; }
194 const gsSparseMatrix<T> & massH() const { return m_massH; }
195
196 const gsMatrix<T> primalL() const { return m_pL; }
197 const gsMatrix<T> primalH() const { return m_pH; }
198
199 const gsMatrix<T> dualL() const { return m_dL; }
200 const gsMatrix<T> dualH() const { return m_dH; }
201
202 void updateMultiPatchL(const gsMatrix<T> & solVector, gsMultiPatch<T> & result);
203 void updateMultiPatchH(const gsMatrix<T> & solVector, gsMultiPatch<T> & result);
204
205 void constructMultiPatchL(const gsMatrix<T> & solVector, gsMultiPatch<T> & result);
206 void constructMultiPatchH(const gsMatrix<T> & solVector, gsMultiPatch<T> & result);
207
208 void constructSolutionL(const gsMatrix<T> & solVector, gsMultiPatch<T> & deformed);
209 void constructSolutionH(const gsMatrix<T> & solVector, gsMultiPatch<T> & deformed);
210 gsMultiPatch<T> constructSolutionL(const gsMatrix<T> & solVector);
211 gsMultiPatch<T> constructSolutionH(const gsMatrix<T> & solVector);
212 gsMultiPatch<T> constructDisplacementL(const gsMatrix<T> & solVector);
213 gsMultiPatch<T> constructDisplacementH(const gsMatrix<T> & solVector);
214 gsVector<T> constructSolutionVectorL(const gsMultiPatch<T> & deformed);
215 gsVector<T> constructSolutionVectorH(const gsMultiPatch<T> & deformed);
216
217 gsMatrix<T> projectL2_L(const gsFunction<T> &fun);
218 gsMatrix<T> projectL2_H(const gsFunction<T> &fun);
219
220 T deformationNorm(const gsMultiPatch<T> & deformed, const gsMultiPatch<T> & original)
221 { return m_assemblerL->deformationNorm(deformed,original); }
222
223 // Linear elasticity
224 T computeError(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, bool withLoads=false, std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
225 std::vector<T> computeErrorElements(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, bool withLoads=false, std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
226 std::vector<T> computeErrorDofs(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, bool withLoads=false, std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
227
228 // Nonlinear elasticity
229 T computeError(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, const gsMultiPatch<T> & deformed, bool withLoads=false, std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
230 std::vector<T> computeErrorElements(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, const gsMultiPatch<T> & deformed, bool withLoads=false, std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
231 std::vector<T> computeErrorDofs(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, const gsMultiPatch<T> & deformed, bool withLoads=false, std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
232
233 // Linear elasticity
234 T computeSquaredError(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, bool withLoads=false, std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
235 std::vector<T> computeSquaredErrorElements(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, bool withLoads=false, std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
236 std::vector<T> computeSquaredErrorDofs(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, bool withLoads=false, std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
237
238 // Nonlinear elasticity
239 T computeSquaredError(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, const gsMultiPatch<T> & deformed, bool withLoads=false, std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
240 std::vector<T> computeSquaredErrorElements(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, const gsMultiPatch<T> & deformed, bool withLoads=false, std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
241 std::vector<T> computeSquaredErrorDofs(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, const gsMultiPatch<T> & deformed, bool withLoads=false, std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
242
243
244 // Eigenvalues
245 T computeErrorEig(const T evPrimalL, const T evDualL, const T evDualH,
246 const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH,
247 const gsMultiPatch<T> & primal,std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
248 std::vector<T> computeErrorEigElements(const T evPrimalL, const T evDualL, const T evDualH,
249 const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH,
250 const gsMultiPatch<T> & primal,std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
251 std::vector<T> computeErrorEigDofs(const T evPrimalL, const T evDualL, const T evDualH,
252 const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH,
253 const gsMultiPatch<T> & primal,std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
254
255 T computeErrorEig(const T evPrimalL, const T evDualL, const T evDualH,
256 const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH,
257 const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed,std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
258 std::vector<T> computeErrorEigElements(const T evPrimalL, const T evDualL, const T evDualH,
259 const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH,
260 const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed,std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
261 std::vector<T> computeErrorEigDofs(const T evPrimalL, const T evDualL, const T evDualH,
262 const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH,
263 const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed,std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
264
265 T computeGoal(const gsMultiPatch<T> & deformed);
266 T computeGoal(const bContainer & bnds, const gsMultiPatch<T> & deformed);
267 T computeGoal(const gsMatrix<T> & points, const gsMultiPatch<T> & deformed);
268
269 T matrixNorm(const gsMultiPatch<T> &dualL, const gsMultiPatch<T> &dualH) const;
270 T matrixNorm(const gsMultiPatch<T> &dualL, const gsMultiPatch<T> &dualH, const gsMultiPatch<T> &deformed) const;
271
272 void setGoal(enum GoalFunction GF, short_t component = 9) { m_goalFunction = GF; m_component = component; }
273
275 void constructStress(const gsFunctionSet<T> & deformed, gsPiecewiseFunction<T> & result, stress_type::type type)
276 {
277 m_assemblerL->constructStress(deformed,result,type);
278 }
279
280 T error() const { return m_error; }
281 std::vector<T> errors() const { return m_errors; }
282 std::vector<T> absErrors() const
283 {
284 std::vector<T> result = m_errors;
285 for (typename std::vector<T>::iterator it = result.begin(); it!=result.end(); it++)
286 *it = std::abs(*it);
287 return result;
288 }
289
290 gsDofMapper getMapperL() { return m_assemblerL->getMapper(); };
291 gsDofMapper getMapperH() { return m_assemblerH->getMapper(); };
292
293protected:
294
295 ThinShellAssemblerStatus _assembleMass(gsThinShellAssemblerBase<T> * assembler, gsSparseMatrix<T> & result, bool lumped = false);
296
297 ThinShellAssemblerStatus _assemble(gsThinShellAssemblerBase<T> * assembler, std::pair<gsSparseMatrix<T>,gsVector<T>> & result);
298
299 ThinShellAssemblerStatus _assembleMatrix(gsThinShellAssemblerBase<T> * assembler, gsSparseMatrix<T> & result);
300 ThinShellAssemblerStatus _assembleMatrix(gsThinShellAssemblerBase<T> * assembler, const gsMultiPatch<T> & deformed, gsSparseMatrix<T> & result);
301 ThinShellAssemblerStatus _assemblePrimal(gsThinShellAssemblerBase<T> * assembler, gsVector<T> & result);
302 ThinShellAssemblerStatus _assemblePrimal(gsThinShellAssemblerBase<T> * assembler, const gsMultiPatch<T> & deformed, gsVector<T> & result);
303
314 ThinShellAssemblerStatus _assembleDual(gsThinShellAssemblerBase<T> * assembler, const gsMultiPatch<T> & primal, gsVector<T> & result)
315 {
316 gsMultiPatch<T> deformed = m_patches;
317 for ( size_t k =0; k!=primal.nPatches(); ++k) // Deform the geometry
318 deformed.patch(k).coefs() += primal.patch(k).coefs(); // Gdef points to mp_def, therefore updated
319
320 return _assembleDual(assembler,primal,deformed,result);
321 }
322 ThinShellAssemblerStatus _assembleDual(gsThinShellAssemblerBase<T> * assembler, const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed, gsVector<T> & result);
323
335 ThinShellAssemblerStatus _assembleDual(const bContainer & bnds, gsThinShellAssemblerBase<T> * assembler, const gsMultiPatch<T> & primal, gsVector<T> & result)
336 {
337 gsMultiPatch<T> deformed = m_patches;
338 for ( size_t k =0; k!=primal.nPatches(); ++k) // Deform the geometry
339 deformed.patch(k).coefs() += primal.patch(k).coefs(); // Gdef points to mp_def, therefore updated
340
341 return _assembleDual(bnds,assembler,primal,deformed,result);
342 }
343 ThinShellAssemblerStatus _assembleDual(const bContainer & bnds, gsThinShellAssemblerBase<T> * assembler, const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed, gsVector<T> & result);
344
356 ThinShellAssemblerStatus _assembleDual(const gsMatrix<T> & points, gsThinShellAssemblerBase<T> * assembler, const gsMultiPatch<T> & primal, gsVector<T> & result)
357 {
358 gsMultiPatch<T> deformed = m_patches;
359 for ( size_t k =0; k!=primal.nPatches(); ++k) // Deform the geometry
360 deformed.patch(k).coefs() += primal.patch(k).coefs(); // Gdef points to mp_def, therefore updated
361
362 return _assembleDual(points,assembler,primal,deformed,result);
363 }
364 ThinShellAssemblerStatus _assembleDual(const gsMatrix<T> & points, gsThinShellAssemblerBase<T> * assembler, const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed, gsVector<T> & result);
365
366 template<int _elWise>
367 void computeError_impl(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, bool withLoads,
368 std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
369
370 template<int _d, bool _bending, int _elWise>
371 typename std::enable_if<(_d==3) && _bending, void>::type
372 computeError_impl(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, const gsMultiPatch<T> & deformed, bool withLoads,
373 std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
374
375 template<int _d, bool _bending, int _elWise>
376 typename std::enable_if<!(_d==3 && _bending), void>::type
377 computeError_impl(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, const gsMultiPatch<T> & deformed, bool withLoads,
378 std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
379
380 template<int _elWise>
381 void computeSquaredError_impl(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, bool withLoads,
382 std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
383
384 template<int _d, bool _bending, int _elWise>
385 typename std::enable_if<(_d==3) && _bending, void>::type
386 computeSquaredError_impl(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, const gsMultiPatch<T> & deformed, bool withLoads,
387 std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
388
389 template<int _d, bool _bending, int _elWise>
390 typename std::enable_if<!(_d==3 && _bending), void>::type
391 computeSquaredError_impl(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, const gsMultiPatch<T> & deformed, bool withLoads,
392 std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
393
394
395 template<int _elWise>
396 void computeErrorEig_impl(const T evPrimalL, const T evDualL, const T evDualH,
397 const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH,
398 const gsMultiPatch<T> & primal,
399 std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
400
401 template<int _elWise>
402 void computeErrorEig_impl(const T evPrimalL, const T evDualL, const T evDualH,
403 const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH,
404 const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed,
405 std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
406
407 void _applyLoadsToElWiseError(const gsMultiPatch<T> &dualL, const gsMultiPatch<T> &dualH, std::vector<T> & errors) const;
408 void _applyLoadsToError(const gsMultiPatch<T> &dualL, const gsMultiPatch<T> &dualH, T & error) const;
409
410
411protected:
412 mutable gsThinShellAssemblerBase<T> * m_assemblerL;
413 mutable gsThinShellAssemblerBase<T> * m_assemblerH;
414
415 gsVector<T> m_pL, m_pH, m_dL, m_dH;
416 gsSparseMatrix<T> m_matrixL;
417 gsSparseMatrix<T> m_matrixH;
418 gsSparseMatrix<T> m_massL;
419 gsSparseMatrix<T> m_massH;
420
421 gsPointLoads<T> m_pLoads;
422
423 T m_error;
424 std::vector<T> m_errors;
425
426 typedef gsExprAssembler<>::geometryMap geometryMap;
427 typedef gsExprAssembler<>::space space;
428 typedef gsExprAssembler<>::solution solution;
429
430 gsBoundaryConditions<T> m_bcs;
431 // using Base::m_bcs;
432 // using Base::m_basis;
433
434 using Base::m_patches;
435 // using Base::m_defpatches;
436 using Base::m_materialMatrices;
437 using Base::m_forceFun;
438 using Base::m_options;
439 using Base::m_foundFun;
440 using Base::m_foundInd;
441 using Base::m_pressFun;
442 using Base::m_pressInd;
443
444 enum GoalFunction m_goalFunction;
445 short_t m_component;
446};
447
448#ifdef GISMO_WITH_PYBIND11
449
453 void pybind11_init_gsThinShellAssemblerDWR2(pybind11::module &m);
454 void pybind11_init_gsThinShellAssemblerDWR3(pybind11::module &m);
455 void pybind11_init_gsThinShellAssemblerDWR3nb(pybind11::module &m);
456
457#endif // GISMO_WITH_PYBIND11
458
466template <class T>
467class gsThinShellAssemblerDWRBase //: public virtual gsThinShellAssemblerBase<T>
468{
469public:
470 typedef gsBoxTopology::bContainer bContainer;
471
474
475 virtual gsOptionList & optionsL() =0;
476 virtual gsOptionList & optionsH() =0;
477
478 virtual gsExprAssembler<T> assemblerL() =0;
479 virtual gsExprAssembler<T> assemblerH() =0;
480
481 virtual void setOptions(gsOptionList & options) =0;
482 virtual void setPointLoads(const gsPointLoads<T> & pLoads) =0;
483
484 virtual void setFoundation(const gsFunction<T> & foundation) =0;
485 virtual void setPressure(const gsFunction<T> & pressure) =0;
486
487 virtual void updateBCs(const gsBoundaryConditions<T> & bconditions) =0;
488
489 virtual void setBasisL(const gsMultiBasis<T> & basis) =0;
490 virtual void setBasisH(const gsMultiBasis<T> & basis) =0;
491
492 virtual void setUndeformed(const gsMultiPatch<T> & patches) =0;
493
494 virtual void homogenizeDirichlet() =0;
495
496 virtual index_t numDofsL() const =0;
497 virtual index_t numDofsH() const =0;
498
499 virtual void setSpaceBasisL(const gsMultiPatch<T> & spaceBasis) =0;
500 virtual void setSpaceBasisH(const gsMultiPatch<T> & spaceBasis) =0;
501
502 virtual ThinShellAssemblerStatus assembleL() =0;
503 virtual ThinShellAssemblerStatus assembleH() =0;
504
505 virtual ThinShellAssemblerStatus assembleMassL(bool lumped = false) =0;
506 virtual ThinShellAssemblerStatus assembleMassH(bool lumped = false) =0;
507 virtual ThinShellAssemblerStatus assembleMatrixL() =0;
508 virtual ThinShellAssemblerStatus assembleMatrixH() =0;
509 virtual ThinShellAssemblerStatus assembleMatrixL(const gsMultiPatch<T> & deformed) =0;
510 virtual ThinShellAssemblerStatus assembleMatrixH(const gsMultiPatch<T> & deformed) =0;
511
512 virtual ThinShellAssemblerStatus assemblePrimalL() =0;
513 virtual ThinShellAssemblerStatus assemblePrimalH() =0;
514 virtual ThinShellAssemblerStatus assemblePrimalL(const gsMultiPatch<T> & deformed) =0;
515 virtual ThinShellAssemblerStatus assemblePrimalH(const gsMultiPatch<T> & deformed) =0;
516
517
518 virtual ThinShellAssemblerStatus assembleDualL(const gsMultiPatch<T> & primal) =0;
519 virtual ThinShellAssemblerStatus assembleDualH(const gsMultiPatch<T> & primal) =0;
520
521 virtual ThinShellAssemblerStatus assembleDualL(const bContainer & bnds, const gsMultiPatch<T> & primal) =0;
522 virtual ThinShellAssemblerStatus assembleDualH(const bContainer & bnds, const gsMultiPatch<T> & primal) =0;
523
524 virtual ThinShellAssemblerStatus assembleDualL(const gsMatrix<T> & points, const gsMultiPatch<T> & primal) =0;
525 virtual ThinShellAssemblerStatus assembleDualH(const gsMatrix<T> & points, const gsMultiPatch<T> & primal) =0;
526
527 virtual ThinShellAssemblerStatus assembleDualL(const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed) =0;
528 virtual ThinShellAssemblerStatus assembleDualH(const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed) =0;
529
530 virtual ThinShellAssemblerStatus assembleDualL(const bContainer & bnds, const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed) =0;
531 virtual ThinShellAssemblerStatus assembleDualH(const bContainer & bnds, const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed) =0;
532
533 virtual ThinShellAssemblerStatus assembleDualL(const gsMatrix<T> & points, const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed) =0;
534 virtual ThinShellAssemblerStatus assembleDualH(const gsMatrix<T> & points, const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed) =0;
535
536 virtual const gsSparseMatrix<T> & matrixL() const =0;
537 virtual const gsSparseMatrix<T> & matrixH() const =0;
538
539 virtual const gsSparseMatrix<T> & massL() const =0;
540 virtual const gsSparseMatrix<T> & massH() const =0;
541
542 virtual const gsMatrix<T> primalL() const =0;
543 virtual const gsMatrix<T> primalH() const =0;
544
545 virtual const gsMatrix<T> dualL() const =0;
546 virtual const gsMatrix<T> dualH() const =0;
547
548 virtual void updateMultiPatchL(const gsMatrix<T> & solVector, gsMultiPatch<T> & result) =0;
549 virtual void updateMultiPatchH(const gsMatrix<T> & solVector, gsMultiPatch<T> & result) =0;
550
551 virtual void constructMultiPatchL(const gsMatrix<T> & solVector, gsMultiPatch<T> & result) =0;
552 virtual void constructMultiPatchH(const gsMatrix<T> & solVector, gsMultiPatch<T> & result) =0;
553
554 virtual void constructSolutionL(const gsMatrix<T> & solVector, gsMultiPatch<T> & deformed) =0;
555 virtual void constructSolutionH(const gsMatrix<T> & solVector, gsMultiPatch<T> & deformed) =0;
556 virtual gsMultiPatch<T> constructSolutionL(const gsMatrix<T> & solVector) =0;
557 virtual gsMultiPatch<T> constructSolutionH(const gsMatrix<T> & solVector) =0;
558 virtual gsMultiPatch<T> constructDisplacementL(const gsMatrix<T> & solVector) =0;
559 virtual gsMultiPatch<T> constructDisplacementH(const gsMatrix<T> & solVector) =0;
560 virtual gsVector<T> constructSolutionVectorL(const gsMultiPatch<T> & deformed) =0;
561 virtual gsVector<T> constructSolutionVectorH(const gsMultiPatch<T> & deformed) =0;
562
563 virtual gsMatrix<T> projectL2_L(const gsFunction<T> &fun) =0;
564 virtual gsMatrix<T> projectL2_H(const gsFunction<T> &fun) =0;
565
566 virtual T deformationNorm(const gsMultiPatch<T> & deformed, const gsMultiPatch<T> & original) =0;
567
568 // Linear elasticity ;
569 virtual T computeError(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, bool withLoads=false, std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false) = 0;
570 virtual std::vector<T> computeErrorElements(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, bool withLoads=false, std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false) = 0;
571 virtual std::vector<T> computeErrorDofs(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, bool withLoads=false, std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false) = 0;
572
573 // Nonlinear elasticity
574 virtual T computeError(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, const gsMultiPatch<T> & deformed, bool withLoads=false,std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false) = 0;
575 virtual std::vector<T> computeErrorElements(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, const gsMultiPatch<T> & deformed, bool withLoads=false, std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false) = 0;
576 virtual std::vector<T> computeErrorDofs(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, const gsMultiPatch<T> & deformed, bool withLoads=false, std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false) = 0;
577
578 // Linear elasticity ;
579 virtual T computeSquaredError(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, bool withLoads=false, std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false) = 0;
580 virtual std::vector<T> computeSquaredErrorElements(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, bool withLoads=false, std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false) = 0;
581 virtual std::vector<T> computeSquaredErrorDofs(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, bool withLoads=false, std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false) = 0;
582
583 // Nonlinear elasticity
584 virtual T computeSquaredError(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, const gsMultiPatch<T> & deformed, bool withLoads=false,std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false) = 0;
585 virtual std::vector<T> computeSquaredErrorElements(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, const gsMultiPatch<T> & deformed, bool withLoads=false, std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false) = 0;
586 virtual std::vector<T> computeSquaredErrorDofs(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, const gsMultiPatch<T> & deformed, bool withLoads=false, std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false) = 0;
587
588
589 // Eigenvalues
590 virtual T computeErrorEig(const T evPrimalL, const T evDualL, const T evDualH,
591 const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH,
592 const gsMultiPatch<T> & primal,std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false) = 0;
593 virtual std::vector<T> computeErrorEigElements(const T evPrimalL, const T evDualL, const T evDualH,
594 const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH,
595 const gsMultiPatch<T> & primal,std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false) = 0;
596 virtual std::vector<T> computeErrorEigDofs(const T evPrimalL, const T evDualL, const T evDualH,
597 const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH,
598 const gsMultiPatch<T> & primal,std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false) = 0;
599
600 virtual T computeErrorEig(const T evPrimalL, const T evDualL, const T evDualH,
601 const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH,
602 const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed,std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false) = 0;
603 virtual std::vector<T> computeErrorEigElements(const T evPrimalL, const T evDualL, const T evDualH,
604 const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH,
605 const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed,std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false) = 0;
606 virtual std::vector<T> computeErrorEigDofs(const T evPrimalL, const T evDualL, const T evDualH,
607 const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH,
608 const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed,std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false) = 0;
609
610 virtual T computeGoal(const gsMultiPatch<T> & deformed) =0;
611 virtual T computeGoal(const bContainer & bnds, const gsMultiPatch<T> & deformed) =0;
612 virtual T computeGoal(const gsMatrix<T> & points, const gsMultiPatch<T> & deformed) =0;
613
614 virtual T matrixNorm(const gsMultiPatch<T> &dualL, const gsMultiPatch<T> &dualH) const = 0;
615 virtual T matrixNorm(const gsMultiPatch<T> &dualL, const gsMultiPatch<T> &dualH, const gsMultiPatch<T> &deformed) const = 0;
616
617 virtual void setGoal(enum GoalFunction GF, short_t component = 9) = 0;
618
620 virtual void constructStress(const gsFunctionSet<T> & deformed,
621 gsPiecewiseFunction<T> & result,
622 stress_type::type type) = 0;
623
624 virtual T error() const =0;
625 virtual std::vector<T> errors() const =0;
626 virtual std::vector<T> absErrors() const =0;
627
628 virtual gsDofMapper getMapperL()=0;
629 virtual gsDofMapper getMapperH()=0;
630
631
632};
633
634#ifdef GISMO_WITH_PYBIND11
635
639 void pybind11_init_gsThinShellAssemblerDWRBase(pybind11::module &m);
640
641#endif // GISMO_WITH_PYBIND11
642
643
644} // namespace gismo
645
646
649
650
651#ifndef GISMO_BUILD_LIB
652#include GISMO_HPP_HEADER(gsThinShellAssemblerDWR.hpp)
653#endif
Class containing a set of boundary conditions.
Definition gsBoundaryConditions.h:342
Maintains a mapping from patch-local dofs to global dof indices and allows the elimination of individ...
Definition gsDofMapper.h:69
Definition gsExprAssembler.h:33
gsExprHelper< T >::geometryMap geometryMap
Geometry map type.
Definition gsExprAssembler.h:65
expr::gsFeSolution< T > solution
Solution type.
Definition gsExprAssembler.h:68
gsExprHelper< T >::space space
Space type.
Definition gsExprAssembler.h:67
Interface for the set of functions defined on a domain (the total number of functions in the set equa...
Definition gsFunctionSet.h:219
A function from a n-dimensional domain to an m-dimensional image.
Definition gsFunction.h:60
memory::unique_ptr< gsMaterialMatrixBase > uPtr
Unique pointer for gsGeometry.
Definition gsMaterialMatrixBase.h:44
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
A function depending on an index i, typically referring to a patch/sub-domain. On each patch a differ...
Definition gsPiecewiseFunction.h:29
Class containing a set of points on a multi-patch isogeometric domain, together with boundary conditi...
Definition gsPointLoads.h:65
Sparse matrix class, based on gsEigen::SparseMatrix.
Definition gsSparseMatrix.h:139
Base class for the gsThinShellAssembler.
Definition gsThinShellAssemblerDWR.h:468
virtual ~gsThinShellAssemblerDWRBase()
Default empty constructor.
Definition gsThinShellAssemblerDWR.h:473
virtual void constructStress(const gsFunctionSet< T > &deformed, gsPiecewiseFunction< T > &result, stress_type::type type)=0
See gsThinShellAssemblerBase for details.
gsExprAssembler< T > assembler()
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.h:158
const gsMultiBasis< T > & basis() const
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.h:360
void constructStress(const gsFunctionSet< T > &original, const gsFunctionSet< T > &deformed, gsPiecewiseFunction< T > &result, stress_type::type type)
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.hpp:2803
gsOptionList & options()
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.h:155
A vector with arbitrary coefficient type and fixed or dynamic size.
Definition gsVector.h:37
#define short_t
Definition gsConfig.h:35
#define index_t
Definition gsConfig.h:32
Provides a base class for material matrices.
Provides linear and nonlinear assemblers for thin shells.
Provides evaluation function for stresses.
The G+Smo namespace, containing all definitions for the library.
ThinShellAssemblerStatus
Definition gsThinShellAssembler.h:54
type
Definition gsThinShellFunctions.h:39