contactHeatFluxSource.H
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
26Class
27 Foam::fa::contactHeatFluxSource
28
29Group
30 grpFaOptionsSources
31
32Description
33 Applies contact heat flux between specified \c faMesh
34 and \c fvMesh within a specified region for compressible flows.
35
36Usage
37 Minimal example by using \c constant/faOptions:
38 \verbatim
39 contactHeatFluxSource1
40 {
41 // Mandatory entries (unmodifiable)
42 type contactHeatFluxSource;
43 Tprimary <TprimaryFieldName>;
44
45 // Optional entries (runtime modifiable)
46 T <Tname>;
47 thicknessLayers (<layer1> <layer2> ... <layerN>);
48
49 // Conditional optional entries (runtime modifiable)
50
51 // when the entry "thicknessLayers" is present
52 kappaLayers (<layer1> <layer2> ... <layerN>);
53
54 // Mandatory/Optional (inherited) entries
55 ...
56 }
57 \endverbatim
58
59 where the entries mean:
60 \table
61 Property | Description | Type | Reqd | Dflt
62 type | Type name: contactHeatFluxSource | word | yes | -
63 Tprimary | Name of primary temperature field | word | yes | -
64 T | Name of operand temperature field | word | no | T
65 thicknessLayers | List of thicknesses of layers | scalarList | no | -
66 kappaLayers | List of conductivities of layers | scalarList | cndtnl | -
67 \endtable
68
69 The inherited entries are elaborated in:
70 - \link faOption.H \endlink
71 - \link faceSetOption.H \endlink
72 - \link temperatureCoupledBase.H \endlink
73
74SourceFiles
75 contactHeatFluxSource.C
76
77\*---------------------------------------------------------------------------*/
78
79#ifndef fa_contactHeatFluxSource_H
80#define fa_contactHeatFluxSource_H
81
82#include "faOption.H"
83#include "Function1.H"
84#include "areaFields.H"
85#include "faceSetOption.H"
87
88// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
89
90namespace Foam
91{
92namespace fa
93{
94
95/*---------------------------------------------------------------------------*\
96 Class contactHeatFluxSource Declaration
97\*---------------------------------------------------------------------------*/
98
99class contactHeatFluxSource
100:
101 public fa::faceSetOption,
102 public temperatureCoupledBase
103{
104 // Private Data
105
106 //- Name of temperature field
107 word TName_;
108
109 //- Name of primary temperature field
110 word TprimaryName_;
111
112 //- Primary region temperature
113 const volScalarField& Tp_;
114
115 //- Temperature - wall [K]
116 areaScalarField Tw1_;
117
118 //- Thickness of layers
119 scalarList thicknessLayers_;
120
121 //- Conductivity of layers
122 scalarList kappaLayers_;
123
124 //- Total contact resistance
125 scalar contactRes_;
126
127 //- Current time index (used for updating)
128 label curTimeIndex_;
129
130
131 // Private Member Functions
132
133 //- Return htc from the primary region
135
136
137public:
138
139 //- Runtime type information
140 TypeName("contactHeatFluxSource");
141
142
143 // Constructors
144
145 //- Construct from explicit source name and mesh
147 (
148 const word& sourceName,
149 const word& modelType,
150 const dictionary& dict,
151 const fvPatch& patch
152 );
153
154 //- No copy construct
156
157 //- No copy assignment
158 void operator=(const contactHeatFluxSource&) = delete;
159
160
161 //- Destructor
162 virtual ~contactHeatFluxSource() = default;
163
164
165 // Member Functions
166
167 // Mapping functions
168
169 //- Map (and resize as needed) from self given a mapping object
170 virtual void autoMap
171 (
172 const fvPatchFieldMapper& mapper
173 )
174 {
176 }
177
178 //- Reverse map the given fvPatchField onto this fvPatchField
179 virtual void rmap
180 (
182 const labelList& map
183 )
184 {
186 }
187
188
189 // Evaluation
191 //- Add explicit contribution to compressible momentum equation
192 virtual void addSup
194 const areaScalarField& h,
195 const areaScalarField& rho,
196 faMatrix<scalar>& eqn,
197 const label fieldi
198 );
199
200
201 // IO
202
203 //- Read source dictionary
204 virtual bool read(const dictionary& dict);
206
207
208// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
209
210} // End namespace fa
211} // End namespace Foam
212
213// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
215
216#endif
217
218// ************************************************************************* //
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.mesh().time().globalPath()/"finiteArea-edges"));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A special matrix type and solver, designed for finite area solutions of scalar equations....
Definition: faMatrix.H:76
Applies contact heat flux between specified faMesh and fvMesh within a specified region for compressi...
void operator=(const contactHeatFluxSource &)=delete
No copy assignment.
virtual void addSup(const areaScalarField &h, const areaScalarField &rho, faMatrix< scalar > &eqn, const label fieldi)
Add explicit contribution to compressible momentum equation.
virtual bool read(const dictionary &dict)
Read source dictionary.
virtual ~contactHeatFluxSource()=default
Destructor.
contactHeatFluxSource(const contactHeatFluxSource &)=delete
No copy construct.
virtual void autoMap(const fvPatchFieldMapper &mapper)
Map (and resize as needed) from self given a mapping object.
TypeName("contactHeatFluxSource")
Runtime type information.
virtual void rmap(const fvPatchField< scalar > &fld, const labelList &map)
Reverse map the given fvPatchField onto this fvPatchField.
Intermediate abstract class for handling face-set options for the derived faOptions.
const fvPatch & patch() const noexcept
Return const access to fvPatch.
Definition: faOptionI.H:42
A FieldMapper for finite-volume patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:82
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
Common functions used in temperature coupled boundaries.
virtual void autoMap(const fvPatchFieldMapper &)=0
Map (and resize as needed) from self given a mapping object.
virtual void rmap(const fvPatchField< scalar > &, const labelList &)=0
Reverse map the given fvPatchField onto this fvPatchField.
A class for managing temporary objects.
Definition: tmp.H:65
A class for handling words, derived from Foam::string.
Definition: word.H:68
Namespace for OpenFOAM.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
List< scalar > scalarList
A List of scalars.
Definition: scalarList.H:64
GeometricField< scalar, faPatchField, areaMesh > areaScalarField
Definition: areaFieldsFwd.H:78
dictionary dict
volScalarField & h
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition: typeInfo.H:73