TSamuraiFDC2Physics.cxx 16.7 KB
Newer Older
Adrien Matta's avatar
Adrien Matta committed
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35
/*****************************************************************************
 * Copyright (C) 2009-2016   this file is part of the NPTool Project         *
 *                                                                           *
 * For the licensing terms see $NPTOOL/Licence/NPTool_Licence                *
 * For the list of contributors see $NPTOOL/Licence/Contributors             *
 *****************************************************************************/

/*****************************************************************************
 * Original Author: Adrien MATTA  contact address: matta@lpccaen.in2p3.fr    *
 *                                                                           *
 * Creation Date  : October 2020                                             *
 * Last update    :                                                          *
 *---------------------------------------------------------------------------*
 * Decription:                                                               *
 *  This class hold SamuraiFDC2 treated data                                 *
 *                                                                           *
 *---------------------------------------------------------------------------*
 * Comment:                                                                  *
 *                                                                           *
 *****************************************************************************/
#include "TSamuraiFDC2Physics.h"

//   STL
#include <sstream>
#include <iostream>
#include <cmath>
#include <stdlib.h>
#include <limits>

//   NPL
#include "RootInput.h"
#include "RootOutput.h"
#include "TAsciiFile.h"
#include "NPOptionManager.h"
#include "NPDetectorFactory.h"
Adrien Matta's avatar
Adrien Matta committed
36
#include "NPSystemOfUnits.h"
Adrien Matta's avatar
Adrien Matta committed
37 38
//   ROOT
#include "TChain.h"
Adrien Matta's avatar
Adrien Matta committed
39 40 41 42
#include "TVector2.h"
#include "Math/Minimizer.h"
#include "Math/Functor.h"
#include "Math/Factory.h"
Adrien Matta's avatar
Adrien Matta committed
43
using namespace NPUNITS;
Adrien Matta's avatar
Adrien Matta committed
44 45 46 47 48 49 50 51 52
///////////////////////////////////////////////////////////////////////////

ClassImp(TSamuraiFDC2Physics)
  ///////////////////////////////////////////////////////////////////////////
  TSamuraiFDC2Physics::TSamuraiFDC2Physics(){
    m_EventData         = new TSamuraiFDC2Data ;
    m_PreTreatedData    = new TSamuraiFDC2Data ;
    m_EventPhysics      = this ;
    //m_Spectra           = NULL;
Adrien Matta's avatar
Adrien Matta committed
53
    ToTThreshold = 100;
Adrien Matta's avatar
Adrien Matta committed
54 55 56 57 58 59 60 61 62 63
  }

///////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Physics::BuildSimplePhysicalEvent(){
  BuildPhysicalEvent();
}

///////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Physics::BuildPhysicalEvent(){
  PreTreat();
Adrien Matta's avatar
Adrien Matta committed
64 65
  RemoveNoise();

Adrien Matta's avatar
Adrien Matta committed
66
  // Map[detector & plane angle, vector of spatial information]
Adrien Matta's avatar
Adrien Matta committed
67 68 69 70
  static map<std::pair<unsigned int,double>, vector<double> > X ; 
  static map<std::pair<unsigned int,double>, vector<double> > Z ; 
  static map<std::pair<unsigned int,double>, vector<double> > R ; 
  X.clear();Z.clear();R.clear();
Adrien Matta's avatar
Adrien Matta committed
71 72
  unsigned int size = Detector.size();
  for(unsigned int i = 0 ; i < size ; i++){
Adrien Matta's avatar
Adrien Matta committed
73 74 75 76 77 78 79 80 81 82 83 84 85
    if(Matched[i]){
      int det = Detector[i];
      int layer = Layer[i];
      int wire = Wire[i]; 
      SamuraiDCIndex idx(det,layer,wire);
      std::pair<unsigned int, double> p(det,Wire_T[idx]);
      X[p].push_back(Wire_X[idx]); 
      Z[p].push_back(Wire_Z[idx]); 
      R[p].push_back(DriftLength[i]); 
    }
  }

  // Reconstruct the vector for each of the plane of each of the detector
