fixedShearStressFvPatchVectorField.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 \*---------------------------------------------------------------------------*/
25 
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "surfaceFields.H"
31 #include "momentumTransportModel.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
36 (
37  const fvPatch& p,
39  const dictionary& dict
40 )
41 :
42  fixedValueFvPatchVectorField(p, iF, dict, false),
43  tau0_(dict.lookupOrDefault<vector>("tau", Zero))
44 {
45  fvPatchField<vector>::operator=(patchInternalField());
46 }
47 
48 
50 (
52  const fvPatch& p,
54  const fvPatchFieldMapper& mapper
55 )
56 :
57  fixedValueFvPatchVectorField(ptf, p, iF, mapper),
58  tau0_(ptf.tau0_)
59 {}
60 
61 
63 (
66 )
67 :
68  fixedValueFvPatchVectorField(ptf, iF),
69  tau0_(ptf.tau0_)
70 {}
71 
72 
73 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
74 
76 {
77  if (updated())
78  {
79  return;
80  }
81 
82  const momentumTransportModel& turbModel =
83  db().lookupType<momentumTransportModel>(internalField().group());
84 
85  const scalarField nuEff(turbModel.nuEff(patch().index()));
86 
87  const vectorField Uc(patchInternalField());
88 
89  vector tauHat = tau0_/(mag(tau0_) + rootVSmall);
90 
91  const scalarField& ry = patch().deltaCoeffs();
92 
93  operator==(tauHat*(tauHat & (tau0_*(1.0/(ry*nuEff)) + Uc)));
94 
95  fixedValueFvPatchVectorField::updateCoeffs();
96 }
97 
98 
100 {
102  writeEntry(os, "tau", tau0_);
103  writeEntry(os, "value", *this);
104 }
105 
106 
107 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
108 
109 namespace Foam
110 {
112  (
115  );
116 }
117 
118 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
Set a constant shear stress as tau0 = -nuEff dU/dn.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
fixedShearStressFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:231
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:251
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
Abstract base class for turbulence models (RAS, LES and laminar).
virtual tmp< volScalarField > nuEff() const =0
Return the effective viscosity.
const scalar nuEff
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensioned< scalar > mag(const dimensioned< Type > &)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dictionary dict
volScalarField & p
Foam::surfaceFields.