Loading...
Searching...
No Matches
steadySolver.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
26while (simple.loop(runTime))
27{
28 Info << "Time = " << runTime.timeName() << nl << endl;
29 // --- Pressure-velocity SIMPLE corrector
30 {
31#include "UEqn.H"
32#include "TEqn.H"
33#include "pEqn.H"
34 rho = rhoRef * (1 - betaTE * (T - Tref));
35 logT = log(max(0.5, (T / TrefXS)));
36 NSF1 = (1 / Keff) * (NSF1_0 + alfa_NSF1 * logT);
37 M11 = NSF1 * (1 - betaTot) - (A1_0 + alfa_A1 * logT);
38 D1 = (D1_0 + alfa_D1 * logT);
39 SP1 = (1 / Keff) * (SP1_0 + alfa_SP1 * logT);
40 alphaEff = turbulence->nu() / Pr + turbulence->nut() / Prt;
41 diffEff = turbulence->nu() / Sc + turbulence->nut() / Sct;
42#include "DiffEqn.H"
43#include "precEqns.H"
44#include "decEqns.H"
45#include "TEqn.H"
46 }
47 laminarTransport.correct();
48 turbulence->correct();
49 runTime.write();
50 Info << "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
51 << " ClockTime = " << runTime.elapsedClockTime() << " s"
52 << nl << endl;
53 powerDens = (1 - decbetaTot) * flux * SP1 + (decLam1 * dec1 + decLam2 * dec2 +
54 decLam3 * dec3);
55 exportSolution(U, name(counter), "./ITHACAoutput/Offline/");
56 exportSolution(p, name(counter), "./ITHACAoutput/Offline/");
57 exportSolution(flux, name(counter), "./ITHACAoutput/Offline/");
58 exportSolution(prec1, name(counter), "./ITHACAoutput/Offline/");
59 exportSolution(prec2, name(counter), "./ITHACAoutput/Offline/");
60 exportSolution(prec3, name(counter), "./ITHACAoutput/Offline/");
61 exportSolution(prec4, name(counter), "./ITHACAoutput/Offline/");
62 exportSolution(prec5, name(counter), "./ITHACAoutput/Offline/");
63 exportSolution(prec6, name(counter), "./ITHACAoutput/Offline/");
64 exportSolution(prec7, name(counter), "./ITHACAoutput/Offline/");
65 exportSolution(prec8, name(counter), "./ITHACAoutput/Offline/");
66 exportSolution(T, name(counter), "./ITHACAoutput/Offline/");
67 exportSolution(dec1, name(counter), "./ITHACAoutput/Offline/");
68 exportSolution(dec2, name(counter), "./ITHACAoutput/Offline/");
69 exportSolution(dec3, name(counter), "./ITHACAoutput/Offline/");
70 exportSolution(powerDens, name(counter), "./ITHACAoutput/Offline/");
71 exportSolution(rho, name(counter), "./ITHACAoutput/Offline/");
72 exportSolution(NSF1, name(counter), "./ITHACAoutput/Offline/");
73 exportSolution(M11, name(counter), "./ITHACAoutput/Offline/");
74 exportSolution(D1, name(counter), "./ITHACAoutput/Offline/");
75 exportSolution(SP1, name(counter), "./ITHACAoutput/Offline/");
76 exportSolution(alphaEff, name(counter), "./ITHACAoutput/Offline/");
77 exportSolution(diffEff, name(counter), "./ITHACAoutput/Offline/");
78 Ufield.append(U);
79 Pfield.append(p);
80 Fluxfield.append(flux);
81 Prec1field.append(prec1);
82 Prec2field.append(prec2);
83 Prec3field.append(prec3);
84 Prec4field.append(prec4);
85 Prec5field.append(prec5);
86 Prec6field.append(prec6);
87 Prec7field.append(prec7);
88 Prec8field.append(prec8);
89 Tfield.append(T);
90 Dec1field.append(dec1);
91 Dec2field.append(dec2);
92 Dec3field.append(dec3);
93 rhofield.append(rho);
94 NSF1field.append(NSF1);
95 M11field.append(M11);
96 D1field.append(D1);
97 SP1field.append(SP1);
98 aEfield.append(alphaEff);
99 dEfield.append(diffEff);
100 counter++;
101}
102
103
104
Foam::Time & runTime
Definition createTime.H:33
volScalarField alphaEff("alphaEff", turbulence->nu()/Pr+alphat)
volScalarField & T
Definition createT.H:46
volScalarField & powerDens
volScalarField & logT
dimensionedScalar & decbetaTot
dimensionedScalar & betaTot
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 & Sc
dimensionedScalar & Tref
dimensionedScalar & rhoRef
dimensionedScalar & betaTE
dimensionedScalar & Prt
dimensionedScalar & Pr
dimensionedScalar & Sct
dimensionedScalar & D1_0
dimensionedScalar & A1_0
dimensionedScalar & Keff
dimensionedScalar & alfa_NSF1
dimensionedScalar & alfa_A1
dimensionedScalar & NSF1_0
dimensionedScalar & TrefXS
dimensionedScalar & alfa_D1
dimensionedScalar & alfa_SP1
dimensionedScalar & SP1_0
dimensionedScalar & decLam1
dimensionedScalar & decLam3
dimensionedScalar & decLam2
simpleControl simple(mesh)
volVectorField & U
volScalarField & p
turbulence
singlePhaseTransportModel & laminarTransport