ParticleDose.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) 2022 OpenCFD Ltd.
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 "ParticleDose.H"
29#include "volFields.H"
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
33template<class CloudType>
35(
36 const dictionary& dict,
37 CloudType& owner,
38 const word& modelName
39)
40:
41 CloudFunctionObject<CloudType>(dict, owner, modelName, typeName),
42 GName_(this->coeffDict().getWord("GName"))
43{}
44
45
46template<class CloudType>
48(
50)
51:
53 GName_(re.GName_)
54{}
55
56
57// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
58
59template<class CloudType>
61(
62 const typename parcelType::trackingData& td
63)
64{
65 auto& c = this->owner();
66
67 if (!c.template foundObject<IOField<scalar>>("D"))
68 {
69 auto* DPtr =
71 (
73 (
74 "D",
75 c.time().timeName(),
76 c,
78 )
79 );
80
81 DPtr->store();
82 }
83
84 auto& D = c.template lookupObjectRef<IOField<scalar>>("D");
85
86 D.resize(c.size(), Zero);
87
88 const fvMesh& mesh = this->owner().mesh();
89
90 const auto& G = mesh.lookupObject<volScalarField>(GName_);
91
92 label parceli = 0;
93 forAllConstIters(c, parcelIter)
94 {
95 const parcelType& p = parcelIter();
96
97 D[parceli] += G[p.cell()]*mesh.time().deltaTValue();
98 parceli++;
99 }
100
101 if (c.size() && c.time().writeTime())
102 {
103 D.write();
104 }
105}
106
107
108// ************************************************************************* //
Templated cloud function object base class.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
A primitive field of type <T> with automated input and output.
Definition: IOField.H:58
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
virtual void postEvolve()
Post-evolve hook.
Calculate the doses absorbed by a particle as the time integral of the particle track along the radia...
Definition: ParticleDose.H:121
scalar deltaTValue() const noexcept
Return time step value.
Definition: TimeStateI.H:43
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 Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
const Type & lookupObject(const word &name, const bool recursive=false) const
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
dynamicFvMesh & mesh
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
dictionary dict
const dimensionedScalar & D
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition: stdFoam.H:278