solidParticleCloud.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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "solidParticleCloud.H"
29#include "fvMesh.H"
30#include "volFields.H"
32
33// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34
36(
37 const fvMesh& mesh,
38 const word& cloudName,
39 bool readFields
40)
41:
43 mesh_(mesh),
44 particleProperties_
45 (
47 (
48 "particleProperties",
49 mesh_.time().constant(),
50 mesh_,
51 IOobject::MUST_READ_IF_MODIFIED,
52 IOobject::NO_WRITE
53 )
54 ),
55 rhop_(dimensionedScalar("rhop", particleProperties_).value()),
56 e_(dimensionedScalar("e", particleProperties_).value()),
57 mu_(dimensionedScalar("mu", particleProperties_).value())
58{
59 if (readFields)
60 {
62 }
63}
64
65
66// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
67
69{
70 const volScalarField& rho = mesh_.lookupObject<const volScalarField>("rho");
71 const volVectorField& U = mesh_.lookupObject<const volVectorField>("U");
72 const volScalarField& nu = mesh_.lookupObject<const volScalarField>("nu");
73
77
79 td(*this, rhoInterp, UInterp, nuInterp, g.value());
80
81 Cloud<solidParticle>::move(*this, td, mesh_.time().deltaTValue());
82}
83
84
85// ************************************************************************* //
const uniformDimensionedVectorField & g
Base cloud calls templated on particle type.
Definition: Cloud.H:68
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const Type & value() const
Return const reference to value.
virtual void move()=0
Reads fields from the time directories and adds them to the mesh database for further post-processing...
Definition: readFields.H:158
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Given cell centre values and point (vertex) values decompose into tetrahedra and linear interpolate w...
constant condensation/saturation model.
A Cloud of solid particles.
Class used to pass tracking data to the trackToFace function.
Definition: solidParticle.H:86
Simple solid spherical particle class with one-way coupling with the continuous phase.
Definition: solidParticle.H:68
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
dynamicFvMesh & mesh
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const wordHashSet &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type.
volScalarField & nu
const word cloudName(propsDict.get< word >("cloud"))