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