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-2022 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 
28 #include "fvcDdt.H"
31 
32 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  namespace fv
37  {
38  defineTypeNameAndDebug(VoFSolidificationMeltingSource, 0);
39 
41  (
42  fvModel,
43  VoFSolidificationMeltingSource,
44  dictionary
45  );
46  }
47 }
48 
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
52 void Foam::fv::VoFSolidificationMeltingSource::readCoeffs()
53 {
54  alphaSolidT_.reset(Function1<scalar>::New("alphaSolidT", coeffs()).ptr());
55  L_ = dimensionedScalar("L", dimEnergy/dimMass, coeffs());
56  relax_ = coeffs().lookupOrDefault<scalar>("relax", 0.9);
57  Cu_ = coeffs().lookupOrDefault<scalar>("Cu", 100000);
58  q_ = coeffs().lookupOrDefault<scalar>("q", 0.001);
59 }
60 
61 
62 Foam::word Foam::fv::VoFSolidificationMeltingSource::alphaSolidName() const
63 {
64  const compressibleTwoPhaseMixture& thermo
65  (
66  mesh().lookupObject<compressibleTwoPhaseMixture>
67  (
68  "phaseProperties"
69  )
70  );
71 
72  const volScalarField& alphaVoF = thermo.alpha1();
73 
74  return IOobject::groupName(alphaVoF.name(), "solid");
75 }
76 
77 
78 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
79 
81 (
82  const word& name,
83  const word& modelType,
84  const dictionary& dict,
85  const fvMesh& mesh
86 )
87 :
88  fvModel(name, modelType, dict, mesh),
89  set_(coeffs(), mesh),
90  alphaSolidT_(),
91  L_("L", dimEnergy/dimMass, NaN),
92  relax_(NaN),
93  Cu_(NaN),
94  q_(NaN),
95  alphaSolid_
96  (
97  IOobject
98  (
99  alphaSolidName(),
100  mesh.time().timeName(),
101  mesh,
102  IOobject::READ_IF_PRESENT,
103  IOobject::AUTO_WRITE
104  ),
105  mesh,
107  zeroGradientFvPatchScalarField::typeName
108  )
109 {
110  readCoeffs();
111 }
112 
113 
114 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
115 
117 {
118  return wordList({"U", "T"});
119 }
120 
121 
123 (
124  const volScalarField& rho,
125  fvMatrix<scalar>& eqn,
126  const word& fieldName
127 ) const
128 {
129  if (debug)
130  {
131  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
132  }
133 
134  const compressibleTwoPhaseMixture& thermo
135  (
136  mesh().lookupObject<compressibleTwoPhaseMixture>
137  (
138  "phaseProperties"
139  )
140  );
141 
142  const volScalarField CpVoF(thermo.thermo1().Cp());
143 
144  if (eqn.psi().dimensions() == dimTemperature)
145  {
146  eqn += L_/CpVoF*(fvc::ddt(rho, alphaSolid_));
147  }
148  else
149  {
150  eqn += L_*(fvc::ddt(rho, alphaSolid_));
151  }
152 }
153 
154 
156 (
157  const volScalarField& rho,
158  fvMatrix<vector>& eqn,
159  const word& fieldName
160 ) const
161 {
162  if (debug)
163  {
164  Info<< type() << ": applying source to " << eqn.psi().name() << endl;
165  }
166 
167  scalarField& Sp = eqn.diag();
168  const scalarField& V = mesh().V();
169 
170  const labelList& cells = set_.cells();
171 
172  forAll(cells, i)
173  {
174  const label celli = cells[i];
175  const scalar Vc = V[celli];
176  const scalar alphaFluid = 1 - alphaSolid_[celli];
177 
178  const scalar S = Cu_*sqr(1 - alphaFluid)/(pow3(alphaFluid) + q_);
179 
180  Sp[celli] -= Vc*rho[celli]*S;
181  }
182 }
183 
184 
186 {
187  if (debug)
188  {
189  Info<< type() << ": " << name()
190  << " - updating solid phase fraction" << endl;
191  }
192 
193  alphaSolid_.oldTime();
194 
195  const compressibleTwoPhaseMixture& thermo
196  (
197  mesh().lookupObject<compressibleTwoPhaseMixture>
198  (
199  "phaseProperties"
200  )
201  );
202 
203  const volScalarField& TVoF = thermo.thermo1().T();
204  const volScalarField CpVoF(thermo.thermo1().Cp());
205  const volScalarField& alphaVoF = thermo.alpha1();
206 
207  const labelList& cells = set_.cells();
208 
209  forAll(cells, i)
210  {
211  const label celli = cells[i];
212 
213  alphaSolid_[celli] = min
214  (
215  relax_*alphaVoF[celli]*alphaSolidT_->value(TVoF[celli])
216  + (1 - relax_)*alphaSolid_[celli],
217  alphaVoF[celli]
218  );
219  }
220 
221  alphaSolid_.correctBoundaryConditions();
222 }
223 
224 
226 (
227  const polyTopoChangeMap& map
228 )
229 {
230  set_.topoChange(map);
231 }
232 
233 
234 void Foam::fv::VoFSolidificationMeltingSource::mapMesh(const polyMeshMap& map)
235 {
236  set_.mapMesh(map);
237 }
238 
239 
241 (
242  const polyDistributionMap& map
243 )
244 {
245  set_.distribute(map);
246 }
247 
248 
250 {
251  set_.movePoints();
252  return true;
253 }
254 
255 
256 bool Foam::fv::VoFSolidificationMeltingSource::read(const dictionary& dict)
257 {
258  if (fvModel::read(dict))
259  {
260  set_.read(coeffs());
261  readCoeffs();
262  return true;
263  }
264  else
265  {
266  return false;
267  }
268 
269  return false;
270 }
271 
272 
273 // ************************************************************************* //
defineTypeNameAndDebug(fixedTemperatureConstraint, 0)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
fluidReactionThermo & thermo
Definition: createFields.H:28
const word & name() const
Return const access to the source name.
Definition: fvModelI.H:28
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvModel.C:164
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
virtual bool movePoints()
Update for mesh motion.
virtual void correct()
Correct the model.
const dimensionSet dimless
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
fvMesh & mesh
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvModelI.H:34
Macros for easy insertion into run-time selection tables.
VoFSolidificationMeltingSource(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from explicit source name and mesh.
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
Calculate the first temporal derivative.
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const GeometricField< Type, fvPatchField, volMesh > &)
static word groupName(Name name, const word &group)
const dictionary & coeffs() const
Return dictionary.
Definition: fvModelI.H:40
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
word timeName
Definition: getTimeIndex.H:3
List< label > labelList
A List of labels.
Definition: labelList.H:56
virtual bool read(const dictionary &dict)
Read source dictionary.
const dimensionSet dimEnergy
const dimensionSet dimMass
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
addToRunTimeSelectionTable(fvConstraint, fixedTemperatureConstraint, dictionary)
List< word > wordList
A List of words.
Definition: fileName.H:54
dimensionedScalar pow3(const dimensionedScalar &ds)
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
virtual void addSup(const volScalarField &rho, fvMatrix< scalar > &eqn, const word &fieldName) const
Add explicit contribution to compressible enthalpy equation.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
messageStream Info
virtual wordList addSupFields() const
Return the list of fields for which the option adds source term.
const dimensionSet dimTemperature
Namespace for OpenFOAM.
static autoPtr< Function1< scalar > > New(const word &name, const dictionary &dict)
Selector.
Definition: Function1New.C:32