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-2024 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 "fieldMapper.H"
29 #include "volFields.H"
30 #include "surfaceFields.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
36 (
37  const fvPatch& p,
39  const dictionary& dict
40 )
41 :
42  mixedFvPatchVectorField(p, iF, dict, false),
43  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
44  rhoName_(dict.lookupOrDefault<word>("rho", "rho"))
45 {
47  (
48  vectorField("value", iF.dimensions(), dict, p.size())
49  );
50  refValue() = *this;
51  refGrad() = Zero;
52  valueFraction() = 0.0;
53 }
54 
55 
58 (
60  const fvPatch& p,
62  const fieldMapper& mapper
63 )
64 :
65  mixedFvPatchVectorField(ptf, p, iF, mapper),
66  phiName_(ptf.phiName_),
67  rhoName_(ptf.rhoName_)
68 {}
69 
70 
73 (
76 )
77 :
78  mixedFvPatchVectorField(pivpvf, iF),
79  phiName_(pivpvf.phiName_),
80  rhoName_(pivpvf.rhoName_)
81 {}
82 
83 
84 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
85 
87 {
88  if (updated())
89  {
90  return;
91  }
92 
93  const surfaceScalarField& phi =
94  db().lookupObject<surfaceScalarField>(phiName_);
95 
96  const fvsPatchField<scalar>& phip =
97  patch().patchField<surfaceScalarField, scalar>(phi);
98 
99  tmp<vectorField> n = patch().nf();
100  const Field<scalar>& magSf = patch().magSf();
101 
102  // Get the tangential component from the internalField (zero-gradient)
103  vectorField Ut(patchInternalField());
104  Ut -= n()*(Ut & n());
105 
106  if (phi.dimensions() == dimVolumetricFlux)
107  {
108  refValue() = Ut + n*phip/magSf;
109  }
110  else if (phi.dimensions() == dimMassFlux)
111  {
112  const fvPatchField<scalar>& rhop =
113  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
114 
115  refValue() = Ut + n*phip/(rhop*magSf);
116  }
117  else
118  {
120  << "dimensions of phi are not correct" << nl
121  << " on patch " << this->patch().name()
122  << " of field " << this->internalField().name()
123  << " in file " << this->internalField().objectPath()
124  << exit(FatalError);
125  }
126 
127  valueFraction() = neg(phip);
128 
129  mixedFvPatchVectorField::updateCoeffs();
130 }
131 
132 
134 (
135  Ostream& os
136 ) const
137 {
139  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
140  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
141  writeEntry(os, "value", *this);
142 }
143 
144 
145 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
146 
147 void Foam::pressureInletOutletParSlipVelocityFvPatchVectorField::operator=
148 (
149  const fvPatchField<vector>& pvf
150 )
151 {
153 }
154 
155 
156 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
157 
158 namespace Foam
159 {
161  (
164  );
165 }
166 
167 // ************************************************************************* //
label n
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...
const dimensionSet & dimensions() const
Return dimensions.
Generic GeometricField class.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:229
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:249
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:83
This velocity inlet/outlet boundary condition for pressure boundary where the pressure is specified....
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
pressureInletOutletParSlipVelocityFvPatchVectorField(const pressureInletOutletParSlipVelocityFvPatchVectorField &, const fvPatch &, const DimensionedField< vector, volMesh > &, const fieldMapper &)
Construct by mapping given.
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
const dimensionSet dimMassFlux
const dimensionSet dimVolumetricFlux
SurfaceField< scalar > surfaceScalarField
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
Field< vector > vectorField
Specialisation of Field<T> for vector.
dimensionedScalar neg(const dimensionedScalar &ds)
error FatalError
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
static const char nl
Definition: Ostream.H:266
dictionary dict
volScalarField & p
Foam::surfaceFields.