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-2021 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 
82  (
83  IOobject
84  (
85  name_ + ":Cp",
86  mesh_.time().timeName(),
87  mesh_,
90  ),
91  mesh_,
93  (
94  "Cp",
96  CpRef
97  ),
99  );
100  }
101  else
102  {
103  return mesh_.lookupObject<volScalarField>(CpName_);
104  }
105 
106  break;
107  }
108  default:
109  {
111  << "Unhandled thermo mode: " << thermoModeTypeNames_[mode_]
112  << abort(FatalError);
113  }
114  }
115 
116  return nullptr;
117 }
118 
119 
120 void Foam::fv::solidificationMeltingSource::update(const volScalarField& Cp)
121 {
122  if (curTimeIndex_ == mesh_.time().timeIndex())
123  {
124  return;
125  }
126 
127  if (debug)
128  {
129  Info<< type() << ": " << name_ << " - updating phase indicator" << endl;
130  }
131 
132  if (mesh_.topoChanging())
133  {
134  deltaT_.resize(cells_.size());
135  }
136 
137  // update old time alpha1 field
138  alpha1_.oldTime();
139 
140  const auto& T = mesh_.lookupObject<volScalarField>(TName_);
141 
142  forAll(cells_, i)
143  {
144  label celli = cells_[i];
145 
146  scalar Tc = T[celli];
147  scalar Cpc = Cp[celli];
148  scalar alpha1New = alpha1_[celli] + relax_*Cpc*(Tc - Tmelt_)/L_;
149 
150  alpha1_[celli] = clamp(alpha1New, zero_one{});
151  deltaT_[i] = Tc - Tmelt_;
152  }
153 
154  alpha1_.correctBoundaryConditions();
155 
156  curTimeIndex_ = mesh_.time().timeIndex();
157 }
158 
159 
160 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
161 
163 (
164  const word& sourceName,
165  const word& modelType,
166  const dictionary& dict,
167  const fvMesh& mesh
168 )
169 :
170  fv::cellSetOption(sourceName, modelType, dict, mesh),
171  Tmelt_(coeffs_.get<scalar>("Tmelt")),
172  L_(coeffs_.get<scalar>("L")),
173  relax_(coeffs_.getOrDefault<scalar>("relax", 0.9)),
174  mode_(thermoModeTypeNames_.get("thermoMode", coeffs_)),
175  rhoRef_(coeffs_.get<scalar>("rhoRef")),
176  TName_(coeffs_.getOrDefault<word>("T", "T")),
177  CpName_(coeffs_.getOrDefault<word>("Cp", "Cp")),
178  UName_(coeffs_.getOrDefault<word>("U", "U")),
179  phiName_(coeffs_.getOrDefault<word>("phi", "phi")),
180  Cu_(coeffs_.getOrDefault<scalar>("Cu", 100000)),
181  q_(coeffs_.getOrDefault<scalar>("q", 0.001)),
182  beta_(coeffs_.get<scalar>("beta")),
183  alpha1_
184  (
185  IOobject
186  (
187  name_ + ":alpha1",
188  mesh.time().timeName(),
189  mesh,
190  IOobject::READ_IF_PRESENT,
191  IOobject::AUTO_WRITE
192  ),
193  mesh,
196  ),
197  curTimeIndex_(-1),
198  deltaT_(cells_.size(), 0)
199 {
200  fieldNames_.resize(2);
201  fieldNames_[0] = UName_;
202 
203  switch (mode_)
204  {
205  case mdThermo:
206  {
207  const auto& thermo =
209 
210  fieldNames_[1] = thermo.he().name();
211  break;
212  }
213  case mdLookup:
214  {
215  fieldNames_[1] = TName_;
216  break;
217  }
218  default:
219  {
221  << "Unhandled thermo mode: " << thermoModeTypeNames_[mode_]
222  << abort(FatalError);
223  }
224  }
225 
227 }
228 
229 
230 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
231 
233 (
234  fvMatrix<scalar>& eqn,
235  const label fieldi
236 )
237 {
238  apply(geometricOneField(), eqn);
239 }
240 
241 
243 (
244  const volScalarField& rho,
245  fvMatrix<scalar>& eqn,
246  const label fieldi
247 )
248 {
249  apply(rho, eqn);
250 }
251 
252 
254 (
255  fvMatrix<vector>& eqn,
256  const label fieldi
257 )
258 {
259  if (debug)
260  {
261  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
262  }
263 
264  const volScalarField Cp(this->Cp());
265 
266  update(Cp);
267 
268  const vector& g = meshObjects::gravity::New(mesh_.time()).value();
269 
270  scalarField& Sp = eqn.diag();
271  vectorField& Su = eqn.source();
272  const scalarField& V = mesh_.V();
273 
274  forAll(cells_, i)
275  {
276  label celli = cells_[i];
277 
278  scalar Vc = V[celli];
279  scalar alpha1c = alpha1_[celli];
280 
281  scalar S = -Cu_*sqr(1.0 - alpha1c)/(pow3(alpha1c) + q_);
282  vector Sb(rhoRef_*g*beta_*deltaT_[i]);
283 
284  Sp[celli] += Vc*S;
285  Su[celli] += Vc*Sb;
286  }
287 }
288 
289 
291 (
292  const volScalarField& rho,
293  fvMatrix<vector>& eqn,
294  const label fieldi
295 )
296 {
297  // Momentum source uses a Boussinesq approximation - redirect
298  addSup(eqn, fieldi);
299 }
300 
301 
303 {
305  {
306  coeffs_.readEntry("Tmelt", Tmelt_);
307  coeffs_.readEntry("L", L_);
308 
309  coeffs_.readIfPresent("relax", relax_);
310 
311  thermoModeTypeNames_.readEntry("thermoMode", coeffs_, mode_);
312 
313  coeffs_.readEntry("rhoRef", rhoRef_);
314  coeffs_.readIfPresent("T", TName_);
315  coeffs_.readIfPresent("U", UName_);
316 
317  coeffs_.readIfPresent("Cu", Cu_);
318  coeffs_.readIfPresent("q", q_);
319 
320  coeffs_.readEntry("beta", beta_);
321 
322  return true;
323  }
324 
325  return false;
326 }
327 
328 
329 // ************************************************************************* //
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
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:160
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:598
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
Ignore writing from objectRegistry::writeObject()
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.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
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:81
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
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
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
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...
addToRunTimeSelectionTable(option, atmAmbientTurbSource, dictionary)
const uniformDimensionedVectorField & g
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
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
Nothing to be read.
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
scalarField & diag()
Definition: lduMatrix.C:197
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:172
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