#include "hparticletree.h"
#include "hades.h"
#include "hdatasource.h"
#include "hrecevent.h"
#include "hpartialevent.h"
#include "hcategory.h"
#include "hlinearcategory.h"
#include "hmatrixcategory.h"
#include "hcategorymanager.h"
#include "hparticletool.h"
#include "hsrckeeper.h"
#include "hgeantmedia.h"
#include "hmdclayer.h"
#include "haddef.h"
#include "hstartdef.h"
#include "richdef.h"
#include "hmdcdef.h"
#include "hmdctrackddef.h"
#include "hmdctrackgdef.h"
#include "rpcdef.h"
#include "tofdef.h"
#include "emcdef.h"
#include "showerdef.h"
#include "walldef.h"
#include "hpiontrackerdef.h"
#include "hgeantdef.h"
#include "hparticlecand.h"
#include "hparticlecandsim.h"
#include "hparticleevtinfo.h"
#include "hparticlemdc.h"
#include "hmetamatch2.h"
#include "hstart2hit.h"
#include "hstart2cal.h"
#include "htboxchan.h"
#include "hrichhit.h"
#include "hrichhitsim.h"
#include "hrichdirclus.h"
#include "hrichcal.h"
#include "hrichcalsim.h"
#include "htofhit.h"
#include "htofhitsim.h"
#include "htofcluster.h"
#include "htofclustersim.h"
#include "hrpccluster.h"
#include "hrpcclustersim.h"
#include "hshowerhit.h"
#include "hshowerhitsim.h"
#include "hemccluster.h"
#include "hemcclustersim.h"
#include "hemccal.h"
#include "hwallhit.h"
#include "hwallhitsim.h"
#include "hwalleventplane.h"
#include "hpiontrackertrack.h"
#include "hpiontrackerhit.h"
#include "hpiontrackercal.h"
#include "hpiontrackerraw.h"
#include "hmdctrkcand.h"
#include "hmdcseg.h"
#include "hmdcsegsim.h"
#include "hmdchit.h"
#include "hmdchitsim.h"
#include "hmdcclus.h"
#include "hmdcclussim.h"
#include "hmdcclusinf.h"
#include "hmdcclusinfsim.h"
#include "hmdcclusfit.h"
#include "hmdcclusfitsim.h"
#include "hmdcwirefit.h"
#include "hmdcwirefitsim.h"
#include "hmdccal1.h"
#include "hmdccal1sim.h"
#include "hgeantkine.h"
#include "hgeantrich.h"
#include "hgeantmdc.h"
#include "hgeanttof.h"
#include "hgeantrpc.h"
#include "hgeantshower.h"
#include "hgeantwall.h"
#include "hgeantstart.h"
#include "hgeantemc.h"
#include "TBranch.h"
#include "TSystem.h"
#include "TObjArray.h"
#include <iostream>
#include <iomanip>
#include <vector>
#include <algorithm>
using namespace std;
ClassImp(HParticleTree)
HParticleTree::HParticleTree(const Text_t *name,const Text_t *title)
: HReconstructor(name,title)
{
fCurrentEvent = 0;
fOutputFile = 0;
fTree = 0;
fCycleNumber = 0;
fOutputFileName = "";
fOutputTitle = "";
fOutputOption = "RECREATE";
fOutputCompression = 0;
fOutputFileSuffix = "_filter_tree";
fOutputDir = "";
fOutFound = kFALSE;
kSkipEmptyEvents = kFALSE;
kSkipTracks = kTRUE;
kDoSorter = kTRUE;
sortType = Particle::kIsBestRKSorter;
keeper = 0;
media = 0;
mdclayer = 0;
Cat_t supported [] = {
catParticleCand,
catParticleMdc,
catParticleEvtInfo,
catStart2Hit,
catStart2Cal,
catTBoxChan,
catWallHit,
catWallEventPlane,
catPionTrackerRaw,
catPionTrackerCal,
catPionTrackerHit,
catPionTrackerTrack,
catTofHit,
catTofCluster,
catRpcCluster,
catShowerHit,
catEmcCluster,
catEmcCal,
catRichHit,
catRichDirClus,
catRichCal,
catMdcSeg,
catMdcHit,
catMdcCal1,
catMdcClus,
catMdcClusInf,
catMdcClusFit,
catMdcWireFit,
catGeantKine,
catMdcGeantRaw,
catTofGeantRaw,
catRpcGeantRaw,
catShowerGeantRaw,
catEmcGeantRaw,
catWallGeantRaw,
catRichGeantRaw,
catRichGeantRaw+1,
catRichGeantRaw+2,
catEmcGeantRaw,
catStartGeantRaw
};
for(UInt_t i = 0 ; i < sizeof(supported) / sizeof(Cat_t); i++){
fmCatNumSupport.push_back(supported[i]);
}
sort(fmCatNumSupport.begin(),fmCatNumSupport.end());
Cat_t supportFullCopy [] ={
catRichHit,
catShowerHit,
catTofHit,
catTofCluster,
catRpcCluster
};
for(UInt_t i = 0 ; i < sizeof(supportFullCopy) / sizeof(Cat_t); i++){
fmCatNumFullCopySupport.push_back(supportFullCopy[i]);
}
sort(fmCatNumFullCopySupport.begin(),fmCatNumFullCopySupport.end());
pUserSelectEvent = NULL;
pUserSelectLeptons = NULL;
pUserSelectHadrons = NULL;
pUserKeepTrack = NULL;
fParamSelectEvent = NULL;
}
HParticleTree::~HParticleTree()
{
}
Bool_t HParticleTree::init()
{
sorter.init();
return kTRUE;
}
Int_t HParticleTree::execute()
{
memset(ctMdcSeg ,0,6*2*sizeof(Int_t));
memset(ctMdcClus,0,6*2*sizeof(Int_t));
memset(ctMdcHit ,0,6*4*sizeof(Int_t));
memset(ctRpcClus,0,6*sizeof(Int_t));
if(!fOutFound)
{
const HRecEvent& event = *((HRecEvent*)gHades->getCurrentEvent());
if(fOutputFileName != ""){
Info("execute()","CREATING OUTPUT FILE FROM SPECIFIED FILENAME= %s !",fOutputFileName.Data());
fOutputFile = new TFile(fOutputFileName.Data(),"RECREATE",fOutputTitle.Data(),fOutputCompression);
fOutFound = kTRUE;
setEvent();
makeTree();
} else {
if( fOutputDir == "" ) fOutputDir = gSystem->WorkingDirectory();
if( fOutputDir.EndsWith("/") ) fOutputDir.Replace(fOutputDir.Length() -1,1,"");
TFile* f = gHades->getOutputFile();
if(f){
fOutputOption = f->GetOption();
fOutputTitle = f->GetTitle();
fOutputCompression = f->GetCompressionLevel();
TString fname = f->GetName();
fname = gSystem->BaseName(fname.Data());
fname.ReplaceAll(".root","");
fOutputFileName = Form("%s/%s%s.root",fOutputDir.Data(),fname.Data(),fOutputFileSuffix.Data());
Info("execute()","CREATING OUTPUT FILE FROM HADES OUTPUT= %s !",fOutputFileName.Data());
fOutputFile = new TFile(fOutputFileName.Data(),"RECREATE",fOutputTitle.Data(),fOutputCompression);
fOutFound = kTRUE;
fCurrentEvent = new HRecEvent(event);
setEvent();
makeTree();
} else {
HDataSource* source = gHades->getDataSource();
if(source){
TString fname = source->getCurrentFileName();
if(fname.EndsWith(".root"))
{
fOutputOption = "RECREATE";
fOutputTitle = "Filter";
fOutputCompression = 2;
fname = gSystem->BaseName(fname.Data());
fname.ReplaceAll(".root","");
fOutputFileName = Form("%s/%s%s.root",fOutputDir.Data(),fname.Data(),fOutputFileSuffix.Data());
Info("","CREATING OUTPUT FILE FROM ROOT INPUT = %s !",fOutputFileName.Data());
fOutputFile = new TFile(fOutputFileName.Data(),"RECREATE",fOutputTitle.Data(),fOutputCompression);
fOutFound = kTRUE;
setEvent();
makeTree();
} else {
fOutputOption = "RECREATE";
fOutputTitle = "Filter";
fOutputCompression = 2;
fname = gSystem->BaseName(fname.Data());
fname.ReplaceAll(".hld","");
fOutputFileName = Form("%s/%s%s.root",fOutputDir.Data(),fname.Data(),fOutputFileSuffix.Data());
Info("","CREATING OUTPUT FILE FROM HLD INPUT = %s !",fOutputFileName.Data());
fOutputFile = new TFile(fOutputFileName.Data(),fOutputOption.Data(),fOutputTitle.Data(),fOutputCompression);
fOutFound = kTRUE;
setEvent();
makeTree();
}
}
}
}
if(gHades->getSrcKeeper()){
HSrcKeeper& sk = *gHades->getSrcKeeper();
keeper = new HSrcKeeper(sk);
} else {
Info("execute()","Could not retrieve HSrcKeeper from gHades!");
}
if(gHades->getGeantMedia()){
HGeantMedia& m = *gHades->getGeantMedia();
media = new HGeantMedia(m);
} else {
Info("execute()","Could not retrieve HGeantMedia from gHades!");
}
const TObjArray* objects = gHades->getObjectsAddedToOutput();
HMdcLayer* ml = (HMdcLayer*)objects ->FindObject("HMdcLayer");
if(ml){
mdclayer = new HMdcLayer(*ml);
} else {
Info("execute()","Could not retrieve HMdcLayer from gHades!");
}
}
if (fOutputFile) {
if (fOutputFile->GetBytesWritten() > gHades->getOutputSizeLimit()) {
recreateOutput();
}
}
if(fTree)
{
Bool_t goodEvent = kTRUE;
fCurrentEvent->Clear();
HEventHeader* header = fCurrentEvent->getHeader();
HEventHeader* hadesheader = gHades->getCurrentEvent()->getHeader();
new (header) HEventHeader(*hadesheader);
if(fCurrentEvent->getCategory(catParticleCand))
{
if(kDoSorter)
{
sorter.backupSetup();
sorter.backupFlags(kTRUE);
sorter.getSetup().copyToStatic();
sorter.cleanUp();
sorter.resetFlags(kTRUE,kTRUE,kTRUE,kTRUE);
if(pUserSelectLeptons){
sorter.fill(pUserSelectLeptons);
} else {
sorter.fill(HParticleTrackSorter::selectLeptons);
}
sorter.selectBest(sortType,Particle::kIsLeptonSorter);
if(pUserSelectHadrons){
sorter.fill(pUserSelectHadrons);
} else {
sorter.fill(HParticleTrackSorter::selectHadrons);
}
sorter.selectBest(sortType,Particle::kIsHadronSorter);
}
goodEvent = kTRUE;
if(pUserSelectEvent){
goodEvent = (*pUserSelectEvent)(fParamSelectEvent);
}
if(goodEvent)
{
Bool_t isSim = fmCatNumToName[catParticleCand] == "HParticleCandSim" ? kTRUE : kFALSE;
if(fCurrentEvent->getCategory(catGeantKine))
{
HCategory* catIn = HCategoryManager::getCategory(catGeantKine);
HCategory* catOut = fCurrentEvent->getCategory(catGeantKine);
Int_t n = catIn->getEntries();
HGeantKine* kine1 = 0;
HGeantKine* kine2 = 0;
Int_t index;
for(Int_t i = 0; i < n; i++){
kine1 = HCategoryManager::getObject (kine1,catIn,i);
kine2 = HCategoryManager::newObjectCopy(kine2,kine1,catOut,index);
}
}
if(fCurrentEvent->getCategory(catMdcGeantRaw))
{
HCategory* catIn = HCategoryManager::getCategory(catMdcGeantRaw);
HLocation loc;
loc.set(4,0,0,0,0);
HMatrixCategory* catOut = (HMatrixCategory*)fCurrentEvent->getCategory(catMdcGeantRaw);
Int_t n = catIn->getEntries();
HGeantMdc* hit1 = 0;
HGeantMdc* hit2 = 0;
Int_t ind[6][4][7] ;
for(Int_t s=0;s<6;s++){
for(Int_t m=0;m<4;m++){
for(Int_t l=0;l<7;l++){
ind[s][m][l]=0;
}
}
}
Int_t index;
for(Int_t i = 0; i < n; i++){
hit1 = HCategoryManager::getObject(hit1,catIn,i);
if(hit1){
loc[0] = hit1->getSector();
loc[1] = hit1->getModule();
loc[2] = hit1->getLayer();
loc[3] = ind[loc[0]][loc[1]][loc[2]];
hit2 =(HGeantMdc*) catOut->getSlot(loc,&index);
new (hit2) HGeantMdc(*hit1);
hit2->setNextHitIndex(hit1->getNextHitIndex());
ind[loc[0]][loc[1]][loc[2]]++;
}
}
}
if(fCurrentEvent->getCategory(catTofGeantRaw))
{
HCategory* catIn = HCategoryManager::getCategory(catTofGeantRaw);
HLocation loc;
loc.set(2,0,0);
HMatrixCategory* catOut = (HMatrixCategory*)fCurrentEvent->getCategory(catTofGeantRaw);
Int_t n = catIn->getEntries();
HGeantTof* hit1 = 0;
HGeantTof* hit2 = 0;
Int_t ind[6] ;
for(Int_t s=0;s<6;s++){
ind[s]=0;
}
Int_t index;
for(Int_t i = 0; i < n; i++){
hit1 = HCategoryManager::getObject(hit1,catIn,i);
if(hit1){
loc[0] = hit1->getSector();
loc[1] = ind[loc[0]];
hit2 =(HGeantTof*) catOut->getSlot(loc,&index);
new (hit2) HGeantTof(*hit1);
hit2->setNextHitIndex(hit1->getNextHitIndex());
ind[loc[0]]++;
}
}
}
if(fCurrentEvent->getCategory(catRpcGeantRaw))
{
HCategory* catIn = HCategoryManager::getCategory(catRpcGeantRaw);
HLocation loc;
loc.set(2,0,0);
HMatrixCategory* catOut = (HMatrixCategory*)fCurrentEvent->getCategory(catRpcGeantRaw);
Int_t n = catIn->getEntries();
HGeantRpc* hit1 = 0;
HGeantRpc* hit2 = 0;
Int_t ind[6] ;
for(Int_t s=0;s<6;s++){
ind[s]=0;
}
Int_t index;
for(Int_t i = 0; i < n; i++){
hit1 = HCategoryManager::getObject(hit1,catIn,i);
if(hit1){
loc[0] = hit1->getSector();
loc[1] = ind[loc[0]];
hit2 =(HGeantRpc*) catOut->getSlot(loc,&index);
new (hit2) HGeantRpc(*hit1);
hit2->setNextHitIndex(hit1->getNextHitIndex());
ind[loc[0]]++;
}
}
}
if(fCurrentEvent->getCategory(catShowerGeantRaw))
{
HCategory* catIn = HCategoryManager::getCategory(catShowerGeantRaw);
HLocation loc;
loc.set(2,0,0);
HMatrixCategory* catOut = (HMatrixCategory*)fCurrentEvent->getCategory(catShowerGeantRaw);
Int_t n = catIn->getEntries();
HGeantShower* hit1 = 0;
HGeantShower* hit2 = 0;
Int_t ind[6] ;
for(Int_t s=0;s<6;s++){
ind[s]=0;
}
Int_t index;
for(Int_t i = 0; i < n; i++){
hit1 = HCategoryManager::getObject(hit1,catIn,i);
if(hit1){
loc[0] = hit1->getSector();
loc[1] = ind[loc[0]];
hit2 =(HGeantShower*) catOut->getSlot(loc,&index);
new (hit2) HGeantShower(*hit1);
hit2->setNextHitIndex(hit1->getNextHitIndex());
ind[loc[0]]++;
}
}
}
if(fCurrentEvent->getCategory(catEmcGeantRaw))
{
HCategory* catIn = HCategoryManager::getCategory(catEmcGeantRaw);
HLocation loc;
loc.set(2,0,0);
HMatrixCategory* catOut = (HMatrixCategory*)fCurrentEvent->getCategory(catEmcGeantRaw);
Int_t n = catIn->getEntries();
HGeantEmc* hit1 = 0;
HGeantEmc* hit2 = 0;
Int_t ind[6] ;
for(Int_t s=0;s<6;s++){
ind[s]=0;
}
Int_t index;
for(Int_t i = 0; i < n; i++){
hit1 = HCategoryManager::getObject(hit1,catIn,i);
if(hit1){
loc[0] = hit1->getSector();
loc[1] = ind[loc[0]];
hit2 =(HGeantEmc*) catOut->getSlot(loc,&index);
new (hit2) HGeantEmc(*hit1);
hit2->setNextHitIndex(hit1->getNextHitIndex());
ind[loc[0]]++;
}
}
}
if(fCurrentEvent->getCategory(catWallGeantRaw))
{
HCategory* catIn = HCategoryManager::getCategory(catWallGeantRaw);
HCategory* catOut = fCurrentEvent->getCategory(catWallGeantRaw);
Int_t n = catIn->getEntries();
HGeantWall* wall1 = 0;
HGeantWall* wall2 = 0;
Int_t index;
for(Int_t i = 0; i < n; i++){
wall1 = HCategoryManager::getObject (wall1,catIn,i);
wall2 = HCategoryManager::newObjectCopy(wall2,wall1,catOut,index);
wall2->setNextHitIndex(wall1->getNextHitIndex());
}
}
if(fCurrentEvent->getCategory(catStartGeantRaw))
{
HCategory* catIn = HCategoryManager::getCategory(catStartGeantRaw);
HCategory* catOut = fCurrentEvent->getCategory(catStartGeantRaw);
Int_t n = catIn->getEntries();
HGeantStart* start1 = 0;
HGeantStart* start2 = 0;
Int_t index;
for(Int_t i = 0; i < n; i++){
start1 = HCategoryManager::getObject (start1,catIn,i);
start2 = HCategoryManager::newObjectCopy(start2,start1,catOut,index);
start2->setNextHitIndex(start1->getNextHitIndex());
}
}
if(fCurrentEvent->getCategory(catRichGeantRaw))
{
HCategory* catIn = HCategoryManager::getCategory(catRichGeantRaw);
HLocation loc;
loc.set(2,0,0);
HMatrixCategory* catOut = (HMatrixCategory*)fCurrentEvent->getCategory(catRichGeantRaw);
Int_t n = catIn->getEntries();
HGeantRichPhoton* hit1 = 0;
HGeantRichPhoton* hit2 = 0;
Int_t ind[6] ;
for(Int_t s=0;s<6;s++){
ind[s]=0;
}
Int_t index;
for(Int_t i = 0; i < n; i++){
hit1 = HCategoryManager::getObject(hit1,catIn,i);
if(hit1){
loc[0] = hit1->getSector();
loc[1] = ind[loc[0]];
hit2 =(HGeantRichPhoton*) catOut->getSlot(loc,&index);
new (hit2) HGeantRichPhoton(*hit1);
hit2->setNextHitIndex(hit1->getNextHitIndex());
ind[loc[0]]++;
}
}
}
if(fCurrentEvent->getCategory(catRichGeantRaw+1))
{
HCategory* catIn = HCategoryManager::getCategory(catRichGeantRaw+1);
HLocation loc;
loc.set(2,0,0);
HMatrixCategory* catOut = (HMatrixCategory*)fCurrentEvent->getCategory(catRichGeantRaw+1);
Int_t n = catIn->getEntries();
HGeantRichDirect* hit1 = 0;
HGeantRichDirect* hit2 = 0;
Int_t ind[6] ;
for(Int_t s=0;s<6;s++){
ind[s]=0;
}
Int_t index;
for(Int_t i = 0; i < n; i++){
hit1 = HCategoryManager::getObject(hit1,catIn,i);
if(hit1){
loc[0] = hit1->getSector();
loc[1] = ind[loc[0]];
hit2 =(HGeantRichDirect*) catOut->getSlot(loc,&index);
new (hit2) HGeantRichDirect(*hit1);
hit2->setNextHitIndex(hit1->getNextHitIndex());
ind[loc[0]]++;
}
}
}
if(fCurrentEvent->getCategory(catRichGeantRaw+2))
{
HCategory* catIn = HCategoryManager::getCategory(catRichGeantRaw+2);
HLocation loc;
loc.set(2,0,0);
HMatrixCategory* catOut = (HMatrixCategory*)fCurrentEvent->getCategory(catRichGeantRaw+2);
Int_t n = catIn->getEntries();
HGeantRichMirror* hit1 = 0;
HGeantRichMirror* hit2 = 0;
Int_t ind[6] ;
for(Int_t s=0;s<6;s++){
ind[s]=0;
}
Int_t index;
for(Int_t i = 0; i < n; i++){
hit1 = HCategoryManager::getObject(hit1,catIn,i);
if(hit1){
loc[0] = hit1->getSector();
loc[1] = ind[loc[0]];
hit2 =(HGeantRichMirror*) catOut->getSlot(loc,&index);
new (hit2) HGeantRichMirror(*hit1);
hit2->setNextHitIndex(hit1->getNextHitIndex());
ind[loc[0]]++;
}
}
}
if(fCurrentEvent->getCategory(catParticleEvtInfo))
{
HCategory* catIn = HCategoryManager::getCategory(catParticleEvtInfo);
HCategory* catOut = fCurrentEvent->getCategory(catParticleEvtInfo);
Int_t n = catIn->getEntries();
HParticleEvtInfo* info1 = 0;
HParticleEvtInfo* info2 = 0;
Int_t index;
for(Int_t i = 0; i < n; i++){
info1 = HCategoryManager::getObject (info1,catIn,i);
info2 = HCategoryManager::newObjectCopy(info2,info1,catOut,index);
}
}
if(fCurrentEvent->getCategory(catStart2Hit))
{
HCategory* catIn = HCategoryManager::getCategory(catStart2Hit);
HCategory* catOut = fCurrentEvent->getCategory(catStart2Hit);
Int_t n = catIn->getEntries();
HStart2Hit* hit1 = 0;
HStart2Hit* hit2 = 0;
Int_t index;
for(Int_t i = 0; i < n; i++){
hit1 = HCategoryManager::getObject (hit1,catIn,i);
hit2 = HCategoryManager::newObjectCopy(hit2,hit1,catOut,index);
}
}
if(fCurrentEvent->getCategory(catStart2Cal))
{
HCategory* catIn = HCategoryManager::getCategory(catStart2Cal);
HLocation loc;
loc.set(2,0,0);
HMatrixCategory* catOut = (HMatrixCategory*)fCurrentEvent->getCategory(catStart2Cal);
Int_t n = catIn->getEntries();
HStart2Cal* hit1 = 0;
HStart2Cal* hit2 = 0;
Int_t index;
for(Int_t i = 0; i < n; i++){
hit1 = HCategoryManager::getObject(hit1,catIn,i);
loc[0] = hit1->getModule();
loc[1] = hit1->getStrip();
hit2 =(HStart2Cal*) catOut->getSlot(loc,&index);
new (hit2) HStart2Cal(*hit1);
}
}
if(fCurrentEvent->getCategory(catTBoxChan))
{
HCategory* catIn = HCategoryManager::getCategory(catTBoxChan);
HCategory* catOut = fCurrentEvent->getCategory(catTBoxChan);
Int_t n = catIn->getEntries();
HTBoxChan* hit1 = 0;
HTBoxChan* hit2 = 0;
Int_t index;
for(Int_t i = 0; i < n; i++){
hit1 = HCategoryManager::getObject (hit1,catIn,i);
hit2 = HCategoryManager::newObjectCopy(hit2,hit1,catOut,index);
}
}
if(fCurrentEvent->getCategory(catWallHit))
{
HCategory* catIn = HCategoryManager::getCategory(catWallHit);
HMatrixCategory* catOut = (HMatrixCategory*)fCurrentEvent->getCategory(catWallHit);
HLocation loc;
loc.set(1,0);
Int_t n = catIn->getEntries();
Int_t index;
for(Int_t i = 0; i < n; i++){
if(isSim){
HWallHitSim* hit1 = 0;
HWallHitSim* hit2 = 0;
hit1 = HCategoryManager::getObject (hit1,catIn,i);
loc[0] = hit1->getCell();
hit2 =(HWallHitSim*) catOut->getSlot(loc,&index);
new (hit2) HWallHitSim(*hit1);
} else {
HWallHit* hit1 = 0;
HWallHit* hit2 = 0;
hit1 = HCategoryManager::getObject (hit1,catIn,i);
loc[0] = hit1->getCell();
hit2 =(HWallHit*) catOut->getSlot(loc,&index);
new (hit2) HWallHit(*hit1);
}
}
}
if(fCurrentEvent->getCategory(catWallEventPlane))
{
HCategory* catIn = HCategoryManager::getCategory(catWallEventPlane);
HCategory* catOut = fCurrentEvent->getCategory(catWallEventPlane);
Int_t n = catIn->getEntries();
HWallEventPlane* info1 = 0;
HWallEventPlane* info2 = 0;
Int_t index;
for(Int_t i = 0; i < n; i++){
info1 = HCategoryManager::getObject (info1,catIn,i);
info2 = HCategoryManager::newObjectCopy(info2,info1,catOut,index);
}
}
if(fCurrentEvent->getCategory(catRichCal))
{
HCategory* catIn = HCategoryManager::getCategory(catRichCal);
HMatrixCategory* catOut = (HMatrixCategory*)fCurrentEvent->getCategory(catRichCal);
HLocation loc;
loc.set(3,0,0,0);
Int_t n = catIn->getEntries();
Int_t index;
for(Int_t i = 0; i < n; i++){
if(isSim){
HRichCalSim* hit1 = 0;
HRichCalSim* hit2 = 0;
hit1 = HCategoryManager::getObject(hit1,catIn,i);
loc[0] = hit1->getSector();
loc[1] = hit1->getRow();
loc[2] = hit1->getCol();
hit2 =(HRichCalSim*) catOut->getSlot(loc,&index);
new (hit2) HRichCalSim(*hit1);
} else {
HRichCal* hit1 = 0;
HRichCal* hit2 = 0;
hit1 = HCategoryManager::getObject(hit1,catIn,i);
loc[0] = hit1->getSector();
loc[1] = hit1->getRow();
loc[2] = hit1->getCol();
hit2 =(HRichCal*) catOut->getSlot(loc,&index);
new (hit2) HRichCal(*hit1);
}
}
}
if(fCurrentEvent->getCategory(catRichDirClus))
{
HCategory* catIn = HCategoryManager::getCategory(catRichDirClus);
HMatrixCategory* catOut = (HMatrixCategory*)fCurrentEvent->getCategory(catRichDirClus);
HLocation loc;
loc.set(1,0);
Int_t n = catIn->getEntries();
HRichDirClus* hit1 = 0;
HRichDirClus* hit2 = 0;
Int_t index;
for(Int_t i = 0; i < n; i++){
hit1 = HCategoryManager::getObject (hit1,catIn,i);
hit2 =(HRichDirClus*) catOut->getNewSlot(loc,&index);
new (hit2) HRichDirClus(*hit1);
}
}
if(fCurrentEvent->getCategory(catPionTrackerTrack))
{
HCategory* catIn = HCategoryManager::getCategory(catPionTrackerTrack);
HCategory* catOut = fCurrentEvent->getCategory(catPionTrackerTrack);
Int_t n = catIn->getEntries();
HPionTrackerTrack* hit1 = 0;
HPionTrackerTrack* hit2 = 0;
Int_t index;
for(Int_t i = 0; i < n; i++){
hit1 = HCategoryManager::getObject (hit1,catIn,i);
hit2 = HCategoryManager::newObjectCopy(hit2,hit1,catOut,index);
}
}
if(fCurrentEvent->getCategory(catPionTrackerCal))
{
HCategory* catIn = HCategoryManager::getCategory(catPionTrackerCal);
HMatrixCategory* catOut = (HMatrixCategory*)fCurrentEvent->getCategory(catPionTrackerCal);
HLocation loc;
loc.set(2,0,0);
Int_t n = catIn->getEntries();
HPionTrackerCal* hit1 = 0;
HPionTrackerCal* hit2 = 0;
Int_t index;
for(Int_t i = 0; i < n; i++){
hit1 = HCategoryManager::getObject (hit1,catIn,i);
hit2 =(HPionTrackerCal*) catOut->getNewSlot(loc,&index);
new (hit2) HPionTrackerCal(*hit1);
}
}
if(fCurrentEvent->getCategory(catPionTrackerHit))
{
HCategory* catIn = HCategoryManager::getCategory(catPionTrackerHit);
HMatrixCategory* catOut = (HMatrixCategory*)fCurrentEvent->getCategory(catPionTrackerHit);
HLocation loc;
loc.set(2,0,0);
Int_t n = catIn->getEntries();
HPionTrackerHit* hit1 = 0;
HPionTrackerHit* hit2 = 0;
Int_t index;
for(Int_t i = 0; i < n; i++){
hit1 = HCategoryManager::getObject (hit1,catIn,i);
hit2 =(HPionTrackerHit*) catOut->getNewSlot(loc,&index);
new (hit2) HPionTrackerHit(*hit1);
}
}
if(fCurrentEvent->getCategory(catEmcCluster))
{
HCategory* catIn = HCategoryManager::getCategory(catEmcCluster);
HLocation loc;
loc.set(2,0,0);
HMatrixCategory* catOut = (HMatrixCategory*)fCurrentEvent->getCategory(catEmcCluster);
Int_t n = catIn->getEntries();
HEmcCluster* hit1 = 0;
HEmcCluster* hit2 = 0;
Int_t index;
for(Int_t i = 0; i < n; i++){
hit1 = HCategoryManager::getObject(hit1,catIn,i);
loc[0] = hit1->getSector();
loc[1] = hit1->getCell();
hit2 =(HEmcCluster*) catOut->getSlot(loc,&index);
new (hit2) HEmcCluster(*hit1);
}
}
if(fCurrentEvent->getCategory(catEmcCal))
{
HCategory* catIn = HCategoryManager::getCategory(catEmcCal);
HLocation loc;
loc.set(2,0,0);
HMatrixCategory* catOut = (HMatrixCategory*)fCurrentEvent->getCategory(catEmcCal);
Int_t n = catIn->getEntries();
HEmcCal* hit1 = 0;
HEmcCal* hit2 = 0;
Int_t index;
for(Int_t i = 0; i < n; i++){
hit1 = HCategoryManager::getObject(hit1,catIn,i);
loc[0] = hit1->getSector();
loc[1] = hit1->getCell();
hit2 =(HEmcCal*) catOut->getSlot(loc,&index);
new (hit2) HEmcCal(*hit1);
}
}
if(fCurrentEvent->getCategory(catRichHit) && doFullCopy(catRichHit)){
HLinearCategory* catOut = (HLinearCategory*)fCurrentEvent->getCategory(catRichHit);
HCategory* catIn = HCategoryManager::getCategory(catRichHit);
Int_t n = catIn->getEntries();
Int_t index;
for(Int_t i = 0; i < n; i++){
if(isSim)
{
HRichHitSim* hit1 = 0;
HRichHitSim* hit2 = 0;
hit1 = HCategoryManager::getObject(hit1,catIn,i);
hit2 = HCategoryManager::newObjectCopy(hit2,hit1,catOut,index);
} else {
HRichHit* hit1 = 0;
HRichHit* hit2 = 0;
hit1 = HCategoryManager::getObject(hit1,catIn,i);
hit2 = HCategoryManager::newObjectCopy(hit2,hit1,catOut,index);
}
}
}
if(fCurrentEvent->getCategory(catTofHit) && doFullCopy(catTofHit)){
HMatrixCategory* catOut = (HMatrixCategory*)fCurrentEvent->getCategory(catTofHit);
HCategory* catIn = HCategoryManager::getCategory(catTofHit);
HLocation loc;
loc.set(3,0,0,0);
Int_t n = catIn->getEntries();
Int_t index;
for(Int_t i = 0; i < n; i++){
if(isSim)
{
HTofHitSim* hit1 = 0;
HTofHitSim* hit2 = 0;
hit1 = HCategoryManager::getObject(hit1,catTofHit,i);
loc[0] = hit1->getSector();
loc[1] = hit1->getModule();
loc[2] = hit1->getCell();
hit2 =(HTofHitSim*) catOut->getSlot(loc,&index);
new (hit2) HTofHitSim(*hit1);
} else {
HTofHit* hit1 = 0;
HTofHit* hit2 = 0;
hit1 = HCategoryManager::getObject(hit1,catTofHit,i);
loc[0] = hit1->getSector();
loc[1] = hit1->getModule();
loc[2] = hit1->getCell();
hit2 =(HTofHit*) catOut->getSlot(loc,&index);
new (hit2) HTofHit(*hit1);
}
}
}
if(fCurrentEvent->getCategory(catTofCluster) && doFullCopy(catTofCluster)){
HMatrixCategory* catOut = (HMatrixCategory*)fCurrentEvent->getCategory(catTofCluster);
HCategory* catIn = HCategoryManager::getCategory(catTofCluster);
HLocation loc;
loc.set(3,0,0,0);
Int_t n = catIn->getEntries();
Int_t index;
for(Int_t i = 0; i < n; i++){
if(isSim)
{
HTofClusterSim* hit1 = 0;
HTofClusterSim* hit2 = 0;
hit1 = HCategoryManager::getObject(hit1,catTofCluster,i);
loc[0] = hit1->getSector();
loc[1] = hit1->getModule();
loc[2] = hit1->getCell();
hit2 =(HTofClusterSim*) catOut->getSlot(loc,&index);
new (hit2) HTofClusterSim(*hit1);
} else {
HTofCluster* hit1 = 0;
HTofCluster* hit2 = 0;
hit1 = HCategoryManager::getObject(hit1,catTofCluster,i);
loc[0] = hit1->getSector();
loc[1] = hit1->getModule();
loc[2] = hit1->getCell();
hit2 =(HTofCluster*) catOut->getSlot(loc,&index);
new (hit2) HTofCluster(*hit1);
}
}
}
if(fCurrentEvent->getCategory(catRpcCluster) && doFullCopy(catRpcCluster)){
HLocation loc;
loc.set(2,0,0);
HMatrixCategory* catOut = (HMatrixCategory*)fCurrentEvent->getCategory(catRpcCluster);
HCategory* catIn = HCategoryManager::getCategory(catRpcCluster);
Int_t n = catIn->getEntries();
Int_t index;
for(Int_t i = 0; i < n; i++){
if(isSim)
{
HRpcClusterSim* hit1 = 0;
HRpcClusterSim* hit2 = 0;
hit1 = HCategoryManager::getObject(hit1,catRpcCluster,i);
loc[0] = hit1->getSector();
loc[1] = ctRpcClus[hit1->getSector()];
hit2 =(HRpcClusterSim*) catOut->getSlot(loc,&index);
new (hit2) HRpcClusterSim(*hit1);
hit2->setAddress(loc[0],index);
ctRpcClus[hit1->getSector()] ++;
} else {
HRpcCluster* hit1 = 0;
HRpcCluster* hit2 = 0;
hit1 = HCategoryManager::getObject(hit1,catRpcCluster,i);
loc[0] = hit1->getSector();
loc[1] = ctRpcClus[hit1->getSector()];
hit2 =(HRpcCluster*) catOut->getSlot(loc,&index);
new (hit2) HRpcCluster(*hit1);
hit2->setAddress(loc[0],index);
ctRpcClus[hit1->getSector()] ++;
}
}
}
if(fCurrentEvent->getCategory(catShowerHit) && doFullCopy(catShowerHit)){
HCategory* catOut = fCurrentEvent->getCategory(catShowerHit);
HCategory* catIn = HCategoryManager::getCategory(catShowerHit);
Int_t n = catIn->getEntries();
Int_t index;
for(Int_t i = 0; i < n; i++){
if(isSim)
{
HShowerHitSim* hit1 = 0;
HShowerHitSim* hit2 = 0;
hit1 = HCategoryManager::getObject(hit1,catShowerHit,i);
hit2 = HCategoryManager::newObjectCopy(hit2,hit1,catOut,index);
} else {
HShowerHit* hit1 = 0;
HShowerHit* hit2 = 0;
hit1 = HCategoryManager::getObject(hit1,catShowerHit,i);
hit2 = HCategoryManager::newObjectCopy(hit2,hit1,catOut,index);
}
}
}
HCategory* catIn = HCategoryManager::getCategory(catParticleCand);
HCategory* catOut = fCurrentEvent->getCategory(catParticleCand);
Int_t n = catIn->getEntries();
HParticleCand *cand, *cand1, *cand2; cand=cand1=cand2=0;
HParticleCandSim* cands1, *cands2; cands1=cands2=0;
Int_t indexCand2;
Int_t particleMdcInd1;
Int_t richInd1;
Int_t seg0Ind1;
Int_t seg1Ind1;
Int_t rpcInd1;
Int_t showerInd1;
Int_t tofInd1;
Int_t tofClstInd1;
Int_t particleMdcInd2;
Int_t richInd2;
Int_t seg0Ind2;
Int_t seg1Ind2;
Int_t rpcInd2;
Int_t showerInd2;
Int_t tofInd2;
Int_t tofClstInd2;
map<Int_t,Int_t>mPartMdc;
map<Int_t,Int_t>mSeg;
map<Int_t,Int_t>mRich;
map<Int_t,Int_t>mTofHit;
map<Int_t,Int_t>mTofClst;
map<Int_t,Int_t>mRpcClst;
map<Int_t,Int_t>mShrHit;
for(Int_t i = 0; i < n; i++){
if(isSim){
cands1 = HCategoryManager::getObject (cands1,catIn,i);
if(kSkipTracks ){
if(pUserKeepTrack){
if(!(*pUserKeepTrack)(cands1)) continue;
} else {
if(!cands1->isFlagBit(Particle::kIsUsed)) continue;
}
}
cands2 = HCategoryManager::newObjectCopy(cands2,cands1,catOut,indexCand2);
cands2->setIndex(indexCand2);
cand = cands2;
particleMdcInd1 = cands1->getIndex();
richInd1 = cands1->getRichInd();
seg0Ind1 = cands1->getInnerSegInd();
seg1Ind1 = cands1->getOuterSegInd();
rpcInd1 = cands1->getRpcInd();
showerInd1 = cands1->getShowerInd();
tofInd1 = cands1->getTofHitInd();
tofClstInd1 = cands1->getTofClstInd();
} else {
cand1 = HCategoryManager::getObject (cand1,catIn,i);
if(kSkipTracks){
if(pUserKeepTrack){
if(!(*pUserKeepTrack)(cand1)) continue;
} else {
if(!cand1->isFlagBit(Particle::kIsUsed)) continue;
}
}
cand2 = HCategoryManager::newObjectCopy(cand2,cand1,catOut,indexCand2);
cand2->setIndex(indexCand2);
cand = cand2;
particleMdcInd1 = cand1->getIndex();
richInd1 = cand1->getRichInd();
seg0Ind1 = cand1->getInnerSegInd();
seg1Ind1 = cand1->getOuterSegInd();
rpcInd1 = cand1->getRpcInd();
showerInd1 = cand1->getShowerInd();
tofInd1 = cand1->getTofHitInd();
tofClstInd1 = cand1->getTofClstInd();
}
particleMdcInd2 = -1;
richInd2 = -1;
seg0Ind2 = -1;
seg1Ind2 = -1;
rpcInd2 = -1;
showerInd2 = -1;
tofInd2 = -1;
tofClstInd2 = -1;
if(catParticleMdc)
{
if(particleMdcInd1 !=-1 && fCurrentEvent->getCategory(catParticleMdc))
{
if(mPartMdc.find(particleMdcInd1) == mPartMdc.end())
{
HLinearCategory* catOut = (HLinearCategory*)fCurrentEvent->getCategory(catParticleMdc);
HParticleMdc* particleMdc1 = 0;
HParticleMdc* particleMdc2 = 0;
particleMdc1 = HCategoryManager::getObject(particleMdc1,catParticleMdc,particleMdcInd1);
particleMdc2 = HCategoryManager::newObjectCopy(particleMdc2,particleMdc1,catOut,particleMdcInd2);
particleMdcInd2 = indexCand2;
particleMdc2->setIndex(particleMdcInd2);
mPartMdc[particleMdcInd1] = particleMdcInd2;
}
}
}
if(!doFullCopy(catRichHit))
{
if(richInd1 !=-1 && fCurrentEvent->getCategory(catRichHit))
{
if(mRich.find(richInd1) == mRich.end())
{
HLinearCategory* catOut = (HLinearCategory*)fCurrentEvent->getCategory(catRichHit);
if(isSim)
{
HRichHitSim* hit1 = 0;
HRichHitSim* hit2 = 0;
hit1 = HCategoryManager::getObject(hit1,catRichHit,richInd1);
hit2 = HCategoryManager::newObjectCopy(hit2,hit1,catOut,richInd2);
} else {
HRichHit* hit1 = 0;
HRichHit* hit2 = 0;
hit1 = HCategoryManager::getObject(hit1,catRichHit,richInd1);
hit2 = HCategoryManager::newObjectCopy(hit2,hit1,catOut,richInd2);
}
mRich[richInd1] = richInd2;
} else {
richInd2 = mRich[richInd1];
}
cand->setRichInd(richInd2);
}
} else { cand->setRichInd(richInd1); }
if(!doFullCopy(catTofHit))
{
if(tofInd1 !=-1 && fCurrentEvent->getCategory(catTofHit))
{
if(mTofHit.find(tofInd1) == mTofHit.end())
{
HLocation loc;
loc.set(3,0,0,0);
HMatrixCategory* catOut = (HMatrixCategory*)fCurrentEvent->getCategory(catTofHit);
if(isSim)
{
HTofHitSim* hit1 = 0;
HTofHitSim* hit2 = 0;
hit1 = HCategoryManager::getObject(hit1,catTofHit,tofInd1);
loc[0] = hit1->getSector();
loc[1] = hit1->getModule();
loc[2] = hit1->getCell();
hit2 =(HTofHitSim*) catOut->getSlot(loc,&tofInd2);
new (hit2) HTofHitSim(*hit1);
} else {
HTofHit* hit1 = 0;
HTofHit* hit2 = 0;
hit1 = HCategoryManager::getObject(hit1,catTofHit,tofInd1);
loc[0] = hit1->getSector();
loc[1] = hit1->getModule();
loc[2] = hit1->getCell();
hit2 =(HTofHit*) catOut->getSlot(loc,&tofInd2);
new (hit2) HTofHit(*hit1);
}
mTofHit[tofInd1] = tofInd2;
} else {
tofInd2 = mTofHit[tofInd1];
}
cand->setTofHitInd(tofInd2);
}
} else { cand->setTofHitInd(tofInd1); }
if(!doFullCopy(catTofCluster))
{
if(tofClstInd1 !=-1 && fCurrentEvent->getCategory(catTofCluster))
{
if(mTofClst.find(tofClstInd1) == mTofClst.end())
{
HLocation loc;
loc.set(3,0,0,0);
HMatrixCategory* catOut = (HMatrixCategory*)fCurrentEvent->getCategory(catTofCluster);
if(isSim)
{
HTofClusterSim* hit1 = 0;
HTofClusterSim* hit2 = 0;
hit1 = HCategoryManager::getObject(hit1,catTofCluster,tofClstInd1);
loc[0] = hit1->getSector();
loc[1] = hit1->getModule();
loc[2] = hit1->getCell();
hit2 =(HTofClusterSim*) catOut->getSlot(loc,&tofClstInd2);
new (hit2) HTofClusterSim(*hit1);
} else {
HTofCluster* hit1 = 0;
HTofCluster* hit2 = 0;
hit1 = HCategoryManager::getObject(hit1,catTofCluster,tofClstInd1);
loc[0] = hit1->getSector();
loc[1] = hit1->getModule();
loc[2] = hit1->getCell();
hit2 =(HTofCluster*) catOut->getSlot(loc,&tofClstInd2);
new (hit2) HTofCluster(*hit1);
}
mTofClst[tofClstInd1] = tofClstInd2;
} else {
tofClstInd2 = mTofClst[tofClstInd1];
}
cand->setTofClstInd(tofClstInd2);
}
} else { cand->setTofClstInd(tofClstInd1); }
if(!doFullCopy(catRpcCluster))
{
if(rpcInd1 !=-1 && fCurrentEvent->getCategory(catRpcCluster))
{
if(mRpcClst.find(rpcInd1) == mRpcClst.end())
{
HLocation loc;
loc.set(2,0,0);
HMatrixCategory* catOut = (HMatrixCategory*)fCurrentEvent->getCategory(catRpcCluster);
if(isSim)
{
HRpcClusterSim* hit1 = 0;
HRpcClusterSim* hit2 = 0;
hit1 = HCategoryManager::getObject(hit1,catRpcCluster,rpcInd1);
loc[0] = hit1->getSector();
loc[1] = ctRpcClus[hit1->getSector()];
hit2 =(HRpcClusterSim*) catOut->getSlot(loc,&rpcInd2);
new (hit2) HRpcClusterSim(*hit1);
hit2->setAddress(loc[0],rpcInd2);
ctRpcClus[hit1->getSector()] ++;
} else {
HRpcCluster* hit1 = 0;
HRpcCluster* hit2 = 0;
hit1 = HCategoryManager::getObject(hit1,catRpcCluster,rpcInd1);
loc[0] = hit1->getSector();
loc[1] = ctRpcClus[hit1->getSector()];
hit2 =(HRpcCluster*) catOut->getSlot(loc,&rpcInd2);
new (hit2) HRpcCluster(*hit1);
hit2->setAddress(loc[0],rpcInd2);
ctRpcClus[hit1->getSector()] ++;
}
mRpcClst[rpcInd1] = rpcInd2;
} else {
rpcInd2 = mRpcClst[rpcInd1];
}
cand->setRpcInd(rpcInd2);
}
} else { cand->setRpcInd(rpcInd1); }
if(!doFullCopy(catShowerHit))
{
if(showerInd1 !=-1 && fCurrentEvent->getCategory(catShowerHit))
{
if(mShrHit.find(showerInd1) == mShrHit.end())
{
HCategory* catOut = fCurrentEvent->getCategory(catShowerHit);
if(isSim)
{
HShowerHitSim* hit1 = 0;
HShowerHitSim* hit2 = 0;
hit1 = HCategoryManager::getObject(hit1,catShowerHit,showerInd1);
hit2 = HCategoryManager::newObjectCopy(hit2,hit1,catOut,showerInd2);
} else {
HShowerHit* hit1 = 0;
HShowerHit* hit2 = 0;
hit1 = HCategoryManager::getObject(hit1,catShowerHit,showerInd1);
hit2 = HCategoryManager::newObjectCopy(hit2,hit1,catOut,showerInd2);
}
mShrHit[showerInd1] = showerInd2;
} else {
showerInd2 = mShrHit[showerInd1];
}
cand->setShowerInd(showerInd2);
}
} else { cand->setShowerInd(showerInd1); }
if(seg0Ind1 !=-1 && fCurrentEvent->getCategory(catMdcSeg))
{
if(mSeg.find(seg0Ind1) == mSeg.end()){
extractMdcSeg(cand,isSim,seg0Ind1,seg0Ind2);
mSeg[seg0Ind1] = seg0Ind2;
} else {
cand->setInnerSegInd(mSeg[seg0Ind1]);
}
}
if(seg1Ind1 !=-1 && fCurrentEvent->getCategory(catMdcSeg))
{
if(mSeg.find(seg1Ind1) == mSeg.end()){
extractMdcSeg(cand,isSim,seg1Ind1,seg1Ind2);
mSeg[seg1Ind1] = seg1Ind2;
} else {
cand->setOuterSegInd(mSeg[seg1Ind1]);
}
}
}
}
if(kDoSorter){
sorter.restoreFlags(kTRUE);
sorter.restoreSetup();
}
}
if( (goodEvent && kSkipEmptyEvents) || !kSkipEmptyEvents ) fTree->Fill();
}
return 0;
}
Bool_t HParticleTree::finalize()
{
closeOutput();
return kTRUE;
}
void HParticleTree::setEventStructure(Int_t n,Cat_t PersistentCat[],Bool_t fullCopy)
{
for(Int_t i = 0; i < n; i++){
if(find(fmCatNumSupport.begin(),fmCatNumSupport.end(),PersistentCat[i]) != fmCatNumSupport.end()){
if(find(fCatNums.begin(),fCatNums.end(),PersistentCat[i]) == fCatNums.end()){
fCatNums.push_back(PersistentCat[i]);
} else {
Warning("setEventStructure()","Cat %i already in event structure ! Skipped second time.",PersistentCat[i]);
continue;
}
if(fullCopy ){
if(find(fmCatNumFullCopySupport.begin(),fmCatNumFullCopySupport.end(),PersistentCat[i]) != fmCatNumFullCopySupport.end() ){
fmCatNumToFullCopy [PersistentCat[i]] = 1;
} else {
Warning("setEventStructure()","Cat %i not supported for fullCopy ! added in filtered copy mode.",PersistentCat[i]);
fmCatNumToFullCopy [PersistentCat[i]] = 0;
}
} else {
fmCatNumToFullCopy [PersistentCat[i]] = 0;
}
} else {
Warning("setEventStructure()","Cat %i not supported, will be skipped!",PersistentCat[i]);
}
}
sort(fCatNums.begin(),fCatNums.end());
}
void HParticleTree::setEvent()
{
HRecEvent* event = ((HRecEvent*) gHades->getCurrentEvent());
fCurrentEvent = new HRecEvent();
TObjArray* pevts = event->getPartialEvents();
map<TString , Int_t> mpOut;
cout<<"#########################################################"<<endl;
cout<<"From task : name = "<< this->GetName()<<", title = "<<this->GetTitle()<<endl;
cout<<"Partial events found in input = "<<pevts->GetEntries()<<endl;
cout<<"Creating output event structure :"<<endl;
Bool_t doAll = fCatNums.size() == 0 ? kTRUE : kFALSE;
for(Int_t i = 0;i < pevts->GetSize(); i++ ) {
HPartialEvent* p = (HPartialEvent*)pevts->At(i);
if(p){
if(mpOut.find(p->GetName()) == mpOut.end()){
TObjArray* cats = p->getCategories( );
Int_t catnum;
for(Int_t j = 0;j < cats->GetSize(); j++ ) {
HCategory* cat = (HCategory*)cats->At(j);
if(cat) catnum = cat->getCategory( );
if(cat &&
( doAll || find(fCatNums.begin(),fCatNums.end(),catnum) != fCatNums.end() ) &&
mpOut.find(p->GetName()) == mpOut.end())
{
mpOut[p->GetName()] = i;
}
}
}
}
}
map<TString , Int_t> mp;
for(Int_t i = 0;i < pevts->GetSize(); i++ ) {
HPartialEvent* p = (HPartialEvent*)pevts->At(i);
if(p ){
Cat_t catbasenum = p->getBaseCat();
if(mp .find(p->GetName()) == mp .end() &&
mpOut.find(p->GetName()) != mpOut.end()
){
cout<<setw(2)<<i<<" partial "<<setw(15)<<p->GetName()<<" title "<<setw(15)<<p->GetTitle()<<" catbase = "<<setw(5)<<catbasenum<<endl;
fCurrentEvent->addPartialEvent(new HPartialEvent(p->GetName(),p->GetTitle(),p->getBaseCat()));
TObjArray* cats = p->getCategories( );
Int_t catnum = 0;
for(Int_t j = 0;j < cats->GetSize(); j++ ) {
HCategory* cat = (HCategory*) cats->At(j);
if(cat) catnum = cat->getCategory( );
if(cat && (doAll || find(fCatNums.begin(),fCatNums.end(),catnum) != fCatNums.end() ) )
{
TString type = cat->GetName() ;
fmCatNameToNum[cat->getClassName()] = catnum;
fmCatNumToName[catnum] = cat->getClassName();
if (type == "HLinearCategory"){
cout<<"\t"<<setw(2)<<j<<" cat "<<setw(15)<<cat->GetName()<<" "<<setw(20)<<cat->getClassName()<<" catnum "<<setw(5)<<catnum<<endl;
HLinearCategory* catL = ((HLinearCategory*)cat);
const TObjArray* ar = catL->getData();
Int_t size = ar->GetSize();
HLinearCategory* catNew = new HLinearCategory(cat->getClassName(),size);
catNew->setDynamicObjects(catL->getDynamicObjects());
fCurrentEvent->addCategory(cat->getCategory( ),catNew,p->GetName());
fmCatNumToPointer [catnum] = catNew ;
fmCatNameToPointer[cat->getClassName()] = catNew ;
} else if (type == "HMatrixCategory"){
HMatrixCategory* catM = ((HMatrixCategory*)cat);
TArrayI* ar = catM->getSizes();
cout<<"\t"<<setw(2)<<j<<" cat "<<setw(15)<<cat->GetName()<<" "<<setw(20)<<cat->getClassName()<<" catnum "<<setw(5)<<catnum
<<" sizes: "<<ar->GetSize()<<" : "<<flush;
for(Int_t l=0; l < ar->GetSize(); l ++){
cout<<ar->At(l)<<" "<<flush;
}
cout<<endl;
HMatrixCategory* catNew = new HMatrixCategory(cat->getClassName(),ar->GetSize(),ar->GetArray(),1.);
fCurrentEvent->addCategory(cat->getCategory( ),catNew,p->GetName());
fmCatNumToPointer [catnum] = catNew ;
fmCatNameToPointer[cat->getClassName()] = catNew ;
}
}
}
mp[p->GetName()]=i;
}
}
}
cout<<"#########################################################"<<endl;
}
void HParticleTree::setOutputFile (TString fname,TString ftitle,TString fopt,Int_t fcomp)
{
fOutputFileName = fname;
fOutputTitle = ftitle;
fOutputOption = fopt;
fOutputCompression = fcomp;
}
Bool_t HParticleTree::makeTree(void)
{
TBranch *b = NULL;
Bool_t r = kFALSE;
Text_t treeTitle[255];
if (fTree) delete fTree;
sprintf(treeTitle,"T.%i",gHades->getSplitLevel());
if (fOutputFile) { fOutputFile->cd(); }
fTree = new HTree("T",treeTitle);
if (fCurrentEvent && fTree)
{
if (gHades->getSplitLevel() == 0) {
b = fTree->Branch("Event",
fCurrentEvent->ClassName(),
&fCurrentEvent,64000,0);
} else {
b = fTree->Branch("Event.",
fCurrentEvent->ClassName(),
&fCurrentEvent,gHades->getTreeBufferSize(),99);
fCurrentEvent->makeBranch(b,fTree);
}
if (b) {
r = kTRUE;
}
}
return r;
}
void HParticleTree::closeOutput()
{
if (fOutputFile)
{
fOutputFile->cd();
if (fCurrentEvent)
{
if (fCurrentEvent->InheritsFrom("HRecEvent")) {
Bool_t expand = ((HRecEvent *)fCurrentEvent)->hasExpandedStreamer();
( (HRecEvent *)fCurrentEvent)->setExpandedStreamer(kTRUE);
fCurrentEvent->Clear();
fCurrentEvent->Write("Event");
( (HRecEvent *)fCurrentEvent)->setExpandedStreamer(expand);
} else {
fCurrentEvent->Clear();
fCurrentEvent->Write("Event");
}
}
if(keeper) keeper ->Write();
if(media) media ->Write();
if(mdclayer)mdclayer->Write();
Bool_t w = gHades->getWriteEvent();
gHades->setWriteEvent(kFALSE);
gHades->Write();
gHades->setWriteEvent(w);
fOutputFile->Write();
delete fOutputFile;
fOutputFile = NULL;
fTree = 0;
}
}
void HParticleTree::recreateOutput(void)
{
fCycleNumber ++;
Bool_t createTree = (fTree != 0);
closeOutput();
TString fname = fOutputFileName;
fname.ReplaceAll(".root","");
fOutputFileName = Form("%s_%i.root",fname.Data(),fCycleNumber);
fOutputFile = new TFile(fOutputFileName.Data(),fOutputOption.Data(),fOutputTitle.Data(),fOutputCompression);
if(keeper) keeper ->Write();
if(media) media ->Write();
if(mdclayer)mdclayer->Write();
if (createTree) { makeTree(); }
}
void HParticleTree::extractMdcSeg(HParticleCand* cand, Bool_t isSim,Int_t segInd1,Int_t& segInd2)
{
HLocation loc;
loc.set(3,0,0,0);
HMatrixCategory* catOut = (HMatrixCategory*)fCurrentEvent->getCategory(catMdcSeg);
HMdcSeg* seg = 0;
if(isSim)
{
HMdcSegSim* hit1 = 0;
HMdcSegSim* hit2 = 0;
hit1 = HCategoryManager::getObject(hit1,catMdcSeg,segInd1);
loc[0] = hit1->getSec();
loc[1] = hit1->getIOSeg();
loc[2] = ctMdcSeg[hit1->getSec()][hit1->getIOSeg()];
hit2 =(HMdcSegSim*) catOut->getSlot(loc,&segInd2);
new (hit2) HMdcSegSim(*hit1);
seg = hit2;
} else {
HMdcSeg* hit1 = 0;
HMdcSeg* hit2 = 0;
hit1 = HCategoryManager::getObject(hit1,catMdcSeg,segInd1);
loc[0] = hit1->getSec();
loc[1] = hit1->getIOSeg();
loc[2] = ctMdcSeg[hit1->getSec()][hit1->getIOSeg()];
hit2 =(HMdcSeg*) catOut->getSlot(loc,&segInd2);
new (hit2) HMdcSeg(*hit1);
seg = hit2;
}
if(seg->getIOSeg()==0) cand->setInnerSegInd(segInd2);
else cand->setOuterSegInd(segInd2);
ctMdcSeg[seg->getSec()][seg->getIOSeg()]++;
extractMdcCal1 (isSim,segInd1);
extractMdcCal1FromClus(isSim,segInd1);
HMdcClus* clus = HParticleTool::getMdcClus(segInd1);
HMdcHit* mdchit0 = HParticleTool::getMdcHit(segInd1,0);
HMdcHit* mdchit1 = HParticleTool::getMdcHit(segInd1,1);
HMdcClusInf* clusinf0 = HParticleTool::getMdcClusInf(segInd1,0);
HMdcClusInf* clusinf1 = HParticleTool::getMdcClusInf(segInd1,1);
HMdcClusFit* clusfit = HParticleTool::getMdcClusFit(segInd1);
HCategory* wireFitCat = fCurrentEvent->getCategory(catMdcWireFit);
HCategory* clusFitCat = fCurrentEvent->getCategory(catMdcClusFit);
HCategory* clusInfCat = fCurrentEvent->getCategory(catMdcClusInf);
HMatrixCategory* hitCat = (HMatrixCategory*) fCurrentEvent->getCategory(catMdcHit);
HMatrixCategory* clusCat = (HMatrixCategory*) fCurrentEvent->getCategory(catMdcClus);
HMdcClus* clus2 = 0;
HMdcClusFit* clusfit2 = 0;
HMdcClusInf* clusinf0_2 = 0;
HMdcClusInf* clusinf1_2 = 0;
HMdcHit* mdchit0_2 = 0;
HMdcHit* mdchit1_2 = 0;
Int_t clusInd2 = -1;
Int_t clusfitInd2 = -1;
Int_t clusinf0Ind2 = -1;
Int_t clusinf1Ind2 = -1;
Int_t mdchit0Ind2 = -1;
Int_t mdchit1Ind2 = -1;
if(clus && clusCat){
HLocation loc;
loc.set(3,0,0,0);
loc[0] = clus->getSec();
loc[1] = clus->getIOSeg();
loc[2] = ctMdcClus[clus->getSec()][clus->getIOSeg()];
clus2 =(HMdcClusSim*) clusCat->getSlot(loc,&clusInd2);
if(isSim){
new (clus2) HMdcClusSim(*((HMdcClusSim*)clus));
} else {
new (clus2) HMdcClus(*clus);
}
clus2->setOwnIndex(clusInd2);
clus2->setSegIndex(segInd2);
ctMdcClus[clus->getSec()][clus->getIOSeg()]++;
}
if(clusfit && clusFitCat){
if(isSim){
HMdcClusFitSim* clusfits2 = 0;
clusfits2 = HCategoryManager::newObjectCopy(clusfits2,(HMdcClusFitSim*)clusfit,clusFitCat,clusfitInd2);
clusfit2 = clusfits2;
} else {
clusfit2 = HCategoryManager::newObjectCopy(clusfit2,clusfit,clusFitCat,clusfitInd2);
}
clusfit2->setClustCatIndex(clusInd2);
}
if(clusinf0 && clusInfCat){
if(isSim){
HMdcClusInfSim* clusinfs2 = 0;
clusinfs2 = HCategoryManager::newObjectCopy(clusinfs2,(HMdcClusInfSim*)clusinf0,clusInfCat,clusinf0Ind2);
clusinf0_2 = clusinfs2;
} else {
clusinf0_2 = HCategoryManager::newObjectCopy(clusinf0_2,clusinf0,clusInfCat,clusinf0Ind2);
}
clusinf0_2->setClusIndex(clusInd2);
clusinf0_2->setClusFitIndex(clusfitInd2);
}
if(clusinf1 && clusInfCat){
if(isSim){
HMdcClusInfSim* clusinfs2 = 0;
clusinfs2 = HCategoryManager::newObjectCopy(clusinfs2,(HMdcClusInfSim*)clusinf1,clusInfCat,clusinf1Ind2);
clusinf1_2 = clusinfs2;
} else {
clusinf1_2 = HCategoryManager::newObjectCopy(clusinf1_2,clusinf1,clusInfCat,clusinf1Ind2);
}
clusinf1_2->setClusIndex(clusInd2);
clusinf1_2->setClusFitIndex(clusfitInd2);
}
if( mdchit0 && hitCat){
HLocation loc;
loc.set(3,0,0,0);
loc[0] = mdchit0->getSector();
loc[1] = mdchit0->getModule();
loc[2] = ctMdcHit[mdchit0->getSector()][mdchit0->getModule()];
mdchit0_2 =(HMdcHit*) hitCat->getSlot(loc,&mdchit0Ind2);
if(isSim){
new (mdchit0_2) HMdcHitSim(*((HMdcHitSim*)mdchit0));
} else {
new (mdchit0_2) HMdcHit(*mdchit0);
}
mdchit0_2->setClusInfIndex(clusinf0Ind2);
ctMdcHit[mdchit0->getSector()][mdchit0->getModule()]++;
}
if( mdchit1 && hitCat){
HLocation loc;
loc.set(3,0,0,0);
loc[0] = mdchit1->getSector();
loc[1] = mdchit1->getModule();
loc[2] = ctMdcHit[mdchit1->getSector()][mdchit1->getModule()];
mdchit1_2 =(HMdcHit*) hitCat->getSlot(loc,&mdchit1Ind2);
if(isSim){
new (mdchit1_2) HMdcHitSim(*((HMdcHitSim*)mdchit1));
} else {
new (mdchit1_2) HMdcHit(*mdchit1);
}
mdchit1_2->setClusInfIndex(clusinf1Ind2);
ctMdcHit[mdchit1->getSector()][mdchit1->getModule()]++;
}
if(seg) {
seg->setHitInd(0,mdchit0Ind2);
seg->setHitInd(1,mdchit1Ind2);
}
if(wireFitCat)
{
Int_t first = -1;
Int_t last = -1;
if(clusfit) {
first = clusfit->getFirstWireFitInd();
last = clusfit->getLastWireFitInd();
}
TObjArray* awire = HParticleTool::getMdcWireFitSeg(segInd1);
Int_t ind2 = 0;
Int_t first2 = -1;
Int_t last2 = -1;
for(Int_t i = first; i <= last; i ++){
if(isSim){
HMdcWireFitSim* hit1 = 0;
HMdcWireFitSim* hit2 = 0;
hit1 = HCategoryManager::getObject(hit1,catMdcWireFit,i);
hit2 = HCategoryManager::newObjectCopy(hit2,hit1,wireFitCat,ind2);
} else {
HMdcWireFit* hit1 = 0;
HMdcWireFit* hit2 = 0;
hit1 = HCategoryManager::getObject(hit1,catMdcWireFit,i);
hit2 = HCategoryManager::newObjectCopy(hit2,hit1,wireFitCat,ind2);
}
if(i == first) first2 = ind2;
last2 = ind2;
}
if(clusfit2) {
clusfit2->setFirstWireFitInd(first2);
clusfit2->setLastWireFitInd(last2);
}
delete awire;
}
}
void HParticleTree::extractMdcCal1(Bool_t isSim,Int_t segInd)
{
HMatrixCategory* catOut = (HMatrixCategory*) fCurrentEvent->getCategory(catMdcCal1) ;
if(catOut){
TObjArray* ar = HParticleTool::getMdcCal1Seg(segInd);
if(ar){
HLocation loc;
loc.set(4,0,0,0,0);
Int_t n = ar->GetSize();
for(Int_t i = 0; i < n ; i++){
HMdcCal1* hit1 = (HMdcCal1*)ar->At(i);
if(hit1){
loc[0] = hit1 ->getSector();
loc[1] = hit1 ->getModule();
loc[2] = hit1 ->getLayer();
loc[3] = hit1 ->getCell();
if(!catOut->getObject(loc)){
Int_t index=0;
if(isSim){
HMdcCal1Sim* hit2 = 0;
hit2 =(HMdcCal1Sim*) catOut->getSlot(loc,&index);
new (hit2) HMdcCal1Sim(*((HMdcCal1Sim*)hit1));
} else {
HMdcCal1* hit2 = 0;
hit2 =(HMdcCal1*) catOut->getSlot(loc,&index);
new (hit2) HMdcCal1(*hit1);
}
}
}
}
}
delete ar;
}
}
void HParticleTree::extractMdcCal1FromClus(Bool_t isSim,Int_t segInd)
{
HMatrixCategory* catOut = (HMatrixCategory*) fCurrentEvent->getCategory(catMdcCal1) ;
if(catOut){
TObjArray* ar = HParticleTool::getMdcCal1Cluster(segInd);
if(ar){
HLocation loc;
loc.set(4,0,0,0,0);
Int_t n = ar->GetSize();
for(Int_t i = 0; i < n ; i++){
HMdcCal1* hit1 = (HMdcCal1*)ar->At(i);
if(hit1){
loc[0] = hit1 ->getSector();
loc[1] = hit1 ->getModule();
loc[2] = hit1 ->getLayer();
loc[3] = hit1 ->getCell();
if(!catOut->getObject(loc)){
Int_t index=0;
if(isSim){
HMdcCal1Sim* hit2 = 0;
hit2 =(HMdcCal1Sim*) catOut->getSlot(loc,&index);
new (hit2) HMdcCal1Sim(*((HMdcCal1Sim*)hit1));
} else {
HMdcCal1* hit2 = 0;
hit2 =(HMdcCal1*) catOut->getSlot(loc,&index);
new (hit2) HMdcCal1(*hit1);
}
}
}
}
}
delete ar;
}
}
Bool_t HParticleTree::doFullCopy(Cat_t cat)
{
map<Int_t,Int_t>::iterator iter = fmCatNumToFullCopy.find(cat);
if(iter != fmCatNumToFullCopy.end() ) {
return (Bool_t)iter->second;
}
return kFALSE;
}