Loading...
Searching...
No Matches
UnsteadyBB Class Reference

Implementation of a parametrized full order unsteady Boussinesq problem and preparation of the the reduced matrices for the online solve. More...

#include <UnsteadyBB.H>

Inheritance diagram for UnsteadyBB:
unsteadyNS steadyNS UnsteadyProblem reductionProblem tutorial10 tutorial11

Public Member Functions

 UnsteadyBB ()
 Null constructor.
 
 UnsteadyBB (int argc, char *argv[])
 Construct with argc and argv.
 
 ~UnsteadyBB ()
 
void truthSolve (List< scalar > mu_now)
 Perform a truthsolve for parameters mu_now.
 
void truthSolve (fileName folder="./ITHACAOutput/Offline")
 Perform a truthsolve for full order solution.
 
void solvesupremizer (word type="snapshots")
 solve the supremizer either with the use of the pressure snaphots or the pressure modes
 
void liftSolveT ()
 Perform a lift solve for temperature.
 
void liftSolve ()
 Perform a lift solve for velocity field.
 
void projectSUP (fileName folder, label NUmodes, label NPmodes, label NTmodes, label NSUPmodes)
 Project using a supremizer approach.
 
void projectPPE (fileName folder, label NUmodes, label NPrghmodes, label NTmodes, label NSUPmodes)
 Project using a PPE approach.
 
Eigen::MatrixXd pressure_gradient_term (label NUmodes, label NPrghmodes, label NSUPmodes)
 Gradient of pressure.
 
Eigen::MatrixXd diffusive_term_temperature (label NUmodes, label NTmodes, label NSUPmodes)
 Diffusive Term Energy Equation.
 
Eigen::MatrixXd divergence_term (label NUmodes, label NPrghmodes, label NSUPmodes)
 Divergence Term (supremizer approach)
 
List< Eigen::MatrixXd > convective_term_temperature (label NUmodes, label NTmodes, label NSUPmodes)
 Convective Term Energy Equation.
 
Eigen::MatrixXd mass_term_temperature (label NUmodes, label NTmodes, label NSUPmodes)
 Mass Term Energy Equation.
 
Eigen::MatrixXd buoyant_term (label NUmodes, label NTmodes, label NSUPmodes)
 Buoyant Term Momentum Equation.
 
Eigen::MatrixXd buoyant_term_poisson (label NPrghmodes, label NTmodes)
 Buoyant Term PPE Equation.
 
bool checkWrite (Time &timeObject)
 Function to check if the solution must be exported.
 
void change_viscosity (double mu)
 Function to change the viscosity.
 
Eigen::Tensor< double, 3 > convective_term_tens_temperature (label NUmodes, label NTmodes, label NSUPmodes)
 Export convective term energy equation as a tensor.
 
- Public Member Functions inherited from unsteadyNS
 unsteadyNS ()
 Construct Null.
 
 unsteadyNS (int argc, char *argv[])
 Construct with argc and argv.
 
void truthSolve (List< scalar > mu_now, fileName folder="./ITHACAoutput/Offline/")
 Perform a truthsolve.
 
- 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 Member Functions inherited from UnsteadyProblem
 UnsteadyProblem ()
 
void setTimes (Time &timeObject)
 
bool checkWrite (Time &timeObject)
 Function to check if the solution must be exported.
 

Public Attributes

PtrList< volScalarField > Prghfield
 List of pointers used to form the shifted pressure snapshots matrix.
 
PtrList< volScalarField > Tmodes
 List of pointers used to form the temperature modes.
 
PtrList< volScalarField > Prghmodes
 List of pointers used to form the shifted pressure modes.
 
PtrList< volScalarField > Tfield
 List of pointers used to form the temperature snapshots matrix.
 
PtrList< volScalarField > Tfield_on
 List of pointers used to form the temperature snapshots matrix.
 
PtrList< volScalarField > Pfield_on
 List of pointers used to form the temperature snapshots matrix.
 
PtrList< volVectorField > Ufield_on
 List of pointers used to form the temperature snapshots matrix.
 
