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-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 
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 
57 namespace Foam
58 {
59  template<>
61  {"thermo", "lookup"};
62 }
63 
66 
67 
68 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
69 
70 void Foam::fv::solidificationMelting::readCoeffs()
71 {
72  Tsol_ = coeffs().lookup<scalar>("Tsol");
73  Tliq_ = coeffs().lookupOrDefault<scalar>("Tliq", Tsol_);
74  alpha1e_ = coeffs().lookupOrDefault<scalar>("alpha1e", 0.0);
75  L_ = coeffs().lookup<scalar>("L");
76 
77  relax_ = coeffs().lookupOrDefault("relax", 0.9);
78 
79  mode_ = thermoModeTypeNames_.read(coeffs().lookup("thermoMode"));
80 
81  rhoRef_ = coeffs().lookup<scalar>("rhoRef");
82  TName_ = coeffs().lookupOrDefault<word>("T", "T");
83  CpName_ = coeffs().lookupOrDefault<word>("Cp", "Cp");
84  UName_ = coeffs().lookupOrDefault<word>("U", "U");
85  phiName_ = coeffs().lookupOrDefault<word>("phi", "phi");
86 
87  Cu_ = coeffs().lookupOrDefault<scalar>("Cu", 100000);
88  q_ = coeffs().lookupOrDefault("q", 0.001);
89 
90  beta_ = coeffs().lookup<scalar>("beta");
91 }
92 
93 
95 Foam::fv::solidificationMelting::Cp() const
96 {
97  switch (mode_)
98  {
99  case thermoMode::thermo:
100  {
101  const basicThermo& thermo =
102  mesh().lookupObject<basicThermo>(physicalProperties::typeName);
103 
104  return thermo.Cp();
105  break;
106  }
107  case thermoMode::lookup:
108  {
109  if (CpName_ == "CpRef")
110  {
111  scalar CpRef = coeffs().lookup<scalar>("CpRef");
112 
113  return volScalarField::New
114  (
115  name() + ":Cp",
116  mesh(),
118  (
120  CpRef
121  ),
122  extrapolatedCalculatedFvPatchScalarField::typeName
123  );
124  }
125  else
126  {
127  return mesh().lookupObject<volScalarField>(CpName_);
128  }
129 
130  break;
131  }
132  default:
133  {
135  << "Unhandled thermo mode: " << thermoModeTypeNames_[mode_]
136  << abort(FatalError);
137  }
138  }
139 
140  return tmp<volScalarField>(nullptr);
141 }
142 
143 
144 Foam::vector Foam::fv::solidificationMelting::g() const
145 {
146  if (mesh().foundObject<uniformDimensionedVectorField>("g"))
147  {
148  const uniformDimensionedVectorField& value =
149  mesh().lookupObject<uniformDimensionedVectorField>("g");
150  return value.value();
151  }
152  else
153  {
154  return coeffs().lookup("g");
155  }
156 }
157 
158 
159 void Foam::fv::solidificationMelting::update
160 (
161  const volScalarField& Cp
162 ) const
163 {
164  if (curTimeIndex_ == mesh().time().timeIndex())
165  {
166  return;
167  }
168 
169  if (debug)
170  {
171  Info<< type() << ": " << name()
172  << " - updating phase indicator" << endl;
173  }
174 
175  // update old time alpha1 field
176  alpha1_.oldTime();
177 
178  const volScalarField& T = mesh().lookupObject<volScalarField>(TName_);
179 
180  const labelUList cells = set_.cells();
181 
182  forAll(cells, i)
183  {
184  const label celli = cells[i];
185 
186  const scalar Tc = T[celli];
187  const scalar Cpc = Cp[celli];
188  const scalar alpha1New =
189  alpha1_[celli]
190  + relax_*Cpc
191  *(
192  Tc
193  - max
194  (
195  Tsol_,
196  Tsol_
197  + (Tliq_ - Tsol_)*(alpha1_[celli] - alpha1e_)/(1 - alpha1e_)
198  )
199  )/L_;
200 
201  alpha1_[celli] = max(0, min(alpha1New, 1));
202  deltaT_[i] =
203  Tc
204  - max
205  (
206  Tsol_,
207  Tsol_
208  + (Tliq_ - Tsol_)*(alpha1_[celli] - alpha1e_)/(1 - alpha1e_)
209  );
210  }
211 
212  alpha1_.correctBoundaryConditions();
213 
214  curTimeIndex_ = mesh().time().timeIndex();
215 }
216 
217 
218 template<class RhoFieldType>
219 void Foam::fv::solidificationMelting::apply
220 (
221  const RhoFieldType& rho,
222  fvMatrix<scalar>& eqn
223 ) const
224 {
225  if (debug)
226  {
227  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
228  }
229 
230  const volScalarField Cp(this->Cp());
231 
232  update(Cp);
233 
234  dimensionedScalar L("L", dimEnergy/dimMass, L_);
235 
236  // Contributions added to rhs of solver equation
237  if (eqn.psi().dimensions() == dimTemperature)
238  {
239  eqn -= L/Cp*(fvc::ddt(rho, alpha1_));
240  }
241  else
242  {
243  eqn -= L*(fvc::ddt(rho, alpha1_));
244  }
245 }
246 
247 
248 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
249 
251 (
252  const word& name,
253  const word& modelType,
254  const fvMesh& mesh,
255  const dictionary& dict
256 )
257 :
258  fvModel(name, modelType, mesh, dict),
259  set_(mesh, coeffs()),
260  Tsol_(NaN),
261  Tliq_(NaN),
262  alpha1e_(NaN),
263  L_(NaN),
264  relax_(NaN),
265  mode_(thermoMode::thermo),
266  rhoRef_(NaN),
267  TName_(word::null),
268  CpName_(word::null),
269  UName_(word::null),
270  phiName_(word::null),
271  Cu_(NaN),
272  q_(NaN),
273  beta_(NaN),
274  alpha1_
275  (
276  IOobject
277  (
278  this->name() + ":alpha1",
279  mesh.time().name(),
280  mesh,
281  IOobject::READ_IF_PRESENT,
282  IOobject::AUTO_WRITE
283  ),
284  mesh,
286  zeroGradientFvPatchScalarField::typeName
287  ),
288  curTimeIndex_(-1),
289  deltaT_(set_.nCells(), 0)
290 {
291  readCoeffs();
292 }
293 
294 
295 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
296 
298 {
299  switch (mode_)
300  {
301  case thermoMode::thermo:
302  {
303  const basicThermo& thermo =
304  mesh().lookupObject<basicThermo>(physicalProperties::typeName);
305 
306  return wordList({UName_, thermo.he().name()});
307  }
308  case thermoMode::lookup:
309  {
310  return wordList({UName_, TName_});
311  }
312  }
313 
314  return wordList::null();
315 }
316 
317 
319 (
320  const volScalarField& he,
321  fvMatrix<scalar>& eqn
322 ) const
323 {
324  apply(geometricOneField(), eqn);
325 }
326 
327 
329 (
330  const volScalarField& rho,
331  const volScalarField& he,
332  fvMatrix<scalar>& eqn
333 ) const
334 {
335  apply(rho, eqn);
336 }
337 
338 
340 (
341  const volVectorField& U,
342  fvMatrix<vector>& eqn
343 ) const
344 {
345  if (debug)
346  {
347  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
348  }
349 
350  const volScalarField Cp(this->Cp());
351 
352  update(Cp);
353 
354  vector g = this->g();
355 
356  scalarField& Sp = eqn.diag();
357  vectorField& Su = eqn.source();
358  const scalarField& V = mesh().V();
359 
360  const labelUList cells = set_.cells();
361 
362  forAll(cells, i)
363  {
364  const label celli = cells[i];
365 
366  const scalar Vc = V[celli];
367  const scalar alpha1c = alpha1_[celli];
368 
369  const scalar S = -Cu_*sqr(1.0 - alpha1c)/(pow3(alpha1c) + q_);
370  const vector Sb = rhoRef_*g*beta_*deltaT_[i];
371 
372  Sp[celli] += Vc*S;
373  Su[celli] += Vc*Sb;
374  }
375 }
376 
377 
379 (
380  const volScalarField& rho,
381  const volVectorField& U,
382  fvMatrix<vector>& eqn
383 ) const
384 {
385  addSup(U, eqn);
386 }
387 
388 
390 {
391  set_.movePoints();
392  return true;
393 }
394 
395 
397 (
398  const polyTopoChangeMap& map
399 )
400 {
401  set_.topoChange(map);
402 }
403 
404 
406 {
407  set_.mapMesh(map);
408 }
409 
410 
412 (
413  const polyDistributionMap& map
414 )
415 {
416  set_.distribute(map);
417 }
418 
419 
421 {
422  if (fvModel::read(dict))
423  {
424  set_.read(coeffs());
425  readCoeffs();
426  return true;
427  }
428  else
429  {
430  return false;
431  }
432 
433  return false;
434 }
435 
436 
437 // ************************************************************************* //
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:434
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, 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:54
Enum read(Istream &) const
Read a word from Istream and return the corresponding.
Definition: NamedEnum.C:61
Type & value()
Return a reference to the value.
Base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:78
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.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
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: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
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
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
#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
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
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet dimless
messageStream Info
const dimensionSet dimTemperature
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
UniformDimensionedField< vector > uniformDimensionedVectorField
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimMass
error FatalError
UList< label > labelUList
Definition: UList.H:65
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
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