G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsVisitorLinearElasticity.h
Go to the documentation of this file.
1
16#pragma once
17
20
22#include <gsCore/gsFuncData.h>
23
24namespace gismo
25{
26
27template <class T>
28class gsVisitorLinearElasticity
29{
30public:
31
32 gsVisitorLinearElasticity(const gsPde<T> & pde_, gsSparseMatrix<T> * elimMatrix = nullptr)
33 : pde_ptr(static_cast<const gsBasePde<T>*>(&pde_)),
34 elimMat(elimMatrix)
35 {}
36
37 void initialize(const gsBasisRefs<T> & basisRefs,
38 const index_t patchIndex,
39 const gsOptionList & options,
40 gsQuadRule<T> & rule)
41 {
42 GISMO_UNUSED(patchIndex);
43 // parametric dimension of the first displacement component
44 dim = basisRefs.front().dim();
45 // a quadrature rule is defined by the basis for the first displacement component.
46 rule = gsQuadrature::get(basisRefs.front(), options);
47 // saving necessary info
48 T E = options.getReal("YoungsModulus");
49 T pr = options.getReal("PoissonsRatio");
50 lambda = E * pr / ( ( 1. + pr ) * ( 1. - 2. * pr ) );
51 mu = E / ( 2. * ( 1. + pr ) );
52 forceScaling = options.getReal("ForceScaling");
53 localStiffening = options.getReal("LocalStiff");
54 // linear elasticity tensor
55 I = gsMatrix<T>::Identity(dim,dim);
56 matrixTraceTensor<T>(C,I,I);
57 C *= lambda;
58 symmetricIdentityTensor<T>(Ctemp,I);
59 C += mu*Ctemp;
60 // resize containers for global indices
61 globalIndices.resize(dim);
62 blockNumbers.resize(dim);
63 }
64
65 inline void evaluate(const gsBasisRefs<T> & basisRefs,
66 const gsGeometry<T> & geo,
67 const gsMatrix<T> & quNodes)
68 {
69 // store quadrature points of the element for geometry evaluation
70 md.points = quNodes;
71 // NEED_VALUE to get points in the physical domain for evaluation of the RHS
72 // NEED_MEASURE to get the Jacobian determinant values for integration
73 // NEED_GRAD_TRANSFORM to get the Jacobian matrix to transform gradient from the parametric to physical domain
75 // Compute image of the quadrature points plus gradient, jacobian and other necessary data
76 geo.computeMap(md);
77 // find local indices of the displacement basis functions active on the element
78 basisRefs.front().active_into(quNodes.col(0),localIndicesDisp);
79 N_D = localIndicesDisp.rows();
80 // Evaluate displacement basis functions and their derivatives on the element
81 basisRefs.front().evalAllDers_into(quNodes,1,basisValuesDisp);
82 // Evaluate right-hand side at the image of the quadrature points
83 pde_ptr->rhs()->eval_into(md.values[0],forceValues);
84 }
85
86 inline void assemble(gsDomainIterator<T> & element,
87 const gsVector<T> & quWeights)
88 {
89 GISMO_UNUSED(element);
90 // initialize local matrix and rhs
91 localMat.setZero(dim*N_D,dim*N_D);
92 localRhs.setZero(dim*N_D,1);
93 // Loop over the quadrature nodes
94 for (index_t q = 0; q < quWeights.rows(); ++q)
95 {
96 // Multiply quadrature weight by the geometry measure
97 const T weightForce = quWeights[q] * md.measure(q);
98 const T weightBody = quWeights[q] * pow(md.measure(q),1-localStiffening);
99 // Compute physical gradients of basis functions at q as a dim x numActiveFunction matrix
100 transformGradients(md,q,basisValuesDisp[1],physGrad);
101 // loop over active basis functions (v_j)
102 for (index_t i = 0; i < N_D; i++)
103 {
104 // stiffness matrix K = B_i^T * C * B_j;
105 setB<T>(B_i,I,physGrad.col(i));
106 tempK = B_i.transpose() * C;
107 // loop over active basis functions (v_j)
108 for (index_t j = 0; j < N_D; j++)
109 {
110 setB<T>(B_j,I,physGrad.col(j));
111 K = tempK * B_j;
112 for (short_t di = 0; di < dim; ++di)
113 for (short_t dj = 0; dj < dim; ++dj)
114 localMat(di*N_D+i,dj*N_D+j) += weightBody * K(di,dj);
115 }
116 }
117 // rhs contribution
118 for (short_t d = 0; d < dim; ++d)
119 localRhs.middleRows(d*N_D,N_D).noalias() += weightForce * forceScaling * forceValues(d,q) * basisValuesDisp[0].col(q) ;
120 }
121 }
122
123 inline void localToGlobal(const int patchIndex,
124 const std::vector<gsMatrix<T> > & eliminatedDofs,
125 gsSparseSystem<T> & system)
126 {
127 // computes global indices for displacement components
128 for (short_t d = 0; d < dim; ++d)
129 {
130 system.mapColIndices(localIndicesDisp, patchIndex, globalIndices[d], d);
131 blockNumbers.at(d) = d;
132 }
133 // push to global system
134 system.pushToRhs(localRhs,globalIndices,blockNumbers);
135 system.pushToMatrix(localMat,globalIndices,eliminatedDofs,blockNumbers,blockNumbers);
136 // push to the elimination system
137 if (elimMat != nullptr)
138 {
139 index_t globalI,globalElimJ;
140 index_t elimSize = 0;
141 for (short_t dJ = 0; dJ < dim; ++dJ)
142 {
143 for (short_t dI = 0; dI < dim; ++dI)
144 for (index_t i = 0; i < N_D; ++i)
145 if (system.colMapper(dI).is_free_index(globalIndices[dI].at(i)))
146 {
147 system.mapToGlobalRowIndex(localIndicesDisp.at(i),patchIndex,globalI,dI);
148 for (index_t j = 0; j < N_D; ++j)
149 if (!system.colMapper(dJ).is_free_index(globalIndices[dJ].at(j)))
150 {
151 globalElimJ = system.colMapper(dJ).global_to_bindex(globalIndices[dJ].at(j));
152 elimMat->coeffRef(globalI,elimSize+globalElimJ) += localMat(N_D*dI+i,N_D*dJ+j);
153 }
154 }
155 elimSize += eliminatedDofs[dJ].rows();
156 }
157 }
158 }
159
160protected:
161 // problem info
162 short_t dim;
163 const gsBasePde<T> * pde_ptr;
164 // Lame coefficients and force scaling factor
165 T lambda, mu, forceScaling, localStiffening;
166 // geometry mapping
167 gsMapData<T> md;
168 // local components of the global linear system
169 gsMatrix<T> localMat;
170 gsMatrix<T> localRhs;
171 // local indices (at the current patch) of the displacement basis functions active at the current element
172 gsMatrix<index_t> localIndicesDisp;
173 // number of displacement basis functions active at the current element
174 index_t N_D;
175 // values and derivatives of displacement basis functions at quadrature points at the current element
176 // values are stored as a N_D x numQuadPoints matrix; not sure about derivatives, must be smth like N_D*dim x numQuadPoints
177 std::vector<gsMatrix<T> > basisValuesDisp;
178 // RHS values at quadrature points at the current element; stored as a dim x numQuadPoints matrix
179 gsMatrix<T> forceValues;
180 // elimination matrix to efficiently change Dirichlet degrees of freedom
181 gsSparseMatrix<T> * elimMat;
182 // all temporary matrices defined here for efficiency
183 gsMatrix<T> C, Ctemp,physGrad, B_i, tempK, B_j, K, I;
184 // containers for global indices
185 std::vector< gsMatrix<index_t> > globalIndices;
186 gsVector<index_t> blockNumbers;
187};
188
189} // namespace gismo
IMHO, a useless class, but it is necessary to use the gsAssembler class. Contains proper information ...
#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.
Tensor operations for elasticity.
The G+Smo namespace, containing all definitions for the library.
@ NEED_VALUE
Value of the object.
Definition gsForwardDeclarations.h:72
@ NEED_GRAD_TRANSFORM
Gradient transformation matrix.
Definition gsForwardDeclarations.h:77
@ 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