solidificationMelting.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) 2014-2025 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 "solidificationMelting.H"
27 #include "fvcDdt.H"
28 #include "fvMatrices.H"
29 #include "basicThermo.H"
34 #include "geometricOneField.H"
35 
36 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace fv
41 {
43 
46  (
47  fvModel,
49  dictionary,
50  solidificationMeltingSource,
51  "solidificationMeltingSource"
52  );
53 }
54 }
55 
56 
59 {
60  "thermo",
61  "lookup"
62 };
63 
64 
65 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
66 
67 void Foam::fv::solidificationMelting::readCoeffs(const dictionary& dict)
68 {
69  Tsol_ = dict.lookup<scalar>("Tsol");
70  Tliq_ = dict.lookupOrDefault<scalar>("Tliq", Tsol_);
71  alpha1e_ = dict.lookupOrDefault<scalar>("alpha1e", 0.0);
72  L_ = dict.lookup<scalar>("L");
73 
74  relax_ = dict.lookupOrDefault("relax", 0.9);
75 
76  mode_ = thermoModeTypeNames_.read(dict.lookup("thermoMode"));
77 
78  rhoRef_ = dict.lookup<scalar>("rhoRef");
79  TName_ = dict.lookupOrDefault<word>("T", "T");
80  CpName_ = dict.lookupOrDefault<word>("Cp", "Cp");
81  UName_ = dict.lookupOrDefault<word>("U", "U");
82  phiName_ = dict.lookupOrDefault<word>("phi", "phi");
83 
84  Cu_ = dict.lookupOrDefault<scalar>("Cu", 100000);
85  q_ = dict.lookupOrDefault("q", 0.001);
86 
87  beta_ = dict.lookup<scalar>("beta");
88 
89  if (mode_ == thermoMode::lookup)
90  {
91  CpRef_ = dict.lookup<scalar>("CpRef");
92  }
93 
95  {
96  g_ = dict.lookup("g");
97  }
98 }
99 
100 
102 Foam::fv::solidificationMelting::Cp() const
103 {
104  switch (mode_)
105  {
106  case thermoMode::thermo:
107  {
108  const basicThermo& thermo =
109  mesh().lookupObject<basicThermo>(physicalProperties::typeName);
110 
111  return thermo.Cp();
112  break;
113  }
114  case thermoMode::lookup:
115  {
116  if (CpName_ == "CpRef")
117  {
118  return volScalarField::New
119  (
120  name() + ":Cp",
121  mesh(),
123  (
125  CpRef_
126  ),
127  extrapolatedCalculatedFvPatchScalarField::typeName
128  );
129  }
130  else
131  {
132  return mesh().lookupObject<volScalarField>(CpName_);
133  }
134 
135  break;
136  }
137  default:
138  {
140  << "Unhandled thermo mode: " << thermoModeTypeNames_[mode_]
141  << abort(FatalError);
142  }
143  }
144 
145  return tmp<volScalarField>(nullptr);
146 }
147 
148 
149 Foam::vector Foam::fv::solidificationMelting::g() const
150 {
151  if (mesh().foundObject<uniformDimensionedVectorField>("g"))
152  {
153  const uniformDimensionedVectorField& value =
155  return value.value();
156  }
157  else
158  {
159  return g_;
160  }
161 }
162 
163 
164 void Foam::fv::solidificationMelting::update
165 (
166  const volScalarField& Cp
167 ) const
168 {
169  if (curTimeIndex_ == mesh().time().timeIndex())
170  {
171  return;
172  }
173 
174  if (debug)
175  {
176  Info<< type() << ": " << name()
177  << " - updating phase indicator" << endl;
178  }
179 
180  // update old time alpha1 field
181  alpha1_.oldTime();
182 
183  const volScalarField& T = mesh().lookupObject<volScalarField>(TName_);
184 
185  const labelList& cells = zone_.zone();
186 
187  forAll(cells, i)
188  {
189  const label celli = cells[i];
190 
191  const scalar Tc = T[celli];
192  const scalar Cpc = Cp[celli];
193  const scalar alpha1New =
194  alpha1_[celli]
195  + relax_*Cpc
196  *(
197  Tc
198  - max
199  (
200  Tsol_,
201  Tsol_
202  + (Tliq_ - Tsol_)*(alpha1_[celli] - alpha1e_)/(1 - alpha1e_)
203  )
204  )/L_;
205 
206  alpha1_[celli] = max(0, min(alpha1New, 1));
207  deltaT_[i] =
208  Tc
209  - max
210  (
211  Tsol_,
212  Tsol_
213  + (Tliq_ - Tsol_)*(alpha1_[celli] - alpha1e_)/(1 - alpha1e_)
214  );
215  }
216 
217  alpha1_.correctBoundaryConditions();
218 
219  curTimeIndex_ = mesh().time().timeIndex();
220 }
221 
222 
223 template<class RhoFieldType>
224 void Foam::fv::solidificationMelting::apply
225 (
226  const RhoFieldType& rho,
227  fvMatrix<scalar>& eqn
228 ) const
229 {
230  if (debug)
231  {
232  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
233  }
234 
235  const volScalarField Cp(this->Cp());
236 
237  update(Cp);
238 
239  dimensionedScalar L("L", dimEnergy/dimMass, L_);
240 
241  // Contributions added to rhs of solver equation
242  if (eqn.psi().dimensions() == dimTemperature)
243  {
244  eqn -= L/Cp*(fvc::ddt(rho, alpha1_));
245  }
246  else
247  {
248  eqn -= L*(fvc::ddt(rho, alpha1_));
249  }
250 }
251 
252 
253 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
254 
256 (
257  const word& name,
258  const word& modelType,
259  const fvMesh& mesh,
260  const dictionary& dict
261 )
262 :
263  fvModel(name, modelType, mesh, dict),
264  zone_(mesh, coeffs(dict)),
265  Tsol_(NaN),
266  Tliq_(NaN),
267  alpha1e_(NaN),
268  L_(NaN),
269  relax_(NaN),
270  mode_(thermoMode::thermo),
271  rhoRef_(NaN),
272  TName_(word::null),
273  CpName_(word::null),
274  UName_(word::null),
275  phiName_(word::null),
276  Cu_(NaN),
277  q_(NaN),
278  beta_(NaN),
279  alpha1_
280  (
281  IOobject
282  (
283  this->name() + ":alpha1",
284  mesh.time().name(),
285  mesh,
286  IOobject::READ_IF_PRESENT,
287  IOobject::AUTO_WRITE
288  ),
289  mesh,
291  zeroGradientFvPatchScalarField::typeName
292  ),
293  curTimeIndex_(-1),
294  deltaT_(zone_.nCells(), 0)
295 {
296  readCoeffs(coeffs(dict));
297 }
298 
299 
300 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
301 
303 {
304  switch (mode_)
305  {
306  case thermoMode::thermo:
307  {
308  const basicThermo& thermo =
309  mesh().lookupObject<basicThermo>(physicalProperties::typeName);
310 
311  return wordList({UName_, thermo.he().name()});
312  }
313  case thermoMode::lookup:
314  {
315  return wordList({UName_, TName_});
316  }
317  }
318 
319  return wordList::null();
320 }
321 
322 
324 (
325  const volScalarField& he,
326  fvMatrix<scalar>& eqn
327 ) const
328 {
329  apply(geometricOneField(), eqn);
330 }
331 
332 
334 (
335  const volScalarField& rho,
336  const volScalarField& he,
337  fvMatrix<scalar>& eqn
338 ) const
339 {
340  apply(rho, eqn);
341 }
342 
343 
345 (
346  const volVectorField& U,
347  fvMatrix<vector>& eqn
348 ) const
349 {
350  if (debug)
351  {
352  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
353  }
354 
355  const volScalarField Cp(this->Cp());
356 
357  update(Cp);
358 
359  const vector g = this->g();
360 
361  scalarField& Sp = eqn.diag();
362  vectorField& Su = eqn.source();
363  const scalarField& V = mesh().V();
364 
365  const labelList& cells = zone_.zone();
366 
367  forAll(cells, i)
368  {
369  const label celli = cells[i];
370 
371  const scalar Vc = V[celli];
372  const scalar alpha1c = alpha1_[celli];
373 
374  const scalar S = -Cu_*sqr(1.0 - alpha1c)/(pow3(alpha1c) + q_);
375  const vector Sb = rhoRef_*g*beta_*deltaT_[i];
376 
377  Sp[celli] += Vc*S;
378  Su[celli] += Vc*Sb;
379  }
380 }
381 
382 
384 (
385  const volScalarField& rho,
386  const volVectorField& U,
387  fvMatrix<vector>& eqn
388 ) const
389 {
390  addSup(U, eqn);
391 }
392 
393 
395 {
396  zone_.movePoints();
397  return true;
398 }
399 
400 
402 (
403  const polyTopoChangeMap& map
404 )
405 {
406  zone_.topoChange(map);
407 }
408 
409 
411 {
412  zone_.mapMesh(map);
413 }
414 
415 
417 (
418  const polyDistributionMap& map
419 )
420 {
421  zone_.distribute(map);
422 }
423 
424 
426 {
427  if (fvModel::read(dict))
428  {
429  zone_.read(coeffs(dict));
430  readCoeffs(coeffs(dict));
431  return true;
432  }
433  else
434  {
435  return false;
436  }
437 
438  return false;
439 }
440 
441 
442 // ************************************************************************* //
scalar Cp(const scalar p, const scalar T) const
Definition: EtoHthermo.H:2
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
static const List< word > & null()
Return a null List.
Definition: ListI.H:118
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:55
Enum read(Istream &) const
Read a word from Istream and return the corresponding.
Definition: NamedEnum.C:54
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:28
Type & value()
Return a reference to the value.
Base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:78
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Field< Type > & source()
Definition: fvMatrix.H:307
VolField< Type > & psi()
Definition: fvMatrix.H:289
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:420
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
Finite volume model abstract base class.
Definition: fvModel.H:60
static const dictionary & coeffs(const word &modelType, const dictionary &)
Return the coefficients sub-dictionary for a given model type.
Definition: fvModelI.H:31
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:200
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:69
This source is designed to model the effect of solidification and melting processes,...
virtual bool movePoints()
Update for mesh motion.
virtual wordList addSupFields() const
Return the list of fields for which the fvModel adds source term.
virtual void addSup(const volScalarField &he, fvMatrix< scalar > &eqn) const
Add explicit contribution to enthalpy equation.
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.
static const NamedEnum< thermoMode, 2 > thermoModeTypeNames_
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
solidificationMelting(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from explicit source name and mesh.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
scalarField & diag()
Definition: lduMatrix.C:186
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
bool foundObject(const word &name) const
Is the named Type in registry.
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 managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Calculate the first temporal derivative.
const cellShapeList & cells
U
Definition: pEqn.H:72
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 > > Su(const VolField< Type > &su, const VolField< Type > &vf)
Definition: fvcSup.C:44
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
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet dimless
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
messageStream Info
const dimensionSet dimTemperature
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
UniformDimensionedField< vector > uniformDimensionedVectorField
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:62
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void pow3(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
const dimensionSet dimMass
error FatalError
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
label timeIndex
Definition: getTimeIndex.H:4
thermo he()
labelList fv(nPoints)
dictionary dict
fluidMulticomponentThermo & thermo
Definition: createFields.H:31