Implementation of a parametrized full order steady NS problem and preparation of the the reduced matrices for the online solve. More...
#include <SteadyNSSimple.H>
Public Member Functions | |
SteadyNSSimple () | |
Null constructor. | |
SteadyNSSimple (int argc, char *argv[]) | |
Construct with argc and argv. | |
~SteadyNSSimple () | |
void | getTurbRBF (label NNutModes) |
Function to calculate RBF weights for turbulence. | |
void | truthSolve2 (List< scalar > mu_now, word Folder="./ITHACAoutput/Offline/") |
Offline solver for the whole Navier-Stokes problem. | |
Public Member Functions inherited from steadyNS | |
steadyNS () | |
Null constructor. | |
steadyNS (int argc, char *argv[]) | |
Construct with argc and argv. | |
~steadyNS () | |
void | truthSolve (List< scalar > mu_now) |
Perform a truthsolve. | |
void | solvesupremizer (word type="snapshots") |
solve the supremizer either with the use of the pressure snaphots or the pressure modes | |
void | liftSolve () |
Perform a lift solve. | |
void | projectPPE (fileName folder, label NUmodes, label NPmodes, label NSUPmodes=0) |
Project using the Poisson Equation for pressure. | |
void | projectSUP (fileName folder, label NUmodes, label NPmodes, label NSUPmodes) |
Project using a supremizer approach. | |
void | discretizeThenProject (fileName folder, label NUmodes, label NPmodes, label NSUPmodes=0) |
Project using the Discretize-then-project approach. | |
Eigen::MatrixXd | diffusive_term (label NUmodes, label NPmodes, label NSUPmodes) |
Diffusive Term. | |
Eigen::MatrixXd | diffusive_term_sym (label NUmodes, label NPmodes, label NSUPmodes) |
Symetric diffusive Term. | |
Eigen::MatrixXd | pressure_gradient_term (label NUmodes, label NPmodes, label NSUPmodes) |
Gradient of pressure. | |
List< Eigen::MatrixXd > | convective_term (label NUmodes, label NPmodes, label NSUPmodes) |
Convective Term. | |
Eigen::MatrixXd | mass_term (label NUmodes, label NPmodes, label NSUPmodes) |
Mass Term. | |
Eigen::MatrixXd | divergence_term (label NUmodes, label NPmodes, label NSUPmodes) |
Divergence Term (supremizer approach) | |
List< Eigen::MatrixXd > | div_momentum (label NUmodes, label NPmodes) |
Divergence of convective term (PPE approach) | |
Eigen::Tensor< double, 3 > | divMomentum (label NUmodes, label NPmodes) |
Divergence of convective term (PPE approach) | |
Eigen::MatrixXd | laplacian_pressure (label NPmodes) |
Laplacian of pressure term (PPE approach) | |
Eigen::MatrixXd | pressure_BC1 (label NPmodes, label NUmodes) |
Term N° 1 given by the additional boundary condition using a PPE approach. | |
List< Eigen::MatrixXd > | pressure_BC2 (label NPmodes, label NUmodes) |
Term N° 2 given by the additional boundary condition using a PPE approach. | |
Eigen::Tensor< double, 3 > | pressureBC2 (label NPmodes, label NUmodes) |
Term N° 2 given by the additional boundary condition using a PPE approach. | |
Eigen::MatrixXd | pressure_BC3 (label NPmodes, label NUmodes) |
Term N° 3 given by the additional boundary condition using a PPE approach. | |
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. | |
List< Eigen::MatrixXd > | bcVelocityVec (label NUmodes, label NSUPmodes) |
Boundary integral modes on boundary used by the penaly method. | |
List< Eigen::MatrixXd > | bcVelocityMat (label NUmodes, label NSUPmodes) |
Boundary integral modes on boundary used by the penaly method. | |
Eigen::MatrixXd | diffusive_term_flux_method (label NUmodes, label NPmodes, label NSUPmodes) |
Diffusive Flux Method. | |
List< Eigen::MatrixXd > | boundary_vector_diffusion (label NUmodes, label NPmodes, label NSUPmodes) |
Boundary vector diffusion term. | |
List< Eigen::MatrixXd > | boundary_vector_convection (label NUmodes, label NPmodes, label NSUPmodes) |
Boundary vector convection term. | |
Eigen::Tensor< double, 3 > | convective_term_flux_tens (label NUmodes, label NPmodes, label NSUPmodes) |
Convective Term. | |
List< Eigen::MatrixXd > | pressure_gradient_term_linsys_div (label NPmodes) |
Laplacian of pressure Linear System - Divergence term. | |
List< Eigen::MatrixXd > | pressure_gradient_term_linsys_diff (label NPmodes) |
Laplacian of pressure Linear System - Diffusion term. | |
List< Eigen::MatrixXd > | pressure_gradient_term_linsys_conv (label NPmodes) |
Laplacian of pressure Linear System - Convection term. | |
Eigen::MatrixXd | diffusive_term_consistent (label NUmodes, label NPmodes, label NSUPmodes) |
Diffusion Term (consistent flux method) | |
List< Eigen::MatrixXd > | boundary_vector_diffusion_consistent (label NUmodes, label NSUPmodes) |
Boundary vector diffusion term (consistent flux method) | |
List< Eigen::MatrixXd > | boundary_vector_convection_consistent (label NUmodes, label NSUPmodes) |
Boundary vector convection term - Consistent Flux Method. | |
Eigen::MatrixXd | mass_matrix_newtime_consistent (label NUmodes, label NPmodes, label NSUPmodes) |
Mass Matrix new time step (consistent flux method) | |
Eigen::MatrixXd | mass_matrix_oldtime_consistent (label NUmodes, label NPmodes, label NSUPmodes) |
Mass Matrix old time step (consistent flux method) | |
Eigen::MatrixXd | pressure_gradient_term_consistent (label NUmodes, label NPmodes, label NSUPmodes) |
Pressure Gradient Term (consistent flux method) | |
Eigen::Tensor< double, 3 > | convective_term_consistent_tens (label NUmodes, label NPmodes, label NSUPmodes) |
Convective Term (consistent flux method) | |
void | change_viscosity (double mu) |
Function to change the viscosity. | |
void | forcesMatrices (label NUmodes, label NPmodes, label NSUPmodes) |
Compute lift and drag matrices. | |
void | forcesMatrices (label nModes) |
Compute lift and drag matrices offline matrices for the case of same number of velocity and pressure modes. | |
void | reconstructLiftAndDrag (const Eigen::MatrixXd &velCoeffs, const Eigen::MatrixXd &pressureCoeffs, fileName folder) |
Method to reconstruct the forces using velocity and pressure coefficients. | |
Eigen::Tensor< double, 3 > | convective_term_tens (label NUmodes, label NPmodes, label NSUPmodes) |
Export convective term as a tensor. | |
void | restart () |
set U and P back to the values into the 0 folder | |
Public Member Functions inherited from reductionProblem | |
reductionProblem () | |
Construct Null. | |
~reductionProblem () | |
void | setParameters () |
Set Parameters Problems. | |
void | genRandPar () |
Generate Random Numbers. | |
void | genRandPar (label tsize) |
Generate Random Numbers given the dimension of the training set. | |
void | genEquiPar () |
Generate Equidistributed Numbers. | |
void | truthSolve () |
Perform a TruthSolve. | |
void | assignBC (volVectorField &s, label BC_ind, Vector< double > &value) |
Assign Boundary Condition to a volVectorField. | |
void | assignBC (volScalarField &s, label BC_ind, double &value) |
Assign Boundary Condition to a volScalarField. | |
void | reconstructFromMatrix (PtrList< volVectorField > &rec_field2, PtrList< volVectorField > &modes, label Nmodes, Eigen::MatrixXd coeff_matrix) |
Exact reconstruction using a certain number of modes for vector list of fields and the projection coefficients (volVectorField) | |
void | reconstructFromMatrix (PtrList< volScalarField > &rec_field2, PtrList< volScalarField > &modes, label Nmodes, Eigen::MatrixXd coeff_matrix) |
Exact reconstruction using a certain number of modes for vector list of fields and the projection coefficients (volScalarField) | |
template<typename T , typename G > | |
void | assignIF (T &s, G &value) |
Assign internal field condition. | |
template<typename T > | |
void | computeLift (T &Lfield, T &liftfield, T &omfield) |
Homogenize the snapshot matrix, it works with PtrList of volVectorField and volScalarField. | |
template<typename T > | |
void | computeLiftT (T &Lfield, T &liftfield, T &omfield) |
Virtual function to compute the lifting function. | |
void | liftSolve () |
Virtual function to compute the lifting function for scalar field. | |
void | liftSolveT () |
void | project () |
General projection operation. | |
void | writeMu (List< scalar > mu_now) |
Write out a list of scalar corresponding to the parameters used in the truthSolve. | |
std::vector< SPLINTER::RBFSpline > | getCoeffManifoldRBF (PtrList< volVectorField > snapshots, PtrList< volVectorField > &modes, word rbfBasis="GAUSSIAN") |
Constructs the parameters-coefficients manifold for vector fields, based on RBF-spline model. | |
std::vector< SPLINTER::RBFSpline > | getCoeffManifoldRBF (PtrList< volScalarField > snapshots, PtrList< volScalarField > &modes, word rbfBasis="GAUSSIAN") |
Constructs the parameters-coefficients manifold for scalar fields, based on RBF-spline model. | |
std::vector< SPLINTER::BSpline > | getCoeffManifoldSPL (PtrList< volVectorField > snapshots, PtrList< volVectorField > &modes, label splDeg=3) |
Constructs the parameters-coefficients manifold for vector fields, based on the B-spline model. | |
std::vector< SPLINTER::BSpline > | getCoeffManifoldSPL (PtrList< volScalarField > snapshots, PtrList< volScalarField > &modes, label splDeg=3) |
Constructs the parameters-coefficients manifold for scalar fields, based on the B-spline model. | |
Public Attributes | |
PtrList< volScalarField > | nutFields |
List of snapshots for the solution for eddy viscosity. | |
volScalarModes | nutModes |
List of POD modes for eddy viscosity. | |
std::vector< SPLINTER::DataTable * > | samples |
Create a samples for interpolation. | |
std::vector< SPLINTER::RBFSpline * > | rbfSplines |
Create a RBF splines for interpolation. | |
Eigen::MatrixXd | coeffL2 |
The matrix of L2 projection coefficients for the eddy viscosity. | |
label | middleStep |
Distancing between intermediate steps (for turbulent case only) | |
bool | middleExport |
Export also intermediate fields. | |
label | saver |
Counter to check if the middleStep has been reached or not (for turbulent case only) | |
label | folderN |
Counter to save intermediate steps in the correct folder (for turbulent case only) | |
fvVectorMatrix * | Ueqn_global |
Initialization for the full velocity linear system. | |
fvScalarMatrix * | Peqn_global |
Initialization for the full pressure linear system. | |
Public Attributes inherited from steadyNS | |
ITHACAparameters * | para |
PtrList< volScalarField > | Pfield |
List of pointers used to form the pressure snapshots matrix. | |
PtrList< volVectorField > | Ufield |
List of pointers used to form the velocity snapshots matrix. | |
PtrList< volVectorField > | supfield |
List of pointers used to form the supremizer snapshots matrix. | |
PtrList< surfaceScalarField > | Phifield |
List of pointers used to form the flux snapshots matrix. | |
volScalarModes | Pmodes |
List of pointers used to form the pressure modes. | |
volVectorModes | Umodes |
List of pointers used to form the velocity modes. | |
volVectorModes | supmodes |
List of pointers used to form the supremizer modes. | |
surfaceScalarModes | Phimodes |
List of pointers used to form the flux modes. | |
PtrList< volVectorField > | liftfield |
List of pointers used to form the list of lifting functions. | |
PtrList< volVectorField > | Uomfield |
List of pointers used to form the homogeneous velocity snapshots. | |
volVectorModes | L_U_SUPmodes |
List of pointers containing the total number of lift, supremizer and velocity modes. | |
surfaceScalarModes | L_PHImodes |
List of pointers containing the total number of flux modes. | |
scalar | tolerance |
Tolerance for the residual of the stationary problems, there is the same tolerance for velocity and pressure. | |
scalar | maxIter |
Number of maximum iterations to be done for the computation of the truth solution. | |
label | NUmodesOut |
Number of velocity modes to be calculated. | |
label | NPmodesOut |
Number of pressure modes to be calculated. | |
label | NSUPmodesOut |
Number of supremizer modes to be calculated. | |
label | NNutModesOut |
Number of nut modes to be calculated. | |
label | NUmodes |
Number of velocity modes used for the projection. | |
label | NPmodes |
Number of pressure modes used for the projection. | |
label | NSUPmodes |
Number of supremizer modes used for the projection. | |
label | NNutModes |
Number of nut modes used for the projection. | |
Eigen::MatrixXd | tauMatrix |
Viscous forces. | |
Eigen::MatrixXd | nMatrix |
Pressure forces. | |
List< Eigen::MatrixXd > | bcVelVec |
Boundary term for penalty method - vector. | |
List< Eigen::MatrixXd > | bcVelMat |
Boundary term for penalty method - matrix. | |
Eigen::MatrixXd | BP_matrix |
Diffusion term for flux method PPE. | |
List< Eigen::MatrixXd > | RD_matrix |
Boundary term for diffusion term. | |
List< Eigen::MatrixXd > | RC_matrix |
Boundary vector for convection term. | |
List< Eigen::MatrixXd > | SD_matrix |
Boundary term for diffusion term - Consistent Flux Method. | |
List< Eigen::MatrixXd > | SC_matrix |
Boundary term for convection term - Consistent Flux Method. | |
Eigen::Tensor< double, 3 > | Cf_tensor |
Convection term for flux method. | |
Eigen::Tensor< double, 3 > | Ci_tensor |
Convection term - Consistent Flux Method. | |
List< Eigen::MatrixXd > | LinSysDiv |
Projection Peqn onto Pressure modes - Divergence term. | |
List< Eigen::MatrixXd > | LinSysDiff |
Projection Peqn onto Pressure modes - Diffusion term. | |
List< Eigen::MatrixXd > | LinSysConv |
Projection Peqn onto Pressure modes - Convection term. | |
bool | supex |
Boolean variable to check the existence of the supremizer modes. | |
autoPtr< volScalarField > | _p |
Pressure field. | |
autoPtr< volVectorField > | _U |
Velocity field. | |
autoPtr< volScalarField > | _p0 |
Initial Pressure field (for restart purposes) | |
autoPtr< volVectorField > | _U0 |
Initial Velocity field (for restart purposes) | |
autoPtr< volVectorField > | Uinl |
Initial dummy field with all Dirichlet boundary conditions. | |
autoPtr< dimensionedScalar > | dt_dummy |
Dummy time step including unit. | |
autoPtr< dimensionedScalar > | nu_dummy |
Dummy viscocity including unit. | |
autoPtr< fvMesh > | _mesh |
Mesh. | |
autoPtr< simpleControl > | _simple |
simpleControl | |
autoPtr< fv::options > | _fvOptions |
fvOptions | |
autoPtr< Time > | _runTime |
Time. | |
autoPtr< surfaceScalarField > | _phi |
Flux. | |
autoPtr< surfaceScalarField > | _phi0 |
Initial Flux (for restart purposes) | |
autoPtr< incompressible::turbulenceModel > | turbulence |
Turbulence model. | |
autoPtr< singlePhaseTransportModel > | _laminarTransport |
Laminar transport (used by turbulence model) | |
autoPtr< IOMRFZoneList > | _MRF |
MRF variable. | |
label | pRefCell |
Reference pressure cell. | |
scalar | pRefValue |
Reference pressure value. | |
scalar | cumulativeContErr = 0 |
continuity error | |
word | bcMethod |
Boundary Method. | |
word | fluxMethod |
Flux Method. | |
Eigen::MatrixXd | B_matrix |
Diffusion term. | |
Eigen::MatrixXd | M_matrix |
Mass Matrix. | |
Eigen::MatrixXd | K_matrix |
Gradient of pressure matrix. | |
List< Eigen::MatrixXd > | C_matrix |
Non linear term. | |
Eigen::Tensor< double, 3 > | C_tensor |
Diffusion term. | |
Eigen::MatrixXd | P_matrix |
Div of velocity. | |
Eigen::MatrixXd | D_matrix |
Laplacian term PPE. | |
List< Eigen::MatrixXd > | G_matrix |
Divergence of momentum PPE. | |
Eigen::Tensor< double, 3 > | gTensor |
Divergence of momentum PPE. | |
Eigen::MatrixXd | BC1_matrix |
PPE BC1. | |
List< Eigen::MatrixXd > | BC2_matrix |
PPE BC2. | |
Eigen::Tensor< double, 3 > | bc2Tensor |
PPE BC2. | |
Eigen::MatrixXd | BC3_matrix |
PPE BC3. | |
Eigen::MatrixXd | BC4_matrix |
PPE BC4. | |
Eigen::MatrixXd | W_matrix |
Mass Matrix New Time Step - Consistent Flux Method. | |
Eigen::MatrixXd | I_matrix |
Mass Matrix Old Time Step - Consistent Flux Method. | |
Eigen::MatrixXd | DF_matrix |
Diffusion Term - Consistent Flux Method. | |
Eigen::MatrixXd | KF_matrix |
Pressure Gradient Term - Consistent Flux Method. | |
Public Attributes inherited from reductionProblem | |
label | Pnumber |
Number of parameters. | |
label | Tnumber |
Dimension of the training set (used only when gerating parameters without input) | |
Eigen::MatrixXd | mu |
Row matrix of parameters. | |
Eigen::MatrixXd | mu_range |
Range of the parameter spaces. | |
Eigen::MatrixXd | mu_samples |
Matrix of parameters to be used for PODI, where each row corresponds to a sample point. In this matrix the time dimension is regarded as a parameter for unsteady problems. | |
double | mu_cur |
Current value of the parameter. | |
bool | podex |
Boolean variable, it is 1 if the POD has already been computed, else 0. | |
bool | offline |
Boolean variable, it is 1 if the Offline phase has already been computed, else 0. | |
IOdictionary * | ITHACAdict |
dictionary to store input output infos | |
autoPtr< argList > | _args |
argList | |
ITHACAparallel * | paral |
parallel handling | |
label | folderN = 1 |
Counter to save intermediate steps in the correct folder, for unsteady and some stationary cases. | |
label | counter = 1 |
Counter used for the output of the full order solutions. | |
Eigen::MatrixXi | inletIndex |
Matrix that contains informations about the inlet boundaries. | |
Eigen::MatrixXi | inletPatch |
Matrix that contains informations about the inlet boundaries without specifing the direction Rows = Number of parametrized boundary conditions Cols = 1 Example: example.inletIndex.resize(2, 1); example.inletIndex(0, 0) = 0; example.inletIndex(1, 0) = 1; Means that there are two parametrized boundary conditions of which the first row is of patch 0 and the second row of patch 1. | |
Eigen::MatrixXi | inletIndexT |
Implementation of a parametrized full order steady NS problem and preparation of the the reduced matrices for the online solve.
In this class are implemented the methods for the offline solve of a steady NS problem and the for the generation of the reduced matrices for subsequent online solve, this class is a son of the reduction problem class
Definition at line 69 of file SteadyNSSimple.H.
SteadyNSSimple::SteadyNSSimple | ( | ) |
Null constructor.
Definition at line 41 of file SteadyNSSimple.C.
SteadyNSSimple::SteadyNSSimple | ( | int | argc, |
char * | argv[] ) |
Construct with argc and argv.
Number of velocity modes to be calculated
Number of pressure modes to be calculated
Number of nut modes to be calculated
Number of velocity modes used for the projection
Number of pressure modes used for the projection
Number of nut modes used for the projection
Definition at line 43 of file SteadyNSSimple.C.
|
inline |
Definition at line 79 of file SteadyNSSimple.H.
void SteadyNSSimple::getTurbRBF | ( | label | NNutModes | ) |
Function to calculate RBF weights for turbulence.
[in] | NNutModes | Number of modes to be used for nut. |
Definition at line 62 of file SteadyNSSimple.C.
void SteadyNSSimple::truthSolve2 | ( | List< scalar > | mu_now, |
word | Folder = "./ITHACAoutput/Offline/" ) |
Offline solver for the whole Navier-Stokes problem.
[in] | mu_now | Viscosity (parametrized) we want to solve the problem with. |
Definition at line 114 of file SteadyNSSimple.C.
Eigen::MatrixXd SteadyNSSimple::coeffL2 |
The matrix of L2 projection coefficients for the eddy viscosity.
Definition at line 94 of file SteadyNSSimple.H.
label SteadyNSSimple::folderN |
Counter to save intermediate steps in the correct folder (for turbulent case only)
Definition at line 106 of file SteadyNSSimple.H.
bool SteadyNSSimple::middleExport |
Export also intermediate fields.
Definition at line 100 of file SteadyNSSimple.H.
label SteadyNSSimple::middleStep |
Distancing between intermediate steps (for turbulent case only)
Definition at line 97 of file SteadyNSSimple.H.
PtrList<volScalarField> SteadyNSSimple::nutFields |
List of snapshots for the solution for eddy viscosity.
Definition at line 82 of file SteadyNSSimple.H.
volScalarModes SteadyNSSimple::nutModes |
List of POD modes for eddy viscosity.
Definition at line 85 of file SteadyNSSimple.H.
fvScalarMatrix* SteadyNSSimple::Peqn_global |
Initialization for the full pressure linear system.
Definition at line 112 of file SteadyNSSimple.H.
std::vector<SPLINTER::RBFSpline*> SteadyNSSimple::rbfSplines |
Create a RBF splines for interpolation.
Definition at line 91 of file SteadyNSSimple.H.
std::vector<SPLINTER::DataTable*> SteadyNSSimple::samples |
Create a samples for interpolation.
Definition at line 88 of file SteadyNSSimple.H.
label SteadyNSSimple::saver |
Counter to check if the middleStep has been reached or not (for turbulent case only)
Definition at line 103 of file SteadyNSSimple.H.
fvVectorMatrix* SteadyNSSimple::Ueqn_global |
Initialization for the full velocity linear system.
Definition at line 110 of file SteadyNSSimple.H.