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 
70 {
71  switch (mode_)
72  {
73  case mdThermo:
74  {
75  const basicThermo& thermo =
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  {
113  }
114 
115  break;
116  }
117  default:
118  {
120  << "Unhandled thermo mode: " << thermoModeTypeNames_[mode_]
121  << abort(FatalError);
122  }
123  }
124 
125  return tmp<volScalarField>(NULL);
126 }
127 
128 
130 {
131  if (mesh_.foundObject<uniformDimensionedVectorField>("g"))
132  {
133  const uniformDimensionedVectorField& value =
134  mesh_.lookupObject<uniformDimensionedVectorField>("g");
135  return value.value();
136  }
137  else
138  {
139  return coeffs_.lookup("g");
140  }
141 }
142 
143 
145 {
146  if (curTimeIndex_ == mesh_.time().timeIndex())
147  {
148  return;
149  }
150 
151  if (debug)
152  {
153  Info<< type() << ": " << name_ << " - updating phase indicator" << endl;
154  }
155 
156  // update old time alpha1 field
157  alpha1_.oldTime();
158 
159  const volScalarField& T = mesh_.lookupObject<volScalarField>(TName_);
160 
161  forAll(cells_, i)
162  {
163  label cellI = cells_[i];
164 
165  scalar Tc = T[cellI];
166  scalar Cpc = Cp[cellI];
167  scalar alpha1New = alpha1_[cellI] + relax_*Cpc*(Tc - Tmelt_)/L_;
168 
169  alpha1_[cellI] = max(0, min(alpha1New, 1));
170  deltaT_[i] = Tc - Tmelt_;
171  }
172 
173  alpha1_.correctBoundaryConditions();
174 
175  curTimeIndex_ = mesh_.time().timeIndex();
176 }
177 
178 
179 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
180 
182 (
183  const word& sourceName,
184  const word& modelType,
185  const dictionary& dict,
186  const fvMesh& mesh
187 )
188 :
189  cellSetOption(sourceName, modelType, dict, mesh),
190  Tmelt_(readScalar(coeffs_.lookup("Tmelt"))),
191  L_(readScalar(coeffs_.lookup("L"))),
192  relax_(coeffs_.lookupOrDefault("relax", 0.9)),
193  mode_(thermoModeTypeNames_.read(coeffs_.lookup("thermoMode"))),
194  rhoRef_(readScalar(coeffs_.lookup("rhoRef"))),
195  TName_(coeffs_.lookupOrDefault<word>("TName", "T")),
196  CpName_(coeffs_.lookupOrDefault<word>("CpName", "Cp")),
197  UName_(coeffs_.lookupOrDefault<word>("UName", "U")),
198  phiName_(coeffs_.lookupOrDefault<word>("phiName", "phi")),
199  Cu_(coeffs_.lookupOrDefault<scalar>("Cu", 100000)),
200  q_(coeffs_.lookupOrDefault("q", 0.001)),
201  beta_(readScalar(coeffs_.lookup("beta"))),
202  alpha1_
203  (
204  IOobject
205  (
206  name_ + ":alpha1",
207  mesh.time().timeName(),
208  mesh,
211  ),
212  mesh,
213  dimensionedScalar("alpha1", dimless, 0.0),
214  zeroGradientFvPatchScalarField::typeName
215  ),
216  curTimeIndex_(-1),
217  deltaT_(cells_.size(), 0)
218 {
219  fieldNames_.setSize(2);
220  fieldNames_[0] = UName_;
221 
222  switch (mode_)
223  {
224  case mdThermo:
225  {
226  const basicThermo& thermo =
227  mesh_.lookupObject<basicThermo>(basicThermo::dictName);
228 
229  fieldNames_[1] = thermo.he().name();
230  break;
231  }
232  case mdLookup:
233  {
234  fieldNames_[1] = TName_;
235  break;
236  }
237  default:
238  {
240  << "Unhandled thermo mode: " << thermoModeTypeNames_[mode_]
241  << abort(FatalError);
242  }
243  }
244 
245  applied_.setSize(2, false);
246 }
247 
248 
249 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
250 
252 (
253  fvMatrix<scalar>& eqn,
254  const label fieldI
255 )
256 {
257  apply(geometricOneField(), eqn);
258 }
259 
260 
262 (
263  const volScalarField& rho,
264  fvMatrix<scalar>& eqn,
265  const label fieldI
266 )
267 {
268  apply(rho, eqn);
269 }
270 
271 
273 (
274  fvMatrix<vector>& eqn,
275  const label fieldI
276 )
277 {
278  if (debug)
279  {
280  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
281  }
282 
283  const volScalarField Cp(this->Cp());
284 
285  update(Cp);
286 
287  vector g = this->g();
288 
289  scalarField& Sp = eqn.diag();
290  vectorField& Su = eqn.source();
291  const scalarField& V = mesh_.V();
292 
293  forAll(cells_, i)
294  {
295  label cellI = cells_[i];
296 
297  scalar Vc = V[cellI];
298  scalar alpha1c = alpha1_[cellI];
299 
300  scalar S = -Cu_*sqr(1.0 - alpha1c)/(pow3(alpha1c) + q_);
301  vector Sb = rhoRef_*g*beta_*deltaT_[i];
302 
303  Sp[cellI] += Vc*S;
304  Su[cellI] += Vc*Sb;
305  }
306 }
307 
308 
310 (
311  const volScalarField& rho,
312  fvMatrix<vector>& eqn,
313  const label fieldI
314 )
315 {
316  // Momentum source uses a Boussinesq approximation - redirect
317  addSup(eqn, fieldI);
318 }
319 
320 
321 // ************************************************************************* //
solidificationMeltingSource.H
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
basicThermo.H
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::fv::cellSetOption
Cell-set options abtract base class. Provides a base set of controls, e.g.
Definition: cellSetOption.H:71
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::dimEnergy
const dimensionSet dimEnergy
Foam::basicThermo::dictName
static const word dictName
Definition: basicThermo.H:176
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::geometricOneField
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Definition: geometricOneField.H:52
Foam::fv::solidificationMeltingSource::g
vector g() const
Return the gravity vector.
Definition: solidificationMeltingSource.C:129
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::basicThermo
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:52
g
const dimensionedVector & g
Definition: setRegionFluidFields.H:33
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:261
Foam::fv::solidificationMeltingSource::update
void update(const volScalarField &Cp)
Update the model.
Definition: solidificationMeltingSource.C:144
Foam::fv::option::name_
const word name_
Source name.
Definition: fvOption.H:72
Foam::fv::solidificationMeltingSource::thermoModeTypeNames_
static const NamedEnum< thermoMode, 2 > thermoModeTypeNames_
Definition: solidificationMeltingSource.H:200
Foam::fvMatrix::psi
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:281
Foam::fv::solidificationMeltingSource::thermoMode
thermoMode
Definition: solidificationMeltingSource.H:194
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::fv::option::mesh_
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:78
fvMatrices.H
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Foam::fv::solidificationMeltingSource::CpName_
word CpName_
Name of specific heat capacity field - default = "Cp" (optional)
Definition: solidificationMeltingSource.H:226
Foam::UniformDimensionedField< vector >
Foam::label
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
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::fvc::Su
tmp< GeometricField< Type, fvPatchField, volMesh > > Su(const GeometricField< Type, fvPatchField, volMesh > &su, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:44
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:87
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::Info
messageStream Info
Foam::fv::option::coeffs_
dictionary coeffs_
Dictionary containing source coefficients.
Definition: fvOption.H:84
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
Foam::fv::solidificationMeltingSource::Cp
tmp< volScalarField > Cp() const
Return the specific heat capacity field.
Definition: solidificationMeltingSource.C:69
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::fv::solidificationMeltingSource::solidificationMeltingSource
solidificationMeltingSource(const solidificationMeltingSource &)
Disallow default bitwise copy construct.
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::fv::solidificationMeltingSource::mode_
thermoMode mode_
Thermodynamics mode.
Definition: solidificationMeltingSource.H:217
rho
rho
Definition: pEqn.H:3
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
fv
labelList fv(nPoints)
Foam::fv::defineTypeNameAndDebug
defineTypeNameAndDebug(option, 0)
uniformDimensionedFields.H
Foam::fv::solidificationMeltingSource::mdThermo
@ mdThermo
Definition: solidificationMeltingSource.H:196
Foam::fv::solidificationMeltingSource::addSup
virtual void addSup(fvMatrix< scalar > &eqn, const label fieldI)
Add explicit contribution to enthalpy equation.
Definition: solidificationMeltingSource.C:252
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::fv::addToRunTimeSelectionTable
addToRunTimeSelectionTable(option, fixedTemperatureConstraint, dictionary)
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
geometricOneField.H
Foam::readScalar
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
Foam::Vector< scalar >
Foam::Time::timeName
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
Foam::fv::solidificationMeltingSource::mdLookup
@ mdLookup
Definition: solidificationMeltingSource.H:197
Foam::fvMatrix::source
Field< Type > & source()
Definition: fvMatrix.H:291
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::dimTemperature
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Definition: objectRegistryTemplates.C:165
Foam::fvc::Sp
tmp< GeometricField< Type, fvPatchField, volMesh > > Sp(const volScalarField &sp, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:67
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::IOobject::READ_IF_PRESENT
@ READ_IF_PRESENT
Definition: IOobject.H:110
zeroGradientFvPatchFields.H
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::NamedEnum
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:52