Adrien Matta's avatar
Adrien Matta committed
86 87 88 89
  static double X0,X100;
  static map<std::pair<unsigned int,double>, TVector3 > VX0 ;  
  static map<std::pair<unsigned int,double>, TVector3 > VX100 ;  
  VX0.clear();VX100.clear();
Adrien Matta's avatar
Adrien Matta committed
90
  for(auto it = X.begin();it!=X.end();++it){
Adrien Matta's avatar
Adrien Matta committed
91
      Track2D(X[it->first],Z[it->first],R[it->first],X0,X100); 
Adrien Matta's avatar
Adrien Matta committed
92 93
      Mult=X[it->first].size();
      // Position at z=0
Adrien Matta's avatar
Adrien Matta committed
94
      TVector3 P(X0,0,0);
Adrien Matta's avatar
Adrien Matta committed
95
      P.RotateZ(-it->first.second);
Adrien Matta's avatar
Adrien Matta committed
96
      VX0[it->first]=P;
Adrien Matta's avatar
Adrien Matta committed
97 98

      // Direction of the vector in the plane
Adrien Matta's avatar
Adrien Matta committed
99
      TVector3 D = TVector3(X100,0,0);
Adrien Matta's avatar
Adrien Matta committed
100
      D.RotateZ(-it->first.second);
Adrien Matta's avatar
Adrien Matta committed
101
      VX100[it->first]=D;
Adrien Matta's avatar
Adrien Matta committed
102 103
    }

Adrien Matta's avatar
Adrien Matta committed
104
  // Reconstruct the central position (z=0) for each detector
Adrien Matta's avatar
Adrien Matta committed
105 106
  static map<unsigned int,vector<TVector3> > C ;  
  C.clear();
Adrien Matta's avatar
Adrien Matta committed
107
  TVector3 P;
Adrien Matta's avatar
Adrien Matta committed
108 109
  for(auto it1 = VX0.begin();it1!=VX0.end();++it1){
    for(auto it2 = it1;it2!=VX0.end();++it2){
Adrien Matta's avatar
Adrien Matta committed
110 111 112 113 114 115
      if(it1!=it2 && it1->first.first==it2->first.first){// different plane, same detector
        ResolvePlane(it1->second,it1->first.second,it2->second,it2->first.second,P);
        C[it1->first.first].push_back(P);
        }
      }
  }
Adrien Matta's avatar
Adrien Matta committed
116 117 118 119 120 121 122 123 124 125 126
   // Reconstruct the central position (z=0) for each detector
  static map<unsigned int,vector<TVector3> > C100 ;  
  C100.clear();
  for(auto it1 = VX100.begin();it1!=VX100.end();++it1){
    for(auto it2 = it1;it2!=VX100.end();++it2){
      if(it1!=it2 && it1->first.first==it2->first.first){// different plane, same detector
        ResolvePlane(it1->second,it1->first.second,it2->second,it2->first.second,P);
        C100[it1->first.first].push_back(P);
        }
      }
  }
Adrien Matta's avatar
Adrien Matta committed
127
  
Adrien Matta's avatar
Adrien Matta committed
128 129
  // Build the Reference position by averaging all possible pair 
  size = C[2].size();
Adrien Matta's avatar
Adrien Matta committed
130
  double PosX100,PosY100;
Adrien Matta's avatar
Adrien Matta committed
131 132
  if(size){
    PosX=0;
Adrien Matta's avatar
Adrien Matta committed
133
    PosY=0;
Adrien Matta's avatar
Adrien Matta committed
134 135
    PosX100=0;
    PosY100=0;
Adrien Matta's avatar
Adrien Matta committed
136 137
    for(unsigned int i = 0 ; i < size ; i++){
      PosX+= C[2][i].X(); 
Adrien Matta's avatar
Adrien Matta committed
138
      PosY+= C[2][i].Y(); 
Adrien Matta's avatar
Adrien Matta committed
139 140
      PosX100+= C100[2][i].X(); 
      PosY100+= C100[2][i].Y(); 
Adrien Matta's avatar
Adrien Matta committed
141
    } 
Adrien Matta's avatar
Adrien Matta committed
142 143

    // Mean position at Z=0
Adrien Matta's avatar
Adrien Matta committed
144
    PosX/=size; 
Adrien Matta's avatar
Adrien Matta committed
145
    PosY/=size;
Adrien Matta's avatar
Adrien Matta committed
146 147 148 149 150 151 152 153 154 155 156 157 158
    // Mean position at Z=100
    PosX100/=size; 
    PosY100/=size;
    // Compute ThetaX, angle between the Direction vector projection in XZ with
    // the Z axis
    TVector3 Proj=TVector3(PosX100-PosX,0,100);
    ThetaX=atan((PosX100-PosX)/100.);
    // Compute PhiY, angle between the Direction vector projection in YZ with
    // the Z axis
    Proj.SetXYZ(0,PosY100-PosY,100);
    Proj=Proj.Unit();
    PhiY=asin(Proj.X());
    Dir=TVector3(PosX100-PosX,PosY100-PosY,100).Unit();
Adrien Matta's avatar
Adrien Matta committed
159
  }
