G+Smo  24.08.0
Geometry + Simulation Modules
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
gsVisitorElPoisson.h
Go to the documentation of this file.
1 
15 #pragma once
16 
18 #include <gsCore/gsFuncData.h>
19 
20 namespace gismo
21 {
22 
23 template <class T>
24 class gsVisitorElPoisson
25 {
26 public:
27 
28  gsVisitorElPoisson(const gsPde<T> & pde_, gsSparseMatrix<T> * elimMatrix = nullptr)
29  : pde_ptr(static_cast<const gsPoissonPde<T>*>(&pde_)),
30  elimMat(elimMatrix)
31  {}
32 
33  void initialize(const gsBasisRefs<T> & basisRefs,
34  const index_t patchIndex,
35  const gsOptionList & options,
36  gsQuadRule<T> & rule)
37  {
38  // a quadrature rule is defined by the basis for the first displacement component.
39  rule = gsQuadrature::get(basisRefs.front(), options);
40  // saving necessary info
41  localStiffening = options.getReal("LocalStiff");
42  // resize containers for global indices
43  globalIndices.resize(1);
44  blockNumbers.resize(1);
45  }
46 
47  inline void evaluate(const gsBasisRefs<T> & basisRefs,
48  const gsGeometry<T> & geo,
49  const gsMatrix<T> & quNodes)
50  {
51  // store quadrature points of the element for geometry evaluation
52  md.points = quNodes;
53  // NEED_MEASURE to get the Jacobian determinant values for integration
55  // Compute the geometry mapping at the quadrature points
56  geo.computeMap(md);
57  // Evaluate displacement basis functions on the element
58  basisRefs.front().evalAllDers_into(quNodes,1,basisValues);
59  // find local indices of the displacement basis functions active on the element
60  basisRefs.front().active_into(quNodes.col(0),localIndices);
61  N = localIndices.rows();
62  pde_ptr->rhs()->eval_into(md.values[0],forceValues);
63  }
64 
65  inline void assemble(gsDomainIterator<T> & element,
66  const gsVector<T> & quWeights)
67  {
68  // initialize local matrix and rhs
69  localMat.setZero(N,N);
70  localRhs.setZero(N,pde_ptr->numRhs());
71  for (index_t q = 0; q < quWeights.rows(); ++q)
72  {
73  // Multiply quadrature weight by the geometry measure
74  const T weightMatrix = quWeights[q] * pow(md.measure(q),1-localStiffening);
75  const T weightRHS = quWeights[q] * md.measure(q);
76  transformGradients(md,q,basisValues[1],physGrad);
77  localMat.noalias() += weightMatrix * (physGrad.transpose() * physGrad);
78  localRhs.noalias() += weightRHS * basisValues[0].col(q) * forceValues.col(q).transpose();
79  }
80  }
81 
82  inline void localToGlobal(const int patchIndex,
83  const std::vector<gsMatrix<T> > & eliminatedDofs,
84  gsSparseSystem<T> & system)
85  {
86  // computes global indices for displacement components
87  system.mapColIndices(localIndices, patchIndex, globalIndices[0], 0);
88  blockNumbers.at(0) = 0;
89  // push to global system
90  system.pushToMatrix(localMat,globalIndices,eliminatedDofs,blockNumbers,blockNumbers);
91  system.pushToRhs(localRhs,globalIndices,blockNumbers);
92 
93  // push to the elimination matrix
94  if (elimMat != nullptr)
95  {
96  index_t globalI, globalElimJ;
97  for (index_t i = 0; i < N; ++i)
98  if (system.colMapper(0).is_free_index(globalIndices[0].at(i)))
99  {
100  system.mapToGlobalRowIndex(localIndices.at(i),patchIndex,globalI,0);
101  for (index_t j = 0; j < N; ++j)
102  if (!system.colMapper(0).is_free_index(globalIndices[0].at(j)))
103  {
104  globalElimJ = system.colMapper(0).global_to_bindex(globalIndices[0].at(j));
105  elimMat->coeffRef(globalI,globalElimJ) += localMat(i,j);
106  }
107  }
108  }
109  }
110 
111 protected:
112  // geometry mapping
113  gsMapData<T> md;
114  const gsPoissonPde<T> * pde_ptr;
115  // local components of the global linear system
116  gsMatrix<T> localMat;
117  gsMatrix<T> localRhs;
118  // local indices (at the current patch) of the displacement basis functions active at the current element
119  gsMatrix<index_t> localIndices;
120  // number of displacement basis functions active at the current element
121  index_t N;
122  // values of displacement basis functions at quadrature points at the current element stored as a N_D x numQuadPoints matrix;
123  std::vector<gsMatrix<T> >basisValues;
124  // RHS values at quadrature points at the current element; stored as a dim x numQuadPoints matrix
125  gsMatrix<T> forceValues;
126  // elimination matrix to efficiently change Dirichlet degrees of freedom
127  gsSparseMatrix<T> * elimMat;
128 
129  // all temporary matrices defined here for efficiency
130  gsMatrix<T> physGrad;
131  real_t localStiffening;
132  // containers for global indices
133  std::vector< gsMatrix<index_t> > globalIndices;
134  gsVector<index_t> blockNumbers;
135 };
136 
137 } // namespace gismo
138 
#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.
The density of the measure pull back.
Definition: gsForwardDeclarations.h:76
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.