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 | Copyright (C) 2011-2016 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 "setRootCase.H"
52 
53  #include "createTime.H"
54  #include "createMesh.H"
55  #include "createFields.H"
56  #include "createFvOptions.H"
57  #include "interrogateWallPatches.H"
58 
59  turbulence->validate();
60 
61  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
62 
63  Info<< "\nStarting time loop\n" << endl;
64 
65  while (runTime.loop())
66  {
67  Info<< "Time = " << runTime.timeName() << nl << endl;
68 
69  fvVectorMatrix divR(turbulence->divDevReff(U));
70  divR.source() = flowMask & divR.source();
71 
73  (
74  divR == gradP + fvOptions(U)
75  );
76 
77  UEqn.relax();
78 
79  fvOptions.constrain(UEqn);
80 
81  UEqn.solve();
82 
83  fvOptions.correct(U);
84 
85 
86  // Correct driving force for a constant volume flow rate
87  dimensionedVector UbarStar = flowMask & U.weightedAverage(mesh.V());
88 
89  U += (Ubar - UbarStar);
90  gradP += (Ubar - UbarStar)/(1.0/UEqn.A())().weightedAverage(mesh.V());
91 
92  laminarTransport.correct();
93  turbulence->correct();
94 
95  Info<< "Uncorrected Ubar = " << (flowDirection & UbarStar.value())
96  << ", pressure gradient = " << (flowDirection & gradP.value())
97  << endl;
98 
99  #include "evaluateNearWall.H"
100 
101  if (runTime.writeTime())
102  {
103  #include "makeGraphs.H"
104  }
105 
106  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
107  << " ClockTime = " << runTime.elapsedClockTime() << " s"
108  << nl << endl;
109  }
110 
111  Info<< "End\n" << endl;
112 
113  return 0;
114 }
115 
116 
117 // ************************************************************************* //
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
U
Definition: pEqn.H:83
vector flowDirection
Definition: createFields.H:41
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
dimensionedVector gradP("gradP", dimensionSet(0, 1,-2, 0, 0), Zero)
fv::options & fvOptions
singlePhaseTransportModel laminarTransport(U, phi)
dynamicFvMesh & mesh
tensor flowMask
Definition: createFields.H:42
static const char nl
Definition: Ostream.H:262
dimensionedVector Ubar("Ubar", dimVelocity, laminarTransport)
fvVectorMatrix & UEqn
Definition: UEqn.H:13
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
messageStream Info