viewFactor.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2019 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::radiationModels::viewFactor
26 
27 Description
28  View factor radiation model. The system solved is: C q = b
29  where:
30  Cij = deltaij/Ej - (1/Ej - 1)Fij
31  q = heat flux
32  b = A eb - Ho
33  and:
34  eb = sigma*T^4
35  Ej = emissivity
36  Aij = deltaij - Fij
37  Fij = view factor matrix
38 
39 
40 SourceFiles
41  viewFactor.C
42 
43 \*---------------------------------------------------------------------------*/
44 
45 #ifndef radiation_viewFactor_H
46 #define radiation_viewFactor_H
47 
48 #include "radiationModel.H"
49 #include "singleCellFvMesh.H"
50 #include "scalarMatrices.H"
51 #include "globalIndex.H"
52 #include "scalarListIOList.H"
53 #include "mapDistribute.H"
54 #include "volFields.H"
55 
56 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57 
58 namespace Foam
59 {
60 namespace radiationModels
61 {
62 
63 /*---------------------------------------------------------------------------*\
64  Class viewFactor Declaration
65 \*---------------------------------------------------------------------------*/
66 
67 class viewFactor
68 :
69  public radiationModel
70 {
71  // Private Data
72 
73  //- Agglomeration List
74  labelListIOList finalAgglom_;
75 
76  //- Map distributed
78 
79  //- Coarse mesh
80  singleCellFvMesh coarseMesh_;
81 
82  //- Net radiative heat flux [W/m^2]
83  volScalarField qr_;
84 
85  //- View factor matrix
87 
88  //- Inverse of C matrix
90 
91  //- Selected patches
92  labelList selectedPatches_;
93 
94  //- Total global coarse faces
95  label totalNCoarseFaces_;
96 
97  //- Total local coarse faces
98  label nLocalCoarseFaces_;
99 
100  //- Constant emissivity
101  bool constEmissivity_;
102 
103  //- Iterations Counter
104  label iterCounter_;
105 
106  //- Pivot Indices for LU decomposition
107  labelList pivotIndices_;
108 
109 
110  // Private Member Functions
111 
112  //- Initialise
113  void initialise();
114 
115  //- Insert view factors into main matrix
116  void insertMatrixElements
117  (
118  const globalIndex& index,
119  const label fromProci,
120  const labelListList& globalFaceFaces,
121  const scalarListList& viewFactors,
122  scalarSquareMatrix& matrix
123  );
124 
125 
126 public:
127 
128  //- Runtime type information
129  TypeName("viewFactor");
130 
131 
132  // Constructors
133 
134  //- Construct from components
135  viewFactor(const volScalarField& T);
136 
137  //- Construct from components
138  viewFactor(const dictionary& dict, const volScalarField& T);
139 
140  //- Disallow default bitwise copy construction
141  viewFactor(const viewFactor&) = delete;
142 
143 
144  //- Destructor
145  virtual ~viewFactor();
146 
147 
148  // Member Functions
149 
150  // Edit
151 
152  //- Solve system of equation(s)
153  void calculate();
154 
155  //- Read radiation properties dictionary
156  bool read();
157 
158  //- Source term component (for power of T^4)
159  virtual tmp<volScalarField> Rp() const;
160 
161  //- Source term component (constant)
162  virtual tmp<volScalarField::Internal> Ru() const;
163 
164 
165  // Access
166 
167  //- Const access to total radiative heat flux field
168  inline const volScalarField& qr() const;
169 
170 
171  // Member Operators
172 
173  //- Disallow default bitwise assignment
174  void operator=(const viewFactor&) = delete;
175 };
176 
177 
178 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
179 
180 #include "viewFactorI.H"
181 
182 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
183 
184 } // End namespace radiationModels
185 } // End namespace Foam
186 
187 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
188 
189 #endif
190 
191 // ************************************************************************* //
dictionary dict
fvMesh as subset of other mesh. Consists of one cell and all original bounday faces. Useful when manipulating boundary data. Single internal cell only needed to be able to manipulate in a standard way.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
TypeName("viewFactor")
Runtime type information.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
const volScalarField & qr() const
Const access to total radiative heat flux field.
Definition: viewFactor.H:27
void calculate()
Solve system of equation(s)
Definition: viewFactor.C:395
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
Top level model for radiation modelling.
virtual tmp< volScalarField::Internal > Ru() const
Source term component (constant)
Definition: viewFactor.C:694
virtual ~viewFactor()
Destructor.
Definition: viewFactor.C:353
bool read()
Read radiation properties dictionary.
Definition: viewFactor.C:359
void operator=(const viewFactor &)=delete
Disallow default bitwise assignment.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
viewFactor(const volScalarField &T)
Construct from components.
Definition: viewFactor.C:239
virtual tmp< volScalarField > Rp() const
Source term component (for power of T^4)
Definition: viewFactor.C:678
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
A class for managing temporary objects.
Definition: PtrList.H:53
View factor radiation model. The system solved is: C q = b where: Cij = deltaij/Ej - (1/Ej - 1)Fij q ...
Definition: viewFactor.H:66
Namespace for OpenFOAM.