incompressibleDenseParticleFluid.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) 2023 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 "surfaceInterpolate.H"
28 #include "fvMatrix.H"
29 #include "fvcFlux.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace solvers
38 {
41  (
42  solver,
44  fvMesh
45  );
46 }
47 }
48 
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
53 {
54  fluidSolver::correctCoNum(phic);
55 }
56 
57 
59 {
61 }
62 
63 
64 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
65 
68 (
69  fvMesh& mesh
70 )
71 :
72  fluidSolver(mesh),
73 
74  continuousPhaseName
75  (
77  (
78  IOobject
79  (
80  "physicalProperties",
81  runTime.constant(),
82  mesh,
83  IOobject::MUST_READ
84  )
85  ).lookup("continuousPhaseName")
86  ),
87 
88  p_
89  (
90  IOobject
91  (
92  "p",
93  runTime.name(),
94  mesh,
95  IOobject::MUST_READ,
96  IOobject::AUTO_WRITE
97  ),
98  mesh
99  ),
100 
102 
103  g
104  (
105  IOobject
106  (
107  "g",
108  runTime.constant(),
109  mesh,
110  IOobject::MUST_READ,
111  IOobject::NO_WRITE
112  )
113  ),
114 
115  Uc_
116  (
117  IOobject
118  (
119  IOobject::groupName("U", continuousPhaseName),
120  runTime.name(),
121  mesh,
122  IOobject::MUST_READ,
123  IOobject::AUTO_WRITE
124  ),
125  mesh
126  ),
127 
128  phic_
129  (
130  IOobject
131  (
132  IOobject::groupName("phi", continuousPhaseName),
133  runTime.name(),
134  mesh,
135  IOobject::READ_IF_PRESENT,
136  IOobject::AUTO_WRITE
137  ),
138  linearInterpolate(Uc_) & mesh.Sf()
139  ),
140 
141  viscosity(viscosityModel::New(mesh)),
142 
143  rhoc
144  (
145  IOobject
146  (
147  IOobject::groupName("rho", continuousPhaseName),
148  runTime.name(),
149  mesh
150  ),
151  mesh,
153  (
154  IOobject::groupName("rho", continuousPhaseName),
155  dimDensity,
156  viscosity->lookup
157  (
158  IOobject::groupName("rho", continuousPhaseName)
159  )
160  )
161  ),
162 
163  muc
164  (
165  IOobject
166  (
167  IOobject::groupName("mu", continuousPhaseName),
168  runTime.name(),
169  mesh
170  ),
171  rhoc*viscosity->nu()
172  ),
173 
174  alphac_
175  (
176  IOobject
177  (
178  IOobject::groupName("alpha", continuousPhaseName),
179  runTime.name(),
180  mesh,
181  IOobject::READ_IF_PRESENT,
182  IOobject::AUTO_WRITE
183  ),
184  mesh,
186  ),
187 
188  alphacMin
189  (
190  1 - mesh.solution().solverDict(alphac_.name()).lookup<scalar>("max")
191  ),
192 
193  alphacf("alphacf", fvc::interpolate(alphac_)),
194 
195  alphaPhic
196  (
197  IOobject::groupName("alphaPhi", continuousPhaseName),
198  alphacf*phic_
199  ),
200 
201  momentumTransport
202  (
203  phaseIncompressible::momentumTransportModel::New
204  (
205  alphac_,
206  Uc_,
207  alphaPhic,
208  phic_,
209  viscosity
210  )
211  ),
212 
213  clouds
214  (
215  parcelClouds::New(mesh, rhoc, Uc_, muc, g)
216  ),
217 
218  p(p_),
219  Uc(Uc_),
220  phic(phic_),
221  alphac(alphac_)
222 {
224 
225  momentumTransport->validate();
226 
227  // Update alphac from the particle locations
228  alphac_ = max(1 - clouds.theta(), alphacMin);
232 
233  correctCoNum();
234 }
235 
236 
237 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
238 
241 {}
242 
243 
244 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
245 
247 {
248  // Read the controls
249  readControls();
250 
251  if (mesh.dynamic() && !Ucf.valid())
252  {
253  Info<< "Constructing face momentum Ucf" << endl;
254 
255  // Ensure the U BCs are up-to-date before constructing Ucf
256  Uc_.correctBoundaryConditions();
257 
258  Ucf = new surfaceVectorField
259  (
260  IOobject
261  (
262  "Ucf",
263  runTime.name(),
264  mesh,
267  ),
268  fvc::interpolate(Uc)
269  );
270  }
271 
272  // Store the particle positions
273  if (mesh.topoChanging())
274  {
275  clouds.storeGlobalPositions();
276  }
277 
279 
280  correctCoNum();
281 
282  // Update the mesh for topology change, mesh to mesh mapping
283  mesh_.update();
284 }
285 
286 
288 {
289  if (pimple.firstIter())
290  {
291  clouds.evolve();
292 
293  // Update continuous phase volume fraction field
294  alphac_ = max(1 - clouds.theta(), alphacMin);
295  alphac_.correctBoundaryConditions();
296  alphacf = fvc::interpolate(alphac);
297 
298  // ... and continuous phase volumetric flux
299  alphaPhic = alphacf*phic;
300 
301  // Update the continuous phase dynamic viscosity
302  muc = rhoc*viscosity->nu();
303 
304  Fd = new volVectorField
305  (
306  IOobject
307  (
308  "Fd",
309  runTime.name(),
310  mesh
311  ),
312  mesh,
314  zeroGradientFvPatchVectorField::typeName
315  );
316 
317  Dc = new volScalarField
318  (
319  IOobject
320  (
321  "Dc",
322  runTime.name(),
323  mesh
324  ),
325  mesh,
327  zeroGradientFvPatchVectorField::typeName
328  );
329 
330  const fvVectorMatrix cloudSU(clouds.SU(Uc));
331 
332  Fd().primitiveFieldRef() = -cloudSU.source()/mesh.V()/rhoc;
333  Fd().correctBoundaryConditions();
334 
335  Dc().primitiveFieldRef() = -cloudSU.diag()/mesh.V()/rhoc;
336  Dc().correctBoundaryConditions();
337 
338  Dcf = fvc::interpolate(Dc()).ptr();
339 
340  phid =
341  (
342  fvc::flux(Fd())
343  /(Dcf() + dimensionedScalar(Dc().dimensions(), small))
344  ).ptr();
345  }
346 
347  if (pimple.predictTransport())
348  {
349  momentumTransport->predict();
350  }
351 }
352 
353 
355 {}
356 
357 
359 {
360  while (pimple.correct())
361  {
362  correctPressure();
363  }
364 
365  tUcEqn.clear();
366 }
367 
368 
370 {
371  if (pimple.correctTransport())
372  {
373  viscosity->correct();
374  momentumTransport->correct();
375  }
376 }
377 
378 
380 {
381  Fd.clear();
382  Dc.clear();
383  Dcf.clear();
384  phid.clear();
385 }
386 
387 
388 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
void correctBoundaryConditions()
Correct boundary field.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const word & name() const
Return name.
Definition: IOobject.H:310
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Field< Type > & source()
Definition: fvMatrix.H:307
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const fvSchemes & schemes() const
Return the fvSchemes.
Definition: fvMesh.C:1748
virtual void preUpdateMesh()
Prepare for mesh update.
Definition: fvModels.C:260
void setFluxRequired(const word &name) const
Definition: fvSchemes.C:503
scalarField & diag()
Definition: lduMatrix.C:186
Abstract base class for turbulence models (RAS, LES and laminar).
const tmp< volScalarField > theta() const
Return the particle volume fraction field.
List of parcel clouds, with the same interface as an individual parcel cloud. Is a mesh object,...
Definition: parcelClouds.H:60
Provides controls for the pressure reference in closed-volume simulations.
Selector class for relaxation factors, solver type and solution.
Definition: solution.H:51
Abstract base class for run-time selectable region solvers.
Definition: solver.H:55
const fvMesh & mesh
Region mesh.
Definition: solver.H:94
Base solver module for fluid solvers.
Definition: fluidSolver.H:62
void continuityErrors(const surfaceScalarField &phi)
Calculate and print the continuity errors.
Definition: fluidSolver.C:131
Solver module for transient flow of incompressible isothermal fluids coupled with particle clouds inc...
virtual void thermophysicalPredictor()
Construct and solve the energy equation,.
autoPtr< phaseIncompressible::momentumTransportModel > momentumTransport
Pointer to the momentum transport model.
surfaceScalarField alphacf
Interpolated continuous phase-fraction.
virtual void prePredictor()
Called at the start of the PIMPLE loop.
const volScalarField & alphac
Reference continuous phase-fraction.
virtual void postSolve()
Called after the PIMPLE loop at the end of the time-step.
surfaceScalarField alphaPhic
Continuous phase volumetric-flux field.
virtual void pressureCorrector()
Construct and solve the pressure equation in the PISO loop.
virtual void postCorrector()
Correct the momentum and thermophysical transport modelling.
void continuityErrors()
Calculate and print the continuity errors.
void correctCoNum()
Correct the cached Courant numbers.
incompressibleDenseParticleFluid(fvMesh &mesh)
Construct from region mesh.
virtual void preSolve()
Called at the start of the time-step, before the PIMPLE loop.
const surfaceScalarField & phic
Reference to the continuous phase volumetric-flux field.
const volScalarField & p
Reference to the pressure field.
An abstract base class for Newtonian viscosity models.
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
pimpleControl pimple(mesh)
Calculate the face-flux of the given field.
tmp< SurfaceField< typename innerProduct< vector, Type >::type > > flux(const VolField< Type > &vf)
Return the face-flux field obtained from the given volVectorField.
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
addToRunTimeSelectionTable(solver, compressibleMultiphaseVoF, fvMesh)
defineTypeNameAndDebug(compressibleMultiphaseVoF, 0)
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
static const zero Zero
Definition: zero.H:97
VolField< vector > volVectorField
Definition: volFieldsFwd.H:62
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const dimensionSet dimless
messageStream Info
const dimensionSet dimAcceleration
tmp< SurfaceField< Type > > linearInterpolate(const VolField< Type > &vf)
Definition: linear.H:108
const dimensionSet dimTime
const dimensionSet dimDensity
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
SurfaceField< vector > surfaceVectorField
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
dictionary dict
volScalarField & p