VoFSolidificationMelting.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) 2017-2024 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 
28 #include "fvcDdt.H"
31 
32 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  namespace fv
37  {
39 
41  (
42  fvModel,
45  );
47  (
48  fvModel,
50  dictionary,
51  VoFSolidificationMeltingSource,
52  "VoFSolidificationMeltingSource"
53  );
54  }
55 }
56 
57 
58 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
59 
60 void Foam::fv::VoFSolidificationMelting::readCoeffs()
61 {
62  alphaSolidT_.reset
63  (
65  (
66  "alphaSolidT",
69  coeffs()
70  ).ptr()
71  );
73  relax_ = coeffs().lookupOrDefault<scalar>("relax", dimless, 0.9);
74  Cu_ = coeffs().lookupOrDefault<scalar>("Cu", dimless/dimTime, 100000);
75  q_ = coeffs().lookupOrDefault<scalar>("q", dimless, 0.001);
76 }
77 
78 
79 Foam::word Foam::fv::VoFSolidificationMelting::alphaSolidName() const
80 {
81  const compressibleTwoPhaseVoFMixture& thermo
82  (
83  mesh().lookupObject<compressibleTwoPhaseVoFMixture>
84  (
85  "phaseProperties"
86  )
87  );
88 
89  const volScalarField& alphaVoF = thermo.alpha1();
90 
91  return IOobject::groupName(alphaVoF.name(), "solid");
92 }
93 
94 
95 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
96 
98 (
99  const word& name,
100  const word& modelType,
101  const fvMesh& mesh,
102  const dictionary& dict
103 )
104 :
105  fvModel(name, modelType, mesh, dict),
106  set_(mesh, coeffs()),
107  alphaSolidT_(),
108  L_("L", dimEnergy/dimMass, NaN),
109  relax_(NaN),
110  Cu_(NaN),
111  q_(NaN),
112 
113  thermo_
114  (
115  mesh().lookupObject<compressibleTwoPhaseVoFMixture>
116  (
117  "phaseProperties"
118  )
119  ),
120 
121  alphaSolid_
122  (
123  IOobject
124  (
125  alphaSolidName(),
126  mesh.time().name(),
127  mesh,
128  IOobject::READ_IF_PRESENT,
129  IOobject::AUTO_WRITE
130  ),
131  mesh,
133  zeroGradientFvPatchScalarField::typeName
134  )
135 {
136  readCoeffs();
137 }
138 
139 
140 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
141 
143 {
144  return wordList({"U", thermo_.thermo1().he().name()});
145 }
146 
147 
149 (
150  const volScalarField& alpha,
151  const volScalarField& rho,
152  const volScalarField& he,
153  fvMatrix<scalar>& eqn
154 ) const
155 {
156  if (debug)
157  {
158  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
159  }
160 
161  eqn += L_*(fvc::ddt(rho, alphaSolid_));
162 }
163 
164 
166 (
167  const volScalarField& rho,
168  const volVectorField& U,
169  fvMatrix<vector>& eqn
170 ) const
171 {
172  if (debug)
173  {
174  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
175  }
176 
177  scalarField& Sp = eqn.diag();
178  const scalarField& V = mesh().V();
179 
180  const labelUList cells = set_.cells();
181 
182  forAll(cells, i)
183  {
184  const label celli = cells[i];
185  const scalar Vc = V[celli];
186  const scalar alphaFluid = 1 - alphaSolid_[celli];
187 
188  const scalar S = Cu_*sqr(1 - alphaFluid)/(pow3(alphaFluid) + q_);
189 
190  Sp[celli] -= Vc*rho[celli]*S;
191  }
192 }
193 
194 
196 {
197  if (debug)
198  {
199  Info<< type() << ": " << name()
200  << " - updating solid phase fraction" << endl;
201  }
202 
203  alphaSolid_.oldTime();
204 
206  (
207  mesh().lookupObject<compressibleTwoPhaseVoFMixture>
208  (
209  "phaseProperties"
210  )
211  );
212 
213  const volScalarField& TVoF = thermo.thermo1().T();
214  const volScalarField& alphaVoF = thermo.alpha1();
215 
216  const labelUList cells = set_.cells();
217 
218  forAll(cells, i)
219  {
220  const label celli = cells[i];
221 
222  alphaSolid_[celli] = min
223  (
224  relax_*alphaVoF[celli]*alphaSolidT_->value(TVoF[celli])
225  + (1 - relax_)*alphaSolid_[celli],
226  alphaVoF[celli]
227  );
228  }
229 
230  alphaSolid_.correctBoundaryConditions();
231 }
232 
233 
235 (
236  const polyTopoChangeMap& map
237 )
238 {
239  set_.topoChange(map);
240 }
241 
242 
244 {
245  set_.mapMesh(map);
246 }
247 
248 
250 (
251  const polyDistributionMap& map
252 )
253 {
254  set_.distribute(map);
255 }
256 
257 
259 {
260  set_.movePoints();
261  return true;
262 }
263 
264 
266 {
267  if (fvModel::read(dict))
268  {
269  set_.read(coeffs());
270  readCoeffs();
271  return true;
272  }
273  else
274  {
275  return false;
276  }
277 
278  return false;
279 }
280 
281 
282 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
static autoPtr< Function1< scalar > > New(const word &name, const Function1s::unitConversions &units, const dictionary &dict)
Select from dictionary.
Definition: Function1New.C:32
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
static word groupName(Name name, const word &group)
Class to represent a mixture of two rhoFluidThermo-based phases.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T, if not found return the given default.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
VolField< Type > & psi()
Definition: fvMatrix.H:289
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Finite volume model abstract base class.
Definition: fvModel.H:59
const dictionary & coeffs() const
Return dictionary.
Definition: fvModelI.H:59
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:199
Solidification and melting model for VoF simulations.
virtual bool movePoints()
Update for mesh motion.
virtual wordList addSupFields() const
Return the list of fields for which the option adds source term.
virtual void addSup(const volScalarField &alpha, const volScalarField &rho, const volScalarField &he, fvMatrix< scalar > &eqn) const
Add explicit contribution to phase energy equation.
virtual void correct()
Correct the model.
VoFSolidificationMelting(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual bool read(const dictionary &dict)
Read source dictionary.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
scalarField & diag()
Definition: lduMatrix.C:186
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for handling words, derived from string.
Definition: word.H:62
Calculate the first temporal derivative.
const cellShapeList & cells
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
addBackwardCompatibleToRunTimeSelectionTable(fvConstraint, fixedTemperature, dictionary, fixedTemperatureConstraint, "fixedTemperatureConstraint")
tmp< VolField< Type > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
tmp< VolField< Type > > Sp(const volScalarField &sp, const VolField< Type > &vf)
Definition: fvcSup.C:67
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const VolField< Type > &)
Namespace for OpenFOAM.
List< word > wordList
A List of words.
Definition: fileName.H:54
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const dimensionSet dimEnergy
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
dimensionedScalar pow3(const dimensionedScalar &ds)
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimless
messageStream Info
const dimensionSet dimTemperature
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
const dimensionSet dimMass
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
const unitConversion unitFraction
thermo he()
labelList fv(nPoints)
dictionary dict
fluidMulticomponentThermo & thermo
Definition: createFields.H:31