/[PAMELA software]/PamVMC/include/PamVMCFieldMgr.h
ViewVC logotype

Annotation of /PamVMC/include/PamVMCFieldMgr.h

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.6 - (hide annotations) (download)
Wed Sep 15 07:05:42 2010 UTC (14 years, 5 months ago) by pizzolot
Branch: MAIN
CVS Tags: HEAD
Changes since 1.5: +16 -17 lines
File MIME type: text/plain
new particle generation methods

1 nikolas 1.1 #ifndef PAMVMCFIELDMGR_H
2     #define PAMVMCFIELDMGR_H
3     #include <iostream>
4 pam-rm2 1.5 #include <RVersion.h>
5 nikolas 1.1
6     #include "TString.h"
7    
8     #define PAMFIELD 1
9    
10    
11     #define ZCAL 13.05 //cm
12     #define CALS3 10.639 //cm
13     #define ZSPEBP 2.97 //cm
14     #define HORMG 22.57 //cm
15    
16 pizzolot 1.6 //#ifdef PAMFIELD
17 nikolas 1.1 #include "TrkParams.h"
18 pizzolot 1.6 //#endif
19 nikolas 1.1
20 pam-rm2 1.5 #if ROOT_VERSION_CODE >= 333572
21     #include <TVirtualMagField.h>
22     class PamVMCFieldMgr: public TVirtualMagField
23     #else
24     class PamVMCFieldMgr: public TObject
25     #endif
26     {
27 nikolas 1.1 private:
28    
29     static PamVMCFieldMgr * fm;
30    
31     protected:
32     PamVMCFieldMgr() {};
33    
34     public:
35    
36     static PamVMCFieldMgr * Instance();
37    
38    
39     void LoadField(){
40     #ifdef PAMFIELD
41     // Load the PAMELA magnetic field
42    
43     std::string pamcal=getenv("PAM_CALIB"); pamcal+="/trk-param/field_param-0/";
44     std::cout << "PAMELA env: " << pamcal << std::endl;
45     TrkParams::SetTrackingMode();
46     TrkParams::SetPrecisionFactor();
47     TrkParams::SetStepMin();
48     TrkParams::Set(pamcal.c_str(),1);
49     TrkParams::Load(1);
50     #endif
51     }
52    
53    
54 pam-rm2 1.5 void Field(const Double_t* x, Double_t* B) {
55 nikolas 1.1
56     #ifdef PAMFIELD
57     float xm[3]={x[0],x[1],x[2]-ZCAL-CALS3-ZSPEBP-HORMG};
58    
59 pam-rm2 1.5 B[0] = TrkParams::GetBX((float *)xm)*10;
60     B[1] = TrkParams::GetBY((float *)xm)*10;
61     B[2] = TrkParams::GetBZ((float *)xm)*10;
62 nikolas 1.1
63     #else
64     if(x) {
65 pam-rm2 1.5 for (Int_t i=0; i<3; i++) B[i] = 0.0;
66 nikolas 1.1 } else {
67 pam-rm2 1.5 B[0]=B[1]=B[2]=0.0;
68 nikolas 1.1 }
69     #endif
70    
71     }
72    
73 pam-rm2 1.5
74     void Field(const Float_t* x, Float_t* b) {
75 nikolas 1.1
76    
77 pizzolot 1.6 #ifdef PAMFIELD
78     float xm[3]={x[0],x[1],x[2]-ZCAL-CALS3-ZSPEBP-HORMG};
79    
80     b[0] = TrkParams::GetBX((float *)xm)*10;
81     b[1] = TrkParams::GetBY((float *)xm)*10;
82     b[2] = TrkParams::GetBZ((float *)xm)*10;
83    
84     #else
85     if(x) {
86     for (Int_t i=0; i<3; i++) b[i] = 0.0;
87     } else {
88     b[0]=b[1]=b[2]=0.0;
89     }
90     #endif
91 pam-rm2 1.5
92 nikolas 1.1 }
93 pam-rm2 1.5
94 nikolas 1.1 };
95    
96     #endif

  ViewVC Help
Powered by ViewVC 1.1.23