Loading...
Searching...
No Matches
thermalShellFvPatchScalarField.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) 2019-2025 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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
29#include "dictionaryContent.H"
30#include "regionProperties.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
36{
37namespace compressible
38{
39
40// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
41
42void thermalShellFvPatchScalarField::create_baffle()
43{
44 if (!baffle_)
45 {
46 baffle_.reset
47 (
48 baffleType::New(this->patch().boundaryMesh().mesh(), dict_)
49 );
50 }
51}
52
53
54// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55
57(
58 const fvPatch& p,
61:
62 parent_bctype(p, iF),
63 dict_()
64{}
65
66
68(
70 const fvPatch& p,
72 const fvPatchFieldMapper& mapper
74:
75 parent_bctype(ptf, p, iF, mapper),
76 dict_(ptf.dict_)
77{}
78
79
81(
82 const fvPatch& p,
84 const dictionary& dict
85)
86:
87 parent_bctype(p, iF, dict),
88 dict_
89 (
90 // Copy dictionary, but without "heavy" data chunks
91 dictionaryContent::copyDict
92 (
93 dict,
94 wordList(), // allow
95 wordList // deny
96 ({
97 "type", // redundant
98 "value"
99 })
100 )
101 )
102{
103 // Create baffle now.
104 // Lazy evaluation has issues with loading the finite-area fields
106 {
107 create_baffle();
108 }
109}
110
111
113(
114 const thermalShellFvPatchScalarField& ptf,
116)
117:
118 parent_bctype(ptf, iF),
119 dict_(ptf.dict_)
120{}
121
122
123// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
124
126{
127 if (this->updated())
128 {
129 return;
130 }
131
132 scalarField& pfld = *this;
133
134 // Create baffle if needed, but ignore if regionFaModels are disabled
135 if (!baffle_)
136 {
138 {
139 create_baffle();
140 }
141 else
142 {
143 static bool warned = false;
144 if (!warned)
145 {
146 warned = true;
148 << "Ignoring, regionFaModels are disabled" << endl;
149 }
150 return;
151 }
152 }
153 auto& baffle = baffle_();
154
155 baffle.evolve();
157 baffle.vsm().mapToVolumePatch(baffle.T(), pfld, patch().index());
158
160}
161
162
164{
166 dict_.write(os, false);
167}
168
169
170// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
171
173(
175 thermalShellFvPatchScalarField
176);
177
178// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
179
180} // End namespace compressible
181} // End namespace Foam
182
183
184// ************************************************************************* //
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:59
This boundary condition provides a coupled temperature condition between a primary region (3D mesh) a...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
thermalShellFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
A wrapper for dictionary content, without operators that could affect inheritance patterns.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
void write(Ostream &os, const bool subDict=true) const
Write dictionary, normally with sub-dictionary formatting.
virtual void write(Ostream &) const
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.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
static autoPtr< thermalShellModel > New(const fvMesh &mesh, const dictionary &dict)
Return a reference to the selected model using dictionary.
volScalarField & p
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
#define InfoInFunction
Report an information message using Foam::Info.
bool allowFaModels() noexcept
The enable/disable state for regionFaModel (default: true).
Namespace for OpenFOAM.
List< word > wordList
List of word.
Definition fileName.H:60
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
fvPatchField< scalar > fvPatchScalarField
dictionary dict