Loading...
Searching...
No Matches
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-2023 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\*---------------------------------------------------------------------------*/
28#include "ParticleDose.H"
29#include "volFields.H"
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
33template<class CloudType>
35(
36 const dictionary& dict,
38 const word& modelName
39)
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 auto* resultPtr = c.template getObjectPtr<IOField<scalar>>("D");
68
69 if (!resultPtr)
70 {
71 resultPtr = new IOField<scalar>
72 (
74 (
75 "D",
76 c.time().timeName(),
77 c,
81 )
82 );
83
84 resultPtr->store();
85 }
86 auto& D = *resultPtr;
87
88 D.resize(c.size(), Zero);
89
90 const fvMesh& mesh = this->owner().mesh();
91
92 const auto& G = mesh.lookupObject<volScalarField>(GName_);
93
94 label parceli = 0;
95 for (const parcelType& p : c)
96 {
97 D[parceli] += G[p.cell()]*mesh.time().deltaTValue();
98 parceli++;
99 }
100
101 const bool haveParticles = c.size();
102 if (c.time().writeTime() && returnReduceOr(haveParticles))
103 {
104 D.write(haveParticles);
105 }
106}
107
108
109// ************************************************************************* //
Templated cloud function object base class.
CloudFunctionObject(CloudType &owner)
Construct null from owner.
const CloudType & owner() const
Return const access to the owner cloud.
particle::trackingData trackingData
Definition DSMCParcel.H:155
A primitive field of type <T> with automated input and output.
Definition IOField.H:53
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
Calculate the doses absorbed by a particle as the time integral of the particle track along the radia...
ParticleDose(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
virtual void postEvolve(const typename parcelType::trackingData &td)
Post-evolve hook.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
const dictionary & dict() const
Return const access to the cloud dictionary.
const word & modelName() const
Return const access to the name of the sub-model.
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
dynamicFvMesh & mesh
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
const dimensionedScalar c
Speed of light in a vacuum.
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
DSMCCloud< dsmcParcel > CloudType
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
const dimensionedScalar & D