VoFSolidificationMeltingSource.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) 2017-2018 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 "twoPhaseMixtureThermo.H"
30 
31 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  namespace fv
36  {
37  defineTypeNameAndDebug(VoFSolidificationMeltingSource, 0);
38 
40  (
41  option,
42  VoFSolidificationMeltingSource,
43  dictionary
44  );
45  }
46 }
47 
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
51 void Foam::fv::VoFSolidificationMeltingSource::update()
52 {
53  if (curTimeIndex_ == mesh_.time().timeIndex())
54  {
55  return;
56  }
57 
58  if (debug)
59  {
60  Info<< type() << ": " << name_
61  << " - updating solid phase fraction" << endl;
62  }
63 
64  alphaSolid_.oldTime();
65 
66  const twoPhaseMixtureThermo& thermo
67  (
68  mesh_.lookupObject<twoPhaseMixtureThermo>
69  (
71  )
72  );
73 
74  const volScalarField& TVoF = thermo.thermo1().T();
75  const volScalarField CpVoF(thermo.thermo1().Cp());
76  const volScalarField& alphaVoF = thermo.alpha1();
77 
78  forAll(cells_, i)
79  {
80  const label celli = cells_[i];
81 
82  alphaSolid_[celli] = min
83  (
84  relax_*alphaVoF[celli]*alphaSolidT_->value(TVoF[celli])
85  + (1 - relax_)*alphaSolid_[celli],
86  alphaVoF[celli]
87  );
88  }
89 
90  alphaSolid_.correctBoundaryConditions();
91 
92  curTimeIndex_ = mesh_.time().timeIndex();
93 }
94 
95 
96 Foam::word Foam::fv::VoFSolidificationMeltingSource::alphaSolidName() const
97 {
98  const twoPhaseMixtureThermo& thermo
99  (
100  mesh_.lookupObject<twoPhaseMixtureThermo>
101  (
103  )
104  );
105 
106  const volScalarField& alphaVoF = thermo.alpha1();
107 
108  return IOobject::groupName(alphaVoF.name(), "solid");
109 }
110 
111 
112 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
113 
115 (
116  const word& sourceName,
117  const word& modelType,
118  const dictionary& dict,
119  const fvMesh& mesh
120 )
121 :
122  cellSetOption(sourceName, modelType, dict, mesh),
123  alphaSolidT_(Function1<scalar>::New("alphaSolidT", coeffs_)),
124  L_("L", dimEnergy/dimMass, coeffs_),
125  relax_(coeffs_.lookupOrDefault("relax", 0.9)),
126  Cu_(coeffs_.lookupOrDefault<scalar>("Cu", 100000)),
127  q_(coeffs_.lookupOrDefault("q", 0.001)),
128  alphaSolid_
129  (
130  IOobject
131  (
132  alphaSolidName(),
133  mesh.time().timeName(),
134  mesh,
135  IOobject::READ_IF_PRESENT,
136  IOobject::AUTO_WRITE
137  ),
138  mesh,
140  zeroGradientFvPatchScalarField::typeName
141  ),
142  curTimeIndex_(-1)
143 {
144  fieldNames_.setSize(2);
145  fieldNames_[0] = "U";
146  fieldNames_[1] = "T";
147  applied_.setSize(fieldNames_.size(), false);
148 }
149 
150 
151 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
152 
154 (
155  fvMatrix<scalar>& eqn,
156  const label fieldi
157 )
158 {
159  apply(geometricOneField(), eqn);
160 }
161 
162 
164 (
165  const volScalarField& rho,
166  fvMatrix<scalar>& eqn,
167  const label fieldi
168 )
169 {
170  apply(rho, eqn);
171 }
172 
173 
175 (
176  fvMatrix<vector>& eqn,
177  const label fieldi
178 )
179 {
180  if (debug)
181  {
182  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
183  }
184 
185  update();
186 
187  scalarField& Sp = eqn.diag();
188  const scalarField& V = mesh_.V();
189 
190  forAll(cells_, i)
191  {
192  const label celli = cells_[i];
193  const scalar Vc = V[celli];
194  const scalar alphaFluid = 1 - alphaSolid_[celli];
195 
196  const scalar S = Cu_*sqr(1 - alphaFluid)/(pow3(alphaFluid) + q_);
197 
198  Sp[celli] -= Vc*S;
199  }
200 }
201 
202 
204 (
205  const volScalarField& rho,
206  fvMatrix<vector>& eqn,
207  const label fieldi
208 )
209 {
210  // Momentum source uses a Boussinesq approximation - redirect
211  addSup(eqn, fieldi);
212 }
213 
214 
215 // ************************************************************************* //
defineTypeNameAndDebug(option, 0)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
addToRunTimeSelectionTable(option, fixedTemperatureConstraint, dictionary)
rhoReactionThermo & thermo
Definition: createFields.H:28
Macros for easy insertion into run-time selection tables.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
const word dictName() const
Return the local dictionary name (final part of scoped name)
Definition: dictionary.H:123
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
static word groupName(Name name, const word &group)
word timeName
Definition: getTimeIndex.H:3
virtual void addSup(fvMatrix< scalar > &eqn, const label fieldi)
Add explicit contribution to enthalpy equation.
VoFSolidificationMeltingSource(const word &sourceName, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
const dimensionSet dimEnergy
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
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
messageStream Info
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
Namespace for OpenFOAM.