Loading...
Searching...
No Matches
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
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_,
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
74 interpolationCellPoint<scalar> rhoInterp(rho);
75 interpolationCellPoint<vector> UInterp(U);
76 interpolationCellPoint<scalar> nuInterp(nu);
77
78 solidParticle::trackingData
79 td(*this, rhoInterp, UInterp, nuInterp, g.value());
80
81 Cloud<solidParticle>::move(*this, td, mesh_.time().deltaTValue());
82}
83
84
85// ************************************************************************* //
const word cloudName(propsDict.get< word >("cloud"))
const uniformDimensionedVectorField & g
Cloud(const polyMesh &mesh, const Foam::zero, const word &cloudName)
Definition Cloud.C:57
void move(TrackCloudType &cloud, typename ParticleType::trackingData &td, const scalar trackTime)
Move the particles.
Definition Cloud.C:137
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Given cell centre values and point (vertex) values decompose into tetrahedra and linear interpolate w...
const Time & time() const noexcept
Return time registry.
void move(const dimensionedVector &g)
Move the particles under the influence of the given gravitational acceleration.
const fvMesh & mesh() const
solidParticleCloud(const solidParticleCloud &)=delete
No copy construct.
Class used to pass tracking data to the trackToFace function.
Simple solid spherical particle class with one-way coupling with the continuous phase.
static void readFields(Cloud< solidParticle > &c)
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
Different types of constants.
GeometricField< vector, fvPatchField, volMesh > volVectorField
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
volScalarField & nu