dnsFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 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 "argList.H"
33 #include "timeSelector.H"
34 #include "Kmesh.H"
35 #include "UOprocess.H"
36 #include "fft.H"
37 #include "writeEk.H"
38 #include "writeFile.H"
39 
40 #include "pisoControl.H"
41 #include "constrainPressure.H"
42 #include "constrainHbyA.H"
43 
44 #include "fvcDdt.H"
45 #include "fvcGrad.H"
46 #include "fvcFlux.H"
47 
48 #include "fvmDdt.H"
49 #include "fvmDiv.H"
50 #include "fvmLaplacian.H"
51 
52 using namespace Foam;
53 
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 
56 int main(int argc, char *argv[])
57 {
58  #include "postProcess.H"
59 
60  #include "setRootCase.H"
61  #include "createTime.H"
62  #include "createMeshNoClear.H"
63  #include "createControl.H"
64  #include "createFields.H"
65  #include "initContinuityErrs.H"
66 
67  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
68 
69  Info<< nl << "Starting time loop" << endl;
70 
71  while (runTime.loop())
72  {
73  Info<< "Time = " << runTime.userTimeName() << nl << endl;
74 
75  force.primitiveFieldRef() = ReImSum
76  (
78  (
79  K/(mag(K) + 1.0e-6) ^ forceGen.newField(), K.nn()
80  )
81  );
82 
83  #include "globalProperties.H"
84 
86  (
87  fvm::ddt(U)
88  + fvm::div(phi, U)
89  - fvm::laplacian(nu, U)
90  ==
91  force
92  );
93 
94  solve(UEqn == -fvc::grad(p));
95 
96 
97  // --- PISO loop
98  while (piso.correct())
99  {
100  volScalarField rAU(1.0/UEqn.A());
101  surfaceScalarField rAUf("rAUf", fvc::interpolate(rAU));
104  (
105  "phiHbyA",
106  fvc::flux(HbyA)
107  + rAUf*fvc::ddtCorr(U, phi)
108  );
109 
110  // Update the pressure BCs to ensure flux consistency
111  constrainPressure(p, U, phiHbyA, rAUf);
112 
113  fvScalarMatrix pEqn
114  (
115  fvm::laplacian(rAUf, p) == fvc::div(phiHbyA)
116  );
117 
118  pEqn.solve();
119 
120  phi = phiHbyA - pEqn.flux();
121 
122  #include "continuityErrs.H"
123 
124  U = HbyA - rAU*fvc::grad(p);
125  U.correctBoundaryConditions();
126  }
127 
128  runTime.write();
129 
130  if (runTime.writeTime())
131  {
132  writeEk(U, K);
133  }
134 
135  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
136  << " ClockTime = " << runTime.elapsedClockTime() << " s"
137  << nl << endl;
138  }
139 
140  Info<< "End\n" << endl;
141 
142  return 0;
143 }
144 
145 
146 // ************************************************************************* //
Generic GeometricField class.
static tmp< complexField > reverseTransform(const tmp< complexField > &field, const labelList &nn)
Definition: fft.C:203
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
SolverPerformance< Type > solve(const dictionary &)
Solve segregated or coupled returning the solution statistics.
Definition: fvMatrixSolve.C:58
tmp< volScalarField > A() const
Return the central coefficient.
Definition: fvMatrix.C:809
tmp< VolField< Type > > H() const
Return the H operation source.
Definition: fvMatrix.C:867
tmp< SurfaceField< Type > > flux() const
Return the face-flux field from the matrix.
Definition: fvMatrix.C:963
bool correct(const bool finalIter=false)
Piso loop within outer loop.
Definition: pisoControl.C:77
fvVectorMatrix & UEqn
Definition: UEqn.H:11
Calculates and prints the continuity errors.
pisoControl piso(mesh)
int main(int argc, char *argv[])
Definition: dnsFoam.C:53
Calculate the first temporal derivative.
Calculate the face-flux of the given field.
Calculate the gradient of the given field.
Calculate the matrix for the first temporal derivative.
Calculate the matrix for the divergence of the given field and flux.
Calculate the matrix for the laplacian of the field.
Declare and initialise the cumulative continuity error.
K
Definition: pEqn.H:75
volScalarField rAU(1.0/UEqn.A())
U
Definition: pEqn.H:72
volVectorField & HbyA
Definition: pEqn.H:13
surfaceScalarField phiHbyA("phiHbyA", fvc::interpolate(rho) *fvc::flux(HbyA))
tmp< SurfaceField< typename innerProduct< vector, Type >::type > > flux(const VolField< Type > &vf)
Return the face-flux field obtained from the given volVectorField.
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< SurfaceField< typename Foam::flux< Type >::type > > ddtCorr(const VolField< Type > &U, const SurfaceField< Type > &Uf)
Definition: fvcDdt.C:196
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
Namespace for OpenFOAM.
void constrainPressure(volScalarField &p, const RhoType &rho, const volVectorField &U, const surfaceScalarField &phiHbyA, const RAUType &rhorAU, const MRFType &MRF)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:34
void writeEk(const volVectorField &U, const Kmesh &K)
Definition: writeEk.C:41
static const char nl
Definition: Ostream.H:260
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
scalarField ReImSum(const UList< complex > &cf)
Definition: complexFields.C:84
Execute application functionObjects to post-process existing results.
volScalarField & p
This function evaluates q(k) (McComb, p61) by summing all wavevectors in a k-shell....