pressureInletOutletParSlipVelocityFvPatchVectorField.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-2018 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 "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "surfaceFields.H"
31 
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
37 (
38  const fvPatch& p,
40 )
41 :
42  mixedFvPatchVectorField(p, iF),
43  phiName_("phi"),
44  rhoName_("rho")
45 {
46  refValue() = *this;
47  refGrad() = Zero;
48  valueFraction() = 0.0;
49 }
50 
51 
54 (
56  const fvPatch& p,
58  const fvPatchFieldMapper& mapper
59 )
60 :
61  mixedFvPatchVectorField(ptf, p, iF, mapper),
62  phiName_(ptf.phiName_),
63  rhoName_(ptf.rhoName_)
64 {}
65 
66 
69 (
70  const fvPatch& p,
72  const dictionary& dict
73 )
74 :
75  mixedFvPatchVectorField(p, iF),
76  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
77  rhoName_(dict.lookupOrDefault<word>("rho", "rho"))
78 {
79  fvPatchVectorField::operator=(vectorField("value", dict, p.size()));
80  refValue() = *this;
81  refGrad() = Zero;
82  valueFraction() = 0.0;
83 }
84 
85 
88 (
90 )
91 :
92  mixedFvPatchVectorField(pivpvf),
93  phiName_(pivpvf.phiName_),
94  rhoName_(pivpvf.rhoName_)
95 {}
96 
97 
100 (
103 )
104 :
105  mixedFvPatchVectorField(pivpvf, iF),
106  phiName_(pivpvf.phiName_),
107  rhoName_(pivpvf.rhoName_)
108 {}
109 
110 
111 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
112 
114 {
115  if (updated())
116  {
117  return;
118  }
119 
120  const surfaceScalarField& phi =
121  db().lookupObject<surfaceScalarField>(phiName_);
122 
123  const fvsPatchField<scalar>& phip =
124  patch().patchField<surfaceScalarField, scalar>(phi);
125 
126  tmp<vectorField> n = patch().nf();
127  const Field<scalar>& magSf = patch().magSf();
128 
129  // Get the tangential component from the internalField (zero-gradient)
130  vectorField Ut(patchInternalField());
131  Ut -= n()*(Ut & n());
132 
133  if (phi.dimensions() == dimVelocity*dimArea)
134  {
135  refValue() = Ut + n*phip/magSf;
136  }
137  else if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
138  {
139  const fvPatchField<scalar>& rhop =
140  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
141 
142  refValue() = Ut + n*phip/(rhop*magSf);
143  }
144  else
145  {
147  << "dimensions of phi are not correct" << nl
148  << " on patch " << this->patch().name()
149  << " of field " << this->internalField().name()
150  << " in file " << this->internalField().objectPath()
151  << exit(FatalError);
152  }
153 
154  valueFraction() = 1.0 - pos0(phip);
155 
156  mixedFvPatchVectorField::updateCoeffs();
157 }
158 
159 
161 (
162  Ostream& os
163 ) const
164 {
166  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
167  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
168  writeEntry("value", os);
169 }
170 
171 
172 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
173 
174 void Foam::pressureInletOutletParSlipVelocityFvPatchVectorField::operator=
175 (
176  const fvPatchField<vector>& pvf
177 )
178 {
180 }
181 
182 
183 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
184 
185 namespace Foam
186 {
188  (
191  );
192 }
193 
194 // ************************************************************************* //
Foam::surfaceFields.
pressureInletOutletParSlipVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
surfaceScalarField & phi
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
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
volVectorField vectorField(fieldObject, mesh)
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:362
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
const dimensionSet & dimensions() const
Return dimensions.
A class for handling words, derived from string.
Definition: word.H:59
Foam::fvPatchFieldMapper.
static const zero Zero
Definition: zero.H:97
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)
This velocity inlet/outlet boundary condition for pressure boundary where the pressure is specified...
static const char nl
Definition: Ostream.H:265
const dimensionSet dimDensity
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:395
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
label n
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.
const dimensionSet dimVelocity