VoidFraction.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 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "VoidFraction.H"
30
31// * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
32
33template<class CloudType>
35{
36 if (thetaPtr_)
37 {
38 thetaPtr_->write();
39 }
40 else
41 {
43 << "thetaPtr not valid" << abort(FatalError);
44 }
45}
46
47
48// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49
50template<class CloudType>
52(
53 const dictionary& dict,
54 CloudType& owner,
55 const word& modelName
56)
57:
58 CloudFunctionObject<CloudType>(dict, owner, modelName, typeName),
59 thetaPtr_(nullptr)
60{}
61
62
63template<class CloudType>
65(
67)
68:
70 thetaPtr_(nullptr)
71{}
72
73
74// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
75
76template<class CloudType>
78(
79 const typename parcelType::trackingData& td
80)
81{
82 if (thetaPtr_)
83 {
84 thetaPtr_->primitiveFieldRef() = 0.0;
85 }
86 else
87 {
88 const fvMesh& mesh = this->owner().mesh();
89
90 thetaPtr_.reset
91 (
93 (
95 (
96 this->owner().name() + "Theta",
97 mesh.time().timeName(),
98 mesh,
101 ),
102 mesh,
104 )
105 );
106 }
107}
108
109
110template<class CloudType>
112(
113 const typename parcelType::trackingData& td
114)
115{
116 volScalarField& theta = thetaPtr_();
117
118 const fvMesh& mesh = this->owner().mesh();
119
120 theta.primitiveFieldRef() /= mesh.time().deltaTValue()*mesh.V();
121
123}
124
125
126template<class CloudType>
128(
129 parcelType& p,
130 const scalar dt,
131 const point&,
132 bool&
133)
134{
135 volScalarField& theta = thetaPtr_();
136
137 theta[p.cell()] += dt*p.nParticle()*p.volume();
138}
139
140
141// ************************************************************************* //
Templated cloud function object base class.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
virtual void postEvolve()
Post-evolve hook.
scalar deltaTValue() const noexcept
Return time step value.
Definition: TimeStateI.H:43
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Creates particle void fraction field on carrier phase.
Definition: VoidFraction.H:59
virtual void postMove(parcelType &p, const scalar dt, const point &position0, bool &keepParticle)
Post-move hook.
Definition: VoidFraction.C:128
virtual void write()
Write post-processing info.
Definition: VoidFraction.C:34
virtual void preEvolve(const typename parcelType::trackingData &td)
Pre-evolve hook.
Definition: VoidFraction.C:78
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Class used to pass data into container.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const DimensionedField< scalar, volMesh > & V() const
Return cell volumes.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
void reset(const label nPoints, const label nInternalFaces, const label nFaces, const label nCells)
Reset this primitiveMesh given the primitive array sizes.
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const dimensionSet dimless
Dimensionless.
errorManip< error > abort(error &err)
Definition: errorManip.H:144
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
dictionary dict