Introduction
The problem consists of steady Navier-Stokes problem with parametrized viscosity. The physical problem is the backward facing step depicted in the following image:

At the inlet a uniform and constant velocity equal to 1 m/s is prescribed.
A detailed look into the code
In this section are explained the main steps necessary to construct the tutorial N°3.
The necessary header files
First of all let's have a look to the header files that need to be included and what they are responsible for.
The header files of ITHACA-FV necessary for this tutorial are the following:
#include "forces.H"
#include "IOmanip.H"
Header file of the ITHACAPOD class.
Header file of the ITHACAstream class, it contains the implementation of several methods for input ou...
Header file of the reducedSteadyNS class.
Header file of the steadyNS class.
Implementation of the tutorial03 class
We can define the tutorial03 class as a child of the steadyNS class
Implementation of a parametrized full order steady NS problem and preparation of the the reduced ma...
The members of the class are the fields that need to be manipulated during the resolution of the problem.
volVectorField& U;
volScalarField& p;
Inside the class it is defined the offlineSolve method according to the specific parametrized problem that needs to be solved.
void offlineSolve()
Perform an Offline solve.
If the offline solve has already been performed, then read the existing snapshots:
{
}
Eigen::MatrixXd mu_samples
Matrix of parameters to be used for PODI, where each row corresponds to a sample point....
bool offline
Boolean variable, it is 1 if the Offline phase has already been computed, else 0.
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.
List< Eigen::MatrixXd > readMatrix(word folder, word mat_name)
Read a three dimensional matrix from a txt file in Eigen format.
void read_fields(PtrList< GeometricField< Type, PatchField, GeoMesh > > &Lfield, word Name, fileName casename, int first_snap, int n_snap)
Function to read a list of fields from the name of the field and casename.
Else, perform the offline solve, which consists of a loop over all the parameters:
Vector<double>
Uinl(0, 0, 0);
for (label i = 0; i <
mu.cols(); i++)
{
}
Eigen::MatrixXd mu
Row matrix of parameters.
void truthSolve()
Perform a TruthSolve.
void change_viscosity(double mu)
Function to change the viscosity.
void restart()
set U and P back to the values into the 0 folder
autoPtr< volVectorField > Uinl
Initial dummy field with all Dirichlet boundary conditions.
See also the steadyNS class for the definition of the methods.
Definition of the main function
This tutorial is divided into Offline and Online stage, as can be seen from the main:
if (example._args().get("stage").match("offline"))
{
offline_stage(example);
}
else if (example._args().get("stage").match("online"))
{
offline_stage(example);
online_stage(example);
}
else
{
std::cout << "Pass '-stage offline', '-stage online'" << std::endl;
}
Pass "Offline" or "Online" as arguments from the command line to run either the offline or the online stage:
An example of type tutorial03 is constructed here:
Offline stage
The inlet boundary is initially set:
example.inletIndex.resize(1, 2);
example.inletIndex(0, 0) = 0;
example.inletIndex(0, 1) = 0;
and the offline stage is performed:
Then, the supremizer problem is solved:
example.solvesupremizer();
In order to show the functionality of reading fields, in this case the lifting function is read from a precomputed simulation with a unitary inlet velocity:
void normalizeFields(PtrList< GeometricField< Type, fvPatchField, volMesh > > &fields)
Normalize list of Geometric fields.
Then the snapshots matrix is homogenized:
example.computeLift(example.Ufield, example.liftfield, example.Uomfield);
and the modes for velocity, pressure and supremizers are obtained:
example.podex, 0, 0, NmodesUout);
example.podex, 0, 0,
NmodesPout);
example.podex,
example.supex, 1, NmodesSUPout);
void getModes(PtrList< GeometricField< Type, PatchField, GeoMesh > > &snapshots, PtrList< GeometricField< Type, PatchField, GeoMesh > > &modes, word fieldName, bool podex, bool supex, bool sup, label nmodes, bool correctBC)
Computes the bases or reads them for a field.
Then, the projection onto the POD modes is performed with:
example.projectSUP("./Matrices", NmodesUproj, NmodesPproj, NmodesSUPproj);
Online stage
The modes and reduced matrices computed during the offline stage are loaded
The reduced object is constructed:
Class where it is implemented a reduced problem for the steady Navier-stokes problem.
and the online solve is performed for some values of the viscosity:
word filename("./parOnline");
The vel_now matrix in this case is not used since there are no parametrized boundary conditions.
The viscosity is set with the command:
reduced.nu = example.mu(k, 0)
Finally, the online solution stored during the online solve is exported to file in three different formats with the lines:
"./ITHACAoutput/red_coeff");
"./ITHACAoutput/red_coeff");
"./ITHACAoutput/red_coeff");
void exportMatrix(Eigen::Matrix< T, -1, dim > &matrix, word Name, word type, word folder)
Export the reduced matrices in numpy (type=python), matlab (type=matlab) and txt (type=eigen) format ...
and the online solution is reconstructed and exported to file:
reduced.reconstruct(true, "./ITHACAoutput/Reconstruction/");
Parallel run
To speed up the offline stage, ITHACA-FV employs OpenFOAM facilities to run applications in parallel on distributed processors: the method of parallel computing used by OpenFOAM is known as domain decomposition.
First, the domain is decomposed as indicated in the directory system/decomposeParDict, with the command
In the parallel run also the lift field needs to be split with
decomposePar -case lift -time 0,1
Then, the offline solve is performed in parallel, evaluating the modes and the reduced matrices on the whole domain
mpirun -np 4 -quiet 03steadyNS -parallel -stage offline
The online stage is performed analogously
mpirun -np 4 -quiet 03steadyNS -parallel -stage online
To run the tutorial in parallel, the users can simple run the script Allrun_parallel.
In the case of OF1812, the parallel run is not supported.
The plain code can be found here.