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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Application
27  dnsFoam
28 
29 Group
30  grpDNSSolvers
31 
32 Description
33  Direct numerical simulation solver for boxes of isotropic turbulence.
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #include "fvCFD.H"
38 #include "Kmesh.H"
39 #include "UOprocess.H"
40 #include "fft.H"
41 #include "calcEk.H"
42 #include "graph.H"
43 #include "pisoControl.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 int main(int argc, char *argv[])
48 {
49  argList::addNote
50  (
51  "Direct numerical simulation for boxes of isotropic turbulence."
52  );
53 
54  #include "postProcess.H"
55 
56  #include "setRootCaseLists.H"
57  #include "createTime.H"
58  #include "createMeshNoClear.H"
59  #include "createControl.H"
60  #include "createFields.H"
61  #include "initContinuityErrs.H"
62 
63  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
64 
65  label ntot = 1;
66  forAll(K.nn(), idim)
67  {
68  ntot *= K.nn()[idim];
69  }
70  const scalar recRootN = 1.0/Foam::sqrt(scalar(ntot));
71 
72  Info<< nl << "Starting time loop" << endl;
73 
74  while (runTime.loop())
75  {
76  Info<< "Time = " << runTime.timeName() << nl << endl;
77 
78  force.primitiveFieldRef() = ReImSum
79  (
80  fft::reverseTransform
81  (
82  K/(mag(K) + 1.0e-6) ^ forceGen.newField(), K.nn()
83  )*recRootN
84  );
85 
86  #include "globalProperties.H"
87 
89  (
90  fvm::ddt(U)
91  + fvm::div(phi, U)
92  - fvm::laplacian(nu, U)
93  ==
94  force
95  );
96 
97  solve(UEqn == -fvc::grad(p));
98 
99 
100  // --- PISO loop
101  while (piso.correct())
102  {
103  volScalarField rAU(1.0/UEqn.A());
107  (
108  "phiHbyA",
109  fvc::flux(HbyA)
110  + rAUf*fvc::ddtCorr(U, phi)
111  );
112 
113  // Update the pressure BCs to ensure flux consistency
115 
116  fvScalarMatrix pEqn
117  (
119  );
120 
121  pEqn.solve(mesh.solver(p.select(piso.finalInnerIter())));
122 
123  phi = phiHbyA - pEqn.flux();
124 
125  #include "continuityErrs.H"
126 
127  U = HbyA - rAU*fvc::grad(p);
128  U.correctBoundaryConditions();
129  }
130 
131  runTime.write();
132 
133  if (runTime.writeTime())
134  {
135  calcEk(U, K).write
136  (
137  runTime.path()/"graphs"/runTime.timeName(),
138  "Ek",
139  runTime.graphFormat()
140  );
141  }
142 
143  runTime.printExecutionTime(Info);
144  }
145 
146  Info<< "End\n" << endl;
147 
148  return 0;
149 }
150 
151 
152 // ************************************************************************* //
Foam::constrainPressure
void constrainPressure(volScalarField &p, const RhoType &rho, const volVectorField &U, const surfaceScalarField &phiHbyA, const RAUType &rhorAU, const MRFType &MRF)
Definition: constrainPressure.C:39
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::constrainHbyA
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:35
Kmesh.H
pisoControl.H
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::fvc::flux
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
rAU
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
globalProperties.H
createMeshNoClear.H
Required Variables.
Foam::fac::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh >> grad(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facGrad.C:56
Foam::calcEk
graph calcEk(const volVectorField &U, const Kmesh &K)
Definition: calcEk.C:27
HbyA
HbyA
Definition: pcEqn.H:74
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::fac::div
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:50
solve
CEqn solve()
UOprocess.H
calcEk.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
K
CGAL::Exact_predicates_exact_constructions_kernel K
Definition: CGALTriangulation3DKernel.H:58
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
piso
pisoControl piso(mesh)
setRootCaseLists.H
Foam::fvVectorMatrix
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:47
phiHbyA
phiHbyA
Definition: pcEqn.H:73
Foam::fvScalarMatrix
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:44
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
ddtCorr
ddtCorr
Definition: readControls.H:9
Foam::graph::write
void write(Ostream &, const word &format) const
Write graph to stream in given format.
Definition: graph.C:267
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
postProcess.H
Execute application functionObjects to post-process existing results.
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:62
Foam::ReImSum
scalarField ReImSum(const UList< complex > &cf)
Sum real and imag components.
Definition: complexField.C:146
U
U
Definition: pEqn.H:72
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::fac::ddt
tmp< GeometricField< Type, faPatchField, areaMesh > > ddt(const dimensioned< Type > dt, const faMesh &mesh)
Definition: facDdt.C:47
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
createTime.H
graph.H
UEqn
fvVectorMatrix & UEqn
Definition: UEqn.H:13
fvCFD.H
Foam::fac::laplacian
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: facLaplacian.C:47
fft.H
Foam::fac::interpolate
static tmp< GeometricField< Type, faePatchField, edgeMesh > > interpolate(const GeometricField< Type, faPatchField, areaMesh > &tvf, const edgeScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
rAUf
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))