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-2022 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 (
55  const fvPatch& p,
57  const dictionary& dict
58 )
59 :
60  mixedFvPatchVectorField(p, iF),
61  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
62  rhoName_(dict.lookupOrDefault<word>("rho", "rho"))
63 {
64  fvPatchVectorField::operator=(vectorField("value", dict, p.size()));
65  refValue() = *this;
66  refGrad() = Zero;
67  valueFraction() = 0.0;
68 }
69 
70 
73 (
75  const fvPatch& p,
77  const fvPatchFieldMapper& mapper
78 )
79 :
80  mixedFvPatchVectorField(ptf, p, iF, mapper),
81  phiName_(ptf.phiName_),
82  rhoName_(ptf.rhoName_)
83 {}
84 
85 
88 (
91 )
92 :
93  mixedFvPatchVectorField(pivpvf, iF),
94  phiName_(pivpvf.phiName_),
95  rhoName_(pivpvf.rhoName_)
96 {}
97 
98 
99 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
100 
102 {
103  if (updated())
104  {
105  return;
106  }
107 
108  const surfaceScalarField& phi =
109  db().lookupObject<surfaceScalarField>(phiName_);
110 
111  const fvsPatchField<scalar>& phip =
112  patch().patchField<surfaceScalarField, scalar>(phi);
113 
114  tmp<vectorField> n = patch().nf();
115  const Field<scalar>& magSf = patch().magSf();
116 
117  // Get the tangential component from the internalField (zero-gradient)
118  vectorField Ut(patchInternalField());
119  Ut -= n()*(Ut & n());
120 
121  if (phi.dimensions() == dimFlux)
122  {
123  refValue() = Ut + n*phip/magSf;
124  }
125  else if (phi.dimensions() == dimMassFlux)
126  {
127  const fvPatchField<scalar>& rhop =
128  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
129 
130  refValue() = Ut + n*phip/(rhop*magSf);
131  }
132  else
133  {
135  << "dimensions of phi are not correct" << nl
136  << " on patch " << this->patch().name()
137  << " of field " << this->internalField().name()
138  << " in file " << this->internalField().objectPath()
139  << exit(FatalError);
140  }
141 
142  valueFraction() = neg(phip);
143 
144  mixedFvPatchVectorField::updateCoeffs();
145 }
146 
147 
149 (
150  Ostream& os
151 ) const
152 {
154  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
155  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
156  writeEntry(os, "value", *this);
157 }
158 
159 
160 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
161 
162 void Foam::pressureInletOutletParSlipVelocityFvPatchVectorField::operator=
163 (
164  const fvPatchField<vector>& pvf
165 )
166 {
168 }
169 
170 
171 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
172 
173 namespace Foam
174 {
176  (
179  );
180 }
181 
182 // ************************************************************************* //
Foam::surfaceFields.
pressureInletOutletParSlipVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
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:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
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
volVectorField vectorField(fieldObject, mesh)
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:243
dimensionedScalar neg(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
const dimensionSet & dimensions() const
Return dimensions.
A class for handling words, derived from string.
Definition: word.H:59
Foam::fvPatchFieldMapper.
const dimensionSet dimFlux
phi
Definition: correctPhi.H:3
static const zero Zero
Definition: zero.H:97
virtual label size() const
Return size.
Definition: fvPatch.H:157
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
This velocity inlet/outlet boundary condition for pressure boundary where the pressure is specified...
static const char nl
Definition: Ostream.H:260
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
const dimensionSet dimMassFlux
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:258
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
Namespace for OpenFOAM.