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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
31 #include "basicThermo.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
36 (
37  const fvPatch& p,
39 )
40 :
42 {}
43 
44 
46 (
48  const fvPatch& p,
50  const fvPatchFieldMapper& mapper
51 )
52 :
53  fixedJumpAMIFvPatchField<scalar>(ptf, p, iF, mapper)
54 {}
55 
56 
58 (
59  const fvPatch& p,
61  const dictionary& dict
62 )
63 :
65 {
66  if (dict.found("value"))
67  {
68  fvPatchScalarField::operator=
69  (
70  scalarField("value", dict, p.size())
71  );
72  }
73  else
74  {
75  evaluate(Pstream::commsTypes::blocking);
76  }
77 }
78 
79 
81 (
83 )
84 :
86 {}
87 
88 
90 (
93 )
94 :
96 {}
97 
98 
99 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
100 
102 {
103  if (this->updated())
104  {
105  return;
106  }
107 
108  if (this->cyclicAMIPatch().owner())
109  {
111  label patchID = patch().index();
112 
113  const scalarField& pp = thermo.p().boundaryField()[patchID];
114  const fixedJumpAMIFvPatchScalarField& TbPatch =
115  refCast<const fixedJumpAMIFvPatchScalarField>
116  (
117  thermo.T().boundaryField()[patchID]
118  );
119 
120  fixedJumpAMIFvPatchScalarField& Tbp =
121  const_cast<fixedJumpAMIFvPatchScalarField&>(TbPatch);
122 
123  // force update of jump
124  Tbp.evaluate(Pstream::commsTypes::blocking);
125 
126  const labelUList& faceCells = this->patch().faceCells();
127 
128  jump_ = thermo.he(pp, Tbp.jump(), faceCells);
129  }
130 
132 }
133 
134 
136 {
138  this->writeEntry("value", os);
139 }
140 
141 
142 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
143 
144 namespace Foam
145 {
147  (
150  );
151 }
152 
153 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:51
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::energyJumpAMIFvPatchScalarField::energyJumpAMIFvPatchScalarField
energyJumpAMIFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: energyJumpAMIFvPatchScalarField.C:36
p
volScalarField & p
Definition: createFieldRefs.H:8
basicThermo.H
Foam::energyJumpAMIFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: energyJumpAMIFvPatchScalarField.C:135
Foam::fixedJumpAMIFvPatchField< scalar >::jump_
Field< scalar > jump_
"jump" field
Definition: fixedJumpAMIFvPatchField.H:108
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::energyJumpAMIFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients.
Definition: energyJumpAMIFvPatchScalarField.C:101
Foam::basicThermo
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:63
fixedJumpAMIFvPatchFields.H
Foam::fixedJumpAMIFvPatchField< scalar >
Foam::basicThermo::lookupThermo
static const basicThermo & lookupThermo(const fvPatchScalarField &pf)
Definition: basicThermo.C:446
Foam::fvPatchField< scalar >::updated
bool updated() const
Return true if the boundary condition has already been updated.
Definition: fvPatchField.H:387
Foam::Field< scalar >
patchID
label patchID
Definition: boundaryProcessorFaPatchPoints.H:5
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::fvPatch::index
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:197
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::cyclicAMIFvPatchField< scalar >::cyclicAMIPatch
const cyclicAMIFvPatch & cyclicAMIPatch() const
Return local reference cast into the cyclic AMI patch.
Definition: cyclicAMIFvPatchField.H:186
energyJumpAMIFvPatchScalarField.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::fvPatchField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:321
Foam::UList< label >
Foam::fixedJumpAMIFvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fixedJumpAMIFvPatchField.C:171
Foam::fvPatch::faceCells
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:113
Foam::energyJumpAMIFvPatchScalarField
This boundary condition provides an energy jump condition across a pair of coupled patches with an ar...
Definition: energyJumpAMIFvPatchScalarField.H:59
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::fvPatchField< scalar >::patch
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:357
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::basicThermo::he
virtual volScalarField & he()=0
Enthalpy/Internal energy [J/kg].
Foam::faceCells
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:56
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::UPstream::commsTypes::blocking
Foam::stringOps::evaluate
string evaluate(label fieldWidth, const std::string &s, size_t pos=0, size_t len=std::string::npos)
String evaluation with specified (positive, non-zero) field width.
Definition: stringOpsEvaluate.C:37