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-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  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 "argList.H"
39 #include "viscosityModel.H"
41 #include "fvModels.H"
42 #include "fvConstraints.H"
43 #include "wallFvPatch.H"
44 #include "setWriter.H"
45 #include "writeFile.H"
46 
47 using namespace Foam;
48 
49 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51 int main(int argc, char *argv[])
52 {
54 
55  #include "setRootCase.H"
56 
57  #include "createTime.H"
58  #include "createMesh.H"
59  #include "createFields.H"
60  #include "interrogateWallPatches.H"
61 
62  turbulence->validate();
63 
64  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
65 
66  Info<< "\nStarting time loop\n" << endl;
67 
68  while (runTime.loop())
69  {
70  Info<< "Time = " << runTime.userTimeName() << nl << endl;
71 
72  fvModels.correct();
73 
74  turbulence->predict();
75 
76  fvVectorMatrix divR(turbulence->divDevSigma(U));
77  divR.source() = flowMask & divR.source();
78 
80  (
81  divR == gradP + fvModels.source(U)
82  );
83 
84  UEqn.relax();
85 
87 
88  UEqn.solve();
89 
91 
92 
93  // Correct driving force for a constant volume flow rate
94  const dimensionedVector UbarStar
95  (
96  flowMask & U.weightedAverage(mesh.V())
97  );
98 
99  U += (Ubar - UbarStar);
100  gradP += (Ubar - UbarStar)/(1.0/UEqn.A())().weightedAverage(mesh.V());
101 
102  viscosity->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  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
117  << " ClockTime = " << runTime.elapsedClockTime() << " s"
118  << nl << endl;
119  }
120 
121  Info<< "End\n" << endl;
122 
123  return 0;
124 }
125 
126 
127 // ************************************************************************* //
static void noParallel()
Remove the parallel options.
Definition: argList.C:175
const Type & value() const
Return const reference to value.
Finite volume constraints.
Definition: fvConstraints.H:62
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:604
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
Finite volume models.
Definition: fvModels.H:65
virtual void correct()
Correct the fvModels.
Definition: fvModels.C:358
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
fvVectorMatrix & UEqn
Definition: UEqn.H:11
int main(int argc, char *argv[])
Definition: financialFoam.C:44
U
Definition: pEqn.H:72
Namespace for OpenFOAM.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
static const char nl
Definition: Ostream.H:260
vector flowDirection
Definition: createFields.H:41
dimensionedVector gradP("gradP", dimensionSet(0, 1, -2, 0, 0), Zero)
tensor flowMask
Definition: createFields.H:42
autoPtr< incompressible::momentumTransportModel > turbulence(incompressible::momentumTransportModel::New(U, phi, viscosity))
dimensionedVector Ubar("Ubar", dimVelocity, viscosity)