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 | Copyright (C) 2011 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 
26 #include "solidParticleCloud.H"
27 #include "fvMesh.H"
28 #include "volFields.H"
29 #include "interpolationCellPoint.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
33 Foam::solidParticleCloud::solidParticleCloud
34 (
35  const fvMesh& mesh,
36  const word& cloudName,
37  bool readFields
38 )
39 :
41  mesh_(mesh),
42  particleProperties_
43  (
44  IOobject
45  (
46  "particleProperties",
47  mesh_.time().constant(),
48  mesh_,
49  IOobject::MUST_READ_IF_MODIFIED,
50  IOobject::NO_WRITE
51  )
52  ),
53  rhop_(dimensionedScalar(particleProperties_.lookup("rhop")).value()),
54  e_(dimensionedScalar(particleProperties_.lookup("e")).value()),
55  mu_(dimensionedScalar(particleProperties_.lookup("mu")).value())
56 {
57  if (readFields)
58  {
60  }
61 }
62 
63 
64 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
65 
67 {
68  return true;
69 }
70 
71 
73 {
74  const volScalarField& rho = mesh_.lookupObject<const volScalarField>("rho");
75  const volVectorField& U = mesh_.lookupObject<const volVectorField>("U");
76  const volScalarField& nu = mesh_.lookupObject<const volScalarField>("nu");
77 
78  interpolationCellPoint<scalar> rhoInterp(rho);
80  interpolationCellPoint<scalar> nuInterp(nu);
81 
83  td(*this, rhoInterp, UInterp, nuInterp, g.value());
84 
86 }
87 
88 
89 // ************************************************************************* //
virtual bool hasWallImpactDistance() const
Switch to specify if particles of the cloud can return.
U
Definition: pEqn.H:83
const word cloudName(propsDict.lookup("cloudName"))
const Type & value() const
Return const reference to value.
void move(const dimensionedVector &g)
Move the particles under the influence of the given.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const HashSet< word > &selectedFields, LIFOStack< regIOobject * > &storedObjects)
Read the selected GeometricFields of the specified type.
Definition: ReadFields.C:244
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
void move(TrackData &td, const scalar trackTime)
Move the particles.
Definition: Cloud.C:187
Class used to pass tracking data to the trackToFace function.
Definition: solidParticle.H:86
const dimensionedVector & g
scalar deltaTValue() const
Return time step value.
Definition: TimeStateI.H:41
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
volScalarField & nu
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243