Loading...
Searching...
No Matches
23burgers.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 a heat transfer Reduction Problem
26SourceFiles
27 02thermalBlock.C
28\*---------------------------------------------------------------------------*/
29
30#include <iostream>
31#include "fvCFD.H"
32#include "IOmanip.H"
33#include "Time.H"
34#include "Burgers.H"
35#include "ITHACAPOD.H"
36#include <Eigen/Dense>
37#define _USE_MATH_DEFINES
38#include <cmath>
39
43class tutorial23: public Burgers
44{
45 public:
46 explicit tutorial23(int argc, char* argv[])
47 :
48 Burgers(argc, argv),
49 U(_U())
50 {}
51
53 volVectorField& U;
54
55 void offlineSolve(word folder = "./ITHACAoutput/Offline/")
56 {
57 if (offline)
58 {
59 ITHACAstream::readMiddleFields(Ufield, U, "./ITHACAoutput/Offline/");
60 }
61 else
62 {
63 truthSolve(folder);
64 }
65 }
66
67
68};
69
70
71int main(int argc, char* argv[])
72{
73 // Create the train object of the tutorial02 type
74 tutorial23 train(argc, argv);
75 // Read some parameters from file
77 train._runTime());
78 int NmodesUout = para->ITHACAdict->lookupOrDefault<int>("NmodesUout", 15);
79 train.offlineSolve();
80 ITHACAPOD::getModes(train.Ufield, train.Umodes, train._U().name(),
81 train.podex, 0, 0,
82 NmodesUout);
83}
int main(int argc, char *argv[])
Definition 23burgers.C:71
Header file of the Burgers class.
Header file of the ITHACAPOD class.
Implementation of a parametrized full order Burgers and preparation of the the reduced matrices for...
Definition Burgers.H:58
autoPtr< fvMesh > _mesh
Mesh.
Definition Burgers.H:85
autoPtr< volVectorField > _U
Velocity field.
Definition Burgers.H:79
PtrList< volVectorField > Ufield
List of pointers used to form the velocity snapshots matrix.
Definition Burgers.H:73
autoPtr< Time > _runTime
Time.
Definition Burgers.H:88
volVectorModes Umodes
List of pointers used to form the velocity modes.
Definition Burgers.H:76
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.
bool offline
Boolean variable, it is 1 if the Offline phase has already been computed, else 0.
bool podex
Boolean variable, it is 1 if the POD has already been computed, else 0.
void truthSolve()
Perform a TruthSolve.
Class where the tutorial number 2 is implemented.
Definition 23burgers.C:44
void offlineSolve(word folder="./ITHACAoutput/Offline/")
Definition 23burgers.C:55
volVectorField & U
Velocity field.
Definition 23burgers.C:53
tutorial23(int argc, char *argv[])
Definition 23burgers.C:46
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
void readMiddleFields(PtrList< GeometricField< Type, PatchField, GeoMesh > > &Lfield, GeometricField< Type, PatchField, GeoMesh > &field, fileName casename)
Funtion to read a list of volVectorField from name of the field including all the intermediate snapsh...