Loading...
Searching...
No Matches
14DMDexample.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 an unsteady NS problem reproduced with DMD
26SourceFiles
27 14DMDexample.C
28\*---------------------------------------------------------------------------*/
29
30#include "unsteadyNS.H"
31#include "ITHACAPOD.H"
32#include "ITHACADMD.H"
33#include "ReducedUnsteadyNS.H"
34#include "ITHACAstream.H"
35#include <chrono>
36#include<math.h>
37#include<iomanip>
38#include "redsvd"
39
41{
42 public:
43 explicit tutorial14(int argc, char* argv[])
44 :
45 unsteadyNS(argc, argv),
46 U(_U()),
47 p(_p())
48 {}
49
50 // Fields To Perform
51 volVectorField& U;
52 volScalarField& p;
53
55 {
56 List<scalar> mu_now(1);
57
58 if (offline)
59 {
60 ITHACAstream::read_fields(Ufield, U, "./ITHACAoutput/Offline/");
61 ITHACAstream::read_fields(Pfield, p, "./ITHACAoutput/Offline/");
62 }
63 else
64 {
65 for (label i = 0; i < mu.cols(); i++)
66 {
67 mu_now[0] = mu(0, i);
68 change_viscosity( mu(0, i));
69 truthSolve(mu_now);
70 }
71 }
72 }
73};
74
75/*---------------------------------------------------------------------------*\
76 Starting the MAIN
77\*---------------------------------------------------------------------------*/
78int main(int argc, char* argv[])
79{
80 // Construct the tutorial04 object
81 tutorial14 example(argc, argv);
82 // Read parameters from ITHACAdict file
84 example._runTime());
85 // Read some input parameters for DMD
86 double dtDMD = readScalar(para->ITHACAdict->lookup("dtDMD"));
87 double finalTimeDMD = readScalar(para->ITHACAdict->lookup("finalTimeDMD"));
88 double startTimeDMD = readScalar(para->ITHACAdict->lookup("startTimeDMD"));
89 int numberOfModesDMD = readInt(para->ITHACAdict->lookup("numberOfModesDMD"));
90 bool exactDMD = readBool(para->ITHACAdict->lookup("exactDMD"));
91 bool exportDMDModes = readBool(para->ITHACAdict->lookup("exportDMDModes"));
92 word exportFolder = string(para->ITHACAdict->lookup("exportFolder"));
94 example.Pnumber = 1;
96 example.Tnumber = 1;
98 example.setParameters();
99 // Set the parameter ranges
100 example.mu_range(0, 0) = 0.01;
101 example.mu_range(0, 1) = 0.01;
102 // Generate equispaced samples inside the parameter range
103 example.genEquiPar();
104 // Time parameters
105 example.startTime = 60;
106 example.finalTime = 85;
107 example.timeStep = 0.01;
108 example.writeEvery = 0.1;
109 // Read the snapshots from Offline folder or compute new snapshots
110 example.offlineSolve();
111 // Create The DMD class for the velocity field
112 word exportFieldNameU = example.Ufield[0].name() + "_" + name(numberOfModesDMD);
113 ITHACADMDvolVector DMDv(example.Ufield, example.writeEvery);
114 DMDv.getModes(numberOfModesDMD, exactDMD, exportDMDModes);
115 DMDv.getDynamics(startTimeDMD, finalTimeDMD, dtDMD);
116 DMDv.reconstruct(exportFolder, exportFieldNameU);
117 // Create The DMD class for the pressure field
118 word exportFieldNameP = example.Pfield[0].name() + "_" + name(numberOfModesDMD);
119 ITHACADMDvolScalar DMDp(example.Pfield, example.writeEvery);
120 DMDp.getModes(numberOfModesDMD, exactDMD, exportDMDModes);
121 DMDp.getDynamics(startTimeDMD, finalTimeDMD, dtDMD);
122 DMDp.reconstruct(exportFolder, exportFieldNameP);
123 return 0;
124}
125
int main(int argc, char *argv[])
Header file of the ITHACADMD class.
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 reducedUnsteadyNS class.
Class of the computation of the DMD, it exploits the SVD methods.
Definition ITHACADMD.H:70
void getDynamics(double tStart, double tFinal, double dt)
Export the dynamics of DMD given an initial time step, a final one and a time step.
Definition ITHACADMD.C:165
void reconstruct(word exportFolder, word fieldName)
Reconstruct and export the solution using the computed dynamics.
Definition ITHACADMD.C:194
void getModes(label SVD_rank=-1, bool exact=true, bool exportDMDmodes=false)
Get the DMD modes.
Definition ITHACADMD.C:50
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 startTime
Start Time (initial time to start storing the snapshots)
scalar writeEvery
Time step of the writing procedure.
scalar timeStep
Time step of the simulation.
scalar finalTime
Final time (final time of the simulation and consequently of the acquisition of the snapshots)
label Pnumber
Number of parameters.
label Tnumber
Dimension of the training set (used only when gerating parameters without input)
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::MatrixXd mu_range
Range of the parameter spaces.
void setParameters()
Set Parameters Problems.
void truthSolve()
Perform a TruthSolve.
void genEquiPar()
Generate Equidistributed Numbers.
void change_viscosity(double mu)
Function to change the viscosity.
Definition steadyNS.C:2014
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
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
autoPtr< volVectorField > _U
Velocity field.
Definition steadyNS.H:263
autoPtr< volScalarField > _p
Pressure field.
Definition steadyNS.H:260
void offlineSolve()
volScalarField & p
volVectorField & U
tutorial14(int argc, char *argv[])
Implementation of a parametrized full order unsteady NS problem and preparation of the the reduced ...
Definition unsteadyNS.H:62
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 unsteadyNS class.