boundaryFoam.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-2017 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  boundaryFoam
28 
29 Group
30  grpIncompressibleSolvers
31 
32 Description
33  Steady-state solver for incompressible, 1D turbulent flow,
34  typically to generate boundary layer conditions at an inlet.
35 
36  Boundary layer code to calculate the U, k and epsilon distributions.
37  Used to create inlet boundary conditions for experimental comparisons
38  for which U and k have not been measured.
39  Turbulence model is runtime selectable.
40 
41 \*---------------------------------------------------------------------------*/
42 
43 #include "fvCFD.H"
46 #include "fvOptions.H"
47 #include "wallFvPatch.H"
48 #include "makeGraph.H"
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 int main(int argc, char *argv[])
53 {
54  argList::addNote
55  (
56  "Steady-state solver for incompressible, 1D turbulent flow,"
57  " typically to generate boundary layer conditions at an inlet."
58  );
59 
60  argList::noParallel();
61  #include "addCheckCaseOptions.H"
62 
63  #include "setRootCaseLists.H"
64  #include "createTime.H"
65  #include "createMesh.H"
66  #include "createFields.H"
67  #include "interrogateWallPatches.H"
68 
69  turbulence->validate();
70 
71  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
72 
73  Info<< "\nStarting time loop\n" << endl;
74 
75  while (runTime.loop())
76  {
77  Info<< "Time = " << runTime.timeName() << nl << endl;
78 
79  fvVectorMatrix divR(turbulence->divDevReff(U));
80  divR.source() = flowMask & divR.source();
81 
83  (
84  divR == gradP + fvOptions(U)
85  );
86 
87  UEqn.relax();
88 
89  fvOptions.constrain(UEqn);
90 
91  UEqn.solve();
92 
93  fvOptions.correct(U);
94 
95 
96  // Correct driving force for a constant volume flow rate
97  dimensionedVector UbarStar = flowMask & U.weightedAverage(mesh.V());
98 
99  U += (Ubar - UbarStar);
100  gradP += (Ubar - UbarStar)/(1.0/UEqn.A())().weightedAverage(mesh.V());
101 
102  laminarTransport.correct();
103  turbulence->correct();
104 
105  Info<< "Uncorrected Ubar = " << (flowDirection & UbarStar.value())
106  << ", pressure gradient = " << (flowDirection & gradP.value())
107  << endl;
108 
109  #include "evaluateNearWall.H"
110 
111  if (runTime.writeTime())
112  {
113  #include "makeGraphs.H"
114  }
115 
116  runTime.printExecutionTime(Info);
117  }
118 
119  Info<< "End\n" << endl;
120 
121  return 0;
122 }
123 
124 
125 // ************************************************************************* //
runTime
engineTime & runTime
Definition: createEngineTime.H:13
gradP
volVectorField gradP(fvc::grad(p))
fvOptions.H
turbulence
Info<< "Reading field U\n"<< 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("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\n"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
singlePhaseTransportModel.H
turbulentTransportModel.H
wallFvPatch.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:50
setRootCaseLists.H
addCheckCaseOptions.H
Required Classes.
Foam::fvVectorMatrix
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:47
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:23
flowMask
tensor flowMask
Definition: createFields.H:42
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
makeGraph.H
flowDirection
vector flowDirection
Definition: createFields.H:41
evaluateNearWall.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
U
U
Definition: pEqn.H:72
Foam::nl
constexpr char nl
Definition: Ostream.H:404
makeGraphs.H
createMesh.H
Required Variables.
Ubar
dimensionedVector Ubar("Ubar", dimVelocity, laminarTransport)
createTime.H
UEqn
fvVectorMatrix & UEqn
Definition: UEqn.H:13
fvCFD.H
interrogateWallPatches.H
laminarTransport
singlePhaseTransportModel laminarTransport(U, phi)