Loading...
Searching...
No Matches
UnsteadyNSExplicit.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) 2020 by the ITHACA-FV authors
11-------------------------------------------------------------------------------
12
13 License
14 This file is part of ITHACA-FV
15
16 ITHACA-FV is free software: you can redistribute it and/or modify
17 it under the terms of the GNU Lesser General Public License as published by
18 the Free Software Foundation, either version 3 of the License, or
19 (at your option) any later version.
20
21 ITHACA-FV is distributed in the hope that it will be useful,
22 but WITHOUT ANY WARRANTY; without even the implied warranty of
23 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
24 GNU Lesser General Public License for more details.
25
26 You should have received a copy of the GNU Lesser General Public License
27 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
28
29\*---------------------------------------------------------------------------*/
30
31#include "UnsteadyNSExplicit.H"
32#include "fvCFD.H"
33
36
37// * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * * * //
38
39// Construct Null
41
42// Construct from zero
44{
45 _args = autoPtr<argList>
46 (
47 new argList(argc, argv)
48 );
49
50 if (!_args->checkRootCase())
51 {
52 Foam::FatalError.exit();
53 }
54
55 argList& args = _args();
56#include "createTime.H"
57#include "createMesh.H"
58 ITHACAdict = new IOdictionary
59 (
60 IOobject
61 (
62 "ITHACAdict",
63 runTime.system(),
64 mesh,
65 IOobject::MUST_READ,
66 IOobject::NO_WRITE
67 )
68 );
69#include "createFields.H"
71 bcMethod = ITHACAdict->lookupOrDefault<word>("bcMethod", "none");
72 M_Assert(bcMethod == "lift" || bcMethod == "penalty" || bcMethod == "none",
73 "The BC method must be set to lift or penalty or none in ITHACAdict");
74 fluxMethod = ITHACAdict->lookupOrDefault<word>("fluxMethod", "inconsistent");
75 M_Assert(fluxMethod == "inconsistent" || fluxMethod == "consistent",
76 "The flux method must be set to inconsistent or consistent in ITHACAdict");
80}
81
82// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
83
84void UnsteadyNSExplicit::truthSolve(List<scalar> mu_now, fileName folder)
85{
86 Time& runTime = _runTime();
87 fvMesh& mesh = _mesh();
88 volScalarField& p = _p();
89 volVectorField& U = _U();
90 surfaceScalarField& phi = _phi();
91
92 if (fluxMethod == "inconsistent")
93 {
94 phi = fvc::flux(U);
95 }
96
97#include "initContinuityErrs.H"
98 dimensionedScalar nu = _nu() * mu_now[0];
99 dimensionedScalar dt = timeStep * _dt();
100 IOMRFZoneList& MRF = _MRF();
101 singlePhaseTransportModel& laminarTransport = _laminarTransport();
102 instantList Times = runTime.times();
103 runTime.setEndTime(finalTime);
104 // Perform a TruthSolve
105 runTime.setTime(Times[1], 1);
106 runTime.setDeltaT(timeStep);
108 // Export and store the initial conditions for velocity, pressure and flux
109 ITHACAstream::exportSolution(U, name(counter), folder);
110 ITHACAstream::exportSolution(p, name(counter), folder);
112 std::ofstream of(folder + name(counter) + "/" +
113 runTime.timeName());
114 Ufield.append(U.clone());
115 Pfield.append(p.clone());
116 Phifield.append(phi.clone());
117 counter++;
119
120 // Start the time loop
121 while (runTime.run())
122 {
123#include "readTimeControls.H"
124#include "CourantNo.H"
125#include "setDeltaT.H"
126 runTime.setEndTime(finalTime);
127 runTime++;
128 Info << "Time = " << runTime.timeName() << nl << endl;
129
130 if (fluxMethod == "inconsistent")
131 {
132#include "IFM.H"
133 }
134 else if (fluxMethod == "consistent")
135 {
136#include "CFM.H"
137 }
138
139#include "initContinuityErrs.H"
140 Info << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
141 << " ClockTime = " << runTime.elapsedClockTime() << " s"
142 << nl << endl;
143
144 if (checkWrite(runTime))
145 {
146 ITHACAstream::exportSolution(U, name(counter), folder);
147 ITHACAstream::exportSolution(p, name(counter), folder);
149 std::ofstream of(folder + name(counter) + "/" +
150 runTime.timeName());
151 Ufield.append(U.clone());
152 Pfield.append(p.clone());
153 Phifield.append(phi.clone());
154 counter++;
156 }
157 }
158}
159
160
161
162
163
164
Foam::fvMesh & mesh
Definition createMesh.H:47
Foam::Time & runTime
Definition createTime.H:33
#define M_Assert(Expr, Msg)
Header file of the UnsteadyNSExplicit class.
static ITHACAparameters * getInstance()
Gets an instance of ITHACAparameters, to be used if the instance is already existing.
UnsteadyNSExplicit()
Construct Null.
autoPtr< dimensionedScalar > _dt
dimensionedScalar dt;
autoPtr< dimensionedScalar > _nu
dimensionedScalar nu;
bool checkWrite(Time &timeObject)
Function to check if the solution must be exported.
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 nextWrite
Auxiliary variable to store the next writing instant.
scalar finalTime
Final time (final time of the simulation and consequently of the acquisition of the snapshots)
label counter
Counter used for the output of the full order solutions.
bool offline
Boolean variable, it is 1 if the Offline phase has already been computed, else 0.
IOdictionary * ITHACAdict
dictionary to store input output infos
autoPtr< argList > _args
argList
bool podex
Boolean variable, it is 1 if the POD has already been computed, else 0.
void truthSolve()
Perform a TruthSolve.
bool supex
Boolean variable to check the existence of the supremizer modes.
Definition steadyNS.H:256
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
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< singlePhaseTransportModel > _laminarTransport
Laminar transport (used by turbulence model)
Definition steadyNS.H:302
PtrList< surfaceScalarField > Phifield
List of pointers used to form the flux snapshots matrix.
Definition steadyNS.H:95
word fluxMethod
Flux Method.
Definition steadyNS.H:320
ITHACAparameters * para
Definition steadyNS.H:82
autoPtr< IOMRFZoneList > _MRF
MRF variable.
Definition steadyNS.H:305
autoPtr< volVectorField > _U
Velocity field.
Definition steadyNS.H:263
word bcMethod
Boundary Method.
Definition steadyNS.H:317
autoPtr< volScalarField > _p
Pressure field.
Definition steadyNS.H:260
dimensionedScalar & nu
void exportSolution(GeometricField< Type, PatchField, GeoMesh > &s, fileName subfolder, fileName folder, word fieldName)
Export a field to file in a certain folder and subfolder.
bool check_pod()
Check if the POD data folder "./ITHACAoutput/POD" exists.
bool check_off()
Check if the offline data folder "./ITHACAoutput/Offline" exists.
bool check_sup()
Check if the supremizer folder exists.
surfaceScalarField & phi
volVectorField & U
volScalarField & p
singlePhaseTransportModel & laminarTransport