Adrien Matta's avatar
Adrien Matta committed
160

Adrien Matta's avatar
Adrien Matta committed
161 162
  return;
}
Adrien Matta's avatar
Adrien Matta committed
163 164 165 166 167 168 169 170 171
////////////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Physics::ResolvePlane(const TVector3& H,const double& ThetaU ,const TVector3& L, const double& ThetaV, TVector3& PosXY){
  // direction of U and V wire 
  TVector3 u = TVector3(0,1,0);
  u.RotateZ(ThetaU); 
  
  TVector3 v = TVector3(0,1,0);
  v.RotateZ(ThetaV); 

Adrien Matta's avatar
Adrien Matta committed
172

Adrien Matta's avatar
Adrien Matta committed
173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189
  // Compute the coeff of the two line of vecotr u (v) going through H (L)
  // dv : y = av*x+bv
  TVector3 HH = H+v;
  double av = (H.Y() - HH.Y())/(H.X()-HH.X());
  double bv =  H.Y() - av*H.X();   
   
  // du : y = au*x+bv
  TVector3 LL = L+u;
  double au = (L.Y() - LL.Y())/(L.X()-LL.X());
  double bu =  L.Y() - au*L.X();   

  // We look for M(xM, yM) that intersect du and dv:
  double xM,yM;
  if(!isinf(au) && !isinf(av)){ // au and av are not inf, i.e. not vertical line
    xM = (bu-bv)/(av-au); 
    yM = au*xM+bu;
   } 
Adrien Matta's avatar
Adrien Matta committed
190
  else if(isinf(av)){// av is inf, so v is along Y axis, H is direct measure of X
Adrien Matta's avatar
Adrien Matta committed
191 192 193
    xM = H.X(); 
    yM = au*xM+bu;
    }
Adrien Matta's avatar
Adrien Matta committed
194
  else if (isinf(au)){//au is inf, so u is along Y axis, L is direct measure of X
Adrien Matta's avatar
Adrien Matta committed
195 196 197 198 199 200 201 202 203 204
    xM = L.X(); 
    yM = av*xM+bv;
    }
  else{ // all is lost
    xM=-10000;
    yM=-10000;
    }
  PosXY=TVector3(xM,yM,0);

  };
