G+Smo  24.08.0
Geometry + Simulation Modules
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
gsThinShellAssemblerDWR.h
Go to the documentation of this file.
1 
16 #pragma once
17 
21 
22 namespace gismo
23 {
24 
25 enum 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 
51 template<class T> class gsThinShellAssemblerDWRBase;
52 
53 template <short_t d, class T, bool bending>
54 class gsThinShellAssemblerDWR : public gsThinShellAssembler<d,T,bending>,
56 
57 {
58 public:
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 
72  gsThinShellAssemblerDWR(
73  const gsMultiPatch<T> & patches,
74  const gsMultiBasis<T> & basisL,
75  const gsMultiBasis<T> & basisH,
76  const gsBoundaryConditions<T> & bconditions,
77  const gsFunction<T> & surface_force,
78  gsMaterialMatrixBase<T> * materialmatrix
79  );
80 
82  gsOptionList & optionsL() {return m_assemblerL->options();}
84  gsOptionList & optionsH() {return m_assemblerH->options();}
85 
87  gsExprAssembler<T> assemblerL() {return m_assemblerL->assembler(); }
89  gsExprAssembler<T> assemblerH() {return m_assemblerH->assembler(); }
90 
92  void setOptions(gsOptionList & options) { m_assemblerL->setOptions(options); m_assemblerH->setOptions(options); }
93 
95  void setPointLoads(const gsPointLoads<T> & pLoads) { m_pLoads = pLoads; m_assemblerL->setPointLoads(pLoads); m_assemblerH->setPointLoads(pLoads); }
96 
98  void setFoundation(const gsFunction<T> & foundation) { m_assemblerL->setFoundation(foundation); m_assemblerH->setFoundation(foundation); }
99 
101  void setPressure(const gsFunction<T> & pressure) { m_assemblerL->setPressure(pressure); m_assemblerH->setPressure(pressure); }
102 
104  void updateBCs(const gsBoundaryConditions<T> & bconditions) { m_assemblerL->updateBCs(bconditions); m_assemblerH->updateBCs(bconditions); }
105 
107  void setBasisL(const gsMultiBasis<T> & basis) { m_assemblerL->setBasis(basis); }
108  void setBasisH(const gsMultiBasis<T> & basis) { m_assemblerH->setBasis(basis); }
109 
111  void setUndeformed(const gsMultiPatch<T> & patches) { m_patches = patches; m_assemblerL->setUndeformed(patches); m_assemblerH->setUndeformed(patches); }
112 
114  void homogenizeDirichlet() { m_assemblerL->homogenizeDirichlet(); m_assemblerH->homogenizeDirichlet(); }
115 
117  index_t numDofsL() const { return m_assemblerL->numDofs(); }
118  index_t numDofsH() const { return m_assemblerH->numDofs(); }
119 
120  void setSpaceBasisL(const gsMultiPatch<T> & spaceBasis) { m_assemblerL->setSpaceBasis(spaceBasis); }
121  void setSpaceBasisH(const gsMultiPatch<T> & spaceBasis) { m_assemblerH->setSpaceBasis(spaceBasis); }
122 
123 
125  //--------------------- SYSTEM ASSEMBLY ----------------------------------//
127  ThinShellAssemblerStatus assembleL();
128  ThinShellAssemblerStatus assembleH();
129 
130  ThinShellAssemblerStatus assembleMassL(bool lumped = false)
131  { return _assembleMass(m_assemblerL, m_massL, lumped); }
132  ThinShellAssemblerStatus assembleMassH(bool lumped = false)
133  { return _assembleMass(m_assemblerH, m_massH, lumped); }
134  ThinShellAssemblerStatus assembleMatrixL()
135  { return _assembleMatrix(m_assemblerL, m_matrixL); }
136  ThinShellAssemblerStatus assembleMatrixH()
137  { return _assembleMatrix(m_assemblerH, m_matrixH); }
138  ThinShellAssemblerStatus assembleMatrixL(const gsMultiPatch<T> & deformed)
139  { return _assembleMatrix(m_assemblerL,deformed, m_matrixL); }
140  ThinShellAssemblerStatus assembleMatrixH(const gsMultiPatch<T> & deformed)
141  { return _assembleMatrix(m_assemblerH,deformed, m_matrixH); }
142 
143  ThinShellAssemblerStatus assemblePrimalL()
144  { return _assemblePrimal(m_assemblerL, m_pL); }
145  ThinShellAssemblerStatus assemblePrimalH()
146  { return _assemblePrimal(m_assemblerH, m_pH); }
147  ThinShellAssemblerStatus assemblePrimalL(const gsMultiPatch<T> & deformed)
148  { return _assemblePrimal(m_assemblerL,deformed, m_pL); }
149  ThinShellAssemblerStatus assemblePrimalH(const gsMultiPatch<T> & deformed)
150  { return _assemblePrimal(m_assemblerH,deformed, m_pH); }
151 
152  ThinShellAssemblerStatus assembleDualL(const gsMultiPatch<T> & primal)
153  { return _assembleDual(m_assemblerL,primal, m_dL); }
154  ThinShellAssemblerStatus assembleDualH(const gsMultiPatch<T> & primal)
155  { return _assembleDual(m_assemblerH,primal, m_dH); }
156 
157  ThinShellAssemblerStatus assembleDualL(const bContainer & bnds, const gsMultiPatch<T> & primal)
158  { return _assembleDual(bnds,m_assemblerL,primal, m_dL); }
159  ThinShellAssemblerStatus assembleDualH(const bContainer & bnds, const gsMultiPatch<T> & primal)
160  { return _assembleDual(bnds,m_assemblerH,primal, m_dH); }
161 
162  ThinShellAssemblerStatus assembleDualL(const gsMatrix<T> & points, const gsMultiPatch<T> & primal)
163  { return _assembleDual(points,m_assemblerL,primal, m_dL); }
164  ThinShellAssemblerStatus assembleDualH(const gsMatrix<T> & points, const gsMultiPatch<T> & primal)
165  { return _assembleDual(points,m_assemblerH,primal, m_dH); }
166 
167  ThinShellAssemblerStatus assembleDualL(const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed)
168  { return _assembleDual(m_assemblerL,primal,deformed, m_dL); }
169  ThinShellAssemblerStatus assembleDualH(const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed)
170  { return _assembleDual(m_assemblerH,primal,deformed, m_dH); }
171 
172  ThinShellAssemblerStatus assembleDualL(const bContainer & bnds, const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed)
173  { return _assembleDual(bnds,m_assemblerL,primal,deformed, m_dL); }
174  ThinShellAssemblerStatus assembleDualH(const bContainer & bnds, const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed)
175  { return _assembleDual(bnds,m_assemblerH,primal,deformed, m_dH); }
176 
177  ThinShellAssemblerStatus assembleDualL(const gsMatrix<T> & points, const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed)
178  { return _assembleDual(points,m_assemblerL,primal,deformed, m_dL); }
179  ThinShellAssemblerStatus assembleDualH(const gsMatrix<T> & points, const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed)
180  { return _assembleDual(points,m_assemblerH,primal,deformed, m_dH); }
181 
182  const gsSparseMatrix<T> & matrixL() const { return m_matrixL; }
183  const gsSparseMatrix<T> & matrixH() const { return m_matrixH; }
184 
185  const gsSparseMatrix<T> & massL() const { return m_massL; }
186  const gsSparseMatrix<T> & massH() const { return m_massH; }
187 
188  const gsMatrix<T> primalL() const { return m_pL; }
189  const gsMatrix<T> primalH() const { return m_pH; }
190 
191  const gsMatrix<T> dualL() const { return m_dL; }
192  const gsMatrix<T> dualH() const { return m_dH; }
193 
194  void updateMultiPatchL(const gsMatrix<T> & solVector, gsMultiPatch<T> & result);
195  void updateMultiPatchH(const gsMatrix<T> & solVector, gsMultiPatch<T> & result);
196 
197  void constructMultiPatchL(const gsMatrix<T> & solVector, gsMultiPatch<T> & result);
198  void constructMultiPatchH(const gsMatrix<T> & solVector, gsMultiPatch<T> & result);
199 
200  void constructSolutionL(const gsMatrix<T> & solVector, gsMultiPatch<T> & deformed);
201  void constructSolutionH(const gsMatrix<T> & solVector, gsMultiPatch<T> & deformed);
202  gsMultiPatch<T> constructSolutionL(const gsMatrix<T> & solVector);
203  gsMultiPatch<T> constructSolutionH(const gsMatrix<T> & solVector);
204  gsMultiPatch<T> constructDisplacementL(const gsMatrix<T> & solVector);
205  gsMultiPatch<T> constructDisplacementH(const gsMatrix<T> & solVector);
206  gsVector<T> constructSolutionVectorL(const gsMultiPatch<T> & deformed);
207  gsVector<T> constructSolutionVectorH(const gsMultiPatch<T> & deformed);
208 
209  gsMatrix<T> projectL2_L(const gsFunction<T> &fun);
210  gsMatrix<T> projectL2_H(const gsFunction<T> &fun);
211 
212  T deformationNorm(const gsMultiPatch<T> & deformed, const gsMultiPatch<T> & original)
213  { return m_assemblerL->deformationNorm(deformed,original); }
214 
215  // Linear elasticity
216  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);
217  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);
218  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);
219 
220  // Nonlinear elasticity
221  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);
222  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);
223  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);
224 
225  // Linear elasticity
226  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);
227  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);
228  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);
229 
230  // Nonlinear elasticity
231  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);
232  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);
233  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);
234 
235 
236  // Eigenvalues
237  T computeErrorEig(const T evPrimalL, const T evDualL, const T evDualH,
238  const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH,
239  const gsMultiPatch<T> & primal,std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
240  std::vector<T> computeErrorEigElements(const T evPrimalL, const T evDualL, const T evDualH,
241  const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH,
242  const gsMultiPatch<T> & primal,std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
243  std::vector<T> computeErrorEigDofs(const T evPrimalL, const T evDualL, const T evDualH,
244  const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH,
245  const gsMultiPatch<T> & primal,std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
246 
247  T computeErrorEig(const T evPrimalL, const T evDualL, const T evDualH,
248  const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH,
249  const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed,std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
250  std::vector<T> computeErrorEigElements(const T evPrimalL, const T evDualL, const T evDualH,
251  const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH,
252  const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed,std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
253  std::vector<T> computeErrorEigDofs(const T evPrimalL, const T evDualL, const T evDualH,
254  const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH,
255  const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed,std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
256 
257  T computeGoal(const gsMultiPatch<T> & deformed);
258  T computeGoal(const bContainer & bnds, const gsMultiPatch<T> & deformed);
259  T computeGoal(const gsMatrix<T> & points, const gsMultiPatch<T> & deformed);
260 
261  T matrixNorm(const gsMultiPatch<T> &dualL, const gsMultiPatch<T> &dualH) const;
262  T matrixNorm(const gsMultiPatch<T> &dualL, const gsMultiPatch<T> &dualH, const gsMultiPatch<T> &deformed) const;
263 
264  void setGoal(enum GoalFunction GF, short_t component = 9) { m_goalFunction = GF; m_component = component; }
265 
267  void constructStress(const gsFunctionSet<T> & deformed, gsPiecewiseFunction<T> & result, stress_type::type type)
268  {
269  m_assemblerL->constructStress(deformed,result,type);
270  }
271 
272  T error() const { return m_error; }
273  std::vector<T> errors() const { return m_errors; }
274  std::vector<T> absErrors() const
275  {
276  std::vector<T> result = m_errors;
277  for (typename std::vector<T>::iterator it = result.begin(); it!=result.end(); it++)
278  *it = std::abs(*it);
279  return result;
280  }
281 
282  gsDofMapper getMapperL() { return m_assemblerL->getMapper(); };
283  gsDofMapper getMapperH() { return m_assemblerH->getMapper(); };
284 
285 protected:
286 
287  ThinShellAssemblerStatus _assembleMass(gsThinShellAssemblerBase<T> * assembler, gsSparseMatrix<T> & result, bool lumped = false);
288 
289  ThinShellAssemblerStatus _assemble(gsThinShellAssemblerBase<T> * assembler, std::pair<gsSparseMatrix<T>,gsVector<T>> & result);
290 
291  ThinShellAssemblerStatus _assembleMatrix(gsThinShellAssemblerBase<T> * assembler, gsSparseMatrix<T> & result);
292  ThinShellAssemblerStatus _assembleMatrix(gsThinShellAssemblerBase<T> * assembler, const gsMultiPatch<T> & deformed, gsSparseMatrix<T> & result);
293  ThinShellAssemblerStatus _assemblePrimal(gsThinShellAssemblerBase<T> * assembler, gsVector<T> & result);
294  ThinShellAssemblerStatus _assemblePrimal(gsThinShellAssemblerBase<T> * assembler, const gsMultiPatch<T> & deformed, gsVector<T> & result);
295 
306  ThinShellAssemblerStatus _assembleDual(gsThinShellAssemblerBase<T> * assembler, const gsMultiPatch<T> & primal, gsVector<T> & result)
307  {
308  gsMultiPatch<T> deformed = m_patches;
309  for ( size_t k =0; k!=primal.nPatches(); ++k) // Deform the geometry
310  deformed.patch(k).coefs() += primal.patch(k).coefs(); // Gdef points to mp_def, therefore updated
311 
312  return _assembleDual(assembler,primal,deformed,result);
313  }
314  ThinShellAssemblerStatus _assembleDual(gsThinShellAssemblerBase<T> * assembler, const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed, gsVector<T> & result);
315 
327  ThinShellAssemblerStatus _assembleDual(const bContainer & bnds, gsThinShellAssemblerBase<T> * assembler, const gsMultiPatch<T> & primal, gsVector<T> & result)
328  {
329  gsMultiPatch<T> deformed = m_patches;
330  for ( size_t k =0; k!=primal.nPatches(); ++k) // Deform the geometry
331  deformed.patch(k).coefs() += primal.patch(k).coefs(); // Gdef points to mp_def, therefore updated
332 
333  return _assembleDual(bnds,assembler,primal,deformed,result);
334  }
335  ThinShellAssemblerStatus _assembleDual(const bContainer & bnds, gsThinShellAssemblerBase<T> * assembler, const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed, gsVector<T> & result);
336 
348  ThinShellAssemblerStatus _assembleDual(const gsMatrix<T> & points, gsThinShellAssemblerBase<T> * assembler, const gsMultiPatch<T> & primal, gsVector<T> & result)
349  {
350  gsMultiPatch<T> deformed = m_patches;
351  for ( size_t k =0; k!=primal.nPatches(); ++k) // Deform the geometry
352  deformed.patch(k).coefs() += primal.patch(k).coefs(); // Gdef points to mp_def, therefore updated
353 
354  return _assembleDual(points,assembler,primal,deformed,result);
355  }
356  ThinShellAssemblerStatus _assembleDual(const gsMatrix<T> & points, gsThinShellAssemblerBase<T> * assembler, const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed, gsVector<T> & result);
357 
358  template<int _elWise>
359  void computeError_impl(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, bool withLoads,
360  std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
361 
362  template<int _d, bool _bending, int _elWise>
363  typename std::enable_if<(_d==3) && _bending, void>::type
364  computeError_impl(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, const gsMultiPatch<T> & deformed, bool withLoads,
365  std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
366 
367  template<int _d, bool _bending, int _elWise>
368  typename std::enable_if<!(_d==3 && _bending), void>::type
369  computeError_impl(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, const gsMultiPatch<T> & deformed, bool withLoads,
370  std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
371 
372  template<int _elWise>
373  void computeSquaredError_impl(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, bool withLoads,
374  std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
375 
376  template<int _d, bool _bending, int _elWise>
377  typename std::enable_if<(_d==3) && _bending, void>::type
378  computeSquaredError_impl(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, const gsMultiPatch<T> & deformed, bool withLoads,
379  std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
380 
381  template<int _d, bool _bending, int _elWise>
382  typename std::enable_if<!(_d==3 && _bending), void>::type
383  computeSquaredError_impl(const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH, const gsMultiPatch<T> & deformed, bool withLoads,
384  std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
385 
386 
387  template<int _elWise>
388  void computeErrorEig_impl(const T evPrimalL, const T evDualL, const T evDualH,
389  const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH,
390  const gsMultiPatch<T> & primal,
391  std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
392 
393  template<int _elWise>
394  void computeErrorEig_impl(const T evPrimalL, const T evDualL, const T evDualH,
395  const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH,
396  const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed,
397  std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false);
398 
399  void _applyLoadsToElWiseError(const gsMultiPatch<T> &dualL, const gsMultiPatch<T> &dualH, std::vector<T> & errors) const;
400  void _applyLoadsToError(const gsMultiPatch<T> &dualL, const gsMultiPatch<T> &dualH, T & error) const;
401 
402 
403 protected:
404  mutable gsThinShellAssemblerBase<T> * m_assemblerL;
405  mutable gsThinShellAssemblerBase<T> * m_assemblerH;
406 
407  gsVector<T> m_pL, m_pH, m_dL, m_dH;
408  gsSparseMatrix<T> m_matrixL;
409  gsSparseMatrix<T> m_matrixH;
410  gsSparseMatrix<T> m_massL;
411  gsSparseMatrix<T> m_massH;
412 
413  gsPointLoads<T> m_pLoads;
414 
415  T m_error;
416  std::vector<T> m_errors;
417 
418  typedef gsExprAssembler<>::geometryMap geometryMap;
419  typedef gsExprAssembler<>::space space;
420  typedef gsExprAssembler<>::solution solution;
421 
422  gsBoundaryConditions<T> m_bcs;
423  // using Base::m_bcs;
424  // using Base::m_basis;
425 
426  using Base::m_patches;
427  // using Base::m_defpatches;
428  using Base::m_materialMatrices;
429  using Base::m_forceFun;
430  using Base::m_options;
431  using Base::m_foundFun;
432  using Base::m_foundInd;
433  using Base::m_pressFun;
434  using Base::m_pressInd;
435 
436  enum GoalFunction m_goalFunction;
437  short_t m_component;
438 };
439 
440 #ifdef GISMO_WITH_PYBIND11
441 
445  void pybind11_init_gsThinShellAssemblerDWR2(pybind11::module &m);
446  void pybind11_init_gsThinShellAssemblerDWR3(pybind11::module &m);
447  void pybind11_init_gsThinShellAssemblerDWR3nb(pybind11::module &m);
448 
449 #endif // GISMO_WITH_PYBIND11
450 
458 template <class T>
459 class gsThinShellAssemblerDWRBase //: public virtual gsThinShellAssemblerBase<T>
460 {
461 public:
462  typedef gsBoxTopology::bContainer bContainer;
463 
466 
467  virtual gsOptionList & optionsL() =0;
468  virtual gsOptionList & optionsH() =0;
469 
470  virtual gsExprAssembler<T> assemblerL() =0;
471  virtual gsExprAssembler<T> assemblerH() =0;
472 
473  virtual void setOptions(gsOptionList & options) =0;
474  virtual void setPointLoads(const gsPointLoads<T> & pLoads) =0;
475 
476  virtual void setFoundation(const gsFunction<T> & foundation) =0;
477  virtual void setPressure(const gsFunction<T> & pressure) =0;
478 
479  virtual void updateBCs(const gsBoundaryConditions<T> & bconditions) =0;
480 
481  virtual void setBasisL(const gsMultiBasis<T> & basis) =0;
482  virtual void setBasisH(const gsMultiBasis<T> & basis) =0;
483 
484  virtual void setUndeformed(const gsMultiPatch<T> & patches) =0;
485 
486  virtual void homogenizeDirichlet() =0;
487 
488  virtual index_t numDofsL() const =0;
489  virtual index_t numDofsH() const =0;
490 
491  virtual void setSpaceBasisL(const gsMultiPatch<T> & spaceBasis) =0;
492  virtual void setSpaceBasisH(const gsMultiPatch<T> & spaceBasis) =0;
493 
494  virtual ThinShellAssemblerStatus assembleL() =0;
495  virtual ThinShellAssemblerStatus assembleH() =0;
496 
497  virtual ThinShellAssemblerStatus assembleMassL(bool lumped = false) =0;
498  virtual ThinShellAssemblerStatus assembleMassH(bool lumped = false) =0;
499  virtual ThinShellAssemblerStatus assembleMatrixL() =0;
500  virtual ThinShellAssemblerStatus assembleMatrixH() =0;
501  virtual ThinShellAssemblerStatus assembleMatrixL(const gsMultiPatch<T> & deformed) =0;
502  virtual ThinShellAssemblerStatus assembleMatrixH(const gsMultiPatch<T> & deformed) =0;
503 
504  virtual ThinShellAssemblerStatus assemblePrimalL() =0;
505  virtual ThinShellAssemblerStatus assemblePrimalH() =0;
506  virtual ThinShellAssemblerStatus assemblePrimalL(const gsMultiPatch<T> & deformed) =0;
507  virtual ThinShellAssemblerStatus assemblePrimalH(const gsMultiPatch<T> & deformed) =0;
508 
509 
510  virtual ThinShellAssemblerStatus assembleDualL(const gsMultiPatch<T> & primal) =0;
511  virtual ThinShellAssemblerStatus assembleDualH(const gsMultiPatch<T> & primal) =0;
512 
513  virtual ThinShellAssemblerStatus assembleDualL(const bContainer & bnds, const gsMultiPatch<T> & primal) =0;
514  virtual ThinShellAssemblerStatus assembleDualH(const bContainer & bnds, const gsMultiPatch<T> & primal) =0;
515 
516  virtual ThinShellAssemblerStatus assembleDualL(const gsMatrix<T> & points, const gsMultiPatch<T> & primal) =0;
517  virtual ThinShellAssemblerStatus assembleDualH(const gsMatrix<T> & points, const gsMultiPatch<T> & primal) =0;
518 
519  virtual ThinShellAssemblerStatus assembleDualL(const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed) =0;
520  virtual ThinShellAssemblerStatus assembleDualH(const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed) =0;
521 
522  virtual ThinShellAssemblerStatus assembleDualL(const bContainer & bnds, const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed) =0;
523  virtual ThinShellAssemblerStatus assembleDualH(const bContainer & bnds, const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed) =0;
524 
525  virtual ThinShellAssemblerStatus assembleDualL(const gsMatrix<T> & points, const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed) =0;
526  virtual ThinShellAssemblerStatus assembleDualH(const gsMatrix<T> & points, const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed) =0;
527 
528  virtual const gsSparseMatrix<T> & matrixL() const =0;
529  virtual const gsSparseMatrix<T> & matrixH() const =0;
530 
531  virtual const gsSparseMatrix<T> & massL() const =0;
532  virtual const gsSparseMatrix<T> & massH() const =0;
533 
534  virtual const gsMatrix<T> primalL() const =0;
535  virtual const gsMatrix<T> primalH() const =0;
536 
537  virtual const gsMatrix<T> dualL() const =0;
538  virtual const gsMatrix<T> dualH() const =0;
539 
540  virtual void updateMultiPatchL(const gsMatrix<T> & solVector, gsMultiPatch<T> & result) =0;
541  virtual void updateMultiPatchH(const gsMatrix<T> & solVector, gsMultiPatch<T> & result) =0;
542 
543  virtual void constructMultiPatchL(const gsMatrix<T> & solVector, gsMultiPatch<T> & result) =0;
544  virtual void constructMultiPatchH(const gsMatrix<T> & solVector, gsMultiPatch<T> & result) =0;
545 
546  virtual void constructSolutionL(const gsMatrix<T> & solVector, gsMultiPatch<T> & deformed) =0;
547  virtual void constructSolutionH(const gsMatrix<T> & solVector, gsMultiPatch<T> & deformed) =0;
548  virtual gsMultiPatch<T> constructSolutionL(const gsMatrix<T> & solVector) =0;
549  virtual gsMultiPatch<T> constructSolutionH(const gsMatrix<T> & solVector) =0;
550  virtual gsMultiPatch<T> constructDisplacementL(const gsMatrix<T> & solVector) =0;
551  virtual gsMultiPatch<T> constructDisplacementH(const gsMatrix<T> & solVector) =0;
552  virtual gsVector<T> constructSolutionVectorL(const gsMultiPatch<T> & deformed) =0;
553  virtual gsVector<T> constructSolutionVectorH(const gsMultiPatch<T> & deformed) =0;
554 
555  virtual gsMatrix<T> projectL2_L(const gsFunction<T> &fun) =0;
556  virtual gsMatrix<T> projectL2_H(const gsFunction<T> &fun) =0;
557 
558  virtual T deformationNorm(const gsMultiPatch<T> & deformed, const gsMultiPatch<T> & original) =0;
559 
560  // Linear elasticity ;
561  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;
562  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;
563  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;
564 
565  // Nonlinear elasticity
566  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;
567  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;
568  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;
569 
570  // Linear elasticity ;
571  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;
572  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;
573  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;
574 
575  // Nonlinear elasticity
576  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;
577  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;
578  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;
579 
580 
581  // Eigenvalues
582  virtual T computeErrorEig(const T evPrimalL, const T evDualL, const T evDualH,
583  const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH,
584  const gsMultiPatch<T> & primal,std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false) = 0;
585  virtual std::vector<T> computeErrorEigElements(const T evPrimalL, const T evDualL, const T evDualH,
586  const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH,
587  const gsMultiPatch<T> & primal,std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false) = 0;
588  virtual std::vector<T> computeErrorEigDofs(const T evPrimalL, const T evDualL, const T evDualH,
589  const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH,
590  const gsMultiPatch<T> & primal,std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false) = 0;
591 
592  virtual T computeErrorEig(const T evPrimalL, const T evDualL, const T evDualH,
593  const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH,
594  const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed,std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false) = 0;
595  virtual std::vector<T> computeErrorEigElements(const T evPrimalL, const T evDualL, const T evDualH,
596  const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH,
597  const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed,std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false) = 0;
598  virtual std::vector<T> computeErrorEigDofs(const T evPrimalL, const T evDualL, const T evDualH,
599  const gsMultiPatch<T> & dualL, const gsMultiPatch<T> & dualH,
600  const gsMultiPatch<T> & primal, const gsMultiPatch<T> & deformed,std::string filename = std::string(), unsigned np=1000, bool parametric=false, bool mesh=false) = 0;
601 
602  virtual T computeGoal(const gsMultiPatch<T> & deformed) =0;
603  virtual T computeGoal(const bContainer & bnds, const gsMultiPatch<T> & deformed) =0;
604  virtual T computeGoal(const gsMatrix<T> & points, const gsMultiPatch<T> & deformed) =0;
605 
606  virtual T matrixNorm(const gsMultiPatch<T> &dualL, const gsMultiPatch<T> &dualH) const = 0;
607  virtual T matrixNorm(const gsMultiPatch<T> &dualL, const gsMultiPatch<T> &dualH, const gsMultiPatch<T> &deformed) const = 0;
608 
609  virtual void setGoal(enum GoalFunction GF, short_t component = 9) = 0;
610 
612  virtual void constructStress(const gsFunctionSet<T> & deformed,
613  gsPiecewiseFunction<T> & result,
614  stress_type::type type) = 0;
615 
616  virtual T error() const =0;
617  virtual std::vector<T> errors() const =0;
618  virtual std::vector<T> absErrors() const =0;
619 
620  virtual gsDofMapper getMapperL()=0;
621  virtual gsDofMapper getMapperH()=0;
622 
623 
624 };
625 
626 #ifdef GISMO_WITH_PYBIND11
627 
631  void pybind11_init_gsThinShellAssemblerDWRBase(pybind11::module &m);
632 
633 #endif // GISMO_WITH_PYBIND11
634 
635 
636 } // namespace gismo
637 
638 
641 
642 
643 #ifndef GISMO_BUILD_LIB
644 #include GISMO_HPP_HEADER(gsThinShellAssemblerDWR.hpp)
645 #endif
Definition: gsExprAssembler.h:30
ThinShellAssemblerStatus
Definition: gsThinShellAssembler.h:53
#define short_t
Definition: gsConfig.h:35
virtual ~gsThinShellAssemblerDWRBase()
Default empty constructor.
Definition: gsThinShellAssemblerDWR.h:465
virtual void constructStress(const gsFunctionSet< T > &deformed, gsPiecewiseFunction< T > &result, stress_type::type type)=0
See gsThinShellAssemblerBase for details.
#define index_t
Definition: gsConfig.h:32
A function from a n-dimensional domain to an m-dimensional image.
Definition: gsFunction.h:59
Maintains a mapping from patch-local dofs to global dof indices and allows the elimination of individ...
Definition: gsDofMapper.h:68
Holds a set of patch-wise bases and their topology information.
Definition: gsMultiBasis.h:36
Assembles the system matrix and vectors for 2D and 3D shell problems, including geometric nonlinearit...
Definition: gsThinShellAssembler.h:76
Class containing a set of points on a multi-patch isogeometric domain, together with boundary conditi...
Definition: gsPointLoads.h:64
Interface for the set of functions defined on a domain (the total number of functions in the set equa...
Definition: gsFuncData.h:23
Container class for a set of geometry patches and their topology, that is, the interface connections ...
Definition: gsMultiPatch.h:33
type
Definition: gsThinShellFunctions.h:40
Class containing a set of boundary conditions.
Definition: gsBoundaryConditions.h:341
void constructStress(const gsFunctionSet< T > &original, const gsFunctionSet< T > &deformed, gsPiecewiseFunction< T > &result, stress_type::type type)
See gsThinShellAssemblerBase for details.
Definition: gsThinShellAssembler.hpp:2760
Provides a base class for material matrices.
EIGEN_STRONG_INLINE abs_expr< E > abs(const E &u)
Absolute value.
Definition: gsExpressions.h:4486
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
Provides linear and nonlinear assemblers for thin shells.
Provides evaluation function for stresses.
Base class for the gsThinShellAssembler.
Definition: gsThinShellAssemblerDWR.h:51