filmSurfaceVelocityFvPatchVectorField.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) 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 "mappedPatchBase.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
35 (
36  const fvPatch& p,
38  const dictionary& dict
39 )
40 :
41  mixedFvPatchField<vector>(p, iF, dict, false),
42  Cs_(dict.lookupOrDefault<scalar>("Cs", 0))
43 {
44  refValue() = Zero;
45  refGrad() = Zero;
46  valueFraction() = 0;
47 
48  if (dict.found("value"))
49  {
51  (
52  vectorField("value", dict, p.size())
53  );
54  }
55  else
56  {
57  // If the value entry is not present initialise to zero-gradient
58  fvPatchVectorField::operator=(patchInternalField());
59  }
60 }
61 
62 
65 (
67  const fvPatch& p,
69  const fvPatchFieldMapper& mapper
70 )
71 :
72  mixedFvPatchField<vector>(ptf, p, iF, mapper),
73  Cs_(ptf.Cs_)
74 {}
75 
76 
79 (
82 )
83 :
84  mixedFvPatchField<vector>(ptf, iF),
85  Cs_(ptf.Cs_)
86 {}
87 
88 
89 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
90 
92 {
93  if (updated())
94  {
95  return;
96  }
97 
98  // Since we're inside initEvaluate/evaluate there might be processor
99  // comms underway. Change the tag we use.
100  const int oldTag = UPstream::msgType();
101  UPstream::msgType() = oldTag + 1;
102 
103  // Get the coupling information from the mappedPatchBase
104  const mappedPatchBase& mpp = mappedPatchBase::getMap(patch().patch());
105  const label patchiNbr = mpp.nbrPolyPatch().index();
106  const fvPatch& patchNbr =
107  refCast<const fvMesh>(mpp.nbrMesh()).boundary()[patchiNbr];
108 
109  // Neighbour patch face-cell velocity
110  const vectorField UpNbr
111  (
112  patchNbr.lookupPatchField<volVectorField, scalar>("U")
113  .patchInternalField()
114  );
115 
116  // Set the reference value to the neighbouring fluid internal velocity
117  refValue() = mpp.fromNeighbour(UpNbr);
118 
119  // Remove the normal component of the surface vel
120  const vectorField n(patch().nf());
121  refValue() -= n*(n & refValue());
122 
123  // Lookup the momentum transport model
124  const compressibleMomentumTransportModel& transportModel =
125  db().lookupType<compressibleMomentumTransportModel>();
126 
127  // Lookup the neighbour momentum transport model
128  const compressibleMomentumTransportModel& transportModelNbr =
130 
131  // Patch laminar dynamic viscosity divided by delta
132  const tmp<scalarField> muEffByDelta
133  (
134  transportModel.rho().boundaryField()[patch().index()]
135  *transportModel.nuEff(patch().index())
136  *patch().deltaCoeffs()
137  );
138 
139  if (Cs_ > 0)
140  {
141  // Get the neighbour patch density
142  const tmp<scalarField> rhopNbr
143  (
144  mpp.fromNeighbour
145  (
146  transportModelNbr.rho().boundaryField()[patchiNbr]
147  )
148  );
149 
150  // Calculate the drag coefficient from the drag constant
151  // and the magnitude of the velocity difference
152  const scalarField Ds(Cs_*rhopNbr*mag(refValue() - *this));
153 
154  // Calculate the value-fraction from the balance between the
155  // external fluid drag and internal film stress
156  valueFraction() = Ds/(muEffByDelta + Ds);
157  }
158  else
159  {
160  // Get the neighbour patch laminar dynamic viscosity divided by delta
161  const tmp<scalarField> muEffByDeltaNbr
162  (
163  mpp.fromNeighbour
164  (
165  transportModelNbr.rho().boundaryField()[patchiNbr]
166  *transportModelNbr.nuEff(patchiNbr)
167  *patchNbr.deltaCoeffs()
168  )
169  );
170 
171  // Calculate the value-fraction from the balance between the
172  // external fluid and internal film stresses
173  valueFraction() = muEffByDeltaNbr()/(muEffByDelta + muEffByDeltaNbr());
174  }
175 
177 
178  // Restore tag
179  UPstream::msgType() = oldTag;
180 }
181 
182 
184 (
185  Ostream& os
186 ) const
187 {
189 
190  if (Cs_ > 0)
191  {
192  writeEntry(os, "Cs", Cs_);
193  }
194 
195  writeEntry(os, "value", *this);
196 }
197 
198 
199 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
200 
201 namespace Foam
202 {
204  (
207  );
208 }
209 
210 
211 // ************************************************************************* //
label n
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...
Generic GeometricField class.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:476
Base class for single-phase compressible turbulence models.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
filmSurfaceVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
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
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:251
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:204
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
virtual const scalarField & deltaCoeffs() const
Return the face - cell distance coefficient.
Definition: fvPatch.C:170
const GeometricField::Patch & lookupPatchField(const word &name) const
Lookup and return the patchField of the named field from the.
Engine which provides mapping between two patches.
const polyPatch & nbrPolyPatch() const
Get the patch to map from.
static const mappedPatchBase & getMap(const polyPatch &patch)
Cast the given polyPatch to a mappedPatchBase. Handle errors.
const polyMesh & nbrMesh() const
Get the mesh for the region to map from.
tmp< Field< Type > > fromNeighbour(const Field< Type > &nbrFld) const
Map/interpolate the neighbour patch field to this patch.
This boundary condition provides a base class for 'mixed' type boundary conditions,...
const Type & lookupType(const word &group=word::null) const
Lookup and return the object of the given Type.
label index() const
Return the index of this patch in the boundaryMesh.
A class for managing temporary objects.
Definition: tmp.H:55
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensioned< scalar > mag(const dimensioned< Type > &)
Field< vector > vectorField
Specialisation of Field<T> for vector.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dictionary dict
volScalarField & p