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-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 "levelSet.H"
29 #include "volFields.H"
30 #include "fvMeshSubset.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
35 (
36  const fvPatch& p,
38  const dictionary& dict
39 )
40 :
42 {
43  if (dict.found("value"))
44  {
46  (
47  vectorField("value", dict, p.size())
48  );
49  }
50  else
51  {
53  (
55  );
56  }
57 }
58 
59 
61 (
63  const fvPatch& p,
65  const fvPatchFieldMapper& mapper
66 )
67 :
68  fixedValueInletOutletFvPatchField<vector>(ptf, p, iF, mapper)
69 {}
70 
71 
73 (
76 )
77 :
79 {}
80 
81 
82 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
83 
84 const Foam::fvMeshSubset&
86 {
87  const fvMesh& mesh = patch().boundaryMesh().mesh();
88  const label timeIndex = mesh.time().timeIndex();
89 
90  if
91  (
92  !faceCellSubset_.valid()
93  || (mesh.changing() && faceCellSubsetTimeIndex_ != timeIndex)
94  )
95  {
96  faceCellSubset_.reset(new fvMeshSubset(mesh));
97  faceCellSubset_->setCellSubset(patch().faceCells());
98  faceCellSubsetTimeIndex_ = timeIndex;
99 
100  // Ask for the tetBasePtIs to trigger all processors to build them.
101  // Without this, processors that do not contain this patch will
102  // generate a comms mismatch.
103  faceCellSubset_->subMesh().tetBasePtIs();
104  }
105 
106  return faceCellSubset_();
107 }
108 
109 
112 {
113  const waveSuperposition& waves = waveSuperposition::New(db());
114 
115  return
117  (
118  patch(),
119  waves.height(t, patch().Cf()),
120  waves.height(t, patch().patch().localPoints()),
121  waves.UGas(t, patch().Cf())(),
122  waves.UGas(t, patch().patch().localPoints())(),
123  waves.ULiquid(t, patch().Cf())(),
124  waves.ULiquid(t, patch().patch().localPoints())()
125  );
126 }
127 
128 
131 {
132  const waveSuperposition& waves = waveSuperposition::New(db());
133 
134  const fvMeshSubset& subset = faceCellSubset();
135  const fvMesh& meshs = subset.subMesh();
136  const label patchis = findIndex(subset.patchMap(), patch().index());
137 
138  const vectorField Us
139  (
141  (
142  meshs,
143  waves.height(t, meshs.cellCentres())(),
144  waves.height(t, meshs.points())(),
145  waves.UGas(t, meshs.cellCentres())(),
146  waves.UGas(t, meshs.points())(),
147  waves.ULiquid(t, meshs.cellCentres())(),
148  waves.ULiquid(t, meshs.points())()
149  )
150  );
151 
152  tmp<vectorField> tResult(new vectorField(patch().size()));
153  vectorField& result = tResult.ref();
154 
155  if (patchis != -1)
156  {
157  forAll(meshs.boundary()[patchis], is)
158  {
159  const label fs = is + meshs.boundary()[patchis].patch().start();
160  const label cs = meshs.boundary()[patchis].faceCells()[is];
161  const label f = subset.faceMap()[fs];
162  const label i = patch().patch().whichFace(f);
163  result[i] = Us[cs];
164  }
165  }
166 
167  return tResult;
168 }
169 
170 
172 {
173  if (updated())
174  {
175  return;
176  }
177 
178  const scalar t = db().time().userTimeValue();
179 
180  operator==(U(t));
181 
183 }
184 
185 
186 // * * * * * * * * * * * * * * Build Macro Function * * * * * * * * * * * * //
187 
188 namespace Foam
189 {
191  (
194  );
195 }
196 
197 // ************************************************************************* //
#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...
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:28
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
This boundary condition sets a fixed value. When the flow direction is inwards this acts exactly like...
Post-processing mesh subset tool. Given the original mesh and the list of selected cells,...
Definition: fvMeshSubset.H:74
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:406
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:893
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 tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:172
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
const polyMesh & mesh() const
Return the mesh reference.
bool changing() const
Is mesh changing.
Definition: polyMesh.H:488
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:403
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1361
const vectorField & cellCentres() const
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
A wrapper around a list of wave models. Superimposes the modelled values of elevation and velocity....
virtual tmp< scalarField > height(const scalar t, const vectorField &p) const
Get the height above the waves at a given time and global positions.
virtual tmp< vectorField > ULiquid(const scalar t, const vectorField &p) const
Get the liquid velocity at a given time and global positions.
virtual tmp< vectorField > UGas(const scalar t, const vectorField &p) const
Get the gas velocity at a given time and global positions.
static const waveSuperposition & New(const objectRegistry &db)
Return a reference to the wave model on the given database,.
This boundary condition provides a waveVelocity condition. This sets the velocity to that specified b...
waveVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
tmp< vectorField > U(const scalar t) const
Return the current modelled velocity field on the patch faces.
tmp< vectorField > Un(const scalar t) const
Return the current modelled velocity field in the neighbour.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
const fvMeshSubset & faceCellSubset() const
Access the face-cell subset.
U
Definition: pEqn.H:72
Namespace for OpenFOAM.
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
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Field< vector > vectorField
Specialisation of Field<T> for vector.
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
ListType subset(const UList< T > &select, const T &value, const ListType &)
Extract elements of List when select is a certain value.
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.
label timeIndex
Definition: getTimeIndex.H:4
labelList f(nPoints)
dictionary dict
volScalarField & p