waveVelocityFvPatchVectorField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2017 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 "surfaceFields.H"
30 #include "volFields.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
35 (
36  const fvPatch& p,
38 )
39 :
40  directionMixedFvPatchVectorField(p, iF),
41  waves_(db())
42 {
43  refValue() = Zero;
44  refGrad() = Zero;
45  valueFraction() = Zero;
46 }
47 
48 
50 (
51  const fvPatch& p,
53  const dictionary& dict
54 )
55 :
56  directionMixedFvPatchVectorField(p, iF),
57  waves_(db(), dict)
58 {
59  if (dict.found("value"))
60  {
61  fvPatchVectorField::operator=(vectorField("value", dict, p.size()));
62  }
63  else
64  {
65  fvPatchVectorField::operator=(patchInternalField());
66  }
67 
68  refValue() = *this;
69  refGrad() = Zero;
70  valueFraction() = Zero;
71 }
72 
73 
75 (
77  const fvPatch& p,
79  const fvPatchFieldMapper& mapper
80 )
81 :
82  directionMixedFvPatchVectorField(ptf, p, iF, mapper),
83  waves_(ptf.waves_)
84 {}
85 
86 
88 (
90 )
91 :
92  directionMixedFvPatchVectorField(ptf),
93  waves_(ptf.waves_)
94 {}
95 
96 
98 (
101 )
102 :
103  directionMixedFvPatchVectorField(ptf, iF),
104  waves_(ptf.waves_)
105 {}
106 
107 
108 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
109 
111 {
112  const scalar t = db().time().timeOutputValue();
113 
114  return
116  (
117  patch(),
118  waves_.height(t, patch().Cf()),
119  waves_.height(t, patch().patch().localPoints()),
120  waves_.UGas(t, patch().Cf())(),
121  waves_.UGas(t, patch().patch().localPoints())(),
122  waves_.ULiquid(t, patch().Cf())(),
123  waves_.ULiquid(t, patch().patch().localPoints())()
124  );
125 }
126 
127 
129 {
130  if (updated())
131  {
132  return;
133  }
134 
135  refValue() = U();
136 
137  const scalarField& phip =
138  patch().lookupPatchField<surfaceScalarField, scalar>("phi");
139 
140  const symmTensorField nnp(sqr(patch().nf()));
141 
142  // Fix all components if inflow, just normal if outflow
143  valueFraction() = (1 - pos0(phip))*I + pos0(phip)*nnp;
144 
145  directionMixedFvPatchVectorField::updateCoeffs();
146  directionMixedFvPatchVectorField::evaluate();
147 }
148 
149 
151 (
152  Ostream& os
153 ) const
154 {
156  waves_.write(os);
157 }
158 
159 
160 // * * * * * * * * * * * * * * Build Macro Function * * * * * * * * * * * * //
161 
162 namespace Foam
163 {
165  (
168  );
169 }
170 
171 // ************************************************************************* //
Foam::surfaceFields.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:431
tmp< vectorField > UGas(const scalar t, const vectorField &p) const
Get the gas velocity at a given time and global positions.
void write(Ostream &) const
Write.
waveVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual void write(Ostream &) const
Write.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
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)
Macros for easy insertion into run-time selection tables.
tmp< DimensionedField< Type, volMesh > > levelSetAverage(const fvMesh &mesh, const scalarField &levelC, const scalarField &levelP, const DimensionedField< Type, volMesh > &positiveC, const DimensionedField< Type, pointMesh > &positiveP, const DimensionedField< Type, volMesh > &negativeC, const DimensionedField< Type, pointMesh > &negativeP)
Calculate the average value of two fields, one on each side of a level set.
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
void write(Ostream &, const label, const dictionary &)
Write with dictionary lookup.
tmp< scalarField > height(const scalar t, const vectorField &p) const
Get the height above the waves at a given time and global positions.
static const Identity< scalar > I
Definition: Identity.H:93
Foam::fvPatchFieldMapper.
static const zero Zero
Definition: zero.H:91
This boundary condition provides a waveVelocity condition. This sets the velocity to that specified b...
virtual label size() const
Return size.
Definition: fvPatch.H:161
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
dimensionedScalar pos0(const dimensionedScalar &ds)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Namespace for OpenFOAM.
tmp< vectorField > U() const
Return the current modelled velocity field.