waveVelocityFvPatchVectorField.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) 2017-2021 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 "levelSet.H"
29 #include "volFields.H"
30 #include "fvMeshSubset.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
35 (
36  const fvPatch& p,
38 )
39 :
41 {}
42 
43 
45 (
46  const fvPatch& p,
48  const dictionary& dict
49 )
50 :
52 {
53  if (dict.found("value"))
54  {
56  (
57  vectorField("value", dict, p.size())
58  );
59  }
60  else
61  {
63  (
64  patchInternalField()
65  );
66  }
67 }
68 
69 
71 (
73  const fvPatch& p,
75  const fvPatchFieldMapper& mapper
76 )
77 :
79 {}
80 
81 
83 (
86 )
87 :
89 {}
90 
91 
92 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
93 
94 const Foam::fvMeshSubset&
96 {
97  const fvMesh& mesh = patch().boundaryMesh().mesh();
98  const label timeIndex = mesh.time().timeIndex();
99 
100  if
101  (
102  !faceCellSubset_.valid()
103  || (mesh.changing() && faceCellSubsetTimeIndex_ != timeIndex)
104  )
105  {
106  faceCellSubset_.reset(new fvMeshSubset(mesh));
107  faceCellSubset_->setCellSubset(patch().faceCells());
108  faceCellSubsetTimeIndex_ = timeIndex;
109 
110  // Ask for the tetBasePtIs to trigger all processors to build them.
111  // Without this, processors that do not contain this patch will
112  // generate a comms mismatch.
113  faceCellSubset_->subMesh().tetBasePtIs();
114  }
115 
116  return faceCellSubset_();
117 }
118 
119 
122 {
123  const waveSuperposition& waves = waveSuperposition::New(db());
124 
125  return
127  (
128  patch(),
129  waves.height(t, patch().Cf()),
130  waves.height(t, patch().patch().localPoints()),
131  waves.UGas(t, patch().Cf())(),
132  waves.UGas(t, patch().patch().localPoints())(),
133  waves.ULiquid(t, patch().Cf())(),
134  waves.ULiquid(t, patch().patch().localPoints())()
135  );
136 }
137 
138 
141 {
142  const waveSuperposition& waves = waveSuperposition::New(db());
143 
145  const fvMesh& meshs = subset.subMesh();
146  const label patchis = findIndex(subset.patchMap(), patch().index());
147 
148  const vectorField Us
149  (
151  (
152  meshs,
153  waves.height(t, meshs.cellCentres())(),
154  waves.height(t, meshs.points())(),
155  waves.UGas(t, meshs.cellCentres())(),
156  waves.UGas(t, meshs.points())(),
157  waves.ULiquid(t, meshs.cellCentres())(),
158  waves.ULiquid(t, meshs.points())()
159  )
160  );
161 
162  tmp<vectorField> tResult(new vectorField(patch().size()));
163  vectorField& result = tResult.ref();
164 
165  if (patchis != -1)
166  {
167  forAll(meshs.boundary()[patchis], is)
168  {
169  const label fs = is + meshs.boundary()[patchis].patch().start();
170  const label cs = meshs.boundary()[patchis].faceCells()[is];
171  const label f = subset.faceMap()[fs];
172  const label i = patch().patch().whichFace(f);
173  result[i] = Us[cs];
174  }
175  }
176 
177  return tResult;
178 }
179 
180 
182 {
183  if (updated())
184  {
185  return;
186  }
187 
188  const scalar t = db().time().userTimeValue();
189 
190  operator==(U(t));
191 
193 }
194 
195 
196 // * * * * * * * * * * * * * * Build Macro Function * * * * * * * * * * * * //
197 
198 namespace Foam
199 {
201  (
204  );
205 }
206 
207 // ************************************************************************* //
const fvBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: fvPatch.H:181
bool changing() const
Is mesh changing (topology changing and/or moving)
Definition: polyMesh.H:537
A wrapper around a list of wave models. Superimposes the modelled values of elevation and velocity...
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual tmp< vectorField > UGas(const scalar t, const vectorField &p) const
Get the gas velocity at a given time and global positions.
waveVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
const fvMeshSubset & faceCellSubset() const
Access the face-cell subset.
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:175
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
bool updated() const
Return true if the boundary condition has already been updated.
Definition: fvPatchField.H:369
tmp< vectorField > Un(const scalar t) const
Return the current modelled velocity field in the neighbour.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
tmp< vectorField > U(const scalar t) const
Return the current modelled velocity field on the patch faces.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
virtual tmp< vectorField > ULiquid(const scalar t, const vectorField &p) const
Get the liquid velocity at a given time and global positions.
volVectorField vectorField(fieldObject, mesh)
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:372
fvMesh & mesh
Macros for easy insertion into run-time selection tables.
const labelList & faceMap() const
Return face map.
virtual void operator==(const fvPatchField< Type > &)
Definition: fvPatchField.C:424
scalar userTimeValue() const
Return current user time value.
Definition: Time.C:845
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1211
const labelList & patchMap() const
Return patch map.
const fvMesh & subMesh() const
Return reference to subset mesh.
virtual tmp< scalarField > height(const scalar t, const vectorField &p) const
Get the height above the waves at a given time and global positions.
const fvMesh & mesh() const
Return the mesh reference.
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:139
Foam::fvPatchFieldMapper.
This boundary condition provides a waveVelocity condition. This sets the velocity to that specified b...
const vectorField & cellCentres() const
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:351
static const waveSuperposition & New(const objectRegistry &db)
Return a reference to the wave model on the given database,.
const Time & time() const
Return time.
tmp< Field< Type > > levelSetAverage(const fvMesh &mesh, const scalarField &levelC, const scalarField &levelP, const Field< Type > &positiveC, const Field< Type > &positiveP, const Field< Type > &negativeC, const Field< Type > &negativeP)
Calculate the average value of two fields, one on each side of a level set.
labelList f(nPoints)
Post-processing mesh subset tool. Given the original mesh and the list of selected cells...
Definition: fvMeshSubset.H:73
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:28
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:216
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
label timeIndex
Definition: getTimeIndex.H:4
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:147
ListType subset(const UList< T > &select, const T &value, const ListType &)
Extract elements of List when select is a certain value.
Field< vector > vectorField
Specialisation of Field<T> for vector.
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
label size() const
Return the number of elements in the UList.
Definition: ListI.H:171
Namespace for OpenFOAM.
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:386
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:800