G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsVisitorMixedLinearElasticity.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 gsVisitorMixedLinearElasticity
29{
30public:
31 gsVisitorMixedLinearElasticity(const gsPde<T> & pde_)
32 : pde_ptr(static_cast<const gsBasePde<T>*>(&pde_)) {}
33
34 void initialize(const gsBasisRefs<T> & basisRefs,
35 const index_t patchIndex,
36 const gsOptionList & options,
37 gsQuadRule<T> & rule)
38 {
39 GISMO_UNUSED(patchIndex);
40 // parametric dimension of the first displacement component
41 dim = basisRefs.front().dim();
42 // a quadrature rule is defined by the basis for the first velocity component.
43 // the same rule is used for the presure
44 rule = gsQuadrature::get(basisRefs.front(), options);
45 // saving necessary info
46 T E = options.getReal("YoungsModulus");
47 T pr = options.getReal("PoissonsRatio");
48 lambda_inv = ( 1. + pr ) * ( 1. - 2. * pr ) / E / pr ;
49 mu = E / ( 2. * ( 1. + pr ) );
50 forceScaling = options.getReal("ForceScaling");
51 I = gsMatrix<T>::Identity(dim,dim);
52 // resize containers for global indices
53 globalIndices.resize(dim+1);
54 blockNumbers.resize(dim+1);
55 }
56
57 inline void evaluate(const gsBasisRefs<T> & basisRefs,
58 const gsGeometry<T> & geo,
59 const gsMatrix<T> & quNodes)
60 {
61 // store quadrature points of the element for geometry evaluation
62 md.points = quNodes;
63 // NEED_VALUE to get points in the physical domain for evaluation of the RHS
64 // NEED_MEASURE to get the Jacobian determinant values for integration
65 // NEED_GRAD_TRANSFORM to get the Jacobian matrix to transform gradient from the parametric to physical domain
67 // Compute image of the quadrature points plus gradient, jacobian and other necessary data
68 geo.computeMap(md);
69 // find local indices of the displacement and pressure basis functions active on the element
70 basisRefs.front().active_into(quNodes.col(0),localIndicesDisp);
71 N_D = localIndicesDisp.rows();
72 basisRefs.back().active_into(quNodes.col(0), localIndicesPres);
73 N_P = localIndicesPres.rows();
74 // Evaluate displacement basis functions and their derivatives on the element
75 basisRefs.front().evalAllDers_into(quNodes,1,basisValuesDisp);
76 // Evaluate pressure basis functions on the element
77 basisRefs.back().eval_into(quNodes,basisValuesPres);
78 // Evaluate right-hand side at the image of the quadrature points
79 pde_ptr->rhs()->eval_into(md.values[0],forceValues);
80 }
81
82 inline void assemble(gsDomainIterator<T> & element,
83 const gsVector<T> & quWeights)
84 {
85 GISMO_UNUSED(element);
86 // Initialize local matrix/rhs // A | B^T
87 localMat.setZero(dim*N_D + N_P, dim*N_D + N_P); // --|-- matrix structure
88 localRhs.setZero(dim*N_D + N_P,1); // B | C
89 // elasticity tensor
90 symmetricIdentityTensor<T>(C,I);
91 C *= mu;
92 // Loop over the quadrature nodes
93 for (index_t q = 0; q < quWeights.rows(); ++q)
94 {
95 // Multiply quadrature weight by the geometry measure
96 const T weight = quWeights[q] * md.measure(q);
97
98 // Compute physical gradients of basis functions at q as a dim x numActiveFunction matrix
99 transformGradients(md, q, basisValuesDisp[1], physGradDisp);
100 // A-matrix: Loop over displacement basis functions
101 for (index_t i = 0; i < N_D; i++)
102 {
103 setB<T>(B_i,I,physGradDisp.col(i));
104 tempK = B_i.transpose() * C;
105 // Loop over displacement basis functions
106 for (index_t j = 0; j < N_D; j++)
107 {
108 setB<T>(B_j,I,physGradDisp.col(j));
109 K = tempK * B_j;
110
111 for (short_t di = 0; di < dim; ++di)
112 for (short_t dj = 0; dj < dim; ++dj)
113 localMat(di*N_D+i, dj*N_D+j) += weight * K(di,dj);
114 }
115 }
116 // B-matrix
117 for (short_t d = 0; d < dim; ++d)
118 {
119 block = weight*basisValuesPres.col(q)*physGradDisp.row(d);
120 localMat.block(dim*N_D,d*N_D,N_P,N_D) += block.block(0,0,N_P,N_D);
121 localMat.block(d*N_D,dim*N_D,N_D,N_P) += block.transpose().block(0,0,N_D,N_P);
122 }
123
124 // C-matrix
125 if (abs(lambda_inv) > 0)
126 localMat.block(dim*N_D,dim*N_D,N_P,N_P) -=
127 (weight*lambda_inv*basisValuesPres.col(q)*basisValuesPres.col(q).transpose()).block(0,0,N_P,N_P);
128
129 // rhs contribution
130 for (short_t d = 0; d < dim; ++d)
131 localRhs.middleRows(d*N_D,N_D).noalias() += weight * forceScaling * forceValues(d,q) * basisValuesDisp[0].col(q) ;
132 }
133 }
134
135 inline void localToGlobal(const int patchIndex,
136 const std::vector<gsMatrix<T> > & eliminatedDofs,
137 gsSparseSystem<T> & system)
138 {
139 // computes global indices for displacement components
140 for (short_t d = 0; d < dim; ++d)
141 {
142 system.mapColIndices(localIndicesDisp,patchIndex,globalIndices[d],d);
143 blockNumbers.at(d) = d;
144 }
145 // computes global indices for pressure
146 system.mapColIndices(localIndicesPres, patchIndex, globalIndices[dim], dim);
147 blockNumbers.at(dim) = dim;
148 // push to global system
149 system.pushToRhs(localRhs,globalIndices,blockNumbers);
150 system.pushToMatrix(localMat,globalIndices,eliminatedDofs,blockNumbers,blockNumbers);
151 }
152
153protected:
154 // problem info
155 short_t dim;
156 const gsBasePde<T> * pde_ptr;
157 // Lame coefficients and force scaling factor
158 T lambda_inv, mu, forceScaling;
159 // geometry mapping
160 gsMapData<T> md;
161 // local components of the global linear system
162 gsMatrix<T> localMat;
163 gsMatrix<T> localRhs;
164 // local indices (at the current patch) of basis functions active at the current element
165 gsMatrix<index_t> localIndicesDisp;
166 gsMatrix<index_t> localIndicesPres;
167 // number of displacement and pressure basis functions active at the current element
168 index_t N_D, N_P;
169 // values and derivatives of displacement basis functions at quadrature points at the current element
170 // values are stored as a N_D x numQuadPoints matrix; not sure about derivatives, must be smth like N_D*dim x numQuadPoints
171 std::vector<gsMatrix<T> > basisValuesDisp;
172 // values of pressure basis functions active at the current element;
173 // stores as a N_P x numQuadPoints matrix
174 gsMatrix<T> basisValuesPres;
175 // RHS values at quadrature points at the current element; stored as a dim x numQuadPoints matrix
176 gsMatrix<T> forceValues;
177
178 // all temporary matrices defined here for efficiency
179 gsMatrix<T> C, physGradDisp, B_i, tempK, B_j, K, block, I;
180 // containers for global indices
181 std::vector< gsMatrix<index_t> > globalIndices;
182 gsVector<index_t> blockNumbers;
183};
184
185} // 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.
EIGEN_STRONG_INLINE abs_expr< E > abs(const E &u)
Absolute value.
Definition gsExpressions.h:4486
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