void MeiLift::setForce() const { const volScalarField& nufField = forceSubM(0).nuField(); const volScalarField& rhoField = forceSubM(0).rhoField(); vector position(0,0,0); vector lift(0,0,0); vector Us(0,0,0); vector Ur(0,0,0); scalar magUr(0); scalar magVorticity(0); scalar ds(0); scalar dParcel(0); scalar nuf(0); scalar rho(0); scalar voidfraction(1); scalar Rep(0); scalar Rew(0); scalar Cl(0); scalar Cl_star(0); scalar J_star(0); scalar Omega_eq(0); scalar alphaStar(0); scalar epsilon(0); scalar omega_star(0); vector vorticity(0,0,0); volVectorField vorticity_ = fvc::curl(U_); #include "resetVorticityInterpolator.H" #include "resetUInterpolator.H" #include "setupProbeModel.H" for(int index = 0;index < particleCloud_.numberOfParticles(); index++) { //if(mask[index][0]) //{ lift = vector::zero; label cellI = particleCloud_.cellIDs()[index][0]; if (cellI > -1) // particle Found { Us = particleCloud_.velocity(index); if( forceSubM(0).interpolation() ) { position = particleCloud_.position(index); Ur = UInterpolator_().interpolate(position,cellI) - Us; vorticity = vorticityInterpolator_().interpolate(position,cellI); } else { Ur = U_[cellI] - Us; vorticity=vorticity_[cellI]; } magUr = mag(Ur); magVorticity = mag(vorticity); if (magUr > 0 && magVorticity > 0) { ds = 2*particleCloud_.radius(index); dParcel = ds; forceSubM(0).scaleDia(ds,index); //caution: this fct will scale ds! nuf = nufField[cellI]; rho = rhoField[cellI]; //Update any scalar or vector quantity for (int iFSub=0;iFSub<nrForceSubModels();iFSub++) forceSubM(iFSub).update( index, cellI, ds, nuf, rho, forceSubM(0).verbose() ); // calc dimensionless properties Rep = ds*magUr/nuf; Rew = magVorticity*ds*ds/nuf; alphaStar = magVorticity*ds/magUr/2.0; epsilon = sqrt(2.0*alphaStar /Rep ); omega_star=2.0*alphaStar; //Basic model for the correction to the Saffman lift //Based on McLaughlin (1991) if(epsilon < 0.1) { J_star = -140 *epsilon*epsilon*epsilon*epsilon*epsilon *log( 1./(epsilon*epsilon+SMALL) ); } else if(epsilon > 20) { J_star = 1.0-0.287/(epsilon*epsilon+SMALL); } else { J_star = 0.3 *( 1.0 +tanh( 2.5 * log10(epsilon+0.191) ) ) *( 0.667 +tanh( 6.0 * (epsilon-0.32) ) ); } Cl=J_star*4.11*epsilon; //multiply McLaughlin's correction to the basic Saffman model //Second order terms given by Loth and Dorgan 2009 if(useSecondOrderTerms_) { Omega_eq = omega_star/2.0*(1.0-0.0075*Rew)*(1.0-0.062*sqrt(Rep)-0.001*Rep); Cl_star=1.0-(0.675+0.15*(1.0+tanh(0.28*(omega_star/2.0-2.0))))*tanh(0.18*sqrt(Rep)); Cl += Omega_eq*Cl_star; } lift = 0.125*M_PI *rho *Cl *magUr*Ur^vorticity/magVorticity *ds*ds; //total force on all particles in parcel forceSubM(0).scaleForce(lift,dParcel,index); if (modelType_=="B") { voidfraction = particleCloud_.voidfraction(index); lift /= voidfraction; } } //********************************** //SAMPLING AND VERBOSE OUTOUT if( forceSubM(0).verbose() ) { Pout << "index = " << index << endl; Pout << "Us = " << Us << endl; Pout << "Ur = " << Ur << endl; Pout << "vorticity = " << vorticity << endl; Pout << "dprim = " << ds << endl; Pout << "rho = " << rho << endl; Pout << "nuf = " << nuf << endl; Pout << "Rep = " << Rep << endl; Pout << "Rew = " << Rew << endl; Pout << "alphaStar = " << alphaStar << endl; Pout << "epsilon = " << epsilon << endl; Pout << "J_star = " << J_star << endl; Pout << "lift = " << lift << endl; } //Set value fields and write the probe if(probeIt_) { #include "setupProbeModelfields.H" // Note: for other than ext one could use vValues.append(x) // instead of setSize vValues.setSize(vValues.size()+1, lift); //first entry must the be the force vValues.setSize(vValues.size()+1, Ur); vValues.setSize(vValues.size()+1, vorticity); sValues.setSize(sValues.size()+1, Rep); sValues.setSize(sValues.size()+1, Rew); sValues.setSize(sValues.size()+1, J_star); particleCloud_.probeM().writeProbe(index, sValues, vValues); } // END OF SAMPLING AND VERBOSE OUTOUT //********************************** } // write particle based data to global array forceSubM(0).partToArray(index,lift,vector::zero); //} } }
void KochHillDrag::setForce() const { const volScalarField& nufField = forceSubM(0).nuField(); const volScalarField& rhoField = forceSubM(0).rhoField(); //update force submodels to prepare for loop for (int iFSub=0;iFSub<nrForceSubModels();iFSub++) forceSubM(iFSub).preParticleLoop(forceSubM(iFSub).verbose()); vector position(0,0,0); scalar voidfraction(1); vector Ufluid(0,0,0); vector drag(0,0,0); vector dragExplicit(0,0,0); scalar dragCoefficient(0); label cellI=0; vector Us(0,0,0); vector Ur(0,0,0); scalar ds(0); scalar dParcel(0); scalar nuf(0); scalar rho(0); scalar magUr(0); scalar Rep(0); scalar Vs(0); scalar volumefraction(0); scalar betaP(0); scalar piBySix(M_PI/6); int couplingInterval(particleCloud_.dataExchangeM().couplingInterval()); #include "resetVoidfractionInterpolator.H" #include "resetUInterpolator.H" #include "setupProbeModel.H" for(int index = 0;index < particleCloud_.numberOfParticles(); index++) { cellI = particleCloud_.cellIDs()[index][0]; drag = vector(0,0,0); dragExplicit = vector(0,0,0); dragCoefficient=0; betaP = 0; Vs = 0; Ufluid =vector(0,0,0); voidfraction=0; if (cellI > -1) // particle Found { if(forceSubM(0).interpolation()) { position = particleCloud_.position(index); voidfraction = voidfractionInterpolator_().interpolate(position,cellI); Ufluid = UInterpolator_().interpolate(position,cellI); //Ensure interpolated void fraction to be meaningful // Info << " --> voidfraction: " << voidfraction << endl; if(voidfraction>1.00) voidfraction = 1.00; if(voidfraction<0.40) voidfraction = 0.40; }else { voidfraction = voidfraction_[cellI]; Ufluid = U_[cellI]; } ds = particleCloud_.d(index); dParcel = ds; forceSubM(0).scaleDia(ds); //caution: this fct will scale ds! nuf = nufField[cellI]; rho = rhoField[cellI]; Us = particleCloud_.velocity(index); //Update any scalar or vector quantity for (int iFSub=0;iFSub<nrForceSubModels();iFSub++) forceSubM(iFSub).update( index, cellI, ds, Ufluid, Us, nuf, rho, forceSubM(0).verbose() ); Ur = Ufluid-Us; magUr = mag(Ur); Rep = 0; Vs = ds*ds*ds*piBySix; volumefraction = max(SMALL,min(1-SMALL,1-voidfraction)); if (magUr > 0) { // calc particle Re Nr Rep = ds*voidfraction*magUr/(nuf+SMALL); // calc model coefficient F0 scalar F0=0.; if(volumefraction < 0.4) { F0 = (1. + 3.*sqrt((volumefraction)/2.) + (135./64.)*volumefraction*log(volumefraction) + 16.14*volumefraction )/ (1+0.681*volumefraction-8.48*sqr(volumefraction) +8.16*volumefraction*volumefraction*volumefraction ); } else { F0 = 10*volumefraction/(voidfraction*voidfraction*voidfraction); } // calc model coefficient F3 scalar F3 = 0.0673+0.212*volumefraction+0.0232/pow(voidfraction,5); //Calculate F (the factor 0.5 is introduced, since Koch and Hill, ARFM 33:619–47, use the radius //to define Rep, and we use the particle diameter, see vanBuijtenen et al., CES 66:2368–2376. scalar F = voidfraction * (F0 + 0.5*F3*Rep); // calc drag model coefficient betaP betaP = 18.*nuf*rho/(ds*ds)*voidfraction*F; // calc particle's drag dragCoefficient = Vs*betaP; if (modelType_=="B") dragCoefficient /= voidfraction; forceSubM(0).scaleCoeff(dragCoefficient,dParcel); if(forceSubM(0).switches()[7]) // implForceDEMaccumulated=true { //get drag from the particle itself for (int j=0 ; j<3 ; j++) drag[j] = particleCloud_.fAccs()[index][j]/couplingInterval; }else { drag = dragCoefficient * Ur; // explicitCorr for (int iFSub=0;iFSub<nrForceSubModels();iFSub++) forceSubM(iFSub).explicitCorr( drag, dragExplicit, dragCoefficient, Ufluid, U_[cellI], Us, UsField_[cellI], forceSubM(iFSub).verbose() ); } } if(forceSubM(0).verbose() && index >=0 && index <2) { Pout << "cellI = " << cellI << endl; Pout << "index = " << index << endl; Pout << "Us = " << Us << endl; Pout << "Ur = " << Ur << endl; Pout << "dprim = " << ds << endl; Pout << "rho = " << rho << endl; Pout << "nuf = " << nuf << endl; Pout << "voidfraction = " << voidfraction << endl; Pout << "Rep = " << Rep << endl; Pout << "betaP = " << betaP << endl; Pout << "drag = " << drag << endl; } //Set value fields and write the probe if(probeIt_) { #include "setupProbeModelfields.H" // Note: for other than ext one could use vValues.append(x) // instead of setSize vValues.setSize(vValues.size()+1, drag); //first entry must the be the force vValues.setSize(vValues.size()+1, Ur); sValues.setSize(sValues.size()+1, Rep); sValues.setSize(sValues.size()+1, betaP); sValues.setSize(sValues.size()+1, voidfraction); particleCloud_.probeM().writeProbe(index, sValues, vValues); } } // write particle based data to global array forceSubM(0).partToArray(index,drag,dragExplicit,Ufluid,dragCoefficient); } }