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 -------------------------------------------------------------------------------
10 License
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 
33 template<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 
46 template<class CloudType>
48 (
50 )
51 :
53  GName_(re.GName_)
54 {}
55 
56 
57 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
58 
59 template<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  (
73  IOobject
74  (
75  "D",
76  c.time().timeName(),
77  c,
78  IOobject::NO_READ,
79  IOobject::NO_WRITE,
80  IOobject::REGISTER
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 // ************************************************************************* //
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
ParticleDose(const dictionary &dict, CloudType &owner, const word &modelName)
Construct from dictionary.
Definition: ParticleDose.C:28
const dimensionedScalar G
Newtonian constant of gravitation.
const dimensionedScalar re
Classical electron radius: default SI units: [m].
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
dynamicFvMesh & mesh
A class for handling words, derived from Foam::string.
Definition: word.H:63
Calculate the doses absorbed by a particle as the time integral of the particle track along the radia...
Definition: ParticleDose.H:114
const dimensionedScalar c
Speed of light in a vacuum.
virtual void postEvolve(const typename parcelType::trackingData &td)
Post-evolve hook.
Definition: ParticleDose.C:54
const dimensionedScalar & D
volScalarField & p
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:67
Templated cloud function object base class.
A primitive field of type <T> with automated input and output.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127