// Construct from components virtualMassForce::virtualMassForce ( const dictionary& dict, cfdemCloud& sm ) : forceModel(dict,sm), propsDict_(dict.subDict(typeName + "Props")), velFieldName_(propsDict_.lookup("velFieldName")), U_(sm.mesh().lookupObject<volVectorField> (velFieldName_)), phiFieldName_(propsDict_.lookup("phiFieldName")), phi_(sm.mesh().lookupObject<surfaceScalarField> (phiFieldName_)), UrelOld_(NULL), splitUrelCalculation_(false), Cadd_(0.5) { if (particleCloud_.dataExchangeM().maxNumberOfParticles() > 0) { // get memory for 2d array particleCloud_.dataExchangeM().allocateArray(UrelOld_,NOTONCPU,3); } // init force sub model setForceSubModels(propsDict_); // define switches which can be read from dict forceSubM(0).setSwitchesList(0,true); // activate treatExplicit switch forceSubM(0).setSwitchesList(4,true); // activate search for interpolate switch forceSubM(0).readSwitches(); //Extra switches/settings if(propsDict_.found("splitUrelCalculation")) { splitUrelCalculation_ = readBool(propsDict_.lookup("splitUrelCalculation")); if(splitUrelCalculation_) Info << "Virtual mass model: will split the Urel calculation\n"; Info << "WARNING: be sure that LIGGGHTS integration takes ddtv_p implicitly into account! \n"; } if(propsDict_.found("Cadd")) { Cadd_ = readScalar(propsDict_.lookup("Cadd")); Info << "Virtual mass model: using non-standard Cadd = " << Cadd_ << endl; } particleCloud_.checkCG(true); //Append the field names to be probed particleCloud_.probeM().initialize(typeName, "virtualMass.logDat"); particleCloud_.probeM().vectorFields_.append("virtualMassForce"); //first entry must the be the force particleCloud_.probeM().vectorFields_.append("Urel"); particleCloud_.probeM().vectorFields_.append("UrelOld"); particleCloud_.probeM().vectorFields_.append("ddtUrel"); particleCloud_.probeM().scalarFields_.append("Vs"); particleCloud_.probeM().scalarFields_.append("rho"); particleCloud_.probeM().writeHeader(); }
// Construct from components noDrag::noDrag ( const dictionary& dict, cfdemCloud& sm ) : forceModel(dict,sm), propsDict_(dict.subDict(typeName + "Props")), noDEMForce_(propsDict_.lookupOrDefault("noDEMForce",false)), keepCFDForce_(propsDict_.lookupOrDefault("keepCFDForce",false)) { // init force sub model setForceSubModels(propsDict_); // define switches which can be read from dict forceSubM(0).setSwitchesList(3,true); // activate search for verbose switch // read those switches defined above, if provided in dict forceSubM(0).readSwitches(); particleCloud_.checkCG(true); }
// Construct from components interface::interface ( const dictionary& dict, cfdemCloud& sm ) : forceModel(dict,sm), propsDict_(dict.subDict(typeName + "Props")), VOFvoidfractionFieldName_(propsDict_.lookup("VOFvoidfractionFieldName")), alpha_(sm.mesh().lookupObject<volScalarField> (VOFvoidfractionFieldName_)), gradAlphaName_(propsDict_.lookup("gradAlphaName")), gradAlpha_(sm.mesh().lookupObject<volVectorField> (gradAlphaName_)), sigma_(readScalar(propsDict_.lookup("sigma"))), theta_(readScalar(propsDict_.lookup("theta"))), alphaThreshold_(readScalar(propsDict_.lookup("alphaThreshold"))), deltaAlphaIn_(readScalar(propsDict_.lookup("deltaAlphaIn"))), deltaAlphaOut_(readScalar(propsDict_.lookup("deltaAlphaOut"))), C_(1.0), interpolation_(false) { if (propsDict_.found("C")) C_=readScalar(propsDict_.lookup("C")); // init force sub model setForceSubModels(propsDict_); // define switches which can be read from dict forceSubM(0).setSwitchesList(0,true); // activate treatExplicit switch forceSubM(0).setSwitchesList(4,true); // activate search for interpolate switch // read those switches defined above, if provided in dict forceSubM(0).readSwitches(); //for (int iFSub=0;iFSub<nrForceSubModels();iFSub++) // forceSubM(iFSub).readSwitches(); particleCloud_.checkCG(false); }
// Construct from components KochHillDrag::KochHillDrag ( const dictionary& dict, cfdemCloud& sm ) : forceModel(dict,sm), propsDict_(dict.subDict(typeName + "Props")), velFieldName_(propsDict_.lookup("velFieldName")), U_(sm.mesh().lookupObject<volVectorField> (velFieldName_)), voidfractionFieldName_(propsDict_.lookup("voidfractionFieldName")), voidfraction_(sm.mesh().lookupObject<volScalarField> (voidfractionFieldName_)), UsFieldName_(propsDict_.lookupOrDefault("granVelFieldName",word("Us"))), UsField_(sm.mesh().lookupObject<volVectorField> (UsFieldName_)) { // suppress particle probe if (probeIt_ && propsDict_.found("suppressProbe")) probeIt_=!Switch(propsDict_.lookup("suppressProbe")); if(probeIt_) { particleCloud_.probeM().initialize(typeName, typeName+".logDat"); particleCloud_.probeM().vectorFields_.append("dragForce"); //first entry must the be the force particleCloud_.probeM().vectorFields_.append("Urel"); //other are debug particleCloud_.probeM().scalarFields_.append("Rep"); //other are debug particleCloud_.probeM().scalarFields_.append("beta"); //other are debug particleCloud_.probeM().scalarFields_.append("voidfraction"); //other are debug particleCloud_.probeM().writeHeader(); } // init force sub model setForceSubModels(propsDict_); // define switches which can be read from dict forceSubM(0).setSwitchesList(0,true); // activate search for treatExplicit switch forceSubM(0).setSwitchesList(2,true); // activate search for implDEM switch forceSubM(0).setSwitchesList(3,true); // activate search for verbose switch forceSubM(0).setSwitchesList(4,true); // activate search for interpolate switch forceSubM(0).setSwitchesList(7,true); // activate implForceDEMacc switch forceSubM(0).setSwitchesList(8,true); // activate scalarViscosity switch // read those switches defined above, if provided in dict for (int iFSub=0;iFSub<nrForceSubModels();iFSub++) forceSubM(iFSub).readSwitches(); particleCloud_.checkCG(true); }
void noDrag::setForce() const { if(forceSubM(0).verbose()) { Info << "noDrag::setForce:" << endl; Info << "noDEMForce=" << noDEMForce_ << endl; Info << "keepCFDForce=" << keepCFDForce_ << endl; } label cellI=0; for(int index = 0;index < particleCloud_.numberOfParticles(); ++index) { cellI = particleCloud_.cellIDs()[index][0]; if (cellI > -1) // particle Found { //========================== // set force on particle (proposed new code) // write particle based data to global array //forceSubM(0).partToArray(index,drag,dragExplicit); //========================== // set force on particle (old code) if(!keepCFDForce_) { for(int j=0;j<3;j++) { expForces()[index][j] = 0.; impForces()[index][j] = 0.; } } if(noDEMForce_) { for(int j=0;j<3;j++) DEMForces()[index][j] = 0.; if(particleCloud_.impDEMdrag()) { Cds()[index][0] = 0.; for(int j=0;j<3;j++) fluidVel()[index][j] = 0.; } } //========================== } } }
// Construct from components DiFeliceDrag::DiFeliceDrag ( const dictionary& dict, cfdemCloud& sm ) : forceModel(dict,sm), propsDict_(dict.subDict(typeName + "Props")), velFieldName_(propsDict_.lookup("velFieldName")), U_(sm.mesh().lookupObject<volVectorField> (velFieldName_)), voidfractionFieldName_(propsDict_.lookup("voidfractionFieldName")), voidfraction_(sm.mesh().lookupObject<volScalarField> (voidfractionFieldName_)), UsFieldName_(propsDict_.lookup("granVelFieldName")), UsField_(sm.mesh().lookupObject<volVectorField> (UsFieldName_)), scaleDia_(1.), scaleDrag_(1.) { //Append the field names to be probed particleCloud_.probeM().initialize(typeName, "diFeliceDrag.logDat"); particleCloud_.probeM().vectorFields_.append("dragForce"); //first entry must the be the force particleCloud_.probeM().vectorFields_.append("Urel"); //other are debug particleCloud_.probeM().scalarFields_.append("Rep"); //other are debug particleCloud_.probeM().scalarFields_.append("Cd"); //other are debug particleCloud_.probeM().scalarFields_.append("voidfraction"); //other are debug particleCloud_.probeM().writeHeader(); particleCloud_.checkCG(true); if (propsDict_.found("scale")) scaleDia_=scalar(readScalar(propsDict_.lookup("scale"))); if (propsDict_.found("scaleDrag")) scaleDrag_=scalar(readScalar(propsDict_.lookup("scaleDrag"))); // init force sub model setForceSubModels(propsDict_); // define switches which can be read from dict forceSubM(0).setSwitchesList(0,true); // activate treatExplicit switch forceSubM(0).setSwitchesList(2,true); // activate implDEM switch forceSubM(0).setSwitchesList(3,true); // activate search for verbose switch forceSubM(0).setSwitchesList(4,true); // activate search for interpolate switch forceSubM(0).setSwitchesList(8,true); // activate scalarViscosity switch // read those switches defined above, if provided in dict forceSubM(0).readSwitches(); }
// Construct from components MeiLift::MeiLift ( const dictionary& dict, cfdemCloud& sm ) : forceModel(dict,sm), propsDict_(dict.subDict(typeName + "Props")), velFieldName_(propsDict_.lookup("velFieldName")), U_(sm.mesh().lookupObject<volVectorField> (velFieldName_)), useSecondOrderTerms_(false) { if (propsDict_.found("useSecondOrderTerms")) useSecondOrderTerms_=true; // init force sub model setForceSubModels(propsDict_); // define switches which can be read from dict forceSubM(0).setSwitchesList(0,true); // activate treatExplicit switch forceSubM(0).setSwitchesList(3,true); // activate search for verbose switch forceSubM(0).setSwitchesList(4,true); // activate search for interpolate switch forceSubM(0).setSwitchesList(8,true); // activate scalarViscosity switch //set default switches (hard-coded default = false) forceSubM(0).setSwitches(0,true); // enable treatExplicit, otherwise this force would be implicit in slip vel! - IMPORTANT! for (int iFSub=0;iFSub<nrForceSubModels();iFSub++) forceSubM(iFSub).readSwitches(); particleCloud_.checkCG(false); //Append the field names to be probed particleCloud_.probeM().initialize(typeName, typeName+".logDat"); particleCloud_.probeM().vectorFields_.append("liftForce"); //first entry must the be the force particleCloud_.probeM().vectorFields_.append("Urel"); //other are debug particleCloud_.probeM().vectorFields_.append("vorticity"); //other are debug particleCloud_.probeM().scalarFields_.append("Rep"); //other are debug particleCloud_.probeM().scalarFields_.append("Rew"); //other are debug particleCloud_.probeM().scalarFields_.append("J_star"); //other are debug particleCloud_.probeM().writeHeader(); }
void DiFeliceDrag::setForce() const { if (scaleDia_ > 1) Info << typeName << " using scale = " << scaleDia_ << endl; else if (particleCloud_.cg() > 1){ scaleDia_=particleCloud_.cg(); Info << typeName << " using scale from liggghts cg = " << scaleDia_ << endl; } const volScalarField& nufField = forceSubM(0).nuField(); const volScalarField& rhoField = forceSubM(0).rhoField(); 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 nuf(0); scalar rho(0); scalar magUr(0); scalar Rep(0); scalar Cd(0); interpolationCellPoint<scalar> voidfractionInterpolator_(voidfraction_); interpolationCellPoint<vector> UInterpolator_(U_); #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); Ufluid =vector(0,0,0); if (cellI > -1) // particle Found { if(forceSubM(0).interpolation()) { position = particleCloud_.position(index); voidfraction = voidfractionInterpolator_.interpolate(position,cellI); Ufluid = UInterpolator_.interpolate(position,cellI); }else { voidfraction = voidfraction_[cellI]; Ufluid = U_[cellI]; } Us = particleCloud_.velocity(index); Ur = Ufluid-Us; ds = 2*particleCloud_.radius(index); nuf = nufField[cellI]; rho = rhoField[cellI]; magUr = mag(Ur); Rep = 0; Cd = 0; dragCoefficient = 0; if (magUr > 0) { // calc particle Re Nr Rep = ds/scaleDia_*voidfraction*magUr/(nuf+SMALL); // calc fluid drag Coeff Cd = sqr(0.63 + 4.8/sqrt(Rep)); // calc model coefficient Xi scalar Xi = 3.7 - 0.65 * exp(-sqr(1.5-log10(Rep))/2); // calc particle's drag dragCoefficient = 0.125*Cd*rho *M_PI *ds*ds *scaleDia_ *pow(voidfraction,(2-Xi))*magUr *scaleDrag_; if (modelType_=="B") dragCoefficient /= voidfraction; drag = dragCoefficient*Ur; //total drag force! forceSubM(0).explicitCorr(drag,dragExplicit,dragCoefficient,Ufluid,U_[cellI],Us,UsField_[cellI],forceSubM(0).verbose(),index); } if(forceSubM(0).verbose() && index >-1 && index <102) { Pout << "index = " << index << endl; Pout << "scaleDrag_ = " << scaleDrag_ << endl; Pout << "Us = " << Us << endl; Pout << "Ur = " << Ur << endl; Pout << "ds/scale = " << ds/scaleDia_ << endl; Pout << "rho = " << rho << endl; Pout << "nuf = " << nuf << endl; Pout << "voidfraction = " << voidfraction << endl; Pout << "Rep = " << Rep << endl; Pout << "Cd = " << Cd << endl; Pout << "drag (total) = " << drag << endl; } //Set value fields and write the probe if(probeIt_) { #include "setupProbeModelfields.H" vValues.append(drag); //first entry must the be the force vValues.append(Ur); sValues.append(Rep); sValues.append(Cd); sValues.append(voidfraction); particleCloud_.probeM().writeProbe(index, sValues, vValues); } } // write particle based data to global array forceSubM(0).partToArray(index,drag,dragExplicit,Ufluid,dragCoefficient); } }
void virtualMassForce::setForce() const { reAllocArrays(); scalar dt = U_.mesh().time().deltaT().value(); vector position(0,0,0); vector Ufluid(0,0,0); vector Ur(0,0,0); vector DDtU(0,0,0); //Compute extra vfields in case it is needed volVectorField DDtU_(0.0*U_/U_.mesh().time().deltaT()); if(splitUrelCalculation_) DDtU_ = fvc::ddt(U_) + fvc::div(phi_, U_); //Total Derivative of fluid velocity interpolationCellPoint<vector> UInterpolator_( U_); interpolationCellPoint<vector> DDtUInterpolator_(DDtU_); #include "setupProbeModel.H" bool haveUrelOld_(false); for(int index = 0; index < particleCloud_.numberOfParticles(); index++) { vector virtualMassForce(0,0,0); label cellI = particleCloud_.cellIDs()[index][0]; if (cellI > -1) // particle Found { if(forceSubM(0).interpolation()) { position = particleCloud_.position(index); Ufluid = UInterpolator_.interpolate(position,cellI); } else { Ufluid = U_[cellI]; } if(splitUrelCalculation_) //if split, just use total derivative of fluid velocity if(forceSubM(0).interpolation()) { DDtU = DDtUInterpolator_.interpolate(position,cellI); } else { DDtU = DDtU_[cellI]; } else { vector Us = particleCloud_.velocity(index); Ur = Ufluid - Us; } //Check of particle was on this CPU the last step if(UrelOld_[index][0]==NOTONCPU) //use 1. element to indicate that particle was on this CPU the last time step haveUrelOld_ = false; else haveUrelOld_ = true; vector UrelOld(0.,0.,0.); vector ddtUrel(0.,0.,0.); for(int j=0; j<3; j++) { UrelOld[j] = UrelOld_[index][j]; UrelOld_[index][j] = Ur[j]; } if(haveUrelOld_ ) //only compute force if we have old (relative) velocity ddtUrel = (Ur-UrelOld)/dt; if(splitUrelCalculation_) //we can always compute the total derivative in case we split ddtUrel = DDtU; scalar ds = 2*particleCloud_.radius(index); scalar Vs = ds*ds*ds*M_PI/6; scalar rho = forceSubM(0).rhoField()[cellI]; virtualMassForce = Cadd_ * rho * Vs * ddtUrel; if( forceSubM(0).verbose() ) //&& index>100 && index < 105) { Pout << "index / cellI = " << index << tab << cellI << endl; Pout << "position = " << particleCloud_.position(index) << endl; } //Set value fields and write the probe if(probeIt_) { #include "setupProbeModelfields.H" vValues.append(virtualMassForce); //first entry must the be the force vValues.append(Ur); vValues.append(UrelOld); vValues.append(ddtUrel); sValues.append(Vs); sValues.append(rho); particleCloud_.probeM().writeProbe(index, sValues, vValues); } } else //particle not on this CPU UrelOld_[index][0]=NOTONCPU; // write particle based data to global array forceSubM(0).partToArray(index,virtualMassForce,vector::zero); } }
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 interface::setForce() const { #include "resetAlphaInterpolator.H" #include "resetGradAlphaInterpolator.H" for(int index = 0;index < particleCloud_.numberOfParticles(); ++index) { //if(mask[index][0]) //{ // definition of spherical particle scalar dp = 2*particleCloud_.radius(index); vector position = particleCloud_.position(index); label cellI = particleCloud_.cellIDs()[index][0]; if(cellI >-1.0) // particle found on proc domain { scalar alphap; vector magGradAlphap; if(forceSubM(0).interpolation()) // use intepolated values for alpha (normally off!!!) { // make interpolation object for alpha alphap = alphaInterpolator_().interpolate(position,cellI); // make interpolation object for grad(alpha)/|grad(alpha)| vector gradAlphap = gradAlphaInterpolator_().interpolate(position,cellI); magGradAlphap = gradAlphap/max(mag(gradAlphap),SMALL); } else // use cell centered values for alpha { //// for any reason fvc::grad(alpha_) cannot be executed here!? //volVectorField gradAlpha=fvc::grad(alpha_); //volVectorField a = gradAlpha/ // max(mag(gradAlpha),dimensionedScalar("a",dimensionSet(0,-1,0,0,0), SMALL)); //magGradAlphap = a[cellI]; alphap = alpha_[cellI]; volVectorField a = gradAlpha_/ max(mag(gradAlpha_),dimensionedScalar("a",dimensionSet(0,-1,0,0,0), SMALL)); magGradAlphap = a[cellI]; } // Initialize an interfaceForce vector vector interfaceForce = Foam::vector(0,0,0); // Calculate the interfaceForce (range of alphap needed for stability) if ((alphaThreshold_-deltaAlphaIn_) < alphap && alphap < (alphaThreshold_+deltaAlphaOut_)) { Info << "within threshold limits" << endl; // Calculate estimate attachment force as // |6*sigma*sin(pi-theta/2)*sin(pi+theta/2)|*2*pi*dp scalar Fatt = mag( 6 * sigma_ * sin(M_PI - theta_/2) * sin(M_PI + theta_/2) ) * M_PI * dp; interfaceForce = - magGradAlphap * tanh(alphap-alphaThreshold_) * Fatt * C_; } if(true && mag(interfaceForce) > 0) { Info << "dp = " << dp << endl; Info << "position = " << position << endl; Info << "cellI = " << cellI << endl; Info << "alpha cell = " << alpha_[cellI] << endl; Info << "alphap = " << alphap << endl; Info << "magGradAlphap = " << magGradAlphap << endl; Info << "interfaceForce = " << interfaceForce << endl; Info << "mag(interfaceForce) = " << mag(interfaceForce) << endl; } // limit interface force /*scalar rhoP=3000; scalar mP=dp*dp*dp*3.1415/4*rhoP; scalar fMax=5*mP*9.81; if(mag(interfaceForce)>fMax){ interfaceForce /= mag(interfaceForce)/fMax; Info << "interface force is limited to " << interfaceForce << endl; }*/ // write particle based data to global array forceSubM(0).partToArray(index,interfaceForce,vector::zero); } // end if particle found on proc domain //}// end if in mask }// end loop particles Info << "interface::setForce - done" << endl; }
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // void scalarGeneralExchange::manipulateScalarField(volScalarField& explicitEulerSource, volScalarField& implicitEulerSource, int speciesID) const { // reset Scalar field explicitEulerSource.internalField() = 0.0; implicitEulerSource.internalField() = 0.0; if(speciesID>=0 && particleSpeciesValue_[speciesID]<0.0) //skip if species is not active return; //Set the names of the exchange fields word fieldName; word partDatName; word partFluxName; word partTransCoeffName; word partFluidName; scalar transportParameter; // realloc the arrays to pull particle data if(speciesID<0) //this is the temperature - always pull from LIGGGHTS { fieldName = tempFieldName_; partDatName = partTempName_; partFluxName = partHeatFluxName_; partTransCoeffName = partHeatTransCoeffName_; partFluidName = partHeatFluidName_; transportParameter = lambda_; allocateMyArrays(0.0); particleCloud_.dataExchangeM().getData(partDatName,"scalar-atom", partDat_); } else { fieldName = eulerianFieldNames_[speciesID]; partDatName = partSpeciesNames_[speciesID]; partFluxName = partSpeciesFluxNames_[speciesID]; partTransCoeffName = partSpeciesTransCoeffNames_[speciesID]; partFluidName = partSpeciesFluidNames_[speciesID]; transportParameter = DMolecular_[speciesID]; allocateMyArrays(0.0); if(particleSpeciesValue_[speciesID]>ALARGECONCENTRATION) particleCloud_.dataExchangeM().getData(partDatName,"scalar-atom", partDat_); } if (scaleDia_ > 1) Info << typeName << " using scale = " << scaleDia_ << endl; else if (particleCloud_.cg() > 1) { scaleDia_=particleCloud_.cg(); Info << typeName << " using scale from liggghts cg = " << scaleDia_ << endl; } //============================== // get references const volScalarField& voidfraction_(particleCloud_.mesh().lookupObject<volScalarField> (voidfractionFieldName_)); // ref to voidfraction field const volVectorField& U_(particleCloud_.mesh().lookupObject<volVectorField> (velFieldName_)); const volScalarField& fluidScalarField_(particleCloud_.mesh().lookupObject<volScalarField> (fieldName)); // ref to scalar field const volScalarField& nufField = forceSubM(0).nuField(); //============================== // calc La based heat flux vector position(0,0,0); scalar voidfraction(1); vector Ufluid(0,0,0); scalar fluidValue(0); label cellI=0; vector Us(0,0,0); vector Ur(0,0,0); scalar dscaled(0); scalar nuf(0); scalar magUr(0); scalar As(0); scalar Rep(0); scalar Pr(0); scalar sDth(scaleDia_*scaleDia_*scaleDia_); interpolationCellPoint<scalar> voidfractionInterpolator_(voidfraction_); interpolationCellPoint<vector> UInterpolator_(U_); interpolationCellPoint<scalar> fluidScalarFieldInterpolator_(fluidScalarField_); #include "setupProbeModel.H" for(int index = 0;index < particleCloud_.numberOfParticles(); ++index) { cellI = particleCloud_.cellIDs()[index][0]; if(cellI >= 0) { if(forceSubM(0).interpolation()) { position = particleCloud_.position(index); voidfraction = voidfractionInterpolator_.interpolate(position,cellI); Ufluid = UInterpolator_.interpolate(position,cellI); fluidValue = fluidScalarFieldInterpolator_.interpolate(position,cellI); }else { voidfraction = voidfraction_[cellI]; Ufluid = U_[cellI]; fluidValue = fluidScalarField_[cellI]; } // calc relative velocity Us = particleCloud_.velocity(index); Ur = Ufluid-Us; magUr = mag(Ur); dscaled = 2*particleCloud_.radius(index)/scaleDia_; As = dscaled*dscaled*M_PI*sDth; nuf = nufField[cellI]; Rep = dscaled*magUr/nuf; if(speciesID<0) //have temperature Pr = Prandtl_; else Pr = max(SMALL,nuf/transportParameter); //This is Sc for species scalar alpha = transportParameter*(this->*Nusselt)(Rep,Pr,voidfraction)/(dscaled); // calc convective heat flux [W] scalar areaTimesTransferCoefficient = alpha * As; scalar tmpPartFlux = areaTimesTransferCoefficient * (fluidValue - partDat_[index][0]); partDatFlux_[index][0] = tmpPartFlux; // split implicit/explicit contribution forceSubM(0).explicitCorrScalar( partDatTmpImpl_[index][0], partDatTmpExpl_[index][0], areaTimesTransferCoefficient, fluidValue, fluidScalarField_[cellI], partDat_[index][0], forceSubM(0).verbose() ); if(validPartTransCoeff_) partDatTransCoeff_[index][0] = alpha; if(validPartFluid_) partDatFluid_[index][0] = fluidValue; if( forceSubM(0).verbose()) { Pout << "fieldName = " << fieldName << endl; Pout << "partTransCoeffName = " << partTransCoeffName << endl; Pout << "index = " <<index << endl; Pout << "partFlux = " << tmpPartFlux << endl; Pout << "magUr = " << magUr << endl; Pout << "As = " << As << endl; Pout << "r = " << particleCloud_.radius(index) << endl; Pout << "dscaled = " << dscaled << endl; Pout << "nuf = " << nuf << endl; Pout << "Rep = " << Rep << endl; Pout << "Pr/Sc = " << Pr << endl; Pout << "Nup/Shp = " << (this->*Nusselt)(Rep,Pr,voidfraction) << endl; Pout << "partDatTransCoeff: " << partDatTransCoeff_[index][0] << endl; Pout << "voidfraction = " << voidfraction << endl; Pout << "partDat_[index][0] = " << partDat_[index][0] << endl ; Pout << "fluidValue = " << fluidValue << endl ; } //Set value fields and write the probe if(probeIt_) { #include "setupProbeModelfields.H" vValues.append(Ur); sValues.append((this->*Nusselt)(Rep,Pr,voidfraction)); sValues.append(Rep); particleCloud_.probeM().writeProbe(index, sValues, vValues); } } } //Handle explicit and implicit source terms on the Euler side //these are simple summations! particleCloud_.averagingM().setScalarSum ( explicitEulerSource, partDatTmpExpl_, particleCloud_.particleWeights(), NULL ); particleCloud_.averagingM().setScalarSum ( implicitEulerSource, partDatTmpImpl_, particleCloud_.particleWeights(), NULL ); // scale with the cell volume to get (total) volume-specific source explicitEulerSource.internalField() /= -explicitEulerSource.mesh().V(); implicitEulerSource.internalField() /= -implicitEulerSource.mesh().V(); // limit explicit source term scalar explicitEulerSourceInCell; forAll(explicitEulerSource,cellI) { explicitEulerSourceInCell = explicitEulerSource[cellI]; if(mag(explicitEulerSourceInCell) > maxSource_ ) { explicitEulerSource[cellI] = sign(explicitEulerSourceInCell) * maxSource_; } }
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); } }