coalCloudListI.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) 2012-2016 OpenFOAM Foundation
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
26\*---------------------------------------------------------------------------*/
27
28#include "fvMatrices.H"
29#include "volFields.H"
30#include "DimensionedField.H"
31
34{
36 (
38 (
40 (
41 "UTransEff",
42 mesh_.time().timeName(),
43 mesh_,
46 ),
47 mesh_,
49 )
50 );
51
53
54 forAll(*this, i)
55 {
56 fld += operator[](i).UTrans();
57 }
58
59 return tfld;
60}
61
62
64(
66) const
67{
69 fvVectorMatrix& fvm = tfvm.ref();
70
71 forAll(*this, i)
72 {
73 fvm += operator[](i).SU(U);
74 }
75
76 return tfvm;
77}
78
79
82{
84 (
86 (
88 (
89 "hsTransEff",
90 mesh_.time().timeName(),
91 mesh_,
94 ),
95 mesh_,
97 )
98 );
99
101
102 forAll(*this, i)
103 {
104 fld += operator[](i).hsTrans();
105 }
106
107 return tfld;
108}
109
110
112(
114) const
115{
117 fvScalarMatrix& fvm = tfvm.ref();
118
119 forAll(*this, i)
120 {
121 fvm += operator[](i).Sh(hs);
122 }
123
124 return tfvm;
125}
126
127
129(
130 const label ii,
132) const
133{
135 fvScalarMatrix& fvm = tfvm.ref();
136
137 forAll(*this, i)
138 {
139 fvm += operator[](i).SYi(ii, Yi);
140 }
141
142 return tfvm;
143}
144
145
148{
150 (
152 (
154 (
155 "rhoTransEff",
156 mesh_.time().timeName(),
157 mesh_,
160 ),
161 mesh_,
163 )
164 );
165
167
168 forAll(*this, i)
169 {
170 forAll(operator[](i).rhoTrans(), j)
171 {
172 fld += operator[](i).rhoTrans()[j];
173 }
174 }
175
176 return tfld;
177}
178
179
180
181
184{
186 (
188 (
190 (
191 "rhoTransEff",
192 mesh_.time().timeName(),
193 mesh_,
196 ),
197 mesh_,
199 )
200 );
201
203
204 forAll(*this, i)
205 {
206 fld += operator[](i).Srho();
207 }
208
209 return tfld;
210}
211
212
215(
216 const label i
217) const
218{
220 (
222 (
224 (
225 "rhoTransEff",
226 mesh_.time().timeName(),
227 mesh_,
230 ),
231 mesh_,
233 )
234 );
235
237
238 forAll(*this, j)
239 {
240 fld += operator[](j).Srho(i);
241 }
242
243 return tfld;
244}
245
246
248(
250) const
251{
253 fvScalarMatrix& fvm = tfvm.ref();
254
255 forAll(*this, i)
256 {
257 fvm += operator[](i).Srho(rho);
258 }
259
260 return tfvm;
261}
262
263
264// ************************************************************************* //
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))
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
scalar Sh() const
Sherwood number.
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
const T & operator[](const label i) const
Return const reference to the element.
Definition: UPtrListI.H:234
tmp< fvVectorMatrix > SU(volVectorField &U) const
Return tmp momentum source term.
tmp< volVectorField::Internal > UTrans() const
Return const reference to momentum source.
tmp< fvScalarMatrix > SYi(const label i, volScalarField &Yi) const
Return mass source term for specie i - specie eqn.
tmp< volScalarField::Internal > Srho() const
Return tmp total mass source for carrier phase.
tmp< volScalarField::Internal > hsTrans() const
Sensible enthalpy transfer [J/kg].
tmp< volScalarField::Internal > rhoTrans() const
Return total mass transfer [kg/m3].
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
U
Definition: pEqn.H:72
A special matrix type and solver, designed for finite volume solutions of scalar equations.
const dimensionSet dimEnergy
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:45
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
const dimensionSet dimVelocity
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:46
const dimensionSet dimForce
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
const dimensionSet dimDensity
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333