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 BeetstraDrag::setForce() const
{
    scale_=cg();
    Info << "BeetstraDrag using scale from liggghts cg = " << scale_ << endl;

    // get viscosity field
    #ifdef comp
        const volScalarField nufField = particleCloud_.turbulence().mu()/rho_;
    #else
        const volScalarField& nufField = particleCloud_.turbulence().nu();
    #endif

    vector position(0,0,0);
    scalar voidfraction(1);
    vector Ufluid(0,0,0);
    vector drag(0,0,0);
    label cellI=0;
    scalar beta(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 Vs(0);
	scalar localPhiP(0);
	scalar filterDragPrefactor(1.0);
	scalar cCorrParcelSize_(1.0) ;
	scalar vCell(0);
	scalar FfFilterPrime(1);
	
    scalar F0=0.;
    scalar G0=0.;

    interpolationCellPoint<scalar> voidfractionInterpolator_(voidfraction_);
    interpolationCellPoint<vector> UInterpolator_(U_);

    #include "setupProbeModel.H"

    for(int index = 0;index <  particleCloud_.numberOfParticles(); index++)
    {
        //if(mask[index][0])
        //{
            cellI = particleCloud_.cellIDs()[index][0];
            drag = vector(0,0,0);
            Ufluid= vector(0,0,0);
            
            if (cellI > -1) // particle Found
            {
                if(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.10) voidfraction = 0.10;
                }
                else
                {
					voidfraction = voidfraction_[cellI];
                    Ufluid = U_[cellI];
                }
           
                Us = particleCloud_.velocity(index);
                Ur = Ufluid-Us;
                ds = 2*particleCloud_.radius(index);
                dPrim_ = ds/scale_;
                nuf = nufField[cellI];
                rho = rho_[cellI];

                magUr = mag(Ur);
				Rep = 0;
                Vs = ds*ds*ds*M_PI/6;
                localPhiP = 1-voidfraction+SMALL;
                vCell = U_.mesh().V()[cellI];
                
                 if (magUr > 0)
                {
                    // calc particle Re Nr
                    Rep = dPrim_*voidfraction*magUr/(nuf+SMALL);


// 3 - C - 1 - In this section we could insert parameters
//             or functions
//             that come from a database

                    // calc model coefficient F0
                    F0 = 10.f * localPhiP / voidfraction / voidfraction
                       + voidfraction * voidfraction * ( 1.0+1.5*sqrt(localPhiP) );

                    // calc model coefficient G0
                    G0 =  0.01720833 * Rep / voidfraction / voidfraction  //0.0172083 = 0.413/24
                              *  ( 1.0 / voidfraction + 3.0 * localPhiP  * voidfraction + 8.4 
                                                            * powf(Rep,-0.343) ) 
                              /  ( 1.0 + powf( 10., 3.0* localPhiP ) 
                                       * powf( Rep,-(1.0+4.0*localPhiP)/2.0 ) );

                    // calc model coefficient beta
                    beta =  18.*nuf*rho/(dPrim_*dPrim_)
                              *voidfraction
                              *(F0 + G0);

                    //Apply filtered drag correction
                    if(useFilteredDragModel_)
                    {
                        FfFilterPrime =  FfFilterFunc(
                                                    filtDragParamsLChar2_,
                                                    vCell,
                                                    0
                                                   );
                 
                        filterDragPrefactor = 1 - fFuncFilteredDrag(FfFilterPrime, localPhiP)
                                                * hFuncFilteredDrag(localPhiP);
                                           
                        beta *= filterDragPrefactor;
                    }
                    if(useParcelSizeDependentFilteredDrag_) //Apply filtered drag correction
                    {
                        scalar dParceldPrim = scale_;
                        cCorrParcelSize_ = 
                               cCorrFunctionFilteredDrag(
                                          filtDragParamsK_, 
                                          filtDragParamsALimit_,
                                          filtDragParamsAExponent_,
                                          localPhiP,
                                          dParceldPrim
                                       );                                          
                        beta *= cCorrParcelSize_;
                    }

                    // calc particle's drag
                    drag = Vs * beta * Ur;

                    if (modelType_=="B")
                        drag /= voidfraction;
                }

// 3 - C - 2 - This is a standardized debug and report section

                if( verbose_ ) //&& index>100 && index < 105)
                {
                    Info << "index / cellI = " << index << tab << cellI << endl;
                    Info << "position = " << particleCloud_.position(index) << endl;
                    Info << "Us = " << Us << endl;
                    Info << "Ur = " << Ur << endl;
                    Info << "ds = " << ds << endl;
                    Info << "rho = " << rho << endl;
                    Info << "nuf = " << nuf << endl;
                    Info << "voidfraction = " << voidfraction << endl;
                    Info << "voidfraction_[cellI]: " << voidfraction_[cellI] << endl;
                    Info << "Rep = " << Rep << endl;
                    Info << "filterDragPrefactor = " << filterDragPrefactor << endl;
                    Info << "fFuncFilteredDrag: " << fFuncFilteredDrag(FfFilterPrime, localPhiP) << endl;
                    Info << "hFuncFilteredDrag: " << hFuncFilteredDrag(localPhiP) << endl;
                    Info << "cCorrParcelSize:   " << cCorrParcelSize_ << endl;
                    Info << "F0 / G0: " << F0 << tab << G0 << endl;
                    Info << "beta:   " << beta << endl;
                    Info << "drag = " << drag << endl;

                    Info << "\nBeetstra drag model settings: treatExplicit " << treatExplicit_ << tab
                                    << "verbose: " << verbose_ << tab
                                    << "dPrim: " << dPrim_  << tab
                                    << "interpolation: " << interpolation_  << tab
                                    << "filteredDragM: " << useFilteredDragModel_   << tab
                                    << "parcelSizeDepDrag: " << useParcelSizeDependentFilteredDrag_
                                    << endl << 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(beta);
                    sValues.append(voidfraction);
                    sValues.append(filterDragPrefactor);
                    particleCloud_.probeM().writeProbe(index, sValues, vValues);
                }
            }

            // set force on particle
            if(treatExplicit_) for(int j=0;j<3;j++) expForces()[index][j] += drag[j];
            else  for(int j=0;j<3;j++) impForces()[index][j] += drag[j];

            // set Cd
            if(implDEM_)
            {
                for(int j=0;j<3;j++) fluidVel()[index][j]=Ufluid[j];

                if (modelType_=="B")
                    Cds()[index][0] = Vs*beta/voidfraction;
                else
                    Cds()[index][0] = Vs*beta;

            }else{
                for(int j=0;j<3;j++) DEMForces()[index][j] += drag[j];
            }
        //}
    }
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
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_;
        }
    }