PtrList< volScalarField > liftfieldT
 List of pointers used to form the list of lifting functions.
 
PtrList< volScalarField > Tomfield
 List of pointers used to form the homogeneous velocity snapshots.
 
PtrList< volScalarField > L_T_modes
 List of pointers containing the lift for temperature and the temperature field.
 
autoPtr< volScalarField > _UliftBC
 
label NTmodes
 Number of temperature modes used for the projection.
 
label NPrghmodes
 Number of pressure modes used for the projection.
 
Eigen::MatrixXd H_matrix
 Buoyancy term - momentum equation.
 
Eigen::MatrixXd HP_matrix
 Buoyancy term - PPE equation.
 
Eigen::MatrixXd Y_matrix
 Diffusive term - energy equation.
 
List< Eigen::MatrixXd > Q_matrix
 Non linear convective term - energy equation.
 
Eigen::Tensor< double, 3 > Q_tensor
 
Eigen::MatrixXd W_matrix
 Mass Matrix - energy equation.
 
autoPtr< fvMesh > _mesh
 Mesh.
 
autoPtr< volScalarField > _p_rgh
 Shifted Pressure field.
 
autoPtr< volScalarField > _T
 Temperature field.
 
autoPtr< dimensionedScalar > _beta
 dimensionedScalar beta;
 
autoPtr< dimensionedScalar > _TRef
 dimensionedScalar Tref;
 
autoPtr< dimensionedScalar > _Pr
 dimensionedScalar Pr;
 
autoPtr< dimensionedScalar > _Prt
 dimensionedScalar Prt;
 
autoPtr< volScalarField > _alphat
 dimensionedScalar alphat;
 
autoPtr< volScalarField > _nut
 dimensionedScalar nut;
 
autoPtr< volScalarField > _rhok
 dimensionedScalar rhok;
 
autoPtr< dimensionedScalar > _nu
 dimensionedScalar nu;
 
autoPtr< volScalarField > _gh
 List of pointers used to form the gravitational acceleration.
 
autoPtr< surfaceScalarField > _ghf
 List of pointers used to form the gravitational acceleration.
 
autoPtr< dimensionedVector > _g
 
autoPtr< dimensionedScalar > _hRef
 dimensionedScalar hRef;
 
autoPtr< dimensionedScalar > _ghRef
 dimensionedScalar ghRef;
 
autoPtr< volScalarField > _S
 Source Term Heat.
 
PtrList< volVectorField > UModesWeighted
 Pointers for perfoming Nested-POD Method List of pointers used to form the Weighted nested-POD modes for velocity field.
 
PtrList< volScalarField > PModesWeighted
 List of pointers used to form the Weighted nested-POD modes for pressure field.
 
PtrList< volScalarField > TModesWeighted
 List of pointers used to form the Weighted nested-POD modes for temperature field.
 
PtrList< volScalarField > NUTModesWeighted
 List of pointers used to form the Weighted nested-POD modes for nut-field.
 
- Public Attributes inherited from unsteadyNS
autoPtr< pimpleControl > _pimple
 pimpleControl
 
autoPtr< incompressible::turbulenceModel > turbulence
 Turbulence model.
 
bool adjustTimeStep
 adjustTimeStep
 
scalar maxCo
 maxCourant
 
scalar maxDeltaT
 maxDeltaT
 
label counter2 = 1
 
word method
 
word timedepbcMethod
 Time-dependent Boundary Method.
 
Eigen::MatrixXd timeBCoff
 
word timeDerivativeSchemeOrder
 
- Public Attributes inherited from steadyNS
ITHACAparameterspara
 
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
 
ITHACAparallelparal
 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
 
- Public Attributes inherited from UnsteadyProblem
scalar startTime
 Start Time (initial time to start storing the snapshots)
 
scalar finalTime
 Final time (final time of the simulation and consequently of the acquisition of the snapshots)
 
scalar timeStep
 Time step of the simulation.
 
scalar writeEvery = timeStep
 Time step of the writing procedure.
 
scalar nextWrite
 Auxiliary variable to store the next writing instant.
 

