surfaceFilm.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 #include "surfaceFilm.H"
27 #include "basicSpecieMixture.H"
28 #include "fvMatrix.H"
30 
31 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  namespace fv
36  {
37  defineTypeNameAndDebug(surfaceFilm, 0);
38 
40  (
41  fvModel,
42  surfaceFilm,
43  dictionary
44  );
45  }
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
52 (
53  const word& sourceName,
54  const word& modelType,
55  const dictionary& dict,
56  const fvMesh& mesh
57 )
58 :
59  fvModel(sourceName, modelType, dict, mesh),
60  primaryThermo_
61  (
63  ),
64  surfaceFilm_
65  (
67  (
68  mesh,
70  )
71  ),
72  curTimeIndex_(-1)
73 {}
74 
75 
76 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
77 
79 {
80  wordList fieldNames({"rho", "U", primaryThermo_.he().name()});
81 
82  if (isA<basicSpecieMixture>(primaryThermo_))
83  {
85  refCast<const basicSpecieMixture>(primaryThermo_);
86 
87  const PtrList<volScalarField>& Y = composition.Y();
88 
89  forAll(Y, i)
90  {
91  if (composition.solve(i))
92  {
93  fieldNames.append(Y[i].name());
94  }
95  }
96  }
97 
98  return fieldNames;
99 }
100 
101 
103 {
104  if (curTimeIndex_ == mesh().time().timeIndex())
105  {
106  return;
107  }
108 
109  surfaceFilm_->evolve();
110 
111  curTimeIndex_ = mesh().time().timeIndex();
112 }
113 
114 
116 (
117  fvMatrix<scalar>& eqn,
118  const word& fieldName
119 ) const
120 {
121  if (debug)
122  {
123  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
124  }
125 
126  if (fieldName == "rho")
127  {
128  eqn += surfaceFilm_->Srho();
129  }
130  else
131  {
133  << "Support for field " << fieldName << " is not implemented"
134  << exit(FatalError);
135  }
136 }
137 
138 
140 (
141  const volScalarField& rho,
142  fvMatrix<scalar>& eqn,
143  const word& fieldName
144 ) const
145 {
146  if (debug)
147  {
148  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
149  }
150 
151  if (fieldName == "rho")
152  {
153  eqn += surfaceFilm_->Srho();
154  }
155  else if (fieldName == primaryThermo_.he().name())
156  {
157  eqn += surfaceFilm_->Sh();
158  }
159  else if
160  (
161  isA<basicSpecieMixture>(primaryThermo_)
162  && refCast<const basicSpecieMixture>(primaryThermo_).contains
163  (
164  eqn.psi().name()
165  )
166  )
167  {
168  eqn += surfaceFilm_->SYi
169  (
170  refCast<const basicSpecieMixture>(primaryThermo_).index(eqn.psi())
171  );
172  }
173  else
174  {
176  << "Support for field " << fieldName << " is not implemented"
177  << exit(FatalError);
178  }
179 }
180 
181 
183 (
184  const volScalarField& rho,
185  fvMatrix<vector>& eqn,
186  const word& fieldName
187 ) const
188 {
189  if (debug)
190  {
191  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
192  }
193 
194  if (fieldName == "U")
195  {
196  eqn += surfaceFilm_->SU();
197  }
198  else
199  {
201  << "Support for field " << fieldName << " is not implemented"
202  << exit(FatalError);
203  }
204 }
205 
206 
207 // ************************************************************************* //
virtual void correct()
Solve the Lagrangian surfaceFilm and update the sources.
Definition: surfaceFilm.C:102
bool solve(label speciei) const
Return true if the specie should be solved for.
defineTypeNameAndDebug(fixedTemperatureConstraint, 0)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const word & name() const
Return name.
Definition: IOobject.H:303
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
basicSpecieMixture & composition
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:282
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Finite volume model abstract base class.
Definition: fvModel.H:55
virtual wordList addSupFields() const
Return the list of fields for which the option adds source term.
Definition: surfaceFilm.C:78
Specialisation of basicMixture for a mixture consisting of a number for molecular species...
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Macros for easy insertion into run-time selection tables.
static List< word > fieldNames
Definition: globalFoam.H:46
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
Base-class for fluid thermodynamic properties.
Definition: fluidThermo.H:53
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
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
static const word dictName
Name of the thermophysical properties dictionary.
Definition: basicThermo.H:129
static autoPtr< surfaceFilmModel > New(const fvMesh &mesh, const dimensionedVector &g, const word &regionType="surfaceFilm")
Return a reference to the selected surface film model.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
addToRunTimeSelectionTable(fvConstraint, fixedTemperatureConstraint, dictionary)
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
PtrList< volScalarField > & Y
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
label timeIndex
Definition: getTimeIndex.H:4
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
surfaceFilm(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
Definition: surfaceFilm.C:52
messageStream Info
PtrList< volScalarField > & Y()
Return the mass-fraction fields.
virtual void addSup(fvMatrix< scalar > &eqn, const word &fieldName) const
Add source to continuity equation.
Definition: surfaceFilm.C:116
Namespace for OpenFOAM.