All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
scalarTransportQHDFoam.C
Go to the documentation of this file.
1 
2 
3 /*---------------------------------------------------------------------------*\
4  ========= |
5  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
6  \\ / O peration |
7  \\ / A nd | www.openfoam.com
8  \\/ M anipulation |
9 -------------------------------------------------------------------------------
10  Copyright (C) 2011-2016 OpenFOAM Foundation
11  Copyright (C) 2019 OpenCFD Ltd.
12  Copyright (C) 2016-2019 ISP RAS (www.ispras.ru) UniCFD Group (www.unicfd.ru)
13 -------------------------------------------------------------------------------
14 License
15  This file is part of QGDsolver library, based on OpenFOAM+.
16 
17  OpenFOAM is free software: you can redistribute it and/or modify it
18  under the terms of the GNU General Public License as published by
19  the Free Software Foundation, either version 3 of the License, or
20  (at your option) any later version.
21 
22  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
23  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
24  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
25  for more details.
26 
27  You should have received a copy of the GNU General Public License
28  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
29 
30 Application
31  scalarTransportQHDFoam
32 
33 Description
34  Evolves a passive scalar transport equation with regularized terms.
35 
36  Developed by UniCFD group (www.unicfd.ru) of ISP RAS (www.ispras.ru).
37 
38 
39 \*---------------------------------------------------------------------------*/
40 
41 #include "fvCFD.H"
42 #include "QHD.H"
43 
44 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45 
46 int main(int argc, char *argv[])
47 {
48  #define NO_CONTROL
49  #include "postProcess.H"
50 
51  #include "setRootCase.H"
52  #include "createTime.H"
53  #include "createMesh.H"
54  #include "createFields.H"
55  #include "createFaceFields.H"
56  #include "createFaceFluxes.H"
57  #include "createTimeControls.H"
58 
59  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60 
61  // Courant numbers used to adjust the time-step
62  scalar CoNum = 0.0;
63  scalar meanCoNum = 0.0;
64 
65  //solve for pressure to get zero divergence flux field
66  #include "updateFluxes.H"
67 
68  Info<< "\nStarting time loop\n" << endl;
69 
70  while (runTime.run())
71  {
72  /*
73  *
74  * Update fields
75  *
76  */
77  #include "updateFields.H"
78 
79 
80  /*
81  *
82  * Update time step
83  *
84  */
85  #include "readTimeControls.H"
86  if(runTime.controlDict().lookupOrDefault<bool>("adjustTimeStep", false))
87  {
88  surfaceScalarField Cof
89  (
90  "Cof",
91  runTime.deltaT()*
92  (
93  mag(Uf)/hQGDf
94  )
95  );
96  CoNum = max(Cof).value();
97  Info << "Courant Number = " << CoNum << endl;
98  }
99 
100  #include "setDeltaT-QGDQHD.H"
101 
102  runTime++;
103 
104  Info<< "Time = " << runTime.timeName() << nl << endl;
105 
106  // --- Store old time values
107  T.oldTime();
108 
109  {
110  phiTf = qgdFlux(phiu,T,Tf);
111  surfaceScalarField phiTauTReg = tauQGDf*phiu*(Uf & gradTf);
112 
113  // --- Solve T
114  if (implicitDiffusion)
115  {
116  solve
117  (
118  fvm::ddt(T)
119  + fvc::div(phiTf) - fvc::Sp(fvc::div(phiu),T)
120  - fvm::laplacian(Hif,T)
121  - fvc::div(phiTauTReg)
122  ==
123  TSu
124  );
125  }
126 
127  Info << "max/min of T: " << max(T).value() << "/" << min(T).value() << endl;
128  }
129 
130  runTime.write();
131 
132  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
133  << " ClockTime = " << runTime.elapsedClockTime() << " s"
134  << nl << endl;
135  }
136 
137 
138 
139  Info<< "End\n" << endl;
140 
141  return 0;
142 }
143 
144 // ************************************************************************* //
tmp< GeometricField< T, Foam::fvsPatchField, Foam::surfaceMesh > > qgdFlux(const GeometricField< scalar, Foam::fvsPatchField, Foam::surfaceMesh > &flux, const GeometricField< T, Foam::fvPatchField, Foam::volMesh > &psi, const GeometricField< T, Foam::fvsPatchField, Foam::surfaceMesh > &psif, const word fluxName)
const surfaceScalarField & tauQGDf
Definition: createFields.H:55
const surfaceScalarField & hQGDf
Definition: createFields.H:54
surfaceScalarField phiTauTReg
Definition: QHDTEqn.H:57
int main(int argc, char *argv[])
Switch implicitDiffusion(thermo.implicitDiffusion())
EEqn solve()
surfaceScalarField phiTf("phiTf", phi *Tf)
gradTf
Definition: updateFields.H:10
Reset the timestep to maintain a constant maximum courant Number. Reduction of time-step is immediate...
Tf
Definition: updateFields.H:33
tmp< surfaceScalarField > div(const volVectorField &vF)
Hif
Definition: updateFields.H:39
Uf
Definition: updateFields.H:30
phiu
——–Start———
Definition: updateFluxes.H:33
volScalarField & T
Definition: createFields.H:53
Includation for QHD solver. Equations of state for QHD based on density. Class for fvsc namespace...
fvScalarMatrix TSu(T, T.dimensions()*dimVolume/dimTime)