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-2024 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 "mappedFvPatchBaseBase.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", dimless, 0))
43 {
44  refValue() = Zero;
45  refGrad() = Zero;
46  valueFraction() = 0;
47 
48  if (dict.found("value"))
49  {
51  (
52  vectorField("value", iF.dimensions(), 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 fieldMapper& 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 mappedFvPatchBaseBase& mapper =
106  const fvPatch& patchNbr = mapper.nbrFvPatch();
107 
108  // Neighbour patch face-cell velocity
109  const vectorField UpNbr
110  (
111  patchNbr.lookupPatchField<volVectorField, scalar>("U")
112  .patchInternalField()
113  );
114 
115  // Set the reference value to the neighbouring fluid internal velocity
116  refValue() = mapper.fromNeighbour(UpNbr);
117 
118  // Remove the normal component of the surface vel
119  const vectorField n(patch().nf());
120  refValue() -= n*(n & refValue());
121 
122  // Lookup the momentum transport model
123  const compressibleMomentumTransportModel& transportModel =
124  db().lookupType<compressibleMomentumTransportModel>();
125 
126  // Lookup the neighbour momentum transport model
127  const compressibleMomentumTransportModel& transportModelNbr =
129 
130  // Patch laminar dynamic viscosity divided by delta
131  const tmp<scalarField> muEffByDelta
132  (
133  transportModel.rho().boundaryField()[patch().index()]
134  *transportModel.nuEff(patch().index())
135  *patch().deltaCoeffs()
136  );
137 
138  if (Cs_ > 0)
139  {
140  // Get the neighbour patch density
141  const tmp<scalarField> rhopNbr
142  (
143  mapper.fromNeighbour
144  (
145  transportModelNbr.rho().boundaryField()[patchNbr.index()]
146  )
147  );
148 
149  // Calculate the drag coefficient from the drag constant
150  // and the magnitude of the velocity difference
151  const scalarField Ds(Cs_*rhopNbr*mag(refValue() - *this));
152 
153  // Calculate the value-fraction from the balance between the
154  // external fluid drag and internal film stress
155  valueFraction() = Ds/(muEffByDelta + Ds);
156  }
157  else
158  {
159  // Get the neighbour patch laminar dynamic viscosity divided by delta
160  const tmp<scalarField> muEffByDeltaNbr
161  (
162  mapper.fromNeighbour
163  (
164  transportModelNbr.rho().boundaryField()[patchNbr.index()]
165  *transportModelNbr.nuEff(patchNbr.index())
166  *patchNbr.deltaCoeffs()
167  )
168  );
169 
170  // Calculate the value-fraction from the balance between the
171  // external fluid and internal film stresses
172  valueFraction() = muEffByDeltaNbr()/(muEffByDelta + muEffByDeltaNbr());
173  }
174 
176 
177  // Restore tag
178  UPstream::msgType() = oldTag;
179 }
180 
181 
183 (
184  Ostream& os
185 ) const
186 {
188 
189  if (Cs_ > 0)
190  {
191  writeEntry(os, "Cs", Cs_);
192  }
193 
194  writeEntry(os, "value", *this);
195 }
196 
197 
198 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
199 
200 namespace Foam
201 {
203  (
206  );
207 }
208 
209 
210 // ************************************************************************* //
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...
const dimensionSet & dimensions() const
Return dimensions.
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:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
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.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:229
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:249
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:202
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:156
virtual const scalarField & deltaCoeffs() const
Return the face - cell distance coefficient.
Definition: fvPatch.C:176
const GeometricField::Patch & lookupPatchField(const word &name) const
Lookup and return the patchField of the named field from the.
Base class for fv patches that provide mapping between two fv patches.
const fvPatch & nbrFvPatch() const
Get the patch to map from.
static const mappedFvPatchBaseBase & getMap(const fvPatch &patch)
Cast the given fvPatch to a mappedFvPatchBaseBase. Handle errors.
const fvMesh & nbrMesh() const
Get the mesh for the region to map from.
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.
A class for managing temporary objects.
Definition: tmp.H:55
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
const dimensionSet dimless
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