Loading...
Searching...
No Matches
12simpleSteadyNS.C
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/>.
24Description
25 Example of steady NS Reduction Problem solved by the use of the SIMPLE algorithm
26SourceFiles
27 12simpleSteadyNS.C
28\*---------------------------------------------------------------------------*/
29
30#include "SteadyNSSimple.H"
31#include "ITHACAstream.H"
32#include "ITHACAPOD.H"
34#include "forces.H"
35#include "IOmanip.H"
36
37
39{
40 public:
42 explicit tutorial12(int argc, char* argv[])
43 :
44 SteadyNSSimple(argc, argv),
45 U(_U()),
46 p(_p()),
47 phi(_phi())
48 {}
49
51 volVectorField& U;
53 volScalarField& p;
55 surfaceScalarField& phi;
56
59 {
60 Vector<double> inl(0, 0, 0);
61 List<scalar> mu_now(1);
62
63 // if the offline solution is already performed read the fields
64 if (offline)
65 {
66 ITHACAstream::read_fields(Ufield, U, "./ITHACAoutput/Offline/");
67 ITHACAstream::read_fields(Pfield, p, "./ITHACAoutput/Offline/");
69 ITHACAstream::readMatrix("./ITHACAoutput/Offline/mu_samples_mat.txt");
70 }
71 else
72 {
73 Vector<double> Uinl(1, 0, 0);
74 label BCind = 0;
75
76 for (label i = 0; i < mu.cols(); i++)
77 {
78 mu_now[0] = mu(0, i);
80 assignIF(U, Uinl);
81 truthSolve2(mu_now);
82 }
83 }
84 }
85
86};
87
88int main(int argc, char* argv[])
89{
90 // Construct the tutorial object
91 tutorial12 example(argc, argv);
92 // Read some parameters from file
94 example._runTime());
95 int NmodesUout = para->ITHACAdict->lookupOrDefault<int>("NmodesUout", 15);
96 int NmodesPout = para->ITHACAdict->lookupOrDefault<int>("NmodesPout", 15);
97 int NmodesUproj = para->ITHACAdict->lookupOrDefault<int>("NmodesUproj", 10);
98 int NmodesPproj = para->ITHACAdict->lookupOrDefault<int>("NmodesPproj", 10);
99 // Read the par file where the parameters are stored
100 word filename("./par");
101 example.mu = ITHACAstream::readMatrix(filename);
102 // Set the inlet boundaries patch 0 directions x and y
103 example.inletIndex.resize(1, 2);
104 example.inletIndex(0, 0) = 0;
105 example.inletIndex(0, 1) = 0;
106 // Perform the offline solve
107 example.offlineSolve();
108 ITHACAstream::read_fields(example.liftfield, example.U, "./lift/");
109 // Homogenize the snapshots
110 example.computeLift(example.Ufield, example.liftfield, example.Uomfield);
111 // Perform POD on velocity and pressure and store the first 10 modes
112 ITHACAPOD::getModes(example.Uomfield, example.Umodes, example._U().name(),
113 example.podex, 0, 0,
114 NmodesUout);
115 ITHACAPOD::getModes(example.Pfield, example.Pmodes, example._p().name(),
116 example.podex, 0, 0,
117 NmodesPout);
118 // Create the reduced object
119 reducedSimpleSteadyNS reduced(example);
120 PtrList<volVectorField> U_rec_list;
121 PtrList<volScalarField> P_rec_list;
122 // Reads inlet volocities boundary conditions.
123 word vel_file(para->ITHACAdict->lookup("online_velocities"));
124 Eigen::MatrixXd vel = ITHACAstream::readMatrix(vel_file);
125
126 //Perform the online solutions
127 for (label k = 0; k < (example.mu).size(); k++)
128 {
129 scalar mu_now = example.mu(0, k);
130 example.change_viscosity(mu_now);
131 reduced.setOnlineVelocity(vel);
132 reduced.solveOnline_Simple(mu_now, NmodesUproj, NmodesPproj);
133 }
134
135 exit(0);
136}
int main(int argc, char *argv[])
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.
Class for the definition of some general parameters, the parameters must be defined from the file ITH...
IOdictionary * ITHACAdict
Dictionary for input objects from file.
static ITHACAparameters * getInstance()
Gets an instance of ITHACAparameters, to be used if the instance is already existing.
Implementation of a parametrized full order steady NS problem and preparation of the the reduced ma...
void truthSolve2(List< scalar > mu_now, word Folder="./ITHACAoutput/Offline/")
Offline solver for the whole Navier-Stokes problem.
Class where it is implemented a reduced problem for the steady Navier-stokes problem.
void solveOnline_Simple(scalar mu_now, int NmodesUproj, int NmodesPproj, int NmodesNut=0, int NmodesSup=0, word Folder="./ITHACAoutput/Reconstruct/")
Method to perform an online solve using a PPE stabilisation method.
void setOnlineVelocity(Eigen::MatrixXd vel)
It checks if the number of imposed boundary conditions is correct and set the inlet velocity equal to...
Eigen::MatrixXd mu_samples
Matrix of parameters to be used for PODI, where each row corresponds to a sample point....
void assignIF(T &s, G &value)
Assign internal field condition.
void computeLift(T &Lfield, T &liftfield, T &omfield)
Homogenize the snapshot matrix, it works with PtrList of volVectorField and volScalarField.
bool offline
Boolean variable, it is 1 if the Offline phase has already been computed, else 0.
Eigen::MatrixXd mu
Row matrix of parameters.
Eigen::MatrixXi inletIndex
Matrix that contains informations about the inlet boundaries.
bool podex
Boolean variable, it is 1 if the POD has already been computed, else 0.
void change_viscosity(double mu)
Function to change the viscosity.
Definition steadyNS.C:2014
autoPtr< surfaceScalarField > _phi
Flux.
Definition steadyNS.H:293
PtrList< volScalarField > Pfield
List of pointers used to form the pressure snapshots matrix.
Definition steadyNS.H:86
autoPtr< Time > _runTime
Time.
Definition steadyNS.H:290
volVectorModes Umodes
List of pointers used to form the velocity modes.
Definition steadyNS.H:101
PtrList< volVectorField > Ufield
List of pointers used to form the velocity snapshots matrix.
Definition steadyNS.H:89
autoPtr< fvMesh > _mesh
Mesh.
Definition steadyNS.H:281
PtrList< volVectorField > liftfield
List of pointers used to form the list of lifting functions.
Definition steadyNS.H:110
PtrList< volVectorField > Uomfield
List of pointers used to form the homogeneous velocity snapshots.
Definition steadyNS.H:113
autoPtr< volVectorField > Uinl
Initial dummy field with all Dirichlet boundary conditions.
Definition steadyNS.H:272
autoPtr< volVectorField > _U
Velocity field.
Definition steadyNS.H:263
volScalarModes Pmodes
List of pointers used to form the pressure modes.
Definition steadyNS.H:98
autoPtr< volScalarField > _p
Pressure field.
Definition steadyNS.H:260
volScalarField & p
Pressure field.
void offlineSolve()
Perform an Offline solve.
tutorial12(int argc, char *argv[])
Constructor.
surfaceScalarField & phi
volVectorField & U
Velocity field.
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.
Definition ITHACAPOD.C:93
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.
label i
Definition pEqn.H:46