Detailed Description

Implementation of a parametrized full order unsteady Boussinesq 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 unsteady Boussinesq problem (two way coupling between momentum and energy equation) and the for the generation of the reduced matrices for subsequent online solve, this class is a son of the reduction problem class

Examples
10UnsteadyBBEnclosed.C, and 11UnsteadyBBOpen.C.

Definition at line 59 of file UnsteadyBB.H.

Constructor & Destructor Documentation

◆ UnsteadyBB() [1/2]

UnsteadyBB::UnsteadyBB ( )

Null constructor.

Examples
10UnsteadyBBEnclosed.C, and 11UnsteadyBBOpen.C.

Definition at line 41 of file UnsteadyBB.C.

◆ UnsteadyBB() [2/2]

UnsteadyBB::UnsteadyBB ( int argc,
char * argv[] )

Construct with argc and argv.

Definition at line 42 of file UnsteadyBB.C.

◆ ~UnsteadyBB()

UnsteadyBB::~UnsteadyBB ( )
inline

Definition at line 68 of file UnsteadyBB.H.

Member Function Documentation

◆ buoyant_term()

Eigen::MatrixXd UnsteadyBB::buoyant_term ( label NUmodes,
label NTmodes,
label NSUPmodes )

Buoyant Term Momentum Equation.

Parameters
[in]NUmodesThe number of velocity modes.
[in]NTmodesThe number of temperature modes.
[in]NSUPmodesThe number of supremizer modes.
Returns
reduced third order tensor in List <Eigen::MatrixXd> format for the convective term.

Definition at line 906 of file UnsteadyBB.C.

◆ buoyant_term_poisson()

Eigen::MatrixXd UnsteadyBB::buoyant_term_poisson ( label NPrghmodes,
label NTmodes )

Buoyant Term PPE Equation.

Parameters
[in]NPrghmodesThe number of Prgh Pressure modes.
[in]NTmodesThe number of temperature modes.
Returns
reduced matrix in List <Eigen::MatrixXd> format for the buoyancy-related PPE term.

Definition at line 947 of file UnsteadyBB.C.

◆ change_viscosity()

void UnsteadyBB::change_viscosity ( double mu)

Function to change the viscosity.

Parameters
[in]muviscosity (scalar)

Definition at line 1200 of file UnsteadyBB.C.

◆ checkWrite()

bool UnsteadyBB::checkWrite ( Time & timeObject)

Function to check if the solution must be exported.

Parameters
timeObjectThe time object of OpenFOAM.
Returns
1 if we must write 0 elsewhere.

Definition at line 264 of file UnsteadyBB.C.

◆ convective_term_temperature()

List< Eigen::MatrixXd > UnsteadyBB::convective_term_temperature ( label NUmodes,
label NTmodes,
label NSUPmodes )

Convective Term Energy Equation.

Parameters
[in]NUmodesThe number of velocity modes.
[in]NTmodesThe number of temperature modes.
[in]NSUPmodesThe number of supremizer modes.
Returns
reduced third order tensor in List <Eigen::MatrixXd> format for the convective term.

Definition at line 991 of file UnsteadyBB.C.

◆ convective_term_tens_temperature()

Eigen::Tensor< double, 3 > UnsteadyBB::convective_term_tens_temperature ( label NUmodes,
label NTmodes,
label NSUPmodes )

Export convective term energy equation as a tensor.

Parameters
[in]NUmodesThe N of velocity modes
[in]NTmodesThe N of temperature modes
[in]NSUPmodesThe N of supremizer modes
Returns
tensor_Q

◆ diffusive_term_temperature()

Eigen::MatrixXd UnsteadyBB::diffusive_term_temperature ( label NUmodes,
label NTmodes,
label NSUPmodes )

Diffusive Term Energy Equation.

Parameters
[in]NUmodesThe number of velocity modes.
[in]NTmodesThe number of temperature modes.
[in]NSUPmodesThe number of supremizer modes.
Returns
reduced matrix in Eigen::MatrixXd format for the diffusion term.

