All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
VoFClouds.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 "VoFClouds.H"
28 #include "fvmSup.H"
31 
32 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  namespace fv
37  {
38  defineTypeNameAndDebug(VoFClouds, 0);
39 
41  (
42  fvModel,
43  VoFClouds,
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  phaseName_(dict.lookup("phase")),
62  carrierPhaseName_(dict.lookup("carrierPhase")),
63  thermo_
64  (
65  mesh.lookupObject<fluidThermo>
66  (
67  IOobject::groupName(physicalProperties::typeName, phaseName_)
68  )
69  ),
70  carrierThermo_
71  (
72  mesh.lookupObject<fluidThermo>
73  (
74  IOobject::groupName(physicalProperties::typeName, carrierPhaseName_)
75  )
76  ),
77  clouds_
78  (
79  carrierThermo_.rho(),
80  mesh.lookupObject<volVectorField>("U"),
81  mesh.lookupObject<uniformDimensionedVectorField>("g"),
82  carrierThermo_
83  ),
84  curTimeIndex_(-1)
85 {}
86 
87 
88 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
89 
91 {
92  return wordList({thermo_.rho()().name(), "U", "T"});
93 }
94 
95 
97 {
98  if (curTimeIndex_ == mesh().time().timeIndex())
99  {
100  return;
101  }
102 
103  clouds_.evolve();
104 
105  curTimeIndex_ = mesh().time().timeIndex();
106 }
107 
108 
110 (
111  const volScalarField& rho,
112  fvMatrix<scalar>& eqn,
113  const word& fieldName
114 ) const
115 {
116  if (debug)
117  {
118  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
119  }
120 
121  if (fieldName == thermo_.rho()().name())
122  {
123  eqn += clouds_.Srho();
124  }
125  else if (fieldName == "T")
126  {
127  const volScalarField::Internal Cv(thermo_.Cv());
128 
129  eqn +=
130  clouds_.Sh(eqn.psi())()/Cv
131  + clouds_.Srho()*(eqn.psi() - thermo_.he()/Cv);
132  }
133  else
134  {
136  << "Support for field " << fieldName << " is not implemented"
137  << exit(FatalError);
138  }
139 }
140 
141 
143 (
144  const volScalarField& rho,
145  fvMatrix<vector>& eqn,
146  const word& fieldName
147 ) const
148 {
149  if (debug)
150  {
151  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
152  }
153 
154  eqn += clouds_.SU(eqn.psi());
155 }
156 
157 
159 {
160  // Store the particle positions
161  clouds_.storeGlobalPositions();
162 }
163 
164 
165 void Foam::fv::VoFClouds::topoChange(const polyTopoChangeMap&)
166 {}
167 
168 
169 void Foam::fv::VoFClouds::mapMesh(const polyMeshMap& map)
170 {}
171 
172 
173 void Foam::fv::VoFClouds::distribute(const polyDistributionMap& map)
174 {
175  clouds_.distribute(map);
176 }
177 
178 
180 {
181  return true;
182 }
183 
184 
185 // ************************************************************************* //
scalar Cv(const scalar p, const scalar T) const
Definition: HtoEthermo.H:2
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
virtual bool movePoints()
Update for mesh motion.
virtual wordList addSupFields() const
Return the list of fields for which the option adds source term.
dictionary dict
defineTypeNameAndDebug(fixedTemperatureConstraint, 0)
const word & name() const
Return name.
Definition: IOobject.H:315
virtual void preUpdateMesh()
Prepare for mesh update.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
UniformDimensionedField< vector > uniformDimensionedVectorField
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
VoFClouds(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:59
fvMesh & mesh
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
Macros for easy insertion into run-time selection tables.
virtual void correct()
Solve the Lagrangian clouds and update the sources.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
stressControl lookup("compactNormalStress") >> compactNormalStress
labelList fv(nPoints)
virtual void addSup(const volScalarField &rho, fvMatrix< scalar > &eqn, const word &fieldName) const
Add explicit contribution to compressible enthalpy equation.
addToRunTimeSelectionTable(fvConstraint, fixedTemperatureConstraint, dictionary)
List< word > wordList
A List of words.
Definition: fileName.H:54
const Time & time() const
Return time.
Definition: IOobject.C:318
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
messageStream Info
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Calculate the matrix for implicit and explicit sources.
Namespace for OpenFOAM.