Loading...
Searching...
No Matches
07NonIsothermalMixingElbow.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2Copyright (C) 2017 by the ITHACA-FV authors
3
4License
5 This file is part of ITHACA-FV
6
7 ITHACA-FV is free software: you can redistribute it and/or modify
8 it under the terms of the GNU Lesser General Public License as published by
9 the Free Software Foundation, either version 3 of the License, or
10 (at your option) any later version.
11
12 ITHACA-FV is distributed in the hope that it will be useful,
13 but WITHOUT ANY WARRANTY; without even the implied warranty of
14 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 GNU Lesser General Public License for more details.
16
17 You should have received a copy of the GNU Lesser General Public License
18 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
19
20Description
21 Example of NS-Stokes and heat transport equation Reduction Problem
22SourceFiles
23 06NonIsothermalMixingElbow.C
24
25
26\*---------------------------------------------------------------------------*/
27
28#include "unsteadyNST.H"
29#include "ITHACAPOD.H"
30#include "ReducedUnsteadyNS.H"
31#include "ReducedUnsteadyNST.H"
32#include "ITHACAstream.H"
33#include <chrono>
34#include<math.h>
35#include<iomanip>
36
38{
39 public:
40 explicit tutorial07(int argc, char* argv[])
41 :
42 unsteadyNST(argc, argv),
43 U(_U()),
44 p(_p()),
45 T(_T())
46 {}
47
48 // Fields To Perform
49 volVectorField& U;
50 volScalarField& p;
51 volScalarField& T;
52
54 {
55 Vector<double> inl(0, 0, 0);
56 List<scalar> mu_now(1);
57 Info << "here" << endl;
58
59 if (offline)
60 {
61 ITHACAstream::read_fields(Ufield, U, "./ITHACAoutput/Offline/");
62 ITHACAstream::read_fields(Pfield, p, "./ITHACAoutput/Offline/");
63 ITHACAstream::read_fields(Tfield, T, "./ITHACAoutput/Offline/");
65 ITHACAstream::readMatrix("./ITHACAoutput/Offline/mu_samples_mat.txt");
66 }
67 else
68 {
69 for (label i = 0; i < mu.cols(); i++)
70 {
71 mu_now[0] = mu(0, i);
72 truthSolve(mu_now);
73 }
74 }
75 }
76};
77
78
79/*---------------------------------------------------------------------------*\
80 Starting the MAIN
81\*---------------------------------------------------------------------------*/
82int main(int argc, char* argv[])
83{
84 // Construct the tutorial object
85 tutorial07 example(argc, argv);
86 // Read some parameters from file
88 example._runTime());
89 int NmodesUout = para->ITHACAdict->lookupOrDefault<int>("NmodesUout", 5);
90 int NmodesPout = para->ITHACAdict->lookupOrDefault<int>("NmodesPout", 5);
91 int NmodesTout = para->ITHACAdict->lookupOrDefault<int>("NmodesTout", 5);
92 int NmodesSUPout = para->ITHACAdict->lookupOrDefault<int>("NmodesSUPout", 5);
93 int NmodesUproj = para->ITHACAdict->lookupOrDefault<int>("NmodesUproj", 5);
94 int NmodesPproj = para->ITHACAdict->lookupOrDefault<int>("NmodesPproj", 15);
95 int NmodesTproj = para->ITHACAdict->lookupOrDefault<int>("NmodesTproj", 5);
96 int NmodesSUPproj = para->ITHACAdict->lookupOrDefault<int>("NmodesSUPproj", 5);
98 example.Pnumber = 1;
100 example.Tnumber = 1;
102 example.setParameters();
103 // Set the parameter ranges
104 example.mu_range(0, 0) = 0.1;
105 example.mu_range(0, 1) = 0.1;
106 //Generate equispaced samples inside the parameter range
107 example.genEquiPar();
108 // Set the inlet boundaries where you have non homogeneous boundary conditions
109 example.inletIndex.resize(2, 2);
110 example.inletIndex << 3, 0, 2, 1;
111 example.inletIndexT.resize(3, 1);
112 example.inletIndexT << 3, 2, 0;
113 // Time parameters
114 example.startTime = 0;
115 example.finalTime = 50;
116 example.timeStep = 0.05;
117 example.writeEvery = 0.1;
118 // Perform The Offline Solve;
119 example.offlineSolve();
120 // Solve the supremizer problem
121 example.solvesupremizer();
122 // Search the lift function for the velocity
123 example.liftSolve();
124 // Search the lift function for the temperature
125 example.liftSolveT();
126 // Create homogeneous basis functions for velocity
127 example.computeLift(example.Ufield, example.liftfield, example.Uomfield);
128 // Create homogeneous basis functions for temperature
129 example.computeLiftT(example.Tfield, example.liftfieldT, example.Tomfield);
130 // Perform a POD decomposition for velocity temperature and pressure fields
131 ITHACAPOD::getModes(example.Uomfield, example.Umodes, example._U().name(),
132 example.podex, 0, 0,
133 NmodesUout);
134 ITHACAPOD::getModes(example.Pfield, example.Pmodes, example._p().name(),
135 example.podex, 0, 0,
136 NmodesPout);
137 ITHACAPOD::getModes(example.Tomfield, example.Tmodes, example._T().name(),
138 example.podex, 0, 0,
139 NmodesTout);
140 ITHACAPOD::getModes(example.supfield, example.supmodes, example._U().name(),
141 example.podex,
142 example.supex, 1, NmodesSUPout);
143 example.projectSUP("./Matrices", NmodesUproj, NmodesPproj, NmodesTproj,
144 NmodesSUPproj);
145 reducedUnsteadyNST reduced(example);
146 // Set values of the ridotto stuff
147 reduced.nu = 0.1;
148 reduced.tstart = 0;
149 reduced.finalTime = 50;
150 reduced.dt = 0.05;
151 reduced.DT = 1e-06;
152 // Set the online velocity
153 Eigen::MatrixXd vel_now;
154 vel_now.resize(2, 1);
155 vel_now << 0.6, 1.2;
156 // Set the online temperature and the value of the internal field
157 Eigen::MatrixXd temp_now;
158 temp_now.resize(3, 1);
159 temp_now << 60, 70, 60;
160 reduced.solveOnline_sup(vel_now, temp_now);
161 // Reconstruct the solution and export it
162 reduced.reconstruct_sup("./ITHACAoutput/ReconstructionSUP/", 2);
163 reduced.reconstruct_supt("./ITHACAoutput/ReconstructionSUP/", 2);
164 exit(0);
165}
166
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 reducedUnsteadyNST class.
Header file of the reducedUnsteadyNS 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.
scalar nu
Reduced viscosity in case of parametrized viscosity.
Class where it is implemented a reduced problem for the unsteady Navier-stokes weakly coupled with t...
scalar finalTime
Scalar to store the final time if the online simulation.
void reconstruct_supt(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
Method to reconstruct a solution from an online solve with a supremizer stabilisation technique for t...
scalar tstart
Scalar to store the final time if the online simulation.
void reconstruct_sup(fileName folder="./ITHACAOutput/online_rec", int printevery=1)
Method to reconstruct a solution from an online solve with a supremizer stabilisation technique.
scalar dt
Scalar to store the time increment.
void solveOnline_sup(Eigen::MatrixXd &vel_now, Eigen::MatrixXd &temp_now, int startSnap=0)
Method to perform an online solve using a supremizer stabilisation method.
Eigen::MatrixXi inletIndexT
label Pnumber
Number of parameters.
Eigen::MatrixXd mu_samples
Matrix of parameters to be used for PODI, where each row corresponds to a sample point....
label Tnumber
Dimension of the training set (used only when gerating parameters without input)
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.
void computeLiftT(T &Lfield, T &liftfield, T &omfield)
Virtual function to compute the lifting function.
Eigen::MatrixXd mu
Row matrix of parameters.
Eigen::MatrixXd mu_range
Range of the parameter spaces.
void setParameters()
Set Parameters Problems.
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 truthSolve()
Perform a TruthSolve.
void genEquiPar()
Generate Equidistributed Numbers.
bool supex
Boolean variable to check the existence of the supremizer modes.
Definition steadyNS.H:256
PtrList< volScalarField > Pfield
List of pointers used to form the pressure snapshots matrix.
Definition steadyNS.H:86
volVectorModes supmodes
List of pointers used to form the supremizer modes.
Definition steadyNS.H:104
void solvesupremizer(word type="snapshots")
solve the supremizer either with the use of the pressure snaphots or the pressure modes
Definition steadyNS.C:136
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
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
PtrList< volVectorField > supfield
List of pointers used to form the supremizer snapshots matrix.
Definition steadyNS.H:92
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
tutorial07(int argc, char *argv[])
Implementation of a parametrized full order unsteady NS problem weakly coupled with the energy equat...
Definition unsteadyNST.H:61
PtrList< volScalarField > Tfield
List of pointers used to form the temperature snapshots matrix.
Definition unsteadyNST.H:99
PtrList< volScalarField > liftfieldT
List of pointers used to form the list of the temperature lifting functions.
PtrList< volScalarField > Tomfield
List of pointers used to form the homogeneous temperature snapshots.
scalar timeStep
Time step of the simulation.
Definition unsteadyNST.H:81
void projectSUP(fileName folder, label NUmodes, label NPmodes, label NTmodes, label NSUPmodes)
Specific variable for the unstationary case.
autoPtr< Time > _runTime
Time.
scalar finalTime
Final time (final time of the simulation and consequently of the acquisition of the snapshots)
Definition unsteadyNST.H:78
scalar startTime
Start Time (initial time to start storing the snapshots)
Definition unsteadyNST.H:75
autoPtr< volScalarField > _T
Temperature field.
scalar writeEvery
Time step of the writing procedure.
Definition unsteadyNST.H:84
void liftSolve()
Perform a lift solve for velocity field.
autoPtr< fvMesh > _mesh
Mesh.
PtrList< volScalarField > Tmodes
List of pointers used to form the temperature modes.
void liftSolveT()
Perform a lift solve for temperature 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
Header file of the unsteadyNST class.