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-2023 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 "fvModels.H"
27 
28 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
29 
30 template<class Type, class ... AlphaRhoFieldTypes>
31 Foam::tmp<Foam::fvMatrix<Type>> Foam::fvModels::sourceTerm
32 (
33  const VolField<Type>& eqnField,
34  const dimensionSet& ds,
35  const AlphaRhoFieldTypes& ... alphaRhoFields
36 ) const
37 {
38  checkApplied();
39 
40  tmp<fvMatrix<Type>> tmtx
41  (
42  new fvMatrix<Type>
43  (
44  eqnField,
45  fvModel::sourceDims(ds, alphaRhoFields ...)
46  )
47  );
48  fvMatrix<Type>& mtx = tmtx.ref();
49 
50  const PtrListDictionary<fvModel>& modelList(*this);
51 
52  const word& fieldName = fvModel::fieldName(alphaRhoFields ...);
53 
54  forAll(modelList, i)
55  {
56  const fvModel& model = modelList[i];
57 
58  if (model.addsSupToField(fieldName))
59  {
60  addSupFields_[i].insert(fieldName);
61 
62  if (debug)
63  {
64  Info<< "Applying model " << model.name() << " to field "
65  << fieldName << endl;
66  }
67 
68  model.addSup(alphaRhoFields ..., mtx);
69  }
70  }
71 
72  return tmtx;
73 }
74 
75 
76 // * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * * //
77 
78 template<class Type>
80 (
81  const VolField<Type>& eqnField
82 ) const
83 {
84  return sourceTerm(eqnField, dimVolume/dimTime);
85 }
86 
87 
88 template<class Type>
90 (
91  const VolField<Type>& field
92 ) const
93 {
94  return sourceTerm(field, dimVolume/dimTime, field);
95 }
96 
97 
98 template<class Type>
100 (
101  const VolField<Type>& field,
102  const VolField<Type>& eqnField
103 ) const
104 {
105  return sourceTerm(eqnField, dimVolume/dimTime, field);
106 }
107 
108 
109 template<class Type>
111 (
112  const volScalarField& rho,
113  const VolField<Type>& field
114 ) const
115 {
116  return sourceTerm(field, dimVolume/dimTime, rho, field);
117 }
118 
119 
120 template<class Type>
122 (
123  const volScalarField& rho,
124  const VolField<Type>& field,
125  const VolField<Type>& eqnField
126 ) const
127 {
128  return sourceTerm(eqnField, dimVolume/dimTime, rho, field);
129 }
130 
131 
132 template<class Type>
134 (
135  const volScalarField& alpha,
136  const volScalarField& rho,
137  const VolField<Type>& field
138 ) const
139 {
140  return sourceTerm(field, dimVolume/dimTime, alpha, rho, field);
141 }
142 
143 
144 template<class Type>
146 (
147  const volScalarField& alpha,
148  const volScalarField& rho,
149  const VolField<Type>& field,
150  const VolField<Type>& eqnField
151 ) const
152 {
153  return sourceTerm(eqnField, dimVolume/dimTime, alpha, rho, field);
154 }
155 
156 
157 template<class Type>
159 (
160  const geometricOneField& alpha,
161  const geometricOneField& rho,
162  const VolField<Type>& field
163 ) const
164 {
165  return this->source(field);
166 }
167 
168 
169 template<class Type>
171 (
172  const volScalarField& alpha,
173  const geometricOneField& rho,
174  const VolField<Type>& field
175 ) const
176 {
177  return this->source(alpha, field);
178 }
179 
180 
181 template<class Type>
183 (
184  const geometricOneField& alpha,
185  const volScalarField& rho,
186  const VolField<Type>& field
187 ) const
188 {
189  return this->source(rho, field);
190 }
191 
192 
193 template<class Type>
195 (
196  const VolField<Type>& field
197 ) const
198 {
199  return sourceTerm(field, dimVolume/sqr(dimTime), field);
200 }
201 
202 
203 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Generic GeometricField class.
Dimension set for the base types.
Definition: dimensionSet.H:125
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Finite volume model abstract base class.
Definition: fvModel.H:59
static dimensionSet sourceDims(const dimensionSet &ds, const AlphaRhoFieldType &alphaRhoField, const AlphaRhoFieldTypes &... alphaRhoFields)
Return the dimensions of the matrix of a source term.
virtual void addSup(fvMatrix< scalar > &eqn) const
Add a source term to a field-less proxy equation.
Definition: fvModel.C:178
const word & name() const
Return const access to the source name.
Definition: fvModelI.H:47
virtual bool addsSupToField(const word &fieldName) const
Return true if the fvModel adds a source term to the given.
Definition: fvModel.C:166
static const word & fieldName()
Return the name of the field associated with a source term. Special.
Definition: fvModelI.H:39
tmp< fvMatrix< Type > > d2dt2(const VolField< Type > &field) const
Return source for an equation with a second time derivative.
tmp< fvMatrix< Type > > source(const VolField< Type > &field) const
Return source for an equation.
tmp< fvMatrix< Type > > sourceProxy(const VolField< Type > &eqnField) const
Return source for an equation.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
A class for handling words, derived from string.
Definition: word.H:62
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
messageStream Info
const dimensionSet dimTime
const dimensionSet dimVolume