Loading...
Searching...
No Matches
SuppressionCollision.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) 2013-2017 OpenFOAM Foundation
9 Copyright (C) 2019-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\*---------------------------------------------------------------------------*/
30#include "kinematicCloud.H"
31
32// * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
33
34template<class CloudType>
36(
37 typename CloudType::parcelType::trackingData& td,
38 const scalar dt
39)
40{
41 const kinematicCloud& sc =
42 this->owner().mesh().template
43 lookupObject<kinematicCloud>(suppressionCloud_);
44
45 volScalarField vDotSweep(sc.vDotSweep());
46 dimensionedScalar Dt("dt", dimTime, dt);
47
48 auto tP = volScalarField::New
49 (
51 1.0 - exp(-vDotSweep*Dt)
52 );
53 const auto& P = tP();
54
55 for (typename CloudType::parcelType& p : this->owner())
56 {
57 const label celli = p.cell();
58
59 scalar xx = this->owner().rndGen().template sample01<scalar>();
60
61 if (xx < P[celli])
62 {
63 p.canCombust() = -1;
64 p.typeId() = max(p.typeId(), suppressedParcelType_);
65 }
66 }
67}
68
69
70// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
71
72template<class CloudType>
74(
75 const dictionary& dict,
76 CloudType& owner
77)
78:
79 StochasticCollisionModel<CloudType>(dict, owner, typeName),
80 suppressionCloud_(this->coeffDict().lookup("suppressionCloud")),
81 suppressedParcelType_
82 (
83 this->coeffDict().getOrDefault("suppressedParcelType", -1)
84 )
85{}
86
87
88template<class CloudType>
90(
92)
93:
95 suppressionCloud_(cm.suppressionCloud_),
98
99
100// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
101
102template<class CloudType>
104{}
105
106
107// ************************************************************************* //
const CloudType & owner() const
Return const access to the owner cloud.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
Templated stochastic collision model class.
StochasticCollisionModel(CloudType &owner)
Construct null from owner.
Inter-cloud collision model, whereby the canReact flag can be used to inhibit devolatilisation and su...
SuppressionCollision(const dictionary &dict, CloudType &owner)
Construct from dictionary.
const word suppressionCloud_
Name of cloud used for suppression.
const label suppressedParcelType_
Suppressed parcel type - optional.
virtual ~SuppressionCollision()
Destructor.
virtual void collide(typename CloudType::parcelType::trackingData &td, const scalar dt)
Update the model.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Virtual abstract base class for templated KinematicCloud.
virtual const tmp< volScalarField > vDotSweep() const =0
Volume swept rate of parcels per cell.
Lookup type of boundary radiation properties.
Definition lookup.H:60
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
const dictionary & dict() const
Return const access to the cloud dictionary.
volScalarField & p
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
DSMCCloud< dsmcParcel > CloudType
dimensionedScalar exp(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict