]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STRUCT/AliFieldReader.cxx
Introducing TGeoMatrix where needed
[u/mrichter/AliRoot.git] / STRUCT / AliFieldReader.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 #include <TCanvas.h>
19 #include <TFile.h>
20 #include <TGraph.h>
21 #include <TH1F.h>
22 #include <TMath.h>
23 #include <TMultiGraph.h>
24 #include <TNtuple.h>
25 #include <TObjString.h>
26 #include <TString.h>
27 #include <TStyle.h>
28
29 #include "AliFieldReader.h"
30 #include "AliMagFMaps.h"
31
32 ClassImp(AliFieldReader)
33
34
35
36 //_______________________________________________________________________
37 AliFieldReader::AliFieldReader():
38     fField(0),
39     fMap(0),
40     fCatalogue(0),
41     fHtmlMain(0),
42     fStepSize(0.),
43     fZStart(1383.),
44     fDd(0.08),
45     fDz(0.064),
46     fPolarity(1.),
47     fCatalogueName("goodfiles.list")
48 {
49 //
50 //  Constructor
51 //
52 }
53
54 AliFieldReader::AliFieldReader(const AliFieldReader& reader):
55     TObject(reader),
56     fField(0),
57     fMap(0),
58     fCatalogue(0),
59     fHtmlMain(0),
60     fStepSize(0.),
61     fZStart(0.),
62     fDd(0.),
63     fDz(0.),
64     fPolarity(0.),
65     fCatalogueName(0)
66 {
67 //
68 //  Constructor
69 //
70     reader.Copy(*this);
71 }
72
73 //_______________________________________________________________________
74 AliFieldReader::~AliFieldReader()
75 {
76   //
77   // Destructor
78   //
79 }
80
81 //_______________________________________________________________________
82 void AliFieldReader::Init()
83 {   
84 //
85 // Initialize the reader
86 //
87     // Calculated map
88     fField = new AliMagFMaps("Maps","Maps", 2, 1., 10., 2);
89     // Catalogue
90     fCatalogue = fopen(fCatalogueName, "r");
91     
92     // HTML
93     fHtmlMain = fopen("bmap.html", "w");
94     MakeHtmlHeaderMain(fHtmlMain);
95     
96 }
97
98
99 void AliFieldReader::ReadMap()
100 {   
101 //
102 // Read the measured dipole field map 
103 //    
104     
105     Float_t zA[450], bxzA[200], byzA[200], bzzA[200], bxzcA[200], byzcA[200], bzzcA[200];
106     Float_t yA[450], bxyA[200], byyA[200], bzyA[200], bxycA[200], byycA[200], bzycA[200];
107     Float_t xA[450], bxxA[200], byxA[200], bzxA[200], bxxcA[200], byxcA[200], bzxcA[200];
108     
109     Char_t sLine[255];
110     Char_t fLine[255];
111     
112     Float_t xpos, ypos, zpos;
113     Int_t iboxpos;
114     
115     Float_t h[3], x[3], b[3];
116     Float_t temp;
117     Float_t xnt[17];
118     Int_t ires;
119     
120     Int_t iret, ipos;
121     Int_t ic115 = 0;
122     Float_t  dx = 0.;
123     Float_t  dy = 0.;
124     
125     Init();
126     // Read register
127     ReadRegisterMap();
128     // n-tuple
129     fMap = new TNtuple("Field Map", "Map", 
130                        "x:y:z:ix:iy:iz:bx:by:bz:bxc:byc:bzc:ifile:ireg:irot:temp:cal", 4000);
131     
132     //
133     // Loop over files
134     // 
135     Int_t ifiles = 0;
136     // File catalogue
137     while ((fgets(fLine, 255, fCatalogue)) != NULL && ifiles <= 2000) {
138
139         if (strncmp(fLine,"#HEADER", 7) == 0) {
140           iret = sscanf(&fLine[7],"%f", &fZStart);
141           continue;
142         }       
143         ifiles++;
144         
145 //      if (ifiles != 87) continue;
146         char fileName[32];
147         iret = sscanf(fLine, "%s", fileName);
148         printf("Reading File %s\n", fileName);
149         
150 //  Get run name
151         TString* tsFile          = new TString(fileName);
152         TObjArray * tokens       =  tsFile->Tokenize(".");
153         const char* runName      = (((TObjString*) tokens->At(0))->GetString()).Data();
154         FILE* file               =  fopen(fileName, "r");
155         
156         
157         Float_t bdl   = 0.;
158         Float_t bdlc  = 0.;
159         Int_t iz      = 0;
160         Int_t izz     = 0;
161         Int_t iy      = 0;
162         Int_t ix      = 0;
163         Int_t iheader = 0;
164         Float_t stepsz = 0.;
165     // Graphs go here
166         TMultiGraph* bxzmg  = new TMultiGraph("bxzmg", "B_{x}");
167         TMultiGraph* byzmg  = new TMultiGraph("byzmg", "B_{y}");
168         TMultiGraph* bzzmg  = new TMultiGraph("bzzmg", "B_{z}");
169         while ((fgets(sLine, 255, file)) != NULL) {
170 // read z-position
171             
172             if (strncmp(sLine," Current", 8) == 0) {
173                 Float_t current;
174                 iret = sscanf(&sLine[9],"%f", &current);
175                 printf("Current %f \n", current);
176                 
177                 fPolarity = current / 6000.;
178             }
179             if (strncmp(sLine," z", 2) == 0) {
180                 Float_t zmin, zmax;
181                 Int_t nsteps;
182                 iret = sscanf(&sLine[3],"%f %f %d", &zmin, &zmax, &nsteps);
183                 printf("zmin zmax %13.3f %13.3f %13.3f\n", 
184                        zmin, zmax, TMath::Abs(zmax - zmin)/Float_t(nsteps));
185                 zmax = TMath::Max(zmin, zmax);
186                 stepsz  = TMath::Abs(zmax - zmin)/Float_t(nsteps);
187             }
188             
189             if (strncmp(sLine," X position", 11) == 0)
190             {
191                 TString string;
192                 TString* tsLine = new TString(sLine);
193                 TObjArray * tokens =  tsLine->Tokenize("=");    
194                 string = ((TObjString*) tokens->At(1))->GetString();
195                 iret = sscanf(string.Data(), "%f", &xpos);
196                 string = ((TObjString*) tokens->At(2))->GetString();
197                 iret = sscanf(string.Data(), "%f", &ypos);
198                 string = ((TObjString*) tokens->At(3))->GetString();
199                 iret = sscanf(string.Data(), "%d", &iboxpos);
200                 
201                 
202                 printf("This file is for  x = %13.3f y = %13.3f Box Position %1d\n", xpos, ypos, iboxpos);
203                 iheader++;
204                 
205             }
206             
207             if (strncmp(sLine," R       Z-pos",8) == 0)
208             {
209                 iret = sscanf(&sLine[8],"%e", &zpos);
210                 ipos  = 0;
211                 ic115 = 0;
212                 izz++;
213             }
214             
215             
216             if (strncmp(sLine,"C",1) == 0) {
217                 ipos++;
218                 Int_t ireg;
219                 iret = sscanf(&sLine[2],"%d %f %f %f %f %d", &ireg, &h[0], &h[1], &h[2], &temp, &ires);
220                 //
221                 // fix for address 115
222                 //
223                 if (ireg == 115) ic115++;
224                 if (ic115 == 2 && ireg == 115) ireg = 119;
225                 
226                 if (iheader == 1) {
227
228                     Float_t bx = 0., by = 0., bz = 0.;
229                     Int_t jx = fRegMap[ireg][0];
230                     Int_t jy = fRegMap[ireg][1];
231                     Int_t jz = fRegMap[ireg][2];
232
233                     switch (iboxpos) {
234                     case 0:
235                         bx =  h[1];
236                         by =  h[2];
237                         dx = -0.36 + jx * fDd;
238                         dy =  jy * fDd;
239                         break;
240                     case 1:
241                         bx = -h[2];
242                         by =  h[1];
243                         dx =  0.36 + jy * fDd;
244                         dy = -(-0.36 + jx * fDd);
245                         break;
246                     case 2:
247                         bx = -h[1];
248                         by = -h[2];
249                         dx =  0.36 - jx * fDd;
250                         dy = -jy * fDd;
251                         break;
252                     case 3:
253                         bx =   h[2];
254                         by =  -h[1];
255                         dx =  -0.36 - jy * fDd;
256                         dy =   -(0.36 - jx * fDd);
257                         break;
258                     } // switch
259
260                     bz =  h[0];
261                     if (jz == 0) {
262                         bz = -bz;
263                         if (iboxpos == 1 || iboxpos == 3) {
264                             by = -by;
265                         } else {
266                             bx = -bx;
267                         }
268                     }
269                     
270                     
271                     Float_t dz = (jz == 0)?  fDz/2. : - fDz/2.;
272                     dz *= 100.;
273                     
274
275                     Float_t xc = xpos + dx;
276                     Float_t yc = ypos + dy;
277
278                     
279                     x[0] =  - xc * 100.;
280                     x[1] =  + yc * 100.;
281                     x[2] =  - (-zpos * 100. + fZStart + dz);
282
283                     fField->Field(x, b);
284                     b[0] *= fPolarity;
285                     b[1] *= fPolarity;
286                     b[2] *= fPolarity;
287                     
288                     xnt[ 0] = xc;
289                     xnt[ 1] = yc;
290                     xnt[ 2] = -x[2] / 100.;
291                     xnt[ 3] = Float_t (jx);
292                     xnt[ 4] = Float_t (jy);
293                     xnt[ 5] = Float_t (jz);
294                     xnt[ 6] = bx;
295                     xnt[ 7] = by;
296                     xnt[ 8] = bz;
297                     xnt[ 9] = b[0]/10.;
298                     xnt[10] = b[1]/10.;
299                     xnt[11] = b[2]/10.;             
300                     xnt[12] = Float_t (ifiles);
301                     xnt[13] = Float_t (ireg);
302                     xnt[14] = Float_t(iboxpos);
303                     xnt[15] = temp;
304                     xnt[16] = Float_t(ires);
305                     
306                     fMap->Fill(xnt);
307                     
308
309 //
310 // Calculated field             
311
312                     if (jy != -1 && jz == 1 && jy == 0 && izz == 40){
313                         x[1] = x[1] + 10.;
314                         fField->Field(x, b);
315                         yA  [jx]  = yc;
316                         bxyA[jx]  = bx;
317                         byyA[jx]  = by;
318                         bzyA[jx]  = bz;
319                         bxycA[jx] =   b[0] / 10.;
320                         byycA[jx] =   b[1] / 10.;
321                         bzycA[jx] =   b[2] / 10.;
322                         iy++;
323                     } // if
324                     
325                     if (jy != -1 && jz == 1 && jy == 0 && izz == 1){
326                         x[1] = x[1] + 10.;
327                         fField->Field(x, b);
328                         xA  [jx]  = xc;
329                         bxxA[jx]  = bx;
330                         byxA[jx]  = by;
331                         bzxA[jx]  = bz;
332                         bxxcA[jx] =   b[0] / 10.;
333                         byxcA[jx] =   b[1] / 10.;
334                         bzxcA[jx] =   b[2] / 10.;
335                         ix++;
336                     } // if
337
338
339                     if (ireg == 181) {
340                         fField->Field(x, b);
341 //                      printf("Field %f %f %f %f %f %f \n", x[0], x[1], x[2], b[0], b[1] , b[2]);
342                         
343                         bdl  +=  stepsz * bx;
344                         bdlc +=  stepsz * b[0];
345                         
346                         zA [iz]  = x[2];
347                         bxzA[iz]  = bx;
348                         byzA[iz]  = by;
349                         bzzA[iz]  = bz;
350                         bxzcA[iz] =   b[0] / 10.;
351                         byzcA[iz] =   b[1] / 10.;
352                         bzzcA[iz] =   b[2] / 10.;
353                         iz++;
354                     }
355                 } // if 1st header 
356             } // if C
357         } // next line
358         gStyle->SetOptStat(0);
359         char title[128];
360         sprintf(title, "File#: %5d, X = %13.2f m , Y = %13.2f m, Box Orientation: %2d", ifiles, xpos, ypos, iboxpos);
361         TCanvas* c2 = new TCanvas("c2", title, 1200, 800);
362         c2->Divide(2,2);
363         c2->cd(1);
364         TGraph* bxg =  new TGraph(iz, zA, bxzA);
365         TGraph* bxcg = new TGraph(iz, zA, bxzcA);    
366         bxcg->SetLineColor(2);
367         bxzmg->Add(bxg);
368         bxzmg->Add(bxcg);
369         bxzmg->Draw("lA");
370         bxzmg->GetHistogram()->SetXTitle("z[m]");
371         bxzmg->GetHistogram()->SetYTitle("B_{x} [T]");
372         bxzmg->Draw("lA");
373         
374         c2->cd(2);
375         TGraph* byg =  new TGraph(iz, zA, byzA);
376         TGraph* bycg = new TGraph(iz, zA, byzcA);    
377         bycg->SetLineColor(2);
378         byzmg->Add(byg);
379         byzmg->Add(bycg);
380         byzmg->Draw("lA");
381         byzmg->GetHistogram()->SetXTitle("z[m]");
382         byzmg->GetHistogram()->SetYTitle("B_{y} [T]");
383         byzmg->Draw("lA");
384         
385         c2->cd(3);
386         TGraph* bzg =  new TGraph(iz, zA, bzzA);
387         TGraph* bzcg = new TGraph(iz, zA, bzzcA);    
388         bzcg->SetLineColor(2);
389         bzzmg->Add(bzg);
390         bzzmg->Add(bzcg);
391         bzzmg->Draw("lA");
392         bzzmg->GetHistogram()->SetXTitle("z[m]");
393         bzzmg->GetHistogram()->SetYTitle("B_{z} [T]");
394         bzzmg->Draw("lA");
395         
396         c2->Update();
397         char pictFile[64];
398         sprintf(pictFile, "%s.gif", runName);
399         c2->SaveAs(pictFile);
400         //
401         // Html generation
402         //
403         char htmlFile[64];
404         sprintf(htmlFile, "%s.html", runName);
405         FILE* chtml = fopen(htmlFile, "w");
406         MakeHtmlHeaderPict(chtml);
407         MakeHtmlPict(chtml, pictFile);
408         MakeHtmlTableEntry(fHtmlMain, fileName, htmlFile, xpos, ypos, iboxpos, bdl, ifiles);
409
410         //
411         //
412         printf("Bdl [Tm] %f %f \n", 2. * bdl, 2 * bdlc / 10.);
413     } // files
414     MakeHtmlTrailor(fHtmlMain);
415     TFile* out = new TFile("fmap.root", "recreate");
416     fMap->Write();
417     out->Close();
418 }
419
420 void AliFieldReader::ReadMapSolenoid(){
421 //
422 //  Read map for solenoid measurement
423 // 
424     Float_t phiA[450], bzPhiA[200], brPhiA[200], btPhiA[200], bbPhiA[200];
425     Float_t bzcPhiA[200], brcPhiA[200], btcPhiA[200], bbcPhiA[200];    
426     Char_t sLine[255];
427     Char_t fLine[255];
428     
429     Float_t zpos, phipos, skewing, temp;
430     Int_t ical;
431     Float_t zmin, zmax;
432     Float_t h[3], x[3], b[3];
433     Int_t iret;
434     Int_t ipos = 0;
435     
436     Init();
437     ReadRegisterMapSolenoid();
438     fMap = new TNtuple("Field Map", "Map", 
439                        "r:phi:z:br:bt:bz:brc:btc:bzc:ifile:ireg:temp:cal:arm", 4000);
440     
441     //
442     // Loop over files
443     // 
444     Int_t ifiles = 0;
445     // File catalogue
446     while ((fgets(fLine, 255, fCatalogue)) != NULL && ifiles <= 2000) {
447
448         if (strncmp(fLine,"#HEADER", 7) == 0) {
449           iret = sscanf(&fLine[7],"%f", &fZStart);
450           continue;
451         }       
452         ifiles++;
453         
454         char fileName[32];
455         iret = sscanf(fLine, "%s", fileName);
456         printf("Reading File %s\n", fileName);
457         
458 //  Get run name
459         TString* tsFile          = new TString(fileName);
460         TObjArray * tokens       =  tsFile->Tokenize(".");
461         Int_t n = tokens->GetEntries();
462         char* runName = new  char[256]; 
463         sprintf(runName, "%s", (((TObjString*) tokens->At(0))->GetString()).Data());
464         if (n > 2) {
465             for (Int_t i = 1; i < n-1; i++)
466             {
467                 sprintf(runName, "%s.%s", 
468                         runName, (((TObjString*) tokens->At(i))->GetString()).Data());
469             }
470         }
471         FILE* file               =  fopen(fileName, "r");
472         
473         
474         Float_t bdl  = 0.;
475         Float_t bdlc = 0.;
476         Int_t iA   = 0;
477         Int_t izz  = 0;
478         Int_t iphi = 0;
479     // Graphs go here
480         TMultiGraph* bxzmg  = new TMultiGraph("bxzmg", "B_{z}");
481         TMultiGraph* byzmg  = new TMultiGraph("byzmg", "B_{r}");
482         TMultiGraph* bzzmg  = new TMultiGraph("bzzmg", "B_{t}");
483         TMultiGraph* bbmg   = new TMultiGraph("bbmg", "|B|");
484
485         while ((fgets(sLine, 255, file)) != NULL) {
486             if (strncmp(sLine," z", 2) == 0) {
487                 Int_t nsteps;
488                 iret = sscanf(&sLine[3],"%f %f %d", &zmin, &zmax, &nsteps);
489                 printf("zmin zmax %13.3f %13.3f %13.3f\n", 
490                        zmin, zmax, TMath::Abs(zmax - zmin)/Float_t(nsteps));
491             }
492             if (strncmp(sLine," R\tPOSITION NUMBER", 18) == 0)
493             {                 
494                 //
495                 // Current z-position
496                 TString string;
497                 TString* tsLine = new TString(sLine);
498                 TObjArray * tokens =  tsLine->Tokenize("=");    
499                 string = ((TObjString*) tokens->At(1))->GetString();
500                 iret = sscanf(string.Data(), "%f", &zpos);
501                 printf("POSITION NUMBER Z: %f\n", zpos);
502                 izz ++;
503                 delete tsLine;
504             }
505
506
507             if (strncmp(sLine," SKEWING ON Z:", 14) == 0)
508             {                 
509                 //
510                 // Skew in z
511               iret = sscanf(&sLine[14],"%f", &skewing);
512               printf("SKEWING ON Z: %f\n", skewing);
513             }
514
515
516             if (strncmp(sLine,"Phi", 3) == 0)
517             {
518               Float_t phiStart, phiStop;
519               
520               iret = sscanf(&sLine[3],"%f %f", &phiStart, &phiStop);
521
522               printf("phiStart phiStop %f %f\n", phiStart, phiStop);
523
524             }
525             
526
527             if (strncmp(sLine," R\tPhi-Angle",12) == 0)
528             {
529                 iret = sscanf(&sLine[12],"%e", &phipos);
530                 ipos = 0;
531                 iphi++;
532             }
533             
534
535             if (strncmp(sLine,"C",1) == 0) {
536         
537                 Int_t ireg;
538                 iret = sscanf(&sLine[2],"%d %f %f %f %f %d", &ireg, &h[0], &h[1], &h[2], &temp, &ical);
539                 ipos++;
540
541                 Int_t ir = fRegMap[ireg][0];
542                 Int_t ia = fRegMap[ireg][1];
543                 Float_t rpos = 0.;
544
545                 
546                 if (ia == 0) {
547                     rpos    = 0.2295 + ir * 0.16;
548                 } else {
549                     if (ireg ==  81) rpos = 0.2295;
550                     if (ireg ==  59) rpos = 1.0295;
551                     if (ireg == 142) rpos = 2.1495;                 
552                     if (ireg == 180) rpos = 3.1095;                 
553                     if (ireg ==  69) rpos = 4.2295;                 
554
555                     //              if (ireg ==  55) rpos = 0.2295;
556                     //if (ireg == 195) rpos = 1.0295;
557                     //if (ireg == 129) rpos = 2.1495;               
558                     //if (ireg == 167) rpos = 3.1095;               
559                     //if (ireg == 142) rpos = 4.2295;               
560                 }
561                 
562                 Float_t phi = phipos;
563                 if (ia == 1) {
564                     phi += 180.;
565                     if (phi > 360.) phi -= 360.;                            
566                 }
567                 
568                 phi = - phi * TMath::Pi() / 180.;
569                 Float_t xpos = rpos * TMath::Cos(phi);
570                 Float_t ypos = rpos * TMath::Sin(phi);
571                 x[0] = - xpos * 100.;
572                 x[1] = ypos * 100.;
573                 x[2] = -400. + zpos * 100.;
574                 
575                 fField->Field(x, b);
576                 Float_t phi0 = TMath::Pi() - phi;
577                 
578                 Float_t brc =   b[0] * TMath::Cos(phi0) +  b[1] * TMath::Sin(phi0);
579                 Float_t btc = - b[0] * TMath::Sin(phi0) +  b[1] * TMath::Cos(phi0);     
580                 Float_t bzc =   b[2];
581                 
582                 brc /= 10.;
583                 btc /= 10.;
584                 bzc /= 10.;
585         
586                 fMap->Fill(rpos, -phi, -(615.5 - fZStart) / 100. + zpos, h[2], -h[1], h[0], brc, btc, bzc, ifiles, ireg, temp, Float_t(ical), Int_t(ia));
587         
588                 if (ireg  == 174) {
589                     printf("Field (Bx, By, Bz) at position %d: %13.3f %13.3f %13.3f %13.3f %13.3f %13.3f %5d\n", 
590                            ipos, zpos, rpos, phipos, h[0], h[1], h[2], ia);     
591                     if (izz == 1) {
592                         phiA[iA]   = phi;
593                         bzPhiA[iA] = -h[0];
594                         brPhiA[iA] =  h[2];
595                         btPhiA[iA] = -h[1];
596                         bbPhiA[iA] = TMath::Sqrt(h[0] * h[0] + h[1] * h[1] + h[2] * h[2]);
597                         bzcPhiA[iA] = bzc;
598                         brcPhiA[iA] = brc;
599                         btcPhiA[iA] = btc;
600                         bbcPhiA[iA] = TMath::Sqrt(brc * brc + bzc * bzc + btc * btc);
601                         iA++;
602                     }
603                 }
604             } // if R
605         } // next line
606
607         gStyle->SetOptStat(0);
608         char title[128];
609         sprintf(title, "Z = %13.2f m , Phi = %13.2f m", zpos, phipos);
610
611
612         TCanvas* c2 = new TCanvas("c2", title, 1200, 800);
613                 
614         c2->Divide(2,2);
615         c2->cd(1);
616         TGraph* bzg  =  new TGraph(iA, phiA, bzPhiA);
617         TGraph* bzcg =  new TGraph(iA, phiA, bzcPhiA);
618         bzcg->SetLineColor(2);
619         bxzmg->Add(bzg);
620         bxzmg->Add(bzcg);
621         bxzmg->Draw("lA");
622         bxzmg->GetHistogram()->SetXTitle("#phi[rad]");
623         bxzmg->GetHistogram()->SetYTitle("B_{z} [T]");
624         bxzmg->Draw("lA");
625
626         c2->cd(2);
627         TGraph* brg  =  new TGraph(iA, phiA, brPhiA);
628         TGraph* brcg =  new TGraph(iA, phiA, brcPhiA);
629         brcg->SetLineColor(2);
630         byzmg->Add(brcg);
631         byzmg->Add(brg);
632
633         byzmg->SetMaximum(0.03);
634         byzmg->SetMinimum(-0.03);       
635         byzmg->Draw("lA");
636         byzmg->GetHistogram()->SetXTitle("#phi[rad]");
637         byzmg->GetHistogram()->SetYTitle("B_{r} [T]");
638         byzmg->Draw("lA");
639
640         c2->cd(3);
641         TGraph* btg  =  new TGraph(iA, phiA, btPhiA);
642         TGraph* btcg =  new TGraph(iA, phiA, btcPhiA);
643         btcg->SetLineColor(2);
644         bzzmg->Add(btcg);
645         bzzmg->Add(btg);
646         bzzmg->Draw("lA");
647         bzzmg->SetMaximum(0.03);
648         bzzmg->SetMinimum(-0.03);       
649         bzzmg->GetHistogram()->SetXTitle("#phi[rad]");
650         bzzmg->GetHistogram()->SetYTitle("B_{t} [T]");
651         bzzmg->Draw("lA");
652         
653
654         c2->cd(4);
655         TGraph* bg  =  new TGraph(iA, phiA, bbPhiA);
656         TGraph* bcg =  new TGraph(iA, phiA, bbcPhiA);
657         bcg->SetLineColor(2);
658         bbmg->Add(bg);
659         bbmg->Add(bcg);
660         bbmg->Draw("lA");
661         bbmg->GetHistogram()->SetXTitle("#phi[rad]");
662         bbmg->GetHistogram()->SetYTitle("|B| [T]");
663         bbmg->Draw("lA");
664         
665
666
667
668         char pictFile[64];
669         sprintf(pictFile, "%s.gif", runName);
670         c2->SaveAs(pictFile);
671
672         //
673         // Html generation
674         //
675         char htmlFile[64];
676         sprintf(htmlFile, "%s.html", runName);
677         FILE* chtml = fopen(htmlFile, "w");
678         MakeHtmlHeaderPict(chtml);
679         MakeHtmlPict(chtml, pictFile);
680         MakeHtmlTableEntry(fHtmlMain, fileName, htmlFile, zmin, zmax, ifiles, 0., 0);
681
682         //
683         //
684         printf("Bdl [Tm] %f %f \n", 2. * bdl, 2 * bdlc / 10.);
685     } // files
686     MakeHtmlTrailor(fHtmlMain);
687     TFile* out = new TFile("fmap.root", "recreate");
688     fMap->Write();
689     out->Close();
690 }
691
692
693
694 void AliFieldReader::MakeHtmlHeaderMain(FILE* file)
695 {
696 //
697 //  Write the header of the heml output
698 //
699     fprintf(file,"<!DOCTYPE html PUBLIC \"-//W3C//DTD HTML 4.01 Transitional//EN\">\n");
700     fprintf(file, "<html>\n");
701     fprintf(file, "<head>\n");
702     fprintf(file, "<meta http-equiv=\"content-type\"\n");
703     fprintf(file, "content=\"text/html; charset=ISO-8859-1\"\n");
704     fprintf(file, "<title>main.html</title>\n");
705     fprintf(file, "<\\head>\n");
706     fprintf(file, "<body>\n");
707     fprintf(file, "<table cellpadding=\"1\" cellspacing=\"1\" border=\"1\"\n");
708     fprintf(file, "style=\"text-align: left; width: 80;\">\n");
709     fprintf(file, "<tbody>\n");
710     fprintf(file, "<td style=\"vertical-align: top;\">File#  <br></td>\n");
711     fprintf(file, "<td style=\"vertical-align: top;\">File Name  <br></td>\n");
712     fprintf(file, "<td style=\"vertical-align: top;\">X-Position <br></td>\n");
713     fprintf(file, "<td style=\"vertical-align: top;\">Y-Position <br></td>\n");
714     fprintf(file, "<td style=\"vertical-align: top;\">Box Orientation <br></td>\n");
715     fprintf(file, "<td style=\"vertical-align: top;\">B.dl <br></td>\n");
716     fprintf(file, "<td style=\"vertical-align: top;\">Link to plots <br></td>\n");
717
718 }
719
720 void  AliFieldReader::MakeHtmlHeaderPict(FILE* file)
721 {
722 //
723 //  Write header for picture
724 //
725     fprintf(file,"<!DOCTYPE html PUBLIC \"-//W3C//DTD HTML 4.01 Transitional//EN\">\n");
726     fprintf(file, "<html>\n");
727     fprintf(file, "<head>\n");
728     fprintf(file, "<meta http-equiv=\"content-type\"\n");
729     fprintf(file, "content=\"text/html; charset=ISO-8859-1\"\n");
730     fprintf(file, "<title>main.html</title>\n");
731     fprintf(file, "</head>\n");
732     fprintf(file, "<body>\n");
733 }
734
735 void AliFieldReader:: MakeHtmlPict(FILE* chtml, char* pictFile)
736 {
737 //
738 //  Write html for including picture
739 //
740     fprintf(chtml, "<img src=\"./%s\" alt=\"%s\" style=\"width: 1196px; height: 772px;\">\n", 
741             pictFile, pictFile);
742     
743     fprintf(chtml, "<br> <br> <br>\n"); 
744     fprintf(chtml, "<a href=\"./bmap.html\">Back to main page</a><br>\n");
745     
746     fprintf(chtml, "</body>\n");
747     fprintf(chtml, "</header>\n");      
748     fclose(chtml);
749 }
750
751 void AliFieldReader::MakeHtmlTableEntry(FILE* htmlmain, char* fileName, char* htmlFile, Float_t x, Float_t y, Int_t i, Float_t bdl, Int_t ifile)
752 {       
753     fprintf(htmlmain, "<tr>\n");
754 //
755     fprintf(htmlmain, "<td style=\"vertical-align: top;\">%5d</td>\n", ifile);
756     fprintf(htmlmain, "<td style=\"vertical-align: top;\">%s</td>\n", fileName);
757
758 //
759     fprintf(htmlmain, "<td style=\"vertical-align: top;\">%13.2f</td>\n",x);
760     fprintf(htmlmain, "<td style=\"vertical-align: top;\">%13.2f</td>\n",y);
761     fprintf(htmlmain, "<td style=\"vertical-align: top;\">%3d   </td>\n",i);
762     fprintf(htmlmain, "<td style=\"vertical-align: top;\">%13.3f</td>\n",bdl);
763 //
764     fprintf(htmlmain, "<td style=\"vertical-align: top;\">\n");
765     fprintf(htmlmain, "<span style=\"text-decoration: underline;\"><a href=\"./%s\">gif</a></span></td>\n", htmlFile);
766     fprintf(htmlmain, "</tr>\n");
767 }
768
769
770 void  AliFieldReader::MakeHtmlTrailor(FILE* htmlmain)
771 {
772 //
773 //  Write the html trailor
774 //
775     fprintf(htmlmain, "</tbody>\n");
776     fprintf(htmlmain, "</table>\n");    
777     fprintf(htmlmain, "</body>\n");
778     fprintf(htmlmain, "</header>\n");   
779     fclose(htmlmain);
780 }
781
782 void AliFieldReader::ReadRegisterMap()
783 {
784 //
785 //  Read the register map
786 //
787     FILE* regmap = fopen("register.map", "r");
788     Int_t ireg;
789     for (ireg = 0; ireg < 200; ireg++) {
790         fRegMap[ireg][0] = -1;
791         fRegMap[ireg][1] = -1;
792         fRegMap[ireg][2] = -1;
793     }
794     for (Int_t iz = 0; iz < 2; iz++) {
795         for (Int_t iy = 2; iy >= 0; iy--) {
796             for (Int_t ix = 0; ix < 10; ix++) {
797                 fscanf(regmap, "%d\n", &ireg);
798                 printf("Address %5d %5d %5d %5d \n", iz, iy, ix, ireg);
799                     fRegMap[ireg][1] = iy;
800                     fRegMap[ireg][2] = iz;
801                 if (iz == 1) {
802                     fRegMap[ireg][0] = ix;
803                 } else {
804                     fRegMap[ireg][0] = 9 - ix;
805                 }
806             } // ix
807         } // iy
808     } // iz
809     fclose(regmap);
810     printf("-> ReadRegisterMap()\n\n");
811 }
812
813
814 void AliFieldReader::ReadRegisterMapSolenoid()
815 {
816 //
817 //  Read the register map
818 //
819     FILE* regmap = fopen("register.map", "r");
820     Int_t ireg;
821     
822 // Initialize
823     for (ireg = 0; ireg < 200; ireg++) {
824         fRegMap[ireg][0] = -1;
825         fRegMap[ireg][1] = -1;
826         fRegMap[ireg][2] = -1;
827     }
828
829 // Main arm 
830     for (Int_t ir = 0; ir < 33; ir++) {
831         fscanf(regmap, "%d\n", &ireg);
832         fRegMap[ireg][0] = ir;
833         fRegMap[ireg][1] = 0;   
834     }
835 // Opposite arm 
836     for (Int_t ir = 0; ir < 5; ir++) {
837         fscanf(regmap, "%d\n", &ireg);
838         fRegMap[ireg][0] = ir;
839         fRegMap[ireg][1] = 1;   
840     }
841     
842     fclose(regmap);
843     
844 }
845
846
847 AliFieldReader& AliFieldReader::operator=(const  AliFieldReader& rhs)
848 {
849 // Assignment operator
850     rhs.Copy(*this);
851     return *this;
852 }
853
854
855 void AliFieldReader::Copy( TObject&) const
856 {
857     //
858     // Copy 
859     //
860     Fatal("Copy","Not implemented!\n");
861 }