FairRoot/PandaRoot
PndMultiField.cxx
Go to the documentation of this file.
1 // -------------------------------------------------------------------------
2 // ----- PndMultiField source file -----
3 // ----- Created 29/01/07 by M. Al/Turany -----
4 // -------------------------------------------------------------------------
5 
6 
7 
8 #include "PndMultiField.h"
9 #include "PndRegion.h"
10 #include "PndConstField.h"
11 #include "PndFieldMap.h"
12 #include "PndMapPar.h"
13 #include "PndMultiFieldPar.h"
14 #include "PndSolenoidMap.h"
15 #include "PndTransMap.h"
16 #include "PndDipoleMap.h"
17 
18 
19 #include "FairRunSim.h"
20 #include "FairRuntimeDb.h"
21 
22 
23 #include "TObjArray.h"
24 
25 
26 
27 #include <iomanip>
28 #include <iostream>
29 #include <fstream>
30 
31 
32 
33 using namespace std;
34 
35 // ------------- Default constructor ----------------------------------
37  : fMaps(new TObjArray(10)), fNoOfMaps(0), fFieldMaps(), fMapIter(), fBeamMom(0.)
38 {
39  fType = 5;
40 }
41 
42 
43 // ------------- Default constructor ----------------------------------
45  : fMaps(new TObjArray(10)), fNoOfMaps(0), fFieldMaps(), fMapIter(), fBeamMom(BeamMom)
46 {
47  fType = 5;
48  if (BeamMom<0.)
49  {
50  FairRunSim *fRun= FairRunSim::Instance();
51  if(fRun) fBeamMom= fRun->GetBeamMom();
52  }
53 
54  Map.ToUpper();
55 
56  if (Map=="FULL") {
57 
58  PndTransMap *map_t= new PndTransMap("TransMap", "R", fBeamMom);
59  PndDipoleMap *map_d1= new PndDipoleMap("DipoleMap1", "R", fBeamMom);
60  PndDipoleMap *map_d2= new PndDipoleMap("DipoleMap2", "R", fBeamMom);
61  PndSolenoidMap *map_s1= new PndSolenoidMap("SolenoidMap1", "R");
62  PndSolenoidMap *map_s2= new PndSolenoidMap("SolenoidMap2", "R");
63  PndSolenoidMap *map_s3= new PndSolenoidMap("SolenoidMap3", "R");
64  PndSolenoidMap *map_s4= new PndSolenoidMap("SolenoidMap4", "R");
65 
66  AddField(map_t);
67  AddField(map_d1);
68  AddField(map_d2);
69  AddField(map_s1);
70  AddField(map_s2);
71  AddField(map_s3);
72  AddField(map_s4);
73 
74  }
75 
76  else if (Map=="AUTO") {
77  PndDipoleMap *map_d1= new PndDipoleMap("DipoleMap1", "R", fBeamMom);
78  PndDipoleMap *map_d2= new PndDipoleMap("DipoleMap2", "R", fBeamMom);
79 
85 
86  if(fBeamMom >= 3.0) {
87  map_t= new PndTransMap("TransMap", "R", fBeamMom);
88  map_s1= new PndSolenoidMap("SolenoidMap1", "R");
89  map_s2= new PndSolenoidMap("SolenoidMap2", "R");
90  map_s3= new PndSolenoidMap("SolenoidMap3", "R");
91  map_s4= new PndSolenoidMap("SolenoidMap4", "R");
92  }
93  else {
94  map_t= new PndTransMap("TransMap_Half", "R", fBeamMom);
95  map_s1= new PndSolenoidMap("SolenoidMap_Half1", "R");
96  map_s2= new PndSolenoidMap("SolenoidMap_Half2", "R");
97  map_s3= new PndSolenoidMap("SolenoidMap_Half3", "R");
98  map_s4= new PndSolenoidMap("SolenoidMap_Half4", "R");
99  }
100 
101  AddField(map_t);
102  AddField(map_d1);
103  AddField(map_d2);
104  AddField(map_s1);
105  AddField(map_s2);
106  AddField(map_s3);
107  AddField(map_s4);
108  }
109 
110  else if (Map=="HALF") {
111  PndTransMap *map_t= new PndTransMap("TransMap_Half", "R", fBeamMom);
112  PndDipoleMap *map_d1= new PndDipoleMap("DipoleMap1", "R", fBeamMom);
113  PndDipoleMap *map_d2= new PndDipoleMap("DipoleMap2", "R", fBeamMom);
114  PndSolenoidMap *map_s1= new PndSolenoidMap("SolenoidMap_Half1", "R");
115  PndSolenoidMap *map_s2= new PndSolenoidMap("SolenoidMap_Half2", "R");
116  PndSolenoidMap *map_s3= new PndSolenoidMap("SolenoidMap_Half3", "R");
117  PndSolenoidMap *map_s4= new PndSolenoidMap("SolenoidMap_Half4", "R");
118 
119  AddField(map_t);
120  AddField(map_d1);
121  AddField(map_d2);
122  AddField(map_s1);
123  AddField(map_s2);
124  AddField(map_s3);
125  AddField(map_s4);
126 
127  }
128 
129 
130  else if (Map=="DIPOLE") {
131  PndDipoleMap *map_d1= new PndDipoleMap("DipoleMap1", "R", fBeamMom);
132  PndDipoleMap *map_d2= new PndDipoleMap("DipoleMap2", "R", fBeamMom);
133 
134  AddField(map_d1);
135  AddField(map_d2);
136 
137  }else if (Map=="SOLENOID") {
138 
139  PndSolenoidMap *map_s1= new PndSolenoidMap("SolenoidMap1", "R");
140  PndSolenoidMap *map_s2= new PndSolenoidMap("SolenoidMap2", "R");
141  PndSolenoidMap *map_s3= new PndSolenoidMap("SolenoidMap3", "R");
142  PndSolenoidMap *map_s4= new PndSolenoidMap("SolenoidMap4", "R");
143 
144  AddField(map_s1);
145  AddField(map_s2);
146  AddField(map_s3);
147  AddField(map_s4);
148  }
149  else if (Map=="SOLENOID_HALF") {
150 
151 
152  PndSolenoidMap *map_s1= new PndSolenoidMap("SolenoidMap_Half1", "R");
153  PndSolenoidMap *map_s2= new PndSolenoidMap("SolenoidMap_Half2", "R");
154  PndSolenoidMap *map_s3= new PndSolenoidMap("SolenoidMap_Half3", "R");
155  PndSolenoidMap *map_s4= new PndSolenoidMap("SolenoidMap_Half4", "R");
156 
157  AddField(map_s1);
158  AddField(map_s2);
159  AddField(map_s3);
160  AddField(map_s4);
161  }
162 
163 
164 
165 }
166 
167 
168 
169 
170 // ------------------------------------------------------------------------
171 
172 // ------------ Constructor from PndFieldPar --------------------------
174  : fMaps( new TObjArray(10)), fNoOfMaps(0), fFieldMaps(), fMapIter(), fBeamMom(0.)
175 {
176  fType = 5;
177  TObjArray *fArray= fieldPar->GetParArray();
178  if(fArray->IsEmpty()) fType=-1;
179 
180 }
181 
182 // ------------ Copy Constructor ----------------------------------------
183 PndMultiField::PndMultiField(const PndMultiField& field) : FairField(), fMaps(field.fMaps), fNoOfMaps(field.fNoOfMaps), fFieldMaps(field.fFieldMaps), fMapIter(field.fMapIter), fBeamMom(field.fBeamMom)
184 {
185  fType = field.fType;
186 }
187 
188 
189 // ------------ Destructor --------------------------------------------
191 
192  //printf("PndMultiField::~PndMultiField() \n");
193  fMaps->Delete();
194  delete fMaps;
195 
196 }
197 
198 // ----------- Adding fields ------------------------------------------
199 void PndMultiField::AddField(FairField *field){
200 
201  if(field){
202  fMaps->AddLast(field);
203  fNoOfMaps++;
204  }
205 
206 }
207 
208 // ----------- Intialisation ------------------------------------------
210  PndConstField *field=0;
211  PndFieldMap *fieldMap=0;
212  for (Int_t n=0; n<=fNoOfMaps; n++){
213  fieldMap = dynamic_cast<PndFieldMap *>(fMaps->At(n));
214  field = dynamic_cast<PndConstField *>(fMaps->At(n));
215  if(fieldMap){
216  fieldMap->Init();
217  fFieldMaps.insert( pair<PndRegion*, FairField*>(new PndRegion(fieldMap->GetZmin(),fieldMap->GetZmax()) , fieldMap ));
218  }else if(field){
219  field->Init();
220  fFieldMaps.insert( pair<PndRegion*, FairField*>(new PndRegion(field->GetZmin(),field->GetZmax() ), field ));
221  }
222  }
223 
224 }
225 
226 // --------- Screen output --------------------------------------------
228  for (Int_t n=0; n<=fNoOfMaps; n++){
229  FairField *fieldMap = dynamic_cast<FairField *>(fMaps->At(n));
230  if(fieldMap) fieldMap->Print("");
231  }
232 }
233 
234 
235 
236 
237 
238 // --------- Fill the parameters --------------------------------------------
240 // for (Int_t n=0; n<=fNoOfMaps; n++){
241 // FairField *fieldMap = dynamic_cast<FairField *>(fMaps->At(n));
242 // if(fieldMap) fieldMap->FillParContainer();
243 // }
244  FairRun *fRun=FairRun::Instance();
245  FairRuntimeDb *rtdb=fRun->GetRuntimeDb();
246  //Bool_t kParameterMerged=kTRUE;
247  PndMultiFieldPar* Par = (PndMultiFieldPar*) rtdb->getContainer("PndMultiFieldPar");
248  Par->SetParameters(this);
249  Par->setInputVersion(fRun->GetRunId(),1);
250  Par->setChanged();
251 
252 
253 
254 }
255 
256 // -------------------------------------------------------------------------
257 void PndMultiField::GetFieldValue(const Double_t point[3], Double_t* bField)
258 {
259 
260  PndRegion *fReg=0;
261  FairField *fField=0;
262  for (fMapIter=fFieldMaps.begin(); fMapIter!= fFieldMaps.end();fMapIter++ ){
263  fReg=fMapIter->first;
264  if(fReg->IsInside(point[2])){
265  fField=fMapIter->second;
266  break;
267  }
268  }
269  if(fField){
270  fField->GetBxyz(point, bField);
271  }else{
272  bField[0] = 0;
273  bField[1] = 0;
274  bField[2] = 0;
275 
276  }
277 }
278 
279 
virtual ~PndMultiField()
Double_t GetZmin() const
Definition: PndFieldMap.h:90
PndMultiField * fField
Definition: sim_emc_apd.C:97
TObjArray * fMaps
Definition: PndMultiField.h:61
virtual void Init()
void SetParameters(FairField *field)
std::map< PndRegion *, FairField * >::iterator fMapIter
Definition: PndMultiField.h:68
int n
virtual void Print()
Double_t fBeamMom
Definition: PndMultiField.h:69
TObjArray * GetParArray()
PndDipoleMap * map_d2
Definition: sim_hit_emc.C:111
PndSolenoidMap * map_s1
Definition: sim_hit_emc.C:112
Double_t GetZmax() const
Definition: PndConstField.h:91
FairRunAna * fRun
Definition: hit_dirc.C:58
Double_t GetZmin() const
Definition: PndConstField.h:90
Double_t
PndSolenoidMap * map_s3
Definition: sim_hit_emc.C:114
PndMultiFieldPar * fieldPar
Definition: sim_ftof.C:102
TClonesArray * point
Definition: anaLmdDigi.C:29
std::map< PndRegion *, FairField * > fFieldMaps
Definition: PndMultiField.h:67
FairRuntimeDb * rtdb
Definition: hit_dirc.C:66
PndSolenoidMap * map_s2
Definition: sim_hit_emc.C:113
Bool_t IsInside(Double_t Z)
Definition: PndRegion.cxx:20
void FillParContainer()
PndMultiFieldPar * Par
Definition: sim_emc_apd.C:115
ClassImp(PndAnaContFact)
void AddField(FairField *field)
PndDipoleMap * map_d1
Definition: sim_hit_emc.C:110
PndTransMap * map_t
Definition: sim_hit_emc.C:109
Double_t GetZmax() const
Definition: PndFieldMap.h:93
PndSolenoidMap * map_s4
Definition: sim_hit_emc.C:115