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-2022 OpenCFD Ltd.
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 
30 #include "dictionaryContent.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace compressible
37 {
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
42 (
43  const fvPatch& p,
45 )
46 :
48  baffle_(nullptr),
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_(nullptr),
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_(nullptr),
82  dict_
83  (
84  // Copy dictionary, but without "heavy" data chunks
85  dictionaryContent::copyDict
86  (
87  dict,
88  wordList(), // allow
89  wordList // deny
90  ({
91  "type", // redundant
92  "value"
93  })
94  )
95  )
96 {
97  if (!baffle_)
98  {
99  baffle_.reset(baffleType::New(p.boundaryMesh().mesh(), dict_));
100  }
101 }
102 
103 
105 (
106  const thermalShellFvPatchScalarField& ptf,
108 )
109 :
110  fixedValueFvPatchField<scalar>(ptf, iF),
111  baffle_(nullptr),
112  dict_(ptf.dict_)
113 {}
114 
115 
116 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
117 
119 {
120  if (this->updated())
121  {
122  return;
123  }
124 
125  baffle_->evolve();
126 
127  scalarField& pfld = *this;
129  baffle_->vsm().mapToVolumePatch(baffle_->T(), pfld, patch().index());
130 
132 }
133 
134 
136 {
138  dict_.write(os, false);
139 }
140 
141 
142 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
143 
145 (
147  thermalShellFvPatchScalarField
148 );
149 
150 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
151 
152 } // End namespace compressible
153 } // End namespace Foam
154 
155 
156 // ************************************************************************* //
dictionary dict
virtual void write(Ostream &) const
Write.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const fvPatch & patch() const noexcept
Return the patch.
Definition: fvPatchField.H:269
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Macros for easy insertion into run-time selection tables.
makePatchTypeField(fvPatchScalarField, alphatFixedDmdtWallBoilingWallFunctionFvPatchScalarField)
bool updated() const noexcept
True if the boundary condition has already been updated.
Definition: fvPatchField.H:296
A FieldMapper for finite-volume patch fields.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A wrapper for dictionary content, without operators that could affect inheritance patterns...
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
bool compressible
Definition: pEqn.H:2
OBJstream os(runTime.globalPath()/outputName)
This boundary condition provides a coupled temperature condition between a primary region (3D mesh) a...
void write(Ostream &os, const bool subDict=true) const
Write dictionary, normally with sub-dictionary formatting.
Definition: dictionaryIO.C:204
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:309
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
thermalShellFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
volScalarField & p
label index() const noexcept
The index of this patch in the boundary mesh.
Definition: fvPatch.H:218
Namespace for OpenFOAM.