fvModelsTemplates.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) 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 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
27 
28 template<class Type, class ... AlphaRhoFieldTypes>
29 Foam::tmp<Foam::fvMatrix<Type>> Foam::fvModels::source
30 (
32  const word& fieldName,
33  const dimensionSet& ds,
34  const AlphaRhoFieldTypes& ... alphaRhos
35 ) const
36 {
37  checkApplied();
38 
39  tmp<fvMatrix<Type>> tmtx
40  (
41  new fvMatrix<Type>
42  (
43  field,
44  fvModel::sourceDims(field, ds, alphaRhos ...)
45  )
46  );
47  fvMatrix<Type>& mtx = tmtx.ref();
48 
49  const PtrListDictionary<fvModel>& modelList(*this);
50 
51  forAll(modelList, i)
52  {
53  const fvModel& model = modelList[i];
54 
55  if (model.addsSupToField(fieldName))
56  {
57  addSupFields_[i].insert(fieldName);
58 
59  if (debug)
60  {
61  Info<< "Applying model " << model.name() << " to field "
62  << fieldName << endl;
63  }
64 
65  model.addSup(alphaRhos ..., mtx, fieldName);
66  }
67  }
68 
69  return tmtx;
70 }
71 
72 
73 // * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * * //
74 
75 template<class Type>
76 Foam::tmp<Foam::fvMatrix<Type>> Foam::fvModels::source
77 (
79 ) const
80 {
81  return this->source(field, field.name());
82 }
83 
84 
85 template<class Type>
86 Foam::tmp<Foam::fvMatrix<Type>> Foam::fvModels::source
87 (
89  const word& fieldName
90 ) const
91 {
92  return source(field, fieldName, dimVolume/dimTime);
93 }
94 
95 
96 template<class Type>
97 Foam::tmp<Foam::fvMatrix<Type>> Foam::fvModels::source
98 (
99  const volScalarField& rho,
101 ) const
102 {
103  return this->source(rho, field, field.name());
104 }
105 
106 
107 template<class Type>
108 Foam::tmp<Foam::fvMatrix<Type>> Foam::fvModels::source
109 (
110  const volScalarField& rho,
112  const word& fieldName
113 ) const
114 {
115  return source(field, fieldName, dimVolume/dimTime, rho);
116 }
117 
118 
119 template<class Type>
120 Foam::tmp<Foam::fvMatrix<Type>> Foam::fvModels::source
121 (
122  const volScalarField& alpha,
123  const volScalarField& rho,
125 ) const
126 {
127  return this->source(alpha, rho, field, field.name());
128 }
129 
130 
131 template<class Type>
132 Foam::tmp<Foam::fvMatrix<Type>> Foam::fvModels::source
133 (
134  const volScalarField& alpha,
135  const volScalarField& rho,
137  const word& fieldName
138 ) const
139 {
140  return source(field, fieldName, dimVolume/dimTime, alpha, rho);
141 }
142 
143 
144 template<class Type>
145 Foam::tmp<Foam::fvMatrix<Type>> Foam::fvModels::source
146 (
147  const geometricOneField& alpha,
148  const geometricOneField& rho,
150 ) const
151 {
152  return this->source(field, field.name());
153 }
154 
155 
156 template<class Type>
157 Foam::tmp<Foam::fvMatrix<Type>> Foam::fvModels::source
158 (
159  const volScalarField& alpha,
160  const geometricOneField& rho,
162 ) const
163 {
165  (
166  IOobject
167  (
168  "one",
169  this->mesh_.time().timeName(),
170  this->mesh_,
171  IOobject::NO_READ,
172  IOobject::NO_WRITE,
173  false
174  ),
175  this->mesh_,
177  );
178 
179  return this->source(alpha, one, field, field.name());
180 }
181 
182 
183 template<class Type>
184 Foam::tmp<Foam::fvMatrix<Type>> Foam::fvModels::source
185 (
186  const geometricOneField& alpha,
187  const volScalarField& rho,
189 ) const
190 {
191  return this->source(rho, field, field.name());
192 }
193 
194 
195 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
196 
197 template<class Type>
199 (
201 ) const
202 {
203  return this->d2dt2(field, field.name());
204 }
205 
206 
207 template<class Type>
209 (
211  const word& fieldName
212 ) const
213 {
214  return source(field, fieldName, dimVolume/sqr(dimTime));
215 }
216 
217 
218 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
tmp< GeometricField< Type, fvPatchField, volMesh > > d2dt2(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcD2dt2.C:45
const word & name() const
Return name.
Definition: IOobject.H:315
const word & name() const
Return const access to the source name.
Definition: fvModelI.H:28
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Finite volume model abstract base class.
Definition: fvModel.H:57
Generic GeometricField class.
const dimensionSet dimless
void insert(const word &, T *)
Add at head of dictionary.
const dimensionSet dimTime
Dimension set for the base types.
Definition: dimensionSet.H:121
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
A class for handling words, derived from string.
Definition: word.H:59
virtual bool addsSupToField(const word &fieldName) const
Return true if the fvModel adds a source term to the given.
Definition: fvModel.C:143
tmp< fvMatrix< Type > > d2dt2(const GeometricField< Type, fvPatchField, volMesh > &field) const
Return source for an equation with a second time derivative.
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
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimVolume
messageStream Info
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
A class representing the concept of 1 (scalar(1)) used to avoid unnecessary manipulations for objects...
Definition: one.H:50