Definition at line 1035 of file UnsteadyBB.C.

◆ divergence_term()

Eigen::MatrixXd UnsteadyBB::divergence_term ( label NUmodes,
label NPrghmodes,
label NSUPmodes )

Divergence Term (supremizer approach)

Parameters
[in]NUmodesThe number of velocity modes.
[in]NPmodesThe number of pressure modes.
[in]NSUPmodesThe number of supremizer modes.
Returns
reduced matrix in Eigen::MatrixXd format for the divergence term.

Definition at line 873 of file UnsteadyBB.C.

◆ liftSolve()

void UnsteadyBB::liftSolve ( )

Perform a lift solve for velocity field.

Definition at line 1097 of file UnsteadyBB.C.

◆ liftSolveT()

void UnsteadyBB::liftSolveT ( )

Perform a lift solve for temperature.

Definition at line 422 of file UnsteadyBB.C.

◆ mass_term_temperature()

Eigen::MatrixXd UnsteadyBB::mass_term_temperature ( label NUmodes,
label NTmodes,
label NSUPmodes )

Mass Term Energy Equation.

Parameters
[in]NUmodesThe number of velocity modes.
[in]NTmodesThe number of temperature modes.
[in]NSUPmodesThe number of supremizer modes.
Returns
reduced matrix in Eigen::MatrixXd format for the mass matrix.

Definition at line 1066 of file UnsteadyBB.C.

◆ pressure_gradient_term()

Eigen::MatrixXd UnsteadyBB::pressure_gradient_term ( label NUmodes,
label NPrghmodes,
label NSUPmodes )

Gradient of pressure.

Parameters
[in]NUmodesThe number of velocity modes.
[in]NPmodesThe number of pressure modes.
[in]NSUPmodesThe number of supremizer modes.
Returns
reduced matrix in Eigen::MatrixXd format for the Gradient of pressure term.

Definition at line 835 of file UnsteadyBB.C.

◆ projectPPE()

void UnsteadyBB::projectPPE ( fileName folder,
label NUmodes,
label NPrghmodes,
label NTmodes,
label NSUPmodes )

Project using a PPE approach.

Parameters
[in]folderThe folder used to save the reduced matrices.
[in]NUmodesThe number of velocity modes.
[in]NPmodesThe number of pressure modes.
[in]NSUPmodesThe number of supremizer modes.

Definition at line 652 of file UnsteadyBB.C.

◆ projectSUP()

void UnsteadyBB::projectSUP ( fileName folder,
label NUmodes,
label NPmodes,
label NTmodes,
label NSUPmodes )

Project using a supremizer approach.

Parameters
[in]folderThe folder used to save the reduced matrices.
[in]NUmodesThe number of velocity modes.
[in]NPmodesThe number of pressure modes.
[in]NSUPmodesThe number of supremizer modes.

Definition at line 488 of file UnsteadyBB.C.

◆ solvesupremizer()

void UnsteadyBB::solvesupremizer ( word type = "snapshots")

solve the supremizer either with the use of the pressure snaphots or the pressure modes

Parameters
[in]typeThe type of the supremizer approach, either done on the pressure snapshots or on pressure modes.

Definition at line 281 of file UnsteadyBB.C.

◆ truthSolve() [1/2]

void UnsteadyBB::truthSolve ( fileName folder = "./ITHACAOutput/Offline")

Perform a truthsolve for full order solution.

Definition at line 173 of file UnsteadyBB.C.

◆ truthSolve() [2/2]

void UnsteadyBB::truthSolve ( List< scalar > mu_now)

Perform a truthsolve for parameters mu_now.

Definition at line 79 of file UnsteadyBB.C.

Member Data Documentation

◆ _alphat

autoPtr<volScalarField> UnsteadyBB::_alphat

dimensionedScalar alphat;

Examples
10UnsteadyBBEnclosed.C, and 11UnsteadyBBOpen.C.

Definition at line 149 of file UnsteadyBB.H.

◆ _beta

autoPtr<dimensionedScalar> UnsteadyBB::_beta

dimensionedScalar beta;