Beispiel #5
0
void GidaspowDrag::setForce() const
{
    if (scaleDia_ > 1)
        Info << "Gidaspow using scale = " << scaleDia_ << endl;
    else if (particleCloud_.cg() > 1){
        scaleDia_=particleCloud_.cg();
        Info << "Gidaspow using scale from liggghts cg = " << scaleDia_ << endl;
    }

    // get viscosity field
    #ifdef comp
        const volScalarField nufField = particleCloud_.turbulence().mu() / rho_;
    #else
        const volScalarField& nufField = particleCloud_.turbulence().nu();
    #endif

    vector position(0,0,0);
    scalar voidfraction(1);
    vector Ufluid(0,0,0);
    vector drag(0,0,0);
    label cellI=0;

    vector Us(0,0,0);
    vector Uturb(0,0,0);
    vector Ur(0,0,0);
    scalar ds(0);
    scalar nuf(0);
    scalar rho(0);
    scalar magUr(0);
    scalar Rep(0);
    scalar Vs(0);
    scalar localPhiP(0);

    scalar CdMagUrLag(0);       //Cd of the very particle
    scalar KslLag(0);           //momentum exchange of the very particle (per unit volume)
    scalar betaP(0);             //momentum exchange of the very particle

    interpolationCellPoint<scalar> voidfractionInterpolator_(voidfraction_);
    interpolationCellPoint<vector> UInterpolator_(U_);

    #include "setupProbeModel.H"

    for(int index = 0;index <  particleCloud_.numberOfParticles(); ++index)
    {
        //if(mask[index][0])
        //{
            cellI = particleCloud_.cellIDs()[index][0];
            drag = vector(0,0,0);
            betaP = 0;
            Vs = 0;
            Ufluid =vector(0,0,0);
            voidfraction=0;

            if (cellI > -1) // particle Found
            {
                position     = particleCloud_.position(index);
                if(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.0;
                    if(voidfraction<0.10) voidfraction = 0.10;
                }
                else
                {
                    voidfraction = voidfraction_[cellI];
                    Ufluid = U_[cellI];                    
                }

                Us = particleCloud_.velocity(index);                    
                Ur = Ufluid-Us;
                
                //accounting for turbulent dispersion
                if(particleCloud_.dispersionM().isActive())
                {
                    Uturb=particleCloud_.fluidTurbVel(index);
                    Ur=(Ufluid+Uturb)-Us;                                        
                }
                                        
                magUr = mag(Ur);
                ds = 2*particleCloud_.radius(index)*phi_;
                rho = rho_[cellI];
                nuf = nufField[cellI];

                Rep=0.0;
                localPhiP = 1.0f-voidfraction+SMALL;
                Vs = ds*ds*ds*M_PI/6;

                //Compute specific drag coefficient (i.e., Force per unit slip velocity and per m³ SUSPENSION)
                //Wen and Yu, 1966
                if(voidfraction > 0.8) //dilute
                {
                    Rep=ds/scaleDia_*voidfraction*magUr/nuf;
                    CdMagUrLag = (24.0*nuf/(ds/scaleDia_*voidfraction)) //1/magUr missing here, but compensated in expression for KslLag!
                                 *(scalar(1)+0.15*Foam::pow(Rep, 0.687));

                    KslLag = 0.75*(
                                            rho*localPhiP*voidfraction*CdMagUrLag
                                          /
                                            (ds/scaleDia_*Foam::pow(voidfraction,2.65))
                                          );
                }
                //Ergun, 1952
                else  //dense
                {
                    KslLag = (150*Foam::pow(localPhiP,2)*nuf*rho)/
                             (voidfraction*ds/scaleDia_*ds/scaleDia_+SMALL)
                            +
                             (1.75*(localPhiP) * magUr * rho)/
                             ((ds/scaleDia_));
                }

                // calc particle's drag coefficient (i.e., Force per unit slip velocity and per m³ PARTICLE)
                betaP = KslLag / localPhiP;

                // calc particle's drag
                drag = Vs * betaP * Ur * scaleDrag_;

                if (modelType_=="B")
                    drag /= voidfraction;

                if(verbose_ && index >=0 && index <1)
                {
                    Pout << " "<< endl;
                    Pout << "Gidaspow drag force verbose: "  << endl;
                    
                    
                    Pout << "magUr = " << magUr << endl;
                    Pout << "localPhiP = " << localPhiP << endl;
                    Pout << "CdMagUrLag = " << CdMagUrLag << endl;
                    
                    Pout << "treatExplicit_ = " << treatExplicit_ << endl;
                    Pout << "implDEM_ = " << implDEM_ << endl;
                    Pout << "modelType_ = " << modelType_ << endl;                    
                    
                    Pout << "KslLag = " << KslLag << endl;                    
                                        
                    Pout << "cellI = " << cellI << endl;
                    Pout << "index = " << index << endl;                    
                    
                    Pout << "Ufluid = " << Ufluid << endl; 
                    Pout << "Uturb = " << Uturb << endl; 
                    Pout << "Us = " << Us << endl; 
                    
                    Pout << "Ur = " << Ur << endl;
                    Pout << "Vs = " << Vs << endl;                    
                    Pout << "ds = " << ds << endl;
                    Pout << "ds/scale = " << ds/scaleDia_ << endl;
                    Pout << "phi = " << phi_ << endl;
                    Pout << "rho = " << rho << endl;
                    Pout << "nuf = " << nuf << endl;
                    Pout << "voidfraction = " << voidfraction << endl;
                    Pout << "Rep = " << Rep << endl;
                    
                    Pout << "localPhiP = " << localPhiP << endl;
                    Pout << "betaP = " << betaP << endl;
                    
                    Pout << "drag = " << drag << endl;                    
                    Pout << "position = " << position << endl;
                    Pout << " "<< 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(betaP);
                    sValues.append(voidfraction);
                    particleCloud_.probeM().writeProbe(index, sValues, vValues);
                }
            }

            // set force on particle
            //treatExplicit_ = false by default
            if(treatExplicit_) for(int j=0;j<3;j++) expForces()[index][j] += drag[j];
            else  for(int j=0;j<3;j++) impForces()[index][j] += drag[j];

            // set Cd
            //implDEM_ = false by default
            if(implDEM_)
            {
                for(int j=0;j<3;j++) fluidVel()[index][j]=Ufluid[j];

                if (modelType_=="B" && cellI > -1)
                    Cds()[index][0] = Vs*betaP/voidfraction*scaleDrag_;
                else
                    Cds()[index][0] = Vs*betaP*scaleDrag_;

            }else{
                for(int j=0;j<3;j++) DEMForces()[index][j] += drag[j];
            }

        //}// end if mask
    }// end loop particles
}
void KochHillDrag::setForce
(
    double** const& mask,
    double**& impForces,
    double**& expForces,
    double**& DEMForces
) const
{
    // get viscosity field
    #ifdef comp
        const volScalarField nufField = particleCloud_.turbulence().mu()/rho_;
    #else
        const volScalarField& nufField = particleCloud_.turbulence().nu();
    #endif

    vector position(0,0,0);
    scalar voidfraction(1);
    vector Ufluid(0,0,0);
    vector drag(0,0,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 Vs(0);
	scalar volumefraction(0);

    interpolationCellPoint<scalar> voidfractionInterpolator_(voidfraction_);
    interpolationCellPoint<vector> UInterpolator_(U_);

    for(int index = 0;index <  particleCloud_.numberOfParticles(); index++)
    {
        if(mask[index][0])
        {

            cellI = particleCloud_.cellIDs()[index][0];
            drag = vector(0,0,0);

            if (cellI > -1) // particle Found
            {
                if(interpolation_)
                {
	                position = particleCloud_.position(index);
                    voidfraction = voidfractionInterpolator_.interpolate(position,cellI);
                    Ufluid = UInterpolator_.interpolate(position,cellI);
                }else
                {
					voidfraction = particleCloud_.voidfraction(index);
                    Ufluid = U_[cellI];
                }

                Us = particleCloud_.velocity(index);
                Ur = Ufluid-Us;
                ds = 2*particleCloud_.radius(index);
                nuf = nufField[cellI];
                rho = rho_[cellI];
                magUr = mag(Ur);
				Rep = 0;
                Vs = ds*ds*ds*M_PI/6;
                volumefraction = 1-voidfraction+SMALL;

                if (magUr > 0)
                {
                    // calc particle Re Nr
                    Rep = ds/scale_*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);

                    // calc model coefficient beta
                    scalar beta = 18*nuf*rho*voidfraction*voidfraction*volumefraction/(ds/scale_*ds/scale_)*
                                  (F0 + 0.5*F3*Rep);

                    // calc particle's drag
                    drag = Vs*beta/volumefraction*Ur;

                    if (modelType_=="B")
                        drag /= voidfraction;
                }

                if(verbose_ && index >=0 && index <2)
                {
                    Info << "index = " << index << endl;
                    Info << "Us = " << Us << endl;
                    Info << "Ur = " << Ur << endl;
                    Info << "ds = " << ds << endl;
                    Info << "ds/scale = " << ds/scale_ << endl;
                    Info << "rho = " << rho << endl;
                    Info << "nuf = " << nuf << endl;
                    Info << "voidfraction = " << voidfraction << endl;
                    Info << "Rep = " << Rep << endl;
                    Info << "drag = " << drag << endl;
                }
            }
            // set force on particle
            if(treatExplicit_) for(int j=0;j<3;j++) expForces[index][j] += drag[j];
            else  for(int j=0;j<3;j++) impForces[index][j] += drag[j];
        }
    }
}
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);
    }
}
void DiFeliceDrag::setForce() const
{
    if (scaleDia_ > 1)
        Info << "DiFeliceDrag using scale = " << scaleDia_ << endl;
    else if (particleCloud_.cg() > 1){
        scaleDia_=particleCloud_.cg();
        Info << "DiFeliceDrag using scale from liggghts cg = " << scaleDia_ << endl;
    }
    
    // get viscosity field
    #ifdef comp
        const volScalarField nufField = particleCloud_.turbulence().mu() / rho_;
    #else
        const volScalarField& nufField = particleCloud_.turbulence().nu();
    #endif

    vector position(0,0,0);
    scalar voidfraction(1);
    vector Ufluid(0,0,0);
    vector drag(0,0,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);

	vector UfluidFluct(0,0,0);
    vector UsFluct(0,0,0);
    vector dragExplicit(0,0,0);
  	scalar dragCoefficient(0);

    interpolationCellPoint<scalar> voidfractionInterpolator_(voidfraction_);
    interpolationCellPoint<vector> UInterpolator_(U_);

    #include "setupProbeModel.H"

    for(int index = 0;index <  particleCloud_.numberOfParticles(); index++)
    {
        //if(mask[index][0])
        //{

            cellI = particleCloud_.cellIDs()[index][0];
            drag = vector(0,0,0);

            if (cellI > -1) // particle Found
            {
                if(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 = rho_[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!

                    //Split forces
                    if(splitImplicitExplicit_)
                    {
                        UfluidFluct  = Ufluid - U_[cellI];
                        UsFluct      = Us     - UsField_[cellI];
                        dragExplicit = dragCoefficient*(UfluidFluct - UsFluct); //explicit part of force
                    }
                }

                if(verbose_ && index >-1 && index <102)
                {
                    Pout << "index = " << index << 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;
                    if(splitImplicitExplicit_)
                    {
                        Pout << "UfluidFluct = " << UfluidFluct << endl;
                        Pout << "UsFluct = " << UsFluct << endl;
                        Pout << "dragExplicit = " << dragExplicit << 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);
                }
            }
            // set force on particle
            if(treatExplicit_) for(int j=0;j<3;j++) expForces()[index][j] += drag[j];
            else   //implicit treatment, taking explicit force contribution into account
            {
               for(int j=0;j<3;j++) 
               { 
                    impForces()[index][j] += drag[j] - dragExplicit[j]; //only consider implicit part!
                    expForces()[index][j] += dragExplicit[j];
               }
            }
            
            for(int j=0;j<3;j++) DEMForces()[index][j] += drag[j];
        }
    //}
}
Beispiel #9
0
	void CalculateDragForce
	(
		cfdemCloud& sm,
		const volScalarField& alpf_,
		const volVectorField& Uf_,
		const volScalarField& rho_,
		const bool& verbose_,
		vectorField& DragForce_,
		const labelListList& particleList_
	)
	{
		
		// get viscosity field
		#ifdef comp
		    const volScalarField nufField = sm.turbulence().mu()/rho_;
		#else
		    const volScalarField& nufField = sm.turbulence().nu();
		#endif

		// Local variables	
		label  cellI(-1);
		vector drag(0,0,0);
		vector Ufluid(0,0,0);
		
		vector position(0,0,0);
		scalar voidfraction(1);
		
		vector Up(0,0,0);
		vector Ur(0,0,0);
		scalar ds(0);
		
		scalar nuf(0);
		scalar rhof(0);
		
		vector WenYuDrag(0,0,0);
		
		interpolationCellPoint<scalar> voidfractionInterpolator_(alpf_);
		interpolationCellPoint<vector> UInterpolator_(Uf_);	
				
		// 
		//_AO_Parallel
		DragForce_.resize(particleList_.size());
					
	    	for(int ii =0; ii < particleList_.size(); ii++)
		{
			int index = particleList_[ii][0];
			cellI = sm.cellIDs()[index][0];
			drag = vector(0,0,0);
			Ufluid = vector(0,0,0);
			WenYuDrag = vector(0,0,0);
			DragForce_[ii] = vector(0,0,0);
			    
			if (cellI > -1) // particle Found
			{
				position = sm.position(index);
			
				if ( alpf_[cellI] > 1. ) Pout << " voidfraction > 1 " << alpf_[cellI] << endl;
			
				voidfraction = voidfractionInterpolator_.interpolate(position,cellI);
				Ufluid = UInterpolator_.interpolate(position,cellI);
			
				if ( voidfraction > 1. ) 
				{
						Pout << " Int. voidfraction > 1 " << " value= " << voidfraction;
						voidfraction = alpf_[cellI];
						Pout << " mod. value = " << voidfraction << endl;
				}		
				Up = sm.velocity(index);
			
				Ur = Ufluid-Up;				
				ds = 2*sm.radius(index);		
				rhof = rho_[cellI];
				nuf = nufField[cellI];
			
				// Drag force
				WenYuDragForce(Ur,ds,rhof,nuf,voidfraction,WenYuDrag);
						
				if(verbose_ && index <= 1)
				{
					Info << "" << endl;
					Pout << " index = " << index << endl;
					Pout << " position = " << position << endl; 
					Pout << " Up = " << Up << endl;
					Pout << " Ur = " << Ur << endl;
					Pout << " dp = " << ds << endl;
					Pout << " rho = " << rhof << endl;
					Pout << " nuf = " << nuf << endl;
					Pout << " voidfraction = " << voidfraction << endl;
					Pout << " drag = " << WenYuDrag << endl;
					Info << " " << endl;
				}
			}	
			for(int j=0;j<3;j++) DragForce_[ii][j] = WenYuDrag[j];
		}	
	}
Beispiel #10
0
	void EulerianParticleVelocityForce
	(
		cfdemCloud& sm,			
		const fvMesh& mesh,
		volVectorField& Uf_,
		volVectorField&	Up_,
		volScalarField& rho_,
		volScalarField& alpf_,
		volScalarField& Pg_,
		volVectorField& MappedDragForce_,
		const labelListList& particleList_,
		const bool& weighting_
	)
	{		
		// Neighbouring cells
		CPCCellToCellStencil neighbourCells(mesh);
				
		// get viscosity field
		#ifdef comp
		    const volScalarField nufField = sm.turbulence().mu()/rho_;
		#else
		    const volScalarField& nufField = sm.turbulence().nu();
		#endif

		// Gas pressure gradient
		volVectorField gradPg_ = fvc::grad(Pg_);
		interpolationCellPoint<vector> gradPgInterpolator_(gradPg_);

		// Local variables	
		label  cellID(-1);
		vector drag(0,0,0);
		vector Ufluid(0,0,0);
		
		vector position(0,0,0);
		scalar voidfraction(1);
		
		vector Up(0,0,0);
		vector Ur(0,0,0);
		scalar ds(0);
		
		scalar nuf(0);
		scalar rhof(0);
		
		vector WenYuDrag(0,0,0);
		
		interpolationCellPoint<scalar> voidfractionInterpolator_(alpf_);
		interpolationCellPoint<vector> UInterpolator_(Uf_);	
		
		scalar dist_s(0);
		scalar sumWeights(0);
		
		scalarField               weightScalar(27,scalar(0.0));
		Field <Field <scalar> >   particleWeights(particleList_.size(),weightScalar);
		
		//Info << " particle size " << particleList_.size() << endl;
		
		// Number of particle in a cell
		scalarField np(mesh.cells().size(),scalar(0));
		
		// Particle volume
		scalar Volp(0);
		vector gradPg_int(0,0,0);
		
		for(int ii = 0; ii < particleList_.size(); ii++)
		{
			int index = particleList_[ii][0];
			
			cellID = sm.cellIDs()[index][0];
			position = sm.position(index);			    

                        Ufluid = UInterpolator_.interpolate(position,cellID); 
			Up = sm.velocity(index);
                        Ur = Ufluid-Up;

                        ds = 2*sm.radius(index);

                        // Calculate WenYu Drag 
                        voidfraction = voidfractionInterpolator_.interpolate(position,cellID);
                        nuf = nufField[cellID];
                        rhof = rho_[cellID];	
                        WenYuDragForce(Ur,ds,rhof,nuf,voidfraction,WenYuDrag);	
    
        		Volp = ds*ds*ds*M_PI/6;
			gradPg_int = gradPgInterpolator_.interpolate(position,cellID);
			
			//if (cellID > -1)  // particle centre is in domain
            		//{
				if(weighting_)
				{
					labelList& cellsNeigh = neighbourCells[cellID];
					sumWeights = 0;
					dist_s = 0;

					//Info << " index = " << index << " ii = " << ii << " cellID = " << cellID << endl;

					forAll(cellsNeigh,jj)
					{
						// Find distances between particle and neighbouring cells					
						dist_s = mag(sm.mesh().C()[cellsNeigh[jj]]-position)/pow(sm.mesh().V()[cellsNeigh[jj]],1./3.);

						if(dist_s <= 0.5)
						{		
							particleWeights[ii][jj] =  1./4.*pow(dist_s,4)-5./8.*pow(dist_s,2)+115./192.;
						}
						else if (dist_s > 0.5 && dist_s <= 1.5)
						{		
							particleWeights[ii][jj] = -1./6.*pow(dist_s,4)+5./6.*pow(dist_s,3)-5./4.*pow(dist_s,2)+5./24.*dist_s+55./96.;
						}
						else if (dist_s > 1.5 && dist_s <= 2.5)
						{		
							particleWeights[ii][jj] =  pow(2.5-dist_s,4)/24.;
						}
						else
						{		
							particleWeights[ii][jj] = 0;
						}

						sumWeights += particleWeights[ii][jj];

					}	

					forAll(cellsNeigh,jj)
					{	
						if ( sumWeights != 0 )
						{
							Up_[cellID] 	         +=  Up*particleWeights[ii][jj]/sumWeights;
							MappedDragForce_[cellID] += (WenYuDrag + Volp * gradPg_int) * particleWeights[ii][jj]/sumWeights;
						}
						else
						{
							Up_[cellID] 		 = vector(0,0,0);
							MappedDragForce_[cellID] = vector(0,0,0);	
						}
					}
				}
				else
				{
void DiFeliceDrag::setForce() const
{
    // get viscosity field
    #ifdef comp
        const volScalarField nufField = particleCloud_.turbulence().mu() / rho_;
    #else
        const volScalarField& nufField = particleCloud_.turbulence().nu();
    #endif

    vector position(0,0,0);
    scalar voidfraction(1);
    vector Ufluid(0,0,0);
    vector drag(0,0,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++)
    {
        //if(mask[index][0])
        //{

            cellI = particleCloud_.cellIDs()[index][0];
            drag = vector(0,0,0);

            if (cellI > -1) // particle Found
            {
                if(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 = rho_[cellI];
                magUr = mag(Ur);
                Rep = 0;
                Cd = 0;

                if (magUr > 0)
                {
                    // calc particle Re Nr
                    Rep = ds*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
                    drag = 0.125*Cd*rho*M_PI*ds*ds*pow(voidfraction,(2-Xi))*magUr*Ur;

                    if (modelType_=="B")
                        drag /= voidfraction;
                }

                if(verbose_ && index >100 && index <102)
                {
                    Pout << "index = " << index << endl;
                    Pout << "Us = " << Us << endl;
                    Pout << "Ur = " << Ur << endl;
                    Pout << "ds = " << ds << endl;
                    Pout << "rho = " << rho << endl;
                    Pout << "nuf = " << nuf << endl;
                    Pout << "voidfraction = " << voidfraction << endl;
                    Pout << "Rep = " << Rep << endl;
                    Pout << "Cd = " << Cd << endl;
                    Pout << "drag = " << 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);
                }
            }
            // set force on particle
            if(treatExplicit_) for(int j=0;j<3;j++) expForces()[index][j] += drag[j];
            else  for(int j=0;j<3;j++) impForces()[index][j] += drag[j];
            for(int j=0;j<3;j++) DEMForces()[index][j] += drag[j];
        }
    //}
}
void KochHillDrag::setForce() const
{
    if (scaleDia_ > 1)
        Info << "KochHill using scale = " << scaleDia_ << endl;
    else if (particleCloud_.cg() > 1){
        scaleDia_=particleCloud_.cg();
        Info << "KochHill using scale from liggghts cg = " << scaleDia_ << endl;
    }

    // get viscosity field
    #ifdef comp
        const volScalarField nufField = particleCloud_.turbulence().mu()/rho_;
    #else
        const volScalarField& nufField = particleCloud_.turbulence().nu();
    #endif

    vector position(0,0,0);
    scalar voidfraction(1);
    vector Ufluid(0,0,0);
    vector drag(0,0,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 Vs(0);
	scalar volumefraction(0);
    scalar betaP(0);

    interpolationCellPoint<scalar> voidfractionInterpolator_(voidfraction_);
    interpolationCellPoint<vector> UInterpolator_(U_);

    #include "setupProbeModel.H"

    for(int index = 0;index <  particleCloud_.numberOfParticles(); index++)
    {
        //if(mask[index][0])
        //{
            cellI = particleCloud_.cellIDs()[index][0];
            drag = vector(0,0,0);
            betaP = 0;
            Vs = 0;
            Ufluid =vector(0,0,0);
            voidfraction=0;

            if (cellI > -1) // particle Found
            {
                if(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];
                }

                Us = particleCloud_.velocity(index);
                Ur = Ufluid-Us;
                ds = particleCloud_.d(index);
                nuf = nufField[cellI];
                rho = rho_[cellI];
                magUr = mag(Ur);
				Rep = 0;
                Vs = ds*ds*ds*M_PI/6;
                volumefraction = 1-voidfraction+SMALL;

                if (magUr > 0)
                {
                    // calc particle Re Nr
                    Rep = ds/scaleDia_*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
                    scalar F = voidfraction * (F0 + 0.5*F3*Rep);

                    // calc drag model coefficient betaP
                    betaP = 18.*nuf*rho/(ds/scaleDia_*ds/scaleDia_)*voidfraction*F;

                    // calc particle's drag
                    drag = Vs*betaP*Ur*scaleDrag_;

                    if (modelType_=="B")
                        drag /= voidfraction;
                }

                if(verbose_ && index >=0 && index <2)
                {
                    Pout << "cellI = " << cellI << endl;
                    Pout << "index = " << index << endl;
                    Pout << "Us = " << Us << endl;
                    Pout << "Ur = " << Ur << endl;
                    Pout << "ds = " << ds << endl;
                    Pout << "ds/scale = " << ds/scaleDia_ << 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"
                    vValues.append(drag);           //first entry must the be the force
                    vValues.append(Ur);
                    sValues.append(Rep);
                    sValues.append(betaP);
                    sValues.append(voidfraction);
                    particleCloud_.probeM().writeProbe(index, sValues, vValues);
                }    
            }
            // set force on particle
            if(treatExplicit_) for(int j=0;j<3;j++) expForces()[index][j] += drag[j];
            else  for(int j=0;j<3;j++) impForces()[index][j] += drag[j];

            // set Cd
            if(implDEM_)
            {
                for(int j=0;j<3;j++) fluidVel()[index][j]=Ufluid[j];

                if (modelType_=="B" && cellI > -1)
                    Cds()[index][0] = Vs*betaP/voidfraction*scaleDrag_;
                else
                    Cds()[index][0] = Vs*betaP*scaleDrag_;

            }else{
                for(int j=0;j<3;j++) DEMForces()[index][j] += drag[j];
            }

        //}
    }
}