61template<
class T>
class gsThinShellAssemblerBase;
75template <
short_t d,
class T,
bool bending>
79 typedef gsBoxTopology::ifContainer ifContainer;
181 this->_assembleDirichlet();
183 m_ddofs = m_space.fixedPart();
184 m_mapper = m_space.mapper();
211 template<
short_t _d,
bool _bending>
215 template<
short_t _d,
bool _bending>
257 template<
short_t _d,
bool _bending>
262 template<
short_t _d,
bool _bending>
267 template<
short_t _d,
bool _bending>
272 template<
short_t _d,
bool _bending>
307 template<
short_t _d,
bool _bending>
312 template<
short_t _d,
bool _bending>
322 std::vector<patchSide> patchSides(1);
330 template<
short_t _d,
bool _bending>
331 typename std::enable_if<(_d==3) && _bending,
gsMatrix<T> >::type
335 template<
short_t _d,
bool _bending>
336 typename std::enable_if<!(_d==3 && _bending),
gsMatrix<T> >::type
345 void setGeometry(
const gsMultiPatch<T> & patches)
355 this->setGeometry(patches);
363 gsMultiBasis<T> &
basis() {
return m_basis;}
376 m_spaceBasis = &spaceBasis;
405 void addStrongC0(
const gsBoxTopology::ifContainer & interfaces);
406 void addStrongC1(
const gsBoxTopology::ifContainer & interfaces);
407 void addWeakC0(
const gsBoxTopology::ifContainer & interfaces);
408 void addWeakC1(
const gsBoxTopology::ifContainer & interfaces);
409 void addUncoupled(
const gsBoxTopology::ifContainer & interfaces);
410 void initInterfaces();
502 void _defaultOptions();
505 void _assembleNeumann();
507 template <
bool _matrix>
509 template <
bool _matrix>
512 template <
bool _matrix>
514 template <
bool _matrix>
517 template <
bool _matrix>
518 void _assembleWeakBCs();
519 template <
bool _matrix>
522 template <
bool _matrix>
523 void _assembleWeakIfc();
524 template <
bool _matrix>
527 void _assembleDirichlet();
532 void _ifcTest(
const T tol = 1e-2);
537 typename std::enable_if<(_d==3),
void>::type
538 _assembleNeumann_impl();
541 typename std::enable_if<!(_d==3),
void>::type
542 _assembleNeumann_impl();
544 template<
short_t _d,
bool _matrix>
545 typename std::enable_if<(_d==3) && _matrix,
void>::type
548 template<
short_t _d,
bool _matrix>
549 typename std::enable_if<(_d==3) && !_matrix,
void>::type
552 template<
short_t _d,
bool _matrix>
553 typename std::enable_if<!(_d==3),
void>::type
556 template<
short_t _d,
bool _matrix>
557 typename std::enable_if<(_d==3) && _matrix,
void>::type
560 template<
short_t _d,
bool _matrix>
561 typename std::enable_if<(_d==3) && !_matrix,
void>::type
564 template<
short_t _d,
bool _matrix>
565 typename std::enable_if<!(_d==3),
void>::type
568 template<
short_t _d,
bool _matrix>
569 typename std::enable_if<(_d==3) && _matrix,
void>::type
572 template<
short_t _d,
bool _matrix>
573 typename std::enable_if<(_d==3) && !_matrix,
void>::type
576 template<
short_t _d,
bool _matrix>
577 typename std::enable_if<!(_d==3),
void>::type
580 template<
short_t _d,
bool _matrix>
581 typename std::enable_if<(_d==3) && _matrix,
void>::type
584 template<
short_t _d,
bool _matrix>
585 typename std::enable_if<(_d==3) && !_matrix,
void>::type
588 template<
short_t _d,
bool _matrix>
589 typename std::enable_if<!(_d==3),
void>::type
592 template<
short_t _d,
bool _matrix>
593 typename std::enable_if<(_d==3) && _matrix,
void>::type
594 _assembleWeakBCs_impl();
596 template<
short_t _d,
bool _matrix>
597 typename std::enable_if<(_d==3) && !_matrix,
void>::type
598 _assembleWeakBCs_impl();
600 template<
short_t _d,
bool _matrix>
601 typename std::enable_if<!(_d==3) && _matrix,
void>::type
602 _assembleWeakBCs_impl();
604 template<
short_t _d,
bool _matrix>
605 typename std::enable_if<!(_d==3) && !_matrix,
void>::type
606 _assembleWeakBCs_impl();
608 template<
short_t _d,
bool _matrix>
609 typename std::enable_if<(_d==3) && _matrix,
void>::type
612 template<
short_t _d,
bool _matrix>
613 typename std::enable_if<(_d==3) && !_matrix,
void>::type
616 template<
short_t _d,
bool _matrix>
617 typename std::enable_if<!(_d==3) && _matrix,
void>::type
620 template<
short_t _d,
bool _matrix>
621 typename std::enable_if<!(_d==3) && !_matrix,
void>::type
624 template<
short_t _d,
bool _matrix>
625 typename std::enable_if<(_d==3) && _matrix,
void>::type
626 _assembleWeakIfc_impl();
628 template<
short_t _d,
bool _matrix>
629 typename std::enable_if<(_d==3) && !_matrix,
void>::type
630 _assembleWeakIfc_impl();
632 template<
short_t _d,
bool _matrix>
633 typename std::enable_if<!(_d==3) && _matrix,
void>::type
634 _assembleWeakIfc_impl();
636 template<
short_t _d,
bool _matrix>
637 typename std::enable_if<!(_d==3) && !_matrix,
void>::type
638 _assembleWeakIfc_impl();
640 template<
short_t _d,
bool _matrix>
641 typename std::enable_if<(_d==3) && _matrix,
void>::type
644 template<
short_t _d,
bool _matrix>
645 typename std::enable_if<(_d==3) && !_matrix,
void>::type
648 template<
short_t _d,
bool _matrix>
649 typename std::enable_if<!(_d==3) && _matrix,
void>::type
652 template<
short_t _d,
bool _matrix>
653 typename std::enable_if<!(_d==3) && !_matrix,
void>::type
679 bool m_parametricForce;
693 mutable bool m_foundInd;
694 mutable bool m_pressInd;
698 mutable T m_alpha_d_bc,m_alpha_r_bc,m_alpha_d_ifc,m_alpha_r_ifc;
701 mutable ifContainer m_inPlane, m_outPlane, m_uncoupled, m_strongC0, m_weakC0, m_strongC1, m_weakC1, m_unassigned;
706#ifdef GISMO_WITH_PYBIND11
711 void pybind11_init_gsThinShellAssembler2(pybind11::module &m);
712 void pybind11_init_gsThinShellAssembler3(pybind11::module &m);
713 void pybind11_init_gsThinShellAssembler3nb(pybind11::module &m);
1013 virtual void addStrongC0(
const gsBoxTopology::ifContainer & interfaces) = 0;
1014 virtual void addStrongC1(
const gsBoxTopology::ifContainer & interfaces) = 0;
1015 virtual void addWeakC0(
const gsBoxTopology::ifContainer & interfaces) = 0;
1016 virtual void addWeakC1(
const gsBoxTopology::ifContainer & interfaces) = 0;
1017 virtual void addUncoupled(
const gsBoxTopology::ifContainer & interfaces) = 0;
1018 virtual void initInterfaces() = 0;
1083 virtual void plotSolution(std::string
string,
const gsMatrix<T> & solVector) = 0;
1092#ifdef GISMO_WITH_PYBIND11
1097 void pybind11_init_gsThinShellAssemblerBase(pybind11::module &m);
1098 void pybind11_enum_gsThinShellAssemblerStatus(pybind11::module &m);
1109#ifndef GISMO_BUILD_LIB
1110#include GISMO_HPP_HEADER(gsThinShellAssembler.hpp)
Definition gsExpressions.h:973
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
space trialSpace(const index_t id) const
Definition gsExprAssembler.h:257
gsExprHelper< T >::geometryMap geometryMap
Geometry map type.
Definition gsExprAssembler.h:65
index_t numDofs() const
Returns the number of degrees of freedom (after initialization)
Definition gsExprAssembler.h:92
const gsMatrix< T > & rhs() const
Returns the right-hand side vector(s)
Definition gsExprAssembler.h:154
gsExprHelper< T >::element element
Current element.
Definition gsExprAssembler.h:64
expr::gsFeSolution< T > solution
Solution type.
Definition gsExprAssembler.h:68
const gsSparseMatrix< T > & matrix() const
Returns the left-hand global matrix.
Definition gsExprAssembler.h:127
Generic evaluator of isogeometric expressions.
Definition gsExprEvaluator.h:39
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
This class defines the base class for material matrices.
Definition gsMaterialMatrixBase.h:33
memory::unique_ptr< gsMaterialMatrixBase > uPtr
Unique pointer for gsGeometry.
Definition gsMaterialMatrixBase.h:44
This class serves as the evaluator of material matrices, based on gsMaterialMatrixBase.
Definition gsMaterialMatrixContainer.h:34
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 gsThinShellAssembler.h:726
virtual gsMultiPatch< T > constructDisplacement(const gsMatrix< T > &solVector) const =0
Construct displacement field from computed solution vector solVector and returns a multipatch.
virtual void projectL2_into(const gsFunction< T > &fun, gsMatrix< T > &result)=0
Projects function fun on the basis and geometry stored in the class and returns the coefficients in r...
gsThinShellAssemblerBase()
Default deconstructor.
Definition gsThinShellAssembler.h:734
virtual gsMatrix< T > fullSolutionVector(const gsMatrix< T > &vector) const =0
Reconstruct the solution vector based on the currently stored boundary conditions (thus the mapper).
virtual gsVector< T > constructSolutionVector(const gsMultiPatch< T > &deformed) const =0
Reconstruct the solution vector based on the currently stored boundary conditions (thus the mapper).
virtual ThinShellAssemblerStatus assemblePressureMatrix(const gsFunction< T > &pressFun, const gsFunctionSet< T > &deformed)=0
Assembles the pressure contribution in the system matrix (non-linear)
virtual gsMultiPatch< T > constructMultiPatch(const gsMatrix< T > &solVector) const =0
Construct solution field from computed solution vector solVector and returns a multipatch.
gsBoxTopology::ifContainer ifContainer
Default deconstructor.
Definition gsThinShellAssembler.h:729
virtual gsMaterialMatrixContainer< T > materials() const =0
Returns the material matrices used in the class.
virtual T getElasticEnergy(const gsFunctionSet< T > &deformed)=0
Returns the elastic energy norm, i.e. norm = 0.5 * u'*F_int.
virtual ThinShellAssemblerStatus assembleMatrix(const gsFunctionSet< T > &deformed)=0
Assembles the tangential stiffness matrix (nonlinear)
virtual T getDisplacementNorm(const gsFunctionSet< T > &deformed)=0
Returns the displacement norm, i.e. norm = sqrt(u'*u/area)
virtual void setSpaceBasis(const gsFunctionSet< T > &spaceBasis)=0
Set the basis that is used for assembly (but not for quadrature!)
virtual gsMatrix< T > computePrincipalStresses(const gsMatrix< T > &points, const gsFunctionSet< T > &deformed, const T z=0)=0
Compute the principal stretches in points given a deformed geometry. Optionally, the stretches can be...
virtual ThinShellAssemblerStatus assemble()=0
Assembles the linear system and corresponding right-hand side.
virtual ThinShellAssemblerStatus status() const =0
Returns the assembler status.
virtual T getArea(const gsFunctionSet< T > &geometry)=0
Returns the area of geometry.
virtual ThinShellAssemblerStatus assemblePressureVector(const gsFunction< T > &pressFun)=0
Assembles the pressure contribution in the system vector (linear)
virtual void setPointLoads(const gsPointLoads< T > &pLoads)=0
Registers a gsPointLoads object for point loads acting on the shell.
virtual ThinShellAssemblerStatus assemblePressureVector(const T pressure, const gsFunctionSet< T > &deformed)=0
Assembles the pressure contribution in the system vector (non-linear)
virtual void projectL2_into(const gsFunction< T > &fun, gsMultiPatch< T > &result)=0
Projects function fun on the basis and geometry stored in the class and returns a multipatch in resul...
virtual ThinShellAssemblerStatus assemblePressureMatrix(const T pressure)=0
Assembles the pressure contribution in the system matrix (linear)
virtual void updateBCs(const gsBoundaryConditions< T > &bconditions)=0
Overwrites the boundary conditions.
virtual T interfaceErrorG1(const gsFunctionSet< T > &deformed)=0
Returns the G1 error over the interface.
virtual ThinShellAssemblerStatus assembleMatrix(const gsMatrix< T > &solVector)=0
Assembles the tangential stiffness matrix (nonlinear)
virtual const gsPointLoads< T > & getPointMass()=0
Gets the registered point masses.
virtual const gsMatrix< T > & rhs() const =0
Returns a reference to the right-hand side vector that is assembled.
virtual void constructDisplacement(const gsMatrix< T > &solVector, gsMultiPatch< T > &deformed) const =0
Construct displacement field from computed solution vector solVector and returns the result in deform...
virtual T interfaceErrorNormal(const gsFunctionSet< T > &deformed)=0
Returns the normal vector error over the interface.
virtual ThinShellAssemblerStatus assembleMass(const bool lumped=false)=0
Assembles the mass matrix (including density and thickness!); if lumped=true, a lumped mass matrix wi...
virtual const gsMultiBasis< T > & getBasis()=0
Gets the basis.
virtual const gsSparseMatrix< T > & matrix() const =0
Returns a reference to the system matrix that is assembled.
virtual gsSparseMatrix< T > & massMatrix()=0
Returns a reference to the mass matrix that is assembled.
virtual ThinShellAssemblerStatus assemblePressureVector(const T pressure)=0
Assembles the pressure contribution in the system vector (linear)
virtual ~gsThinShellAssemblerBase()
Default empty constructor.
Definition gsThinShellAssembler.h:737
virtual ThinShellAssemblerStatus assembleVector(const gsFunctionSet< T > &deformed, const bool homogenize=true)=0
Assembles the residual vector.
virtual const gsMultiPatch< T > & geometry() const =0
Returns the undeformed geometry.
virtual void setOptions(gsOptionList &options)=0
Sets the options of the assembler.
virtual const gsFunctionSet< T > & getSpaceBasis()=0
Get the basis that is used for assembly (but not for quadrature!)
virtual void setUndeformed(const gsMultiPatch< T > &patches)=0
Overwrites the undeformed geometry.
virtual void setBasis(const gsMultiBasis< T > &basis)=0
Overwrites the basis.
virtual gsExprAssembler< T > assembler()=0
Returns the internal expression assembler.
virtual gsMatrix< T > computePrincipalStretches(const gsMatrix< T > &points, const gsFunctionSet< T > &deformed, const T z=0)=0
Compute the principal stretches in points given a deformed geometry. Optionally, the stretches can be...
virtual ThinShellAssemblerStatus assemble(const gsFunctionSet< T > &deformed, const bool matrix=true, const bool homogenize=true)=0
Assembles the tangential stiffness matrix and the residual for an iteration of Newton's method.
virtual T interfaceErrorMeanCurvature(const gsFunctionSet< T > &deformed)=0
Returns the mean curvature error over the interface.
virtual T interfaceErrorC0(const gsFunctionSet< T > &deformed)=0
Returns the C1 error over the interface.
virtual gsMatrix< T > boundaryForce(const gsFunctionSet< T > &deformed, const std::vector< patchSide > &patchSides) const =0
Computes the force on a set of boundaries.
virtual ThinShellAssemblerStatus assembleVector(const gsMatrix< T > &solVector, const bool homogenize=true)=0
Assembles the residual vector.
virtual ThinShellAssemblerStatus assembleMatrix(const gsMatrix< T > &solVector, const gsMatrix< T > &prevVector)=0
Assembles the tangential stiffness matrix (nonlinear) using the Mixed Integration Point (MIP) method.
virtual T interfaceErrorGaussCurvature(const gsFunctionSet< T > &deformed)=0
Returns the Gaussian curvature error over the interface.
virtual void setPointMass(const gsPointLoads< T > &pLoads)=0
Registers a gsPointLoads object for a point mass acting on the shell. The point masss must be 1-dimen...
virtual gsDofMapper getMapper()=0
Returns the gsDofMapper.
virtual gsMaterialMatrixBase< T > * material(const index_t p) const =0
Returns the material matrix on patch p used in the class.
virtual gsOptionList & options()=0
Returns the options of the assembler.
virtual void homogenizeDirichlet()=0
Sets the Dirichlet BCs to zero.
virtual ThinShellAssemblerStatus assembleMatrix(const gsFunctionSet< T > &deformed, const gsFunctionSet< T > &previous, gsMatrix< T > &update)=0
Assembles the tangential stiffness matrix (nonlinear) using the Mixed Integration Point (MIP) method.
virtual void constructSolution(const gsMatrix< T > &solVector, gsMultiPatch< T > &deformed) const =0
Construct deformed shell geometry from computed solution vector solVector and returns the result in d...
virtual void setPressure(const gsFunction< T > &pressure)=0
Registers a scalar function acting as pressure (in normal direction) in normal direction.
virtual ThinShellAssemblerStatus assemblePressureMatrix(const T pressure, const gsFunctionSet< T > &deformed)=0
Assembles the pressure contribution in the system matrix (non-linear)
virtual gsMatrix< T > projectL2(const gsFunction< T > &fun)=0
Projects function fun on the basis and geometry stored in the class and returns the coefficients as a...
virtual index_t numDofs() const =0
Returns the number of degrees of freedom in the assembler.
virtual ThinShellAssemblerStatus assemblePressureVector(const gsFunction< T > &pressFun, const gsFunctionSet< T > &deformed)=0
Assembles the pressure contribution in the system vector (non-linear)
virtual ThinShellAssemblerStatus assemble(const gsMatrix< T > &solVector, const bool matrix=true, const bool homogenize=true)=0
Assembles the tangential stiffness matrix and the residual for an iteration of Newton's method.
virtual gsMultiPatch< T > constructSolution(const gsMatrix< T > &solVector) const =0
Construct deformed shell geometry from computed solution vector solVector and returns a multipatch.
virtual ThinShellAssemblerStatus assemblePressureMatrix(const gsFunction< T > &pressFun)=0
Assembles the pressure contribution in the system matrix (linear)
virtual ThinShellAssemblerStatus assembleFoundation()=0
Assembles the elastic foundation matrix.
virtual void constructStress(const gsFunctionSet< T > &deformed, gsPiecewiseFunction< T > &result, stress_type::type type)=0
Construct Cauchy stress tensor for visualization.
virtual void setFoundation(const gsFunction< T > &foundation)=0
Registers a stiffness function to be used for handling an elastic foundation, only relevant for 3D sh...
virtual const gsPointLoads< T > & getPointLoads()=0
Gets the registered point loads.
Assembles the system matrix and vectors for 2D and 3D shell problems, including geometric nonlinearit...
Definition gsThinShellAssembler.h:77
void setPointLoads(const gsPointLoads< T > &pLoads)
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.h:165
gsMultiBasis< T > & getBasis()
Gets the basis.
Definition gsThinShellAssembler.h:364
gsExprAssembler< T > assembler()
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.h:158
gsThinShellAssembler & operator=(const gsThinShellAssembler &other)
Assignment operator.
Definition gsThinShellAssembler.hpp:99
index_t numDofs() const
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.h:191
gsMatrix< T > boundaryForce(const gsFunctionSet< T > &deformed, const std::vector< patchSide > &patchSides) const
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.hpp:2300
ThinShellAssemblerStatus assembleFoundationVector(const gsFunctionSet< T > &deformed)
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.hpp:2270
void setOptions(gsOptionList &options)
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.hpp:229
T getDisplacementNorm(const gsFunctionSet< T > &deformed)
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.hpp:2520
const gsMatrix< T > & rhs() const
Returns a reference to the right-hand side vector that is assembled.
Definition gsThinShellAssembler.h:401
gsMaterialMatrixContainer< T > materials() const
Returns the material matrices used in the class.
Definition gsThinShellAssembler.h:394
gsMultiPatch< T > constructDisplacement(const gsMatrix< T > &solVector) const
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.hpp:2679
std::enable_if<!(_d==3 &&_bending), ThinShellAssemblerStatus >::type assembleMatrix_impl(const gsFunctionSet< T > &, const gsFunctionSet< T > &, gsMatrix< T > &)
Implementation of assembleMatrix for planar geometries (2D)
Definition gsThinShellAssembler.h:274
gsSparseMatrix< T > & massMatrix()
Returns a reference to the mass matrix that is assembled.
Definition gsThinShellAssembler.h:399
void _initialize()
Initializes the method.
Definition gsThinShellAssembler.hpp:243
const gsMultiBasis< T > & basis() const
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.h:360
T deformationNorm(const gsMultiPatch< T > &deformed, const gsMultiPatch< T > &original)
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.hpp:2895
gsThinShellAssembler()
Constructor for te shell assembler.
Definition gsThinShellAssembler.h:134
gsMatrix< T > computePrincipalStresses(const gsMatrix< T > &points, const gsFunctionSet< T > &deformed, const T z=0)
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.hpp:2763
gsMatrix< T > boundaryForce(const gsFunctionSet< T > &deformed, const patchSide &ps) const
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.h:320
const gsPointLoads< T > & getPointMass()
Gets the registered point masses.
Definition gsThinShellAssembler.h:168
gsDofMapper getMapper() const
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.h:207
std::enable_if<(_d==3)&&_bending, ThinShellAssemblerStatus >::type assemble_impl()
Specialisation of assemble() for surfaces (3D)
Definition gsThinShellAssembler.hpp:1540
ThinShellAssemblerStatus assembleVector(const gsFunctionSet< T > &deformed, const bool homogenize=true)
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.hpp:1980
void plotSolution(std::string string, const gsMatrix< T > &solVector)
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.hpp:2563
T interfaceErrorNormal(const gsFunctionSet< T > &deformed)
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.h:485
T getArea(const gsFunctionSet< T > &geometry)
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.hpp:2506
std::enable_if<(_d==3)&&_bending, ThinShellAssemblerStatus >::type assembleMatrix_impl(const gsFunctionSet< T > &deformed)
Implementation of assembleMatrix for surfaces (3D)
Definition gsThinShellAssembler.hpp:1711
T interfaceErrorG1(const gsFunctionSet< T > &deformed)
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.h:480
ThinShellAssemblerStatus assemblePressureVector(const gsFunction< T > &pressFun)
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.hpp:2193
T getElasticEnergy(const gsFunctionSet< T > &deformed)
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.hpp:2539
gsThinShellAssembler(const gsThinShellAssembler &other)
Copy constructor (makes deep copy)
Definition gsThinShellAssembler.h:137
ThinShellAssemblerStatus assemble()
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.hpp:1525
std::enable_if<(_d==3)&&_bending, gsMatrix< T > >::type boundaryForce_impl(const gsFunctionSet< T > &deformed, const std::vector< patchSide > &patchSides) const
Implementation of the boundary force vector for surfaces (3D)
Definition gsThinShellAssembler.hpp:2308
ThinShellAssemblerStatus assembleMass(const bool lumped=false)
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.hpp:1436
void constructStress(const gsFunctionSet< T > &deformed, gsPiecewiseFunction< T > &result, stress_type::type type)
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.hpp:2795
std::enable_if<(_d==3)&&_bending, ThinShellAssemblerStatus >::type assembleVector_impl(const gsFunctionSet< T > &deformed, const bool homogenize)
Implementation of assembleVector for surfaces (3D)
Definition gsThinShellAssembler.hpp:1988
const gsMultiPatch< T > & geometry() const
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.h:343
gsMatrix< T > fullSolutionVector(const gsMatrix< T > &vector) const
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.hpp:2685
void updateBCs(const gsBoundaryConditions< T > &bconditions)
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.h:177
void setUndeformed(const gsMultiPatch< T > &patches)
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.h:353
void setFoundation(const gsFunction< T > &foundation)
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.h:171
ThinShellAssemblerStatus assembleMatrix(const gsFunctionSet< T > &deformed)
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.hpp:1703
void setPressure(const gsFunction< T > &pressure)
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.h:174
gsOptionList & options()
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.h:155
T interfaceErrorGaussCurvature(const gsFunctionSet< T > &deformed)
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.h:490
gsMultiPatch< T > constructSolution(const gsMatrix< T > &solVector) const
Construct deformed shell geometry from computed solution vector solVector and returns a multipatch.
Definition gsThinShellAssembler.hpp:2488
void projectL2_into(const gsFunction< T > &fun, gsMatrix< T > &result)
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.hpp:2832
ThinShellAssemblerStatus status() const
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.h:194
void setPointMass(const gsPointLoads< T > &pMass)
Registers a gsPointLoads object for a point mass acting on the shell. The point masss must be 1-dimen...
Definition gsThinShellAssembler.h:167
const gsPointLoads< T > & getPointLoads()
Gets the registered point loads.
Definition gsThinShellAssembler.h:166
void setBasis(const gsMultiBasis< T > &basis)
Overwrites the basis.
Definition gsThinShellAssembler.h:366
gsMultiPatch< T > constructMultiPatch(const gsMatrix< T > &solVector) const
Construct solution field from computed solution vector solVector and returns a multipatch.
Definition gsThinShellAssembler.hpp:2623
ThinShellAssemblerStatus assembleFoundation()
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.hpp:1496
T interfaceErrorMeanCurvature(const gsFunctionSet< T > &deformed)
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.h:495
gsVector< T > constructSolutionVector(const gsMultiPatch< T > &deformed) const
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.hpp:2697
T interfaceErrorC0(const gsFunctionSet< T > &deformed)
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.h:475
void setSpaceBasis(const gsFunctionSet< T > &spaceBasis)
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.h:374
gsMatrix< T > projectL2(const gsFunction< T > &fun)
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.hpp:2885
const gsFunctionSet< T > & getSpaceBasis()
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.h:204
gsMultiPatch< T > _constructSolution(const gsMatrix< T > &solVector, const gsMultiPatch< T > &undeformed) const
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.hpp:2477
gsMaterialMatrixBase< T > * material(const index_t p) const
Returns the material matrix on patch p used in the class.
Definition gsThinShellAssembler.h:395
void homogenizeDirichlet()
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.hpp:1314
ThinShellAssemblerStatus assemblePressureMatrix(const gsFunction< T > &pressFun)
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.hpp:2131
const gsSparseMatrix< T > & matrix() const
Returns a reference to the system matrix that is assembled.
Definition gsThinShellAssembler.h:398
gsMatrix< T > computePrincipalStretches(const gsMatrix< T > &points, const gsFunctionSet< T > &deformed, const T z=0)
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.hpp:2733
gsDofMapper getMapper()
See gsThinShellAssemblerBase for details.
Definition gsThinShellAssembler.h:472
gsThinShellAssembler(gsThinShellAssembler &&other)
Move constructor.
Definition gsThinShellAssembler.h:143
A vector with arbitrary coefficient type and fixed or dynamic size.
Definition gsVector.h:37
#define index_t
Definition gsConfig.h:32
#define GISMO_NO_IMPLEMENTATION
Definition gsDebug.h:129
Generic expressions matrix assembly.
Generic expressions evaluator.
Provides a base class for material matrices.
Provides a container for material matrices for thin shells.
Provides a simple container for point loads on multi-patch domains.
Provides evaluation function for stresses.
The G+Smo namespace, containing all definitions for the library.
ThinShellAssemblerStatus
Definition gsThinShellAssembler.h:54
@ DimensionError
Assembly failed due to a dimension error.
@ Success
Assembly is successful.
@ AssemblyError
Assembly failed due to an error in the expression (e.g. overflow)
S give(S &x)
Definition gsMemory.h:266
Struct which represents an interface between two patches.
Definition gsBoundary.h:650
Struct which represents a certain side of a patch.
Definition gsBoundary.h:232
Defines the coupling type over interfaces.
Definition gsThinShellAssembler.h:34
type
Definition gsThinShellFunctions.h:39