Loading...
Searching...
No Matches
NLmsrProblem.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\*---------------------------------------------------------------------------*/
25
26scalar residual = 1;
27scalar uresidual = 1;
28Vector<double> uresidual_v(0, 0, 0);
29scalar presidual = 1;
30scalar fluxresidual = 1;
31scalar prec1residual = 1;
32scalar prec2residual = 1;
33scalar prec3residual = 1;
34scalar prec4residual = 1;
35scalar prec5residual = 1;
36scalar prec6residual = 1;
37scalar prec7residual = 1;
38scalar prec8residual = 1;
39scalar Tresidual = 1;
40scalar dec1residual = 1;
41scalar dec2residual = 1;
42scalar dec3residual = 1;
43scalar csolve = 0;
44
45
46// Variable that can be changed
47turbulence->read();
48std::ofstream res_os;
49res_os.open("./ITHACAoutput/Offline/residuals", std::ios_base::app);
50std::string equation;
51
52#if defined(OFVER) && (OFVER == 6)
53
54while (simple.loop(runTime) && residual > tolerance && csolve < maxIter )
55#else
56while (simple.loop() && residual > tolerance && csolve < maxIter )
57#endif
58{
59 Info << "Time = " << runTime.timeName() << nl << endl;
60 {
61#include "UEqn.H"
62#include "TEqn.H"
63 TEqn.relax();
64 Tresidual = TEqn.solve().initialResidual();
65#include "pEqn.H"
66#include "updateConsts.H"
67#include "DiffEqn.H"
68#include "precEqns.H"
69#include "decEqns.H"
70#include "updateK.H"
71 flux_old = flux;
72 TEqn.relax();
73 Tresidual = TEqn.solve().initialResidual();
74 scalar C = 0;
75
76 for (label i = 0; i < 3; i++)
77 {
78 if (C < uresidual_v[i])
79 {
80 C = uresidual_v[i];
81 }
82 }
83
84 uresidual = C;
99 Info << "\nResidual: " << residual << endl << endl;
100 }
101 laminarTransport.correct();
102 turbulence->correct();
103 powerDens = (1 - decbetaTot) * flux * SP1_0 + (decLam1 * dec1 + decLam2 * dec2 +
104 decLam3 * dec3);
105 csolve = csolve + 1;
106 Info << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
107 << " ClockTime = " << runTime.elapsedClockTime() << " s"
108 << nl << endl;
109 Info << "\nIteration number: " << csolve << endl << endl;
110}
111
112res_os << residual << std::endl;
113res_os.close();
114runTime.setTime(runTime.startTime(), 0);
115ITHACAstream::exportSolution(U, name(counter), "./ITHACAoutput/Offline/");
116ITHACAstream::exportSolution(p, name(counter), "./ITHACAoutput/Offline/");
117ITHACAstream::exportSolution(flux, name(counter), "./ITHACAoutput/Offline/");
118ITHACAstream::exportSolution(prec1, name(counter), "./ITHACAoutput/Offline/");
119ITHACAstream::exportSolution(prec2, name(counter), "./ITHACAoutput/Offline/");
120ITHACAstream::exportSolution(prec3, name(counter), "./ITHACAoutput/Offline/");
121ITHACAstream::exportSolution(prec4, name(counter), "./ITHACAoutput/Offline/");
122ITHACAstream::exportSolution(prec5, name(counter), "./ITHACAoutput/Offline/");
123ITHACAstream::exportSolution(prec6, name(counter), "./ITHACAoutput/Offline/");
124ITHACAstream::exportSolution(prec7, name(counter), "./ITHACAoutput/Offline/");
125ITHACAstream::exportSolution(prec8, name(counter), "./ITHACAoutput/Offline/");
126ITHACAstream::exportSolution(T, name(counter), "./ITHACAoutput/Offline/");
127ITHACAstream::exportSolution(dec1, name(counter), "./ITHACAoutput/Offline/");
128ITHACAstream::exportSolution(dec2, name(counter), "./ITHACAoutput/Offline/");
129ITHACAstream::exportSolution(dec3, name(counter), "./ITHACAoutput/Offline/");
131 "./ITHACAoutput/Offline/");
132ITHACAstream::exportSolution(v, name(counter), "./ITHACAoutput/Offline/");
133ITHACAstream::exportSolution(D, name(counter), "./ITHACAoutput/Offline/");
134ITHACAstream::exportSolution(NSF, name(counter), "./ITHACAoutput/Offline/");
135ITHACAstream::exportSolution(A, name(counter), "./ITHACAoutput/Offline/");
136ITHACAstream::exportSolution(SP, name(counter), "./ITHACAoutput/Offline/");
137Ufield.append(U.clone());
138Pfield.append(p.clone());
139Fluxfield.append(flux.clone());
140Prec1field.append(prec1.clone());
141Prec2field.append(prec2.clone());
142Prec3field.append(prec3.clone());
143Prec4field.append(prec4.clone());
144Prec5field.append(prec5.clone());
145Prec6field.append(prec6.clone());
146Prec7field.append(prec7.clone());
147Prec8field.append(prec8.clone());
148Tfield.append(T.clone());
149Dec1field.append(dec1.clone());
150Dec2field.append(dec2.clone());
151Dec3field.append(dec3.clone());
152vFields.append(v.clone());
153DFields.append(D.clone());
154NSFFields.append(NSF.clone());
155AFields.append(A.clone());
156SPFields.append(SP.clone());
157
158
159
160
161
162
163
164
Foam::Time & runTime
Definition createTime.H:33
scalar prec2residual
scalar prec3residual
scalar fluxresidual
scalar dec1residual
scalar prec1residual
scalar prec4residual
scalar Tresidual
std::string equation
scalar prec7residual
scalar csolve
std::ofstream res_os
scalar prec6residual
scalar dec2residual
scalar dec3residual
scalar prec8residual
Vector< double > uresidual_v(0, 0, 0)
scalar uresidual
scalar prec5residual
scalar presidual
scalar residual
volScalarField & T
Definition createT.H:46
volScalarField & powerDens
fvScalarMatrix & TEqn
Definition TEqn.H:15
volScalarField & SP
dimensionedScalar & decbetaTot
volScalarField & NSF
volScalarField & A
volScalarField & D
volScalarField & v
volScalarField & prec5
volScalarField & prec6
volScalarField & prec8
volScalarField & flux
volScalarField & prec3
volScalarField & prec2
volScalarField & prec7
volScalarField & prec4
volScalarField & prec1
volScalarField & dec3
volScalarField & dec1
volScalarField & dec2
dimensionedScalar & SP1_0
dimensionedScalar & decLam1
dimensionedScalar & decLam3
dimensionedScalar & decLam2
void exportSolution(GeometricField< Type, PatchField, GeoMesh > &s, fileName subfolder, fileName folder, word fieldName)
Export a field to file in a certain folder and subfolder.
simpleControl simple(mesh)
volVectorField & U
volScalarField & p
turbulence
singlePhaseTransportModel & laminarTransport
label i
Definition pEqn.H:46