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 Description
28  Direct numerical simulation solver for boxes of isotropic turbulence
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #include "fvCFD.H"
33 #include "Kmesh.H"
34 #include "UOprocess.H"
35 #include "fft.H"
36 #include "calcEk.H"
37 #include "graph.H"
38 #include "pisoControl.H"
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 int main(int argc, char *argv[])
43 {
44  #include "setRootCase.H"
45 
46  #include "createTime.H"
47  #include "createMeshNoClear.H"
48 
49  pisoControl piso(mesh);
50 
51  #include "readTransportProperties.H"
52  #include "createFields.H"
54  #include "initContinuityErrs.H"
55 
56  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57 
58  Info<< nl << "Starting time loop" << endl;
59 
60  while (runTime.loop())
61  {
62  Info<< "Time = " << runTime.timeName() << nl << endl;
63 
64  force.internalField() = ReImSum
65  (
66  fft::reverseTransform
67  (
68  K/(mag(K) + 1.0e-6) ^ forceGen.newField(), K.nn()
69  )
70  );
71 
72  #include "globalProperties.H"
73 
75  (
76  fvm::ddt(U)
77  + fvm::div(phi, U)
78  - fvm::laplacian(nu, U)
79  ==
80  force
81  );
82 
83  solve(UEqn == -fvc::grad(p));
84 
85 
86  // --- PISO loop
87  while (piso.correct())
88  {
89  volScalarField rAU(1.0/UEqn.A());
91  volVectorField HbyA("HbyA", U);
92  HbyA = rAU*UEqn.H();
93 
95  (
96  "phiHbyA",
97  (fvc::interpolate(HbyA) & mesh.Sf())
99  );
100 
101  fvScalarMatrix pEqn
102  (
104  );
105 
106  pEqn.solve(mesh.solver(p.select(piso.finalInnerIter())));
107 
108  phi = phiHbyA - pEqn.flux();
109 
110  #include "continuityErrs.H"
111 
112  U = HbyA - rAU*fvc::grad(p);
113  U.correctBoundaryConditions();
114  }
115 
116  runTime.write();
117 
118  if (runTime.outputTime())
119  {
120  calcEk(U, K).write
121  (
122  runTime.path()/"graphs"/runTime.timeName(),
123  "Ek",
124  runTime.graphFormat()
125  );
126  }
127 
128  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
129  << " ClockTime = " << runTime.elapsedClockTime() << " s"
130  << nl << endl;
131  }
132 
133  Info<< "End\n" << endl;
134 
135  return 0;
136 }
137 
138 
139 // ************************************************************************* //
void write(Ostream &, const word &format) const
Write graph to stream in given format.
Definition: graph.C:225
rhoEqn solve()
dimensioned< scalar > mag(const dimensioned< Type > &)
CGAL::Exact_predicates_exact_constructions_kernel K
surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU))
phiHbyA
Definition: pcEqn.H:74
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
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
messageStream Info
tmp< surfaceScalarField > interpolate(const RhoType &rho)
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
dynamicFvMesh & mesh
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
static const char nl
Definition: Ostream.H:260
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
scalarField ReImSum(const UList< complex > &cf)
Definition: complexFields.C:84
fvVectorMatrix UEqn(fvm::ddt(rho, U)+fvm::div(phi, U)+turbulence->divDevRhoReff(U)==fvOptions(rho, U))
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
volScalarField & p
Definition: createFields.H:51
const dimensionedScalar e
Elementary charge.
Definition: doubleFloat.H:78
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
graph calcEk(const volVectorField &U, const Kmesh &K)
Definition: calcEk.C:27
surfaceScalarField & phi
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
volScalarField rAU(1.0/UEqn.A())
volScalarField & nu
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
int main(int argc, char *argv[])
Definition: postCalc.C:54
U
Definition: pEqn.H:82
HbyA
Definition: pEqn.H:7