Loading...
Searching...
No Matches
steadyNS.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ██╗████████╗██╗ ██╗ █████╗ ██████╗ █████╗ ███████╗██╗ ██╗
3 ██║╚══██╔══╝██║ ██║██╔══██╗██╔════╝██╔══██╗ ██╔════╝██║ ██║
4 ██║ ██║ ███████║███████║██║ ███████║█████╗█████╗ ██║ ██║
5 ██║ ██║ ██╔══██║██╔══██║██║ ██╔══██║╚════╝██╔══╝ ╚██╗ ██╔╝
6 ██║ ██║ ██║ ██║██║ ██║╚██████╗██║ ██║ ██║ ╚████╔╝
7 ╚═╝ ╚═╝ ╚═╝ ╚═╝╚═╝ ╚═╝ ╚═════╝╚═╝ ╚═╝ ╚═╝ ╚═══╝
8
9 * In real Time Highly Advanced Computational Applications for Finite Volumes
10 * Copyright (C) 2017 by the ITHACA-FV authors
11-------------------------------------------------------------------------------
12License
13 This file is part of ITHACA-FV
14 ITHACA-FV is free software: you can redistribute it and/or modify
15 it under the terms of the GNU Lesser General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18 ITHACA-FV is distributed in the hope that it will be useful,
19 but WITHOUT ANY WARRANTY; without even the implied warranty of
20 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21 GNU Lesser General Public License for more details.
22 You should have received a copy of the GNU Lesser General Public License
23 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
24Class
25 steadyNS
26Description
27 Reduction class for a Stationary Navier-Stokes problem.
28SourceFiles
29 steadyNS.C
30\*---------------------------------------------------------------------------*/
31
36
37#ifndef steadyNS_H
38#define steadyNS_H
39#include "fvCFD.H"
40#include "singlePhaseTransportModel.H"
41#include "turbulentTransportModel.H"
42#include "simpleControl.H"
43#include "pisoControl.H"
44#include "fvOptions.H"
45#include "reductionProblem.H"
46#include "ITHACAstream.H"
47#include "ITHACAparameters.H"
48#if OPENFOAM >= 1812
49#include "ITHACAforces18.H"
50#else
51#include "ITHACAforces.H"
52#endif
53#include "volFields.H"
54#include <iostream>
55#include "IPstream.H"
56#include "OPstream.H"
57#include "Modes.H"
58
59// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60
61/*---------------------------------------------------------------------------*\
62 Class SteadyNS Declaration
63\*---------------------------------------------------------------------------*/
64
66
70{
71
72
73 public:
74 // Constructors
76 steadyNS();
77
79 steadyNS(int argc, char* argv[]);
81
83
84 // Member Functions
86 PtrList<volScalarField> Pfield;
87
89 PtrList<volVectorField> Ufield;
90
92 PtrList<volVectorField> supfield;
93
95 PtrList<surfaceScalarField> Phifield;
96
99
102
105
108
110 PtrList<volVectorField> liftfield;
111
113 PtrList<volVectorField> Uomfield;
114
117
120
122 scalar tolerance;
123
125 scalar maxIter;
126
129
132
135
138
140 label NUmodes;
141
143 label NPmodes;
144
147
150
155
157 Eigen::MatrixXd B_matrix;
158
160 Eigen::MatrixXd M_matrix;
161
163 Eigen::MatrixXd K_matrix;
164
166 List <Eigen::MatrixXd> C_matrix;
167 Eigen::Tensor<double, 3 > C_tensor;
168
170 Eigen::MatrixXd P_matrix;
171
173 Eigen::MatrixXd D_matrix;
174
176 List <Eigen::MatrixXd> G_matrix;
177
179 Eigen::Tensor<double, 3 > gTensor;
180
182 Eigen::MatrixXd BC1_matrix;
183
185 List <Eigen::MatrixXd> BC2_matrix;
186
188 Eigen::Tensor<double, 3 > bc2Tensor;
189
191 Eigen::MatrixXd BC3_matrix;
192
194 Eigen::MatrixXd BC4_matrix;
195
197 Eigen::MatrixXd W_matrix;
198
200 Eigen::MatrixXd I_matrix;
201
203 Eigen::MatrixXd DF_matrix;
204
206 Eigen::MatrixXd KF_matrix;
207
209 //
210
212 Eigen::MatrixXd tauMatrix;
213
215 Eigen::MatrixXd nMatrix;
216
218 List <Eigen::MatrixXd> bcVelVec;
219
221 List <Eigen::MatrixXd> bcVelMat;
222
224 Eigen::MatrixXd BP_matrix;
225
227 List<Eigen::MatrixXd> RD_matrix;
228
230 List<Eigen::MatrixXd> RC_matrix;
231
233 List<Eigen::MatrixXd> SD_matrix;
234
236 List<Eigen::MatrixXd> SC_matrix;
237
239 Eigen::Tensor<double, 3 > Cf_tensor;
240
242 Eigen::Tensor<double, 3 > Ci_tensor;
243
245 List<Eigen::MatrixXd> LinSysDiv;
246
248 List<Eigen::MatrixXd> LinSysDiff;
249
251 List<Eigen::MatrixXd> LinSysConv;
252
253
254 // Other Variables
256 bool supex;
257
258 // Dummy variables to transform simplefoam into a class
260 autoPtr<volScalarField> _p;
261
263 autoPtr<volVectorField> _U;
264
266 autoPtr<volScalarField> _p0;
267
269 autoPtr<volVectorField> _U0;
270
272 autoPtr<volVectorField> Uinl;
273
275 autoPtr<dimensionedScalar> dt_dummy;
276
278 autoPtr<dimensionedScalar> nu_dummy;
279
281 mutable autoPtr<fvMesh> _mesh;
282
284 autoPtr<simpleControl> _simple;
285
287 autoPtr<fv::options> _fvOptions;
288
290 autoPtr<Time> _runTime;
291
293 autoPtr<surfaceScalarField> _phi;
294
296 autoPtr<surfaceScalarField> _phi0;
297
299 autoPtr<incompressible::turbulenceModel> turbulence;
300
302 autoPtr<singlePhaseTransportModel> _laminarTransport;
303
305 autoPtr<IOMRFZoneList> _MRF;
306
308 label pRefCell;
309
311 scalar pRefValue;
312
315
318
321
322 // Functions
323
324 //--------------------------------------------------------------------------
330 void truthSolve(List<scalar> mu_now);
331
337 void solvesupremizer(word type = "snapshots");
338
340 void liftSolve();
341
342 // Wrapped Proj. Methods;
343
344 //--------------------------------------------------------------------------
352 void projectPPE(fileName folder, label NUmodes, label NPmodes,
353 label NSUPmodes = 0);
354
355 //--------------------------------------------------------------------------
363 void projectSUP(fileName folder, label NUmodes, label NPmodes, label NSUPmodes);
364
365 //--------------------------------------------------------------------------
373 void discretizeThenProject(fileName folder, label NUmodes, label NPmodes,
374 label NSUPmodes = 0);
375
376 //--------------------------------------------------------------------------
377 // Projection Methods Momentum Equation
386 Eigen::MatrixXd diffusive_term(label NUmodes, label NPmodes, label NSUPmodes);
387
388 //--------------------------------------------------------------------------
389 // Projection Methods Momentum Equation
398 Eigen::MatrixXd diffusive_term_sym(label NUmodes, label NPmodes, label NSUPmodes);
399
400 //--------------------------------------------------------------------------
409 Eigen::MatrixXd pressure_gradient_term(label NUmodes, label NPmodes,
410 label NSUPmodes);
411
412 //--------------------------------------------------------------------------
421 List < Eigen::MatrixXd > convective_term(label NUmodes, label NPmodes,
422 label NSUPmodes);
423
424 //--------------------------------------------------------------------------
433 Eigen::MatrixXd mass_term(label NUmodes, label NPmodes, label NSUPmodes);
434
435 // Projection Methods Continuity Equation
436
437 //--------------------------------------------------------------------------
446 Eigen::MatrixXd divergence_term(label NUmodes, label NPmodes, label NSUPmodes);
447
448 //--------------------------------------------------------------------------
456 List < Eigen::MatrixXd > div_momentum(label NUmodes, label NPmodes);
457
466 Eigen::Tensor<double, 3 > divMomentum(label NUmodes, label NPmodes);
467
468 //--------------------------------------------------------------------------
475 Eigen::MatrixXd laplacian_pressure(label NPmodes);
476
477 //--------------------------------------------------------------------------
485 Eigen::MatrixXd pressure_BC1(label NPmodes, label NUmodes);
486
487 //--------------------------------------------------------------------------
495 List < Eigen::MatrixXd > pressure_BC2(label NPmodes, label NUmodes);
496
497 //--------------------------------------------------------------------------
505 Eigen::Tensor<double, 3 > pressureBC2(label NPmodes, label NUmodes);
506
507 //--------------------------------------------------------------------------
515 Eigen::MatrixXd pressure_BC3(label NPmodes, label NUmodes);
516
517 //--------------------------------------------------------------------------
526 Eigen::MatrixXd pressure_BC4(label NPmodes, label NUmodes);
527
528 //--------------------------------------------------------------------------
536 List< Eigen::MatrixXd > bcVelocityVec(label NUmodes, label NSUPmodes);
537
538 //--------------------------------------------------------------------------
546 List< Eigen::MatrixXd > bcVelocityMat(label NUmodes, label NSUPmodes);
547
548 //--------------------------------------------------------------------------
549 // Projection Methods Momentum Equation
558 Eigen::MatrixXd diffusive_term_flux_method(label NUmodes, label NPmodes,
559 label NSUPmodes);
560
561 //--------------------------------------------------------------------------
570 List<Eigen::MatrixXd> boundary_vector_diffusion(label NUmodes, label NPmodes,
571 label NSUPmodes);
572
573 //--------------------------------------------------------------------------
582 List<Eigen::MatrixXd> boundary_vector_convection(label NUmodes, label NPmodes,
583 label NSUPmodes);
584
585 //--------------------------------------------------------------------------
594 Eigen::Tensor<double, 3 > convective_term_flux_tens(label NUmodes,
595 label NPmodes,
596 label NSUPmodes);
597
598 //--------------------------------------------------------------------------
607 List < Eigen::MatrixXd > pressure_gradient_term_linsys_div(label NPmodes);
608
609 //--------------------------------------------------------------------------
618 List < Eigen::MatrixXd > pressure_gradient_term_linsys_diff(label NPmodes);
619
620 //--------------------------------------------------------------------------
629 List < Eigen::MatrixXd > pressure_gradient_term_linsys_conv(label NPmodes);
630
631 //--------------------------------------------------------------------------
632 // Projection Methods Flux Equation
641 Eigen::MatrixXd diffusive_term_consistent(label NUmodes, label NPmodes,
642 label NSUPmodes);
643
644 //--------------------------------------------------------------------------
652 List < Eigen::MatrixXd > boundary_vector_diffusion_consistent(label NUmodes,
653 label NSUPmodes);
654
655 //--------------------------------------------------------------------------
663 List < Eigen::MatrixXd > boundary_vector_convection_consistent(label NUmodes,
664 label NSUPmodes);
665
666 //--------------------------------------------------------------------------
675 Eigen::MatrixXd mass_matrix_newtime_consistent(label NUmodes, label NPmodes,
676 label NSUPmodes);
677
678 //--------------------------------------------------------------------------
687 Eigen::MatrixXd mass_matrix_oldtime_consistent(label NUmodes, label NPmodes,
688 label NSUPmodes);
689
690 //--------------------------------------------------------------------------
699 Eigen::MatrixXd pressure_gradient_term_consistent(label NUmodes, label NPmodes,
700 label NSUPmodes);
701 //--------------------------------------------------------------------------
710 Eigen::Tensor<double, 3 > convective_term_consistent_tens(label NUmodes,
711 label NPmodes,
712 label NSUPmodes);
713
714 //--------------------------------------------------------------------------
719 void change_viscosity(double mu);
720
721 //--------------------------------------------------------------------------
728 void forcesMatrices(label NUmodes, label NPmodes, label NSUPmodes);
729
730 //--------------------------------------------------------------------------
736 void forcesMatrices(label nModes);
737
738
739 //--------------------------------------------------------------------------
746 void reconstructLiftAndDrag(const Eigen::MatrixXd& velCoeffs,
747 const Eigen::MatrixXd& pressureCoeffs, fileName folder);
748
749 //--------------------------------------------------------------------------
758 Eigen::Tensor<double, 3 > convective_term_tens(label NUmodes,
759 label NPmodes,
760 label NSUPmodes);
761
763 void restart();
764
765
766};
767
768#endif
769
770
771
772
773
774
775
776
777
778
Header file of the ITHACAstream class, it contains the implementation of several methods for input ou...
Header file of the Modes class.
Class for the definition of some general parameters, the parameters must be defined from the file ITH...
A general class for the implementation of a full order parametrized problem.
Eigen::MatrixXd mu
Row matrix of parameters.
void truthSolve()
Perform a TruthSolve.
Implementation of a parametrized full order steady NS problem and preparation of the the reduced ma...
Definition steadyNS.H:70
void change_viscosity(double mu)
Function to change the viscosity.
Definition steadyNS.C:2014
List< Eigen::MatrixXd > bcVelVec
Boundary term for penalty method - vector.
Definition steadyNS.H:218
scalar maxIter
Number of maximum iterations to be done for the computation of the truth solution.
Definition steadyNS.H:125
void restart()
set U and P back to the values into the 0 folder
Definition steadyNS.C:2253
List< Eigen::MatrixXd > pressure_gradient_term_linsys_div(label NPmodes)
Laplacian of pressure Linear System - Divergence term.
Definition steadyNS.C:1671
label NPmodes
Number of pressure modes used for the projection.
Definition steadyNS.H:143
List< Eigen::MatrixXd > pressure_gradient_term_linsys_conv(label NPmodes)
Laplacian of pressure Linear System - Convection term.
Definition steadyNS.C:1703
bool supex
Boolean variable to check the existence of the supremizer modes.
Definition steadyNS.H:256
void projectPPE(fileName folder, label NUmodes, label NPmodes, label NSUPmodes=0)
Project using the Poisson Equation for pressure.
Definition steadyNS.C:353
Eigen::MatrixXd nMatrix
Pressure forces.
Definition steadyNS.H:215
Eigen::MatrixXd BC1_matrix
PPE BC1.
Definition steadyNS.H:182
Eigen::MatrixXd diffusive_term(label NUmodes, label NPmodes, label NSUPmodes)
Diffusive Term.
Definition steadyNS.C:861
autoPtr< surfaceScalarField > _phi
Flux.
Definition steadyNS.H:293
Eigen::MatrixXd BC4_matrix
PPE BC4.
Definition steadyNS.H:194
Eigen::MatrixXd tauMatrix
Viscous forces.
Definition steadyNS.H:212
Eigen::MatrixXd BC3_matrix
PPE BC3.
Definition steadyNS.H:191
Eigen::Tensor< double, 3 > C_tensor
Diffusion term.
Definition steadyNS.H:167
autoPtr< volVectorField > _U0
Initial Velocity field (for restart purposes)
Definition steadyNS.H:269
void forcesMatrices(label NUmodes, label NPmodes, label NSUPmodes)
Compute lift and drag matrices.
Definition steadyNS.C:2027
autoPtr< simpleControl > _simple
simpleControl
Definition steadyNS.H:284
surfaceScalarModes L_PHImodes
List of pointers containing the total number of flux modes.
Definition steadyNS.H:119
PtrList< volScalarField > Pfield
List of pointers used to form the pressure snapshots matrix.
Definition steadyNS.H:86
volVectorModes supmodes
List of pointers used to form the supremizer modes.
Definition steadyNS.H:104
autoPtr< surfaceScalarField > _phi0
Initial Flux (for restart purposes)
Definition steadyNS.H:296
List< Eigen::MatrixXd > bcVelocityVec(label NUmodes, label NSUPmodes)
Boundary integral modes on boundary used by the penaly method.
Definition steadyNS.C:1418
List< Eigen::MatrixXd > boundary_vector_convection_consistent(label NUmodes, label NSUPmodes)
Boundary vector convection term - Consistent Flux Method.
Definition steadyNS.C:1956
Eigen::MatrixXd W_matrix
Mass Matrix New Time Step - Consistent Flux Method.
Definition steadyNS.H:197
Eigen::MatrixXd mass_matrix_oldtime_consistent(label NUmodes, label NPmodes, label NSUPmodes)
Mass Matrix old time step (consistent flux method)
Definition steadyNS.C:1773
List< Eigen::MatrixXd > pressure_gradient_term_linsys_diff(label NPmodes)
Laplacian of pressure Linear System - Diffusion term.
Definition steadyNS.C:1738
autoPtr< fv::options > _fvOptions
fvOptions
Definition steadyNS.H:287
Eigen::MatrixXd pressure_BC4(label NPmodes, label NUmodes)
Term N° 4 given by the additional boundary condition using a PPE approach for time-dependent BCs.
Definition steadyNS.C:1378
void solvesupremizer(word type="snapshots")
solve the supremizer either with the use of the pressure snaphots or the pressure modes
Definition steadyNS.C:136
autoPtr< Time > _runTime
Time.
Definition steadyNS.H:290
Eigen::MatrixXd I_matrix
Mass Matrix Old Time Step - Consistent Flux Method.
Definition steadyNS.H:200
~steadyNS()
Definition steadyNS.H:80
volVectorModes Umodes
List of pointers used to form the velocity modes.
Definition steadyNS.H:101
Eigen::MatrixXd diffusive_term_sym(label NUmodes, label NPmodes, label NSUPmodes)
Symetric diffusive Term.
Definition steadyNS.C:892
PtrList< volVectorField > Ufield
List of pointers used to form the velocity snapshots matrix.
Definition steadyNS.H:89
List< Eigen::MatrixXd > boundary_vector_diffusion_consistent(label NUmodes, label NSUPmodes)
Boundary vector diffusion term (consistent flux method)
Definition steadyNS.C:1917
List< Eigen::MatrixXd > LinSysDiff
Projection Peqn onto Pressure modes - Diffusion term.
Definition steadyNS.H:248
autoPtr< dimensionedScalar > dt_dummy
Dummy time step including unit.
Definition steadyNS.H:275
Eigen::MatrixXd diffusive_term_flux_method(label NUmodes, label NPmodes, label NSUPmodes)
Diffusive Flux Method.
Definition steadyNS.C:1497
Eigen::MatrixXd diffusive_term_consistent(label NUmodes, label NPmodes, label NSUPmodes)
Diffusion Term (consistent flux method)
Definition steadyNS.C:1807
Eigen::MatrixXd pressure_BC1(label NPmodes, label NUmodes)
Term N° 1 given by the additional boundary condition using a PPE approach.
Definition steadyNS.C:1212
steadyNS()
Null constructor.
Definition steadyNS.C:40
Eigen::MatrixXd divergence_term(label NUmodes, label NPmodes, label NSUPmodes)
Divergence Term (supremizer approach)
Definition steadyNS.C:1073
List< Eigen::MatrixXd > G_matrix
Divergence of momentum PPE.
Definition steadyNS.H:176
scalar tolerance
Tolerance for the residual of the stationary problems, there is the same tolerance for velocity and p...
Definition steadyNS.H:122
Eigen::MatrixXd mass_term(label NUmodes, label NPmodes, label NSUPmodes)
Mass Term.
Definition steadyNS.C:1041
autoPtr< fvMesh > _mesh
Mesh.
Definition steadyNS.H:281
autoPtr< singlePhaseTransportModel > _laminarTransport
Laminar transport (used by turbulence model)
Definition steadyNS.H:302
label NUmodes
Number of velocity modes used for the projection.
Definition steadyNS.H:140
autoPtr< dimensionedScalar > nu_dummy
Dummy viscocity including unit.
Definition steadyNS.H:278
List< Eigen::MatrixXd > boundary_vector_diffusion(label NUmodes, label NPmodes, label NSUPmodes)
Boundary vector diffusion term.
Definition steadyNS.C:1530
List< Eigen::MatrixXd > LinSysDiv
Projection Peqn onto Pressure modes - Divergence term.
Definition steadyNS.H:245
PtrList< surfaceScalarField > Phifield
List of pointers used to form the flux snapshots matrix.
Definition steadyNS.H:95
PtrList< volVectorField > liftfield
List of pointers used to form the list of lifting functions.
Definition steadyNS.H:110
Eigen::Tensor< double, 3 > convective_term_flux_tens(label NUmodes, label NPmodes, label NSUPmodes)
Convective Term.
Definition steadyNS.C:1622
Eigen::Tensor< double, 3 > convective_term_consistent_tens(label NUmodes, label NPmodes, label NSUPmodes)
Convective Term (consistent flux method)
Definition steadyNS.C:1878
scalar pRefValue
Reference pressure value.
Definition steadyNS.H:311
label pRefCell
Reference pressure cell.
Definition steadyNS.H:308
word fluxMethod
Flux Method.
Definition steadyNS.H:320
volVectorModes L_U_SUPmodes
List of pointers containing the total number of lift, supremizer and velocity modes.
Definition steadyNS.H:116
ITHACAparameters * para
Definition steadyNS.H:82
PtrList< volVectorField > Uomfield
List of pointers used to form the homogeneous velocity snapshots.
Definition steadyNS.H:113
Eigen::MatrixXd pressure_BC3(label NPmodes, label NUmodes)
Term N° 3 given by the additional boundary condition using a PPE approach.
Definition steadyNS.C:1338
List< Eigen::MatrixXd > boundary_vector_convection(label NUmodes, label NPmodes, label NSUPmodes)
Boundary vector convection term.
Definition steadyNS.C:1576
label NNutModesOut
Number of nut modes to be calculated.
Definition steadyNS.H:137
Eigen::MatrixXd B_matrix
Diffusion term.
Definition steadyNS.H:157
List< Eigen::MatrixXd > RD_matrix
Boundary term for diffusion term.
Definition steadyNS.H:227
Eigen::MatrixXd D_matrix
Laplacian term PPE.
Definition steadyNS.H:173
autoPtr< IOMRFZoneList > _MRF
MRF variable.
Definition steadyNS.H:305
List< Eigen::MatrixXd > SD_matrix
Boundary term for diffusion term - Consistent Flux Method.
Definition steadyNS.H:233
label NSUPmodes
Number of supremizer modes used for the projection.
Definition steadyNS.H:146
Eigen::Tensor< double, 3 > gTensor
Divergence of momentum PPE.
Definition steadyNS.H:179
label NNutModes
Number of nut modes used for the projection.
Definition steadyNS.H:149
Eigen::MatrixXd KF_matrix
Pressure Gradient Term - Consistent Flux Method.
Definition steadyNS.H:206
Eigen::MatrixXd K_matrix
Gradient of pressure matrix.
Definition steadyNS.H:163
List< Eigen::MatrixXd > C_matrix
Non linear term.
Definition steadyNS.H:166
Eigen::Tensor< double, 3 > convective_term_tens(label NUmodes, label NPmodes, label NSUPmodes)
Export convective term as a tensor.
Definition steadyNS.C:995
autoPtr< incompressible::turbulenceModel > turbulence
Turbulence model.
Definition steadyNS.H:299
List< Eigen::MatrixXd > convective_term(label NUmodes, label NPmodes, label NSUPmodes)
Convective Term.
Definition steadyNS.C:955
Eigen::MatrixXd pressure_gradient_term_consistent(label NUmodes, label NPmodes, label NSUPmodes)
Pressure Gradient Term (consistent flux method)
Definition steadyNS.C:1844
List< Eigen::MatrixXd > pressure_BC2(label NPmodes, label NUmodes)
Term N° 2 given by the additional boundary condition using a PPE approach.
Definition steadyNS.C:1251
label NUmodesOut
Number of velocity modes to be calculated.
Definition steadyNS.H:128
List< Eigen::MatrixXd > LinSysConv
Projection Peqn onto Pressure modes - Convection term.
Definition steadyNS.H:251
PtrList< volVectorField > supfield
List of pointers used to form the supremizer snapshots matrix.
Definition steadyNS.H:92
Eigen::Tensor< double, 3 > divMomentum(label NUmodes, label NPmodes)
Divergence of convective term (PPE approach)
Definition steadyNS.C:1147
List< Eigen::MatrixXd > SC_matrix
Boundary term for convection term - Consistent Flux Method.
Definition steadyNS.H:236
Eigen::MatrixXd pressure_gradient_term(label NUmodes, label NPmodes, label NSUPmodes)
Gradient of pressure.
Definition steadyNS.C:923
Eigen::MatrixXd DF_matrix
Diffusion Term - Consistent Flux Method.
Definition steadyNS.H:203
List< Eigen::MatrixXd > bcVelMat
Boundary term for penalty method - matrix.
Definition steadyNS.H:221
Eigen::MatrixXd P_matrix
Div of velocity.
Definition steadyNS.H:170
List< Eigen::MatrixXd > RC_matrix
Boundary vector for convection term.
Definition steadyNS.H:230
Eigen::MatrixXd BP_matrix
Diffusion term for flux method PPE.
Definition steadyNS.H:224
Eigen::Tensor< double, 3 > Cf_tensor
Convection term for flux method.
Definition steadyNS.H:239
Eigen::MatrixXd M_matrix
Mass Matrix.
Definition steadyNS.H:160
List< Eigen::MatrixXd > BC2_matrix
PPE BC2.
Definition steadyNS.H:185
autoPtr< volVectorField > Uinl
Initial dummy field with all Dirichlet boundary conditions.
Definition steadyNS.H:272
void discretizeThenProject(fileName folder, label NUmodes, label NPmodes, label NSUPmodes=0)
Project using the Discretize-then-project approach.
Definition steadyNS.C:687
void projectSUP(fileName folder, label NUmodes, label NPmodes, label NSUPmodes)
Project using a supremizer approach.
Definition steadyNS.C:543
Eigen::MatrixXd laplacian_pressure(label NPmodes)
Laplacian of pressure term (PPE approach)
Definition steadyNS.C:1183
label NPmodesOut
Number of pressure modes to be calculated.
Definition steadyNS.H:131
surfaceScalarModes Phimodes
List of pointers used to form the flux modes.
Definition steadyNS.H:107
Eigen::MatrixXd mass_matrix_newtime_consistent(label NUmodes, label NPmodes, label NSUPmodes)
Mass Matrix new time step (consistent flux method)
Definition steadyNS.C:1990
autoPtr< volScalarField > _p0
Initial Pressure field (for restart purposes)
Definition steadyNS.H:266
List< Eigen::MatrixXd > bcVelocityMat(label NUmodes, label NSUPmodes)
Boundary integral modes on boundary used by the penaly method.
Definition steadyNS.C:1455
Eigen::Tensor< double, 3 > bc2Tensor
PPE BC2.
Definition steadyNS.H:188
autoPtr< volVectorField > _U
Velocity field.
Definition steadyNS.H:263
Eigen::Tensor< double, 3 > pressureBC2(label NPmodes, label NUmodes)
Term N° 2 given by the additional boundary condition using a PPE approach.
Definition steadyNS.C:1293
Eigen::Tensor< double, 3 > Ci_tensor
Convection term - Consistent Flux Method.
Definition steadyNS.H:242
void liftSolve()
Perform a lift solve.
Definition steadyNS.C:252
List< Eigen::MatrixXd > div_momentum(label NUmodes, label NPmodes)
Divergence of convective term (PPE approach)
Definition steadyNS.C:1106
volScalarModes Pmodes
List of pointers used to form the pressure modes.
Definition steadyNS.H:98
word bcMethod
Boundary Method.
Definition steadyNS.H:317
scalar cumulativeContErr
continuity error
Definition steadyNS.H:314
void reconstructLiftAndDrag(const Eigen::MatrixXd &velCoeffs, const Eigen::MatrixXd &pressureCoeffs, fileName folder)
Method to reconstruct the forces using velocity and pressure coefficients.
Definition steadyNS.C:2203
autoPtr< volScalarField > _p
Pressure field.
Definition steadyNS.H:260
label NSUPmodesOut
Number of supremizer modes to be calculated.
Definition steadyNS.H:134
Header file of the reductionProblem class.