Definition at line 137 of file UnsteadyBB.H.

◆ _g

autoPtr<dimensionedVector> UnsteadyBB::_g

Definition at line 166 of file UnsteadyBB.H.

◆ _gh

autoPtr<volScalarField> UnsteadyBB::_gh

List of pointers used to form the gravitational acceleration.

Definition at line 161 of file UnsteadyBB.H.

◆ _ghf

autoPtr<surfaceScalarField> UnsteadyBB::_ghf

List of pointers used to form the gravitational acceleration.

Definition at line 164 of file UnsteadyBB.H.

◆ _ghRef

autoPtr<dimensionedScalar> UnsteadyBB::_ghRef

dimensionedScalar ghRef;

Definition at line 172 of file UnsteadyBB.H.

◆ _hRef

autoPtr<dimensionedScalar> UnsteadyBB::_hRef

dimensionedScalar hRef;

Definition at line 169 of file UnsteadyBB.H.

◆ _mesh

autoPtr<fvMesh> UnsteadyBB::_mesh
mutable

Mesh.

Examples
10UnsteadyBBEnclosed.C, and 11UnsteadyBBOpen.C.

Definition at line 128 of file UnsteadyBB.H.

◆ _nu

autoPtr<dimensionedScalar> UnsteadyBB::_nu

dimensionedScalar nu;

Definition at line 158 of file UnsteadyBB.H.

◆ _nut

autoPtr<volScalarField> UnsteadyBB::_nut

dimensionedScalar nut;

Definition at line 152 of file UnsteadyBB.H.

◆ _p_rgh

autoPtr<volScalarField> UnsteadyBB::_p_rgh

Shifted Pressure field.

Examples
10UnsteadyBBEnclosed.C, and 11UnsteadyBBOpen.C.

Definition at line 131 of file UnsteadyBB.H.

◆ _Pr

autoPtr<dimensionedScalar> UnsteadyBB::_Pr

dimensionedScalar Pr;

Examples
10UnsteadyBBEnclosed.C, and 11UnsteadyBBOpen.C.

Definition at line 143 of file UnsteadyBB.H.

◆ _Prt

autoPtr<dimensionedScalar> UnsteadyBB::_Prt

dimensionedScalar Prt;

Examples
10UnsteadyBBEnclosed.C, and 11UnsteadyBBOpen.C.

Definition at line 146 of file UnsteadyBB.H.

◆ _rhok

autoPtr<volScalarField> UnsteadyBB::_rhok

dimensionedScalar rhok;

Definition at line 155 of file UnsteadyBB.H.

◆ _S

autoPtr<volScalarField> UnsteadyBB::_S

Source Term Heat.

Definition at line 175 of file UnsteadyBB.H.

◆ _T

autoPtr<volScalarField> UnsteadyBB::_T

Temperature field.

Examples
10UnsteadyBBEnclosed.C, and 11UnsteadyBBOpen.C.

Definition at line 134 of file UnsteadyBB.H.

◆ _TRef

autoPtr<dimensionedScalar> UnsteadyBB::_TRef

dimensionedScalar Tref;

Definition at line 140 of file UnsteadyBB.H.

◆ _UliftBC

autoPtr<volScalarField> UnsteadyBB::_UliftBC

Definition at line 103 of file UnsteadyBB.H.

◆ H_matrix

Eigen::MatrixXd UnsteadyBB::H_matrix

Buoyancy term - momentum equation.

Definition at line 112 of file UnsteadyBB.H.

◆ HP_matrix

Eigen::MatrixXd UnsteadyBB::HP_matrix

Buoyancy term - PPE equation.

Definition at line 115 of file UnsteadyBB.H.

◆ L_T_modes

PtrList<volScalarField> UnsteadyBB::L_T_modes

List of pointers containing the lift for temperature and the temperature field.

Definition at line 100 of file UnsteadyBB.H.

◆ liftfieldT

PtrList<volScalarField> UnsteadyBB::liftfieldT

List of pointers used to form the list of lifting functions.

