37#ifndef ReducedUnsteadyNSTurb_H
38#define ReducedUnsteadyNSTurb_H
46#include <unsupported/Eigen/NonLinearOptimization>
47#include <unsupported/Eigen/NumericalDiff>
64 int operator()(
const Eigen::VectorXd& x, Eigen::VectorXd& fvec)
const;
65 int df(
const Eigen::VectorXd& x, Eigen::MatrixXd& fjac)
const;
79 std::vector<SPLINTER::RBFSpline*>
SPLINES;
98 int operator()(
const Eigen::VectorXd& x, Eigen::VectorXd& fvec)
const;
99 int df(
const Eigen::VectorXd& x, Eigen::MatrixXd& fjac)
const;
131 int operator()(
const Eigen::VectorXd& x, Eigen::VectorXd& fvec)
const;
132 int df(
const Eigen::VectorXd& x, Eigen::MatrixXd& fjac)
const;
166 int operator()(
const Eigen::VectorXd& x, Eigen::VectorXd& fvec)
const;
167 int df(
const Eigen::VectorXd& x, Eigen::MatrixXd& fjac)
const;
301 fileName folder =
"./online_rec");
Header file of the reducedSteadyNS class.
Header file of the reducedUnsteadyNS class.
Header file of the UnsteadyNSTurb class.
Class where it is implemented a reduced problem for the unsteady Navier-stokes problem.
void solveOnlinePPEAve(Eigen::MatrixXd velNow)
Method to perform an online solve using a PPE stabilisation method with the use of the average splitt...
ReducedUnsteadyNSTurb()
Construct Null.
int nphiNut
Number of viscosity modes.
newtonUnsteadyNSTurbSUP newtonObjectSUP
Function object to call the non linear solver sup approach.
newtonUnsteadyNSTurbSUPAve newtonObjectSUPAve
Function object to call the non linear solver sup approach with the splitting of the eddy viscosity.
PtrList< volScalarField > nutRecFields
Reconstructed eddy viscosity fields list.
void solveOnlineSUPAve(Eigen::MatrixXd velNow)
Method to perform an online solve using a supremizer stabilisation method with the use of the average...
Eigen::MatrixXd rbfCoeffMat
The matrix of the eddy viscosity RBF interoplated coefficients.
bool skipLift
Interpolation boolean variable to skip lifting functions.
int dimA
Dimension of the interpolation independent variable.
Eigen::VectorXd muStar
Online parameter value.
Eigen::MatrixXd initCond
The matrix of the initial velocity and pressure reduced coefficients.
Eigen::VectorXd gNutAve
The reduced average vector for viscosity coefficients.
void solveOnlineSUP(Eigen::MatrixXd velNow)
Method to perform an online solve using a supremizer stabilisation method.
int interChoice
Interpolation independent variable choice.
newtonUnsteadyNSTurbPPE newtonObjectPPE
Function object to call the non linear solver PPE approach.
void reconstruct(bool exportFields=false, fileName folder="./online_rec")
Method to reconstruct the solutions from an online solve with any of the two techniques SUP or the PP...
PtrList< volScalarField > nutModes
List of pointers to store the modes for the eddy viscosity.
void solveOnlinePPE(Eigen::MatrixXd velNow)
Method to perform an online solve using a PPE stabilisation method.
UnsteadyNSTurb * problem
Pointer to the FOM problem.
newtonUnsteadyNSTurbPPEAve newtonObjectPPEAve
Function object to call the non linear solver PPE approach with the splitting of the eddy viscosity.
Eigen::MatrixXd setOnlineVelocity(Eigen::MatrixXd vel)
Sets the online velocity.
Eigen::VectorXd nut0
The initial eddy viscosity reduced coefficients.
Implementation of a parametrized full order unsteady NS problem and preparation of the the reduced ...
Template object created to solve non linear problems using the Eigen library.
Class where it is implemented a reduced problem for the unsteady Navier-stokes problem.
int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec) const
int df(const Eigen::VectorXd &x, Eigen::MatrixXd &fjac) const
newtonUnsteadyNSTurbPPEAve()
newtonUnsteadyNSTurbPPEAve(int Nx, int Ny, UnsteadyNSTurb &problem)
std::vector< SPLINTER::RBFSpline * > SPLINES
int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec) const
newtonUnsteadyNSTurbPPE()
newtonUnsteadyNSTurbPPE(int Nx, int Ny, UnsteadyNSTurb &problem)
int df(const Eigen::VectorXd &x, Eigen::MatrixXd &fjac) const
std::vector< SPLINTER::RBFSpline * > SPLINES
int df(const Eigen::VectorXd &x, Eigen::MatrixXd &fjac) const
std::vector< SPLINTER::RBFSpline * > SPLINES
int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec) const
newtonUnsteadyNSTurbSUPAve()
newtonUnsteadyNSTurbSUPAve(int Nx, int Ny, UnsteadyNSTurb &problem)
std::vector< SPLINTER::RBFSpline * > SPLINES
int df(const Eigen::VectorXd &x, Eigen::MatrixXd &fjac) const
newtonUnsteadyNSTurbSUP()
newtonUnsteadyNSTurbSUP(int Nx, int Ny, UnsteadyNSTurb &problem)
int operator()(const Eigen::VectorXd &x, Eigen::VectorXd &fvec) const