Loading...
Searching...
No Matches
NLsolve.H
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-------------------------------------------------------------------------------
12 License
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/>.
24\*---------------------------------------------------------------------------*/
25fv::options& fvOptions = _fvOptions();
26surfaceScalarField& phi = _phi();
27volScalarField& rho = _rho();
28fluidThermo& thermo = pThermo();
29volScalarField& psi = _psi();
31fvMesh& mesh = _mesh();
32dimensionedScalar& initialMass = _initialMass();
33bool closedVolume = false;
34
35scalar residual = 1;
36scalar eresidual = 1;
37Vector<double> uresidual_v(0, 0, 0);
38scalar presidual = 1;
39
41middleStep = para->ITHACAdict->lookupOrDefault<int>("middleStep", 20);
42label csolve = 0;
45label saver = 0;
46label folderN = 0;
47std::clock_t startOff;
48startOff = std::clock();
50
51// files to save data
52std::ofstream res_os;
53std::ofstream snaps_os;
54std::ofstream iters;
55std::ofstream res_U;
56std::ofstream res_P;
57std::ofstream res_E;
58std::ofstream cpuTimes;
59
61res_os.open(folder + "/residuals", std::ios_base::app);
62snaps_os.open(folder + "/snaps", std::ios_base::app);
63iters.open(folder + "/iters", std::ios_base::app);
64cpuTimes.open(folder + "/cpuTimes", std::ios_base::app);
66res_U.open(folder + name(counter) + "/res_U", std::ios_base::app);
67res_P.open(folder + name(counter) + "/res_P", std::ios_base::app);
68res_E.open(folder + name(counter) + "/res_E", std::ios_base::app);
69
70// Variable that can be changed
71// std::ofstream res_os;
72// res_os.open("./ITHACAoutput/Offline/residuals", std::ios_base::app);
73
74tolerance = ITHACAdict->lookupOrDefault<scalar>("tolerance", 1e-5);
75maxIter = ITHACAdict->lookupOrDefault<scalar>("maxIter", 2000);
76
77turbulence->validate();
78
79#if OFVER == 6
80
81while (simple.loop(runTime) && residual > tolerance && csolve < maxIter )
82#else
83while (simple.loop() && residual > tolerance && csolve < maxIter )
84#endif
85{
86 Info << "Time = " << runTime.timeName() << nl << endl;
87 //Momentum equation phase
88 getUmatrix(U);
89
90 if (simple.momentumPredictor())
91 {
92 uresidual_v = solve(Ueqn_global() == - getGradP(p)).initialResidual(); //Working
93 fvOptions.correct(U);
94 }
95
96 //Energy equation phase
97 getEmatrix(U, p);
98 eresidual = Eeqn_global().solve().initialResidual();
99 fvOptions.correct(thermo.he());
100 thermo.correct(); // Here are calculated both temperature and density based on P,U and he.
101 // Pressure equation phase
102 constrainPressure(p, rho, U, getPhiHbyA(Ueqn_global(), U, p),
103 getRhorAUf(Ueqn_global()));// Update the pressure BCs to ensure flux consistency
104 surfaceScalarField phiHbyACalculated = getPhiHbyA(Ueqn_global(), U, p);
105 closedVolume = adjustPhi(phiHbyACalculated, U, p);
106
107 while (simple.correctNonOrthogonal())
108 {
109 getPmatrix(Ueqn_global(), U, p);
110 presidual = Peqn_global().solve().initialResidual();
111
112 if (simple.finalNonOrthogonalIter())
113 {
114 phi = getPhiHbyA(Ueqn_global(), U, p) + Peqn_global().flux();
115 }
116 }
117
118#include "continuityErrs.H"
119 p.relax();// Explicitly relax pressure for momentum corrector
120 U = HbyA() - (1.0 / Ueqn_global().A()) * getGradP(p);
121 U.correctBoundaryConditions();
122 fvOptions.correct(U);
123 bool pLimited = pressureControl.limit(p);
124
125 // For closed-volume cases adjust the pressure and density levels to obey overall mass continuity
126 if (closedVolume)
127 {
128 p += (initialMass - fvc::domainIntegrate(psi * p))
129 / fvc::domainIntegrate(psi);
130 }
131
132 if (pLimited || closedVolume)
133 {
134 p.correctBoundaryConditions();
135 }
136
137 rho = thermo.rho(); // Here rho is calculated as p*psi = p/(R*T)
138 rho.relax();
139 auto uresidual = max(max(uresidual_v[0], uresidual_v[1]), uresidual_v[2]);
141 Info << "\nResidual: " << residual << endl << endl;
142 turbulence->correct();
143 csolve = csolve + 1;
144 Info << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
145 << " ClockTime = " << runTime.elapsedClockTime() << " s"
146 << nl << endl;
147 saver++;
148 std::cout << "saver ================== \t" << saver << std::endl;
149
150 if (middleExport && saver == middleStep)
151 {
152 saver = 0;
153 folderN++;
154 ITHACAstream::exportSolution(U, name(folderN), folder + name(counter));
155 ITHACAstream::exportSolution(p, name(folderN), folder + name(counter));
156 ITHACAstream::exportSolution(E, name(folderN), folder + name(counter));
157 auto nut = _mesh().lookupObject<volScalarField>("nut");
158 ITHACAstream::exportSolution(nut, name(folderN), folder + name(counter));
159 Ufield.append(U.clone());
160 Pfield.append(p.clone());
161 Efield.append(E.clone());
162 nutFields.append(nut.clone());
163 res_U << uresidual << std::endl;
164 res_P << presidual << std::endl;
165 res_E << eresidual << std::endl;
166 }
167}
168
169snaps_os << folderN + 1 << std::endl;
170//snaps_os << folderN << std::endl;
171iters << csolve << std::endl;
172res_os << residual << std::endl;
174cpuTimes << std::clock() - startOff << std::endl;
175res_os.close();
176res_U.close();
177res_P.close();
178snaps_os.close();
179iters.close();
180cpuTimes.close();
181
182if (middleExport)
183{
184 ITHACAstream::exportSolution(U, name(folderN + 1), folder + name(counter));
185 ITHACAstream::exportSolution(p, name(folderN + 1), folder + name(counter));
186 ITHACAstream::exportSolution(E, name(folderN + 1), folder + name(counter));
187 auto nut = _mesh().lookupObject<volScalarField>("nut");
188 ITHACAstream::exportSolution(nut, name(folderN + 1), folder + name(counter));
189}
190else
191{
192 ITHACAstream::exportSolution(U, name(counter), folder);
193 ITHACAstream::exportSolution(p, name(counter), folder);
194 ITHACAstream::exportSolution(E, name(counter), folder);
195 auto nut = _mesh().lookupObject<volScalarField>("nut");
196 ITHACAstream::exportSolution(nut, name(counter), folder);
197}
198
199Ufield.append(U.clone());
200Pfield.append(p.clone());
201Efield.append(E.clone());
202auto nut = _mesh().lookupObject<volScalarField>("nut");
203nutFields.append(nut.clone());
204
205// res_os << residual << std::endl;
206// res_os.close();
207runTime.setTime(runTime.startTime(), 0);
208
209
210
211
212
213
214
215
eresidual
Definition EEqn.H:21
Foam::fvMesh & mesh
Definition createMesh.H:47
Foam::Time & runTime
Definition createTime.H:33
scalar uresidual
std::ofstream res_P
Definition NLsolve.H:56
std::ofstream snaps_os
Definition NLsolve.H:53
fv::options & fvOptions
Definition NLsolve.H:25
tolerance
Definition NLsolve.H:74
middleStep
Definition NLsolve.H:41
label csolve
Definition NLsolve.H:42
dimensionedScalar & initialMass
Definition NLsolve.H:32
maxIter
Definition NLsolve.H:75
std::ofstream cpuTimes
Definition NLsolve.H:58
bool closedVolume
Definition NLsolve.H:33
std::ofstream iters
Definition NLsolve.H:54
std::clock_t startOff
Definition NLsolve.H:47
std::ofstream res_U
Definition NLsolve.H:55
volScalarField & psi
Definition NLsolve.H:29
ITHACAparameters * para
Definition NLsolve.H:40
std::ofstream res_os
Definition NLsolve.H:52
label saver
Counter to check if the middleStep has been reached or not (for turbulent case only)
Definition NLsolve.H:45
auto nut
Definition NLsolve.H:195
Vector< double > uresidual_v(0, 0, 0)
double durationOff
Definition NLsolve.H:49
std::ofstream res_E
Definition NLsolve.H:57
scalar presidual
Definition NLsolve.H:38
label folderN
Definition NLsolve.H:46
scalar residual
Definition NLsolve.H:35
pressureControl & pressureControl
Definition NLsolve.H:30
TEqn solve()
Class for the definition of some general parameters, the parameters must be defined from the file ITH...
static ITHACAparameters * getInstance()
Gets an instance of ITHACAparameters, to be used if the instance is already existing.
void exportSolution(GeometricField< Type, PatchField, GeoMesh > &s, fileName subfolder, fileName folder, word fieldName)
Export a field to file in a certain folder and subfolder.
constrainPressure(p, rho, U, phiHbyA, rhorAtU, MRF)
bool pLimited
Definition pcEqn.H:98
HbyA
Definition pcEqn.H:62
simpleControl simple(mesh)
surfaceScalarField & phi
volVectorField & U
_fvOptions
_psi
volScalarField & p
volScalarField & rho
fluidThermo & thermo
_rho
_pressureControl
_initialMass
_phi
turbulence
adjustPhi(phiHbyA, U, p)