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