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-2015 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 "volFields.H"
54 #include "IOmapDistribute.H"
55 #include "solarLoad.H"
56 
57 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58 
59 namespace Foam
60 {
61 namespace radiation
62 {
63 
64 /*---------------------------------------------------------------------------*\
65  Class viewFactor Declaration
66 \*---------------------------------------------------------------------------*/
67 
68 class viewFactor
69 :
70  public radiationModel
71 {
72 protected:
73 
74 // Static data
75 
76  //- Static name for view factor walls
77  static const word viewFactorWalls;
78 
79 
80  // Private data
81 
82  //- Agglomeration List
84 
85  //- Map distributed
87 
88  //- Coarse mesh
90 
91  //- Net radiative heat flux [W/m2]
93 
94  //- View factor matrix
96 
97  //- Inverse of C matrix
99 
100  //- Selected patches
102 
103  //- Total global coarse faces
105 
106  //- Total local coarse faces
108 
109  //- Constant emissivity
110  bool constEmissivity_;
111 
112  //- Iterations Counter
114 
115  //- Pivot Indices for LU decomposition
117 
118  //- Use Solar Load model
119  bool useSolarLoad_;
120 
121  //- Solar load radiation model
123 
124 
125  // Private Member Functions
126 
127  //- Initialise
128  void initialise();
129 
130  //- Insert view factors into main matrix
132  (
133  const globalIndex& index,
134  const label fromProcI,
135  const labelListList& globalFaceFaces,
136  const scalarListList& viewFactors,
137  scalarSquareMatrix& matrix
138  );
139 
140  //- Disallow default bitwise copy construct
141  viewFactor(const viewFactor&);
142 
143  //- Disallow default bitwise assignment
144  void operator=(const viewFactor&);
145 
146 
147 public:
148 
149  //- Runtime type information
150  TypeName("viewFactor");
151 
152 
153  // Constructors
154 
155  //- Construct from components
156  viewFactor(const volScalarField& T);
157 
158  //- Construct from components
159  viewFactor(const dictionary& dict, const volScalarField& T);
160 
161 
162  //- Destructor
163  virtual ~viewFactor();
164 
165 
166  // Member functions
167 
168  // Edit
169 
170  //- Solve system of equation(s)
171  void calculate();
172 
173  //- Read radiation properties dictionary
174  bool read();
175 
176  //- Source term component (for power of T^4)
177  virtual tmp<volScalarField> Rp() const;
178 
179  //- Source term component (constant)
180  virtual tmp<DimensionedField<scalar, volMesh> > Ru() const;
181 
182 
183  // Access
184 
185  //- Const access to total radiative heat flux field
186  inline const volScalarField& Qr() const;
187 };
188 
189 
190 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
191 
192 #include "viewFactorI.H"
193 
194 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
195 
196 } // End namespace radiation
197 } // End namespace Foam
198 
199 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
200 
201 #endif
202 
203 // ************************************************************************* //
volFields.H
Foam::radiation::viewFactor::insertMatrixElements
void insertMatrixElements(const globalIndex &index, const label fromProcI, const labelListList &globalFaceFaces, const scalarListList &viewFactors, scalarSquareMatrix &matrix)
Insert view factors into main matrix.
Definition: viewFactor.C:365
solarLoad.H
Foam::radiation::viewFactor::CLU_
autoPtr< scalarSquareMatrix > CLU_
Inverse of C matrix.
Definition: viewFactor.H:97
Foam::radiation::viewFactor::Qr
const volScalarField & Qr() const
Const access to total radiative heat flux field.
Definition: viewFactorI.H:26
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::radiation::viewFactor::TypeName
TypeName("viewFactor")
Runtime type information.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::radiation::viewFactor::useSolarLoad_
bool useSolarLoad_
Use Solar Load model.
Definition: viewFactor.H:118
globalIndex.H
Foam::radiation::viewFactor::viewFactor
viewFactor(const viewFactor &)
Disallow default bitwise copy construct.
Foam::radiation::viewFactor::Rp
virtual tmp< volScalarField > Rp() const
Source term component (for power of T^4)
Definition: viewFactor.C:680
Foam::radiation::viewFactor::~viewFactor
virtual ~viewFactor()
Destructor.
Definition: viewFactor.C:345
Foam::radiation::viewFactor::solarLoad_
autoPtr< solarLoad > solarLoad_
Solar load radiation model.
Definition: viewFactor.H:121
Foam::radiation::viewFactor::initialise
void initialise()
Initialise.
Definition: viewFactor.C:52
viewFactorI.H
Foam::singleCellFvMesh
fvMesh as subset of other mesh. Consists of one cell and all original bounday faces....
Definition: singleCellFvMesh.H:53
Foam::label
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
Foam::radiation::viewFactor::selectedPatches_
labelList selectedPatches_
Selected patches.
Definition: viewFactor.H:100
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
Foam::radiation::viewFactor
View factor radiation model. The system solved is: C q = b where: Cij = deltaij/Ej - (1/Ej - 1)Fij q ...
Definition: viewFactor.H:67
Foam::radiation::viewFactor::coarseMesh_
singleCellFvMesh coarseMesh_
Coarse mesh.
Definition: viewFactor.H:88
singleCellFvMesh.H
radiation
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::radiation::viewFactor::calculate
void calculate()
Solve system of equation(s)
Definition: viewFactor.C:387
Foam::radiation::viewFactor::nLocalCoarseFaces_
label nLocalCoarseFaces_
Total local coarse faces.
Definition: viewFactor.H:106
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
Foam::radiation::viewFactor::pivotIndices_
labelList pivotIndices_
Pivot Indices for LU decomposition.
Definition: viewFactor.H:115
Foam::radiation::viewFactor::operator=
void operator=(const viewFactor &)
Disallow default bitwise assignment.
Foam::radiation::viewFactor::map_
autoPtr< IOmapDistribute > map_
Map distributed.
Definition: viewFactor.H:85
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
Foam::SquareMatrix< scalar >
IOmapDistribute.H
Foam::radiation::viewFactor::totalNCoarseFaces_
label totalNCoarseFaces_
Total global coarse faces.
Definition: viewFactor.H:103
Foam::radiation::viewFactor::iterCounter_
label iterCounter_
Iterations Counter.
Definition: viewFactor.H:112
scalarListIOList.H
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::radiation::radiationModel
Top level model for radiation modelling.
Definition: radiationModel.H:73
Foam::radiation::viewFactor::Fmatrix_
autoPtr< scalarSquareMatrix > Fmatrix_
View factor matrix.
Definition: viewFactor.H:94
Foam::IOList< labelList >
Foam::radiation::viewFactor::viewFactorWalls
static const word viewFactorWalls
Static name for view factor walls.
Definition: viewFactor.H:76
scalarMatrices.H
Foam::radiation::viewFactor::Ru
virtual tmp< DimensionedField< scalar, volMesh > > Ru() const
Source term component (constant)
Definition: viewFactor.C:708
Foam::radiation::viewFactor::read
bool read()
Read radiation properties dictionary.
Definition: viewFactor.C:351
Foam::radiation::viewFactor::constEmissivity_
bool constEmissivity_
Constant emissivity.
Definition: viewFactor.H:109
Foam::radiation::viewFactor::finalAgglom_
labelListIOList finalAgglom_
Agglomeration List.
Definition: viewFactor.H:82
Foam::radiation::viewFactor::Qr_
volScalarField Qr_
Net radiative heat flux [W/m2].
Definition: viewFactor.H:91
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
radiationModel.H