dnsFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  dnsFoam
26 
27 Group
28  grpDNSSolvers
29 
30 Description
31  Direct numerical simulation solver for boxes of isotropic turbulence
32 
33 \*---------------------------------------------------------------------------*/
34 
35 #include "fvCFD.H"
36 #include "Kmesh.H"
37 #include "UOprocess.H"
38 #include "fft.H"
39 #include "calcEk.H"
40 #include "graph.H"
41 #include "pisoControl.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 int main(int argc, char *argv[])
46 {
47  #include "setRootCase.H"
48 
49  #include "createTime.H"
50  #include "createMeshNoClear.H"
51 
52  pisoControl piso(mesh);
53 
54  #include "readTransportProperties.H"
55  #include "createFields.H"
57  #include "initContinuityErrs.H"
58 
59  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60 
61  Info<< nl << "Starting time loop" << endl;
62 
63  while (runTime.loop())
64  {
65  Info<< "Time = " << runTime.timeName() << nl << endl;
66 
67  force.internalField() = ReImSum
68  (
69  fft::reverseTransform
70  (
71  K/(mag(K) + 1.0e-6) ^ forceGen.newField(), K.nn()
72  )
73  );
74 
75  #include "globalProperties.H"
76 
78  (
79  fvm::ddt(U)
80  + fvm::div(phi, U)
81  - fvm::laplacian(nu, U)
82  ==
83  force
84  );
85 
86  solve(UEqn == -fvc::grad(p));
87 
88 
89  // --- PISO loop
90  while (piso.correct())
91  {
92  volScalarField rAU(1.0/UEqn.A());
94  volVectorField HbyA("HbyA", U);
95  HbyA = rAU*UEqn.H();
96 
98  (
99  "phiHbyA",
100  (fvc::interpolate(HbyA) & mesh.Sf())
101  + rAUf*fvc::ddtCorr(U, phi)
102  );
103 
104  fvScalarMatrix pEqn
105  (
107  );
108 
109  pEqn.solve(mesh.solver(p.select(piso.finalInnerIter())));
110 
111  phi = phiHbyA - pEqn.flux();
112 
113  #include "continuityErrs.H"
114 
115  U = HbyA - rAU*fvc::grad(p);
116  U.correctBoundaryConditions();
117  }
118 
119  runTime.write();
120 
121  if (runTime.outputTime())
122  {
123  calcEk(U, K).write
124  (
125  runTime.path()/"graphs"/runTime.timeName(),
126  "Ek",
127  runTime.graphFormat()
128  );
129  }
130 
131  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
132  << " ClockTime = " << runTime.elapsedClockTime() << " s"
133  << nl << endl;
134  }
135 
136  Info<< "End\n" << endl;
137 
138  return 0;
139 }
140 
141 
142 // ************************************************************************* //
Kmesh.H
pisoControl.H
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
p
p
Definition: pEqn.H:62
globalProperties.H
createMeshNoClear.H
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::MULES::interpolate
tmp< surfaceScalarField > interpolate(const RhoType &rho)
Definition: IMULESTemplates.C:40
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
UOprocess.H
calcEk.H
U
U
Definition: pEqn.H:46
K
CGAL::Exact_predicates_exact_constructions_kernel K
Definition: CGALTriangulation3DKernel.H:56
nu
volScalarField & nu
Definition: readMechanicalProperties.H:179
Foam::fvc::laplacian
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
solve
rhoEqn solve()
Foam::fvVectorMatrix
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
Foam::fvScalarMatrix
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::fvc::ddtCorr
tmp< GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMesh > > ddtCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
Definition: fvcDdt.C:155
UEqn
tmp< fvVectorMatrix > UEqn(fvm::div(phi, U)+MRF.DDt(U)+turbulence->divDevReff(U)==fvOptions(U))
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Foam::graph::write
void write(Ostream &, const word &format) const
Write graph to stream in given format.
Definition: graph.C:221
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::ReImSum
scalarField ReImSum(const UList< complex > &cf)
Definition: complexFields.C:84
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
rAU
volScalarField rAU("rAU", 1.0/UEqn().A())
setRootCase.H
HbyA
HbyA
Definition: pEqn.H:4
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:52
readTurbulenceProperties.H
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: doubleFloat.H:94
createTime.H
graph.H
Foam::calcEk
graph calcEk(const volVectorField &U, const Kmesh &K)
Definition: calcEk.C:27
fvCFD.H
phiHbyA
phiHbyA
Definition: pEqn.H:21
Foam::fvc::ddt
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
fft.H
rAUf
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))