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 | Copyright (C) 2014-2015 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 
27 #include "fvMatrices.H"
28 #include "basicThermo.H"
32 #include "geometricOneField.H"
33 
34 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  template<>
39  const char* NamedEnum
40  <
42  2
43  >::names[] =
44  {
45  "thermo",
46  "lookup"
47  };
48 
49  namespace fv
50  {
51  defineTypeNameAndDebug(solidificationMeltingSource, 0);
52 
54  (
55  option,
56  solidificationMeltingSource,
58  );
59  }
60 }
61 
64 
65 
66 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
67 
69 Foam::fv::solidificationMeltingSource::Cp() const
70 {
71  switch (mode_)
72  {
73  case mdThermo:
74  {
75  const basicThermo& thermo =
76  mesh_.lookupObject<basicThermo>(basicThermo::dictName);
77 
78  return thermo.Cp();
79  break;
80  }
81  case mdLookup:
82  {
83  if (CpName_ == "CpRef")
84  {
85  scalar CpRef = readScalar(coeffs_.lookup("CpRef"));
86 
87  return tmp<volScalarField>
88  (
89  new volScalarField
90  (
91  IOobject
92  (
93  name_ + ":Cp",
94  mesh_.time().timeName(),
95  mesh_,
98  ),
99  mesh_,
101  (
102  "Cp",
104  CpRef
105  ),
106  zeroGradientFvPatchScalarField::typeName
107  )
108  );
109  }
110  else
111  {
112  return mesh_.lookupObject<volScalarField>(CpName_);
113  }
114 
115  break;
116  }
117  default:
118  {
120  (
121  "Foam::tmp<Foam::volScalarField> "
122  "Foam::fv::solidificationMeltingSource::Cp() const"
123  )
124  << "Unhandled thermo mode: " << thermoModeTypeNames_[mode_]
125  << abort(FatalError);
126  }
127  }
128 
129  return tmp<volScalarField>(NULL);
130 }
131 
132 
133 Foam::vector Foam::fv::solidificationMeltingSource::g() const
134 {
135  if (mesh_.foundObject<uniformDimensionedVectorField>("g"))
136  {
137  const uniformDimensionedVectorField& value =
138  mesh_.lookupObject<uniformDimensionedVectorField>("g");
139  return value.value();
140  }
141  else
142  {
143  return coeffs_.lookup("g");
144  }
145 }
146 
147 
148 void Foam::fv::solidificationMeltingSource::update(const volScalarField& Cp)
149 {
150  if (curTimeIndex_ == mesh_.time().timeIndex())
151  {
152  return;
153  }
154 
155  if (debug)
156  {
157  Info<< type() << ": " << name_ << " - updating phase indicator" << endl;
158  }
159 
160  // update old time alpha1 field
161  alpha1_.oldTime();
162 
163  const volScalarField& T = mesh_.lookupObject<volScalarField>(TName_);
164 
165  forAll(cells_, i)
166  {
167  label cellI = cells_[i];
168 
169  scalar Tc = T[cellI];
170  scalar Cpc = Cp[cellI];
171  scalar alpha1New = alpha1_[cellI] + relax_*Cpc*(Tc - Tmelt_)/L_;
172 
173  alpha1_[cellI] = max(0, min(alpha1New, 1));
174  deltaT_[i] = Tc - Tmelt_;
175  }
176 
177  alpha1_.correctBoundaryConditions();
178 
179  curTimeIndex_ = mesh_.time().timeIndex();
180 }
181 
182 
183 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
184 
185 Foam::fv::solidificationMeltingSource::solidificationMeltingSource
186 (
187  const word& sourceName,
188  const word& modelType,
189  const dictionary& dict,
190  const fvMesh& mesh
191 )
192 :
193  cellSetOption(sourceName, modelType, dict, mesh),
194  Tmelt_(readScalar(coeffs_.lookup("Tmelt"))),
195  L_(readScalar(coeffs_.lookup("L"))),
196  relax_(coeffs_.lookupOrDefault("relax", 0.9)),
197  mode_(thermoModeTypeNames_.read(coeffs_.lookup("thermoMode"))),
198  rhoRef_(readScalar(coeffs_.lookup("rhoRef"))),
199  TName_(coeffs_.lookupOrDefault<word>("TName", "T")),
200  CpName_(coeffs_.lookupOrDefault<word>("CpName", "Cp")),
201  UName_(coeffs_.lookupOrDefault<word>("UName", "U")),
202  phiName_(coeffs_.lookupOrDefault<word>("phiName", "phi")),
203  Cu_(coeffs_.lookupOrDefault<scalar>("Cu", 100000)),
204  q_(coeffs_.lookupOrDefault("q", 0.001)),
205  beta_(readScalar(coeffs_.lookup("beta"))),
206  alpha1_
207  (
208  IOobject
209  (
210  name_ + ":alpha1",
211  mesh.time().timeName(),
212  mesh,
215  ),
216  mesh,
217  dimensionedScalar("alpha1", dimless, 0.0),
218  zeroGradientFvPatchScalarField::typeName
219  ),
220  curTimeIndex_(-1),
221  deltaT_(cells_.size(), 0)
222 {
223  fieldNames_.setSize(2);
224  fieldNames_[0] = UName_;
225 
226  switch (mode_)
227  {
228  case mdThermo:
229  {
230  const basicThermo& thermo =
231  mesh_.lookupObject<basicThermo>(basicThermo::dictName);
232 
233  fieldNames_[1] = thermo.he().name();
234  break;
235  }
236  case mdLookup:
237  {
238  fieldNames_[1] = TName_;
239  break;
240  }
241  default:
242  {
244  (
245  "fv::solidificationMeltingSource::solidificationMeltingSource"
246  ) << "Unhandled thermo mode: " << thermoModeTypeNames_[mode_]
247  << abort(FatalError);
248  }
249  }
250 
251  applied_.setSize(2, false);
252 }
253 
254 
255 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
256 
258 (
259  fvMatrix<scalar>& eqn,
260  const label fieldI
261 )
262 {
263  apply(geometricOneField(), eqn);
264 }
265 
266 
268 (
269  const volScalarField& rho,
270  fvMatrix<scalar>& eqn,
271  const label fieldI
272 )
273 {
274  apply(rho, eqn);
275 }
276 
277 
279 (
280  fvMatrix<vector>& eqn,
281  const label fieldI
282 )
283 {
284  if (debug)
285  {
286  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
287  }
288 
289  const volScalarField Cp(this->Cp());
290 
291  update(Cp);
292 
293  vector g = this->g();
294 
295  scalarField& Sp = eqn.diag();
296  vectorField& Su = eqn.source();
297  const scalarField& V = mesh_.V();
298 
299  forAll(cells_, i)
300  {
301  label cellI = cells_[i];
302 
303  scalar Vc = V[cellI];
304  scalar alpha1c = alpha1_[cellI];
305 
306  scalar S = -Cu_*sqr(1.0 - alpha1c)/(pow3(alpha1c) + q_);
307  vector Sb = rhoRef_*g*beta_*deltaT_[i];
308 
309  Sp[cellI] += Vc*S;
310  Su[cellI] += Vc*Sb;
311  }
312 }
313 
314 
316 (
317  const volScalarField& rho,
318  fvMatrix<vector>& eqn,
319  const label fieldI
320 )
321 {
322  // Momentum source uses a Boussinesq approximation - redirect
323  addSup(eqn, fieldI);
324 }
325 
326 
327 // ************************************************************************* //
defineTypeNameAndDebug(cellSetOption, 0)
Cell-set options abtract base class. Provides a base set of controls, e.g.
Definition: cellSetOption.H:71
dimensionedScalar pow3(const dimensionedScalar &ds)
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:52
A class for handling words, derived from string.
Definition: word.H:59
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
tmp< GeometricField< Type, fvPatchField, volMesh > > Sp(const volScalarField &sp, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:67
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
const dimensionSet dimEnergy
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:68
messageStream Info
dynamicFvMesh & mesh
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
virtual tmp< volScalarField > Cp() const =0
Heat capacity at constant pressure [J/kg/K].
Namespace for OpenFOAM.
tmp< GeometricField< Type, fvPatchField, volMesh > > Su(const GeometricField< Type, fvPatchField, volMesh > &su, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:44
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:281
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
addToRunTimeSelectionTable(option, fixedTemperatureConstraint, dictionary)
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define forAll(list, i)
Definition: UList.H:421
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
Macros for easy insertion into run-time selection tables.
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const word & name() const
Return name.
Definition: IOobject.H:260
scalarField & diag()
Definition: lduMatrix.C:183
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:589
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:52
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
Dimensioned<Type> registered with the database as a registered IOobject which has the functionality o...
const dimensionedVector & g
error FatalError
static const NamedEnum< thermoMode, 2 > thermoModeTypeNames_
labelList fv(nPoints)
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
const word dictName() const
Return the local dictionary name (final part of scoped name)
Definition: dictionary.H:115
Field< Type > & source()
Definition: fvMatrix.H:291
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
virtual void addSup(fvMatrix< scalar > &eqn, const label fieldI)
Add explicit contribution to enthalpy equation.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A class for managing temporary objects.
Definition: PtrList.H:118
const Type & value() const
Return const reference to value.
A special matrix type and solver, designed for finite volume solutions of scalar equations.
virtual volScalarField & he()=0
Enthalpy/Internal energy [J/kg].