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