wallShearStress.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-2015 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  wallShearStress
26 
27 Description
28  Calculates and reports wall shear stress for all patches, for the
29  specified times when using RAS turbulence models.
30 
31  Compressible modes is automatically selected based on the existence of the
32  "thermophysicalProperties" dictionary required to construct the
33  thermodynamics package.
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #include "fvCFD.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 void calcIncompressible
45 (
46  const fvMesh& mesh,
47  const Time& runTime,
48  const volVectorField& U,
49  volVectorField& wallShearStress
50 )
51 {
52  #include "createPhi.H"
53 
54  singlePhaseTransportModel laminarTransport(U, phi);
55 
56  autoPtr<incompressible::RASModel> model
57  (
58  incompressible::New<incompressible::RASModel>(U, phi, laminarTransport)
59  );
60 
61  const volSymmTensorField Reff(model->devReff());
62 
63  forAll(wallShearStress.boundaryField(), patchI)
64  {
65  wallShearStress.boundaryField()[patchI] =
66  (
67  -mesh.Sf().boundaryField()[patchI]
68  /mesh.magSf().boundaryField()[patchI]
69  ) & Reff.boundaryField()[patchI];
70  }
71 }
72 
73 
74 void calcCompressible
75 (
76  const fvMesh& mesh,
77  const Time& runTime,
78  const volVectorField& U,
79  volVectorField& wallShearStress
80 )
81 {
82  IOobject rhoHeader
83  (
84  "rho",
85  runTime.timeName(),
86  mesh,
87  IOobject::MUST_READ,
88  IOobject::NO_WRITE
89  );
90 
91  if (!rhoHeader.headerOk())
92  {
93  Info<< " no rho field" << endl;
94  return;
95  }
96 
97  Info<< "Reading field rho\n" << endl;
98  volScalarField rho(rhoHeader, mesh);
99 
100  #include "compressibleCreatePhi.H"
101 
102  autoPtr<fluidThermo> pThermo(fluidThermo::New(mesh));
103  fluidThermo& thermo = pThermo();
104 
105  autoPtr<compressible::RASModel> model
106  (
107  compressible::New<compressible::RASModel>
108  (
109  rho,
110  U,
111  phi,
112  thermo
113  )
114  );
115 
116  const volSymmTensorField Reff(model->devRhoReff());
117 
118  forAll(wallShearStress.boundaryField(), patchI)
119  {
120  wallShearStress.boundaryField()[patchI] =
121  (
122  -mesh.Sf().boundaryField()[patchI]
123  /mesh.magSf().boundaryField()[patchI]
124  ) & Reff.boundaryField()[patchI];
125  }
126 }
127 
128 
129 int main(int argc, char *argv[])
130 {
131  timeSelector::addOptions();
132  #include "addRegionOption.H"
133  #include "setRootCase.H"
134  #include "createTime.H"
135  instantList timeDirs = timeSelector::select0(runTime, args);
136  #include "createNamedMesh.H"
137 
138  forAll(timeDirs, timeI)
139  {
140  runTime.setTime(timeDirs[timeI], timeI);
141  Info<< "Time = " << runTime.timeName() << endl;
142  mesh.readUpdate();
143 
144  volVectorField wallShearStress
145  (
146  IOobject
147  (
148  "wallShearStress",
149  runTime.timeName(),
150  mesh,
151  IOobject::NO_READ,
152  IOobject::AUTO_WRITE
153  ),
154  mesh,
156  (
157  "wallShearStress",
159  vector::zero
160  )
161  );
162 
163  IOobject UHeader
164  (
165  "U",
166  runTime.timeName(),
167  mesh,
168  IOobject::MUST_READ,
169  IOobject::NO_WRITE
170  );
171 
172  if (UHeader.headerOk())
173  {
174  Info<< "Reading field U\n" << endl;
176 
177  if
178  (
179  IOobject
180  (
182  runTime.constant(),
183  mesh
184  ).headerOk()
185  )
186  {
187  calcCompressible(mesh, runTime, U, wallShearStress);
188  }
189  else
190  {
191  calcIncompressible(mesh, runTime, U, wallShearStress);
192  }
193  }
194  else
195  {
196  Info<< " no U field" << endl;
197  }
198 
199  Info<< "Writing wall shear stress to field " << wallShearStress.name()
200  << nl << endl;
201 
202  wallShearStress.write();
203  }
204 
205  Info<< "End" << endl;
206 
207  return 0;
208 }
209 
210 
211 // ************************************************************************* //
pThermo
Info<< "Reading thermophysical properties\n"<< endl;autoPtr< psiThermo > pThermo(psiThermo::New(mesh))
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
singlePhaseTransportModel.H
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
turbulentTransportModel.H
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::instantList
List< instant > instantList
List of instants.
Definition: instantList.H:42
U
U
Definition: pEqn.H:46
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:48
dictName
const word dictName("particleTrackDict")
UHeader
IOobject UHeader("U", runTime.timeName(), mesh, IOobject::MUST_READ)
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::volSymmTensorField
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:58
Foam::Ostream::write
virtual Ostream & write(const token &)=0
Write next token to stream.
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
addRegionOption.H
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
createNamedMesh.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
compressibleCreatePhi.H
Creates and initialises the face-flux field phi.
rho
rho
Definition: pEqn.H:3
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
setRootCase.H
timeDirs
static instantList timeDirs
Definition: globalFoam.H:44
laminarTransport
singlePhaseTransportModel laminarTransport(U, phi)
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
createTime.H
fvCFD.H
args
Foam::argList args(argc, argv)
turbulentFluidThermoModel.H