energyJumpFvPatchScalarField.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) 2012-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 "fixedJumpFvPatchFields.H"
29 #include "basicThermo.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
34 (
35  const fvPatch& p,
37  const dictionary& dict
38 )
39 :
40  fixedJumpFvPatchField<scalar>(p, iF)
41 {
42  if (dict.found("value"))
43  {
45  (
46  scalarField("value", dict, p.size())
47  );
48  }
49  else
50  {
52  }
53 }
54 
55 
57 (
59  const fvPatch& p,
61  const fvPatchFieldMapper& mapper
62 )
63 :
64  fixedJumpFvPatchField<scalar>(ptf, p, iF, mapper)
65 {}
66 
67 
69 (
72 )
73 :
74  fixedJumpFvPatchField<scalar>(ptf, iF)
75 {}
76 
77 
78 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
79 
81 {
82  if (this->updated())
83  {
84  return;
85  }
86 
87  if (this->cyclicPatch().owner())
88  {
90  label patchID = patch().index();
91 
92  const fixedJumpFvPatchScalarField& TbPatch =
93  refCast<const fixedJumpFvPatchScalarField>
94  (
95  thermo.T().boundaryField()[patchID]
96  );
97 
98  fixedJumpFvPatchScalarField& Tbp =
99  const_cast<fixedJumpFvPatchScalarField&>(TbPatch);
100 
101  // force update of jump
102  Tbp.updateCoeffs();
103 
104  const labelUList& faceCells = this->patch().faceCells();
105 
106  jump_ = thermo.he(Tbp.jump(), faceCells);
107  }
108 
110 }
111 
112 
114 {
116  writeEntry(os, "value", *this);
117 }
118 
119 
120 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
121 
122 namespace Foam
123 {
125  (
128  );
129 }
130 
131 // ************************************************************************* //
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 Boundary & boundaryField() const
Return const-reference to the boundary field.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
Base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:78
static const basicThermo & lookupThermo(const fvPatchScalarField &pf)
Lookup the thermo associated with the given patch field.
Definition: basicThermo.C:81
virtual volScalarField & he()=0
Enthalpy/Internal energy [J/kg].
virtual const volScalarField & T() const =0
Temperature [K].
virtual void evaluate(const Pstream::commsTypes commsType)
Evaluate the patch field.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
This boundary condition provides an energy jump condition across a pair of coupled patches....
virtual void write(Ostream &) const
Write.
virtual void updateCoeffs()
Update the coefficients.
energyJumpFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
This boundary condition provides a jump condition, using the cyclic condition as a base.
virtual void write(Ostream &) const
Write.
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
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
Namespace for OpenFOAM.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dictionary dict
volScalarField & p
fluidMulticomponentThermo & thermo
Definition: createFields.H:31