27 class gsVisitorNavierStokes
31 gsVisitorNavierStokes(
const gsPde<T> & pde_,
const gsMultiPatch<T> & velocity_,
32 const gsMultiPatch<T> & pressure_)
33 : pde_ptr(static_cast<const gsBasePde<T>*>(&pde_)),
35 pressure(pressure_) {}
37 void initialize(
const gsBasisRefs<T> & basisRefs,
39 const gsOptionList & options,
43 dim = basisRefs.front().dim();
48 assemblyType = options.getInt(
"Assembly");
49 viscosity = options.getReal(
"Viscosity");
50 density = options.getReal(
"Density");
52 forceScaling = options.getReal(
"ForceScaling");
54 globalIndices.resize(dim+1);
55 blockNumbers.resize(dim+1);
58 inline void evaluate(
const gsBasisRefs<T> & basisRefs,
59 const gsGeometry<T> & geo,
60 const gsMatrix<T> & quNodes)
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();
77 basisRefs.front().evalAllDers_into(quNodes,1,basisValuesVel);
79 basisRefs.back().eval_into(quNodes,basisValuesPres);
81 basisRefs.back().deriv_into(quNodes,basisGradsPres);
83 pde_ptr->rhs()->eval_into(md.values[0],forceValues);
85 mdVelocity.points = quNodes;
91 velocity.patch(patch).computeMap(mdVelocity);
93 pressure.patch(patch).eval_into(quNodes,pressureValues);
96 inline void assemble(gsDomainIterator<T> & element,
97 const gsVector<T> & quWeights)
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);
107 inline void localToGlobal(
const int patchIndex,
108 const std::vector<gsMatrix<T> > & eliminatedDofs,
109 gsSparseSystem<T> & system)
112 for (
short_t d = 0; d < dim; ++d)
114 system.mapColIndices(localIndicesVel, patchIndex, globalIndices[d], d);
115 blockNumbers.at(d) = d;
118 system.mapColIndices(localIndicesPres, patchIndex, globalIndices[dim], dim);
119 blockNumbers.at(dim) = dim;
121 system.pushToRhs(localRhs,globalIndices,blockNumbers);
122 system.pushToMatrix(localMat,globalIndices,eliminatedDofs,blockNumbers,blockNumbers);
127 void assembleNewtonUpdate(gsDomainIterator<T> & element,
128 const gsVector<T> & quWeights)
131 localMat.setZero(dim*N_V + N_P, dim*N_V + N_P);
132 localRhs.setZero(dim*N_V + N_P,1);
135 for (
index_t q = 0; q < quWeights.rows(); ++q)
138 const T weight = quWeights[q] * md.measure(q);
140 transformGradients(md, q, basisValuesVel[1], physGradVel);
142 physJacCurVel = mdVelocity.jacobian(q)*(md.jacobian(q).cramerInverse());
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);
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);
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);
157 for (
short_t d = 0; d < dim; ++d)
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);
161 localMat.block(d*N_V,dim*N_V,N_V,N_P) -= block.transpose().block(0,0,N_V,N_P);
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);
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();
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);
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();
180 localRhs.middleRows(dim*N_V,N_P).noalias() += weight *
181 basisValuesPres.col(q) * physJacCurVel.trace();
185 void assembleNewtonFull(gsDomainIterator<T> & element,
186 const gsVector<T> & quWeights)
189 localMat.setZero(dim*N_V + N_P, dim*N_V + N_P);
190 localRhs.setZero(dim*N_V + N_P,1);
193 for (
index_t q = 0; q < quWeights.rows(); ++q)
196 const T weight = quWeights[q] * md.measure(q);
198 transformGradients(md, q, basisValuesVel[1], physGradVel);
200 physJacCurVel = mdVelocity.jacobian(q)*(md.jacobian(q).cramerInverse());
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);
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);
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);
215 for (
short_t d = 0; d < dim; ++d)
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);
219 localMat.block(d*N_V,dim*N_V,N_V,N_P) -= block.transpose().block(0,0,N_V,N_P);
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);
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);
232 void assembleOseen(gsDomainIterator<T> & element,
233 const gsVector<T> & quWeights)
236 localMat.setZero(dim*N_V + N_P, dim*N_V + N_P);
237 localRhs.setZero(dim*N_V + N_P,1);
240 for (
index_t q = 0; q < quWeights.rows(); ++q)
243 const T weight = quWeights[q] * md.measure(q);
245 transformGradients(md, q, basisValuesVel[1], physGradVel);
247 physJacCurVel = mdVelocity.jacobian(q)*(md.jacobian(q).cramerInverse());
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);
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);
258 for (
short_t d = 0; d < dim; ++d)
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);
262 localMat.block(d*N_V,dim*N_V,N_V,N_P) -= block.transpose().block(0,0,N_V,N_P);
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);
274 const gsBasePde<T> * pde_ptr;
279 T viscosity, forceScaling, density;
283 gsMatrix<T> localMat;
284 gsMatrix<T> localRhs;
286 gsMatrix<index_t> localIndicesVel;
287 gsMatrix<index_t> localIndicesPres;
293 std::vector<gsMatrix<T> > basisValuesVel;
296 gsMatrix<T> basisValuesPres;
298 gsMatrix<T> basisGradsPres;
300 gsMatrix<T> forceValues;
302 const gsMultiPatch<T> & velocity;
304 gsMapData<T> mdVelocity;
306 const gsMultiPatch<T> & pressure;
308 gsMatrix<T> pressureValues;
310 gsMatrix<T> pressureGrads;
313 gsMatrix<T> block, physGradVel, physJacCurVel;
315 std::vector< gsMatrix<index_t> > globalIndices;
316 gsVector<index_t> blockNumbers;
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.