SND@LHC Software
Loading...
Searching...
No Matches
ShipGoliathField Class Reference

#include <ShipGoliathField.h>

Inheritance diagram for ShipGoliathField:
Collaboration diagram for ShipGoliathField:

Public Member Functions

 ShipGoliathField ()
 
 ShipGoliathField (const char *name)
 
virtual ~ShipGoliathField ()
 
void Init (const char *fieldfile)
 
void getpos (TString vol, TVector3 &bot, TVector3 &top) const
 return value at position
 
void close ()
 
virtual Double_t GetBx (Double_t x, Double_t y, Double_t z)
 
virtual Double_t GetBy (Double_t x, Double_t y, Double_t z)
 
virtual Double_t GetBz (Double_t x, Double_t y, Double_t z)
 
virtual void Print ()
 
void sethistbxyz (TH3D *histbx, TH3D *histby, TH3D *histbz)
 
TH3D * gethistbx () const
 
TH3D * gethistby () const
 
TH3D * gethistbz () const
 

Public Attributes

Double_t coords [13][6]
 
TFile * fieldmap = NULL
 
Float_t kilogauss = 1.
 
Float_t tesla = 10*kilogauss
 
Float_t cm = 1
 
Float_t m = 100*cm
 
Float_t mm = 0.1*cm
 

Private Member Functions

 ClassDef (ShipGoliathField, 2)
 

Private Attributes

double fMiddle
 
double fPeak
 
int fOrient
 
TH3D * fhistbx
 
TH3D * fhistby
 
TH3D * fhistbz
 
Double_t xmin
 
Double_t xmax
 
Double_t ymin
 
Double_t ymax
 
Double_t zmin
 
Double_t zmax
 

Detailed Description

ShipGoliathField.h

Definition at line 21 of file ShipGoliathField.h.

Constructor & Destructor Documentation

◆ ShipGoliathField() [1/2]

ShipGoliathField::ShipGoliathField ( )

Default constructor

Definition at line 26 of file ShipGoliathField.cxx.

27 : FairField()
28{
29}

◆ ShipGoliathField() [2/2]

ShipGoliathField::ShipGoliathField ( const char *  name)

Standard constructor

Parameters
nameObject name

Definition at line 35 of file ShipGoliathField.cxx.

36 : FairField(name)
37{
38}

◆ ~ShipGoliathField()

ShipGoliathField::~ShipGoliathField ( )
virtual

Destructor

Definition at line 43 of file ShipGoliathField.cxx.

43{ }

Member Function Documentation

◆ ClassDef()

ShipGoliathField::ClassDef ( ShipGoliathField  ,
 
)
private

◆ close()

void ShipGoliathField::close ( )

Definition at line 128 of file ShipGoliathField.cxx.

128 {
129 //char* mypath = std::getenv("mySHIPSOFT");
130 //strcat(mypath,"/FairShip/field/GoliathFieldMap.root");
131 //fieldmap->Close(mypath);
132 //gROOT->cd();
133}

◆ GetBx()

Double_t ShipGoliathField::GetBx ( Double_t  x,
Double_t  y,
Double_t  z 
)
virtual

Get components of field at a given point

Parameters
x,y,zPoint coordinates [cm]

Definition at line 136 of file ShipGoliathField.cxx.

136 {
137 Double_t bx=0.;
138 TH3D* hbx;
139
140 if ((x < xmin )|| (x> xmax) || (y < ymin) || (y>ymax) || ((z-350.75)<zmin) || ((z-350.75)>zmax)) {
141 return bx;}
142 // FairShip: 0 after absorber -384.5<target<-240 -239.9<absorber<0<T1T2<121<Goliath<481<T3T4<755 766.6<RPC<966.6
144 Int_t binx = hbx->GetXaxis()->FindBin(x);
145 Int_t biny = hbx->GetYaxis()->FindBin(y);
146 Int_t binz = hbx->GetZaxis()->FindBin(z - 350.75);
147 bx=hbx->GetBinContent(binx,biny,binz)*tesla;
148 //cout << "GetBX " << x << ", binx " << binx <<" y "<< y << " biny "<<biny<<" z "<< z << " binz "<<binz<<" Bx= " << bx << endl;
149 return bx;
150}
TH3D * gethistbx() const

◆ GetBy()

Double_t ShipGoliathField::GetBy ( Double_t  x,
Double_t  y,
Double_t  z 
)
virtual

Definition at line 156 of file ShipGoliathField.cxx.

