solidificationMeltingSource.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2014-2017 OpenFOAM Foundation
9  Copyright (C) 2018-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
30 #include "fvMatrices.H"
31 #include "basicThermo.H"
32 #include "gravityMeshObject.H"
36 #include "geometricOneField.H"
37 
38 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42 namespace fv
43 {
44  defineTypeNameAndDebug(solidificationMeltingSource, 0);
45  addToRunTimeSelectionTable(option, solidificationMeltingSource, dictionary);
46 }
47 }
48 
49 const Foam::Enum
50 <
52 >
54 ({
55  { thermoMode::mdThermo, "thermo" },
56  { thermoMode::mdLookup, "lookup" },
57 });
58 
59 
60 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
61 
63 Foam::fv::solidificationMeltingSource::Cp() const
64 {
65  switch (mode_)
66  {
67  case mdThermo:
68  {
69  const auto& thermo =
71 
72  return thermo.Cp();
73  break;
74  }
75  case mdLookup:
76  {
77  if (CpName_ == "CpRef")
78  {
79  const scalar CpRef = coeffs_.get<scalar>("CpRef");
80 
81  return volScalarField::New
82  (
84  mesh_,
86  (
87  "Cp",
89  CpRef
90  ),
92  );
93  }
94  else
95  {
96  return mesh_.lookupObject<volScalarField>(CpName_);
97  }
98 
99  break;
100  }
101  default:
102  {
104  << "Unhandled thermo mode: " << thermoModeTypeNames_[mode_]
105  << abort(FatalError);
106  }
107  }
108 
109  return nullptr;
110 }
111 
112 
113 void Foam::fv::solidificationMeltingSource::update(const volScalarField& Cp)
114 {
115  if (curTimeIndex_ == mesh_.time().timeIndex())
116  {
117  return;
118  }
119 
120  if (debug)
121  {
122  Info<< type() << ": " << name_ << " - updating phase indicator" << endl;
123  }
124 
125  if (mesh_.topoChanging())
126  {
127  deltaT_.resize(cells_.size());
128  }
129 
130  // update old time alpha1 field
131  alpha1_.oldTime();
132 
133  const auto& T = mesh_.lookupObject<volScalarField>(TName_);
134 
135  forAll(cells_, i)
136  {
137  label celli = cells_[i];
138 
139  scalar Tc = T[celli];
140  scalar Cpc = Cp[celli];
141  scalar alpha1New = alpha1_[celli] + relax_*Cpc*(Tc - Tmelt_)/L_;
142 
143  alpha1_[celli] = clamp(alpha1New, zero_one{});
144  deltaT_[i] = Tc - Tmelt_;
145  }
146 
147  alpha1_.correctBoundaryConditions();
148 
149  curTimeIndex_ = mesh_.time().timeIndex();
150 }
151 
152 
153 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
154 
156 (
157  const word& sourceName,
158  const word& modelType,
159  const dictionary& dict,
160  const fvMesh& mesh
161 )
162 :
163  fv::cellSetOption(sourceName, modelType, dict, mesh),
164  Tmelt_(coeffs_.get<scalar>("Tmelt")),
165  L_(coeffs_.get<scalar>("L")),
166  relax_(coeffs_.getOrDefault<scalar>("relax", 0.9)),
167  mode_(thermoModeTypeNames_.get("thermoMode", coeffs_)),
168  rhoRef_(coeffs_.get<scalar>("rhoRef")),
169  TName_(coeffs_.getOrDefault<word>("T", "T")),
170  CpName_(coeffs_.getOrDefault<word>("Cp", "Cp")),
171  UName_(coeffs_.getOrDefault<word>("U", "U")),
172  phiName_(coeffs_.getOrDefault<word>("phi", "phi")),
173  Cu_(coeffs_.getOrDefault<scalar>("Cu", 100000)),
174  q_(coeffs_.getOrDefault<scalar>("q", 0.001)),
175  beta_(coeffs_.get<scalar>("beta")),
176  alpha1_
177  (
178  IOobject
179  (
180  IOobject::scopedName(name_, "alpha1"),
181  mesh.time().timeName(),
182  mesh,
183  IOobject::READ_IF_PRESENT,
184  IOobject::AUTO_WRITE,
185  IOobject::REGISTER
186  ),
187  mesh,
190  ),
191  curTimeIndex_(-1),
192  deltaT_(cells_.size(), 0)
193 {
194  fieldNames_.resize(2);
195  fieldNames_[0] = UName_;
196 
197  switch (mode_)
198  {
199  case mdThermo:
200  {
201  const auto& thermo =
203 
204  fieldNames_[1] = thermo.he().name();
205  break;
206  }
207  case mdLookup:
208  {
209  fieldNames_[1] = TName_;
210  break;
211  }
212  default:
213  {
215  << "Unhandled thermo mode: " << thermoModeTypeNames_[mode_]
216  << abort(FatalError);
217  }
218  }
219 
221 }
222 
223 
224 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
225 
227 (
228  fvMatrix<scalar>& eqn,
229  const label fieldi
230 )
231 {
232  apply(geometricOneField(), eqn);
233 }
234 
235 
237 (
238  const volScalarField& rho,
239  fvMatrix<scalar>& eqn,
240  const label fieldi
241 )
242 {
243  apply(rho, eqn);
244 }
245 
246 
248 (
249  fvMatrix<vector>& eqn,
250  const label fieldi
251 )
252 {
253  if (debug)
254  {
255  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
256  }
257 
258  const volScalarField Cp(this->Cp());
259 
260  update(Cp);
261 
262  const vector& g = meshObjects::gravity::New(mesh_.time()).value();
263 
264  scalarField& Sp = eqn.diag();
265  vectorField& Su = eqn.source();
266  const scalarField& V = mesh_.V();
267 
268  forAll(cells_, i)
269  {
270  label celli = cells_[i];
271 
272  scalar Vc = V[celli];
273  scalar alpha1c = alpha1_[celli];
274 
275  scalar S = -Cu_*sqr(1.0 - alpha1c)/(pow3(alpha1c) + q_);
276  vector Sb(rhoRef_*g*beta_*deltaT_[i]);
277 
278  Sp[celli] += Vc*S;
279  Su[celli] += Vc*Sb;
280  }
281 }
282 
283 
285 (
286  const volScalarField& rho,
287  fvMatrix<vector>& eqn,
288  const label fieldi
289 )
290 {
291  // Momentum source uses a Boussinesq approximation - redirect
292  addSup(eqn, fieldi);
293 }
294 
295 
297 {
299  {
300  coeffs_.readEntry("Tmelt", Tmelt_);
301  coeffs_.readEntry("L", L_);
302 
303  coeffs_.readIfPresent("relax", relax_);
304 
305  thermoModeTypeNames_.readEntry("thermoMode", coeffs_, mode_);
306 
307  coeffs_.readEntry("rhoRef", rhoRef_);
308  coeffs_.readIfPresent("T", TName_);
309  coeffs_.readIfPresent("U", UName_);
310 
311  coeffs_.readIfPresent("Cu", Cu_);
312  coeffs_.readIfPresent("q", q_);
313 
314  coeffs_.readEntry("beta", beta_);
315 
316  return true;
317  }
318 
319  return false;
320 }
321 
322 
323 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:53
const scalarField & diag() const
Definition: lduMatrix.C:163
dictionary dict
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
thermoMode
Options for the thermo mode specification.
const word zeroGradientType
A zeroGradient patch field type.
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:59
zeroField Su
Definition: alphaSuSp.H:1
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition: fvOption.H:157
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:153
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
defineTypeNameAndDebug(atmAmbientTurbSource, 0)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:142
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static void apply(bitSet &selection, const Detail::parcelSelection::actionType action, const Predicate &accept, const UList< Type > &list, const AccessOp &aop)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
virtual bool read(const dictionary &dict)
Read source dictionary.
Template invariant parts for fvPatchField.
Definition: fvPatchField.H:77
const dimensionSet dimless
Dimensionless.
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
static const Enum< thermoMode > thermoModeTypeNames_
Names for thermoMode.
Macros for easy insertion into run-time selection tables.
solidificationMeltingSource(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
psiReactionThermo & thermo
Definition: createFields.H:28
word timeName
Definition: getTimeIndex.H:3
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:40
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition: fvMatrix.H:485
dynamicFvMesh & mesh
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
const word name_
Source name.
Definition: fvOption.H:132
A class for handling words, derived from Foam::string.
Definition: word.H:63
labelList fv(nPoints)
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
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:64
errorManip< error > abort(error &err)
Definition: errorManip.H:139
virtual bool read(const dictionary &dict)
Read source dictionary.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
addToRunTimeSelectionTable(option, atmAmbientTurbSource, dictionary)
const uniformDimensionedVectorField & g
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
const volScalarField & Cp
Definition: EEqn.H:7
int debug
Static debugging option.
mesh update()
static const gravity & New(const word &name, const Time &runTime)
Return named gravity field cached or construct on Time.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const dimensionSet dimEnergy
virtual void addSup(fvMatrix< scalar > &eqn, const label fieldi)
Add explicit contribution to enthalpy equation.
Field< Type > & source() noexcept
Definition: fvMatrix.H:533
dimensionedScalar pow3(const dimensionedScalar &ds)
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: error.H:64
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition: fvOption.C:41
A special matrix type and solver, designed for finite volume solutions of scalar equations.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Intermediate abstract class for handling cell-set options for the derived fvOptions.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: HashPtrTable.H:50
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
Definition: fvPatchField.H:213
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
dictionary coeffs_
Dictionary containing source coefficients.
Definition: fvOption.H:152
Namespace for OpenFOAM.
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
Definition: dimensionSet.C:271
zeroField Sp
Definition: alphaSuSp.H:2
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127