Adrien Matta's avatar
Adrien Matta committed
205
///////////////////////////////////////////////////////////////////////////
Adrien Matta's avatar
Adrien Matta committed
206
void TSamuraiFDC2Physics::Track2D(const vector<double>& X,const vector<double>& Z,const vector<double>& R,double& X0,double& X100 ){
Adrien Matta's avatar
Adrien Matta committed
207 208 209 210 211 212
  fitX=X;
  fitZ=Z;
  fitR=R;
  // assume all X,Z,R of same size
  unsigned int size = X.size();
  // Define the starting point of the fit: a straight line passing through the 
Adrien Matta's avatar
Adrien Matta committed
213
  // the first and last wire
Adrien Matta's avatar
Adrien Matta committed
214 215 216 217 218 219 220 221 222 223 224
  // z = ax+b -> x=(z-b)/a
  double a = fitZ[size-1]-fitZ[0]/(fitX[size-1]-fitX[0]);
  double b = fitZ[0]-a*fitX[0];
  double parameter[2]={a,b};
  static ROOT::Math::Minimizer* min = ROOT::Math::Factory::CreateMinimizer("Minuit2", "Migrad");
  static ROOT::Math::Functor f(this,&TSamuraiFDC2Physics::SumD,2); 
  min->SetFunction(f);
  min->SetVariable(0,"a",parameter[0], 1);
  min->SetVariable(1,"b",parameter[1], 1);
  min->Minimize(); 
  const double *xs = min->X();
Adrien Matta's avatar
Adrien Matta committed
225 226
  X0=-xs[1]/xs[0];
  X100=(100-xs[1])/xs[0];
Adrien Matta's avatar
Adrien Matta committed
227 228 229 230 231 232
}
////////////////////////////////////////////////////////////////////////////////
double TSamuraiFDC2Physics::SumD(const double* parameter ){
  unsigned int size = fitX.size();
  // Compute the sum P of the distance between the circle and the track
  double P = 0;
Adrien Matta's avatar
Adrien Matta committed
233 234 235 236 237
  double a = parameter[0];
  double b = parameter[1];
  double ab= a*b;
  double a2=a*a;

Adrien Matta's avatar
Adrien Matta committed
238
  for(unsigned int i = 0 ; i < size ; i++){
Adrien Matta's avatar
Adrien Matta committed
239 240 241 242 243
    double c = fitX[i];
    double d = fitZ[i];
    double x = (a*d-ab+c)/(1+a2);
    double z = a*x+b;
    P+= sqrt(abs( (x-c)*(x-c)+(z-d)*(z-d)-fitR[i]*fitR[i]));
Adrien Matta's avatar
Adrien Matta committed
244
    }
Adrien Matta's avatar
Adrien Matta committed
245

Adrien Matta's avatar
Adrien Matta committed
246 247
  return P;
  }
Adrien Matta's avatar
Adrien Matta committed
248 249 250
///////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Physics::PreTreat(){
  ClearPreTreatedData();
Adrien Matta's avatar
Adrien Matta committed
251 252 253 254 255 256 257
  static CalibrationManager* Cal = CalibrationManager::getInstance();
  static string channel;
   // one map per detector
  map<unsigned int, vector<double> > X ; 
  map<unsigned int, vector<double> > Z ; 
  map<unsigned int, vector<double> > R ; 

Adrien Matta's avatar
Adrien Matta committed
258 259
  unsigned int size = m_EventData->Mult();
  for(unsigned int i = 0 ; i < size ; i++){
Adrien Matta's avatar
Adrien Matta committed
260 261 262
    // EDGE=1 is the leading edge, IE, the real time.
    // EDGE=0 is the trailing edge, so it helps build Tot
    if(m_EventData->GetEdge(i)==1){
Adrien Matta's avatar
Adrien Matta committed
263 264 265 266 267 268 269
      int det   = m_EventData->GetDetectorNbr(i); 
      int layer = m_EventData->GetLayerNbr(i); 
      int wire  = m_EventData->GetWireNbr(i); 
      double time = m_EventData->GetTime(i);
      double etime = 0;
      // look for matching trailing edge   
      for(unsigned int j = 0 ; j < size ; j++){
Adrien Matta's avatar
Adrien Matta committed
270
        if(m_EventData->GetEdge(j)==0){
Adrien Matta's avatar
Adrien Matta committed
271 272 273 274 275 276 277 278
          int edet   = m_EventData->GetDetectorNbr(j); 
          int elayer = m_EventData->GetLayerNbr(j); 
          int ewire  = m_EventData->GetWireNbr(j); 
          // same wire
          if(wire==ewire && layer==elayer && det==edet){
            etime = m_EventData->GetTime(j); 
          }    
        }
Adrien Matta's avatar
Adrien Matta committed
279 280 281 282
        if(etime && etime>time)
          break;
        else
          etime=0;
Adrien Matta's avatar
Adrien Matta committed
283 284
      }
      // a valid wire must have an edge
Adrien Matta's avatar
Adrien Matta committed
285
      if(etime && time && etime-time>ToTThreshold){
Adrien Matta's avatar
Adrien Matta committed
286 287 288 289
        Detector.push_back(det);
        Layer.push_back(layer);       
        Wire.push_back(wire);
        Time.push_back(time);
Adrien Matta's avatar
Adrien Matta committed
290
        ToT.push_back(etime-time);
Adrien Matta's avatar
Adrien Matta committed
291 292
        channel="SamuraiFDC2/L" + NPL::itoa(layer);
        DriftLength.push_back(Cal->ApplySigmoid(channel,etime));
Adrien Matta's avatar
Adrien Matta committed
293
      }
Adrien Matta's avatar
Adrien Matta committed
294 295 296 297 298
    }

  }
  return;
}
Adrien Matta's avatar
Adrien Matta committed
299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324
///////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Physics::RemoveNoise(){
  // Remove the noise by looking if a matching wire exist in the adjacent layer
  // this done by looking at the closest plane with the same orientation
  unsigned int size = Detector.size(); 
  for(unsigned int i = 0 ; i < size ; i++){
    bool match=false;
    int det = Detector[i];
    int layer = Layer[i];
    int wire = Wire[i];
    // look for matching adjacent wire   

    for(unsigned int j = 0 ; j < size ; j++){
      int adet = Detector[j];
      int alayer = Layer[j];
      int awire = Wire[j];
      bool blayer = false;
      if(layer%2==0){
        if(layer+1==alayer)
          blayer=true;
      }

      else{
        if(layer-1==alayer)
          blayer=true;
      }
Adrien Matta's avatar
Adrien Matta committed
325

Adrien Matta's avatar
Adrien Matta committed
326 327 328 329 330
      if(det==adet && blayer && abs(wire-awire)<=1){
        match=true;
        break;
      }
    }
Adrien Matta's avatar
Adrien Matta committed
331 332 333 334 335 336

    if(match)
      Matched.push_back(true);
    else
      Matched.push_back(false);
  }