156 {
157 Double_t by=0.;
158 TH3D* hby;
159
160 if ((x < xmin )|| (x> xmax) || (y < ymin) || (y>ymax) || ((z-350.75)<zmin) || ((z-350.75)>zmax)) {
161 return by;}
163 Int_t binx = hby->GetXaxis()->FindBin(x);
164 Int_t biny = hby->GetYaxis()->FindBin(y);
165 Int_t binz = hby->GetZaxis()->FindBin(z- 350.75);
166 by=hby->GetBinContent(binx,biny,binz)*tesla;
167 return by;
168}
TH3D * gethistby() const

◆ GetBz()

Double_t ShipGoliathField::GetBz ( Double_t  x,
Double_t  y,
Double_t  z 
)
virtual

Definition at line 175 of file ShipGoliathField.cxx.

175 {
176 Double_t bz=0.;
177 TH3D* hbz;
178
179 if ((x < xmin )|| (x> xmax) || (y < ymin) || (y>ymax) || ((z-350.75)<zmin) || ((z-350.75)>zmax)) {
180 return bz;}
181
183 Int_t binx = hbz->GetXaxis()->FindBin(x);
184 Int_t biny = hbz->GetYaxis()->FindBin(y);
185 Int_t binz = hbz->GetZaxis()->FindBin(z- 350.75);
186 bz=hbz->GetBinContent(binx,biny,binz)*tesla;
187 return bz;
188}
TH3D * gethistbz() const

◆ gethistbx()

TH3D * ShipGoliathField::gethistbx ( ) const
inline

Definition at line 68 of file ShipGoliathField.h.

68{ return fhistbx;};

◆ gethistby()

TH3D * ShipGoliathField::gethistby ( ) const
inline

Definition at line 69 of file ShipGoliathField.h.

69{ return fhistby;};

◆ gethistbz()

TH3D * ShipGoliathField::gethistbz ( ) const
inline

Definition at line 70 of file ShipGoliathField.h.

70{ return fhistbz;};

◆ getpos()

void ShipGoliathField::getpos ( TString  vol,
TVector3 &  bot,
TVector3 &  top 
) const

return value at position

Definition at line 103 of file ShipGoliathField.cxx.

103 {
104 std::cout<<"ShipGoliathField::getpos: GetCurrentNavigator volname "<<volname<<std::endl;
105 TGeoNavigator* nav;
106 if (gGeoManager) {
107 nav = gGeoManager->GetCurrentNavigator();}
108 else { std::cout<<"No geomanager"<<std::endl;}
109 std::cout<<"ShipGoliathField::getpos: cd to volume "<<volname<<std::endl;
110 Bool_t rc = nav->cd(volname);
111 if (not rc){
112 cout << "ShipGoliathfield::getpos, TGeoNavigator failed "<<volname<<endl;
113 return;
114 }
115 std::cout<<"ShipGoliathField::getpos GetCurrentNode"<<std::endl;
116 TGeoNode* W = nav->GetCurrentNode();
117 TGeoBBox* S = dynamic_cast<TGeoBBox*>(W->GetVolume()->GetShape());
118 Double_t top[3] = {S->GetDX(),S->GetDY(),S->GetDZ()};
119 Double_t bot[3] = {-S->GetDX(),-S->GetDY(),-S->GetDZ()};
120 Double_t Gtop[3],Gbot[3];
121 nav->LocalToMaster(top, Gtop);
122 nav->LocalToMaster(bot, Gbot);
123 vtop.SetXYZ(Gtop[0],Gtop[1],Gtop[2]);
124 vbot.SetXYZ(Gbot[0],Gbot[1],Gbot[2]);
125}
dict S
Definition MufiCTR.py:12

◆ Init()

void ShipGoliathField::Init ( const char *  fieldfile)

Definition at line 65 of file ShipGoliathField.cxx.

