viewFactor.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 | Copyright (C) 2011-2017 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::radiation::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 radiationModelviewFactor_H
46 #define radiationModelviewFactor_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 radiation
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/m2]
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  //- Disallow default bitwise copy construct
126  viewFactor(const viewFactor&);
127 
128  //- Disallow default bitwise assignment
129  void operator=(const viewFactor&);
130 
131 
132 public:
133 
134  //- Runtime type information
135  TypeName("viewFactor");
136 
137 
138  // Constructors
139 
140  //- Construct from components
141  viewFactor(const volScalarField& T);
142 
143  //- Construct from components
144  viewFactor(const dictionary& dict, const volScalarField& T);
145 
146 
147  //- Destructor
148  virtual ~viewFactor();
149 
150 
151  // Member functions
152 
153  // Edit
154 
155  //- Solve system of equation(s)
156  void calculate();
157 
158  //- Read radiation properties dictionary
159  bool read();
160 
161  //- Source term component (for power of T^4)
162  virtual tmp<volScalarField> Rp() const;
163 
164  //- Source term component (constant)
165  virtual tmp<volScalarField::Internal> Ru() const;
166 
167 
168  // Access
169 
170  //- Const access to total radiative heat flux field
171  inline const volScalarField& qr() const;
172 };
173 
174 
175 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
176 
177 #include "viewFactorI.H"
178 
179 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
180 
181 } // End namespace radiation
182 } // End namespace Foam
183 
184 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
185 
186 #endif
187 
188 // ************************************************************************* //
virtual tmp< volScalarField > Rp() const
Source term component (for power of T^4)
Definition: viewFactor.C:677
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
const volScalarField & qr() const
Const access to total radiative heat flux field.
Definition: viewFactor.H:27
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
TypeName("viewFactor")
Runtime type information.
virtual tmp< volScalarField::Internal > Ru() const
Source term component (constant)
Definition: viewFactor.C:705
Top level model for radiation modelling.
void calculate()
Solve system of equation(s)
Definition: viewFactor.C:394
virtual ~viewFactor()
Destructor.
Definition: viewFactor.C:352
bool read()
Read radiation properties dictionary.
Definition: viewFactor.C:358
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
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
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Namespace for OpenFOAM.