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