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 | Copyright (C) 2012-2017 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 (
46  const fvPatch& p,
48  const fvPatchFieldMapper& mapper
49 )
50 :
51  fixedJumpFvPatchField<scalar>(ptf, p, iF, mapper)
52 {}
53 
54 
56 (
57  const fvPatch& p,
59  const dictionary& dict
60 )
61 :
63 {
64  if (dict.found("value"))
65  {
66  fvPatchScalarField::operator=
67  (
68  scalarField("value", dict, p.size())
69  );
70  }
71  else
72  {
73  evaluate(Pstream::commsTypes::blocking);
74  }
75 }
76 
77 
79 (
81 )
82 :
84 {}
85 
86 
88 (
91 )
92 :
94 {}
95 
96 
97 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
98 
100 {
101  if (this->updated())
102  {
103  return;
104  }
105 
106  if (this->cyclicPatch().owner())
107  {
109  label patchID = patch().index();
110 
111  const scalarField& pp = thermo.p().boundaryField()[patchID];
112  const fixedJumpFvPatchScalarField& TbPatch =
113  refCast<const fixedJumpFvPatchScalarField>
114  (
115  thermo.T().boundaryField()[patchID]
116  );
117 
118  fixedJumpFvPatchScalarField& Tbp =
119  const_cast<fixedJumpFvPatchScalarField&>(TbPatch);
120 
121  // force update of jump
122  Tbp.updateCoeffs();
123 
124  const labelUList& faceCells = this->patch().faceCells();
125 
126  jump_ = thermo.he(pp, Tbp.jump(), faceCells);
127  }
128 
130 }
131 
132 
134 {
136  this->writeEntry("value", os);
137 }
138 
139 
140 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
141 
142 namespace Foam
143 {
145  (
148  );
149 }
150 
151 // ************************************************************************* //
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:431
virtual void write(Ostream &) const
Write.
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
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:52
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:179
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
bool updated() const
Return true if the boundary condition has already been updated.
Definition: fvPatchField.H:371
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:163
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:61
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition: Field.C:726
virtual const volScalarField & T() const
Temperature [K].
Definition: basicThermo.C:489
Macros for easy insertion into run-time selection tables.
virtual volScalarField & he()=0
Enthalpy/Internal energy [J/kg].
virtual volScalarField & p()
Pressure [Pa].
Definition: basicThermo.C:477
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
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:93
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:61
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:340
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
virtual void write(Ostream &) const
Write.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
volScalarField scalarField(fieldObject, mesh)
static Table::iterator lookupThermo(const dictionary &thermoTypeDict, Table *tablePtr, const int nCmpt, const char *cmptNames[], const word &thermoTypeName)
Generic lookup for thermodynamics package thermoTypeName.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:311
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
volScalarField & p
Namespace for OpenFOAM.
virtual void updateCoeffs()
Update the coefficients.