Loading...
Searching...
No Matches
porosityModelList.C
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-2014 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 "porosityModelList.H"
29#include "volFields.H"
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
33Foam::porosityModelList::porosityModelList
34(
35 const fvMesh& mesh,
36 const dictionary& dict
37)
38:
40 mesh_(mesh)
41{
43 active(true);
44}
45
46
47// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
48
49bool Foam::porosityModelList::active(const bool warn) const
50{
51 bool anyOk = false;
52 forAll(*this, i)
53 {
54 anyOk = anyOk || this->operator[](i).active();
55 }
56
57 if (warn && this->size() && !anyOk)
58 {
59 Info<< "No porosity models active" << endl;
60 }
61
62 return anyOk;
63}
64
65
67{
68 label count = 0;
69 for (const entry& dEntry : dict)
70 {
71 if (dEntry.isDict())
72 {
73 ++count;
74 }
75 }
76
77 this->resize(count);
78
79 count = 0;
80 for (const entry& dEntry : dict)
81 {
82 if (dEntry.isDict())
83 {
84 const word& name = dEntry.keyword();
85 const dictionary& modelDict = dEntry.dict();
86
87 this->set
88 (
89 count++,
91 );
92 }
93 }
94}
95
96
98{
99 bool allOk = true;
100 forAll(*this, i)
101 {
102 porosityModel& pm = this->operator[](i);
103 bool ok = pm.read(dict.subDict(pm.name()));
104 allOk = (allOk && ok);
105 }
106 return allOk;
107}
108
109
110bool Foam::porosityModelList::writeData(Ostream& os) const
111{
112 forAll(*this, i)
113 {
114 os << nl;
115 this->operator[](i).writeData(os);
116 }
117
118 return os.good();
119}
120
121
123(
125)
126{
127 forAll(*this, i)
128 {
129 this->operator[](i).addResistance(UEqn);
130 }
131}
132
133
135(
137 const volScalarField& rho,
138 const volScalarField& mu
139)
140{
141 forAll(*this, i)
142 {
143 this->operator[](i).addResistance(UEqn, rho, mu);
144 }
145}
146
147
149(
150 const fvVectorMatrix& UEqn,
151 volTensorField& AU,
152 bool correctAUprocBC
153)
154{
155 forAll(*this, i)
156 {
157 this->operator[](i).addResistance(UEqn, AU, correctAUprocBC);
158 }
159}
160
161
162// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
163
164Foam::Ostream& Foam::operator<<
165(
166 Ostream& os,
167 const porosityModelList& models
168)
169{
170 models.writeData(os);
171 return os;
172}
173
174
175// ************************************************************************* //
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
const porosityModel * set(const label i) const
Definition PtrList.H:171
constexpr PtrList() noexcept
Definition PtrListI.H:29
const porosityModel & operator[](const label i) const
Definition UPtrListI.H:289
label size() const noexcept
Definition UPtrListI.H:106
label count() const noexcept
Definition UPtrList.H:890
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A keyword and a list of tokens is an 'entry'.
Definition entry.H:66
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
List container for porosity models.
void addResistance(fvVectorMatrix &UEqn)
Add resistance.
void reset(const dictionary &dict)
Reset the source list.
const fvMesh & mesh_
Reference to the mesh database.
bool writeData(Ostream &os) const
Write data to Ostream.
bool read(const dictionary &dict)
Read dictionary.
bool active(const bool warn=false) const
Return active status.
Top level model for porosity models.
virtual bool read(const dictionary &dict)
Read porosity dictionary.
static autoPtr< porosityModel > New(const word &name, const fvMesh &mesh, const dictionary &dict, const wordRe &cellZoneName=wordRe::null)
Selector.
const word & name() const
Return const access to the porosity model name.
A class for handling words, derived from Foam::string.
Definition word.H:66
fvVectorMatrix & UEqn
Definition UEqn.H:13
patchWriters resize(patchIds.size())
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
void set(List< bool > &bools, const labelUList &locations)
Set the listed locations (assign 'true').
Definition BitOps.C:35
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Definition BitOps.H:73
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
GeometricField< tensor, fvPatchField, volMesh > volTensorField
fvMatrix< vector > fvVectorMatrix
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299