boundaryFoam.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-2018 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  boundaryFoam
26 
27 Description
28  Steady-state solver for incompressible, 1D turbulent flow, typically to
29  generate boundary layer conditions at an inlet, for use in a simulation.
30 
31  Boundary layer code to calculate the U, k and epsilon distributions.
32  Used to create inlet boundary conditions for experimental comparisons
33  for which U and k have not been measured.
34  Turbulence model is runtime selectable.
35 
36 \*---------------------------------------------------------------------------*/
37 
38 #include "fvCFD.H"
41 #include "fvOptions.H"
42 #include "wallFvPatch.H"
43 #include "makeGraph.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 int main(int argc, char *argv[])
48 {
49  argList::noParallel();
50 
51  #include "setRootCaseLists.H"
52 
53  #include "createTime.H"
54  #include "createMesh.H"
55  #include "createFields.H"
56  #include "interrogateWallPatches.H"
57 
58  turbulence->validate();
59 
60  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61 
62  Info<< "\nStarting time loop\n" << endl;
63 
64  while (runTime.loop())
65  {
66  Info<< "Time = " << runTime.timeName() << nl << endl;
67 
68  fvVectorMatrix divR(turbulence->divDevReff(U));
69  divR.source() = flowMask & divR.source();
70 
72  (
73  divR == gradP + fvOptions(U)
74  );
75 
76  UEqn.relax();
77 
78  fvOptions.constrain(UEqn);
79 
80  UEqn.solve();
81 
82  fvOptions.correct(U);
83 
84 
85  // Correct driving force for a constant volume flow rate
86  dimensionedVector UbarStar = flowMask & U.weightedAverage(mesh.V());
87 
88  U += (Ubar - UbarStar);
89  gradP += (Ubar - UbarStar)/(1.0/UEqn.A())().weightedAverage(mesh.V());
90 
91  laminarTransport.correct();
92  turbulence->correct();
93 
94  Info<< "Uncorrected Ubar = " << (flowDirection & UbarStar.value())
95  << ", pressure gradient = " << (flowDirection & gradP.value())
96  << endl;
97 
98  #include "evaluateNearWall.H"
99 
100  if (runTime.writeTime())
101  {
102  #include "makeGraphs.H"
103  }
104 
105  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
106  << " ClockTime = " << runTime.elapsedClockTime() << " s"
107  << nl << endl;
108  }
109 
110  Info<< "End\n" << endl;
111 
112  return 0;
113 }
114 
115 
116 // ************************************************************************* //
fv::options & fvOptions
vector flowDirection
Definition: createFields.H:41
dimensionedVector gradP("gradP", dimensionSet(0, 1, -2, 0, 0), Zero)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
dynamicFvMesh & mesh
tensor flowMask
Definition: createFields.H:42
static const char nl
Definition: Ostream.H:265
U
Definition: pEqn.H:72
dimensionedVector Ubar("Ubar", dimVelocity, laminarTransport)
fvVectorMatrix & UEqn
Definition: UEqn.H:13
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
messageStream Info
singlePhaseTransportModel laminarTransport(U, phi)