#include "hparticlecandfiller.h"
#include "hcategorymanager.h"
#include "hparticletool.h"
#include "hparticleconstants.h"
#include "hgeomvolume.h"
#include "hgeantdef.h"
#include "hparticledef.h"
#include "hmdcdef.h"
#include "hmdctrackddef.h"
#include "hmdctrackgdef.h"
#include "tofdef.h"
#include "showerdef.h"
#include "richdef.h"
#include "rpcdef.h"
#include "emcdef.h"
#include "hades.h"
#include "hiterator.h"
#include "hcategory.h"
#include "hlinearcategory.h"
#include "hlocation.h"
#include "hruntimedb.h"
#include "hspectrometer.h"
#include "hdetector.h"
#include "hevent.h"
#include "hmetamatch2.h"
#include "hmdcseg.h"
#include "hmdcsegsim.h"
#include "hmdctrkcand.h"
#include "hsplinetrack.h"
#include "hrktrackB.h"
#include "hshowerhitsim.h"
#include "htofhit.h"
#include "htofcluster.h"
#include "hrichhit.h"
#include "hrichhitsim.h"
#include "hrpccluster.h"
#include "hrpcclustersim.h"
#include "hparticlecand.h"
#include "hparticlecandsim.h"
#include "hparticlemdc.h"
#include "hparticlecandfillerpar.h"
#include "hrich700digipar.h"
#include "htofwalkpar.h"
#include "hparticleanglecor.h"
#include "hmdclayer.h"
#include "TVector2.h"
#include <iostream>
#include <iomanip>
#include <stdlib.h>
#include <algorithm>
#include <numeric>
using namespace std;
using namespace Particle;
ClassImp(HParticleCandFiller)
HParticleCandFiller::HParticleCandFiller (const Option_t par[]):
    HReconstructor ()
{
    
    
    initVars ();
    setConditions(par);
}
HParticleCandFiller::HParticleCandFiller (const Text_t * name,
					  const Text_t * title,
					  const Option_t par[]):
