Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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
32#include "volFields.H"
33#include "surfaceFields.H"
34#include "turbulenceModel.H"
35
36// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37
39(
40 const fvPatch& p,
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
70:
71 fixedValueFvPatchVectorField(ptf, p, iF, mapper),
72 tau0_(ptf.tau0_)
73{}
74
75
77(
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_);
135
136
137// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
138
139namespace Foam
140{
142 (
144 fixedShearStressFvPatchVectorField
145 );
146}
147
148// ************************************************************************* //
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...
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
This boundary condition sets a user-defined shear stress constant and uniform across a given patch by...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
fixedShearStressFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
A FieldMapper for finite-volume patch fields.
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
A class for managing temporary objects.
Definition tmp.H:75
Abstract base class for turbulence models (RAS, LES and laminar).
static const word propertiesName
Default name of the turbulence properties dictionary.
virtual tmp< volScalarField > nuEff() const =0
Return the effective viscosity.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
Namespace for OpenFOAM.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Vector< scalar > vector
Definition vector.H:57
fvPatchField< vector > fvPatchVectorField
dictionary dict
Foam::surfaceFields.