movingMappedWallVelocityFvPatchVectorField.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-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 
27 #include "volFields.H"
28 #include "surfaceFields.H"
29 #include "fvcMeshPhi.H"
30 #include "mappedFvPatchBaseBase.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
37 (
38  const fvPatch& p,
40  const dictionary& dict
41 )
42 :
43  fixedValueFvPatchVectorField(p, iF, dict)
44 {}
45 
46 
49 (
51  const fvPatch& p,
53  const fieldMapper& mapper
54 )
55 :
56  fixedValueFvPatchVectorField(ptf, p, iF, mapper, false)
57 {
58  mapper(*this, ptf, vector::zero);
59 }
60 
61 
64 (
67 )
68 :
69  fixedValueFvPatchVectorField(mwvpvf, iF)
70 {}
71 
72 
73 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
74 
76 (
77  const fvPatchVectorField& ptf,
78  const fieldMapper& mapper
79 )
80 {
81  mapper(*this, ptf, vector::zero);
82 }
83 
84 
86 {
87  if (updated())
88  {
89  return;
90  }
91 
92  const fvPatch& fvp = patch();
93 
94  const mappedFvPatchBaseBase& mapper =
96  const fvMesh& nbrMesh = mapper.nbrMesh();
97  const fvPatch& nbrFvp = mapper.nbrFvPatch();
98 
99  if (nbrMesh.moving())
100  {
101  // Calculate the normal velocity on this mesh from the mesh flux
102  const volVectorField& U =
103  static_cast<const volVectorField&>(internalField());
104  tmp<scalarField> Un =
105  U.mesh().moving()
106  ? fvc::meshPhi(U, fvp.index())/(fvp.magSf() + vSmall)
107  : tmp<scalarField>(new scalarField(fvp.size(), Zero));
108 
109  // Calculate and map the mesh velocity from the neighbour
110  vectorField nbrCf0(nbrFvp.size());
111  forAll(nbrCf0, nbrPatchFacei)
112  {
113  nbrCf0[nbrPatchFacei] =
114  nbrMesh.faces()
115  [
116  nbrMesh.polyFacesBf()[nbrFvp.index()][nbrPatchFacei]
117  ].centre(nbrMesh.oldPoints());
118  }
119  const vectorField nbrUp
120  (
121  (nbrFvp.Cf() - nbrCf0)/db().time().deltaTValue()
122  );
123  const vectorField Up(mapper.fromNeighbour(nbrUp));
124 
125  // Set the velocity equal to the normal component from this mesh plus
126  // the tangential component from the neighbouring mesh. That way the
127  // flux is exactly preserved but the neighbour can impart shear. The
128  // normal components should be the same anyway, otherwise the mapped
129  // patches are moving apart. Using the normal component from this mesh
130  // just prevents mapping inaccuracies from affecting the flux.
131  const vectorField n(fvp.nf());
132  vectorField::operator=(Up + n*(Un - (n & Up)));
133  }
134 
135  fixedValueFvPatchVectorField::updateCoeffs();
136 }
137 
138 
140 {
142  writeEntry(os, "value", *this);
143 }
144 
145 
146 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
147 
148 namespace Foam
149 {
151  (
154  );
155 }
156 
157 // ************************************************************************* //
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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...
void operator=(const Field< vector > &)
Definition: Field.C:550
Generic GeometricField class.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
static const Form zero
Definition: VectorSpace.H:118
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
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const GeometricBoundaryField< label, fvsPatchField, surfaceMesh > & polyFacesBf() const
Return face-poly-face addressing.
Definition: fvMesh.C:880
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
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
virtual label size() const
Return size.
Definition: fvPatch.H:138
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:156
const scalarField & magSf() const
Return face area magnitudes.
Definition: fvPatch.C:142
tmp< vectorField > nf() const
Return face normals.
Definition: fvPatch.C:130
const vectorField & Cf() const
Return face centres.
Definition: fvPatch.C:105
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 no-slip velocity condition for mapped walls. The wall velocity is ...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
movingMappedWallVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void map(const fvPatchVectorField &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1326
virtual const pointField & oldPoints() const
Return old points for mesh motion.
Definition: polyMesh.C:1351
bool moving() const
Is mesh moving.
Definition: polyMesh.H:476
A class for managing temporary objects.
Definition: tmp.H:55
Calculate the mesh motion flux and convert fluxes from absolute to relative and back.
U
Definition: pEqn.H:72
tmp< surfaceScalarField > meshPhi(const volVectorField &U)
Definition: fvcMeshPhi.C:34
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dictionary dict
volScalarField & p
Foam::surfaceFields.