Examples
10UnsteadyBBEnclosed.C, and 11UnsteadyBBOpen.C.

Definition at line 94 of file UnsteadyBB.H.

◆ NPrghmodes

label UnsteadyBB::NPrghmodes

Number of pressure modes used for the projection.

Definition at line 109 of file UnsteadyBB.H.

◆ NTmodes

label UnsteadyBB::NTmodes

Number of temperature modes used for the projection.

Definition at line 106 of file UnsteadyBB.H.

◆ NUTModesWeighted

PtrList<volScalarField> UnsteadyBB::NUTModesWeighted

List of pointers used to form the Weighted nested-POD modes for nut-field.

Definition at line 188 of file UnsteadyBB.H.

◆ Pfield_on

PtrList<volScalarField> UnsteadyBB::Pfield_on

List of pointers used to form the temperature snapshots matrix.

Definition at line 88 of file UnsteadyBB.H.

◆ PModesWeighted

PtrList<volScalarField> UnsteadyBB::PModesWeighted

List of pointers used to form the Weighted nested-POD modes for pressure field.

Definition at line 182 of file UnsteadyBB.H.

◆ Prghfield

PtrList<volScalarField> UnsteadyBB::Prghfield

List of pointers used to form the shifted pressure snapshots matrix.

Examples
11UnsteadyBBOpen.C.

Definition at line 73 of file UnsteadyBB.H.

◆ Prghmodes

PtrList<volScalarField> UnsteadyBB::Prghmodes

List of pointers used to form the shifted pressure modes.

Definition at line 79 of file UnsteadyBB.H.

◆ Q_matrix

List<Eigen::MatrixXd> UnsteadyBB::Q_matrix

Non linear convective term - energy equation.

Definition at line 121 of file UnsteadyBB.H.

◆ Q_tensor

Eigen::Tensor<double, 3 > UnsteadyBB::Q_tensor

Definition at line 122 of file UnsteadyBB.H.

◆ Tfield

PtrList<volScalarField> UnsteadyBB::Tfield

List of pointers used to form the temperature snapshots matrix.

Examples
10UnsteadyBBEnclosed.C, and 11UnsteadyBBOpen.C.

Definition at line 82 of file UnsteadyBB.H.

◆ Tfield_on

PtrList<volScalarField> UnsteadyBB::Tfield_on

List of pointers used to form the temperature snapshots matrix.

Examples
10UnsteadyBBEnclosed.C, and 11UnsteadyBBOpen.C.

Definition at line 85 of file UnsteadyBB.H.

◆ Tmodes

PtrList<volScalarField> UnsteadyBB::Tmodes

List of pointers used to form the temperature modes.

Definition at line 76 of file UnsteadyBB.H.

◆ TModesWeighted

PtrList<volScalarField> UnsteadyBB::TModesWeighted

List of pointers used to form the Weighted nested-POD modes for temperature field.

Definition at line 185 of file UnsteadyBB.H.

◆ Tomfield

PtrList<volScalarField> UnsteadyBB::Tomfield

List of pointers used to form the homogeneous velocity snapshots.

Definition at line 97 of file UnsteadyBB.H.

◆ Ufield_on

PtrList<volVectorField> UnsteadyBB::Ufield_on

List of pointers used to form the temperature snapshots matrix.

Examples
10UnsteadyBBEnclosed.C, and 11UnsteadyBBOpen.C.

Definition at line 91 of file UnsteadyBB.H.

◆ UModesWeighted

PtrList<volVectorField> UnsteadyBB::UModesWeighted

Pointers for perfoming Nested-POD Method List of pointers used to form the Weighted nested-POD modes for velocity field.

Definition at line 179 of file UnsteadyBB.H.

◆ W_matrix

Eigen::MatrixXd UnsteadyBB::W_matrix

Mass Matrix - energy equation.

Definition at line 125 of file UnsteadyBB.H.

◆ Y_matrix

Eigen::MatrixXd UnsteadyBB::Y_matrix

Diffusive term - energy equation.

Definition at line 118 of file UnsteadyBB.H.


The documentation for this class was generated from the following files: