00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #include "InjectionModel.H"
00027 #include <OpenFOAM/mathematicalConstants.H>
00028 #include <meshTools/meshTools.H>
00029
00030
00031
00032 template<class CloudType>
00033 void Foam::InjectionModel<CloudType>::readProps()
00034 {
00035 IOobject propsDictHeader
00036 (
00037 "injectionProperties",
00038 owner_.db().time().timeName(),
00039 "uniform"/cloud::prefix/owner_.name(),
00040 owner_.db(),
00041 IOobject::MUST_READ,
00042 IOobject::NO_WRITE,
00043 false
00044 );
00045
00046 if (propsDictHeader.headerOk())
00047 {
00048 const IOdictionary propsDict(propsDictHeader);
00049
00050 propsDict.readIfPresent("massInjected", massInjected_);
00051 propsDict.readIfPresent("nInjections", nInjections_);
00052 propsDict.readIfPresent("parcelsAddedTotal", parcelsAddedTotal_);
00053 propsDict.readIfPresent("timeStep0", timeStep0_);
00054 }
00055 }
00056
00057
00058 template<class CloudType>
00059 void Foam::InjectionModel<CloudType>::writeProps()
00060 {
00061 if (owner_.db().time().outputTime())
00062 {
00063 IOdictionary propsDict
00064 (
00065 IOobject
00066 (
00067 "injectionProperties",
00068 owner_.db().time().timeName(),
00069 "uniform"/cloud::prefix/owner_.name(),
00070 owner_.db(),
00071 IOobject::NO_READ,
00072 IOobject::NO_WRITE,
00073 false
00074 )
00075 );
00076
00077 propsDict.add("massInjected", massInjected_);
00078 propsDict.add("nInjections", nInjections_);
00079 propsDict.add("parcelsAddedTotal", parcelsAddedTotal_);
00080 propsDict.add("timeStep0", timeStep0_);
00081
00082 propsDict.regIOobject::write();
00083 }
00084 }
00085
00086
00087 template<class CloudType>
00088 void Foam::InjectionModel<CloudType>::prepareForNextTimeStep
00089 (
00090 const scalar time,
00091 label& newParcels,
00092 scalar& newVolume
00093 )
00094 {
00095
00096 newParcels = 0;
00097 newVolume = 0.0;
00098
00099
00100 if (time < SOI_)
00101 {
00102 timeStep0_ = time;
00103 return;
00104 }
00105
00106
00107 scalar t0 = timeStep0_ - SOI_;
00108 scalar t1 = time - SOI_;
00109
00110
00111 newParcels = parcelsToInject(t0, t1);
00112
00113
00114 newVolume = volumeToInject(t0, t1);
00115
00116
00117 if ((newParcels == 0) && (newVolume > 0.0))
00118 {
00119
00120 }
00121 else
00122 {
00123
00124 timeStep0_ = time;
00125 }
00126 }
00127
00128
00129 template<class CloudType>
00130 void Foam::InjectionModel<CloudType>::findCellAtPosition
00131 (
00132 label& cellI,
00133 vector& position
00134 )
00135 {
00136 const vector p0 = position;
00137
00138 bool foundCell = false;
00139
00140 cellI = owner_.mesh().findCell(position);
00141
00142 if (cellI >= 0)
00143 {
00144 const vector& C = owner_.mesh().C()[cellI];
00145 position += SMALL*(C - position);
00146
00147 foundCell = owner_.mesh().pointInCell(position, cellI);
00148 }
00149 reduce(foundCell, orOp<bool>());
00150
00151
00152
00153 if (!foundCell)
00154 {
00155 cellI = owner_.mesh().findNearestCell(position);
00156
00157 if (cellI >= 0)
00158 {
00159 const vector& C = owner_.mesh().C()[cellI];
00160 position += SMALL*(C - position);
00161
00162 foundCell = owner_.mesh().pointInCell(position, cellI);
00163 }
00164 reduce(foundCell, orOp<bool>());
00165 }
00166
00167 if (!foundCell)
00168 {
00169 FatalErrorIn
00170 (
00171 "Foam::InjectionModel<CloudType>::findCellAtPosition"
00172 "("
00173 "label&, "
00174 "vector&"
00175 ")"
00176 )<< "Cannot find parcel injection cell. "
00177 << "Parcel position = " << p0 << nl
00178 << abort(FatalError);
00179 }
00180 }
00181
00182
00183 template<class CloudType>
00184 Foam::scalar Foam::InjectionModel<CloudType>::setNumberOfParticles
00185 (
00186 const label parcels,
00187 const scalar volume,
00188 const scalar diameter,
00189 const scalar rho
00190 )
00191 {
00192 scalar nP = 0.0;
00193 switch (parcelBasis_)
00194 {
00195 case pbMass:
00196 {
00197 nP =
00198 volume/volumeTotal_
00199 *massTotal_/rho
00200 /(parcels*mathematicalConstant::pi/6.0*pow3(diameter));
00201 break;
00202 }
00203 case pbNumber:
00204 {
00205 nP = massTotal_/(rho*volumeTotal_);
00206 break;
00207 }
00208 default:
00209 {
00210 nP = 0.0;
00211 FatalErrorIn
00212 (
00213 "Foam::scalar "
00214 "Foam::InjectionModel<CloudType>::setNumberOfParticles"
00215 "("
00216 "const label, "
00217 "const scalar, "
00218 "const scalar, "
00219 "const scalar"
00220 ")"
00221 )<< "Unknown parcelBasis type" << nl
00222 << exit(FatalError);
00223 }
00224 }
00225
00226 return nP;
00227 }
00228
00229
00230 template<class CloudType>
00231 void Foam::InjectionModel<CloudType>::postInjectCheck
00232 (
00233 const label parcelsAdded,
00234 const scalar massAdded
00235 )
00236 {
00237 const label allParcelsAdded = returnReduce(parcelsAdded, sumOp<label>());
00238
00239 if (allParcelsAdded > 0)
00240 {
00241 Info<< nl
00242 << "--> Cloud: " << owner_.name() << nl
00243 << " Added " << allParcelsAdded << " new parcels" << nl << endl;
00244 }
00245
00246
00247 parcelsAddedTotal_ += allParcelsAdded;
00248
00249
00250 massInjected_ += returnReduce(massAdded, sumOp<scalar>());
00251
00252
00253 time0_ = owner_.db().time().value();
00254
00255
00256 nInjections_++;
00257
00258
00259 writeProps();
00260 }
00261
00262
00263
00264
00265 template<class CloudType>
00266 Foam::InjectionModel<CloudType>::InjectionModel(CloudType& owner)
00267 :
00268 dict_(dictionary::null),
00269 owner_(owner),
00270 coeffDict_(dictionary::null),
00271 SOI_(0.0),
00272 volumeTotal_(0.0),
00273 massTotal_(0.0),
00274 massInjected_(0.0),
00275 nInjections_(0),
00276 parcelsAddedTotal_(0),
00277 parcelBasis_(pbNumber),
00278 time0_(0.0),
00279 timeStep0_(0.0)
00280 {
00281 readProps();
00282 }
00283
00284
00285 template<class CloudType>
00286 Foam::InjectionModel<CloudType>::InjectionModel
00287 (
00288 const dictionary& dict,
00289 CloudType& owner,
00290 const word& type
00291 )
00292 :
00293 dict_(dict),
00294 owner_(owner),
00295 coeffDict_(dict.subDict(type + "Coeffs")),
00296 SOI_(readScalar(coeffDict_.lookup("SOI"))),
00297 volumeTotal_(0.0),
00298 massTotal_(dimensionedScalar(coeffDict_.lookup("massTotal")).value()),
00299 massInjected_(0.0),
00300 nInjections_(0),
00301 parcelsAddedTotal_(0),
00302 parcelBasis_(pbNumber),
00303 time0_(owner.db().time().value()),
00304 timeStep0_(0.0)
00305 {
00306
00307
00308
00309 Info<< " Constructing " << owner.mesh().nGeometricD() << "-D injection"
00310 << endl;
00311
00312 word parcelBasisType = coeffDict_.lookup("parcelBasisType");
00313 if (parcelBasisType == "mass")
00314 {
00315 parcelBasis_ = pbMass;
00316 }
00317 else if (parcelBasisType == "number")
00318 {
00319 parcelBasis_ = pbNumber;
00320 }
00321 else
00322 {
00323 FatalErrorIn
00324 (
00325 "Foam::InjectionModel<CloudType>::InjectionModel"
00326 "("
00327 "const dictionary&, "
00328 "CloudType&, "
00329 "const word&"
00330 ")"
00331 )<< "parcelBasisType must be either 'number' or 'mass'" << nl
00332 << exit(FatalError);
00333 }
00334
00335 readProps();
00336 }
00337
00338
00339
00340
00341 template<class CloudType>
00342 Foam::InjectionModel<CloudType>::~InjectionModel()
00343 {}
00344
00345
00346
00347
00348 template<class CloudType>
00349 template<class TrackData>
00350 void Foam::InjectionModel<CloudType>::inject(TrackData& td)
00351 {
00352 if (!active())
00353 {
00354 return;
00355 }
00356
00357 const scalar time = owner_.db().time().value();
00358 const scalar carrierDt = owner_.db().time().deltaTValue();
00359 const polyMesh& mesh = owner_.mesh();
00360
00361
00362 label parcelsAdded = 0;
00363 scalar massAdded = 0.0;
00364 label newParcels = 0;
00365 scalar newVolume = 0.0;
00366
00367 prepareForNextTimeStep(time, newParcels, newVolume);
00368
00369
00370 const scalar deltaT =
00371 max(0.0, min(carrierDt, min(time - SOI_, timeEnd() - time0_)));
00372
00373
00374 const scalar padTime = max(0.0, SOI_ - time0_);
00375
00376
00377 for (label parcelI=0; parcelI<newParcels; parcelI++)
00378 {
00379 if (validInjection(parcelI))
00380 {
00381
00382 scalar timeInj = time0_ + padTime + deltaT*parcelI/newParcels;
00383
00384
00385 label cellI = -1;
00386 vector pos = vector::zero;
00387 setPositionAndCell(parcelI, newParcels, timeInj, pos, cellI);
00388
00389 if (cellI > -1)
00390 {
00391
00392 scalar dt = time - timeInj;
00393
00394
00395 meshTools::constrainToMeshCentre(mesh, pos);
00396
00397
00398 parcelType* pPtr = new parcelType(td.cloud(), pos, cellI);
00399
00400
00401 setProperties(parcelI, newParcels, timeInj, *pPtr);
00402
00403
00404 td.cloud().checkParcelProperties(*pPtr, dt, fullyDescribed());
00405
00406
00407 meshTools::constrainDirection
00408 (
00409 mesh,
00410 mesh.solutionD(),
00411 pPtr->U()
00412 );
00413
00414
00415 pPtr->nParticle() =
00416 setNumberOfParticles
00417 (
00418 newParcels,
00419 newVolume,
00420 pPtr->d(),
00421 pPtr->rho()
00422 );
00423
00424
00425 td.cloud().addParticle(pPtr);
00426
00427 massAdded += pPtr->nParticle()*pPtr->mass();
00428 parcelsAdded++;
00429 }
00430 }
00431 }
00432
00433 postInjectCheck(parcelsAdded, massAdded);
00434 }
00435
00436
00437
00438
00439 #include "NewInjectionModel.C"
00440
00441