65 {
66
67 fieldmap = TFile::Open(fieldfile);
68
69
70 TH3D* histbx= (TH3D*)fieldmap->Get("Bx");
71 TH3D* histby= (TH3D*)fieldmap->Get("By");
72 TH3D* histbz= (TH3D*)fieldmap->Get("Bz");
73 xmin = histbx->GetXaxis()->GetXmin();
74 xmax = histbx->GetXaxis()->GetXmax();
75 ymin = histbx->GetYaxis()->GetXmin();
76 ymax = histbx->GetYaxis()->GetXmax();
77 zmin = histbx->GetZaxis()->GetXmin();
78 zmax = histbx->GetZaxis()->GetXmax();
79
80 ShipGoliathField::sethistbxyz(histbx, histby, histbz);
81
82 /*
83 TVector3 bot,top;
84 std::cout<<"ShipGoliathField::setup: making string vector"<<std::endl;
85 std::vector<TString> volume={"/volGoliath_1/VolVacuum_1","/volGoliath_1/volLateralS1_1","/volGoliath_1/volLateralS2_1","/volGoliath_1/volLateralSurface1low_1",\
86 "/volGoliath_1/volLateralSurface2low_1","/volGoliath_1/volLateralS1_b_1","/volGoliath_1/volLateralS2_b_1","/volGoliath_1/volLateralSurface1blow_1",\
87 "/volGoliath_1/volLateralSurface2blow_1","/volGoliath_1/volLateralS1_d_1","/volGoliath_1/volLateralS2_d_1","/volGoliath_1/volLateralS1_c_1",\
88 "/volGoliath_1/volLateralS2_c_1"};
89 for (Int_t i=0;i<13;i++){
90 std::cout<<"ShipGoliathField::setup: getpos volume "<<volume[i]<<std::endl;
91 ShipGoliathField::getpos(volume[i],bot,top);
92 for (Int_t j=0;j<3;j++) {
93 coords[i][j]=top[j];
94 coords[i][j+3]=bot[j];
95 }
96
97 //std::cout<<volume[i]<<" "<<coords[i][3] << " posX " << posX <<" "<< coords[i][0] << ";" << coords[i][4] << " posY "<< posY <<" "<< coords[i][1] << ";" << coords[i][5] << " posZ " << posZ<<" "<< coords[i][2]<<std::endl;
98 }
99 */
100}
void sethistbxyz(TH3D *histbx, TH3D *histby, TH3D *histbz)

◆ Print()

void ShipGoliathField::Print ( )
virtual

Screen output

Definition at line 196 of file ShipGoliathField.cxx.

196 {
197 cout << "======================================================" << endl;
198 cout << "---- " << fTitle << " : " << fName << endl;
199 cout << "----" << endl;
200 cout << "---- Field type : constant" << endl;
201 cout << "----" << endl;
202 cout << "---- Field regions : " << endl;
203 cout.precision(4);
204 cout << "======================================================" << endl;
205}

◆ sethistbxyz()

void ShipGoliathField::sethistbxyz ( TH3D *  histbx,
TH3D *  histby,
TH3D *  histbz 
)
inline

Definition at line 62 of file ShipGoliathField.h.

62 {
63 fhistbx=histbx;
64 fhistby=histby;
65 fhistbz=histbz;
66 };

Member Data Documentation

◆ cm

Float_t ShipGoliathField::cm = 1

Definition at line 75 of file ShipGoliathField.h.

◆ coords

Double_t ShipGoliathField::coords[13][6]

Definition at line 41 of file ShipGoliathField.h.

◆ fhistbx

TH3D* ShipGoliathField::fhistbx
private

Definition at line 84 of file ShipGoliathField.h.

◆ fhistby

TH3D* ShipGoliathField::fhistby
private

Definition at line 85 of file ShipGoliathField.h.

◆ fhistbz

TH3D* ShipGoliathField::fhistbz
private

Definition at line 86 of file ShipGoliathField.h.

◆ fieldmap

TFile* ShipGoliathField::fieldmap = NULL

Definition at line 60 of file ShipGoliathField.h.

◆ fMiddle

double ShipGoliathField::fMiddle
private

Definition at line 80 of file ShipGoliathField.h.

◆ fOrient

int ShipGoliathField::fOrient
private

Definition at line 82 of file ShipGoliathField.h.

◆ fPeak

double ShipGoliathField::fPeak
private

Definition at line 81 of file ShipGoliathField.h.

◆ kilogauss

Float_t ShipGoliathField::kilogauss = 1.

Definition at line 72 of file ShipGoliathField.h.

◆ m

Float_t ShipGoliathField::m = 100*cm

Definition at line 76 of file ShipGoliathField.h.

◆ mm

Float_t ShipGoliathField::mm = 0.1*cm

Definition at line 77 of file ShipGoliathField.h.

◆ tesla

Float_t ShipGoliathField::tesla = 10*kilogauss

Definition at line 73 of file ShipGoliathField.h.

◆ xmax

Double_t ShipGoliathField::xmax
private

Definition at line 88 of file ShipGoliathField.h.

◆ xmin

Double_t ShipGoliathField::xmin
private

Definition at line 88 of file ShipGoliathField.h.

◆ ymax

Double_t ShipGoliathField::ymax
private

Definition at line 88 of file ShipGoliathField.h.

◆ ymin

Double_t ShipGoliathField::ymin
private

Definition at line 88 of file ShipGoliathField.h.

◆ zmax

Double_t ShipGoliathField::zmax
private

Definition at line 88 of file ShipGoliathField.h.

◆ zmin

Double_t ShipGoliathField::zmin
private

Definition at line 88 of file ShipGoliathField.h.


The documentation for this class was generated from the following files: