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-2022 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"
28 #include "basicSpecieMixture.H"
29 #include "fvMatrix.H"
31 
32 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  namespace fv
37  {
38  defineTypeNameAndDebug(surfaceFilm, 0);
39 
41  (
42  fvModel,
43  surfaceFilm,
44  dictionary
45  );
46  }
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
53 (
54  const word& sourceName,
55  const word& modelType,
56  const dictionary& dict,
57  const fvMesh& mesh
58 )
59 :
60  fvModel(sourceName, modelType, dict, mesh),
61  surfaceFilm_
62  (
63  regionModels::surfaceFilmModels::thermoSingleLayer::typeName,
64  mesh,
66  "surfaceFilm"
67  ),
68  fieldNames_
69  (
70  {
72  (
73  IOobject::groupName("rho", surfaceFilm_.phaseName())
74  )
75  ? IOobject::groupName("rho", surfaceFilm_.phaseName())
76  : surfaceFilm_.primaryThermo().rho()().name(),
77  surfaceFilm_.UPrimary().name(),
78  surfaceFilm_.primaryThermo().he().name()
79  }
80  ),
81  curTimeIndex_(-1)
82 {
83  if (isA<basicSpecieMixture>(surfaceFilm_.primaryThermo()))
84  {
86  refCast<const basicSpecieMixture>(surfaceFilm_.primaryThermo());
87 
88  const PtrList<volScalarField>& Y = composition.Y();
89 
90  forAll(Y, i)
91  {
92  if (composition.solve(i))
93  {
94  fieldNames_.append(Y[i].name());
95  }
96  }
97  }
98 }
99 
100 
101 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
102 
104 {
105  return fieldNames_;
106 }
107 
108 
110 {
111  return surfaceFilm_.maxDeltaT();
112 }
113 
114 
116 {
117  if (curTimeIndex_ == mesh().time().timeIndex())
118  {
119  return;
120  }
121 
122  surfaceFilm_.evolve();
123 
124  curTimeIndex_ = mesh().time().timeIndex();
125 }
126 
127 
129 (
130  fvMatrix<scalar>& eqn,
131  const word& fieldName
132 ) const
133 {
134  if (debug)
135  {
136  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
137  }
138 
139  if (fieldName == fieldNames_[0])
140  {
141  eqn += surfaceFilm_.Srho();
142  }
143  else
144  {
146  << "Support for field " << fieldName << " is not implemented"
147  << exit(FatalError);
148  }
149 }
150 
151 
153 (
154  const volScalarField& rho,
155  fvMatrix<scalar>& eqn,
156  const word& fieldName
157 ) const
158 {
159  if (debug)
160  {
161  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
162  }
163 
164  if (fieldName == fieldNames_[0])
165  {
166  eqn += surfaceFilm_.Srho();
167  }
168  else if (fieldName == fieldNames_[2])
169  {
170  eqn += surfaceFilm_.Sh();
171  }
172  else if
173  (
174  isA<basicSpecieMixture>(surfaceFilm_.primaryThermo())
175  && refCast<const basicSpecieMixture>(surfaceFilm_.primaryThermo()).contains
176  (
177  eqn.psi().name()
178  )
179  )
180  {
181  eqn += surfaceFilm_.SYi
182  (
183  refCast<const basicSpecieMixture>(surfaceFilm_.primaryThermo())
184  .index(eqn.psi())
185  );
186  }
187  else
188  {
190  << "Support for field " << fieldName << " is not implemented"
191  << exit(FatalError);
192  }
193 }
194 
195 
197 (
198  const volScalarField& rho,
199  fvMatrix<vector>& eqn,
200  const word& fieldName
201 ) const
202 {
203  if (debug)
204  {
205  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
206  }
207 
208  if (fieldName == fieldNames_[1])
209  {
210  eqn += surfaceFilm_.SU();
211  }
212  else
213  {
215  << "Support for field " << fieldName << " is not implemented"
216  << exit(FatalError);
217  }
218 }
219 
220 
222 {}
223 
224 
226 {}
227 
228 
230 {}
231 
232 
234 {
235  return true;
236 }
237 
238 
239 // ************************************************************************* //
virtual void correct()
Solve the Lagrangian surfaceFilm and update the sources.
Definition: surfaceFilm.C:115
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
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
const word & name() const
Return name.
Definition: IOobject.H:315
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:306
basicSpecieMixture & composition
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:290
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
bool foundObject(const word &name) const
Is the named Type found?
Finite volume model abstract base class.
Definition: fvModel.H:57
virtual wordList addSupFields() const
Return the list of fields for which the option adds source term.
Definition: surfaceFilm.C:103
Specialisation of basicMixture for a mixture consisting of a number for molecular species...
fvMesh & mesh
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.
virtual scalar maxDeltaT() const
Return the maximum time-step for stable operation.
Definition: surfaceFilm.C:109
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
Definition: surfaceFilm.C:221
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
static word groupName(Name name, const word &group)
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
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: surfaceFilm.C:225
virtual bool movePoints()
Update for mesh motion.
Definition: surfaceFilm.C:233
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
addToRunTimeSelectionTable(fvConstraint, fixedTemperatureConstraint, dictionary)
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: surfaceFilm.C:229
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:95
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:53
messageStream Info
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:50
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:129
Namespace for OpenFOAM.