energyJumpAMIFvPatchScalarField.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-------------------------------------------------------------------------------
11License
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
32#include "basicThermo.H"
33
34// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35
37(
38 const fvPatch& p,
40)
41:
42 fixedJumpAMIFvPatchField<scalar>(p, iF)
43{}
44
45
47(
49 const fvPatch& p,
51 const fvPatchFieldMapper& mapper
52)
53:
54 fixedJumpAMIFvPatchField<scalar>(ptf, p, iF, mapper)
55{}
56
57
59(
60 const fvPatch& p,
62 const dictionary& dict
63)
64:
65 fixedJumpAMIFvPatchField<scalar>(p, iF)
66{
67 if (dict.found("value"))
68 {
70 (
71 scalarField("value", dict, p.size())
72 );
73 }
74 else
75 {
77 }
78}
79
80
82(
84)
85:
86 fixedJumpAMIFvPatchField<scalar>(ptf)
87{}
88
89
91(
94)
95:
96 fixedJumpAMIFvPatchField<scalar>(ptf, iF)
97{}
98
99
100// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
101
103{
104 if (this->updated())
105 {
106 return;
107 }
108
109 if (this->cyclicAMIPatch().owner())
110 {
112 label patchID = patch().index();
113
114 const scalarField& pp = thermo.p().boundaryField()[patchID];
115 const fixedJumpAMIFvPatchScalarField& TbPatch =
116 refCast<const fixedJumpAMIFvPatchScalarField>
117 (
118 thermo.T().boundaryField()[patchID]
119 );
120
121 fixedJumpAMIFvPatchScalarField& Tbp =
122 const_cast<fixedJumpAMIFvPatchScalarField&>(TbPatch);
123
124 // force update of jump
125 Tbp.evaluate(Pstream::commsTypes::blocking);
126
127 const labelUList& faceCells = this->patch().faceCells();
128
129 jump_ =
130 thermo.he(pp, Tbp+Tbp.jump(), faceCells)
131 - thermo.he(pp, Tbp, faceCells);
132 }
133
135}
136
137
139{
141 this->writeEntry("value", os);
142}
143
144
145// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
146
147namespace Foam
148{
150 (
153 );
154}
155
156// ************************************************************************* //
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...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:66
static const basicThermo & lookupThermo(const fvPatchScalarField &pf)
Definition: basicThermo.C:454
virtual void evaluate(const Pstream::commsTypes commsType)
Evaluate the patch field.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
This boundary condition provides an energy jump condition across a pair of coupled patches with an ar...
virtual void updateCoeffs()
Update the coefficients.
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:59
This boundary condition provides a jump condition, across non-conformal cyclic path-pairs,...
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
friend Ostream & operator(Ostream &, const fvPatchField< scalar > &)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
Namespace for OpenFOAM.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dictionary dict