Adrien Matta's avatar
Adrien Matta committed
337 338
  return;
}
Adrien Matta's avatar
Adrien Matta committed
339 340
///////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Physics::Clear(){
Adrien Matta's avatar
Adrien Matta committed
341 342
  Mult=0;
  PosX=PosY=-10000;
Adrien Matta's avatar
Adrien Matta committed
343
  DriftLength.clear();
Adrien Matta's avatar
Adrien Matta committed
344 345 346 347 348 349 350
  Detector.clear();
  Layer.clear();
  Wire.clear();
  Time.clear();
  ToT.clear();
  ParticleDirection.clear();
  MiddlePosition.clear();
Adrien Matta's avatar
Adrien Matta committed
351
  Matched.clear();
Adrien Matta's avatar
Adrien Matta committed
352 353 354 355 356 357 358
}
///////////////////////////////////////////////////////////////////////////

////   Innherited from VDetector Class   ////

///////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Physics::ReadConfiguration(NPL::InputParser parser){
Adrien Matta's avatar
Adrien Matta committed
359 360 361 362 363 364 365 366 367 368 369 370 371
  vector<NPL::InputBlock*> blocks = parser.GetAllBlocksWithToken("SAMURAIFDC2");
  if(NPOptionManager::getInstance()->GetVerboseLevel())
    cout << "//// " << blocks.size() << " detector(s) found " << endl; 

  vector<string> token= {"XML"};

  for(unsigned int i = 0 ; i < blocks.size() ; i++){
    cout << endl << "////  Samurai FDC2 (" << i+1 << ")" << endl;
    string xmlpath = blocks[i]->GetString("XML");
    NPL::XmlParser xml;
    xml.LoadFile(xmlpath);
    AddDC("SAMURAIFDC2",xml);
  }
Adrien Matta's avatar
Adrien Matta committed
372
}
Adrien Matta's avatar
Adrien Matta committed
373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399

