StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
StFcsCosmicMaker.cxx
1 /*
2  *
3  * \class StFcsCosmicMaker
4  *
5  */
6 
7 #include "StFcsCosmicMaker.h"
8 
9 #include "StRoot/StEvent/StEvent.h"
10 #include "StRoot/St_base/StMessMgr.h"
11 #include "StRoot/StEvent/StFcsCollection.h"
12 #include "StRoot/StEvent/StFcsHit.h"
13 #include "StRoot/StEvent/StFcsCluster.h"
14 #include "StRoot/StFcsDbMaker/StFcsDb.h"
15 
16 #include "TH1F.h"
17 #include "TH2F.h"
18 #include "TString.h"
19 #include "TFile.h"
20 
21 StFcsCosmicMaker::StFcsCosmicMaker(const Char_t* name) : StMaker(name) {
22  mNTowerThre[0]=10;
23  mSigmaMaxThre[0]=3.0;
24  mSigmaMinThre[0]=0.5;
25  mThetaThre[0]=15;
26  mNTowerThre[1]=8;
27  mSigmaMaxThre[1]=2.0;
28  mSigmaMinThre[1]=0.5;
29  mThetaThre[1]=15;
30 };
31 
32 StFcsCosmicMaker::~StFcsCosmicMaker(){};
33 
34 Int_t StFcsCosmicMaker::Init(){
35  SetAttr(".Privilege",1); //give Privilege to skip a event
36 
37  mFcsDb = static_cast<StFcsDb*>(GetDataSet("fcsDb"));
38  if(!mFcsDb){
39  LOG_FATAL << "Error finding StFcsDb"<< endm;
40  return kStFatal;
41  }
42 
43  if(mSetFile==0){
44  int yday=mRun/1000;
45  sprintf(mFilename,"%d/%d.cosmic.root",yday,mRun);
46  printf("StFcsCosmicMaker::Init - Opening %s\n",mFilename);
47  }else{
48  sprintf(mFilename,"%s",mSetFile);
49  }
50  mFile=new TFile(mFilename,"RECREATE");
51 
52  const char* nameEHP[kFcsEHP] = {"Ecal","Hcal","Pres"};
53  char t[100],t2[100];
54  for(int ehp=0; ehp<2; ehp++){
55  sprintf(t,"%4s_ADC",nameEHP[ehp]);
56  mAdc[ehp] = new TH1F(t,t,100,0.0,100.0);
57  sprintf(t,"%4s_NTower",nameEHP[ehp]);
58  mNTower[ehp] = new TH1F(t,t,100,0.0,100.0);
59  sprintf(t,"%4s_SigmaMax",nameEHP[ehp]);
60  mSigmaMax[ehp] = new TH1F(t,t,100,0.0,10.0);
61  sprintf(t,"%4s_SigmaMin",nameEHP[ehp]);
62  mSigmaMin[ehp] = new TH1F(t,t,100,0.0,3.0);
63  sprintf(t,"%4s_Sigma",nameEHP[ehp]);
64  sprintf(t2,"%4s; SigmaMax; SigmaMin",nameEHP[ehp]);
65  mSigma[ehp] = new TH2F(t,t2,100,0.0,10.0,100,0.0,3.0);
66  sprintf(t,"%4s_SigmaNtow",nameEHP[ehp]);
67  sprintf(t2,"%4s; SigmaMax; NTower",nameEHP[ehp]);
68  mSigmaNtow[ehp] = new TH2F(t,t2,100,0.0,10.0,100,0.0,100.0);
69  }
70  return kStOK;
71 };
72 
74  LOG_INFO << "StFcsCosmicMaker::Make()" << endm;
75  mFcsCollection=0;
76 
77  StEvent* event = (StEvent*)GetInputDS("StEvent");
78  if(!event) {
79  LOG_INFO << "No StEvent found" << endm;
80  }else{
81  mFcsCollection=event->fcsCollection();
82  }
83  if(!mFcsCollection){
84  LOG_INFO << "No StFcsCollection found" << endm;
85  return kStErr;
86  }
87 
88  static int nEvent=0;
89  static int nCosmic=0;
90  int flag=0;
91  for(int det=0; det<4; det++){
92  int ehp=det/2;
93  int nCluster=mFcsCollection->numberOfClusters(det);
94  StSPtrVecFcsCluster& clusters = mFcsCollection->clusters(det);
95  for (int i=0; i<nCluster; i++){
96  StFcsCluster* clu=clusters[i];
97  int nt = clu->nTowers();
98  float smax = clu->sigmaMax();
99  float smin= clu->sigmaMin();
100  float theta = clu->theta() * 180.0 / 3.141592654;
101  while(theta>90.0) theta-=180.0; //reduce to +-90 degree
102  theta = fabs(theta); //reduce to 0~90 degree from x axis
103  mNTower[ehp]->Fill(float(nt));
104  mSigmaMax[ehp]->Fill(smax);
105  mSigmaMin[ehp]->Fill(smin);
106  mSigma[ehp]->Fill(smax,smin);
107  mSigmaNtow[ehp]->Fill(smax,nt);
108  if(nt > mNTowerThre[ehp]){
109  LOG_INFO << Form("StFcsCosmicMaker NTow=%3d SMax=%8.5f SMin=%8.5f",nt,smax,smin)<<endm;
110  if(smax > mSigmaMaxThre[ehp] && smin < mSigmaMinThre[ehp] && theta>mThetaThre[ehp]){
111  //cosmic candidate
112  for(int j=0; j<nt; j++){
113  StFcsHit* hit = clu->hits()[j];
114  int adc = hit->adcSum();
115  mAdc[ehp]->Fill(adc);
116  }
117  flag=1;
118  }
119  }
120  }
121  }
122  nEvent++;
123  if(flag==1){ //found cosmic candidate
124  nCosmic++;
125  LOG_INFO << Form("StFcsCosmicMaker found cosmic candidate %d / %d events",nCosmic,nEvent)<<endm;
126  return kStOK;
127  }else{ // oterwise skip event display
128  LOG_INFO << Form("StFcsCosmicMaker did not found cosmic candidate %d / %d events, skipping",nCosmic,nEvent)<<endm;
129  return kStSkip;
130  }
131 };
132 
134  mFile->Write();
135  mFile->Close();
136  printf("StFcsCosmicMaker::Finish - Closing %s\n",mFilename);
137  return kStOK;
138 };
139 
140 ClassImp(StFcsCosmicMaker);
141 
142 /*
143  * $Id: StFcsCosmicMaker.cxx,v 1.9 2021/03/30 13:29:27 akio Exp $
144  * $Log: StFcsCosmicMaker.cxx,v $
145  *
146  */
virtual Int_t Make()
Definition: Stypes.h:49
Definition: Stypes.h:40
Definition: Stypes.h:44
virtual Int_t Finish()