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-2021 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
30#include "dictionaryContent.H"
31
32// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33
34namespace Foam
35{
36namespace compressible
37{
38
39// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40
42(
43 const fvPatch& p,
45)
46:
47 fixedValueFvPatchField<scalar>(p, iF),
48 baffle_(),
49 dict_()
50{}
51
52
54(
56 const fvPatch& p,
58 const fvPatchFieldMapper& mapper
59)
60:
62 (
63 ptf,
64 p,
65 iF,
66 mapper
67 ),
68 baffle_(),
69 dict_(ptf.dict_)
70{}
71
72
74(
75 const fvPatch& p,
77 const dictionary& dict
78)
79:
80 fixedValueFvPatchField<scalar>(p, iF, dict),
81 baffle_(),
82 dict_
83 (
84 // Copy dictionary, but without "heavy" data chunks
85 dictionaryContent::copyDict
86 (
87 dict,
88 wordRes(), // allow
89 wordRes // deny
90 ({
91 "type", // redundant
92 "value"
93 })
94 )
95 )
96{
98
99 if (!baffle_)
100 {
101 baffle_.reset(baffle::New(p, dict).ptr());
102 }
103}
104
105
107(
110)
111:
112 fixedValueFvPatchField<scalar>(ptf, iF),
113 baffle_(),
114 dict_(ptf.dict_)
115{}
116
117
118// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
119
121{
122 if (this->updated())
123 {
124 return;
125 }
126
127 baffle_->evolve();
128
131 (
132 this->internalField().name()
133 ).boundaryFieldRef();
134
135 baffle_->vsm().mapToVolume(baffle_->T(), vfb);
136
138}
139
140
142{
144 dict_.write(os, false);
145}
146
147
148// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
149
151(
154);
155
156// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
157
158} // End namespace compressible
159} // End namespace Foam
160
161
162// ************************************************************************* //
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
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
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.
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:126
void write(Ostream &os, const bool subDict=true) const
Write dictionary, normally with sub-dictionary formatting.
Definition: dictionaryIO.C:206
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:210
const DimensionedField< Type, volMesh > & internalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:368
bool updated() const
Return true if the boundary condition has already been updated.
Definition: fvPatchField.H:392
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
Type & lookupObjectRef(const word &name, const bool recursive=false) const
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
bool compressible
Definition: pEqn.H:2
Namespace for OpenFOAM.
dictionary dict