G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsVisitorMassElasticity.h
Go to the documentation of this file.
1
16#pragma once
17
19#include <gsCore/gsFuncData.h>
20
21namespace gismo
22{
23
24template <class T>
25class gsVisitorMassElasticity
26{
27public:
28
29 gsVisitorMassElasticity(gsSparseMatrix<T> * elimMatrix = nullptr) :
30 elimMat(elimMatrix) {}
31
32 void initialize(const gsBasisRefs<T> & basisRefs,
33 const index_t patchIndex,
34 const gsOptionList & options,
35 gsQuadRule<T> & rule)
36 {
37 GISMO_UNUSED(patchIndex);
38 // parametric dimension of the first displacement component
39 dim = basisRefs.front().dim();
40 // a quadrature rule is defined by the basis for the first displacement component.
41 rule = gsQuadrature::get(basisRefs.front(), options);
42 // saving necessary info
43 density = options.getReal("Density");
44 // resize containers for global indices
45 globalIndices.resize(dim);
46 blockNumbers.resize(dim);
47 }
48
49 inline void evaluate(const gsBasisRefs<T> & basisRefs,
50 const gsGeometry<T> & geo,
51 const gsMatrix<T> & quNodes)
52 {
53 // store quadrature points of the element for geometry evaluation
54 md.points = quNodes;
55 // NEED_MEASURE to get the Jacobian determinant values for integration
56 md.flags = NEED_MEASURE;
57 // Compute the geometry mapping at the quadrature points
58 geo.computeMap(md);
59 // Evaluate displacement basis functions on the element
60 basisRefs.front().eval_into(quNodes,basisValuesDisp);
61 // find local indices of the displacement basis functions active on the element
62 basisRefs.front().active_into(quNodes.col(0),localIndicesDisp);
63 N_D = localIndicesDisp.rows();
64
65 }
66
67 inline void assemble(gsDomainIterator<T> & element,
68 const gsVector<T> & quWeights)
69 {
70 GISMO_UNUSED(element);
71 // initialize local matrix and rhs
72 localMat.setZero(dim*N_D,dim*N_D);
73 block = density*basisValuesDisp * quWeights.asDiagonal() * md.measures.asDiagonal() * basisValuesDisp.transpose();
74 for (short_t d = 0; d < dim; ++d)
75 localMat.block(d*N_D,d*N_D,N_D,N_D) = block.block(0,0,N_D,N_D);
76 }
77
78 inline void localToGlobal(const int patchIndex,
79 const std::vector<gsMatrix<T> > & eliminatedDofs,
80 gsSparseSystem<T> & system)
81 {
82 // computes global indices for displacement components
83 for (short_t d = 0; d < dim; ++d)
84 {
85 system.mapColIndices(localIndicesDisp, patchIndex, globalIndices[d], d);
86 blockNumbers.at(d) = d;
87 }
88 // push to global system
89 system.pushToMatrix(localMat,globalIndices,eliminatedDofs,blockNumbers,blockNumbers);
90
91 // push to the elimination system
92 if (elimMat != nullptr)
93 {
94 index_t globalI,globalElimJ;
95 index_t elimSize = 0;
96 for (short_t dJ = 0; dJ < dim; ++dJ)
97 {
98 for (short_t dI = 0; dI < dim; ++dI)
99 for (index_t i = 0; i < N_D; ++i)
100 if (system.colMapper(dI).is_free_index(globalIndices[dI].at(i)))
101 {
102 system.mapToGlobalRowIndex(localIndicesDisp.at(i),patchIndex,globalI,dI);
103 for (index_t j = 0; j < N_D; ++j)
104 if (!system.colMapper(dJ).is_free_index(globalIndices[dJ].at(j)))
105 {
106 globalElimJ = system.colMapper(dJ).global_to_bindex(globalIndices[dJ].at(j));
107 elimMat->coeffRef(globalI,elimSize+globalElimJ) += localMat(N_D*dI+i,N_D*dJ+j);
108 }
109 }
110 elimSize += eliminatedDofs[dJ].rows();
111 }
112 }
113 }
114
115protected:
116 // problem info
117 short_t dim;
118 //density
119 T density;
120 // geometry mapping
121 gsMapData<T> md;
122 // local components of the global linear system
123 gsMatrix<T> localMat;
124 // local indices (at the current patch) of the displacement basis functions active at the current element
125 gsMatrix<index_t> localIndicesDisp;
126 // number of displacement basis functions active at the current element
127 index_t N_D;
128 // values of displacement basis functions at quadrature points at the current element stored as a N_D x numQuadPoints matrix;
129 gsMatrix<T> basisValuesDisp;
130 bool assembleMatrix;
131 // elimination matrix to efficiently change Dirichlet degrees of freedom
132 gsSparseMatrix<T> * elimMat;
133
134 // all temporary matrices defined here for efficiency
135 gsMatrix<T> block;
136 // containers for global indices
137 std::vector< gsMatrix<index_t> > globalIndices;
138 gsVector<index_t> blockNumbers;
139};
140
141} // namespace gismo
142
#define short_t
Definition gsConfig.h:35
#define index_t
Definition gsConfig.h:32
#define GISMO_UNUSED(x)
Definition gsDebug.h:112
This object is a cache for computed values from an evaluator.
Creates a variety of quadrature rules.
The G+Smo namespace, containing all definitions for the library.
@ NEED_MEASURE
The density of the measure pull back.
Definition gsForwardDeclarations.h:76
static gsQuadRule< T > get(const gsBasis< T > &basis, const gsOptionList &options, short_t fixDir=-1)
Constructs a quadrature rule based on input options.
Definition gsQuadrature.h:51