///////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Physics::AddDC(string name, NPL::XmlParser& xml){
  std::vector<NPL::block*> b = xml.GetAllBlocksWithName(name);  
  // FDC2 case
  if(name=="SAMURAIFDC2"){
    unsigned int det=2;
    unsigned int size = b.size();
    for(unsigned int i = 0 ; i < size ; i++){
      unsigned int layer = b[i]->AsInt("layer"); 
      unsigned int wire  = b[i]->AsInt("wireid"); 
      double X = b[i]->AsDouble("wirepos");  
      double Z = b[i]->AsDouble("wirez");  
      string sDir = b[i]->AsString("anodedir");
      double T=0;
      if(sDir=="X")
        T=0*deg;
      else if(sDir=="Y")
        T=90*deg;
      else if(sDir=="U")
        T=30*deg;
      else if(sDir=="V")
        T=-30*deg;
      else{
        cout << "ERROR: Unknown layer orientation for Samurai FDC2"<< endl;
        exit(1);
      }
Adrien Matta's avatar
Adrien Matta committed
400 401 402 403
      SamuraiDCIndex idx(det,layer,wire);
      Wire_X[idx]=X;
      Wire_Z[idx]=Z;
      Wire_T[idx]=T;
Adrien Matta's avatar
Adrien Matta committed
404 405 406 407
    }
  }
}

Adrien Matta's avatar
Adrien Matta committed
408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445
///////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Physics::InitSpectra(){  
  //m_Spectra = new TSamuraiFDC2Spectra(m_NumberOfDetector);
}

///////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Physics::FillSpectra(){  
  //  m_Spectra -> FillRawSpectra(m_EventData);
  //  m_Spectra -> FillPreTreatedSpectra(m_PreTreatedData);
  //  m_Spectra -> FillPhysicsSpectra(m_EventPhysics);
}
///////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Physics::CheckSpectra(){  
  //  m_Spectra->CheckSpectra();  
}
///////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Physics::ClearSpectra(){  
  // To be done
}
///////////////////////////////////////////////////////////////////////////
map< string , TH1*> TSamuraiFDC2Physics::GetSpectra() {
  /*  if(m_Spectra)
      return m_Spectra->GetMapHisto();
      else{
      map< string , TH1*> empty;
      return empty;
      }*/
  map< string , TH1*> empty;
  return empty;

} 

///////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Physics::WriteSpectra(){
  // m_Spectra->WriteSpectra();
}
///////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Physics::AddParameterToCalibrationManager(){
Adrien Matta's avatar
Adrien Matta committed
446 447 448 449 450 451 452
  CalibrationManager* Cal = CalibrationManager::getInstance();

  // each layer
  for( int l = 0 ; l < 14 ; ++l){
    Cal->AddParameter("SamuraiFDC2", "L"+ NPL::itoa(l),"FDC2_L"+ NPL::itoa(l));
  }

Adrien Matta's avatar
Adrien Matta committed
453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499
}

///////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Physics::InitializeRootInputRaw(){
  TChain* inputChain = RootInput::getInstance()->GetChain()   ;
  inputChain->SetBranchStatus( "SamuraiFDC2" , true );
  // The following line is necessary only for system were the tree is splitted
  // (older root version). The found argument silenced the Branches not found
  // warning for non splitted tree.
  if(inputChain->FindBranch("fDC_*"))
    inputChain->SetBranchStatus( "fDC_*",true);
  inputChain->SetBranchAddress( "SamuraiFDC2" , &m_EventData );

}

///////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Physics::InitializeRootInputPhysics(){
}

///////////////////////////////////////////////////////////////////////////
void TSamuraiFDC2Physics::InitializeRootOutput(){
  TTree* outputTree = RootOutput::getInstance()->GetTree();
  outputTree->Branch( "SamuraiFDC2" , "TSamuraiFDC2Physics" , &m_EventPhysics );
}

////////////////////////////////////////////////////////////////////////////////
//            Construct Method to be pass to the DetectorFactory              //
////////////////////////////////////////////////////////////////////////////////
NPL::VDetector* TSamuraiFDC2Physics::Construct(){
  return (NPL::VDetector*) new TSamuraiFDC2Physics();
}

////////////////////////////////////////////////////////////////////////////////
//            Registering the construct method to the factory                 //
////////////////////////////////////////////////////////////////////////////////
extern "C"{
  class proxy_samuraiFDC2{
    public:
      proxy_samuraiFDC2(){
        NPL::DetectorFactory::getInstance()->AddToken("SAMURAIFDC2","Samurai");
        NPL::DetectorFactory::getInstance()->AddDetector("SAMURAIFDC2",TSamuraiFDC2Physics::Construct);
      }
  };

  proxy_samuraiFDC2 p_samuraiFDC2;
}