HReconstructor (name, title)
{
    
    
    initVars ();
    setConditions(par);
}
HParticleCandFiller::~HParticleCandFiller (void)
{
    
    if(fMetaMatchIter) delete fMetaMatchIter;
}
void HParticleCandFiller::setConditions(const Option_t par[])
{
    
    
    
    
    
    
    
    
    
    
    
    
    
    
    
    TString s = par;
    s.ToUpper();
    fbgoodSeg0   = (strstr(s.Data(), "GOOGSEG0")  != NULL);
    fbgoodSeg1   = (strstr(s.Data(), "GOODSEG1")  != NULL);
    fbgoodMeta   = (strstr(s.Data(), "GOODMETA")  != NULL);
    fbgoodRK     = (strstr(s.Data(), "GOODRK")    != NULL);
    fbgoodLepton = (strstr(s.Data(), "GOODLEPTON")!= NULL);
    fmomSwitch   = (strstr(s.Data(), "KALMAN")    != NULL) ? Particle::kMomKalman : Particle::kMomRK;
    fbIsDebug    = (strstr(s.Data(), "DEBUG")     != NULL);
    fbdoRichAlign= (strstr(s.Data(), "NORICHALIGN")!= NULL)? kFALSE : kTRUE;
    fbdoRichVertexCorr = (strstr(s.Data(), "NORICHVERTEXCORR")!= NULL)? kFALSE : kTRUE;
    fbdoMETAQANorm     = (strstr(s.Data(), "NOMETAQANORM")!= NULL)? kFALSE : kTRUE;
    fbdoMomentumCorr   = (strstr(s.Data(), "NOMOMENTUMCORR")!= NULL)? kFALSE : kTRUE;
    fbdoPathLengthCorr = (strstr(s.Data(), "NOPATHLENGTHCORR")!= NULL)? kFALSE : kTRUE;
    fbdoGeantAcceptance= (strstr(s.Data(), "NOGEANTACCPETANCE")!= NULL)? kFALSE : kTRUE;
    fbnoFake     = (strstr(s.Data(), "NOFAKE")    != NULL);
    if(strstr(s.Data(), "ACCEPTFAKE")!= NULL)fbnoFake=kFALSE;
}
void HParticleCandFiller::initVars ()
{
    
    fCatMdcTrkCand   = NULL;
    fCatMetaMatch    = NULL;
    fCatMdcSeg       = NULL;
    fCatTofHit       = NULL;
    fCatTofCluster   = NULL;
    fCatEmcCluster   = NULL;
    fCatShowerHit    = NULL;
    fCatRichHit      = NULL;
    fCatRpcCluster   = NULL;
    fCatSpline       = NULL;
    fCatRK           = NULL;
    fCatKalman       = NULL;
    fCatParticleCand = NULL;
    fCatGeantKine    = NULL;
    fMetaMatchIter   = NULL;
    fCatParticleCand = NULL;
    fCatParticleDebug= NULL;
    fCatParticleMdc  = NULL;
    fFillerPar       = NULL;
    fTofWalkPar      = NULL;
    fbIsSimulation   = kFALSE;
    fbIsDebug        = kFALSE;
    fbFillMdc        = kFALSE;
    fbgoodSeg0       = kFALSE;
    fbgoodSeg1       = kFALSE;
    fbgoodMeta       = kFALSE;
    fbgoodRK         = kFALSE;
    fbgoodLepton     = kFALSE;
    fmomSwitch       = Particle::kMomRK;
    fsortSwitch      = 0; 
    fbdoRichAlign    = kTRUE;
    fbdoRichVertexCorr=kTRUE;
    fbdoMETAQANorm    =kTRUE;
    fbdoMomentumCorr  =kTRUE;
    fbdoPathLengthCorr=kTRUE;
    fbdoGeantAcceptance=kTRUE;
    fbnoFake          =kFALSE;
    fMinWireGoodTrack=   5;
    fScaleGhostTrack = 0.1;
    fScaleGoodTrack  = 3.0;
    fAngleCloseTrack = 15.;
    all_candidates.resize( 2000, 0 );  
    all_candidates.clear();            
    fRichCorrectionVersion = -1;
}
Bool_t HParticleCandFiller::init (void)
{
    
    HEvent *ev = gHades->getCurrentEvent ();
    if (ev) {
	fCatGeantKine = ev->getCategory (catGeantKine);
	if (fCatGeantKine) { fbIsSimulation = kTRUE; }
	else               { fbIsSimulation = kFALSE;}
	if(!(fCatMetaMatch  = HCategoryManager::getCategory(catMetaMatch ,kFALSE,"catMetaMatch" ))) return kFALSE;
	if(!(fCatMdcTrkCand = HCategoryManager::getCategory(catMdcTrkCand,kFALSE,"catMdcTrkCand"))) return kFALSE;
	if(!(fCatMdcSeg     = HCategoryManager::getCategory(catMdcSeg    ,kFALSE,"catMdcSeg"    ))) return kFALSE;
	fCatTofHit     = HCategoryManager::getCategory(catTofHit      ,kTRUE,"catTofHit");
	fCatTofCluster = HCategoryManager::getCategory(catTofCluster  ,kTRUE,"catTofCluster");
	fCatRichHit    = HCategoryManager::getCategory(catRichHit     ,kTRUE,"catRichHit");
	fCatRpcCluster = HCategoryManager::getCategory(catRpcCluster  ,kTRUE,"catRpcCluster");
	fCatSpline     = HCategoryManager::getCategory(catSplineTrack ,kTRUE,"catSplineTrack");
	if(fmomSwitch==Particle::kMomRK)    fCatRK         = HCategoryManager::getCategory(catRKTrackB    ,kTRUE,"catRKTrackB");
	if(fmomSwitch==Particle::kMomKalman)fCatKalman     = HCategoryManager::getCategory(catKalTrack    ,kTRUE,"catKalTrack");
	fCatShowerHit = HCategoryManager::getCategory(catShowerHit,kTRUE,"catShowerHit");
	fCatEmcCluster= HCategoryManager::getCategory(catEmcCluster,kTRUE,"catEmcCluster");
	if (fbIsSimulation) {
	    fCatParticleCand = HCategoryManager::addCategory(catParticleCand,"HParticleCandSim",5000,"Particle");
	} else {
	    fCatParticleCand = HCategoryManager::addCategory(catParticleCand,"HParticleCand"   ,5000,"Particle");
	}
	if (!fCatParticleCand) { return kFALSE; }
	if (fbIsDebug) {
	    fCatParticleDebug = HCategoryManager::addCategory(catParticleDebug,"candidate",5000,"Particle",kTRUE);
	    if (!fCatParticleDebug) { return kFALSE; }
	}
	if (fbFillMdc) {
	    fCatParticleMdc = HCategoryManager::addCategory(catParticleMdc,"HParticleMdc",5000,"Particle",kTRUE);
	    if (!fCatParticleMdc) { return kFALSE; }
	}
	Int_t beamID = gHades->getBeamTimeID();
	if(beamID == Particle::kUnknownBeam) {
	    Error ("init()", "gHades->getBeamTimeID() == Particle::kUnknownBeam! Not set properly, use Particle::eBeamTime definition (hparticledef.h), RICH correction will depend on it.");
	    HParticleConstants::printSpace("Particle::eBeam");
	    return kFALSE;
	} else {
	    if(beamID == Particle::kApr12 ||
	       beamID == Particle::kJul14 ||
	       beamID == Particle::kAug14 )  fRichCorrectionVersion = 0;
	    else                             fRichCorrectionVersion = 1;
	    Info ("init()", "gHades->getBeamTimeID() == %i !, use Particle::eBeamTime definition (hparticledef.h) to set RICH correction properly.",beamID);
	    HParticleConstants::printSpace("Particle::eBeam");
	}
	fMetaMatchIter = (HIterator*) fCatMetaMatch->MakeIterator("native");
	if(!fMetaMatchIter){
	    Error ("init()", "Retrieve ZERO pointer for MetaMatch iter!");
	    return kFALSE;
	}
	fFillerPar =(HParticleCandFillerPar*) gHades->getRuntimeDb()->getContainer("ParticleCandFillerPar");
	if(!fFillerPar) {
 	    Error ("init()", "Retrieve ZERO pointer for HParticleCandFillerPar!");
	    return kFALSE;
	}
	fRich700DigiPar =(HRich700DigiPar*) gHades->getRuntimeDb()->getContainer("Rich700DigiPar");
	if(!fRich700DigiPar) {
 	    Error ("init()", "Retrieve ZERO pointer for HRich700DigiPar!");
	    return kFALSE;
	}
	fTofWalkPar =(HTofWalkPar*) gHades->getRuntimeDb()->getContainer("TofWalkPar");
	if(!fTofWalkPar) {
 	    Error ("init()", "Retrieve ZERO pointer for HTofWalkPar!");
	    return kFALSE;
	}
	fSizesCells = (HMdcSizesCells*)HMdcSizesCells::getObject();
	if(!fSizesCells)
	{
 	    Error ("init()", "Retrieve ZERO pointer for HMdcSizesCells!");
            return kFALSE;
	}
    } else {
	Error ("init()", "Retrieve ZERO pointer for fMetaMatchIter!");
	return kFALSE;
    }
    return kTRUE;
}
Bool_t HParticleCandFiller::reinit(void)
{
    Bool_t result = fSizesCells->initContainer();
    cropLay = 0;
    if(result){
	cropLay = new HMdcLayer();
	cropLay->fillRKPlane    (fSizesCells);
	cropLay->fillLayerPoints(fSizesCells,kFALSE);
	gHades->addObjectToOutput(cropLay);
	
    }
    return result;
}
Int_t HParticleCandFiller::execute (void)
{
    
    
    
    
    clearVector();
    HMetaMatch2* metaMatch = 0;
    fMetaMatchIter->Reset();
    while( (metaMatch = (HMetaMatch2*)fMetaMatchIter->Next()) != 0){
	fillCand(metaMatch);
    }
    fillSingleProperties();
    fillCollectiveProperties();
    fillOutput();
    clearVector();
    fillGeantAcceptance();
    return 0;
}
Bool_t HParticleCandFiller::finalize (void)
{
    
    return kTRUE;
}
Int_t  HParticleCandFiller::findBestRich(HMetaMatch2* meta,HMdcSeg* mdcSeg1 )
{
    
    
    
    if(!meta || !mdcSeg1 || !fCatRichHit || !fFillerPar) return -1;
    Int_t nrich = meta->getNCandForRich();
    Int_t slot  = -1;
    if(nrich >= 1){
        Float_t zVertex        = gHades->getCurrentEvent()->getHeader()->getVertexZ();
        HVertex& vertexClust   = gHades->getCurrentEvent()->getHeader()->getVertexCluster();
	Float_t chi2Clust = vertexClust.getChi2(); 
	Double_t qabest   = 100000;
        Double_t qa       = 0;
        HRichHit* richHit  = 0;
	for(Int_t i = 0; i < nrich; i ++){
	    richHit      = HCategoryManager::getObject(richHit   ,fCatRichHit   ,meta->getARichInd(i));
	    Float_t richCorr    = 0.;
            Float_t theta       = richHit->getTheta();
            Float_t phi         = richHit->getPhi();   
	    if(fbdoRichVertexCorr && chi2Clust > -0.5 && zVertex > -999.9F && zVertex < 100 )
	    {
		if(fRichCorrectionVersion == 0) {
		    richCorr = fFillerPar->getRichCorr(zVertex,theta,phi); 
                    theta += richCorr; 
		} else  {
		    Float_t x = (Float_t)richHit->getRich700CircleCenterX();
		    Float_t y = (Float_t)richHit->getRich700CircleCenterY();
		    fRich700DigiPar->getInterpolatedSectorThetaPhi(x,y,zVertex,  theta,phi); 
		}
	    }
	    Float_t  thetaCor = theta;
            Float_t  phiCor   = phi;
	    if(!fbIsSimulation && fbdoRichAlign){ 
		if(fRichCorrectionVersion == 0) {
		    Double_t thcor = theta;
                    Double_t phcor = phi;
		    HParticleAngleCor::alignRichRing(theta,phi,thcor,phcor); 
		    thetaCor = (Float_t)thcor;
                    phiCor   = (Float_t)phcor;
		} else {
		    fRich700DigiPar->getAlignedThetaPhi(theta,phi,thetaCor,phiCor); 
		}
	    }
            qa = HParticleTool::calcRichQA(mdcSeg1, thetaCor, phiCor); 
	    if(qa >= 0 && qa < qabest) {
		slot   = i;
		qabest = qa;
	    }
	} 
    }
    return slot;
}
void HParticleCandFiller::fillCandNoMeta(Bool_t rkSuccess,HMetaMatch2* meta,candidate& cand,Int_t num)
{
    HMdcTrkCand*  mdcTrkCand = 0;
    HMdcSeg*      mdcSeg1    = 0;
    HMdcSeg*      mdcSeg2    = 0;
    HTofCluster*  tofClst    = 0;
    HTofHit*      tofHit1    = 0;
    HTofHit*      tofHit2    = 0;
    HShowerHit*   showerHit  = 0;
    HEmcCluster*  emcClst    = 0;
    HRpcCluster*  rpcClst    = 0;
    HRichHit*     richHit    = 0;
    HSplineTrack* spline     = 0;
    HRKTrackB*    rk         = 0;
    HKalTrack*    kal        = 0;
    cand.rkSuccess = rkSuccess;
    cand.nCand     = num;
    
    
    cand           .fillMeta(meta);
    cand.mdctrk    .fillMeta(meta);
    cand.spline    .fillMeta(meta);
    mdcTrkCand = HCategoryManager::getObject(mdcTrkCand,fCatMdcTrkCand,cand.mdctrk.ind);
    if(mdcTrkCand){
	cand.seg1.ind  = mdcTrkCand->getSeg1Ind();
	cand.seg2.ind  = mdcTrkCand->getSeg2Ind();
    }
    mdcSeg1    = HCategoryManager::getObject(mdcSeg1   ,fCatMdcSeg    ,cand.seg1.ind);
    mdcSeg2    = HCategoryManager::getObject(mdcSeg2   ,fCatMdcSeg    ,cand.seg2.ind);
    spline     = HCategoryManager::getObject(spline    ,fCatSpline    ,cand.spline.ind);
    Int_t slot = findBestRich(meta,mdcSeg1);
    cand.richhit.fillMeta(meta,slot);    
    richHit    = HCategoryManager::getObject(richHit   ,fCatRichHit   ,cand.richhit.ind);
    cand.mdctrk    .fill(mdcTrkCand);
    cand.seg1      .fill(mdcSeg1);
    cand.seg2      .fill(mdcSeg2);
    cand.spline    .fill(spline);
    cand.richhit   .fill(richHit);    
    
    cand.system = -1;
    cand.selectTof    = kNoUse;
    cand.usedMeta     = kNoUse;
    cand.rk.selectTof = kNoUse;
    cand.rk.usedMeta  = kNoUse;
    cand.kal.selectTof = kNoUse;
    cand.kal.usedMeta  = kNoUse;
    cand.objects.reset();
    cand.objects.pMdcTrk    =  mdcTrkCand;
    cand.objects.pSeg1      = (HMdcSegSim* )    mdcSeg1;
    cand.objects.pSeg2      = (HMdcSegSim* )    mdcSeg2;
    cand.objects.pRichHit   = (HRichHitSim*)    richHit;
    cand.objects.pTofClst   = (HTofClusterSim*) tofClst;
    cand.objects.pTofHit1   = (HTofHitSim*)     tofHit1;
    cand.objects.pTofHit2   = (HTofHitSim*)     tofHit2;
    cand.objects.pRpcClst   = (HRpcClusterSim*) rpcClst;
    cand.objects.pShowerHit = (HShowerHitSim*)  showerHit;
    cand.objects.pEmcClst   = (HEmcClusterSim*) emcClst;
    cand.objects.pSpline    = spline;
    if     (fmomSwitch == Particle::kMomRK    ) {
	cand.rk.ind = meta->getRungeKuttaInd();
	rk = HCategoryManager::getObject(rk,fCatRK,cand.rk.ind);
	cand.rk.fill(rk);
	cand.objects.pRk = rk;
    }
    else if(fmomSwitch == Particle::kMomKalman) {
	cand.kal.ind = meta->getKalmanFilterInd();
	kal = HCategoryManager::getObject(kal,fCatKalman,cand.kal.ind);
	cand.kal.fill(kal);
	cand.objects.pKalman = kal;
    }
    else   {Error("fillCandNoMeta","Unknown momentum option");}
}
void  HParticleCandFiller::fillCandTof(Bool_t rkSuccess,HMetaMatch2* meta,candidate& cand,Int_t num, Int_t n)
{
    HMdcTrkCand*  mdcTrkCand = 0;
    HMdcSeg*      mdcSeg1    = 0;
    HMdcSeg*      mdcSeg2    = 0;
    HTofCluster*  tofClst    = 0;
    HTofHit*      tofHit1    = 0;
    HTofHit*      tofHit2    = 0;
    HShowerHit*   showerHit  = 0;
    HEmcCluster*  emcClst    = 0;
    HRpcCluster*  rpcClst    = 0;
    HRichHit*     richHit    = 0;
    HSplineTrack* spline     = 0;
    HRKTrackB*    rk         = 0;
    HKalTrack*    kal        = 0;
    cand.rkSuccess = rkSuccess;
    cand.nCand     = num;
    
    
    cand           .fillMeta(meta);
    cand.mdctrk    .fillMeta(meta);
    cand.spline    .fillMeta(meta);
    mdcTrkCand = HCategoryManager::getObject(mdcTrkCand,fCatMdcTrkCand,cand.mdctrk.ind);
    if(mdcTrkCand){
	cand.seg1.ind  = mdcTrkCand->getSeg1Ind();
	cand.seg2.ind  = mdcTrkCand->getSeg2Ind();
    }
    mdcSeg1    = HCategoryManager::getObject(mdcSeg1   ,fCatMdcSeg    ,cand.seg1.ind);
    mdcSeg2    = HCategoryManager::getObject(mdcSeg2   ,fCatMdcSeg    ,cand.seg2.ind);
    spline     = HCategoryManager::getObject(spline    ,fCatSpline    ,cand.spline.ind);
    Int_t slot = findBestRich(meta,mdcSeg1);
    cand.richhit.fillMeta(meta,slot);    
    richHit    = HCategoryManager::getObject(richHit   ,fCatRichHit   ,cand.richhit.ind);
    cand.mdctrk    .fill(mdcTrkCand);
    cand.seg1      .fill(mdcSeg1);
    cand.seg2      .fill(mdcSeg2);
    cand.spline    .fill(spline);
    cand.richhit   .fill(richHit);    
    
    
    
    cand.tof[kTofClst].fillMeta(meta,n,kTofClst);
    cand.tof[kTofHit1].fillMeta(meta,n,kTofHit1);
    cand.tof[kTofHit2].fillMeta(meta,n,kTofHit2);
    tofClst   = HCategoryManager::getObject(tofClst   ,fCatTofCluster,cand.tof[kTofClst].ind);
    tofHit1   = HCategoryManager::getObject(tofHit1   ,fCatTofHit    ,cand.tof[kTofHit1].ind);
    tofHit2   = HCategoryManager::getObject(tofHit2   ,fCatTofHit    ,cand.tof[kTofHit2].ind);
    cand.tof[kTofClst].fill(tofClst,kTofClst);
    cand.tof[kTofHit1].fill(tofHit1,kTofHit1);
    cand.tof[kTofHit2].fill(tofHit2,kTofHit2);
    
    
    
    cand.system =  1;
    cand.selectTof    = kNoUse;
    cand.usedMeta     = kNoUse;
    cand.rk.usedMeta  = kNoUse;
    cand.rk.selectTof = kNoUse;
    cand.kal.usedMeta  = kNoUse;
    cand.kal.selectTof = kNoUse;
    
    
    
    
    
    
    
    Float_t    tofqua[3] = {-1,-1,-1};
    Float_t minTof   = 1e37;
    Int_t minIndMeta = kNoUse;
    if(tofClst||tofHit1||tofHit2) {
	
	for(Int_t i = 0; i < 3; i ++) {
	    tofqua[i] = cand.tof[i].quality;
	    if (tofqua[i] >= 0 && tofqua[i] < minTof){ minTof = tofqua[i]; minIndMeta = i;}
	}
	cand.selectTof = minIndMeta; 
    }
    cand.usedMeta = minIndMeta;
    
    cand.objects.reset();
    cand.objects.pMdcTrk    =  mdcTrkCand;
    cand.objects.pSeg1      = (HMdcSegSim* )    mdcSeg1;
    cand.objects.pSeg2      = (HMdcSegSim* )    mdcSeg2;
    cand.objects.pRichHit   = (HRichHitSim*)    richHit;
    cand.objects.pTofClst   = (HTofClusterSim*) tofClst;
    cand.objects.pTofHit1   = (HTofHitSim*)     tofHit1;
    cand.objects.pTofHit2   = (HTofHitSim*)     tofHit2;
    cand.objects.pRpcClst   = (HRpcClusterSim*) rpcClst;
    cand.objects.pShowerHit = (HShowerHitSim*)  showerHit;
    cand.objects.pEmcClst   = (HEmcClusterSim*) emcClst;
    cand.objects.pSpline    = spline;
    
    if(!rkSuccess){ 
	cand.rk.ind       = meta->getRungeKuttaInd();
	cand.rk.usedMeta  = kNoUse;
        cand.kal.ind      = meta->getKalmanFilterInd();
	cand.kal.usedMeta = kNoUse;
    } else {
	if(fmomSwitch == Particle::kMomRK    )
	{
	    
	    
	    HRKTrackB* rkhit [3] = { 0, 0, 0};
	    Float_t    rkqua [3] = {-1,-1,-1};
	    Float_t min  = 1e37;
	    Int_t minInd = kNoUse;
	    
	    
	    
	    
	    rkhit[kTofClst] = tofClst ? HCategoryManager::getObject(rkhit[kTofClst],fCatRK,meta->getRungeKuttaIndTofClst(n)): 0 ;
	    rkhit[kTofHit1] = tofHit1 ? HCategoryManager::getObject(rkhit[kTofHit1],fCatRK,meta->getRungeKuttaIndTofHit1(n)): 0 ;
	    rkhit[kTofHit2] = tofHit2 ? HCategoryManager::getObject(rkhit[kTofHit2],fCatRK,meta->getRungeKuttaIndTofHit2(n)): 0 ;
	    if(fsortSwitch==0){
		if(tofClst && rkhit[kTofClst]) { rkqua[kTofClst] = rkhit[kTofClst]->getQualityTof();}
		if(tofHit1 && rkhit[kTofHit1]) { rkqua[kTofHit1] = rkhit[kTofHit1]->getQualityTof();}
		if(tofHit2 && rkhit[kTofHit2]) { rkqua[kTofHit2] = rkhit[kTofHit2]->getQualityTof();}
	    } else if (fsortSwitch==1) {
		if(tofClst && rkhit[kTofClst]) { rkqua[kTofClst] = rkhit[kTofClst]->getMetaRadius();}
		if(tofHit1 && rkhit[kTofHit1]) { rkqua[kTofHit1] = rkhit[kTofHit1]->getMetaRadius();}
		if(tofHit2 && rkhit[kTofHit2]) { rkqua[kTofHit2] = rkhit[kTofHit2]->getMetaRadius();}
	    } else   { Error("fillCandTof","Unknown sort option");}
	    
	    for(Int_t i = 0; i < 3; i ++) {
		if (rkqua[i] > 1e37 ) rkqua[i] = 1e36;
		if (rkqua[i] >= 0 && rkqua[i] < min){ min = rkqua[i]; minInd = i;}
	    }
	    cand.rk.selectTof = minInd; 
	    
	    
	    if     (minInd == kTofClst) { cand.rk.ind = meta->getRungeKuttaIndTofClst(n); cand.rk.usedMeta = kTofClst;}
	    else if(minInd == kTofHit1) { cand.rk.ind = meta->getRungeKuttaIndTofHit1(n); cand.rk.usedMeta = kTofHit1;}
	    else if(minInd == kTofHit2) { cand.rk.ind = meta->getRungeKuttaIndTofHit2(n); cand.rk.usedMeta = kTofHit2;}
	    else {
		Error("fillCand()","(system == 1) No valid rk object!");
	    }
	    
	}
	else if(fmomSwitch == Particle::kMomKalman)
	{
	    
	    
	    HKalTrack* rkhit [3] = { 0, 0, 0};
	    Float_t    rkqua [3] = {-1,-1,-1};
	    Float_t min  = 1e37;
	    Int_t minInd = kNoUse;
	    
	    
	    
	    
	    rkhit[kTofClst] = tofClst ? HCategoryManager::getObject(rkhit[kTofClst],fCatKalman,meta->getKalmanFilterIndTofClst(n)): 0 ;
	    rkhit[kTofHit1] = tofHit1 ? HCategoryManager::getObject(rkhit[kTofHit1],fCatKalman,meta->getKalmanFilterIndTofHit1(n)): 0 ;
	    rkhit[kTofHit2] = tofHit2 ? HCategoryManager::getObject(rkhit[kTofHit2],fCatKalman,meta->getKalmanFilterIndTofHit2(n)): 0 ;
	    if(fsortSwitch==0){
		if(tofClst && rkhit[kTofClst]) { rkqua[kTofClst] = rkhit[kTofClst]->getQualityTof();}
		if(tofHit1 && rkhit[kTofHit1]) { rkqua[kTofHit1] = rkhit[kTofHit1]->getQualityTof();}
		if(tofHit2 && rkhit[kTofHit2]) { rkqua[kTofHit2] = rkhit[kTofHit2]->getQualityTof();}
	    } else if(fsortSwitch==1){
		if(tofClst && rkhit[kTofClst]) { rkqua[kTofClst] = rkhit[kTofClst]->getMetaRadius();}
		if(tofHit1 && rkhit[kTofHit1]) { rkqua[kTofHit1] = rkhit[kTofHit1]->getMetaRadius();}
		if(tofHit2 && rkhit[kTofHit2]) { rkqua[kTofHit2] = rkhit[kTofHit2]->getMetaRadius();}
	    } else   { Error("fillCandTof","Unknown sort option");}
	    
	    for(Int_t i = 0; i < 3; i ++) {
		if (rkqua[i] > 1e37 ) rkqua[i] = 1e36;
		if (rkqua[i] >= 0 && rkqua[i] < min){ min = rkqua[i]; minInd = i;}
	    }
	    cand.kal.selectTof = minInd; 
	    
	    
	    if     (minInd == kTofClst) { cand.kal.ind = meta->getKalmanFilterIndTofClst(n); cand.kal.usedMeta = kTofClst;}
	    else if(minInd == kTofHit1) { cand.kal.ind = meta->getKalmanFilterIndTofHit1(n); cand.kal.usedMeta = kTofHit1;}
	    else if(minInd == kTofHit2) { cand.kal.ind = meta->getKalmanFilterIndTofHit2(n); cand.kal.usedMeta = kTofHit2;}
	    else {
		Error("fillCand()","(system == 1) No valid kaltrack object!");
	    }
	    
	} else   {
	    Error("fillCandTof","Unknown momentum option");
	}
    }
    if(fmomSwitch == Particle::kMomRK)
    {
	
	rk = HCategoryManager::getObject(rk,fCatRK,cand.rk.ind);
	cand.rk.fill(rk);
	cand.objects.pRk    = rk;
	
    }    else if(fmomSwitch == Particle::kMomKalman) {
	
	kal = HCategoryManager::getObject(kal,fCatKalman,cand.kal.ind);
	cand.kal.fill(kal);
	cand.objects.pKalman = kal;
	
    } else   {
	Error("fillCandTof","Unknown momentum option");
    }
    
}
void  HParticleCandFiller::fillCandRpc(Bool_t rkSuccess,HMetaMatch2* meta,candidate& cand,Int_t num, Int_t n)
{
    HMdcTrkCand*  mdcTrkCand = 0;
    HMdcSeg*      mdcSeg1    = 0;
    HMdcSeg*      mdcSeg2    = 0;
    HTofCluster*  tofClst    = 0;
    HTofHit*      tofHit1    = 0;
    HTofHit*      tofHit2    = 0;
    HShowerHit*   showerHit  = 0;
    HEmcCluster*  emcClst    = 0;
    HRpcCluster*  rpcClst    = 0;
    HRichHit*     richHit    = 0;
    HSplineTrack* spline     = 0;
    HRKTrackB*    rk         = 0;
    HKalTrack*    kal        = 0;
    cand.rkSuccess = rkSuccess;
    cand.nCand     = num;
    
    
    cand           .fillMeta(meta);
    cand.mdctrk    .fillMeta(meta);
    cand.spline    .fillMeta(meta);
    mdcTrkCand = HCategoryManager::getObject(mdcTrkCand,fCatMdcTrkCand,cand.mdctrk.ind);
    if(mdcTrkCand){
	cand.seg1.ind  = mdcTrkCand->getSeg1Ind();
	cand.seg2.ind  = mdcTrkCand->getSeg2Ind();
    }
    mdcSeg1    = HCategoryManager::getObject(mdcSeg1   ,fCatMdcSeg    ,cand.seg1.ind);
    mdcSeg2    = HCategoryManager::getObject(mdcSeg2   ,fCatMdcSeg    ,cand.seg2.ind);
    spline     = HCategoryManager::getObject(spline    ,fCatSpline    ,cand.spline.ind);
    Int_t slot = findBestRich(meta,mdcSeg1);
    cand.richhit.fillMeta(meta,slot);    
    richHit    = HCategoryManager::getObject(richHit   ,fCatRichHit   ,cand.richhit.ind);
    cand.mdctrk    .fill(mdcTrkCand);
    cand.seg1      .fill(mdcSeg1);
    cand.seg2      .fill(mdcSeg2);
    cand.spline    .fill(spline);
    cand.richhit   .fill(richHit);    
    
    
    
    cand.rpcclst      .fillMeta(meta,n);   
    rpcClst   = HCategoryManager::getObject(rpcClst   ,fCatRpcCluster,cand.rpcclst.ind);
    cand.rpcclst      .fill(rpcClst);
    
    
    
    cand.system = 0;
    cand.selectTof    = kRpcClst;
    cand.usedMeta     = kRpcClst;
    if(rkSuccess) {
	cand.rk.usedMeta  = kRpcClst;
	cand.rk.selectTof = kRpcClst;
	cand.rk.ind       = meta->getRungeKuttaIndRpcClst(n);
	cand.kal.usedMeta  = kRpcClst;
	cand.kal.selectTof = kRpcClst;
	cand.kal.ind       = meta->getKalmanFilterIndRpcClst(n);
    } else {
	cand.rk.usedMeta  = kNoUse;
	cand.rk.selectTof = kNoUse;
	cand.rk.ind       = meta->getRungeKuttaInd();
	cand.kal.usedMeta  = kNoUse;
	cand.kal.selectTof = kNoUse;
	cand.kal.ind       = meta->getKalmanFilterInd();
    }
    
    cand.objects.reset();
    cand.objects.pMdcTrk    =  mdcTrkCand;
    cand.objects.pSeg1      = (HMdcSegSim* )    mdcSeg1;
    cand.objects.pSeg2      = (HMdcSegSim* )    mdcSeg2;
    cand.objects.pRichHit   = (HRichHitSim*)    richHit;
    cand.objects.pTofClst   = (HTofClusterSim*) tofClst;
    cand.objects.pTofHit1   = (HTofHitSim*)     tofHit1;
    cand.objects.pTofHit2   = (HTofHitSim*)     tofHit2;
    cand.objects.pRpcClst   = (HRpcClusterSim*) rpcClst;
    cand.objects.pShowerHit = (HShowerHitSim*)  showerHit;
    cand.objects.pEmcClst   = (HEmcClusterSim*) emcClst;
    cand.objects.pSpline    = spline;
    if(fmomSwitch == Particle::kMomRK)
    {
	
	rk = HCategoryManager::getObject(rk,fCatRK,cand.rk.ind);
	cand.rk.fill(rk);
	cand.objects.pRk    = rk;
	
    }    else if(fmomSwitch == Particle::kMomKalman) {
	
	kal = HCategoryManager::getObject(kal,fCatKalman,cand.kal.ind);
	cand.kal.fill(kal);
	cand.objects.pKalman = kal;
	
    } else   {
	Error("fillCandRpc","Unknown momentum option");
    }
}
void  HParticleCandFiller::fillCandShower(Bool_t rkSuccess,HMetaMatch2* meta,candidate& cand,Int_t num, Int_t n)
{
    HMdcTrkCand*  mdcTrkCand = 0;
    HMdcSeg*      mdcSeg1    = 0;
    HMdcSeg*      mdcSeg2    = 0;
    HTofCluster*  tofClst    = 0;
    HTofHit*      tofHit1    = 0;
    HTofHit*      tofHit2    = 0;
    HShowerHit*   showerHit  = 0;
    HEmcCluster*  emcClst    = 0;
    HRpcCluster*  rpcClst    = 0;
    HRichHit*     richHit    = 0;
    HSplineTrack* spline     = 0;
    HRKTrackB*    rk         = 0;
    HKalTrack*    kal        = 0;
    cand.rkSuccess = rkSuccess;
    cand.nCand     = num;
    
    
    cand           .fillMeta(meta);
    cand.mdctrk    .fillMeta(meta);
    cand.spline    .fillMeta(meta);
    mdcTrkCand = HCategoryManager::getObject(mdcTrkCand,fCatMdcTrkCand,cand.mdctrk.ind);
    if(mdcTrkCand){
	cand.seg1.ind  = mdcTrkCand->getSeg1Ind();
	cand.seg2.ind  = mdcTrkCand->getSeg2Ind();
    }
    mdcSeg1    = HCategoryManager::getObject(mdcSeg1   ,fCatMdcSeg    ,cand.seg1.ind);
    mdcSeg2    = HCategoryManager::getObject(mdcSeg2   ,fCatMdcSeg    ,cand.seg2.ind);
    spline     = HCategoryManager::getObject(spline    ,fCatSpline    ,cand.spline.ind);
    Int_t slot = findBestRich(meta,mdcSeg1);
    cand.richhit.fillMeta(meta,slot);    
    richHit    = HCategoryManager::getObject(richHit   ,fCatRichHit   ,cand.richhit.ind);
    cand.mdctrk    .fill(mdcTrkCand);
    cand.seg1      .fill(mdcSeg1);
    cand.seg2      .fill(mdcSeg2);
    cand.spline    .fill(spline);
    cand.richhit   .fill(richHit);    
    
    
    
    cand.showerhit    .fillMeta(meta,n);
    showerHit = HCategoryManager::getObject(showerHit ,fCatShowerHit ,cand.showerhit.ind);
    cand.showerhit    .fill(showerHit);
    
    
    
    cand.system = 0;
    cand.selectTof    = kShowerHit;
    cand.usedMeta     = kShowerHit;
    if(rkSuccess) {
	cand.rk.usedMeta  = kShowerHit;
	cand.rk.selectTof = kShowerHit;
	cand.rk.ind       = meta->getRungeKuttaIndShowerHit(n);
	cand.kal.usedMeta  = kShowerHit;
	cand.kal.selectTof = kShowerHit;
	cand.kal.ind       = meta->getKalmanFilterIndShowerHit(n);
    } else {
	cand.rk.usedMeta  = kNoUse;
	cand.rk.selectTof = kNoUse;
	cand.rk.ind       = meta->getRungeKuttaInd();
	cand.kal.usedMeta  = kNoUse;
	cand.kal.selectTof = kNoUse;
	cand.kal.ind       = meta->getKalmanFilterInd();
    }
    
    cand.objects.reset();
    cand.objects.pMdcTrk    =  mdcTrkCand;
    cand.objects.pSeg1      = (HMdcSegSim* )    mdcSeg1;
    cand.objects.pSeg2      = (HMdcSegSim* )    mdcSeg2;
    cand.objects.pRichHit   = (HRichHitSim*)    richHit;
    cand.objects.pTofClst   = (HTofClusterSim*) tofClst;
    cand.objects.pTofHit1   = (HTofHitSim*)     tofHit1;
    cand.objects.pTofHit2   = (HTofHitSim*)     tofHit2;
    cand.objects.pRpcClst   = (HRpcClusterSim*) rpcClst;
    cand.objects.pShowerHit = (HShowerHitSim*)  showerHit;
    cand.objects.pEmcClst   = (HEmcClusterSim*) emcClst;
    cand.objects.pSpline    = spline;
    if(fmomSwitch == Particle::kMomRK)
    {
	
	rk = HCategoryManager::getObject(rk,fCatRK,cand.rk.ind);
	cand.rk.fill(rk);
	cand.objects.pRk    = rk;
	
    }    else if(fmomSwitch == Particle::kMomKalman) {
	
	kal = HCategoryManager::getObject(kal,fCatKalman,cand.kal.ind);
	cand.kal.fill(kal);
	cand.objects.pKalman = kal;
	
    } else   {
	Error("fillCandRpc","Unknown momentum option");
    }
}
void  HParticleCandFiller::fillCandEmc(Bool_t rkSuccess,HMetaMatch2* meta,candidate& cand,Int_t num, Int_t n)
{
    HMdcTrkCand*  mdcTrkCand = 0;
    HMdcSeg*      mdcSeg1    = 0;
    HMdcSeg*      mdcSeg2    = 0;
    HTofCluster*  tofClst    = 0;
    HTofHit*      tofHit1    = 0;
    HTofHit*      tofHit2    = 0;
    HShowerHit*   showerHit  = 0;
    HEmcCluster*  emcClst    = 0;
    HRpcCluster*  rpcClst    = 0;
    HRichHit*     richHit    = 0;
    HSplineTrack* spline     = 0;
    HRKTrackB*    rk         = 0;
    HKalTrack*    kal        = 0;
    cand.rkSuccess = rkSuccess;
    cand.nCand     = num;
    
    
    cand           .fillMeta(meta);
    cand.mdctrk    .fillMeta(meta);
    cand.spline    .fillMeta(meta);
    mdcTrkCand = HCategoryManager::getObject(mdcTrkCand,fCatMdcTrkCand,cand.mdctrk.ind);
    if(mdcTrkCand){
	cand.seg1.ind  = mdcTrkCand->getSeg1Ind();
	cand.seg2.ind  = mdcTrkCand->getSeg2Ind();
    }
    mdcSeg1    = HCategoryManager::getObject(mdcSeg1   ,fCatMdcSeg    ,cand.seg1.ind);
    mdcSeg2    = HCategoryManager::getObject(mdcSeg2   ,fCatMdcSeg    ,cand.seg2.ind);
    spline     = HCategoryManager::getObject(spline    ,fCatSpline    ,cand.spline.ind);
    Int_t slot = findBestRich(meta,mdcSeg1);
    cand.richhit.fillMeta(meta,slot);    
    richHit    = HCategoryManager::getObject(richHit   ,fCatRichHit   ,cand.richhit.ind);
    cand.mdctrk    .fill(mdcTrkCand);
    cand.seg1      .fill(mdcSeg1);
    cand.seg2      .fill(mdcSeg2);
    cand.spline    .fill(spline);
    cand.richhit   .fill(richHit);    
    
    
    
    cand.emcclst    .fillMeta(meta,n);
    emcClst = HCategoryManager::getObject(emcClst ,fCatEmcCluster ,cand.emcclst.ind);
    cand.emcclst    .fill(emcClst);
    
    
    
    cand.system = 0;
    cand.selectTof    = kEmcClst;
    cand.usedMeta     = kEmcClst;
    if(rkSuccess) {
	cand.rk.usedMeta  = kEmcClst;
	cand.rk.selectTof = kEmcClst;
	cand.rk.ind       = meta->getRungeKuttaIndEmcCluster(n);
	cand.kal.usedMeta  = kShowerHit;
	cand.kal.selectTof = kShowerHit;
	cand.kal.ind       = meta->getKalmanFilterIndEmcCluster(n);
    } else {
	cand.rk.usedMeta  = kNoUse;
	cand.rk.selectTof = kNoUse;
	cand.rk.ind       = meta->getRungeKuttaInd();
	cand.kal.usedMeta  = kNoUse;
	cand.kal.selectTof = kNoUse;
	cand.kal.ind       = meta->getKalmanFilterInd();
    }
    
    cand.objects.reset();
    cand.objects.pMdcTrk    =  mdcTrkCand;
    cand.objects.pSeg1      = (HMdcSegSim* )    mdcSeg1;
    cand.objects.pSeg2      = (HMdcSegSim* )    mdcSeg2;
    cand.objects.pRichHit   = (HRichHitSim*)    richHit;
    cand.objects.pTofClst   = (HTofClusterSim*) tofClst;
    cand.objects.pTofHit1   = (HTofHitSim*)     tofHit1;
    cand.objects.pTofHit2   = (HTofHitSim*)     tofHit2;
    cand.objects.pRpcClst   = (HRpcClusterSim*) rpcClst;
    cand.objects.pShowerHit = (HShowerHitSim*)  showerHit;
    cand.objects.pEmcClst   = (HEmcClusterSim*) emcClst;
    cand.objects.pSpline    = spline;
    if(fmomSwitch == Particle::kMomRK)
    {
	
	rk = HCategoryManager::getObject(rk,fCatRK,cand.rk.ind);
	cand.rk.fill(rk);
	cand.objects.pRk    = rk;
	
    }    else if(fmomSwitch == Particle::kMomKalman) {
	
	kal = HCategoryManager::getObject(kal,fCatKalman,cand.kal.ind);
	cand.kal.fill(kal);
	cand.objects.pKalman = kal;
	
    } else   {
	Error("fillCandRpc","Unknown momentum option");
    }
}
void HParticleCandFiller::fillCand(HMetaMatch2* meta){
    
    
    
    
    if(fbnoFake) if(meta->isFake()) { return; }
    Bool_t isEMC = meta->isEmcCluster();
    
    
    
    
    
    HRKTrackB*    rk         = 0;
    HKalTrack*    kal        = 0;
    Bool_t  rkSuccess    = kTRUE;
    Float_t rkchi2       = -1;
    if(fmomSwitch == Particle::kMomRK)
    {
	Int_t   rkReplaceInd = meta->getRungeKuttaInd();
	if ( (rk = HCategoryManager::getObject(rk,fCatRK,rkReplaceInd)) != 0 ) {
	    rkchi2 = rk->getChiq();
            if (rkchi2 < 0)  { rkSuccess = kFALSE; }
	}
    } else if(fmomSwitch == Particle::kMomKalman) {
	Int_t   rkReplaceInd = meta->getKalmanFilterInd();
	if ( (kal = HCategoryManager::getObject(kal,fCatKalman,rkReplaceInd)) != 0 ) {
	    rkchi2 = kal->getChi2();
            if (rkchi2 < 0)  { rkSuccess = kFALSE; }
	}
    } else   {
	Error("fillCand","Unknown momentum option");
    }
    
    Int_t       bestShrInd  = -1;
    HShowerHit* bestShr     = 0;
    HRKTrackB*  bestRKShr   = 0;
    HKalTrack*  bestKalShr  = 0;
    Int_t       bestEmcInd  = -1;
    HEmcCluster* bestEmc    = 0;
    HRKTrackB*  bestRKEmc   = 0;
    HKalTrack*  bestKalEmc  = 0;
 
    if(isEMC){
	
	
	Float_t     emcQA       = 10000.;
	HEmcCluster* hit =0;
	Float_t qa      =-1;
	for(UInt_t i=0; i < meta->getNEmcClusters();i++){
	    Int_t ind = meta->getEmcClusterInd (i);
	    hit = HCategoryManager::getObject(hit,fCatEmcCluster,ind);
	    if(hit){
		if(!rkSuccess){
		    qa=meta->getEmcClusterQuality(i);
		    if(qa>=0 && qa < emcQA) {
			emcQA      = qa;
			bestEmc    = hit;
			bestEmcInd = i;
		    }
		} else {
		    HRKTrackB* rkEmc=0;
		    HKalTrack* kalEmc=0;
		    if     (fmomSwitch == Particle::kMomRK)     {
			rkEmc  = HCategoryManager::getObject(rkEmc ,fCatRK    ,meta->getRungeKuttaIndEmcCluster(i));
			if     (fsortSwitch==0)qa=rkEmc->getQualityEmc();
			else if(fsortSwitch==1)qa=rkEmc->getMetaRadius();
			else   { Error("fillCand","Unknown sort option");}
		    }
		    else if(fmomSwitch == Particle::kMomKalman) {
			kalEmc = HCategoryManager::getObject(kalEmc,fCatKalman,meta->getKalmanFilterIndEmcCluster(i));
			if     (fsortSwitch==0)qa=kalEmc->getQualityEmc();
			else if(fsortSwitch==1)qa=kalEmc->getMetaRadius();
			else   { Error("fillCand","Unknown sort option");}
		    }
		    else   { Error("fillCand","Unknown momentum option");}
		    if(qa>=0 && qa < emcQA) {
			emcQA      = qa;
			bestEmc    = hit;
			bestEmcInd = i;
			bestRKEmc  = rkEmc;
			bestKalEmc = kalEmc;
		    }
		}
	    }
	}
	
    } else {
	
	
	Float_t     shrQA       = 10000.;
	HShowerHit* hit =0;
	Float_t qa      =-1;
	for(UInt_t i=0; i < meta->getNShrHits();i++){
	    Int_t ind = meta->getShowerHitInd (i);
	    hit = HCategoryManager::getObject(hit,fCatShowerHit,ind);
	    if(hit){
		if(!rkSuccess){
		    qa=meta->getShowerHitQuality(i);
		    if(qa>=0 && qa < shrQA) {
			shrQA      = qa;
			bestShr    = hit;
			bestShrInd = i;
		    }
		} else {
		    HRKTrackB* rkShr=0;
		    HKalTrack* kalShr=0;
		    if     (fmomSwitch == Particle::kMomRK)     {
			rkShr  = HCategoryManager::getObject(rkShr ,fCatRK    ,meta->getRungeKuttaIndShowerHit(i));
			if     (fsortSwitch==0)qa=rkShr->getQualityShower();
			else if(fsortSwitch==1)qa=rkShr->getMetaRadius();
			else   { Error("fillCand","Unknown sort option");}
		    }
		    else if(fmomSwitch == Particle::kMomKalman) {
			kalShr = HCategoryManager::getObject(kalShr,fCatKalman,meta->getKalmanFilterIndShowerHit(i));
			if     (fsortSwitch==0)qa=kalShr->getQualityShower();
			else if(fsortSwitch==1)qa=kalShr->getMetaRadius();
			else   { Error("fillCand","Unknown sort option");}
		    }
		    else   { Error("fillCand","Unknown momentum option");}
		    if(qa>=0 && qa < shrQA) {
			shrQA      = qa;
			bestShr    = hit;
			bestShrInd = i;
			bestRKShr  = rkShr;
			bestKalShr = kalShr;
		    }
		}
	    }
	}
	
    }
    HShowerHit*   showerHit  = 0;
    HEmcCluster*  emcClst    = 0;
    
    if(meta->getSystem()==-1)    
    {
	Int_t num = 1;
	candidate& cand = *(new candidate);
	cand.reset();
	cand.isEMC = isEMC;
	fillCandNoMeta(rkSuccess,meta,cand,num);
	
        
	if(fbIsSimulation){
	    fillCandSim(cand);
	}
	
        all_candidates.push_back(&cand);
    }
    
    else if(meta->getSystem()==1)  
    {
	Int_t numTof = meta->getNTofHits();
	
	
	for(Int_t n = 0; n < numTof; n ++)
	{
	    candidate& cand = *(new candidate);
	    cand.reset();
	    cand.isEMC = isEMC;
	    fillCandTof(rkSuccess,meta,cand,numTof,n);
	    
	    
	    if(fbIsSimulation){
		fillCandSim(cand);
	    }
	    
	    all_candidates.push_back(&cand);
	} 
    }
    
    else if(meta->getSystem()==0)  
    {
	Int_t numRpc = meta->getNRpcClusters();
	
	
	for(Int_t n = 0; n < numRpc; n ++)
	{
	    candidate& cand = *(new candidate);
	    cand.reset();
	    cand.isEMC = isEMC;
	    fillCandRpc(rkSuccess,meta,cand,numRpc,n);
	    if(isEMC){
		if(bestEmc){
		    if     (fmomSwitch == Particle::kMomRK)     cand.emcclst     .fillMeta(meta,bestEmcInd,bestRKEmc);
		    else if(fmomSwitch == Particle::kMomKalman) cand.emcclst     .fillMetaKal(meta,bestEmcInd,bestKalEmc);
		    else   { Error("fillCand","Unknown momentum option");}
		    emcClst = HCategoryManager::getObject(emcClst ,fCatEmcCluster ,cand.emcclst.ind);
		    cand.emcclst    .fill(emcClst);
		    cand.objects.pEmcClst = (HEmcClusterSim*)  emcClst;
		}
	    } else {
		if(bestShr){
		    if     (fmomSwitch == Particle::kMomRK)     cand.showerhit    .fillMeta(meta,bestShrInd,bestRKShr);
		    else if(fmomSwitch == Particle::kMomKalman) cand.showerhit    .fillMetaKal(meta,bestShrInd,bestKalShr);
		    else   { Error("fillCand","Unknown momentum option");}
		    showerHit = HCategoryManager::getObject(showerHit ,fCatShowerHit ,cand.showerhit.ind);
		    cand.showerhit    .fill(showerHit);
		    cand.objects.pShowerHit = (HShowerHitSim*)  showerHit;
		}
	    }
	    
	    
	    if(fbIsSimulation){
		fillCandSim(cand);
	    }
	    
	    all_candidates.push_back(&cand);
	} 
	
	
	if(numRpc <= 0 && (bestShr||bestEmc)){ 
	    
	    
	    candidate& cand = *(new candidate);
	    cand.reset();
	    cand.isEMC = isEMC;
	    if(isEMC){
		Int_t numEmc=meta->getNEmcClusters();
		fillCandEmc(rkSuccess,meta,cand,numEmc,bestEmcInd );  
		if     (fmomSwitch == Particle::kMomRK)     cand.emcclst    .fillMeta(meta,bestEmcInd,bestRKEmc);
		else if(fmomSwitch == Particle::kMomKalman) cand.emcclst    .fillMetaKal(meta,bestEmcInd,bestKalEmc);
		else   { Error("fillCand","Unknown momentum option");}
	    } else {
		Int_t numShr=meta->getNShrHits();
		fillCandShower(rkSuccess,meta,cand,numShr,bestShrInd );  
		if     (fmomSwitch == Particle::kMomRK)     cand.showerhit    .fillMeta(meta,bestShrInd,bestRKShr);
		else if(fmomSwitch == Particle::kMomKalman) cand.showerhit    .fillMetaKal(meta,bestShrInd,bestKalShr);
		else   { Error("fillCand","Unknown momentum option");}
	    }
	    
	    
	    if(fbIsSimulation){
		fillCandSim(cand);
	    }
	    
	    all_candidates.push_back(&cand);
	}
	
    }
    
    else if(meta->getSystem()==2)  
    {
	Int_t numTof = meta->getNTofHits();
        Int_t numRpc = meta->getNRpcClusters();
	
	
	for(Int_t n = 0; n < numTof; n ++)
	{
	    candidate& cand = *(new candidate);
	    cand.reset();
	    cand.isEMC = isEMC;
	    fillCandTof(rkSuccess,meta,cand,numTof+numRpc,n);
	    cand.system = 2;
	    if(isEMC){
		if(bestEmc){
		    if     (fmomSwitch == Particle::kMomRK)     cand.emcclst   .fillMeta(meta,bestEmcInd,bestRKEmc);
		    else if(fmomSwitch == Particle::kMomKalman) cand.emcclst   .fillMetaKal(meta,bestEmcInd,bestKalEmc);
		    else   { Error("fillCand","Unknown momentum option");}
		    emcClst = HCategoryManager::getObject(emcClst ,fCatEmcCluster ,cand.emcclst.ind);
		    cand.emcclst    .fill(emcClst);
		    cand.objects.pEmcClst = (HEmcClusterSim*)  emcClst;
		}
	    } else {
		if(bestShr){
		    if     (fmomSwitch == Particle::kMomRK)     cand.showerhit    .fillMeta(meta,bestShrInd,bestRKShr);
		    else if(fmomSwitch == Particle::kMomKalman) cand.showerhit    .fillMetaKal(meta,bestShrInd,bestKalShr);
		    else   { Error("fillCand","Unknown momentum option");}
		    showerHit = HCategoryManager::getObject(showerHit ,fCatShowerHit ,cand.showerhit.ind);
		    cand.showerhit    .fill(showerHit);
		    cand.objects.pShowerHit = (HShowerHitSim*)  showerHit;
		}
	    }
	    
	    
	    if(fbIsSimulation){
		fillCandSim(cand);
	    }
	    
	    all_candidates.push_back(&cand);
	} 
	
	
	for(Int_t n = 0; n < numRpc; n ++)
	{
	    candidate& cand = *(new candidate);
	    cand.reset();
	    cand.isEMC = isEMC;
	    fillCandRpc(rkSuccess,meta,cand,numTof+numRpc,n);
	    cand.system=2;
	    if(isEMC){
		if(bestEmc){
		    if     (fmomSwitch == Particle::kMomRK)     cand.emcclst    .fillMeta(meta,bestEmcInd,bestRKEmc);
		    else if(fmomSwitch == Particle::kMomKalman) cand.emcclst    .fillMetaKal(meta,bestEmcInd,bestKalEmc);
		    else   { Error("fillCand","Unknown momentum option");}
		    emcClst = HCategoryManager::getObject(emcClst ,fCatEmcCluster ,cand.emcclst.ind);
		    cand.emcclst    .fill(emcClst);
		    cand.objects.pEmcClst = (HEmcClusterSim*) emcClst;
		}
	    } else {
		if(bestShr){
		    if     (fmomSwitch == Particle::kMomRK)     cand.showerhit    .fillMeta(meta,bestShrInd,bestRKShr);
		    else if(fmomSwitch == Particle::kMomKalman) cand.showerhit    .fillMetaKal(meta,bestShrInd,bestKalShr);
		    else   { Error("fillCand","Unknown momentum option");}
		    showerHit = HCategoryManager::getObject(showerHit ,fCatShowerHit ,cand.showerhit.ind);
		    cand.showerhit    .fill(showerHit);
		    cand.objects.pShowerHit = (HShowerHitSim*)  showerHit;
		}
	    }
	    
	    
	    if(fbIsSimulation){
		fillCandSim(cand);
	    }
	    
	    all_candidates.push_back(&cand);
	} 
	
    } else {
	Error("fillCand()","metamatch system other than -1 to 2 = %i!",meta->getSystem());
        return;
    }
}
void HParticleCandFiller::fillCandSim(candidate& cand)
{
    
    
    
    
    
    
    
    
    
    
    
    Int_t tr = -1;
    pointers& objects = cand.objects;
    
    
    
    
    vector<Int_t> goodTracks;  
    Int_t nGoodInnerTracks = objects.pSeg1->getNNotFakeTracks();
    Int_t nGoodOuterTracks = 0;
    
    
    Int_t nwInnerCut = objects.pSeg1->getSumWires()/2 + 1;
    Int_t ghostFlag=0;
    if( nGoodInnerTracks == 0) {
	ghostFlag = kIsGhost|kIsInnerGhost ;
    } else {
	if(objects.pSeg1->getTrack(0) == gHades->getEmbeddingRealTrackId() )                          ghostFlag |= kIsGhost|kIsInnerGhost ;
	if(nwInnerCut > objects.pSeg1->getNTimes(0) || objects.pSeg1->getNTimes(0)<fMinWireGoodTrack) ghostFlag |= kIsGhost|kIsInnerGhost;
    }
    Int_t nwOuterCut = 0;
    if(objects.pSeg2){
	nGoodOuterTracks = objects.pSeg2->getNNotFakeTracks();
	if( nGoodOuterTracks == 0) {
	    ghostFlag |= kIsGhost|kIsOuterGhost ;
	} else {
	    nwOuterCut = objects.pSeg2->getSumWires()/2 + 1;
	    if(nwOuterCut > objects.pSeg2->getNTimes(0) || objects.pSeg2->getNTimes(0)<fMinWireGoodTrack) ghostFlag |= kIsGhost|kIsOuterGhost;
	    if(objects.pSeg1->getTrack(0) != objects.pSeg2->getTrack(0))                                  ghostFlag |= kIsGhost;
	    if(objects.pSeg2->getTrack(0) == gHades->getEmbeddingRealTrackId() )                          ghostFlag |= kIsGhost|kIsOuterGhost ;
	}
    }
    
    
    for(Int_t i = 0; i < nGoodInnerTracks; i ++) {
 	tr =  objects.pSeg1->getTrack(i);
	if(objects.pSeg2) tr = objects.pSeg1->getGoodTrack(i,objects.pSeg2,fMinWireGoodTrack); 
	if(tr > 0 ) { 
	    goodTracks.push_back(tr);
	}
    }
    
    Int_t ghost = 0;
    if(objects.pSeg1) { 
	for(Int_t i = 0; i < objects.pSeg1->getNTracks(); i++) { 
	    ghost = kIsInnerGhost;
	    if(i < nGoodInnerTracks && objects.pSeg1->getNTimes(i) >= fMinWireGoodTrack) ghost = 0; 
	    tr = objects.pSeg1->getTrack(i);
	    if(tr < 1 && tr != gHades->getEmbeddingRealTrackId() ) continue; 
	    if(tr > 0 && find(goodTracks.begin(),goodTracks.end(),tr) == goodTracks.end()) { 
		cand.mdc1Tracks.addTrack(tr,objects.pSeg1->getNTimes(i),kIsInInnerMDC|kIsGhost|ghost,fScaleGhostTrack);
	    } else {                                                                         
                if(tr > 0 ) cand.mdc1Tracks.addTrack(tr,objects.pSeg1->getNTimes(i),kIsInInnerMDC,fScaleGoodTrack); 
                else        cand.mdc1Tracks.addTrack(tr,objects.pSeg1->getNTimes(i),kIsInInnerMDC);                 
	    }
	}
	cand.mdc1Tracks.calcWeights();
    }
    if(objects.pSeg2) { 
	for(Int_t i = 0; i < objects.pSeg2->getNTracks(); i++) { 
	    ghost = kIsOuterGhost;
	    if(i < nGoodOuterTracks && objects.pSeg2->getNTimes(i) >= fMinWireGoodTrack) ghost = 0; 
	    tr = objects.pSeg2->getTrack(i);
	    if(tr < 1 && tr != gHades->getEmbeddingRealTrackId() ) continue; 
            if(tr > 0 && find(goodTracks.begin(),goodTracks.end(),tr) == goodTracks.end()) { 
		cand.mdc2Tracks.addTrack(tr,objects.pSeg2->getNTimes(i),kIsInOuterMDC|kIsGhost|ghost,fScaleGhostTrack);
	    } else {                                                                         
                if(tr > 0 ) cand.mdc2Tracks.addTrack(tr,objects.pSeg2->getNTimes(i),kIsInOuterMDC,fScaleGoodTrack); 
                else        cand.mdc2Tracks.addTrack(tr,objects.pSeg2->getNTimes(i),kIsInOuterMDC);                
	    }
	}
	cand.mdc2Tracks.calcWeights();
    }
    
    
    
    
    if(objects.pRichHit){
	if(objects.pRichHit->weigTrack1 > 0.0 && (objects.pRichHit->track1 > 0 || objects.pRichHit->track1 == gHades->getEmbeddingRealTrackId())) {
	    cand.richTracks.addTrack(objects.pRichHit->track1,objects.pRichHit->weigTrack1,kIsInRICH);
	}
	if(objects.pRichHit->weigTrack2 > 0.0 && (objects.pRichHit->track2 > 0 || objects.pRichHit->track2 == gHades->getEmbeddingRealTrackId())) {
	    cand.richTracks.addTrack(objects.pRichHit->track2,objects.pRichHit->weigTrack2,kIsInRICH);
	}
	if(objects.pRichHit->weigTrack3 > 0.0 && (objects.pRichHit->track3 > 0 || objects.pRichHit->track3 == gHades->getEmbeddingRealTrackId())) {
	    cand.richTracks.addTrack(objects.pRichHit->track3,objects.pRichHit->weigTrack3,kIsInRICH);
	}
	cand.richTracks.calcWeights();
    }
    
    
    
    
    
    
    
    Int_t selMetaHit = kNoUse;
    if(cand.rkSuccess) {
	if     (fmomSwitch == Particle::kMomRK)     selMetaHit = cand.rk.usedMeta;
	else if(fmomSwitch == Particle::kMomKalman) selMetaHit = cand.kal.usedMeta;
        else   Error("fillCandSim()","Unknown momswitch !");
    } else { selMetaHit = cand.usedMeta;    }  
    
    
    if       (selMetaHit  == kTofClst)  {
	if(objects.pTofClst) {
	    Int_t nTr = objects.pTofClst ->getNParticipants();
	    for(Int_t i = 0; i < nTr && i < 3; i ++)
	    {
		
		tr = objects.pTofClst->getNTrack1(i);
		tr = HParticleTool::findFirstHitInTof(tr,2);
                Bool_t good = kFALSE;
		if( tr > 0 || tr == gHades->getEmbeddingRealTrackId() ) {
		    cand.tofTracks.addTrack(tr,1,kIsInTOF|kIsInMETA);
		    good = kTRUE;
		}
		if(objects.pTofClst->getNTrack1(i) != objects.pTofClst->getNTrack2(i) ){  
		    tr = objects.pTofClst ->getNTrack2(i);
		    tr = HParticleTool::findFirstHitInTof(tr,2);
		    if( tr > 0 || tr == gHades->getEmbeddingRealTrackId() ) {
			cand.tofTracks.addTrack(tr,1,kIsInTOF|kIsInMETA);
		    }
		} else { if (good) cand.tofTracks.addTrack(tr,1,kIsInTOF|kIsInMETA); } 
		cand.tofTracks.calcWeights();
	    }
	}
	
    } else if(selMetaHit  == kTofHit1)  {
	if(objects.pTofHit1){
	    tr = objects.pTofHit1 ->getNTrack1();
	    tr = HParticleTool::findFirstHitInTof(tr,2);
	    Bool_t good = kFALSE;
	    if( tr > 0 || tr == gHades->getEmbeddingRealTrackId() ) {
		cand.tofTracks.addTrack(tr,1,kIsInTOF|kIsInMETA);
		good = kTRUE;
	    }
	    if(objects.pTofHit1 ->getNTrack1() != objects.pTofHit1 ->getNTrack2() ){  
		tr = objects.pTofHit1 ->getNTrack2();
		tr = HParticleTool::findFirstHitInTof(tr,2);
		if( tr > 0 || tr == gHades->getEmbeddingRealTrackId() ) {
		    cand.tofTracks.addTrack(tr,1,kIsInTOF|kIsInMETA);
		}
	    } else { if (good) cand.tofTracks.addTrack(tr,1,kIsInTOF|kIsInMETA); } 
	    cand.tofTracks.calcWeights();
	}
	
    } else if(selMetaHit  == kTofHit2)  {
	if(objects.pTofHit2){
	    tr = objects.pTofHit2 ->getNTrack1();
	    tr = HParticleTool::findFirstHitInTof(tr,2);
            Bool_t good = kFALSE;
	    if( tr > 0 || tr == gHades->getEmbeddingRealTrackId() ) {
		cand.tofTracks.addTrack(tr,1,kIsInTOF|kIsInMETA);
		good = kTRUE;
	    }
	    if(objects.pTofHit2 ->getNTrack1() != objects.pTofHit2 ->getNTrack2() ){  
		tr = objects.pTofHit2 ->getNTrack2();
		tr = HParticleTool::findFirstHitInTof(tr,2);
		if( tr > 0 || tr == gHades->getEmbeddingRealTrackId() ) {
		    cand.tofTracks.addTrack(tr,1,kIsInTOF|kIsInMETA);
		}
	    } else { if (good) cand.tofTracks.addTrack(tr,1,kIsInTOF|kIsInMETA); } 
	    cand.tofTracks.calcWeights();
	}
	
    } else if(selMetaHit  == kRpcClst)  {
	if(objects.pRpcClst){
	    Int_t tracks[4];
	    Bool_t atbox[4];
	    objects.pRpcClst->getTrackList  (tracks);
	    objects.pRpcClst->getIsAtBoxList(atbox);
	    Int_t nTr = objects.pRpcClst->howManyTracks();
	    if(nTr > 4) {
		Warning("fillCandSim()","RPC number of tracks > 4! Should not happen.");
                nTr = 4;
	    }
	    for(Int_t i = 0; i < nTr; i ++) {
		tr = tracks[i];
		if( tr >= 0 || tr == gHades->getEmbeddingRealTrackId()){
		    if(tr >= 0 ) { 
			if(atbox[i]) { cand.rpcTracks.addTrack(tr,1,kIsInRPC|kIsInMETA);}
			else         {
			    Warning("fillCandSim()","RPC track not at box! Should not happen.");
			}
		    } else {     
			cand.rpcTracks.addTrack(tr,1,kIsInRPC|kIsInMETA);
		    }
		}
	    }
	    cand.rpcTracks.calcWeights();
	}
	
    } else if(selMetaHit  == kShowerHit){
	
    } else if(selMetaHit  == kEmcClst){
	
    } else { ; }  
    
    
    
    if(objects.pEmcClst){
	Int_t nTr = objects.pEmcClst->getNTracks();
	for(Int_t i = 0; i < nTr; i ++) {
	    tr = objects.pEmcClst->getTrack(i);
	    if( tr >= 0 || tr == gHades->getEmbeddingRealTrackId()){
		cand.emcTracks.addTrack(tr,1,kIsInEMC|kIsInMETA);
	    }
	}
	cand.emcTracks.calcWeights();
    }
    if(objects.pShowerHit){
	Int_t nTr = objects.pShowerHit->getNTracks();
	for(Int_t i = 0; i < nTr; i ++) {
	    tr = objects.pShowerHit->getTrack(i);
	    if( tr >= 0 || tr == gHades->getEmbeddingRealTrackId()){
		cand.showerTracks.addTrack(tr,1,kIsInSHOWER|kIsInMETA);
	    }
	}
	cand.showerTracks.calcWeights();
    }
    
    
    
    
    cand.commonTracks.addTrackWeight(cand.mdc1Tracks);
    cand.commonTracks.addTrackWeight(cand.mdc2Tracks);
    cand.commonTracks.addTrackWeight(cand.richTracks);
    cand.commonTracks.addTrackWeight(cand.tofTracks);
    cand.commonTracks.addTrackWeight(cand.rpcTracks);
    if(cand.isEMC)cand.commonTracks.addTrackWeight(cand.emcTracks);
    else          cand.commonTracks.addTrackWeight(cand.showerTracks);
    cand.commonTracks.sortWeightNdet();
    if(cand.commonTracks.tracks.size()>0)cand.commonTracks.tracks[0]->flag |= ghostFlag; 
    
}
void HParticleCandFiller::fillSingleProperties()
{
    
    
    
    
    Float_t zVertex        = gHades->getCurrentEvent()->getHeader()->getVertexZ();
    HVertex& vertexClust   = gHades->getCurrentEvent()->getHeader()->getVertexCluster();
    Float_t chi2Clust = vertexClust.getChi2(); 
    for(UInt_t i = 0; i < all_candidates.size(); i ++){
	candidate& cand = *(all_candidates[i]);
	
	
	
	
	if(cand.richhit.ind >= 0){
	    Float_t richCorr    = 0;
	    if(fbdoRichVertexCorr && chi2Clust > -0.5 && zVertex > -999.9F && zVertex < 100 )
	    {
		if(fRichCorrectionVersion == 0) {
		    richCorr = fFillerPar->getRichCorr(zVertex,cand.richhit.theta,cand.richhit.phi); 
		} else  {
                    Float_t th = cand.richhit.theta;
		    fRich700DigiPar->getInterpolatedSectorThetaPhi(cand.richhit.centerx,cand.richhit.centery,zVertex,  th,cand.richhit.phi); 
		    richCorr = th - cand.richhit.theta;
		}
	    }
	    cand.corrThetaRich  = richCorr;
	    cand.richhit.theta += richCorr;
	    if(!fbIsSimulation && fbdoRichAlign){ 
		Float_t thetaCor = cand.richhit.theta;
		Float_t phiCor   = cand.richhit.phi;
		if(fRichCorrectionVersion == 0) {
		    Double_t thcor = cand.richhit.theta;
                    Double_t phcor = cand.richhit.phi;
		    HParticleAngleCor::alignRichRing(cand.richhit.theta,cand.richhit.phi,thcor,phcor); 
		    thetaCor = (Float_t)thcor;
                    phiCor   = (Float_t)phcor;
		} else {
		    fRich700DigiPar->getAlignedThetaPhi(cand.richhit.theta,cand.richhit.phi,thetaCor,phiCor); 
		}
		cand.alignThetaRich  = thetaCor - cand.richhit.theta;
                cand.alignPhiRich    = phiCor - cand.richhit.phi;
		cand.richhit.theta = thetaCor;
                cand.richhit.phi   = phiCor;
            }
	    Float_t dTheta = cand.deltaThetaSeg1(); 
	    Float_t dPhi   = cand.deltaPhiSeg1();   
	    Float_t mom = cand.spline.p;            
	    if     (fmomSwitch == Particle::kMomRK )    { if( cand.rk.chi2  > 0) mom = cand.rk.p; }      
	    else if(fmomSwitch == Particle::kMomKalman) { if( cand.kal.chi2 > 0) mom = cand.kal.p;}      
	    else   Error("fillSingleProperties()","Unknown momswitch !");
	    cand.richMdcQuality     = sqrt( dTheta * dTheta + dPhi * dPhi);
	    cand.richRkQuality      = -1;
	    cand.hasRingCorrelation =  0;
	    if(fFillerPar->acceptPhiTheta(cand.sector,mom,dPhi,dTheta)){
		cand.hasRingCorrelation = kIsRICHMDC;
	    }
	    if( (fmomSwitch == Particle::kMomRK     && cand.rk.chi2  > 0 ) ||
                (fmomSwitch == Particle::kMomKalman && cand.kal.chi2 > 0 )
	      ){
		if(fmomSwitch == Particle::kMomRK){
		    dTheta = cand.deltaThetaRk(); 
		    dPhi   = cand.deltaPhiRk();   
		} else if(fmomSwitch == Particle::kMomKalman){
		    dTheta = cand.deltaThetaKal(); 
		    dPhi   = cand.deltaPhiKal();   
		} else Error("fillSingleProperties()","Unknown momswitch !");
		cand.richRkQuality = sqrt( dTheta * dTheta + dPhi * dPhi);
		if(fFillerPar->acceptPhiTheta(cand.sector,mom,dPhi,dTheta)){
		    cand.hasRingCorrelation = cand.hasRingCorrelation | kIsRICHRK;  
		}
	    }
	}
	
	
	
	
	if(fbgoodSeg0   && cand.seg1.chi2 < 0)          { cand.used = 0; }
	if(fbgoodSeg1   && cand.seg2.chi2 < 0)          { cand.used = 0; }
	if(fbgoodMeta   && cand.system    < 0)          { cand.used = 0; }
	if(fbgoodLepton && cand.hasRingCorrelation == 0){ cand.used = 0; }
	if(fbgoodRK ){
	    if     (fmomSwitch == Particle::kMomRK     ) { if( cand.rk.chi2   < 0) { cand.used = 0; }}
	    else if(fmomSwitch == Particle::kMomKalman ) { if( cand.kal.chi2  < 0) { cand.used = 0; }}
	    else Error("fillSingleProperties()","Unknown momswitch !");
	}
	
    }
}
void HParticleCandFiller::fillCollectiveProperties()
{
    
    
    
    
    UInt_t  size = all_candidates.size();
    if(size < 2) return;
    Float_t phi1=0.,phi2=0.,theta1=0.,theta2=0.;
    Float_t oAngle = -1;
    Int_t ind      =  0;
    
    for(UInt_t i = 0; i < size ; i ++){                   
	candidate& cand = *(all_candidates[i]);
	if(cand.rkSuccess){
	    if(fmomSwitch == Particle::kMomRK){
		phi1   = cand.rk.phi;
		theta1 = cand.rk.theta;
	    } else if(fmomSwitch == Particle::kMomKalman){
		phi1   = cand.kal.phi;
		theta1 = cand.kal.theta;
	    } else Error("fillCollectiveProperties()","Unknown momswitch !");
	} else {
	    phi1   = cand.seg1.phi;
	    theta1 = cand.seg1.theta;
	}
	
	for(UInt_t j = 0; j < size; j ++){              
	    candidate& cand2 = *(all_candidates[j]);
            if(cand.ind      == cand2.ind)      continue;   
	    if(cand.sector   != cand2.sector)   continue;   
            if(cand.seg1.ind == cand2.seg1.ind) continue;   
	    if(cand2.rkSuccess){
		if(fmomSwitch == Particle::kMomRK){
		    phi2   = cand2.rk.phi;
		    theta2 = cand2.rk.theta;
		} else if(fmomSwitch == Particle::kMomKalman){
		    phi1   = cand.kal.phi;
		    theta1 = cand.kal.theta;
		} else Error("fillCollectiveProperties()","Unknown momswitch !");
	    } else {
		phi2   = cand2.seg1.phi;
		theta2 = cand2.seg1.theta;
	    }
	    oAngle = HParticleTool::getOpeningAngle(phi1,theta1,phi2,theta2);
	    if(oAngle < fAngleCloseTrack){                   
		cand.closeTracks.addTrack(j,oAngle);
	    }
	} 
	
        
	UInt_t ntr  = (cand.closeTracks).tracks.size();
	for(UInt_t l = 0; l < ntr; l ++){
	    ind    = (cand.closeTracks).tracks[l].ind;
	    oAngle = (cand.closeTracks).tracks[l].oAngle;
	    candidate& c = *(all_candidates[ind]);
	    if(c.seg1.chi2 > 0){
		if(oAngle < fabs(cand.oAFitted)){
		    cand.oAFitted = oAngle;
		    if( c.hasRingCorrelation == 0) { cand.oAFitted *= -1; }  
		}
	    } else {
		if(oAngle < fabs(cand.oANoFitted)){
		    cand.oANoFitted = oAngle;
		    if( c.hasRingCorrelation == 0) { cand.oANoFitted *= -1; }  
		}
	    }
	}
        
    } 
    
}
void HParticleCandFiller::fillOutput()
{
    
    
    
    
    
    HParticleCand*    part  =  0;
    HParticleCandSim* partS =  0;
    Int_t index             = -1;
    HVertex vertex = (gHades->getCurrentEvent()->getHeader()->getVertexReco());
    Float_t sweights= vertex.getSumOfWeights();
    HGeomVector vvertex = vertex.getPos();
    Int_t nmin = 8;
    for(UInt_t i = 0; i < all_candidates.size(); i ++){
	candidate* cand = all_candidates[i];
	if(!cand->used) continue;
	pointers& objects = cand->objects;
	if(fbIsSimulation) {
	    if((partS = HCategoryManager::newObject(partS,fCatParticleCand,index)) != 0 )
	    {
		cand->fillParticleCand   (partS,index,fmomSwitch);
		cand->fillParticleCandSim(partS,fCatGeantKine);
                if(fbdoMETAQANorm)     HParticleTool::normDX(partS,fTofWalkPar);
                if(fbdoMomentumCorr)   HParticleTool::setCorrectedMomentum(partS);
		if(fbdoPathLengthCorr&&sweights>=nmin) HParticleTool::correctPathLength(partS,vvertex,cropLay->getMdcPlanes(),cropLay->getTargetMidPoint());
		checkCropLayer(partS);
	    }
	} else {
	    if((part = HCategoryManager::newObject(part,fCatParticleCand,index)) != 0 )
	    {
		cand->fillParticleCand   (part,index,fmomSwitch);
                if(fbdoMETAQANorm)     HParticleTool::normDX(part,fTofWalkPar);
		if(fbdoMomentumCorr)   HParticleTool::setCorrectedMomentum(part);
		if(fbdoPathLengthCorr&&sweights>=nmin) HParticleTool::correctPathLength(part,vvertex,cropLay->getMdcPlanes(),cropLay->getTargetMidPoint());
		checkCropLayer(part);
           }
	}
	if(fbIsDebug){
	    candidate* c;
	    c = HCategoryManager::newSlot(c,fCatParticleDebug,index);
	    if(c){
		c =  new (c) candidate(*cand);
	    }
	}
	if(fbFillMdc) {
	    HParticleMdc* partMdc = 0;
	    partMdc = HCategoryManager::newSlot(partMdc,fCatParticleMdc,index);
	    if(partMdc){
		partMdc =  new (partMdc) HParticleMdc();
                partMdc->setIndex(index);
		partMdc->fill(objects.pSeg1);
                partMdc->fill(objects.pSeg2);
                partMdc->fill(objects.pMdcTrk);
	    }
	}
	
    }
}
void  HParticleCandFiller::fillGeantAcceptance()
{
    if(fCatGeantKine&&fbdoGeantAcceptance){
        Int_t n = fCatGeantKine->getEntries();
        HGeantKine* kine = 0;
	Bool_t cropped [4] ;
	for(Int_t i = 0 ; i < n ; i ++){
	    kine = HCategoryManager::getObject(kine,fCatGeantKine,i);
	    if(kine){
		if(!kine->isAcceptanceFilled()) kine->fillAcceptanceBit();
		if(cropLay){
		    for(Int_t j = 0 ; j < 4 ; j ++) cropped[j] = kFALSE;
		    HParticleTool::checkCropedLayer(kine,cropLay,cropped);
		    for(Int_t j = 0 ; j < 4 ; j ++) {
			kine->unsetAtMdcEdge(j);
			if(cropped[j]) kine->setAtMdcEdge(j);
		    }
                    kine->setCropedFilled();
		}
	    }
	}
    }
}
void  HParticleCandFiller::checkCropLayer(HParticleCand* c)
{
    if(!cropLay) return;
    Int_t sec = c->getSector();
    if(!fSizesCells->modStatus(sec,0) && !fSizesCells->modStatus(sec,1)) return;
    Int_t layers[2] = {0,5};
    HMdcSeg* seg[2] = {0,0} ;
    Double_t x,y,z;
    Double_t xmin,xmax,ymin,ymax;
    for(Int_t s = 0; s < 2; s++){ 
	if(s==0) seg[0] = HParticleTool::getMdcSeg(c->getInnerSegInd());
	else     seg[1] = HParticleTool::getMdcSeg(c->getOuterSegInd());
        if(!seg[s]) continue;
	for(Int_t m = 0; m < 2; m++){ 
	    Int_t mod = 2*s+m;        
	    if(!fSizesCells->modStatus(sec,mod)) continue;
	    Int_t lay = layers[m];
	    (*fSizesCells)[sec][mod][lay].calcSegIntersec(seg[s]->getZ(),seg[s]->getR(),seg[s]->getTheta(),seg[s]->getPhi(),x,y,z);
	    (*fSizesCells)[sec][mod][lay].transTo(x,y,z);
	    TVector2  p(x,y);
	    if(!cropLay->getCropedLayersEdge(p, sec,mod,lay, xmin,xmax,ymin,ymax)){
		c->setAtMdcEdge(s*2+m);
	    }
	}
    }
    
}
void HParticleCandFiller::clearVector(void)
{
    
    UInt_t nt = all_candidates.size();
    for(UInt_t i = 0; i < nt; i ++){
	candidate* cand = all_candidates[i];
	delete cand;
    }
    all_candidates.clear();
}
 hparticlecandfiller.cc:10  hparticlecandfiller.cc:11  hparticlecandfiller.cc:12  hparticlecandfiller.cc:13  hparticlecandfiller.cc:14  hparticlecandfiller.cc:15  hparticlecandfiller.cc:16  hparticlecandfiller.cc:17  hparticlecandfiller.cc:18  hparticlecandfiller.cc:19  hparticlecandfiller.cc:20  hparticlecandfiller.cc:21  hparticlecandfiller.cc:22  hparticlecandfiller.cc:23  hparticlecandfiller.cc:24  hparticlecandfiller.cc:25  hparticlecandfiller.cc:26  hparticlecandfiller.cc:27  hparticlecandfiller.cc:28  hparticlecandfiller.cc:29  hparticlecandfiller.cc:30  hparticlecandfiller.cc:31  hparticlecandfiller.cc:32  hparticlecandfiller.cc:33  hparticlecandfiller.cc:34  hparticlecandfiller.cc:35  hparticlecandfiller.cc:36  hparticlecandfiller.cc:37  hparticlecandfiller.cc:38  hparticlecandfiller.cc:39  hparticlecandfiller.cc:40  hparticlecandfiller.cc:41  hparticlecandfiller.cc:42  hparticlecandfiller.cc:43  hparticlecandfiller.cc:44  hparticlecandfiller.cc:45  hparticlecandfiller.cc:46  hparticlecandfiller.cc:47  hparticlecandfiller.cc:48  hparticlecandfiller.cc:49  hparticlecandfiller.cc:50  hparticlecandfiller.cc:51  hparticlecandfiller.cc:52  hparticlecandfiller.cc:53  hparticlecandfiller.cc:54  hparticlecandfiller.cc:55  hparticlecandfiller.cc:56  hparticlecandfiller.cc:57  hparticlecandfiller.cc:58  hparticlecandfiller.cc:59  hparticlecandfiller.cc:60  hparticlecandfiller.cc:61  hparticlecandfiller.cc:62  hparticlecandfiller.cc:63  hparticlecandfiller.cc:64  hparticlecandfiller.cc:65  hparticlecandfiller.cc:66  hparticlecandfiller.cc:67  hparticlecandfiller.cc:68  hparticlecandfiller.cc:69  hparticlecandfiller.cc:70  hparticlecandfiller.cc:71  hparticlecandfiller.cc:72  hparticlecandfiller.cc:73  hparticlecandfiller.cc:74  hparticlecandfiller.cc:75  hparticlecandfiller.cc:76  hparticlecandfiller.cc:77  hparticlecandfiller.cc:78  hparticlecandfiller.cc:79  hparticlecandfiller.cc:80  hparticlecandfiller.cc:81  hparticlecandfiller.cc:82  hparticlecandfiller.cc:83  hparticlecandfiller.cc:84  hparticlecandfiller.cc:85  hparticlecandfiller.cc:86  hparticlecandfiller.cc:87  hparticlecandfiller.cc:88  hparticlecandfiller.cc:89  hparticlecandfiller.cc:90  hparticlecandfiller.cc:91  hparticlecandfiller.cc:92  hparticlecandfiller.cc:93  hparticlecandfiller.cc:94  hparticlecandfiller.cc:95  hparticlecandfiller.cc:96  hparticlecandfiller.cc:97  hparticlecandfiller.cc:98  hparticlecandfiller.cc:99  hparticlecandfiller.cc:100  hparticlecandfiller.cc:101  hparticlecandfiller.cc:102  hparticlecandfiller.cc:103  hparticlecandfiller.cc:104  hparticlecandfiller.cc:105  hparticlecandfiller.cc:106  hparticlecandfiller.cc:107  hparticlecandfiller.cc:108  hparticlecandfiller.cc:109  hparticlecandfiller.cc:110  hparticlecandfiller.cc:111  hparticlecandfiller.cc:112  hparticlecandfiller.cc:113  hparticlecandfiller.cc:114  hparticlecandfiller.cc:115  hparticlecandfiller.cc:116  hparticlecandfiller.cc:117  hparticlecandfiller.cc:118  hparticlecandfiller.cc:119  hparticlecandfiller.cc:120  hparticlecandfiller.cc:121  hparticlecandfiller.cc:122  hparticlecandfiller.cc:123  hparticlecandfiller.cc:124  hparticlecandfiller.cc:125  hparticlecandfiller.cc:126  hparticlecandfiller.cc:127  hparticlecandfiller.cc:128  hparticlecandfiller.cc:129  hparticlecandfiller.cc:130  hparticlecandfiller.cc:131  hparticlecandfiller.cc:132  hparticlecandfiller.cc:133  hparticlecandfiller.cc:134  hparticlecandfiller.cc:135  hparticlecandfiller.cc:136  hparticlecandfiller.cc:137  hparticlecandfiller.cc:138  hparticlecandfiller.cc:139  hparticlecandfiller.cc:140  hparticlecandfiller.cc:141  hparticlecandfiller.cc:142  hparticlecandfiller.cc:143  hparticlecandfiller.cc:144  hparticlecandfiller.cc:145  hparticlecandfiller.cc:146  hparticlecandfiller.cc:147  hparticlecandfiller.cc:148  hparticlecandfiller.cc:149  hparticlecandfiller.cc:150  hparticlecandfiller.cc:151  hparticlecandfiller.cc:152  hparticlecandfiller.cc:153  hparticlecandfiller.cc:154  hparticlecandfiller.cc:155  hparticlecandfiller.cc:156  hparticlecandfiller.cc:157  hparticlecandfiller.cc:158  hparticlecandfiller.cc:159  hparticlecandfiller.cc:160  hparticlecandfiller.cc:161  hparticlecandfiller.cc:162  hparticlecandfiller.cc:163  hparticlecandfiller.cc:164  hparticlecandfiller.cc:165  hparticlecandfiller.cc:166  hparticlecandfiller.cc:167  hparticlecandfiller.cc:168  hparticlecandfiller.cc:169  hparticlecandfiller.cc:170  hparticlecandfiller.cc:171  hparticlecandfiller.cc:172  hparticlecandfiller.cc:173  hparticlecandfiller.cc:174  hparticlecandfiller.cc:175  hparticlecandfiller.cc:176  hparticlecandfiller.cc:177  hparticlecandfiller.cc:178  hparticlecandfiller.cc:179  hparticlecandfiller.cc:180  hparticlecandfiller.cc:181  hparticlecandfiller.cc:182  hparticlecandfiller.cc:183  hparticlecandfiller.cc:184  hparticlecandfiller.cc:185  hparticlecandfiller.cc:186  hparticlecandfiller.cc:187  hparticlecandfiller.cc:188  hparticlecandfiller.cc:189  hparticlecandfiller.cc:190  hparticlecandfiller.cc:191  hparticlecandfiller.cc:192  hparticlecandfiller.cc:193  hparticlecandfiller.cc:194  hparticlecandfiller.cc:195  hparticlecandfiller.cc:196  hparticlecandfiller.cc:197  hparticlecandfiller.cc:198  hparticlecandfiller.cc:199  hparticlecandfiller.cc:200  hparticlecandfiller.cc:201  hparticlecandfiller.cc:202  hparticlecandfiller.cc:203  hparticlecandfiller.cc:204  hparticlecandfiller.cc:205  hparticlecandfiller.cc:206  hparticlecandfiller.cc:207  hparticlecandfiller.cc:208  hparticlecandfiller.cc:209  hparticlecandfiller.cc:210  hparticlecandfiller.cc:211  hparticlecandfiller.cc:212  hparticlecandfiller.cc:213  hparticlecandfiller.cc:214  hparticlecandfiller.cc:215  hparticlecandfiller.cc:216  hparticlecandfiller.cc:217  hparticlecandfiller.cc:218  hparticlecandfiller.cc:219  hparticlecandfiller.cc:220  hparticlecandfiller.cc:221  hparticlecandfiller.cc:222  hparticlecandfiller.cc:223  hparticlecandfiller.cc:224  hparticlecandfiller.cc:225  hparticlecandfiller.cc:226  hparticlecandfiller.cc:227  hparticlecandfiller.cc:228  hparticlecandfiller.cc:229  hparticlecandfiller.cc:230  hparticlecandfiller.cc:231  hparticlecandfiller.cc:232  hparticlecandfiller.cc:233  hparticlecandfiller.cc:234  hparticlecandfiller.cc:235  hparticlecandfiller.cc:236  hparticlecandfiller.cc:237  hparticlecandfiller.cc:238  hparticlecandfiller.cc:239  hparticlecandfiller.cc:240  hparticlecandfiller.cc:241  hparticlecandfiller.cc:242  hparticlecandfiller.cc:243  hparticlecandfiller.cc:244  hparticlecandfiller.cc:245  hparticlecandfiller.cc:246  hparticlecandfiller.cc:247  hparticlecandfiller.cc:248  hparticlecandfiller.cc:249  hparticlecandfiller.cc:250  hparticlecandfiller.cc:251  hparticlecandfiller.cc:252  hparticlecandfiller.cc:253  hparticlecandfiller.cc:254  hparticlecandfiller.cc:255  hparticlecandfiller.cc:256  hparticlecandfiller.cc:257  hparticlecandfiller.cc:258  hparticlecandfiller.cc:259  hparticlecandfiller.cc:260  hparticlecandfiller.cc:261  hparticlecandfiller.cc:262  hparticlecandfiller.cc:263  hparticlecandfiller.cc:264  hparticlecandfiller.cc:265  hparticlecandfiller.cc:266  hparticlecandfiller.cc:267  hparticlecandfiller.cc:268  hparticlecandfiller.cc:269  hparticlecandfiller.cc:270  hparticlecandfiller.cc:271  hparticlecandfiller.cc:272  hparticlecandfiller.cc:273  hparticlecandfiller.cc:274  hparticlecandfiller.cc:275  hparticlecandfiller.cc:276  hparticlecandfiller.cc:277  hparticlecandfiller.cc:278  hparticlecandfiller.cc:279  hparticlecandfiller.cc:280  hparticlecandfiller.cc:281  hparticlecandfiller.cc:282  hparticlecandfiller.cc:283  hparticlecandfiller.cc:284  hparticlecandfiller.cc:285  hparticlecandfiller.cc:286  hparticlecandfiller.cc:287  hparticlecandfiller.cc:288  hparticlecandfiller.cc:289  hparticlecandfiller.cc:290  hparticlecandfiller.cc:291  hparticlecandfiller.cc:292  hparticlecandfiller.cc:293  hparticlecandfiller.cc:294  hparticlecandfiller.cc:295  hparticlecandfiller.cc:296  hparticlecandfiller.cc:297  hparticlecandfiller.cc:298  hparticlecandfiller.cc:299  hparticlecandfiller.cc:300  hparticlecandfiller.cc:301  hparticlecandfiller.cc:302  hparticlecandfiller.cc:303  hparticlecandfiller.cc:304  hparticlecandfiller.cc:305  hparticlecandfiller.cc:306  hparticlecandfiller.cc:307  hparticlecandfiller.cc:308  hparticlecandfiller.cc:309  hparticlecandfiller.cc:310  hparticlecandfiller.cc:311  hparticlecandfiller.cc:312  hparticlecandfiller.cc:313  hparticlecandfiller.cc:314  hparticlecandfiller.cc:315  hparticlecandfiller.cc:316  hparticlecandfiller.cc:317  hparticlecandfiller.cc:318  hparticlecandfiller.cc:319  hparticlecandfiller.cc:320  hparticlecandfiller.cc:321  hparticlecandfiller.cc:322  hparticlecandfiller.cc:323  hparticlecandfiller.cc:324  hparticlecandfiller.cc:325  hparticlecandfiller.cc:326  hparticlecandfiller.cc:327  hparticlecandfiller.cc:328  hparticlecandfiller.cc:329  hparticlecandfiller.cc:330  hparticlecandfiller.cc:331  hparticlecandfiller.cc:332  hparticlecandfiller.cc:333  hparticlecandfiller.cc:334  hparticlecandfiller.cc:335  hparticlecandfiller.cc:336  hparticlecandfiller.cc:337  hparticlecandfiller.cc:338  hparticlecandfiller.cc:339  hparticlecandfiller.cc:340  hparticlecandfiller.cc:341  hparticlecandfiller.cc:342  hparticlecandfiller.cc:343  hparticlecandfiller.cc:344  hparticlecandfiller.cc:345  hparticlecandfiller.cc:346  hparticlecandfiller.cc:347  hparticlecandfiller.cc:348  hparticlecandfiller.cc:349  hparticlecandfiller.cc:350  hparticlecandfiller.cc:351  hparticlecandfiller.cc:352  hparticlecandfiller.cc:353  hparticlecandfiller.cc:354  hparticlecandfiller.cc:355  hparticlecandfiller.cc:356  hparticlecandfiller.cc:357  hparticlecandfiller.cc:358  hparticlecandfiller.cc:359  hparticlecandfiller.cc:360  hparticlecandfiller.cc:361  hparticlecandfiller.cc:362  hparticlecandfiller.cc:363  hparticlecandfiller.cc:364  hparticlecandfiller.cc:365  hparticlecandfiller.cc:366  hparticlecandfiller.cc:367  hparticlecandfiller.cc:368  hparticlecandfiller.cc:369  hparticlecandfiller.cc:370  hparticlecandfiller.cc:371  hparticlecandfiller.cc:372  hparticlecandfiller.cc:373  hparticlecandfiller.cc:374  hparticlecandfiller.cc:375  hparticlecandfiller.cc:376  hparticlecandfiller.cc:377  hparticlecandfiller.cc:378  hparticlecandfiller.cc:379  hparticlecandfiller.cc:380  hparticlecandfiller.cc:381  hparticlecandfiller.cc:382  hparticlecandfiller.cc:383  hparticlecandfiller.cc:384  hparticlecandfiller.cc:385  hparticlecandfiller.cc:386  hparticlecandfiller.cc:387  hparticlecandfiller.cc:388  hparticlecandfiller.cc:389  hparticlecandfiller.cc:390  hparticlecandfiller.cc:391  hparticlecandfiller.cc:392  hparticlecandfiller.cc:393  hparticlecandfiller.cc:394  hparticlecandfiller.cc:395  hparticlecandfiller.cc:396  hparticlecandfiller.cc:397  hparticlecandfiller.cc:398  hparticlecandfiller.cc:399  hparticlecandfiller.cc:400  hparticlecandfiller.cc:401  hparticlecandfiller.cc:402  hparticlecandfiller.cc:403  hparticlecandfiller.cc:404  hparticlecandfiller.cc:405  hparticlecandfiller.cc:406  hparticlecandfiller.cc:407  hparticlecandfiller.cc:408  hparticlecandfiller.cc:409  hparticlecandfiller.cc:410  hparticlecandfiller.cc:411  hparticlecandfiller.cc:412  hparticlecandfiller.cc:413  hparticlecandfiller.cc:414  hparticlecandfiller.cc:415  hparticlecandfiller.cc:416  hparticlecandfiller.cc:417  hparticlecandfiller.cc:418  hparticlecandfiller.cc:419  hparticlecandfiller.cc:420  hparticlecandfiller.cc:421  hparticlecandfiller.cc:422  hparticlecandfiller.cc:423  hparticlecandfiller.cc:424  hparticlecandfiller.cc:425  hparticlecandfiller.cc:426  hparticlecandfiller.cc:427  hparticlecandfiller.cc:428  hparticlecandfiller.cc:429  hparticlecandfiller.cc:430  hparticlecandfiller.cc:431  hparticlecandfiller.cc:432  hparticlecandfiller.cc:433  hparticlecandfiller.cc:434  hparticlecandfiller.cc:435  hparticlecandfiller.cc:436  hparticlecandfiller.cc:437  hparticlecandfiller.cc:438  hparticlecandfiller.cc:439  hparticlecandfiller.cc:440  hparticlecandfiller.cc:441  hparticlecandfiller.cc:442  hparticlecandfiller.cc:443  hparticlecandfiller.cc:444  hparticlecandfiller.cc:445  hparticlecandfiller.cc:446  hparticlecandfiller.cc:447  hparticlecandfiller.cc:448  hparticlecandfiller.cc:449  hparticlecandfiller.cc:450  hparticlecandfiller.cc:451  hparticlecandfiller.cc:452  hparticlecandfiller.cc:453  hparticlecandfiller.cc:454  hparticlecandfiller.cc:455  hparticlecandfiller.cc:456  hparticlecandfiller.cc:457  hparticlecandfiller.cc:458  hparticlecandfiller.cc:459  hparticlecandfiller.cc:460  hparticlecandfiller.cc:461  hparticlecandfiller.cc:462  hparticlecandfiller.cc:463  hparticlecandfiller.cc:464  hparticlecandfiller.cc:465  hparticlecandfiller.cc:466  hparticlecandfiller.cc:467  hparticlecandfiller.cc:468  hparticlecandfiller.cc:469  hparticlecandfiller.cc:470  hparticlecandfiller.cc:471  hparticlecandfiller.cc:472  hparticlecandfiller.cc:473  hparticlecandfiller.cc:474  hparticlecandfiller.cc:475  hparticlecandfiller.cc:476  hparticlecandfiller.cc:477  hparticlecandfiller.cc:478  hparticlecandfiller.cc:479  hparticlecandfiller.cc:480  hparticlecandfiller.cc:481  hparticlecandfiller.cc:482  hparticlecandfiller.cc:483  hparticlecandfiller.cc:484  hparticlecandfiller.cc:485  hparticlecandfiller.cc:486  hparticlecandfiller.cc:487  hparticlecandfiller.cc:488  hparticlecandfiller.cc:489  hparticlecandfiller.cc:490  hparticlecandfiller.cc:491  hparticlecandfiller.cc:492  hparticlecandfiller.cc:493  hparticlecandfiller.cc:494  hparticlecandfiller.cc:495  hparticlecandfiller.cc:496  hparticlecandfiller.cc:497  hparticlecandfiller.cc:498  hparticlecandfiller.cc:499  hparticlecandfiller.cc:500  hparticlecandfiller.cc:501  hparticlecandfiller.cc:502  hparticlecandfiller.cc:503  hparticlecandfiller.cc:504  hparticlecandfiller.cc:505  hparticlecandfiller.cc:506  hparticlecandfiller.cc:507  hparticlecandfiller.cc:508  hparticlecandfiller.cc:509  hparticlecandfiller.cc:510  hparticlecandfiller.cc:511  hparticlecandfiller.cc:512  hparticlecandfiller.cc:513  hparticlecandfiller.cc:514  hparticlecandfiller.cc:515  hparticlecandfiller.cc:516  hparticlecandfiller.cc:517  hparticlecandfiller.cc:518  hparticlecandfiller.cc:519  hparticlecandfiller.cc:520  hparticlecandfiller.cc:521  hparticlecandfiller.cc:522  hparticlecandfiller.cc:523  hparticlecandfiller.cc:524  hparticlecandfiller.cc:525  hparticlecandfiller.cc:526  hparticlecandfiller.cc:527  hparticlecandfiller.cc:528  hparticlecandfiller.cc:529  hparticlecandfiller.cc:530  hparticlecandfiller.cc:531  hparticlecandfiller.cc:532  hparticlecandfiller.cc:533  hparticlecandfiller.cc:534  hparticlecandfiller.cc:535  hparticlecandfiller.cc:536  hparticlecandfiller.cc:537  hparticlecandfiller.cc:538  hparticlecandfiller.cc:539  hparticlecandfiller.cc:540  hparticlecandfiller.cc:541  hparticlecandfiller.cc:542  hparticlecandfiller.cc:543  hparticlecandfiller.cc:544  hparticlecandfiller.cc:545  hparticlecandfiller.cc:546  hparticlecandfiller.cc:547  hparticlecandfiller.cc:548  hparticlecandfiller.cc:549  hparticlecandfiller.cc:550  hparticlecandfiller.cc:551  hparticlecandfiller.cc:552  hparticlecandfiller.cc:553  hparticlecandfiller.cc:554  hparticlecandfiller.cc:555  hparticlecandfiller.cc:556  hparticlecandfiller.cc:557  hparticlecandfiller.cc:558  hparticlecandfiller.cc:559  hparticlecandfiller.cc:560  hparticlecandfiller.cc:561  hparticlecandfiller.cc:562  hparticlecandfiller.cc:563  hparticlecandfiller.cc:564  hparticlecandfiller.cc:565  hparticlecandfiller.cc:566  hparticlecandfiller.cc:567  hparticlecandfiller.cc:568  hparticlecandfiller.cc:569  hparticlecandfiller.cc:570  hparticlecandfiller.cc:571  hparticlecandfiller.cc:572  hparticlecandfiller.cc:573  hparticlecandfiller.cc:574  hparticlecandfiller.cc:575  hparticlecandfiller.cc:576  hparticlecandfiller.cc:577  hparticlecandfiller.cc:578  hparticlecandfiller.cc:579  hparticlecandfiller.cc:580  hparticlecandfiller.cc:581  hparticlecandfiller.cc:582  hparticlecandfiller.cc:583  hparticlecandfiller.cc:584  hparticlecandfiller.cc:585  hparticlecandfiller.cc:586  hparticlecandfiller.cc:587  hparticlecandfiller.cc:588  hparticlecandfiller.cc:589  hparticlecandfiller.cc:590  hparticlecandfiller.cc:591  hparticlecandfiller.cc:592  hparticlecandfiller.cc:593  hparticlecandfiller.cc:594  hparticlecandfiller.cc:595  hparticlecandfiller.cc:596  hparticlecandfiller.cc:597  hparticlecandfiller.cc:598  hparticlecandfiller.cc:599  hparticlecandfiller.cc:600  hparticlecandfiller.cc:601  hparticlecandfiller.cc:602  hparticlecandfiller.cc:603  hparticlecandfiller.cc:604  hparticlecandfiller.cc:605  hparticlecandfiller.cc:606  hparticlecandfiller.cc:607  hparticlecandfiller.cc:608  hparticlecandfiller.cc:609  hparticlecandfiller.cc:610  hparticlecandfiller.cc:611  hparticlecandfiller.cc:612  hparticlecandfiller.cc:613  hparticlecandfiller.cc:614  hparticlecandfiller.cc:615  hparticlecandfiller.cc:616  hparticlecandfiller.cc:617  hparticlecandfiller.cc:618  hparticlecandfiller.cc:619  hparticlecandfiller.cc:620  hparticlecandfiller.cc:621  hparticlecandfiller.cc:622  hparticlecandfiller.cc:623  hparticlecandfiller.cc:624  hparticlecandfiller.cc:625  hparticlecandfiller.cc:626  hparticlecandfiller.cc:627  hparticlecandfiller.cc:628  hparticlecandfiller.cc:629  hparticlecandfiller.cc:630  hparticlecandfiller.cc:631  hparticlecandfiller.cc:632  hparticlecandfiller.cc:633  hparticlecandfiller.cc:634  hparticlecandfiller.cc:635  hparticlecandfiller.cc:636  hparticlecandfiller.cc:637  hparticlecandfiller.cc:638  hparticlecandfiller.cc:639  hparticlecandfiller.cc:640  hparticlecandfiller.cc:641  hparticlecandfiller.cc:642  hparticlecandfiller.cc:643  hparticlecandfiller.cc:644  hparticlecandfiller.cc:645  hparticlecandfiller.cc:646  hparticlecandfiller.cc:647  hparticlecandfiller.cc:648  hparticlecandfiller.cc:649  hparticlecandfiller.cc:650  hparticlecandfiller.cc:651  hparticlecandfiller.cc:652  hparticlecandfiller.cc:653  hparticlecandfiller.cc:654  hparticlecandfiller.cc:655  hparticlecandfiller.cc:656  hparticlecandfiller.cc:657  hparticlecandfiller.cc:658  hparticlecandfiller.cc:659  hparticlecandfiller.cc:660  hparticlecandfiller.cc:661  hparticlecandfiller.cc:662  hparticlecandfiller.cc:663  hparticlecandfiller.cc:664  hparticlecandfiller.cc:665  hparticlecandfiller.cc:666  hparticlecandfiller.cc:667  hparticlecandfiller.cc:668  hparticlecandfiller.cc:669  hparticlecandfiller.cc:670  hparticlecandfiller.cc:671  hparticlecandfiller.cc:672  hparticlecandfiller.cc:673  hparticlecandfiller.cc:674  hparticlecandfiller.cc:675  hparticlecandfiller.cc:676  hparticlecandfiller.cc:677  hparticlecandfiller.cc:678  hparticlecandfiller.cc:679  hparticlecandfiller.cc:680  hparticlecandfiller.cc:681  hparticlecandfiller.cc:682  hparticlecandfiller.cc:683  hparticlecandfiller.cc:684  hparticlecandfiller.cc:685  hparticlecandfiller.cc:686  hparticlecandfiller.cc:687  hparticlecandfiller.cc:688  hparticlecandfiller.cc:689  hparticlecandfiller.cc:690  hparticlecandfiller.cc:691  hparticlecandfiller.cc:692  hparticlecandfiller.cc:693  hparticlecandfiller.cc:694  hparticlecandfiller.cc:695  hparticlecandfiller.cc:696  hparticlecandfiller.cc:697  hparticlecandfiller.cc:698  hparticlecandfiller.cc:699  hparticlecandfiller.cc:700  hparticlecandfiller.cc:701  hparticlecandfiller.cc:702  hparticlecandfiller.cc:703  hparticlecandfiller.cc:704  hparticlecandfiller.cc:705  hparticlecandfiller.cc:706  hparticlecandfiller.cc:707  hparticlecandfiller.cc:708  hparticlecandfiller.cc:709  hparticlecandfiller.cc:710  hparticlecandfiller.cc:711  hparticlecandfiller.cc:712  hparticlecandfiller.cc:713  hparticlecandfiller.cc:714  hparticlecandfiller.cc:715  hparticlecandfiller.cc:716  hparticlecandfiller.cc:717  hparticlecandfiller.cc:718  hparticlecandfiller.cc:719  hparticlecandfiller.cc:720  hparticlecandfiller.cc:721  hparticlecandfiller.cc:722  hparticlecandfiller.cc:723  hparticlecandfiller.cc:724  hparticlecandfiller.cc:725  hparticlecandfiller.cc:726  hparticlecandfiller.cc:727  hparticlecandfiller.cc:728  hparticlecandfiller.cc:729  hparticlecandfiller.cc:730  hparticlecandfiller.cc:731  hparticlecandfiller.cc:732  hparticlecandfiller.cc:733  hparticlecandfiller.cc:734  hparticlecandfiller.cc:735  hparticlecandfiller.cc:736  hparticlecandfiller.cc:737  hparticlecandfiller.cc:738  hparticlecandfiller.cc:739  hparticlecandfiller.cc:740  hparticlecandfiller.cc:741  hparticlecandfiller.cc:742  hparticlecandfiller.cc:743  hparticlecandfiller.cc:744  hparticlecandfiller.cc:745  hparticlecandfiller.cc:746  hparticlecandfiller.cc:747  hparticlecandfiller.cc:748  hparticlecandfiller.cc:749  hparticlecandfiller.cc:750  hparticlecandfiller.cc:751  hparticlecandfiller.cc:752  hparticlecandfiller.cc:753  hparticlecandfiller.cc:754  hparticlecandfiller.cc:755  hparticlecandfiller.cc:756  hparticlecandfiller.cc:757  hparticlecandfiller.cc:758  hparticlecandfiller.cc:759  hparticlecandfiller.cc:760  hparticlecandfiller.cc:761  hparticlecandfiller.cc:762  hparticlecandfiller.cc:763  hparticlecandfiller.cc:764  hparticlecandfiller.cc:765  hparticlecandfiller.cc:766  hparticlecandfiller.cc:767  hparticlecandfiller.cc:768  hparticlecandfiller.cc:769  hparticlecandfiller.cc:770  hparticlecandfiller.cc:771  hparticlecandfiller.cc:772  hparticlecandfiller.cc:773  hparticlecandfiller.cc:774  hparticlecandfiller.cc:775  hparticlecandfiller.cc:776  hparticlecandfiller.cc:777  hparticlecandfiller.cc:778  hparticlecandfiller.cc:779  hparticlecandfiller.cc:780  hparticlecandfiller.cc:781  hparticlecandfiller.cc:782  hparticlecandfiller.cc:783  hparticlecandfiller.cc:784  hparticlecandfiller.cc:785  hparticlecandfiller.cc:786  hparticlecandfiller.cc:787  hparticlecandfiller.cc:788  hparticlecandfiller.cc:789  hparticlecandfiller.cc:790  hparticlecandfiller.cc:791  hparticlecandfiller.cc:792  hparticlecandfiller.cc:793  hparticlecandfiller.cc:794  hparticlecandfiller.cc:795  hparticlecandfiller.cc:796  hparticlecandfiller.cc:797  hparticlecandfiller.cc:798  hparticlecandfiller.cc:799  hparticlecandfiller.cc:800  hparticlecandfiller.cc:801  hparticlecandfiller.cc:802  hparticlecandfiller.cc:803  hparticlecandfiller.cc:804  hparticlecandfiller.cc:805  hparticlecandfiller.cc:806  hparticlecandfiller.cc:807  hparticlecandfiller.cc:808  hparticlecandfiller.cc:809  hparticlecandfiller.cc:810  hparticlecandfiller.cc:811  hparticlecandfiller.cc:812  hparticlecandfiller.cc:813  hparticlecandfiller.cc:814  hparticlecandfiller.cc:815  hparticlecandfiller.cc:816  hparticlecandfiller.cc:817  hparticlecandfiller.cc:818  hparticlecandfiller.cc:819  hparticlecandfiller.cc:820  hparticlecandfiller.cc:821  hparticlecandfiller.cc:822  hparticlecandfiller.cc:823  hparticlecandfiller.cc:824  hparticlecandfiller.cc:825  hparticlecandfiller.cc:826  hparticlecandfiller.cc:827  hparticlecandfiller.cc:828  hparticlecandfiller.cc:829  hparticlecandfiller.cc:830  hparticlecandfiller.cc:831  hparticlecandfiller.cc:832  hparticlecandfiller.cc:833  hparticlecandfiller.cc:834  hparticlecandfiller.cc:835  hparticlecandfiller.cc:836  hparticlecandfiller.cc:837  hparticlecandfiller.cc:838  hparticlecandfiller.cc:839  hparticlecandfiller.cc:840  hparticlecandfiller.cc:841  hparticlecandfiller.cc:842  hparticlecandfiller.cc:843  hparticlecandfiller.cc:844  hparticlecandfiller.cc:845  hparticlecandfiller.cc:846  hparticlecandfiller.cc:847  hparticlecandfiller.cc:848  hparticlecandfiller.cc:849  hparticlecandfiller.cc:850  hparticlecandfiller.cc:851  hparticlecandfiller.cc:852  hparticlecandfiller.cc:853  hparticlecandfiller.cc:854  hparticlecandfiller.cc:855  hparticlecandfiller.cc:856  hparticlecandfiller.cc:857  hparticlecandfiller.cc:858  hparticlecandfiller.cc:859  hparticlecandfiller.cc:860  hparticlecandfiller.cc:861  hparticlecandfiller.cc:862  hparticlecandfiller.cc:863  hparticlecandfiller.cc:864  hparticlecandfiller.cc:865  hparticlecandfiller.cc:866  hparticlecandfiller.cc:867  hparticlecandfiller.cc:868  hparticlecandfiller.cc:869  hparticlecandfiller.cc:870  hparticlecandfiller.cc:871  hparticlecandfiller.cc:872  hparticlecandfiller.cc:873  hparticlecandfiller.cc:874  hparticlecandfiller.cc:875  hparticlecandfiller.cc:876  hparticlecandfiller.cc:877  hparticlecandfiller.cc:878  hparticlecandfiller.cc:879  hparticlecandfiller.cc:880  hparticlecandfiller.cc:881  hparticlecandfiller.cc:882  hparticlecandfiller.cc:883  hparticlecandfiller.cc:884  hparticlecandfiller.cc:885  hparticlecandfiller.cc:886  hparticlecandfiller.cc:887  hparticlecandfiller.cc:888  hparticlecandfiller.cc:889  hparticlecandfiller.cc:890  hparticlecandfiller.cc:891  hparticlecandfiller.cc:892  hparticlecandfiller.cc:893  hparticlecandfiller.cc:894  hparticlecandfiller.cc:895  hparticlecandfiller.cc:896  hparticlecandfiller.cc:897  hparticlecandfiller.cc:898  hparticlecandfiller.cc:899  hparticlecandfiller.cc:900  hparticlecandfiller.cc:901  hparticlecandfiller.cc:902  hparticlecandfiller.cc:903  hparticlecandfiller.cc:904  hparticlecandfiller.cc:905  hparticlecandfiller.cc:906  hparticlecandfiller.cc:907  hparticlecandfiller.cc:908  hparticlecandfiller.cc:909  hparticlecandfiller.cc:910  hparticlecandfiller.cc:911  hparticlecandfiller.cc:912  hparticlecandfiller.cc:913  hparticlecandfiller.cc:914  hparticlecandfiller.cc:915  hparticlecandfiller.cc:916  hparticlecandfiller.cc:917  hparticlecandfiller.cc:918  hparticlecandfiller.cc:919  hparticlecandfiller.cc:920  hparticlecandfiller.cc:921  hparticlecandfiller.cc:922  hparticlecandfiller.cc:923  hparticlecandfiller.cc:924  hparticlecandfiller.cc:925  hparticlecandfiller.cc:926  hparticlecandfiller.cc:927  hparticlecandfiller.cc:928  hparticlecandfiller.cc:929  hparticlecandfiller.cc:930  hparticlecandfiller.cc:931  hparticlecandfiller.cc:932  hparticlecandfiller.cc:933  hparticlecandfiller.cc:934  hparticlecandfiller.cc:935  hparticlecandfiller.cc:936  hparticlecandfiller.cc:937  hparticlecandfiller.cc:938  hparticlecandfiller.cc:939  hparticlecandfiller.cc:940  hparticlecandfiller.cc:941  hparticlecandfiller.cc:942  hparticlecandfiller.cc:943  hparticlecandfiller.cc:944  hparticlecandfiller.cc:945  hparticlecandfiller.cc:946  hparticlecandfiller.cc:947  hparticlecandfiller.cc:948  hparticlecandfiller.cc:949  hparticlecandfiller.cc:950  hparticlecandfiller.cc:951  hparticlecandfiller.cc:952  hparticlecandfiller.cc:953  hparticlecandfiller.cc:954  hparticlecandfiller.cc:955  hparticlecandfiller.cc:956  hparticlecandfiller.cc:957  hparticlecandfiller.cc:958  hparticlecandfiller.cc:959  hparticlecandfiller.cc:960  hparticlecandfiller.cc:961  hparticlecandfiller.cc:962  hparticlecandfiller.cc:963  hparticlecandfiller.cc:964  hparticlecandfiller.cc:965  hparticlecandfiller.cc:966  hparticlecandfiller.cc:967  hparticlecandfiller.cc:968  hparticlecandfiller.cc:969  hparticlecandfiller.cc:970  hparticlecandfiller.cc:971  hparticlecandfiller.cc:972  hparticlecandfiller.cc:973  hparticlecandfiller.cc:974  hparticlecandfiller.cc:975  hparticlecandfiller.cc:976  hparticlecandfiller.cc:977  hparticlecandfiller.cc:978  hparticlecandfiller.cc:979  hparticlecandfiller.cc:980  hparticlecandfiller.cc:981  hparticlecandfiller.cc:982  hparticlecandfiller.cc:983  hparticlecandfiller.cc:984  hparticlecandfiller.cc:985  hparticlecandfiller.cc:986  hparticlecandfiller.cc:987  hparticlecandfiller.cc:988  hparticlecandfiller.cc:989  hparticlecandfiller.cc:990  hparticlecandfiller.cc:991  hparticlecandfiller.cc:992  hparticlecandfiller.cc:993  hparticlecandfiller.cc:994  hparticlecandfiller.cc:995  hparticlecandfiller.cc:996  hparticlecandfiller.cc:997  hparticlecandfiller.cc:998  hparticlecandfiller.cc:999  hparticlecandfiller.cc:1000  hparticlecandfiller.cc:1001  hparticlecandfiller.cc:1002  hparticlecandfiller.cc:1003  hparticlecandfiller.cc:1004  hparticlecandfiller.cc:1005  hparticlecandfiller.cc:1006  hparticlecandfiller.cc:1007  hparticlecandfiller.cc:1008  hparticlecandfiller.cc:1009  hparticlecandfiller.cc:1010  hparticlecandfiller.cc:1011  hparticlecandfiller.cc:1012  hparticlecandfiller.cc:1013  hparticlecandfiller.cc:1014  hparticlecandfiller.cc:1015  hparticlecandfiller.cc:1016  hparticlecandfiller.cc:1017  hparticlecandfiller.cc:1018  hparticlecandfiller.cc:1019  hparticlecandfiller.cc:1020  hparticlecandfiller.cc:1021  hparticlecandfiller.cc:1022  hparticlecandfiller.cc:1023  hparticlecandfiller.cc:1024  hparticlecandfiller.cc:1025  hparticlecandfiller.cc:1026  hparticlecandfiller.cc:1027  hparticlecandfiller.cc:1028  hparticlecandfiller.cc:1029  hparticlecandfiller.cc:1030  hparticlecandfiller.cc:1031  hparticlecandfiller.cc:1032  hparticlecandfiller.cc:1033  hparticlecandfiller.cc:1034  hparticlecandfiller.cc:1035  hparticlecandfiller.cc:1036  hparticlecandfiller.cc:1037  hparticlecandfiller.cc:1038  hparticlecandfiller.cc:1039  hparticlecandfiller.cc:1040  hparticlecandfiller.cc:1041  hparticlecandfiller.cc:1042  hparticlecandfiller.cc:1043  hparticlecandfiller.cc:1044  hparticlecandfiller.cc:1045  hparticlecandfiller.cc:1046  hparticlecandfiller.cc:1047  hparticlecandfiller.cc:1048  hparticlecandfiller.cc:1049  hparticlecandfiller.cc:1050  hparticlecandfiller.cc:1051  hparticlecandfiller.cc:1052  hparticlecandfiller.cc:1053  hparticlecandfiller.cc:1054  hparticlecandfiller.cc:1055  hparticlecandfiller.cc:1056  hparticlecandfiller.cc:1057  hparticlecandfiller.cc:1058  hparticlecandfiller.cc:1059  hparticlecandfiller.cc:1060  hparticlecandfiller.cc:1061  hparticlecandfiller.cc:1062  hparticlecandfiller.cc:1063  hparticlecandfiller.cc:1064  hparticlecandfiller.cc:1065  hparticlecandfiller.cc:1066  hparticlecandfiller.cc:1067  hparticlecandfiller.cc:1068  hparticlecandfiller.cc:1069  hparticlecandfiller.cc:1070  hparticlecandfiller.cc:1071  hparticlecandfiller.cc:1072  hparticlecandfiller.cc:1073  hparticlecandfiller.cc:1074  hparticlecandfiller.cc:1075  hparticlecandfiller.cc:1076  hparticlecandfiller.cc:1077  hparticlecandfiller.cc:1078  hparticlecandfiller.cc:1079  hparticlecandfiller.cc:1080  hparticlecandfiller.cc:1081  hparticlecandfiller.cc:1082  hparticlecandfiller.cc:1083  hparticlecandfiller.cc:1084  hparticlecandfiller.cc:1085  hparticlecandfiller.cc:1086  hparticlecandfiller.cc:1087  hparticlecandfiller.cc:1088  hparticlecandfiller.cc:1089  hparticlecandfiller.cc:1090  hparticlecandfiller.cc:1091  hparticlecandfiller.cc:1092  hparticlecandfiller.cc:1093  hparticlecandfiller.cc:1094  hparticlecandfiller.cc:1095  hparticlecandfiller.cc:1096  hparticlecandfiller.cc:1097  hparticlecandfiller.cc:1098  hparticlecandfiller.cc:1099  hparticlecandfiller.cc:1100  hparticlecandfiller.cc:1101  hparticlecandfiller.cc:1102  hparticlecandfiller.cc:1103  hparticlecandfiller.cc:1104  hparticlecandfiller.cc:1105  hparticlecandfiller.cc:1106  hparticlecandfiller.cc:1107  hparticlecandfiller.cc:1108  hparticlecandfiller.cc:1109  hparticlecandfiller.cc:1110  hparticlecandfiller.cc:1111  hparticlecandfiller.cc:1112  hparticlecandfiller.cc:1113  hparticlecandfiller.cc:1114  hparticlecandfiller.cc:1115  hparticlecandfiller.cc:1116  hparticlecandfiller.cc:1117  hparticlecandfiller.cc:1118  hparticlecandfiller.cc:1119  hparticlecandfiller.cc:1120  hparticlecandfiller.cc:1121  hparticlecandfiller.cc:1122  hparticlecandfiller.cc:1123  hparticlecandfiller.cc:1124  hparticlecandfiller.cc:1125  hparticlecandfiller.cc:1126  hparticlecandfiller.cc:1127  hparticlecandfiller.cc:1128  hparticlecandfiller.cc:1129  hparticlecandfiller.cc:1130  hparticlecandfiller.cc:1131  hparticlecandfiller.cc:1132  hparticlecandfiller.cc:1133  hparticlecandfiller.cc:1134  hparticlecandfiller.cc:1135  hparticlecandfiller.cc:1136  hparticlecandfiller.cc:1137  hparticlecandfiller.cc:1138  hparticlecandfiller.cc:1139  hparticlecandfiller.cc:1140  hparticlecandfiller.cc:1141  hparticlecandfiller.cc:1142  hparticlecandfiller.cc:1143  hparticlecandfiller.cc:1144  hparticlecandfiller.cc:1145  hparticlecandfiller.cc:1146  hparticlecandfiller.cc:1147  hparticlecandfiller.cc:1148  hparticlecandfiller.cc:1149  hparticlecandfiller.cc:1150  hparticlecandfiller.cc:1151  hparticlecandfiller.cc:1152  hparticlecandfiller.cc:1153  hparticlecandfiller.cc:1154  hparticlecandfiller.cc:1155  hparticlecandfiller.cc:1156  hparticlecandfiller.cc:1157  hparticlecandfiller.cc:1158  hparticlecandfiller.cc:1159  hparticlecandfiller.cc:1160  hparticlecandfiller.cc:1161  hparticlecandfiller.cc:1162  hparticlecandfiller.cc:1163  hparticlecandfiller.cc:1164  hparticlecandfiller.cc:1165  hparticlecandfiller.cc:1166  hparticlecandfiller.cc:1167  hparticlecandfiller.cc:1168  hparticlecandfiller.cc:1169  hparticlecandfiller.cc:1170  hparticlecandfiller.cc:1171  hparticlecandfiller.cc:1172  hparticlecandfiller.cc:1173  hparticlecandfiller.cc:1174  hparticlecandfiller.cc:1175  hparticlecandfiller.cc:1176  hparticlecandfiller.cc:1177  hparticlecandfiller.cc:1178  hparticlecandfiller.cc:1179  hparticlecandfiller.cc:1180  hparticlecandfiller.cc:1181  hparticlecandfiller.cc:1182  hparticlecandfiller.cc:1183  hparticlecandfiller.cc:1184  hparticlecandfiller.cc:1185  hparticlecandfiller.cc:1186  hparticlecandfiller.cc:1187  hparticlecandfiller.cc:1188  hparticlecandfiller.cc:1189  hparticlecandfiller.cc:1190  hparticlecandfiller.cc:1191  hparticlecandfiller.cc:1192  hparticlecandfiller.cc:1193  hparticlecandfiller.cc:1194  hparticlecandfiller.cc:1195  hparticlecandfiller.cc:1196  hparticlecandfiller.cc:1197  hparticlecandfiller.cc:1198  hparticlecandfiller.cc:1199  hparticlecandfiller.cc:1200  hparticlecandfiller.cc:1201  hparticlecandfiller.cc:1202  hparticlecandfiller.cc:1203  hparticlecandfiller.cc:1204  hparticlecandfiller.cc:1205  hparticlecandfiller.cc:1206  hparticlecandfiller.cc:1207  hparticlecandfiller.cc:1208  hparticlecandfiller.cc:1209  hparticlecandfiller.cc:1210  hparticlecandfiller.cc:1211  hparticlecandfiller.cc:1212  hparticlecandfiller.cc:1213  hparticlecandfiller.cc:1214  hparticlecandfiller.cc:1215  hparticlecandfiller.cc:1216  hparticlecandfiller.cc:1217  hparticlecandfiller.cc:1218  hparticlecandfiller.cc:1219  hparticlecandfiller.cc:1220  hparticlecandfiller.cc:1221  hparticlecandfiller.cc:1222  hparticlecandfiller.cc:1223  hparticlecandfiller.cc:1224  hparticlecandfiller.cc:1225  hparticlecandfiller.cc:1226  hparticlecandfiller.cc:1227  hparticlecandfiller.cc:1228  hparticlecandfiller.cc:1229  hparticlecandfiller.cc:1230  hparticlecandfiller.cc:1231  hparticlecandfiller.cc:1232  hparticlecandfiller.cc:1233  hparticlecandfiller.cc:1234  hparticlecandfiller.cc:1235  hparticlecandfiller.cc:1236  hparticlecandfiller.cc:1237  hparticlecandfiller.cc:1238  hparticlecandfiller.cc:1239  hparticlecandfiller.cc:1240  hparticlecandfiller.cc:1241  hparticlecandfiller.cc:1242  hparticlecandfiller.cc:1243  hparticlecandfiller.cc:1244  hparticlecandfiller.cc:1245  hparticlecandfiller.cc:1246  hparticlecandfiller.cc:1247  hparticlecandfiller.cc:1248  hparticlecandfiller.cc:1249  hparticlecandfiller.cc:1250  hparticlecandfiller.cc:1251  hparticlecandfiller.cc:1252  hparticlecandfiller.cc:1253  hparticlecandfiller.cc:1254  hparticlecandfiller.cc:1255  hparticlecandfiller.cc:1256  hparticlecandfiller.cc:1257  hparticlecandfiller.cc:1258  hparticlecandfiller.cc:1259  hparticlecandfiller.cc:1260  hparticlecandfiller.cc:1261  hparticlecandfiller.cc:1262  hparticlecandfiller.cc:1263  hparticlecandfiller.cc:1264  hparticlecandfiller.cc:1265  hparticlecandfiller.cc:1266  hparticlecandfiller.cc:1267  hparticlecandfiller.cc:1268  hparticlecandfiller.cc:1269  hparticlecandfiller.cc:1270  hparticlecandfiller.cc:1271  hparticlecandfiller.cc:1272  hparticlecandfiller.cc:1273  hparticlecandfiller.cc:1274  hparticlecandfiller.cc:1275  hparticlecandfiller.cc:1276  hparticlecandfiller.cc:1277  hparticlecandfiller.cc:1278  hparticlecandfiller.cc:1279  hparticlecandfiller.cc:1280  hparticlecandfiller.cc:1281  hparticlecandfiller.cc:1282  hparticlecandfiller.cc:1283  hparticlecandfiller.cc:1284  hparticlecandfiller.cc:1285  hparticlecandfiller.cc:1286  hparticlecandfiller.cc:1287  hparticlecandfiller.cc:1288  hparticlecandfiller.cc:1289  hparticlecandfiller.cc:1290  hparticlecandfiller.cc:1291  hparticlecandfiller.cc:1292  hparticlecandfiller.cc:1293  hparticlecandfiller.cc:1294  hparticlecandfiller.cc:1295  hparticlecandfiller.cc:1296  hparticlecandfiller.cc:1297  hparticlecandfiller.cc:1298  hparticlecandfiller.cc:1299  hparticlecandfiller.cc:1300  hparticlecandfiller.cc:1301  hparticlecandfiller.cc:1302  hparticlecandfiller.cc:1303  hparticlecandfiller.cc:1304  hparticlecandfiller.cc:1305  hparticlecandfiller.cc:1306  hparticlecandfiller.cc:1307  hparticlecandfiller.cc:1308  hparticlecandfiller.cc:1309  hparticlecandfiller.cc:1310  hparticlecandfiller.cc:1311  hparticlecandfiller.cc:1312  hparticlecandfiller.cc:1313  hparticlecandfiller.cc:1314  hparticlecandfiller.cc:1315  hparticlecandfiller.cc:1316  hparticlecandfiller.cc:1317  hparticlecandfiller.cc:1318  hparticlecandfiller.cc:1319  hparticlecandfiller.cc:1320  hparticlecandfiller.cc:1321  hparticlecandfiller.cc:1322  hparticlecandfiller.cc:1323  hparticlecandfiller.cc:1324  hparticlecandfiller.cc:1325  hparticlecandfiller.cc:1326  hparticlecandfiller.cc:1327  hparticlecandfiller.cc:1328  hparticlecandfiller.cc:1329  hparticlecandfiller.cc:1330  hparticlecandfiller.cc:1331  hparticlecandfiller.cc:1332  hparticlecandfiller.cc:1333  hparticlecandfiller.cc:1334  hparticlecandfiller.cc:1335  hparticlecandfiller.cc:1336  hparticlecandfiller.cc:1337  hparticlecandfiller.cc:1338  hparticlecandfiller.cc:1339  hparticlecandfiller.cc:1340  hparticlecandfiller.cc:1341  hparticlecandfiller.cc:1342  hparticlecandfiller.cc:1343  hparticlecandfiller.cc:1344  hparticlecandfiller.cc:1345  hparticlecandfiller.cc:1346  hparticlecandfiller.cc:1347  hparticlecandfiller.cc:1348  hparticlecandfiller.cc:1349  hparticlecandfiller.cc:1350  hparticlecandfiller.cc:1351  hparticlecandfiller.cc:1352  hparticlecandfiller.cc:1353  hparticlecandfiller.cc:1354  hparticlecandfiller.cc:1355  hparticlecandfiller.cc:1356  hparticlecandfiller.cc:1357  hparticlecandfiller.cc:1358  hparticlecandfiller.cc:1359  hparticlecandfiller.cc:1360  hparticlecandfiller.cc:1361  hparticlecandfiller.cc:1362  hparticlecandfiller.cc:1363  hparticlecandfiller.cc:1364  hparticlecandfiller.cc:1365  hparticlecandfiller.cc:1366  hparticlecandfiller.cc:1367  hparticlecandfiller.cc:1368  hparticlecandfiller.cc:1369  hparticlecandfiller.cc:1370  hparticlecandfiller.cc:1371  hparticlecandfiller.cc:1372  hparticlecandfiller.cc:1373  hparticlecandfiller.cc:1374  hparticlecandfiller.cc:1375  hparticlecandfiller.cc:1376  hparticlecandfiller.cc:1377  hparticlecandfiller.cc:1378  hparticlecandfiller.cc:1379  hparticlecandfiller.cc:1380  hparticlecandfiller.cc:1381  hparticlecandfiller.cc:1382  hparticlecandfiller.cc:1383  hparticlecandfiller.cc:1384  hparticlecandfiller.cc:1385  hparticlecandfiller.cc:1386  hparticlecandfiller.cc:1387  hparticlecandfiller.cc:1388  hparticlecandfiller.cc:1389  hparticlecandfiller.cc:1390  hparticlecandfiller.cc:1391  hparticlecandfiller.cc:1392  hparticlecandfiller.cc:1393  hparticlecandfiller.cc:1394  hparticlecandfiller.cc:1395  hparticlecandfiller.cc:1396  hparticlecandfiller.cc:1397  hparticlecandfiller.cc:1398  hparticlecandfiller.cc:1399  hparticlecandfiller.cc:1400  hparticlecandfiller.cc:1401  hparticlecandfiller.cc:1402  hparticlecandfiller.cc:1403  hparticlecandfiller.cc:1404  hparticlecandfiller.cc:1405  hparticlecandfiller.cc:1406  hparticlecandfiller.cc:1407  hparticlecandfiller.cc:1408  hparticlecandfiller.cc:1409  hparticlecandfiller.cc:1410  hparticlecandfiller.cc:1411  hparticlecandfiller.cc:1412  hparticlecandfiller.cc:1413  hparticlecandfiller.cc:1414  hparticlecandfiller.cc:1415  hparticlecandfiller.cc:1416  hparticlecandfiller.cc:1417  hparticlecandfiller.cc:1418  hparticlecandfiller.cc:1419  hparticlecandfiller.cc:1420  hparticlecandfiller.cc:1421  hparticlecandfiller.cc:1422  hparticlecandfiller.cc:1423  hparticlecandfiller.cc:1424  hparticlecandfiller.cc:1425  hparticlecandfiller.cc:1426  hparticlecandfiller.cc:1427  hparticlecandfiller.cc:1428  hparticlecandfiller.cc:1429  hparticlecandfiller.cc:1430  hparticlecandfiller.cc:1431  hparticlecandfiller.cc:1432  hparticlecandfiller.cc:1433  hparticlecandfiller.cc:1434  hparticlecandfiller.cc:1435  hparticlecandfiller.cc:1436  hparticlecandfiller.cc:1437  hparticlecandfiller.cc:1438  hparticlecandfiller.cc:1439  hparticlecandfiller.cc:1440  hparticlecandfiller.cc:1441  hparticlecandfiller.cc:1442  hparticlecandfiller.cc:1443  hparticlecandfiller.cc:1444  hparticlecandfiller.cc:1445  hparticlecandfiller.cc:1446  hparticlecandfiller.cc:1447  hparticlecandfiller.cc:1448  hparticlecandfiller.cc:1449  hparticlecandfiller.cc:1450  hparticlecandfiller.cc:1451  hparticlecandfiller.cc:1452  hparticlecandfiller.cc:1453  hparticlecandfiller.cc:1454  hparticlecandfiller.cc:1455  hparticlecandfiller.cc:1456  hparticlecandfiller.cc:1457  hparticlecandfiller.cc:1458  hparticlecandfiller.cc:1459  hparticlecandfiller.cc:1460  hparticlecandfiller.cc:1461  hparticlecandfiller.cc:1462  hparticlecandfiller.cc:1463  hparticlecandfiller.cc:1464  hparticlecandfiller.cc:1465  hparticlecandfiller.cc:1466  hparticlecandfiller.cc:1467  hparticlecandfiller.cc:1468  hparticlecandfiller.cc:1469  hparticlecandfiller.cc:1470  hparticlecandfiller.cc:1471  hparticlecandfiller.cc:1472  hparticlecandfiller.cc:1473  hparticlecandfiller.cc:1474  hparticlecandfiller.cc:1475  hparticlecandfiller.cc:1476  hparticlecandfiller.cc:1477  hparticlecandfiller.cc:1478  hparticlecandfiller.cc:1479  hparticlecandfiller.cc:1480  hparticlecandfiller.cc:1481  hparticlecandfiller.cc:1482  hparticlecandfiller.cc:1483  hparticlecandfiller.cc:1484  hparticlecandfiller.cc:1485  hparticlecandfiller.cc:1486  hparticlecandfiller.cc:1487  hparticlecandfiller.cc:1488  hparticlecandfiller.cc:1489  hparticlecandfiller.cc:1490  hparticlecandfiller.cc:1491  hparticlecandfiller.cc:1492  hparticlecandfiller.cc:1493  hparticlecandfiller.cc:1494  hparticlecandfiller.cc:1495  hparticlecandfiller.cc:1496  hparticlecandfiller.cc:1497  hparticlecandfiller.cc:1498  hparticlecandfiller.cc:1499  hparticlecandfiller.cc:1500  hparticlecandfiller.cc:1501  hparticlecandfiller.cc:1502  hparticlecandfiller.cc:1503  hparticlecandfiller.cc:1504  hparticlecandfiller.cc:1505  hparticlecandfiller.cc:1506  hparticlecandfiller.cc:1507  hparticlecandfiller.cc:1508  hparticlecandfiller.cc:1509  hparticlecandfiller.cc:1510  hparticlecandfiller.cc:1511  hparticlecandfiller.cc:1512  hparticlecandfiller.cc:1513  hparticlecandfiller.cc:1514  hparticlecandfiller.cc:1515  hparticlecandfiller.cc:1516  hparticlecandfiller.cc:1517  hparticlecandfiller.cc:1518  hparticlecandfiller.cc:1519  hparticlecandfiller.cc:1520  hparticlecandfiller.cc:1521  hparticlecandfiller.cc:1522  hparticlecandfiller.cc:1523  hparticlecandfiller.cc:1524  hparticlecandfiller.cc:1525  hparticlecandfiller.cc:1526  hparticlecandfiller.cc:1527  hparticlecandfiller.cc:1528  hparticlecandfiller.cc:1529  hparticlecandfiller.cc:1530  hparticlecandfiller.cc:1531  hparticlecandfiller.cc:1532  hparticlecandfiller.cc:1533  hparticlecandfiller.cc:1534  hparticlecandfiller.cc:1535  hparticlecandfiller.cc:1536  hparticlecandfiller.cc:1537  hparticlecandfiller.cc:1538  hparticlecandfiller.cc:1539  hparticlecandfiller.cc:1540  hparticlecandfiller.cc:1541  hparticlecandfiller.cc:1542  hparticlecandfiller.cc:1543  hparticlecandfiller.cc:1544  hparticlecandfiller.cc:1545  hparticlecandfiller.cc:1546  hparticlecandfiller.cc:1547  hparticlecandfiller.cc:1548  hparticlecandfiller.cc:1549  hparticlecandfiller.cc:1550  hparticlecandfiller.cc:1551  hparticlecandfiller.cc:1552  hparticlecandfiller.cc:1553  hparticlecandfiller.cc:1554  hparticlecandfiller.cc:1555  hparticlecandfiller.cc:1556  hparticlecandfiller.cc:1557  hparticlecandfiller.cc:1558  hparticlecandfiller.cc:1559  hparticlecandfiller.cc:1560  hparticlecandfiller.cc:1561  hparticlecandfiller.cc:1562  hparticlecandfiller.cc:1563  hparticlecandfiller.cc:1564  hparticlecandfiller.cc:1565  hparticlecandfiller.cc:1566  hparticlecandfiller.cc:1567  hparticlecandfiller.cc:1568  hparticlecandfiller.cc:1569  hparticlecandfiller.cc:1570  hparticlecandfiller.cc:1571  hparticlecandfiller.cc:1572  hparticlecandfiller.cc:1573  hparticlecandfiller.cc:1574  hparticlecandfiller.cc:1575  hparticlecandfiller.cc:1576  hparticlecandfiller.cc:1577  hparticlecandfiller.cc:1578  hparticlecandfiller.cc:1579  hparticlecandfiller.cc:1580  hparticlecandfiller.cc:1581  hparticlecandfiller.cc:1582  hparticlecandfiller.cc:1583  hparticlecandfiller.cc:1584  hparticlecandfiller.cc:1585  hparticlecandfiller.cc:1586  hparticlecandfiller.cc:1587  hparticlecandfiller.cc:1588  hparticlecandfiller.cc:1589  hparticlecandfiller.cc:1590  hparticlecandfiller.cc:1591  hparticlecandfiller.cc:1592  hparticlecandfiller.cc:1593  hparticlecandfiller.cc:1594  hparticlecandfiller.cc:1595  hparticlecandfiller.cc:1596  hparticlecandfiller.cc:1597  hparticlecandfiller.cc:1598  hparticlecandfiller.cc:1599  hparticlecandfiller.cc:1600  hparticlecandfiller.cc:1601  hparticlecandfiller.cc:1602  hparticlecandfiller.cc:1603  hparticlecandfiller.cc:1604  hparticlecandfiller.cc:1605  hparticlecandfiller.cc:1606  hparticlecandfiller.cc:1607  hparticlecandfiller.cc:1608  hparticlecandfiller.cc:1609  hparticlecandfiller.cc:1610  hparticlecandfiller.cc:1611  hparticlecandfiller.cc:1612  hparticlecandfiller.cc:1613  hparticlecandfiller.cc:1614  hparticlecandfiller.cc:1615  hparticlecandfiller.cc:1616  hparticlecandfiller.cc:1617  hparticlecandfiller.cc:1618  hparticlecandfiller.cc:1619  hparticlecandfiller.cc:1620  hparticlecandfiller.cc:1621  hparticlecandfiller.cc:1622  hparticlecandfiller.cc:1623  hparticlecandfiller.cc:1624  hparticlecandfiller.cc:1625  hparticlecandfiller.cc:1626  hparticlecandfiller.cc:1627  hparticlecandfiller.cc:1628  hparticlecandfiller.cc:1629  hparticlecandfiller.cc:1630  hparticlecandfiller.cc:1631  hparticlecandfiller.cc:1632  hparticlecandfiller.cc:1633  hparticlecandfiller.cc:1634  hparticlecandfiller.cc:1635  hparticlecandfiller.cc:1636  hparticlecandfiller.cc:1637  hparticlecandfiller.cc:1638  hparticlecandfiller.cc:1639  hparticlecandfiller.cc:1640  hparticlecandfiller.cc:1641  hparticlecandfiller.cc:1642  hparticlecandfiller.cc:1643  hparticlecandfiller.cc:1644  hparticlecandfiller.cc:1645  hparticlecandfiller.cc:1646  hparticlecandfiller.cc:1647  hparticlecandfiller.cc:1648  hparticlecandfiller.cc:1649  hparticlecandfiller.cc:1650  hparticlecandfiller.cc:1651  hparticlecandfiller.cc:1652  hparticlecandfiller.cc:1653  hparticlecandfiller.cc:1654  hparticlecandfiller.cc:1655  hparticlecandfiller.cc:1656  hparticlecandfiller.cc:1657  hparticlecandfiller.cc:1658  hparticlecandfiller.cc:1659  hparticlecandfiller.cc:1660  hparticlecandfiller.cc:1661  hparticlecandfiller.cc:1662  hparticlecandfiller.cc:1663  hparticlecandfiller.cc:1664  hparticlecandfiller.cc:1665  hparticlecandfiller.cc:1666  hparticlecandfiller.cc:1667  hparticlecandfiller.cc:1668  hparticlecandfiller.cc:1669  hparticlecandfiller.cc:1670  hparticlecandfiller.cc:1671  hparticlecandfiller.cc:1672  hparticlecandfiller.cc:1673  hparticlecandfiller.cc:1674  hparticlecandfiller.cc:1675  hparticlecandfiller.cc:1676  hparticlecandfiller.cc:1677  hparticlecandfiller.cc:1678  hparticlecandfiller.cc:1679  hparticlecandfiller.cc:1680  hparticlecandfiller.cc:1681  hparticlecandfiller.cc:1682  hparticlecandfiller.cc:1683  hparticlecandfiller.cc:1684  hparticlecandfiller.cc:1685  hparticlecandfiller.cc:1686  hparticlecandfiller.cc:1687  hparticlecandfiller.cc:1688  hparticlecandfiller.cc:1689  hparticlecandfiller.cc:1690  hparticlecandfiller.cc:1691  hparticlecandfiller.cc:1692  hparticlecandfiller.cc:1693  hparticlecandfiller.cc:1694  hparticlecandfiller.cc:1695  hparticlecandfiller.cc:1696  hparticlecandfiller.cc:1697  hparticlecandfiller.cc:1698  hparticlecandfiller.cc:1699  hparticlecandfiller.cc:1700  hparticlecandfiller.cc:1701  hparticlecandfiller.cc:1702  hparticlecandfiller.cc:1703  hparticlecandfiller.cc:1704  hparticlecandfiller.cc:1705  hparticlecandfiller.cc:1706  hparticlecandfiller.cc:1707  hparticlecandfiller.cc:1708  hparticlecandfiller.cc:1709  hparticlecandfiller.cc:1710  hparticlecandfiller.cc:1711  hparticlecandfiller.cc:1712  hparticlecandfiller.cc:1713  hparticlecandfiller.cc:1714  hparticlecandfiller.cc:1715  hparticlecandfiller.cc:1716  hparticlecandfiller.cc:1717  hparticlecandfiller.cc:1718  hparticlecandfiller.cc:1719  hparticlecandfiller.cc:1720  hparticlecandfiller.cc:1721  hparticlecandfiller.cc:1722  hparticlecandfiller.cc:1723  hparticlecandfiller.cc:1724  hparticlecandfiller.cc:1725  hparticlecandfiller.cc:1726  hparticlecandfiller.cc:1727  hparticlecandfiller.cc:1728  hparticlecandfiller.cc:1729  hparticlecandfiller.cc:1730  hparticlecandfiller.cc:1731  hparticlecandfiller.cc:1732  hparticlecandfiller.cc:1733  hparticlecandfiller.cc:1734  hparticlecandfiller.cc:1735  hparticlecandfiller.cc:1736  hparticlecandfiller.cc:1737  hparticlecandfiller.cc:1738  hparticlecandfiller.cc:1739  hparticlecandfiller.cc:1740  hparticlecandfiller.cc:1741  hparticlecandfiller.cc:1742  hparticlecandfiller.cc:1743  hparticlecandfiller.cc:1744  hparticlecandfiller.cc:1745  hparticlecandfiller.cc:1746  hparticlecandfiller.cc:1747  hparticlecandfiller.cc:1748  hparticlecandfiller.cc:1749  hparticlecandfiller.cc:1750  hparticlecandfiller.cc:1751  hparticlecandfiller.cc:1752  hparticlecandfiller.cc:1753  hparticlecandfiller.cc:1754  hparticlecandfiller.cc:1755  hparticlecandfiller.cc:1756  hparticlecandfiller.cc:1757  hparticlecandfiller.cc:1758  hparticlecandfiller.cc:1759  hparticlecandfiller.cc:1760  hparticlecandfiller.cc:1761  hparticlecandfiller.cc:1762  hparticlecandfiller.cc:1763  hparticlecandfiller.cc:1764  hparticlecandfiller.cc:1765  hparticlecandfiller.cc:1766  hparticlecandfiller.cc:1767  hparticlecandfiller.cc:1768  hparticlecandfiller.cc:1769  hparticlecandfiller.cc:1770  hparticlecandfiller.cc:1771  hparticlecandfiller.cc:1772  hparticlecandfiller.cc:1773  hparticlecandfiller.cc:1774  hparticlecandfiller.cc:1775  hparticlecandfiller.cc:1776  hparticlecandfiller.cc:1777  hparticlecandfiller.cc:1778  hparticlecandfiller.cc:1779  hparticlecandfiller.cc:1780  hparticlecandfiller.cc:1781  hparticlecandfiller.cc:1782  hparticlecandfiller.cc:1783  hparticlecandfiller.cc:1784  hparticlecandfiller.cc:1785  hparticlecandfiller.cc:1786  hparticlecandfiller.cc:1787  hparticlecandfiller.cc:1788  hparticlecandfiller.cc:1789  hparticlecandfiller.cc:1790  hparticlecandfiller.cc:1791  hparticlecandfiller.cc:1792  hparticlecandfiller.cc:1793  hparticlecandfiller.cc:1794  hparticlecandfiller.cc:1795  hparticlecandfiller.cc:1796  hparticlecandfiller.cc:1797  hparticlecandfiller.cc:1798  hparticlecandfiller.cc:1799  hparticlecandfiller.cc:1800  hparticlecandfiller.cc:1801  hparticlecandfiller.cc:1802  hparticlecandfiller.cc:1803  hparticlecandfiller.cc:1804  hparticlecandfiller.cc:1805  hparticlecandfiller.cc:1806  hparticlecandfiller.cc:1807  hparticlecandfiller.cc:1808  hparticlecandfiller.cc:1809  hparticlecandfiller.cc:1810  hparticlecandfiller.cc:1811  hparticlecandfiller.cc:1812  hparticlecandfiller.cc:1813  hparticlecandfiller.cc:1814  hparticlecandfiller.cc:1815  hparticlecandfiller.cc:1816  hparticlecandfiller.cc:1817  hparticlecandfiller.cc:1818  hparticlecandfiller.cc:1819  hparticlecandfiller.cc:1820  hparticlecandfiller.cc:1821  hparticlecandfiller.cc:1822  hparticlecandfiller.cc:1823  hparticlecandfiller.cc:1824  hparticlecandfiller.cc:1825  hparticlecandfiller.cc:1826  hparticlecandfiller.cc:1827  hparticlecandfiller.cc:1828  hparticlecandfiller.cc:1829  hparticlecandfiller.cc:1830  hparticlecandfiller.cc:1831  hparticlecandfiller.cc:1832  hparticlecandfiller.cc:1833  hparticlecandfiller.cc:1834  hparticlecandfiller.cc:1835  hparticlecandfiller.cc:1836  hparticlecandfiller.cc:1837  hparticlecandfiller.cc:1838  hparticlecandfiller.cc:1839  hparticlecandfiller.cc:1840  hparticlecandfiller.cc:1841  hparticlecandfiller.cc:1842  hparticlecandfiller.cc:1843  hparticlecandfiller.cc:1844  hparticlecandfiller.cc:1845  hparticlecandfiller.cc:1846  hparticlecandfiller.cc:1847  hparticlecandfiller.cc:1848  hparticlecandfiller.cc:1849  hparticlecandfiller.cc:1850  hparticlecandfiller.cc:1851  hparticlecandfiller.cc:1852  hparticlecandfiller.cc:1853  hparticlecandfiller.cc:1854  hparticlecandfiller.cc:1855  hparticlecandfiller.cc:1856  hparticlecandfiller.cc:1857  hparticlecandfiller.cc:1858  hparticlecandfiller.cc:1859  hparticlecandfiller.cc:1860  hparticlecandfiller.cc:1861  hparticlecandfiller.cc:1862  hparticlecandfiller.cc:1863  hparticlecandfiller.cc:1864  hparticlecandfiller.cc:1865  hparticlecandfiller.cc:1866  hparticlecandfiller.cc:1867  hparticlecandfiller.cc:1868  hparticlecandfiller.cc:1869  hparticlecandfiller.cc:1870  hparticlecandfiller.cc:1871  hparticlecandfiller.cc:1872  hparticlecandfiller.cc:1873  hparticlecandfiller.cc:1874  hparticlecandfiller.cc:1875  hparticlecandfiller.cc:1876  hparticlecandfiller.cc:1877  hparticlecandfiller.cc:1878  hparticlecandfiller.cc:1879  hparticlecandfiller.cc:1880  hparticlecandfiller.cc:1881  hparticlecandfiller.cc:1882  hparticlecandfiller.cc:1883  hparticlecandfiller.cc:1884  hparticlecandfiller.cc:1885  hparticlecandfiller.cc:1886  hparticlecandfiller.cc:1887  hparticlecandfiller.cc:1888  hparticlecandfiller.cc:1889  hparticlecandfiller.cc:1890  hparticlecandfiller.cc:1891  hparticlecandfiller.cc:1892  hparticlecandfiller.cc:1893  hparticlecandfiller.cc:1894  hparticlecandfiller.cc:1895  hparticlecandfiller.cc:1896  hparticlecandfiller.cc:1897  hparticlecandfiller.cc:1898  hparticlecandfiller.cc:1899  hparticlecandfiller.cc:1900  hparticlecandfiller.cc:1901  hparticlecandfiller.cc:1902  hparticlecandfiller.cc:1903  hparticlecandfiller.cc:1904  hparticlecandfiller.cc:1905  hparticlecandfiller.cc:1906  hparticlecandfiller.cc:1907  hparticlecandfiller.cc:1908  hparticlecandfiller.cc:1909  hparticlecandfiller.cc:1910  hparticlecandfiller.cc:1911  hparticlecandfiller.cc:1912  hparticlecandfiller.cc:1913  hparticlecandfiller.cc:1914  hparticlecandfiller.cc:1915  hparticlecandfiller.cc:1916  hparticlecandfiller.cc:1917  hparticlecandfiller.cc:1918  hparticlecandfiller.cc:1919  hparticlecandfiller.cc:1920  hparticlecandfiller.cc:1921  hparticlecandfiller.cc:1922  hparticlecandfiller.cc:1923  hparticlecandfiller.cc:1924  hparticlecandfiller.cc:1925  hparticlecandfiller.cc:1926  hparticlecandfiller.cc:1927  hparticlecandfiller.cc:1928  hparticlecandfiller.cc:1929  hparticlecandfiller.cc:1930  hparticlecandfiller.cc:1931  hparticlecandfiller.cc:1932  hparticlecandfiller.cc:1933  hparticlecandfiller.cc:1934  hparticlecandfiller.cc:1935  hparticlecandfiller.cc:1936  hparticlecandfiller.cc:1937  hparticlecandfiller.cc:1938  hparticlecandfiller.cc:1939  hparticlecandfiller.cc:1940  hparticlecandfiller.cc:1941  hparticlecandfiller.cc:1942  hparticlecandfiller.cc:1943  hparticlecandfiller.cc:1944  hparticlecandfiller.cc:1945  hparticlecandfiller.cc:1946  hparticlecandfiller.cc:1947  hparticlecandfiller.cc:1948  hparticlecandfiller.cc:1949  hparticlecandfiller.cc:1950  hparticlecandfiller.cc:1951  hparticlecandfiller.cc:1952  hparticlecandfiller.cc:1953  hparticlecandfiller.cc:1954  hparticlecandfiller.cc:1955  hparticlecandfiller.cc:1956  hparticlecandfiller.cc:1957  hparticlecandfiller.cc:1958  hparticlecandfiller.cc:1959  hparticlecandfiller.cc:1960  hparticlecandfiller.cc:1961  hparticlecandfiller.cc:1962  hparticlecandfiller.cc:1963  hparticlecandfiller.cc:1964  hparticlecandfiller.cc:1965  hparticlecandfiller.cc:1966  hparticlecandfiller.cc:1967  hparticlecandfiller.cc:1968  hparticlecandfiller.cc:1969  hparticlecandfiller.cc:1970  hparticlecandfiller.cc:1971  hparticlecandfiller.cc:1972  hparticlecandfiller.cc:1973  hparticlecandfiller.cc:1974  hparticlecandfiller.cc:1975  hparticlecandfiller.cc:1976  hparticlecandfiller.cc:1977  hparticlecandfiller.cc:1978  hparticlecandfiller.cc:1979  hparticlecandfiller.cc:1980  hparticlecandfiller.cc:1981  hparticlecandfiller.cc:1982  hparticlecandfiller.cc:1983  hparticlecandfiller.cc:1984  hparticlecandfiller.cc:1985  hparticlecandfiller.cc:1986  hparticlecandfiller.cc:1987  hparticlecandfiller.cc:1988  hparticlecandfiller.cc:1989  hparticlecandfiller.cc:1990  hparticlecandfiller.cc:1991  hparticlecandfiller.cc:1992  hparticlecandfiller.cc:1993  hparticlecandfiller.cc:1994  hparticlecandfiller.cc:1995  hparticlecandfiller.cc:1996  hparticlecandfiller.cc:1997  hparticlecandfiller.cc:1998  hparticlecandfiller.cc:1999  hparticlecandfiller.cc:2000  hparticlecandfiller.cc:2001  hparticlecandfiller.cc:2002  hparticlecandfiller.cc:2003  hparticlecandfiller.cc:2004  hparticlecandfiller.cc:2005  hparticlecandfiller.cc:2006  hparticlecandfiller.cc:2007  hparticlecandfiller.cc:2008  hparticlecandfiller.cc:2009  hparticlecandfiller.cc:2010  hparticlecandfiller.cc:2011  hparticlecandfiller.cc:2012  hparticlecandfiller.cc:2013  hparticlecandfiller.cc:2014  hparticlecandfiller.cc:2015  hparticlecandfiller.cc:2016  hparticlecandfiller.cc:2017  hparticlecandfiller.cc:2018  hparticlecandfiller.cc:2019  hparticlecandfiller.cc:2020  hparticlecandfiller.cc:2021  hparticlecandfiller.cc:2022  hparticlecandfiller.cc:2023  hparticlecandfiller.cc:2024  hparticlecandfiller.cc:2025  hparticlecandfiller.cc:2026  hparticlecandfiller.cc:2027  hparticlecandfiller.cc:2028  hparticlecandfiller.cc:2029  hparticlecandfiller.cc:2030  hparticlecandfiller.cc:2031  hparticlecandfiller.cc:2032  hparticlecandfiller.cc:2033  hparticlecandfiller.cc:2034  hparticlecandfiller.cc:2035  hparticlecandfiller.cc:2036  hparticlecandfiller.cc:2037  hparticlecandfiller.cc:2038  hparticlecandfiller.cc:2039  hparticlecandfiller.cc:2040  hparticlecandfiller.cc:2041  hparticlecandfiller.cc:2042  hparticlecandfiller.cc:2043  hparticlecandfiller.cc:2044  hparticlecandfiller.cc:2045  hparticlecandfiller.cc:2046  hparticlecandfiller.cc:2047  hparticlecandfiller.cc:2048  hparticlecandfiller.cc:2049  hparticlecandfiller.cc:2050  hparticlecandfiller.cc:2051  hparticlecandfiller.cc:2052  hparticlecandfiller.cc:2053  hparticlecandfiller.cc:2054  hparticlecandfiller.cc:2055  hparticlecandfiller.cc:2056  hparticlecandfiller.cc:2057  hparticlecandfiller.cc:2058  hparticlecandfiller.cc:2059  hparticlecandfiller.cc:2060  hparticlecandfiller.cc:2061  hparticlecandfiller.cc:2062  hparticlecandfiller.cc:2063  hparticlecandfiller.cc:2064  hparticlecandfiller.cc:2065  hparticlecandfiller.cc:2066  hparticlecandfiller.cc:2067  hparticlecandfiller.cc:2068  hparticlecandfiller.cc:2069  hparticlecandfiller.cc:2070  hparticlecandfiller.cc:2071  hparticlecandfiller.cc:2072  hparticlecandfiller.cc:2073