fixedPressureCompressibleDensityFvPatchScalarField.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-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 "fvPatchFieldMapper.H"
29 #include "surfaceFields.H"
30 #include "volFields.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
36 (
37  const fvPatch& p,
39  const dictionary& dict
40 )
41 :
42  fixedValueFvPatchField<scalar>(p, iF, dict),
43  pName_(dict.lookupOrDefault<word>("p", "p"))
44 {}
45 
46 
49 (
51  const fvPatch& p,
53  const fvPatchFieldMapper& mapper
54 )
55 :
56  fixedValueFvPatchField<scalar>(ptf, p, iF, mapper),
57  pName_(ptf.pName_)
58 {}
59 
60 
63 (
66 )
67 :
68  fixedValueFvPatchField<scalar>(ptf, iF),
69  pName_(ptf.pName_)
70 {}
71 
72 
73 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
74 
76 {
77  if (updated())
78  {
79  return;
80  }
81 
82  const fvPatchField<scalar>& pp =
83  patch().lookupPatchField<volScalarField, scalar>(pName_);
84 
85  const dictionary& thermoProps =
86  db().lookupObject<IOdictionary>("thermodynamicProperties");
87 
88  const scalar rholSat =
89  dimensionedScalar(thermoProps.lookup("rholSat")).value();
90 
91  const scalar pSat =
92  dimensionedScalar(thermoProps.lookup("pSat")).value();
93 
94  const scalar psil = dimensionedScalar(thermoProps.lookup("psil")).value();
95 
96  operator==(rholSat + psil*(pp - pSat));
97 
99 }
100 
101 
103 (
104  Ostream& os
105 ) const
106 {
108  writeEntryIfDifferent<word>(os, "p", "p", pName_);
109  writeEntry(os, "value", *this);
110 }
111 
112 
113 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
114 
115 namespace Foam
116 {
118  (
121  );
122 }
123 
124 // ************************************************************************* //
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...
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
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:160
const Type & value() const
Return const reference to value.
This boundary condition calculates a (liquid) compressible density as a function of pressure and flui...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
fixedPressureCompressibleDensityFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
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 void write(Ostream &) const
Write.
Definition: fvPatchField.C:231
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
A class for handling words, derived from string.
Definition: word.H:62
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dictionary dict
volScalarField & p
Foam::surfaceFields.