energyJumpFvPatchScalarField.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2012-2017 OpenFOAM Foundation
9  Copyright (C) 2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
31 #include "fixedJumpFvPatchFields.H"
32 #include "basicThermo.H"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
37 (
38  const fvPatch& p,
40 )
41 :
42  fixedJumpFvPatchField<scalar>(p, iF)
43 {}
44 
45 
47 (
49  const fvPatch& p,
51  const fvPatchFieldMapper& mapper
52 )
53 :
54  fixedJumpFvPatchField<scalar>(ptf, p, iF, mapper)
55 {}
56 
57 
59 (
60  const fvPatch& p,
62  const dictionary& dict
63 )
64 :
65  fixedJumpFvPatchField<scalar>(p, iF)
66 {
67  if (!this->readValueEntry(dict))
68  {
70  }
71 }
72 
73 
75 (
76  const energyJumpFvPatchScalarField& ptf
77 )
78 :
79  fixedJumpFvPatchField<scalar>(ptf)
80 {}
81 
82 
84 (
87 )
88 :
89  fixedJumpFvPatchField<scalar>(ptf, iF)
90 {}
91 
92 
93 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
94 
96 {
97  if (this->updated())
98  {
99  return;
100  }
101 
102  if (this->cyclicPatch().owner())
103  {
104  const basicThermo& thermo = basicThermo::lookupThermo(*this);
105  label patchID = patch().index();
106 
107  const scalarField& pp = thermo.p().boundaryField()[patchID];
108  const fixedJumpFvPatchScalarField& TbPatch =
109  refCast<const fixedJumpFvPatchScalarField>
110  (
112  );
113 
114  fixedJumpFvPatchScalarField& Tbp =
115  const_cast<fixedJumpFvPatchScalarField&>(TbPatch);
116 
117  // force update of jump
119 
120  const labelUList& faceCells = this->patch().faceCells();
121 
122  setJump
123  (
124  thermo.he(pp, Tbp+Tbp.jump(), faceCells)
125  - thermo.he(pp, Tbp, faceCells)
126  );
127  }
128 
130 }
131 
132 
134 {
137 }
138 
139 
140 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
141 
142 namespace Foam
143 {
145  (
147  energyJumpFvPatchScalarField
148  );
149 }
150 
151 // ************************************************************************* //
bool readValueEntry(const dictionary &dict, IOobjectOption::readOption readOpt=IOobjectOption::LAZY_READ)
Read the "value" entry into *this.
Definition: fvPatchField.C:32
dictionary dict
virtual void write(Ostream &) const
Write.
"blocking" : (MPI_Bsend, MPI_Recv)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
This boundary condition provides an energy jump condition across a pair of coupled patches...
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
virtual const volScalarField & T() const
Temperature [K].
Definition: basicThermo.C:613
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
Macros for easy insertion into run-time selection tables.
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
virtual volScalarField & he()=0
Enthalpy/Internal energy [J/kg].
virtual volScalarField & p()
Pressure [Pa].
Definition: basicThermo.C:601
psiReactionThermo & thermo
Definition: createFields.H:28
void evaluate()
Evaluate boundary conditions.
energyJumpFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
fvPatchField< scalar > fvPatchScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A FieldMapper for finite-volume patch fields.
virtual void evaluate(const Pstream::commsTypes commsType)
Evaluate the patch field.
virtual void write(Ostream &) const
Write.
OBJstream os(runTime.globalPath()/outputName)
This boundary condition provides a jump condition, using the cyclic condition as a base...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:309
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
static const basicThermo & lookupThermo(const fvPatchScalarField &pf)
Definition: basicThermo.C:448
volScalarField & p
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
virtual void updateCoeffs()
Update the coefficients.