From 1a7a489c7ed061d117000c13893a172083b06505 Mon Sep 17 00:00:00 2001 From: Will Bainbridge Date: Mon, 12 Feb 2018 09:08:52 +0000 Subject: [PATCH] InjectionModel: Removed the limit on the number of particles per parcel A lower limit of one on the number of particles represented by a single parcel has been removed from the injection models. It may be appropriate to simulate the statistical behaviour of a particulate flow with more lagrangian elements than physical particles. A unity lower limit does not permit this. The limit was, in some situations, also causing the large-diameter end of an injected distribution to be clipped. This resolves bug report https://bugs.openfoam.org/view.php?id=2837 --- .../InjectionModel/InjectionModel.C | 34 ++++++------------- .../InjectionModel/InjectionModel.H | 6 +--- 2 files changed, 11 insertions(+), 29 deletions(-) diff --git a/src/lagrangian/intermediate/submodels/Kinematic/InjectionModel/InjectionModel/InjectionModel.C b/src/lagrangian/intermediate/submodels/Kinematic/InjectionModel/InjectionModel/InjectionModel.C index 434fd37b5..7d1846640 100644 --- a/src/lagrangian/intermediate/submodels/Kinematic/InjectionModel/InjectionModel/InjectionModel.C +++ b/src/lagrangian/intermediate/submodels/Kinematic/InjectionModel/InjectionModel/InjectionModel.C @@ -198,7 +198,7 @@ Foam::scalar Foam::InjectionModel::setNumberOfParticles scalar volumep = pi/6.0*pow3(diameter); scalar volumeTot = massTotal_/rho; - nP = (volumeFraction*volumeTot + delayedVolume_)/(parcels*volumep); + nP = volumeFraction*volumeTot/(parcels*volumep); break; } case pbNumber: @@ -274,8 +274,7 @@ Foam::InjectionModel::InjectionModel(CloudType& owner) parcelBasis_(pbNumber), nParticleFixed_(0.0), time0_(0.0), - timeStep0_(this->template getModelProperty("timeStep0")), - delayedVolume_(0.0) + timeStep0_(this->template getModelProperty("timeStep0")) {} @@ -302,8 +301,7 @@ Foam::InjectionModel::InjectionModel parcelBasis_(pbNumber), nParticleFixed_(0.0), time0_(owner.db().time().value()), - timeStep0_(this->template getModelProperty("timeStep0")), - delayedVolume_(0.0) + timeStep0_(this->template getModelProperty("timeStep0")) { // Provide some info // - also serves to initialise mesh dimensions - needed for parallel runs @@ -369,8 +367,7 @@ Foam::InjectionModel::InjectionModel parcelBasis_(im.parcelBasis_), nParticleFixed_(im.nParticleFixed_), time0_(im.time0_), - timeStep0_(im.timeStep0_), - delayedVolume_(im.delayedVolume_) + timeStep0_(im.timeStep0_) {} @@ -428,8 +425,6 @@ void Foam::InjectionModel::inject if (prepareForNextTimeStep(time, newParcels, newVolumeFraction)) { - scalar delayedVolume = 0; - const scalar trackTime = this->owner().solution().trackTime(); const polyMesh& mesh = this->owner().mesh(); @@ -505,30 +500,21 @@ void Foam::InjectionModel::inject pPtr->rho() ); - if (pPtr->nParticle() >= 1.0) - { - parcelsAdded++; - massAdded += pPtr->nParticle()*pPtr->mass(); + parcelsAdded ++; - if (pPtr->move(cloud, td, dt)) - { - cloud.addParticle(pPtr); - } - else - { - delete pPtr; - } + massAdded += pPtr->nParticle()*pPtr->mass(); + + if (pPtr->move(cloud, td, dt)) + { + cloud.addParticle(pPtr); } else { - delayedVolume += pPtr->nParticle()*pPtr->volume(); delete pPtr; } } } } - - delayedVolume_ = delayedVolume; } postInjectCheck(parcelsAdded, massAdded); diff --git a/src/lagrangian/intermediate/submodels/Kinematic/InjectionModel/InjectionModel/InjectionModel.H b/src/lagrangian/intermediate/submodels/Kinematic/InjectionModel/InjectionModel/InjectionModel.H index 5ee312f2c..c8aefcda1 100644 --- a/src/lagrangian/intermediate/submodels/Kinematic/InjectionModel/InjectionModel/InjectionModel.H +++ b/src/lagrangian/intermediate/submodels/Kinematic/InjectionModel/InjectionModel/InjectionModel.H @@ -2,7 +2,7 @@ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | - \\ / A nd | Copyright (C) 2011-2017 OpenFOAM Foundation + \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation \\/ M anipulation | ------------------------------------------------------------------------------- License @@ -134,10 +134,6 @@ protected: //- Time at start of injection time step [s] scalar timeStep0_; - //- Volume that should have been injected, but would lead to - // less than 1 particle per parcel - scalar delayedVolume_; - // Protected Member Functions