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 | 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 \*---------------------------------------------------------------------------*/
25 
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "surfaceFields.H"
31 #include "turbulenceModel.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
36 (
37  const fvPatch& p,
39 )
40 :
41  fixedValueFvPatchVectorField(p, iF),
42  tau0_(vector::zero)
43 {}
44 
45 
47 (
48  const fvPatch& p,
50  const dictionary& dict
51 )
52 :
53  fixedValueFvPatchVectorField(p, iF),
54  tau0_(dict.lookupOrDefault<vector>("tau", vector::zero))
55 {
56  fvPatchField<vector>::operator=(patchInternalField());
57 }
58 
59 
61 (
63  const fvPatch& p,
65  const fvPatchFieldMapper& mapper
66 )
67 :
68  fixedValueFvPatchVectorField(ptf, p, iF, mapper),
69  tau0_(ptf.tau0_)
70 {}
71 
72 
74 (
76 )
77 :
78  fixedValueFvPatchVectorField(ptf),
79  tau0_(ptf.tau0_)
80 {}
81 
82 
84 (
87 )
88 :
89  fixedValueFvPatchVectorField(ptf, iF),
90  tau0_(ptf.tau0_)
91 {}
92 
93 
94 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95 
97 {
98  if (updated())
99  {
100  return;
101  }
102 
103  const turbulenceModel& turbModel = db().lookupObject<turbulenceModel>
104  (
106  (
109  )
110  );
111 
112  tmp<scalarField> nuEff(turbModel.nuEff(patch().index()));
113 
114  const vectorField Uc(patchInternalField());
115 
116  vector tauHat = tau0_/(mag(tau0_) + ROOTVSMALL);
117 
118  const scalarField& ry = patch().deltaCoeffs();
119 
120  operator==(tauHat*(tauHat & (tau0_*(1.0/(ry*nuEff())) + Uc)));
121 
122  fixedValueFvPatchVectorField::updateCoeffs();
123 }
124 
125 
127 {
129  os.writeKeyword("tau") << tau0_ << token::END_STATEMENT << nl;
130  writeEntry("value", os);
131 }
132 
133 
134 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
135 
136 namespace Foam
137 {
139  (
142  );
143 }
144 
145 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
volFields.H
Foam::token::END_STATEMENT
@ END_STATEMENT
Definition: token.H:99
Foam::IOobject::groupName
static word groupName(Name name, const word &group)
p
p
Definition: pEqn.H:62
Foam::fixedShearStressFvPatchVectorField::fixedShearStressFvPatchVectorField
fixedShearStressFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Definition: fixedShearStressFvPatchVectorField.C:36
Foam::fixedShearStressFvPatchVectorField
Set a constant shear stress as tau0 = -nuEff dU/dn.
Definition: fixedShearStressFvPatchVectorField.H:52
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
dimensionedInternalField
rDeltaT dimensionedInternalField()
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:97
Foam::operator==
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Foam::constant::atomic::group
const char *const group
Group name for atomic constants.
Definition: atomicConstants.C:40
surfaceFields.H
Foam::surfaceFields.
fvPatchFieldMapper.H
Foam::fixedShearStressFvPatchVectorField::tau0_
const vector tau0_
Constant shear stress.
Definition: fixedShearStressFvPatchVectorField.H:59
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::turbulenceModel
Abstract base class for turbulence models (RAS, LES and laminar).
Definition: turbulenceModel.H:60
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::makePatchTypeField
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
fixedShearStressFvPatchVectorField.H
Foam::fixedShearStressFvPatchVectorField::write
virtual void write(Ostream &) const
Write.
Definition: fixedShearStressFvPatchVectorField.C:126
Foam::Vector< scalar >
Foam::Ostream::writeKeyword
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
Foam::fixedShearStressFvPatchVectorField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fixedShearStressFvPatchVectorField.C:96
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::turbulenceModel::nuEff
virtual tmp< volScalarField > nuEff() const =0
Return the effective viscosity.
turbulenceModel.H
write
Tcoeff write()
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51