porosityModelList.C
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) 2012-2021 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 \*---------------------------------------------------------------------------*/
25 
26 #include "porosityModelList.H"
27 #include "volFields.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
32 (
33  const fvMesh& mesh,
34  const dictionary& dict
35 )
36 :
38  mesh_(mesh)
39 {
40  reset(dict);
41 }
42 
43 
44 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
45 
47 {}
48 
49 
50 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
51 
53 {
54  label count = 0;
55  forAllConstIter(dictionary, dict, iter)
56  {
57  if (iter().isDict())
58  {
59  count++;
60  }
61  }
62 
63  this->setSize(count);
64  label i = 0;
65  forAllConstIter(dictionary, dict, iter)
66  {
67  if (iter().isDict())
68  {
69  const word& name = iter().keyword();
70  const dictionary& modelDict = iter().dict();
71 
72  this->set
73  (
74  i++,
75  porosityModel::New(name, mesh_, modelDict)
76  );
77  }
78  }
79 }
80 
81 
83 {
84  bool allOk = true;
85  forAll(*this, i)
86  {
87  porosityModel& pm = this->operator[](i);
88  bool ok = pm.read(dict.subDict(pm.name()));
89  allOk = (allOk && ok);
90  }
91  return allOk;
92 }
93 
94 
96 {
97  forAll(*this, i)
98  {
99  os << nl;
100  this->operator[](i).writeData(os);
101  }
102 
103  return os.good();
104 }
105 
106 
108 (
109  fvVectorMatrix& UEqn
110 )
111 {
112  forAll(*this, i)
113  {
114  this->operator[](i).addResistance(UEqn);
115  }
116 }
117 
118 
120 (
121  const fvVectorMatrix& UEqn,
122  volTensorField& AU,
123  bool correctAUprocBC
124 )
125 {
126  forAll(*this, i)
127  {
128  this->operator[](i).addResistance(UEqn, AU, correctAUprocBC);
129  }
130 }
131 
132 
133 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
134 
135 Foam::Ostream& Foam::operator<<
136 (
137  Ostream& os,
138  const porosityModelList& models
139 )
140 {
141  models.writeData(os);
142  return os;
143 }
144 
145 
146 // ************************************************************************* //
bool writeData(Ostream &os) const
Write data to Ostream.
bool read(const dictionary &dict)
Read dictionary.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
const fvMesh & mesh_
Reference to the mesh database.
const T & operator[](const label) const
Return element const reference.
Definition: UPtrListI.H:96
porosityModelList(const fvMesh &mesh, const dictionary &dict)
label count(const ListType &l, typename ListType::const_reference x)
Count the number of occurrences of a value in a list.
~porosityModelList()
Destructor.
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:330
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:1002
A class for handling words, derived from string.
Definition: word.H:59
void reset(const dictionary &dict)
Reset the source list.
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:72
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
static const char nl
Definition: Ostream.H:260
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
const word & name() const
Return const access to the porosity model name.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
List container for porosity models.
static autoPtr< porosityModel > New(const word &name, const fvMesh &mesh, const dictionary &dict, const word &cellZoneName=word::null)
Selector.
void addResistance(fvVectorMatrix &UEqn)
Add resistance.
virtual bool read(const dictionary &dict)
Read porosity dictionary.
Top level model for porosity models.
Definition: porosityModel.H:54