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