G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsVisitorStokes.h
Go to the documentation of this file.
1
15#pragma once
16
18#include <gsCore/gsFuncData.h>
20
21namespace gismo
22{
23
24template <class T>
25class gsVisitorStokes
26{
27public:
28
29 gsVisitorStokes(const gsPde<T> & pde_)
30 : dim(), pde_ptr(static_cast<const gsBasePde<T>*>(&pde_)),
31 viscosity(), density(), N_V(), N_P()
32 {}
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 viscosity = options.getReal("Viscosity");
47 density = options.getReal("Density");
48 // resize containers for global indices
49 globalIndices.resize(dim+1);
50 blockNumbers.resize(dim+1);
51 }
52
53 inline void evaluate(const gsBasisRefs<T> & basisRefs,
54 const gsGeometry<T> & geo,
55 const gsMatrix<T> & quNodes)
56 {
57 // store quadrature points of the element for geometry evaluation
58 md.points = quNodes;
59 // NEED_VALUE to get points in the physical domain for evaluation of the RHS
60 // NEED_MEASURE to get the Jacobian determinant values for integration
61 // NEED_GRAD_TRANSFORM to get the Jacobian matrix to transform gradient from the parametric to physical domain
63 // Compute image of the quadrature points plus gradient, jacobian and other necessary data
64 geo.computeMap(md);
65 // find local indices of the velocity and pressure basis functions active on the element
66 basisRefs.front().active_into(quNodes.col(0),localIndicesVel);
67 N_V = localIndicesVel.rows();
68 basisRefs.back().active_into(quNodes.col(0), localIndicesPres);
69 N_P = localIndicesPres.rows();
70 // Evaluate velocity basis functions and their derivatives on the element
71 basisRefs.front().evalAllDers_into(quNodes,1,basisValuesVel);
72 // Evaluate pressure basis functions on the element
73 basisRefs.back().eval_into(quNodes,basisValuesPres);
74 // Evaluate right-hand side at the image of the quadrature points
75 pde_ptr->rhs()->eval_into(md.values[0],forceValues);
76 }
77
78 inline void assemble(gsDomainIterator<T> & element,
79 const gsVector<T> & quWeights)
80 {
81 GISMO_UNUSED(element);
82 // Initialize local matrix/rhs // A | B^T
83 localMat.setZero(dim*N_V + N_P, dim*N_V + N_P); // --|-- matrix structure
84 localRhs.setZero(dim*N_V + N_P,1); // B | 0
85
86 // Loop over the quadrature nodes
87 for (index_t q = 0; q < quWeights.rows(); ++q)
88 {
89 // Multiply quadrature weight by the geometry measure
90 const T weight = quWeights[q] * md.measure(q);
91 // Compute physical gradients of the velocity basis functions at q as a dim x numActiveFunction matrix
92 transformGradients(md, q, basisValuesVel[1], physGradVel);
93 // matrix A
94 block = weight*viscosity*density * physGradVel.transpose()*physGradVel;
95 for (short_t d = 0; d < dim; ++d)
96 localMat.block(d*N_V,d*N_V,N_V,N_V) += block.block(0,0,N_V,N_V);
97 // matrix B
98 for (short_t d = 0; d < dim; ++d)
99 {
100 block = weight*basisValuesPres.col(q)*physGradVel.row(d);
101 localMat.block(dim*N_V,d*N_V,N_P,N_V) -= block.block(0,0,N_P,N_V);
102 localMat.block(d*N_V,dim*N_V,N_V,N_P) -= block.transpose().block(0,0,N_V,N_P);
103 }
104 // rhs contribution
105 for (short_t d = 0; d < dim; ++d)
106 localRhs.middleRows(d*N_V,N_V).noalias() += weight *density * forceValues(d,q) * basisValuesVel[0].col(q) ;
107 }
108 }
109
110 inline void localToGlobal(const int patchIndex,
111 const std::vector<gsMatrix<T> > & eliminatedDofs,
112 gsSparseSystem<T> & system)
113 {
114 // computes global indices for velocity components
115 for (short_t d = 0; d < dim; ++d)
116 {
117 system.mapColIndices(localIndicesVel, patchIndex, globalIndices[d], d);
118 blockNumbers.at(d) = d;
119 }
120 // computes global indices for pressure
121 system.mapColIndices(localIndicesPres, patchIndex, globalIndices[dim], dim);
122 blockNumbers.at(dim) = dim;
123 // push to global system
124 system.pushToRhs(localRhs,globalIndices,blockNumbers);
125 system.pushToMatrix(localMat,globalIndices,eliminatedDofs,blockNumbers,blockNumbers);
126 }
127
128protected:
129 // problem info
130 short_t dim;
131 const gsBasePde<T> * pde_ptr;
132 T viscosity, density;
133 // geometry mapping
134 gsMapData<T> md;
135 // local components of the global linear system
136 gsMatrix<T> localMat;
137 gsMatrix<T> localRhs;
138 // local indices (at the current patch) of basis functions active at the current element
139 gsMatrix<index_t> localIndicesVel;
140 gsMatrix<index_t> localIndicesPres;
141 // number of velocity and pressure basis functions active at the current element
142 index_t N_V, N_P;
143 // values and derivatives of velocity basis functions at quadrature points at the current element
144 // values are stored as a N_V x numQuadPoints matrix; not sure about derivatives, must be smth like N_V*dim x numQuadPoints
145 std::vector<gsMatrix<T> > basisValuesVel;
146 // values of pressure basis functions active at the current element;
147 // stores as a N_P x numQuadPoints matrix
148 gsMatrix<T> basisValuesPres;
149 // RHS values at quadrature points at the current element; stored as a dim x numQuadPoints matrix
150 gsMatrix<T> forceValues;
151
152 // all temporary matrices defined here for efficiency
153 gsMatrix<T> block, physGradVel;
154 // containers for global indices
155 std::vector< gsMatrix<index_t> > globalIndices;
156 gsVector<index_t> blockNumbers;
157};
158
159} // 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.
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