G+Smo  24.08.0
Geometry + Simulation Modules
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
gsVisitorNavierStokes.h
Go to the documentation of this file.
1 
15 #pragma once
16 
18 #include <gsCore/gsFuncData.h>
19 #include <algorithm>
20 
21 #include <gsElasticity/gsBasePde.h>
22 
23 namespace gismo
24 {
25 
26 template <class T>
27 class gsVisitorNavierStokes
28 {
29 public:
30 
31  gsVisitorNavierStokes(const gsPde<T> & pde_, const gsMultiPatch<T> & velocity_,
32  const gsMultiPatch<T> & pressure_)
33  : pde_ptr(static_cast<const gsBasePde<T>*>(&pde_)),
34  velocity(velocity_),
35  pressure(pressure_) {}
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 velocity component.
45  // the same rule is used for the presure
46  rule = gsQuadrature::get(basisRefs.front(), options);
47  // saving necessary info
48  assemblyType = options.getInt("Assembly");
49  viscosity = options.getReal("Viscosity");
50  density = options.getReal("Density");
51  patch = patchIndex;
52  forceScaling = options.getReal("ForceScaling");
53  // resize containers for global indices
54  globalIndices.resize(dim+1);
55  blockNumbers.resize(dim+1);
56  }
57 
58  inline void evaluate(const gsBasisRefs<T> & basisRefs,
59  const gsGeometry<T> & geo,
60  const gsMatrix<T> & quNodes)
61  {
62  // store quadrature points of the element for geometry evaluation
63  md.points = quNodes;
64  // NEED_VALUE to get points in the physical domain for evaluation of the RHS
65  // NEED_MEASURE to get the Jacobian determinant values for integration
66  // NEED_GRAD_TRANSFORM to get the Jacobian matrix to transform gradient from the parametric to physical domain
67  // NEED_2ND_DER to transform hessians to physical domain
69  // Compute image of the quadrature points plus gradient, jacobian and other necessary data
70  geo.computeMap(md);
71  // find local indices of the velocity and pressure basis functions active on the element
72  basisRefs.front().active_into(quNodes.col(0),localIndicesVel);
73  N_V = localIndicesVel.rows();
74  basisRefs.back().active_into(quNodes.col(0), localIndicesPres);
75  N_P = localIndicesPres.rows();
76  // Evaluate velocity basis functions and their derivatives on the element (and hessians, if SUPG is used)
77  basisRefs.front().evalAllDers_into(quNodes,1,basisValuesVel);
78  // Evaluate pressure basis functions on the element
79  basisRefs.back().eval_into(quNodes,basisValuesPres);
80  // Evaluate gradients of pressure basis functions if SUPG is used
81  basisRefs.back().deriv_into(quNodes,basisGradsPres);
82  // Evaluate right-hand side at the image of the quadrature points
83  pde_ptr->rhs()->eval_into(md.values[0],forceValues);
84  // store quadrature points of the element for velocity evaluation
85  mdVelocity.points = quNodes;
86  // NEED_VALUE to compute velocity values
87  // NEED_DERIV to compute velocity gradient
88  // NEED_2ND_DER to compute velocity hessian
89  mdVelocity.flags = NEED_VALUE | NEED_DERIV;
90  // evaluate velocity
91  velocity.patch(patch).computeMap(mdVelocity);
92  // evaluate pressure values; for some reason, pressure evaluation via gsMapData doesn't work
93  pressure.patch(patch).eval_into(quNodes,pressureValues);
94  }
95 
96  inline void assemble(gsDomainIterator<T> & element,
97  const gsVector<T> & quWeights)
98  {
99  if (assemblyType == 0)
100  assembleOseen(element,quWeights);
101  else if (assemblyType == 1)
102  assembleNewtonUpdate(element,quWeights);
103  else if (assemblyType == 2)
104  assembleNewtonFull(element,quWeights);
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 
127  void assembleNewtonUpdate(gsDomainIterator<T> & element,
128  const gsVector<T> & quWeights)
129  {
130  // Initialize local matrix/rhs // A | D
131  localMat.setZero(dim*N_V + N_P, dim*N_V + N_P); // --|-- matrix structure
132  localRhs.setZero(dim*N_V + N_P,1); // B | 0// roughly estimate h - diameter of the element ( for SUPG)
133  // T h = cellSize(element);
134  // Loop over the quadrature nodes
135  for (index_t q = 0; q < quWeights.rows(); ++q)
136  {
137  // Multiply quadrature weight by the geometry measure
138  const T weight = quWeights[q] * md.measure(q);
139  // Compute physical gradients of the velocity basis functions at q as a dim x numActiveFunction matrix
140  transformGradients(md, q, basisValuesVel[1], physGradVel);
141  // Compute physical Jacobian of the current velocity field
142  physJacCurVel = mdVelocity.jacobian(q)*(md.jacobian(q).cramerInverse());
143  // matrix A: diffusion
144  block = weight*density*viscosity * physGradVel.transpose()*physGradVel;
145  for (short_t d = 0; d < dim; ++d)
146  localMat.block(d*N_V,d*N_V,N_V,N_V) += block.block(0,0,N_V,N_V);
147  // matrix A: advection
148  block = weight*basisValuesVel[0].col(q) * (mdVelocity.values[0].col(q).transpose()*physGradVel);
149  for (short_t d = 0; d < dim; ++d)
150  localMat.block(d*N_V,d*N_V,N_V,N_V) += density*block.block(0,0,N_V,N_V);
151  // matrix A: reaction
152  block = weight*density*basisValuesVel[0].col(q) * basisValuesVel[0].col(q).transpose();
153  for (short_t di = 0; di < dim; ++di)
154  for (short_t dj = 0; dj < dim; ++dj)
155  localMat.block(di*N_V,dj*N_V,N_V,N_V) += physJacCurVel(di,dj)*block.block(0,0,N_V,N_V);
156  // matrices B and D
157  for (short_t d = 0; d < dim; ++d)
158  {
159  block = weight*basisValuesPres.col(q)*physGradVel.row(d);
160  localMat.block(dim*N_V,d*N_V,N_P,N_V) -= block.block(0,0,N_P,N_V); // B
161  localMat.block(d*N_V,dim*N_V,N_V,N_P) -= block.transpose().block(0,0,N_V,N_P); // D
162  }
163  // rhs: force
164  for (short_t d = 0; d < dim; ++d)
165  localRhs.middleRows(d*N_V,N_V).noalias() += weight *density* forceScaling *
166  forceValues(d,q) * basisValuesVel[0].col(q);
167  // rhs: residual diffusion
168  for (short_t d = 0; d < dim; ++d)
169  localRhs.middleRows(d*N_V,N_V).noalias() -= weight * viscosity * density*
170  (physJacCurVel.row(d)*physGradVel).transpose();
171  // rhs: residual nonlinear
172  for (short_t d = 0; d < dim; ++d)
173  localRhs.middleRows(d*N_V,N_V).noalias() -= weight * density*
174  (physJacCurVel.row(d) * mdVelocity.values[0].col(q))(0,0) * basisValuesVel[0].col(q);
175  // rhs: residual pressure
176  for (short_t d = 0; d < dim; ++d)
177  localRhs.middleRows(d*N_V,N_V).noalias() += weight *
178  pressureValues.at(q) * physGradVel.row(d).transpose();
179  // rhs: constraint residual
180  localRhs.middleRows(dim*N_V,N_P).noalias() += weight *
181  basisValuesPres.col(q) * physJacCurVel.trace();
182  }
183  }
184 
185  void assembleNewtonFull(gsDomainIterator<T> & element,
186  const gsVector<T> & quWeights)
187  {
188  // Initialize local matrix/rhs // A | D
189  localMat.setZero(dim*N_V + N_P, dim*N_V + N_P); // --|-- matrix structure
190  localRhs.setZero(dim*N_V + N_P,1); // B | 0// roughly estimate h - diameter of the element ( for SUPG)
191  // T h = cellSize(element);
192  // Loop over the quadrature nodes
193  for (index_t q = 0; q < quWeights.rows(); ++q)
194  {
195  // Multiply quadrature weight by the geometry measure
196  const T weight = quWeights[q] * md.measure(q);
197  // Compute physical gradients of the velocity basis functions at q as a dim x numActiveFunction matrix
198  transformGradients(md, q, basisValuesVel[1], physGradVel);
199  // Compute physical Jacobian of the current velocity field
200  physJacCurVel = mdVelocity.jacobian(q)*(md.jacobian(q).cramerInverse());
201  // matrix A: diffusion
202  block = weight*density*viscosity * physGradVel.transpose()*physGradVel;
203  for (short_t d = 0; d < dim; ++d)
204  localMat.block(d*N_V,d*N_V,N_V,N_V) += block.block(0,0,N_V,N_V);
205  // matrix A: advection
206  block = weight*basisValuesVel[0].col(q) * (mdVelocity.values[0].col(q).transpose()*physGradVel);
207  for (short_t d = 0; d < dim; ++d)
208  localMat.block(d*N_V,d*N_V,N_V,N_V) += density*block.block(0,0,N_V,N_V);
209  // matrix A: reaction
210  block = weight*density*basisValuesVel[0].col(q) * basisValuesVel[0].col(q).transpose();
211  for (short_t di = 0; di < dim; ++di)
212  for (short_t dj = 0; dj < dim; ++dj)
213  localMat.block(di*N_V,dj*N_V,N_V,N_V) += physJacCurVel(di,dj)*block.block(0,0,N_V,N_V);
214  // matrices B and D
215  for (short_t d = 0; d < dim; ++d)
216  {
217  block = weight*basisValuesPres.col(q)*physGradVel.row(d);
218  localMat.block(dim*N_V,d*N_V,N_P,N_V) -= block.block(0,0,N_P,N_V); // B
219  localMat.block(d*N_V,dim*N_V,N_V,N_P) -= block.transpose().block(0,0,N_V,N_P); // D
220  }
221  // rhs: force
222  for (short_t d = 0; d < dim; ++d)
223  localRhs.middleRows(d*N_V,N_V).noalias() += weight *density* forceScaling *
224  forceValues(d,q) * basisValuesVel[0].col(q);
225  // rhs: residual nonlinear
226  for (short_t d = 0; d < dim; ++d)
227  localRhs.middleRows(d*N_V,N_V).noalias() += weight * density*
228  (physJacCurVel.row(d) * mdVelocity.values[0].col(q))(0,0) * basisValuesVel[0].col(q);
229  }
230  }
231 
232  void assembleOseen(gsDomainIterator<T> & element,
233  const gsVector<T> & quWeights)
234  {
235  // Initialize local matrix/rhs // A | D
236  localMat.setZero(dim*N_V + N_P, dim*N_V + N_P); // --|-- matrix structure
237  localRhs.setZero(dim*N_V + N_P,1); // B | 0
238 
239  // Loop over the quadrature nodes
240  for (index_t q = 0; q < quWeights.rows(); ++q)
241  {
242  // Multiply quadrature weight by the geometry measure
243  const T weight = quWeights[q] * md.measure(q);
244  // Compute physical gradients of the velocity basis functions at q as a dim x numActiveFunction matrix
245  transformGradients(md, q, basisValuesVel[1], physGradVel);
246  // Compute physical Jacobian of the current velocity field
247  physJacCurVel = mdVelocity.jacobian(q)*(md.jacobian(q).cramerInverse());
248  // matrix A: diffusion
249  block = weight*viscosity *density* physGradVel.transpose()*physGradVel;
250  for (short_t d = 0; d < dim; ++d)
251  localMat.block(d*N_V,d*N_V,N_V,N_V) += block.block(0,0,N_V,N_V);
252  // matrix A: advection
253  block = weight*basisValuesVel[0].col(q) *
254  (mdVelocity.values[0].col(q).transpose()*physGradVel);
255  for (short_t d = 0; d < dim; ++d)
256  localMat.block(d*N_V,d*N_V,N_V,N_V) += density*block.block(0,0,N_V,N_V);
257  // matrices B and D
258  for (short_t d = 0; d < dim; ++d)
259  {
260  block = weight*basisValuesPres.col(q)*physGradVel.row(d);
261  localMat.block(dim*N_V,d*N_V,N_P,N_V) -= block.block(0,0,N_P,N_V); // B
262  localMat.block(d*N_V,dim*N_V,N_V,N_P) -= block.transpose().block(0,0,N_V,N_P); // D
263  }
264  // rhs: force
265  for (short_t d = 0; d < dim; ++d)
266  localRhs.middleRows(d*N_V,N_V).noalias() += weight *density* forceScaling *
267  forceValues(d,q) * basisValuesVel[0].col(q);
268  }
269  }
270 
271 protected:
272  // problem info
273  short_t dim;
274  const gsBasePde<T> * pde_ptr;
275  // switch between assembling different linear systems (Newton or Oseen iterations)
276  index_t assemblyType;
277  index_t patch; // current patch
278  // constants: viscosity and the force scaling factor
279  T viscosity, forceScaling, density;
280  // geometry mapping
281  gsMapData<T> md;
282  // local components of the global linear system
283  gsMatrix<T> localMat;
284  gsMatrix<T> localRhs;
285  // local indices (at the current patch) of basis functions active at the current element
286  gsMatrix<index_t> localIndicesVel;
287  gsMatrix<index_t> localIndicesPres;
288  // number of velocity and pressure basis functions active at the current element
289  index_t N_V, N_P;
290  // values and derivatives of velocity basis functions at quadrature points at the current element
291  // values are stored as a N_V x numQuadPoints matrix; not sure about derivatives, must be smth like N_V*dim x numQuadPoints
292  // if supg is true, also stores the hessian
293  std::vector<gsMatrix<T> > basisValuesVel;
294  // values of pressure basis functions active at the current element;
295  // stores as a N_P x numQuadPoints matrix
296  gsMatrix<T> basisValuesPres;
297  // if supg is true, stores the derivatives
298  gsMatrix<T> basisGradsPres;
299  // RHS values at quadrature points at the current element; stored as a dim x numQuadPoints matrix
300  gsMatrix<T> forceValues;
301  // current displacement field
302  const gsMultiPatch<T> & velocity;
303  // evaluation data of the current velocity field
304  gsMapData<T> mdVelocity;
305  // current pressure field
306  const gsMultiPatch<T> & pressure;
307  // pressure values at the current element; stored as a 1 x numQuadPoints matrix
308  gsMatrix<T> pressureValues;
309  // pressure gradients at the current element (only for supg); stored as a dim x numQuadPoints matrix
310  gsMatrix<T> pressureGrads;
311 
312  // all temporary matrices defined here for efficiency
313  gsMatrix<T> block, physGradVel, physJacCurVel;
314  // containers for global indices
315  std::vector< gsMatrix<index_t> > globalIndices;
316  gsVector<index_t> blockNumbers;
317 };
318 
319 } // namespace gismo
Gradient of the object.
Definition: gsForwardDeclarations.h:73
#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.