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-2021 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 )
38 :
40 {}
41 
42 
44 (
45  const fvPatch& p,
47  const dictionary& dict
48 )
49 :
51 {
52  if (dict.found("value"))
53  {
54  fvPatchScalarField::operator=
55  (
56  scalarField("value", dict, p.size())
57  );
58  }
59  else
60  {
61  evaluate(Pstream::commsTypes::blocking);
62  }
63 }
64 
65 
67 (
69  const fvPatch& p,
71  const fvPatchFieldMapper& mapper
72 )
73 :
74  fixedJumpFvPatchField<scalar>(ptf, p, iF, mapper)
75 {}
76 
77 
79 (
82 )
83 :
85 {}
86 
87 
88 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
89 
91 {
92  if (this->updated())
93  {
94  return;
95  }
96 
97  if (this->cyclicPatch().owner())
98  {
100  label patchID = patch().index();
101 
102  const fixedJumpFvPatchScalarField& TbPatch =
103  refCast<const fixedJumpFvPatchScalarField>
104  (
105  thermo.T().boundaryField()[patchID]
106  );
107 
108  fixedJumpFvPatchScalarField& Tbp =
109  const_cast<fixedJumpFvPatchScalarField&>(TbPatch);
110 
111  // force update of jump
112  Tbp.updateCoeffs();
113 
114  const labelUList& faceCells = this->patch().faceCells();
115 
116  jump_ = thermo.he(Tbp.jump(), faceCells);
117  }
118 
120 }
121 
122 
124 {
126  writeEntry(os, "value", *this);
127 }
128 
129 
130 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
131 
132 namespace Foam
133 {
135  (
138  );
139 }
140 
141 // ************************************************************************* //
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:663
virtual void write(Ostream &) const
Write.
fluidReactionThermo & thermo
Definition: createFields.H:28
Base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:77
void evaluate(GeometricField< Type, PatchField, GeoMesh > &result, const Function1< Type > &func, const GeometricField< Type, PatchField, GeoMesh > &x)
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:175
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
bool updated() const
Return true if the boundary condition has already been updated.
Definition: fvPatchField.H:369
const Boundary & boundaryField() const
Return const-reference to the boundary field.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
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:63
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
Macros for easy insertion into run-time selection tables.
virtual volScalarField & he()=0
Enthalpy/Internal energy [J/kg].
energyJumpFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
const cyclicFvPatch & cyclicPatch() const
Return local reference cast into the cyclic patch.
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:99
Foam::fvPatchFieldMapper.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:351
virtual void write(Ostream &) const
Write.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
volScalarField scalarField(fieldObject, mesh)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
virtual const volScalarField & T() const =0
Temperature [K].
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:216
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static const basicThermo & lookupThermo(const fvPatchScalarField &pf)
Lookup the thermo associated with the given patch field.
Definition: basicThermo.C:83
volScalarField & p
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.
virtual void updateCoeffs()
Update the coefficients.