diff --git a/src/lagrangian/intermediate/submodels/Kinematic/InjectionModel/InjectionModel/InjectionModel.C b/src/lagrangian/intermediate/submodels/Kinematic/InjectionModel/InjectionModel/InjectionModel.C index 37dd01da999860e04ece9168fbcb2163639b46bf..48286f170db457902df79d5164a57b2a370e8cf5 100644 --- a/src/lagrangian/intermediate/submodels/Kinematic/InjectionModel/InjectionModel/InjectionModel.C +++ b/src/lagrangian/intermediate/submodels/Kinematic/InjectionModel/InjectionModel/InjectionModel.C @@ -3,7 +3,7 @@ \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright (C) 2011-2017 OpenFOAM Foundation - \\/ M anipulation | Copyright (C) 2016 OpenCFD Ltd. + \\/ M anipulation | Copyright (C) 2016-2018 OpenCFD Ltd. ------------------------------------------------------------------------------- License This file is part of OpenFOAM. @@ -275,6 +275,7 @@ Foam::InjectionModel<CloudType>::InjectionModel(CloudType& owner) nParticleFixed_(0.0), time0_(0.0), timeStep0_(this->template getModelProperty<scalar>("timeStep0")), + minParticlesPerParcel_(1), delayedVolume_(0.0), injectorID_(-1) {} @@ -304,6 +305,11 @@ Foam::InjectionModel<CloudType>::InjectionModel nParticleFixed_(0.0), time0_(owner.db().time().value()), timeStep0_(this->template getModelProperty<scalar>("timeStep0")), + minParticlesPerParcel_ + ( + this->coeffDict().template + lookupOrDefault<scalar>("minParticlesPerParcel", 1) + ), delayedVolume_(0.0), injectorID_(this->coeffDict().lookupOrDefault("injectorID", -1)) { @@ -382,6 +388,7 @@ Foam::InjectionModel<CloudType>::InjectionModel nParticleFixed_(im.nParticleFixed_), time0_(im.time0_), timeStep0_(im.timeStep0_), + minParticlesPerParcel_(im.minParticlesPerParcel_), delayedVolume_(im.delayedVolume_), injectorID_(im.injectorID_) {} @@ -517,7 +524,7 @@ void Foam::InjectionModel<CloudType>::inject pPtr->rho() ); - if (pPtr->nParticle() >= 1.0) + if (pPtr->nParticle() >= minParticlesPerParcel_) { parcelsAdded++; massAdded += pPtr->nParticle()*pPtr->mass(); diff --git a/src/lagrangian/intermediate/submodels/Kinematic/InjectionModel/InjectionModel/InjectionModel.H b/src/lagrangian/intermediate/submodels/Kinematic/InjectionModel/InjectionModel/InjectionModel.H index f235314642074e3594f5ef1c8e6b26b6893d5a08..952fa334c5bf044ec4183959a9406f082ca7679e 100644 --- a/src/lagrangian/intermediate/submodels/Kinematic/InjectionModel/InjectionModel/InjectionModel.H +++ b/src/lagrangian/intermediate/submodels/Kinematic/InjectionModel/InjectionModel/InjectionModel.H @@ -3,7 +3,7 @@ \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright (C) 2011-2017 OpenFOAM Foundation - \\/ M anipulation | + \\/ M anipulation | Copyright (C) 2018 OpenCFD Ltd. ------------------------------------------------------------------------------- License This file is part of OpenFOAM. @@ -81,7 +81,7 @@ public: // Enumerations //- Parcel basis representation options - // i.e constant number of particles OR constant mass per parcel + //- i.e constant number of particles OR constant mass per parcel enum parcelBasis { pbNumber, @@ -100,7 +100,7 @@ protected: scalar SOI_; //- Total volume of particles introduced by this injector [m^3] - // - scaled to ensure massTotal is achieved + //- Note: scaled to ensure massTotal is achieved scalar volumeTotal_; //- Total mass to inject [kg] @@ -128,7 +128,7 @@ protected: parcelBasis parcelBasis_; //- nParticle to assign to parcels when the 'fixed' basis - // is selected + //- is selected scalar nParticleFixed_; //- Continuous phase time at start of injection time step [s] @@ -137,8 +137,12 @@ protected: //- Time at start of injection time step [s] scalar timeStep0_; + //- Minimum number of particles used to represent each parcel + //- default = 1 + scalar minParticlesPerParcel_; + //- Volume that should have been injected, but would lead to - // less than 1 particle per parcel + //- less than minParticlesPerParcel_ particle per parcel scalar delayedVolume_; //- Optional injector ID