fixedShearStressFvPatchVectorField.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  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
31 #include "fvPatchFieldMapper.H"
32 #include "volFields.H"
33 #include "surfaceFields.H"
34 #include "turbulenceModel.H"
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
39 (
40  const fvPatch& p,
42 )
43 :
44  fixedValueFvPatchVectorField(p, iF),
45  tau0_(Zero)
46 {}
47 
48 
50 (
51  const fvPatch& p,
53  const dictionary& dict
54 )
55 :
56  fixedValueFvPatchVectorField(p, iF, dict, IOobjectOption::NO_READ),
57  tau0_(dict.getOrDefault<vector>("tau", Zero))
58 {
59  this->extrapolateInternal(); // Zero-gradient patch values
60 }
61 
62 
64 (
66  const fvPatch& p,
68  const fvPatchFieldMapper& mapper
69 )
70 :
71  fixedValueFvPatchVectorField(ptf, p, iF, mapper),
72  tau0_(ptf.tau0_)
73 {}
74 
75 
77 (
79 )
80 :
81  fixedValueFvPatchVectorField(ptf),
82  tau0_(ptf.tau0_)
83 {}
84 
85 
87 (
90 )
91 :
92  fixedValueFvPatchVectorField(ptf, iF),
93  tau0_(ptf.tau0_)
94 {}
95 
96 
97 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
98 
100 {
101  if (updated())
102  {
103  return;
104  }
105 
106  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
107  (
109  (
111  internalField().group()
112  )
113  );
114 
115  tmp<scalarField> nuEff(turbModel.nuEff(patch().index()));
116 
117  const vectorField Uc(patchInternalField());
118 
119  vector tauHat = tau0_/(mag(tau0_) + ROOTVSMALL);
120 
121  const scalarField& ry = patch().deltaCoeffs();
123  operator==(tauHat*(tauHat & (tau0_*(1.0/(ry*nuEff)) + Uc)));
124 
125  fixedValueFvPatchVectorField::updateCoeffs();
126 }
127 
128 
130 {
132  os.writeEntry("tau", tau0_);
134 }
135 
136 
137 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
138 
139 namespace Foam
140 {
142  (
144  fixedShearStressFvPatchVectorField
145  );
146 }
147 
148 // ************************************************************************* //
Foam::surfaceFields.
dictionary dict
fvPatchField< vector > fvPatchVectorField
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
This boundary condition sets a user-defined shear stress constant and uniform across a given patch by...
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
Macros for easy insertion into run-time selection tables.
constexpr const char *const group
Group name for atomic constants.
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
static const word propertiesName
Default name of the turbulence properties dictionary.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A FieldMapper for finite-volume patch fields.
Vector< scalar > vector
Definition: vector.H:57
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
OBJstream os(runTime.globalPath()/outputName)
fixedShearStressFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
Field< vector > vectorField
Specialisation of Field<T> for vector.
volScalarField & p
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127