Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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,
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(
84(
87)
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 =
110 (
111 thermo.T().boundaryField()[patchID]
112 );
113
114 fixedJumpFvPatchScalarField& Tbp =
115 const_cast<fixedJumpFvPatchScalarField&>(TbPatch);
116
117 // force update of jump
118 Tbp.evaluate(Pstream::commsTypes::buffered);
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)
131
132
134{
138
139
140// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
141
142namespace Foam
143{
145 (
147 energyJumpFvPatchScalarField
148 );
149}
150
151// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
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:59
@ buffered
"buffered" : (MPI_Bsend, MPI_Recv)
Definition UPstream.H:82
Abstract base-class for fluid and solid thermodynamic properties.
Definition basicThermo.H:62
static const basicThermo & lookupThermo(const fvPatchScalarField &pf)
virtual void evaluate(const Pstream::commsTypes commsType)
const cyclicFvPatch & cyclicPatch() const
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
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 > &)
Construct from patch and internal field.
Smooth ATC in cells next to a set of patches supplied by type.
Definition faceCells.H:55
This boundary condition provides a jump condition, using the cyclic condition as a base.
virtual void write(Ostream &) const
Write.
virtual void setJump(const Field< scalar > &jump)
fixedJumpFvPatchField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
const fvPatch & patch() const noexcept
Return the patch.
bool updated() const noexcept
True if the boundary condition has already been updated.
A FieldMapper for finite-volume patch fields.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
bool readValueEntry(const dictionary &dict, IOobjectOption::readOption readOpt=IOobjectOption::LAZY_READ)
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)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
UList< label > labelUList
A UList of labels.
Definition UList.H:75
fvPatchField< scalar > fvPatchScalarField
dictionary dict