27class gsVisitorNavierStokes
31 gsVisitorNavierStokes(
const gsPde<T> & pde_,
const gsMultiPatch<T> & velocity_,
32 const gsMultiPatch<T> & pressure_)
33 : dim(), pde_ptr(static_cast<const gsBasePde<T>*>(&pde_)), assemblyType(), patch(),
34 viscosity(), forceScaling(), density(),
35 N_V(), N_P(), velocity(velocity_), pressure(pressure_)
38 void initialize(
const gsBasisRefs<T> & basisRefs,
40 const gsOptionList & options,
45 dim = basisRefs.front().dim();
50 assemblyType = options.getInt(
"Assembly");
51 viscosity = options.getReal(
"Viscosity");
52 density = options.getReal(
"Density");
54 forceScaling = options.getReal(
"ForceScaling");
56 globalIndices.resize(dim+1);
57 blockNumbers.resize(dim+1);
60 inline void evaluate(
const gsBasisRefs<T> & basisRefs,
61 const gsGeometry<T> & geo,
62 const gsMatrix<T> & quNodes)
74 basisRefs.front().active_into(quNodes.col(0),localIndicesVel);
75 N_V = localIndicesVel.rows();
76 basisRefs.back().active_into(quNodes.col(0), localIndicesPres);
77 N_P = localIndicesPres.rows();
79 basisRefs.front().evalAllDers_into(quNodes,1,basisValuesVel);
81 basisRefs.back().eval_into(quNodes,basisValuesPres);
83 basisRefs.back().deriv_into(quNodes,basisGradsPres);
85 pde_ptr->rhs()->eval_into(md.values[0],forceValues);
87 mdVelocity.points = quNodes;
93 velocity.patch(patch).computeMap(mdVelocity);
95 pressure.patch(patch).eval_into(quNodes,pressureValues);
98 inline void assemble(gsDomainIterator<T> & element,
99 const gsVector<T> & quWeights)
102 if (assemblyType == 0)
103 assembleOseen(element,quWeights);
104 else if (assemblyType == 1)
105 assembleNewtonUpdate(element,quWeights);
106 else if (assemblyType == 2)
107 assembleNewtonFull(element,quWeights);
110 inline void localToGlobal(
const int patchIndex,
111 const std::vector<gsMatrix<T> > & eliminatedDofs,
112 gsSparseSystem<T> & system)
115 for (
short_t d = 0; d < dim; ++d)
117 system.mapColIndices(localIndicesVel, patchIndex, globalIndices[d], d);
118 blockNumbers.at(d) = d;
121 system.mapColIndices(localIndicesPres, patchIndex, globalIndices[dim], dim);
122 blockNumbers.at(dim) = dim;
124 system.pushToRhs(localRhs,globalIndices,blockNumbers);
125 system.pushToMatrix(localMat,globalIndices,eliminatedDofs,blockNumbers,blockNumbers);
130 void assembleNewtonUpdate(gsDomainIterator<T> & element,
131 const gsVector<T> & quWeights)
135 localMat.setZero(dim*N_V + N_P, dim*N_V + N_P);
136 localRhs.setZero(dim*N_V + N_P,1);
139 for (
index_t q = 0; q < quWeights.rows(); ++q)
142 const T weight = quWeights[q] * md.measure(q);
144 transformGradients(md, q, basisValuesVel[1], physGradVel);
146 physJacCurVel = mdVelocity.jacobian(q)*(md.jacobian(q).cramerInverse());
148 block = weight*density*viscosity * physGradVel.transpose()*physGradVel;
149 for (
short_t d = 0; d < dim; ++d)
150 localMat.block(d*N_V,d*N_V,N_V,N_V) += block.block(0,0,N_V,N_V);
152 block = weight*basisValuesVel[0].col(q) * (mdVelocity.values[0].col(q).transpose()*physGradVel);
153 for (
short_t d = 0; d < dim; ++d)
154 localMat.block(d*N_V,d*N_V,N_V,N_V) += density*block.block(0,0,N_V,N_V);
156 block = weight*density*basisValuesVel[0].col(q) * basisValuesVel[0].col(q).transpose();
157 for (
short_t di = 0; di < dim; ++di)
158 for (
short_t dj = 0; dj < dim; ++dj)
159 localMat.block(di*N_V,dj*N_V,N_V,N_V) += physJacCurVel(di,dj)*block.block(0,0,N_V,N_V);
161 for (
short_t d = 0; d < dim; ++d)
163 block = weight*basisValuesPres.col(q)*physGradVel.row(d);
164 localMat.block(dim*N_V,d*N_V,N_P,N_V) -= block.block(0,0,N_P,N_V);
165 localMat.block(d*N_V,dim*N_V,N_V,N_P) -= block.transpose().block(0,0,N_V,N_P);
168 for (
short_t d = 0; d < dim; ++d)
169 localRhs.middleRows(d*N_V,N_V).noalias() += weight *density* forceScaling *
170 forceValues(d,q) * basisValuesVel[0].col(q);
172 for (
short_t d = 0; d < dim; ++d)
173 localRhs.middleRows(d*N_V,N_V).noalias() -= weight * viscosity * density*
174 (physJacCurVel.row(d)*physGradVel).transpose();
176 for (
short_t d = 0; d < dim; ++d)
177 localRhs.middleRows(d*N_V,N_V).noalias() -= weight * density*
178 (physJacCurVel.row(d) * mdVelocity.values[0].col(q))(0,0) * basisValuesVel[0].col(q);
180 for (
short_t d = 0; d < dim; ++d)
181 localRhs.middleRows(d*N_V,N_V).noalias() += weight *
182 pressureValues.at(q) * physGradVel.row(d).transpose();
184 localRhs.middleRows(dim*N_V,N_P).noalias() += weight *
185 basisValuesPres.col(q) * physJacCurVel.trace();
189 void assembleNewtonFull(gsDomainIterator<T> & element,
190 const gsVector<T> & quWeights)
194 localMat.setZero(dim*N_V + N_P, dim*N_V + N_P);
195 localRhs.setZero(dim*N_V + N_P,1);
198 for (
index_t q = 0; q < quWeights.rows(); ++q)
201 const T weight = quWeights[q] * md.measure(q);
203 transformGradients(md, q, basisValuesVel[1], physGradVel);
205 physJacCurVel = mdVelocity.jacobian(q)*(md.jacobian(q).cramerInverse());
207 block = weight*density*viscosity * physGradVel.transpose()*physGradVel;
208 for (
short_t d = 0; d < dim; ++d)
209 localMat.block(d*N_V,d*N_V,N_V,N_V) += block.block(0,0,N_V,N_V);
211 block = weight*basisValuesVel[0].col(q) * (mdVelocity.values[0].col(q).transpose()*physGradVel);
212 for (
short_t d = 0; d < dim; ++d)
213 localMat.block(d*N_V,d*N_V,N_V,N_V) += density*block.block(0,0,N_V,N_V);
215 block = weight*density*basisValuesVel[0].col(q) * basisValuesVel[0].col(q).transpose();
216 for (
short_t di = 0; di < dim; ++di)
217 for (
short_t dj = 0; dj < dim; ++dj)
218 localMat.block(di*N_V,dj*N_V,N_V,N_V) += physJacCurVel(di,dj)*block.block(0,0,N_V,N_V);
220 for (
short_t d = 0; d < dim; ++d)
222 block = weight*basisValuesPres.col(q)*physGradVel.row(d);
223 localMat.block(dim*N_V,d*N_V,N_P,N_V) -= block.block(0,0,N_P,N_V);
224 localMat.block(d*N_V,dim*N_V,N_V,N_P) -= block.transpose().block(0,0,N_V,N_P);
227 for (
short_t d = 0; d < dim; ++d)
228 localRhs.middleRows(d*N_V,N_V).noalias() += weight *density* forceScaling *
229 forceValues(d,q) * basisValuesVel[0].col(q);
231 for (
short_t d = 0; d < dim; ++d)
232 localRhs.middleRows(d*N_V,N_V).noalias() += weight * density*
233 (physJacCurVel.row(d) * mdVelocity.values[0].col(q))(0,0) * basisValuesVel[0].col(q);
237 void assembleOseen(gsDomainIterator<T> & element,
238 const gsVector<T> & quWeights)
242 localMat.setZero(dim*N_V + N_P, dim*N_V + N_P);
243 localRhs.setZero(dim*N_V + N_P,1);
246 for (
index_t q = 0; q < quWeights.rows(); ++q)
249 const T weight = quWeights[q] * md.measure(q);
251 transformGradients(md, q, basisValuesVel[1], physGradVel);
253 physJacCurVel = mdVelocity.jacobian(q)*(md.jacobian(q).cramerInverse());
255 block = weight*viscosity *density* physGradVel.transpose()*physGradVel;
256 for (
short_t d = 0; d < dim; ++d)
257 localMat.block(d*N_V,d*N_V,N_V,N_V) += block.block(0,0,N_V,N_V);
259 block = weight*basisValuesVel[0].col(q) *
260 (mdVelocity.values[0].col(q).transpose()*physGradVel);
261 for (
short_t d = 0; d < dim; ++d)
262 localMat.block(d*N_V,d*N_V,N_V,N_V) += density*block.block(0,0,N_V,N_V);
264 for (
short_t d = 0; d < dim; ++d)
266 block = weight*basisValuesPres.col(q)*physGradVel.row(d);
267 localMat.block(dim*N_V,d*N_V,N_P,N_V) -= block.block(0,0,N_P,N_V);
268 localMat.block(d*N_V,dim*N_V,N_V,N_P) -= block.transpose().block(0,0,N_V,N_P);
271 for (
short_t d = 0; d < dim; ++d)
272 localRhs.middleRows(d*N_V,N_V).noalias() += weight *density* forceScaling *
273 forceValues(d,q) * basisValuesVel[0].col(q);
280 const gsBasePde<T> * pde_ptr;
285 T viscosity, forceScaling, density;
289 gsMatrix<T> localMat;
290 gsMatrix<T> localRhs;
292 gsMatrix<index_t> localIndicesVel;
293 gsMatrix<index_t> localIndicesPres;
299 std::vector<gsMatrix<T> > basisValuesVel;
302 gsMatrix<T> basisValuesPres;
304 gsMatrix<T> basisGradsPres;
306 gsMatrix<T> forceValues;
308 const gsMultiPatch<T> & velocity;
310 gsMapData<T> mdVelocity;
312 const gsMultiPatch<T> & pressure;
314 gsMatrix<T> pressureValues;
316 gsMatrix<T> pressureGrads;
319 gsMatrix<T> block, physGradVel, physJacCurVel;
321 std::vector< gsMatrix<index_t> > globalIndices;
322 gsVector<index_t> blockNumbers;
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_DERIV
Gradient of the object.
Definition gsForwardDeclarations.h:73
@ 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