HYDRA_development_version
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
hmdcmille.cc
Go to the documentation of this file.
1 //*--- Author : V.Pechenov
2 //*--- Modified: 26.06.07 V.Pechenov
3 //*--- $Id: hmdcmille.cc,v 1.61 2007/07/09 12:38:36 pechenov Exp $
4 
5 using namespace std;
6 #include "hmdcmille.h"
7 #include "Mille.h"
8 #include "hades.h"
9 #include "hruntimedb.h"
10 #include "hcategory.h"
11 #include "hiterator.h"
12 #include "hmdcsizescells.h"
13 #include "hmdcclussim.h"
14 #include "hmdctrackdset.h"
15 #include "hmdctrackparam.h"
16 #include "hmdcwiredata.h"
17 #include "hgeantkine.h"
18 #include "hgeantmdc.h"
19 #include "hmdcgetcontainers.h"
20 #include "hmdclayergeompar.h"
21 #include "hmdclayercorrpar.h"
22 #include "hmdcgeompar.h"
23 #include "hspecgeompar.h"
24 #include "hmdcgeomstruct.h"
25 #include "hparasciifileio.h"
26 #include "hmdcdetector.h"
27 #include "hmdclistcells.h"
28 #include <TH1.h>
29 #include <TH2.h>
30 #include <TSystem.h>
31 
32 //_HADES_CLASS_DESCRIPTION
33 //////////////////////////////////////////////////////////////////////////////
34 //
35 // HMdcMille
36 //
37 // Interface class for Millipede alignment.
38 //
39 //////////////////////////////////////////////////////////////////////////////
40 
41 
43 
44 HMdcMille::HMdcMille() : HMdc12Fit("mdc.Mille","mdc.Mille") {
45  setDef();
46 }
47 
48 HMdcMille::HMdcMille(const Char_t* milleOutFName,const Char_t* geomPFName,
49  Int_t bmTmID,const Char_t* sumOfShFName) : HMdc12Fit("mdc.Mille","mdc.Mille") {
50 
51  setDef();
52  if(bmTmID<0 || bmTmID>9) printErrorAndExit("HMdcMille","beamTimeId must be >= 0 and <=9 !!!");
53  beamTimeId = bmTmID;
54 
55 //! if(milleOutFName) mille = new Mille(milleOutFName);
56 //! else printErrorAndExit("HMdcMille",
57 //! "Mille output file name is not specifaed !!!");
58  milleFileName = milleOutFName;
60 
61  if(sumOfShFName != 0) sumOfShiftsFName = sumOfShFName;
62  if(geomPFName != 0) geomParFileName = geomPFName;
63 
65  readPedeResFile(pedeResFName + ".res");
66  if(!loadGeometryPar()) printErrorAndExit("HMdcMille","Can't load mdc geometry!!! Stop!");
67  if(!makeShifts()) printErrorAndExit("HMdcMille","Can't change mdc geometry!");
68 
69  char buf[50];
70  sprintf(buf,".st%02i",iteration);
71  stepIter = buf;
72 }
73 
74 void HMdcMille::setDef(void) {
76  nParMax = 9;
77  iteration = 0;
78  mille = 0;
79  shiftType = 1;
80  doCopyGeomFile = kFALSE;
81  doCopySumShFile = kTRUE;
82  doCopyResFile = kTRUE;
83  doCopyLogFile = kFALSE;
84  doCopyHisFile = kFALSE;
85 
86  createPedeInFile = kFALSE; // kTRUE - create par.file
87  isGeomChanged = kFALSE;
88  isGeomFileExist = kFALSE;
89  nTracks = 0;
90  nWiresTot = 0;
91  nLayersTot = 0;
92  parSteps[0] = 0.01; // shift along X
93  parSteps[1] = 0.01; // shift along Y
94  parSteps[2] = 0.01; // shift along Z
95  parSteps[3] = 0.003; // rotation around X (in degree)
96  parSteps[4] = 0.003; // rotation around Y (in degree)
97  parSteps[5] = 0.003; // rotation around Z (in degree)
98  parSteps[6] = 0.004; // cell thickness
99  cellThicknFree = kFALSE;
100  nHists = 0;
101  fixFitTOffset = kTRUE;
102  useMdcShParOnly = kFALSE;
103  doLayP2Align = kTRUE;
104  nWiresMin = 15;
105  thetaCut = -1.;
106 
109 
110  nBinaryFiles = 0;
111  milleFileSize = 0;
112  crPedeTaskFile = kFALSE;
113 
114  printDebugFlag = kFALSE;
115  doHists = kTRUE;
116 
117  nClustersCut = 5; //???
118  filteringFlag = 0;
119  for(Int_t s=0;s<6;s++) {
120  useSector[s] = kTRUE;
121  nWiresCut[s] = 320;
122  }
123 
124  // Clean arrais:
125  for(Int_t im=0;im<24;im++) pSCMod[im] = NULL;
126  for(Int_t il=0;il<144;il++) pSCLayer[il] = NULL;
127  for(Int_t il=0;il<144;il++) nLayerParts[il] = 0;
128  for(Int_t s=0;s<6;s++) for(Int_t m=0;m<4;m++) {
129  isMdcInAlign[s][m] = kFALSE;
130  histInd[s] = -1;
131  mdcMods[s][m] = 0;
132  shiftSysMdc[s][m] = -1;
133  for(Int_t p=0;p<nParMax;p++) { // Mdc shifts
134  mShFlag[s][m][p] = mFixFlag[s][m][p] = kFALSE;
135  shiftsMdc[s][m][p] = 0.;
136  sigmaMdc[s][m][p] = 0.;
137  }
138  for(Int_t l=0;l<6;l++) { // Layer shifts
139  shiftSysLay[s][m][l] = -1;
140  for(Int_t p=0;p<nParMax;p++) {
141  lShFlag[s][m][l][p] = kFALSE;
142  lFixFlag[s][m][l][p] = kFALSE;
143  shiftsLay[s][m][l][p] = 0.;
144  sigmaLay[s][m][l][p] = 0.;
145  }
146  }
147  }
148 
149  shitsInfo[0] = "shift along X axis [mm]";
150  shitsInfo[1] = "shift along Y axis [mm]";
151  shitsInfo[2] = "shift along Z axis [mm]";
152  shitsInfo[3] = "rotation around X axis [deg]";
153  shitsInfo[4] = "rotation around Y axis [deg]";
154  shitsInfo[5] = "rotation around Z axis [deg]";
155  shitsInfo[6] = "cell thickness";
156  shitsInfo[7] = "shift along Y axis [mm]";
157  shitsInfo[8] = "rotation around Z axis [deg]";
158 
159  // Default file names:
160  pedeResFName = "millepede";
161  sumOfShiftsFName = "sumShifts.txt";
162  geomParFileName = "geoParameters.txt";
163 }
164 
166  pGetCont->deleteCont();
167 }
168 
169 Bool_t HMdcMille::init(void) {
170  return HMdc12Fit::init();
171 }
172 
173 Bool_t HMdcMille::reinit(void) {
174  if( !HMdc12Fit::reinit() ) return kFALSE;
176  if(pSizesCells == 0) return kFALSE;
177 
179 
180  for(Int_t sec=0;sec<6;sec++) {
181  nMods[sec] = 0;
182  HMdcSizesCellsSec& pSizesCellsSec = (*pSizesCells)[sec];
183  if(&pSizesCellsSec == 0) continue;
184  if(doHists) createHists(sec);
185  for(module=0;module<4;module++) {
186  HMdcSizesCellsMod& pSizesCellsMod = pSizesCellsSec[module];
187  if(&pSizesCellsMod == 0) continue;
188  nMods[sec]++;
189  setIMod(sec,module); // calculate iMod
190  pSCMod[iMod] = &pSizesCellsMod;
191 
192  for(layer=0;layer<6;layer++) {
193  HMdcSizesCellsLayer& rSizesCellsLay = pSizesCellsMod[layer];
194  setILay(layer); // calculate iLay
195  if(&rSizesCellsLay == 0) continue;
196  pSCLayer[iLay] = &rSizesCellsLay;
197  nLayerParts[iLay] = rSizesCellsLay.getLayerNParts();
198  if (shiftType == 0) setTrShiftInSecSys();
199  else if(shiftType == 1) setTrShiftInModSys();
200  else if(shiftType == 2) setTrShiftInLabSys();
201  }
202  }
203  }
204 
205  for(Int_t s=0;s<6;s++) for(Int_t m=0;m<4;m++) for(Int_t l=-1;l<6;l++) {
206  if(l>=0 && useMdcShParOnly) continue;
207  for(Int_t p=0;p<nParMax;p++) {
208  if(l <0 && (!mShFlag[s][m][p] || mFixFlag[s][m][p] || sigmaMdc[s][m][p]<0.) ) continue;
209  if(l>=0 && (!lShFlag[s][m][l][p] || lFixFlag[s][m][l][p] || sigmaLay[s][m][l][p]<0.)) continue;
210  isMdcInAlign[s][m] = kTRUE;
211  }
212  }
213 
214  return kTRUE;
215 }
216 
217 void HMdcMille::createHists(Int_t sec) {
218  if(fHistograms != 0 && nHists>0) return;
219  if(fHistograms == 0) fHistograms = new TList();
220  char name[100];
221  char title[100];
222  histInd[sec] = nHists;
223 
224  sprintf(name,"pThVsPhS%i",sec+1);
225  sprintf(title,"Sec.%i #Theta vs #phi",sec+1);
226  fHistograms->Add(new TH2F(name,title,60,60.,120.,80,10.,90.));
227  nHists++;
228 
229  sprintf(name,"hChi2S%i",sec+1);
230  sprintf(title,"Sec.%i Chi^{2}/NDF",sec+1);
231  fHistograms->Add(new TH1F(name,title,200,0.,20.));
232  nHists++;
233 
234  sprintf(name,"hNWiresS%i",sec+1);
235  sprintf(title,"Sec.%i Num.wires per track",sec+1);
236  Int_t nBins = clusFitAlg==4 ? 100 : 50;
237  fHistograms->Add(new TH1F(name,title,nBins,0.,nBins));
238  nHists++;
239 
240  nBins = clusFitAlg==4 ? 49 : 25;
241  sprintf(name,"hNLayersS%i",sec+1);
242  sprintf(title,"Sec.%i Num.layers per track",sec+1);
243  fHistograms->Add(new TH1F(name,title,nBins,0.,nBins));
244  nHists++;
245 }
246 
248  Double_t steps[6]={0.,0.,0.,0.,0.,0.};
249  for(Int_t parn=0;parn<6;parn++) {
250  derNorm[parn] = 1./parSteps[parn];
251  steps[parn] = parSteps[parn];
253 
254  steps[parn] = -parSteps[parn];
256  steps[parn] = 0.;
257  }
258  // For layer thickness align. parameter:
259  for(Int_t lay=0;lay<6;lay++) {
260  steps[2] = parSteps[6]*(lay-3+0.5);
261  derLThNorm[lay] = 1./steps[2];
263 
264  steps[2] = -steps[2];
266  steps[2] = 0.;
267  }
268 }
269 
271  // mdcSys - mdc<->sec. transfomation
272  // laySec - rotated_layer<->mdc transfomation
273  const HGeomTransform *mdcSys = pSCMod[iMod]->getSecTrans();
274  const HGeomTransform *laySec = pSCLayer[iLay]->getRotLayP1SysRMod();
275  const HGeomTransform *layP2Mod = NULL;
276  if(nLayerParts[iLay] == 2) layP2Mod = pSCLayer[iLay]->getRotLayP2SysRMod();
277  for(Int_t parn=0;parn<6;parn++) {
278  calcShiftInModSys(mdcSys,laySec,posShifts[parn],layPosSh[iLay][parn]);
279  calcShiftInModSys(mdcSys,laySec,negShifts[parn],layNegSh[iLay][parn]);
280  if(layP2Mod == NULL) continue;
281  calcShiftInModSys(mdcSys,layP2Mod,posShifts[parn],layPosSh[iLay][parn+20]);
282  calcShiftInModSys(mdcSys,layP2Mod,negShifts[parn],layNegSh[iLay][parn+20]);
283  }
284  // Cell thickness:
285  calcShiftInModSys(mdcSys,laySec,posLThSh[layer],layPosSh[iLay][6]);
286  calcShiftInModSys(mdcSys,laySec,negLThSh[layer],layNegSh[iLay][6]);
287  if(layP2Mod != NULL) {
288  calcShiftInModSys(mdcSys,layP2Mod,posLThSh[layer],layPosSh[iLay][6+20]);
289  calcShiftInModSys(mdcSys,layP2Mod,negLThSh[layer],layNegSh[iLay][6+20]);
290  }
291 }
292 
294  // laySec - rotated_layer<->sec transfomation
295  const HGeomTransform* laySec = pSCLayer[iLay]->getRotLayP1SysRSec();
296  const HGeomTransform* layP2Sec = NULL;
297  if(nLayerParts[iLay] == 2) layP2Sec = pSCLayer[iLay]->getRotLayP2SysRSec();
298  for(Int_t parn=0;parn<6;parn++) {
299  calcShiftInSecSys(laySec,posShifts[parn],layPosSh[iLay][parn]);
300  calcShiftInSecSys(laySec,negShifts[parn],layNegSh[iLay][parn]);
301  if(layP2Sec == NULL) continue; // Layer second part:
302  calcShiftInSecSys(layP2Sec,posShifts[parn],layPosSh[iLay][parn+20]);
303  calcShiftInSecSys(layP2Sec,negShifts[parn],layNegSh[iLay][parn+20]);
304  }
305  // Cell thickness:
307  calcShiftInSecSys(laySec,negLThSh[layer],layNegSh[iLay][6]);
308  if(layP2Sec == NULL) return;
309  // Layer second part:
310  calcShiftInSecSys(layP2Sec,posLThSh[layer],layPosSh[iLay][6+20]);
311  calcShiftInSecSys(layP2Sec,negLThSh[layer],layNegSh[iLay][6+20]);
312 }
313 
315  // secSys - sec.<->lab. transfomation
316  // laySec - rotated_layer<->mdc transfomation
318  const HGeomTransform *secSys = (*pSCells)[iSec].getLabTrans();
319  const HGeomTransform *laySec = pSCLayer[iLay]->getRotLayP1SysRSec();
320  const HGeomTransform *layP2Sec = NULL;
321  if(nLayerParts[iLay] == 2) layP2Sec = pSCLayer[iLay]->getRotLayP2SysRSec();
322  for(Int_t parn=0;parn<6;parn++) {
323  calcShiftInLabSys(secSys,laySec,posShifts[parn],layPosSh[iLay][parn]);
324  calcShiftInLabSys(secSys,laySec,negShifts[parn],layNegSh[iLay][parn]);
325  if(layP2Sec == NULL) continue;
326  calcShiftInLabSys(secSys,layP2Sec,posShifts[parn],layPosSh[iLay][parn+20]);
327  calcShiftInLabSys(secSys,layP2Sec,negShifts[parn],layNegSh[iLay][parn+20]);
328  }
329  // Cell thickness:
330  calcShiftInLabSys(secSys,laySec,posLThSh[layer],layPosSh[iLay][6]);
331  calcShiftInLabSys(secSys,laySec,negLThSh[layer],layNegSh[iLay][6]);
332  if(layP2Sec != NULL) {
333  calcShiftInLabSys(secSys,layP2Sec,posLThSh[layer],layPosSh[iLay][6+20]);
334  calcShiftInLabSys(secSys,layP2Sec,negLThSh[layer],layNegSh[iLay][6+20]);
335  }
336 }
337 
338 void HMdcMille::calcShiftInModSys(const HGeomTransform* mdcSys, const HGeomTransform* laySysMod,
339  const HGeomTransform& shift, HGeomTransform& laySh) {
340  // mdcSys - mdc<->sec. transfomation
341  // laySysMod - rotated_layer<->mdc transfomation
342  // shift - transfomation for shift of one global parameter
343  // laySh - new (shifted) rotated_layer<->sec transfomation
344  laySh = *laySysMod;
345  laySh.transFrom(shift);
346  laySh.transFrom(*mdcSys);
347 }
348 
349 void HMdcMille::calcShiftInSecSys( const HGeomTransform* laySysSec, const HGeomTransform& shift,
350  HGeomTransform& laySh) {
351  // laySysSec - rotated_layer<->sec transfomation
352  // shift - transfomation for shift of one global parameter
353  // laySh - new (shifted) rotated_layer<->sec transfomation
354  laySh = *laySysSec;
355  laySh.transFrom(shift); // Do shift of layer<->sec. transfomation
356 }
357 
358 void HMdcMille::calcShiftInLabSys(const HGeomTransform* secSys, const HGeomTransform* laySysSec,
359  const HGeomTransform& shift, HGeomTransform& laySh) {
360  // secSys - lab.<->sec. transfomation
361  // laySysSec - rotated_layer<->sec. transfomation
362  // shift - transfomation for shift of one global parameter
363  // laySh - new (shifted) rotated_layer<->sec transfomation
364  laySh = *laySysSec;
365  laySh.transFrom(*secSys);
366  laySh.transFrom(shift);
367  laySh.transTo(*secSys);
368 }
369 
370 Bool_t HMdcMille::finalize(void) {
371  if(mille) delete mille;
372  mille = 0;
374  if(!writeParAsciiFile()) printErrorAndExit("finalize","Can't write geometry parameters ASCII file!!!");
377  return kTRUE;
378 }
379 
380 Int_t HMdcMille::execute(void) {
381  if(HMdcTrackDSet::fPrint()) {
383  if(fitpar.getPrintFlag()) printf(
384  "============================== Event %i =============================\n",
386  }
387  if(clusFitAlg==3) {
388  if(filteringFlag>0) {
389  if( !eventFilter() ) return kSkipEvent;
390  if( filteringFlag==2 ) return 0;
391  }
393  } else if(clusFitAlg==4) fitAlgForMilleCosmic();
394  return 0;
395 }
396 
398  // For beam data only!
399  HMdcClus* pClst;
400  UInt_t nClusters[6] = {0,0,0,0,0,0};
401  iterClus->Reset();
402  while((pClst=(HMdcClus*)iterClus->Next())) if(pClst->getIndexParent()<0) {
403  UChar_t sec = pClst->getSec();
404  nClusters[sec]++;
405  }
407  Bool_t isGoodSectors = kFALSE;
408  for(Int_t s=0;s<6;s++) {
409  useSector[s] = kFALSE;
410  if(nMods[s] <= 1) continue;
411  if(nClusters[s]>0 && nClusters[s]<=nClustersCut) {
412  if(pEvLCells != NULL) {
413  Int_t nCells = (*pEvLCells)[s].getNCells();
414  Int_t nWiresC = (nWiresMin*nMods[s])/4;
415  if(nCells>=nWiresC && nCells<=nWiresCut[s]) {
416  useSector[s] = kTRUE;
417  isGoodSectors = kTRUE;
418  }
419  } else {
420  useSector[s] = kTRUE;
421  isGoodSectors = kTRUE;
422  }
423  }
424  }
425  return isGoodSectors;
426 }
427 
429  // Sector fit
430  // Magnet off + combined clusters in sector
431  HMdcClus* fClst1;
432  iterClus->Reset();
433  while((fClst1=(HMdcClus*)iterClus->Next())) if(fClst1->getIndexParent()<0) {
434  UChar_t sec = fClst1->getSec();
435  if(!useSector[sec]) continue;
436  //--- Angle theta cut:
437 //!!! if(thetaCut>0. && fClst1->getTheta()*TMath::RadToDeg()<thetaCut) continue;
438 // Double_t theta = fClst1->getTheta()*TMath::RadToDeg();
439 // Double_t phi = fClst1->getPhi()*TMath::RadToDeg();
440 // Int_t sec = fClst1->getSec();
441 /*
442 if(sec==1) {
443  if(theta<35.) {
444  if(phi>=108. && theta<35.) continue;
445  if(phi>=102. && theta<25.) continue;
446  if(phi>=90. && theta<20.) continue;
447  }
448 } else if(sec==2) {
449  if(theta<30.) {
450  if(phi>=96. && theta<30.) continue;
451  if(phi>=78. && theta<25.) continue;
452  }
453 } else if(sec==3) {
454  if(phi>=96. && theta<25.) continue;
455 } else if(sec==4) {
456  if(phi>=108. && theta<25.) continue;
457 }
458 */
459  fittersArr[0].reset();
460  Int_t first,last;
461  fClst1->getIndexRegChilds(first,last);
462  HMdcClus* fClst2 = (first<0) ? 0:(HMdcClus*)fClusCat->getObject(first);
463  // Num.wires test:
464  Int_t nWr = fClst1->getNDrTimes();
465  if(fClst2) nWr += fClst2->getNDrTimes();
466  nWiresMinTr = (nWiresMin*nMods[sec])/4;
467  if(nWr < nWiresMinTr) continue;
468 
469  Int_t nWiresInAlign = 0;
470  Int_t mod = fClst1->getIOSeg()*2;
471  if( isMdcInAlign[sec][mod] ) nWiresInAlign += fClst1->getNDrTimesMod(0);
472  if( isMdcInAlign[sec][mod+1] ) nWiresInAlign += fClst1->getNDrTimesMod(1);
473  mod = fClst2->getIOSeg()*2;
474  if( isMdcInAlign[sec][mod] ) nWiresInAlign += fClst2->getNDrTimesMod(0);
475  if( isMdcInAlign[sec][mod+1] ) nWiresInAlign += fClst2->getNDrTimesMod(1);
476  if(nWiresInAlign < 6) continue;
477 
478  Bool_t flag=kFALSE;
479  if(fClst2==0) {
480  if(fClst1->getNDrTimesMod(0)==0 || fClst1->getNDrTimesMod(1)==0) continue;
481  Int_t typeClFn = fClst1->getTypeClFinder();
482  if(typeClFn==0) flag = fitSeg(fClst1); // Two mdc per segment
483  else if(typeClFn==2) flag = fitMixedClus(fClst1);
484  if(!flag) continue;
485  fitter=fittersArr[0].getFitter(0); // ???
486  if(fitter->getSegIndex()<0) continue;
487  fillTrkCandISeg();
488  } else {
489  flag = fitSec(fClst1,fClst2); //??????????????
490  }
491  if(flag && fitter->getFitStatus()) sendToMille();
492  }
493 }
494 
496  // Two sector fit
497  // Magnet off + combined clusters in sector
499  HMdcClus* fClst = NULL;
500  iterClus->Reset();
501  HMdcClus *fClstArr[6][2]; // [sector][0,1 -IOseg]
502  Int_t nSegs[6];
503  Int_t nMods[6];
504  Int_t secList[4];
505  Int_t nLaySec[6];
506  Int_t &s1 = secList[0];
507  Int_t &s2 = secList[1];
508  Int_t &s3 = secList[2];
509  Int_t &s4 = secList[3];
510  // In cosmic one track per event only! First cluster finding:
511  while((fClst=(HMdcClus*)iterClus->Next())) if(fClst->getIndexParent()<0) {
512  Int_t nSectors = 0;
513  Int_t nSegments = 1;
514  Int_t nModules = 0;
515  Int_t nWires = 0;
516  for(Int_t s=0;s<6;s++) {
517  fClstArr[s][0] = NULL;
518  fClstArr[s][1] = NULL;
519  nSegs[s] = 0;
520  nMods[s] = 0;
521  nLaySec[s] = 0;
522  }
523  Int_t nWiresInAlign = 0;
524  while(kTRUE) {
525  Int_t sec = fClst->getSec();
526  Int_t seg = fClst->getIOSeg();
527  if(fClstArr[sec][seg] != 0) {
528  Error("fitAlgForMilleCosmic","Second time HMdcClus for sec.%i seg.%i",sec,seg);
529  return;
530  }
531  Int_t mod = seg*2;
532  if( isMdcInAlign[sec][mod] ) nWiresInAlign += fClst->getNDrTimesMod(0);
533  if( isMdcInAlign[sec][mod+1] ) nWiresInAlign += fClst->getNDrTimesMod(1);
534 
535  fClstArr[sec][seg] = fClst;
536  if(nSegs[sec] == 0) {
537  secList[nSectors] = sec;
538  nSectors++;
539  }
540  nLaySec[sec] += fClst->getNLayers();
541  nSegs[sec]++;
542  nSegments++;
543  Int_t nm = fClst->getActiveModule()<0;
544  nMods[sec] += nm;
545  nModules += nm;
546  nWires += fClst->getNDrTimes();
547  Int_t first,last;
548  fClst->getIndexRegChilds(first,last);
549  if(first < 0 || nSectors == 4) break;
550  fClst = (HMdcClus*)fClusCat->getObject(first);
551  }
552 
553  if(nWiresInAlign < 6) continue;
554 
555  for(Int_t s=0;s<6;s++) if(fClstArr[s][0]==NULL && fClstArr[s][1]!=NULL) {
556  fClstArr[s][0] = fClstArr[s][1];
557  fClstArr[s][1] = NULL;
558  }
559  if(nModules==1) continue;
560  if(nModules>=4) nWiresMinTr = nWiresMin;
561  else nWiresMinTr = (nWiresMin*nModules)/4;
562 
563  if(nWires < nWiresMinTr) continue; // Num.wires test
564 
565 
566  fittersArr[0].reset();
567  Bool_t flag=kFALSE;
568  if(nSegments == 1) {
569  if(nModules < 2) continue;
570  Int_t typeClFn = fClst->getTypeClFinder();
571  if(typeClFn==0) flag = fitSeg(fClst);
572  else if(typeClFn==2) flag = fitMixedClus(fClst);
573  } else {
574  if (nSectors == 1) flag = fitSec(fClstArr[s1][0],fClstArr[s1][1]);
575  else if(nSectors == 2) flag = fit2Sectors(fClstArr[s1][0],fClstArr[s1][1],
576  fClstArr[s2][0],fClstArr[s2][1]);
577 // else if(nSectors > 2) { //???????????
578 // Int_t sec1 = 0;
579 // for(Int_t s=0;s<6;s++) if(s!=sec1 && nLaySec[s]>nLaySec[sec1]) sec1 = s;
580 // Int_t sec2 = sec1>0 ? 0 : 1;
581 // for(Int_t s=0;s<6;s++) if(s!=sec1 && s!=sec2 && nLaySec[s]>nLaySec[sec2]) sec2 = s;
582 // flag = fit2Sectors(fClstArr[sec1][0],fClstArr[sec1][1],fClstArr[sec2][0],fClstArr[sec2][1]);
583 // }
584  else if(nSectors == 3) flag = fitNSectors(fClstArr[s1][0],fClstArr[s1][1],
585  fClstArr[s2][0],fClstArr[s2][1],
586  fClstArr[s3][0],fClstArr[s3][1],NULL,NULL);
587  else if(nSectors == 4) flag = fitNSectors(fClstArr[s1][0],fClstArr[s1][1],
588  fClstArr[s2][0],fClstArr[s2][1],
589  fClstArr[s3][0],fClstArr[s3][1],
590  fClstArr[s4][0],fClstArr[s4][1]);
591 
592 
593 
594  }
595 
596  if(flag && fitter->getFitStatus()) sendToMille();
597 
598 // Int_t first,last;
599 // Int_t sec = fClst1->getSec();
600 // fClst1->getIndexRegChilds(first,last);
601 // HMdcClus* fClst2 = (first<0) ? NULL : (HMdcClus*)fClusCat->getObject(first);
602 // HMdcClus* fClst3 = NULL;
603 // HMdcClus* fClst4 = NULL;
604 // if(fClst2) {
605 // if(sec == fClst2->getSec()) {
606 // fClst2->getIndexRegChilds(first,last);
607 // fClst3 = (first<0) ? 0:(HMdcClus*)fClusCat->getObject(first);
608 // } else {
609 // fClst3 = fClst2;
610 // fClst2 = NULL;
611 // }
612 // if(fClst3) {
613 // fClst3->getIndexRegChilds(first,last);
614 // fClst4 = (first<0) ? 0:(HMdcClus*)fClusCat->getObject(first);
615 // }
616 // }
617 // Int_t sec2 = -1;
618 // if(fClst3) sec2 = fClst3->getSec();
619 // else if(fClst4) sec2 = fClst4->getSec();
620 // Int_t nModules = 0;
621 // nModules += fClst1->getActiveModule()<0 ? 2:1;
622 // if(fClst2) nModules += fClst2->getActiveModule()<0 ? 2:1;
623 // if(fClst3) nModules += fClst3->getActiveModule()<0 ? 2:1;
624 // if(fClst4) nModules += fClst4->getActiveModule()<0 ? 2:1;
625 // if(nModules==1) continue;
626 // if(nModules>=4) nWiresMinTr = nWiresMin;
627 // else nWiresMinTr = (nWiresMin*nModules)/4;
628 //
629 // // Num.wires test:
630 // Int_t nWr = fClst1->getNDrTimes();
631 // if(fClst2) nWr += fClst2->getNDrTimes();
632 // if(fClst3) nWr += fClst3->getNDrTimes();
633 // if(fClst4) nWr += fClst4->getNDrTimes();
634 // if(nWr < nWiresMinTr) continue;
635 //
636 // Bool_t flag=kFALSE;
637 // if(fClst2==0 && fClst3==0 && fClst4==0) {
638 // if(fClst1->getNDrTimesMod(0)==0 || fClst1->getNDrTimesMod(1)==0) continue;
639 // Int_t typeClFn = fClst1->getTypeClFinder();
640 // if(typeClFn==0) flag = fitSeg(fClst1);
641 // else if(typeClFn==2) flag = fitMixedClus(fClst1);
642 // } else {
643 // if(fClst3==0&&fClst4==0) flag = fitSec(fClst1,fClst2);
644 // else flag = fit2Sectors(fClst1,fClst2,fClst3,fClst4);
645 // }
646 // if(flag && fitter->getFitStatus()) sendToMille();
647  }
648 }
649 
651  if(mille == 0) return;
652 
654  HMdcWiresArr& wiresArr = fitter->getWiresArr();
655  Int_t nWiresData = wiresArr.getNDriftTimes();
656 
657  // test chi2 and ...:
658  // if(fitter->getMaxChi2()>9) return;
659 
660  // TEST TEST ???????????????!!!!!!!!!!!!!!!!!!!!!!!!!!!
661  HMdcList24GroupCells& listCells = wiresArr.getOutputListCells();
662  for(Int_t mod=0;mod<4;mod++) {
663  Int_t nLayersMod = listCells.getNLayersMod(mod);
664  if(nLayersMod == 1) return;
665 //??? if(mod==1 && nLayersMod < 4) return;
666  }
667  if(wiresArr.getSector2()>=0) {
668  HMdcList24GroupCells &listCells2 = wiresArr.getOutputListCells2();
669  for(Int_t mod=0;mod<4;mod++) if(listCells2.getNLayersMod(mod) == 1) return;
670  if(wiresArr.getSector3()>=0) {
671  HMdcList24GroupCells &listCells3 = wiresArr.getOutputListCells3();
672  for(Int_t mod=0;mod<4;mod++) if(listCells3.getNLayersMod(mod) == 1) return;
673  if(wiresArr.getSector4()>=0) {
674  HMdcList24GroupCells &listCells4 = wiresArr.getOutputListCells4();
675  for(Int_t mod=0;mod<4;mod++) if(listCells4.getNLayersMod(mod) == 1) return;
676  }
677  }
678  }
679 
680  for(Int_t wInd=0;wInd<nWiresData;wInd++) {
681  wireData = wiresArr.getWireData(wInd);
682  if(wireData->getHitStatus() != 1) continue;
683  if(!wireData->isInCell()) return;
684  }
685 
686 
689 
690  Int_t nWires = finalParam->getNumOfGoodWires();
691  Int_t nLayers = wiresArr.getOutputNLayers();
692  if(nWires < nWiresMinTr) return;
693 
694  // Filling histograms:
695  if(doHists) {
696  Int_t hstInd = histInd[wiresArr.getSector()];
697  if(hstInd>=0) {
698  ((TH2F*)fHistograms->At(hstInd))->Fill(finalParam->getPhiDeg(),finalParam->getThetaDeg());
699  ((TH1F*)fHistograms->At(hstInd+1))->Fill(finalParam->getChi2());
700  ((TH1F*)fHistograms->At(hstInd+2))->Fill(nWires+0.1);
701  ((TH1F*)fHistograms->At(hstInd+3))->Fill(nLayers+0.1);
702  }
703  }
704 
705  // Derivatives calculation:
706  if(fixFitTOffset) {
707  for(Int_t wInd=0;wInd<nWiresData;wInd++) {
708  wireData = wiresArr.getWireData(wInd);
709  if(!wireData->getAnalytDeriv(locDer,0)) continue;
710  int nGlobDer = getGlobalDer();
712  } // nWires
713  } else {
714 
715  for(Int_t mi=0;mi<24;mi++) for(Int_t p=0;p<nParMax;p++) tofDerCorr[mi][p] = 0.;
716 
717  for(Int_t wInd=0;wInd<nWiresData;wInd++) {
718  wireData = wiresArr.getWireData(wInd);
719  if(wireData->getHitStatus() != 1) continue;
720  setIModILay(); // calculate iMod and iLay + ...
721  for(Int_t p=0;p<6;p++) tofDerCorr[iMod][p] += calcDerTofCorr(p);
722  if(mShFlag[sector][module][6]) tofDerCorr[iMod][6] += calcDerTofCorr(6);
723  }
724 
725  for(Int_t wInd=0;wInd<nWiresData;wInd++) {
726  wireData = wiresArr.getWireData(wInd);
727  if(!wireData->getAnalytDeriv(locDer,finalParam)) continue;
728  int nGlobDer = getGlobalDerWTof();
729  mille->mille(4,locDer,nGlobDer,globDer,globLabel,
731 // if(module==0) printf("%iw.Glob: Der %g %g %g %g %g %g %g %g\n",
732 // wInd,globDer[0],globDer[1],globDer[2],globDer[3],globDer[4],globDer[5],
733 // globDer[6],globDer[7]);
734  } // nWires
735  }
736  nWiresTot += nWires;
737  nLayersTot += nLayers;
738  nTracks++;
739 
740  milleFileSize += mille->end();
741  if(milleFileSize > 1800000000) openNewBinaryFile();
742 }
743 
745  // New binary file:
746  if(mille) delete mille;
747  TString newFile(milleFileName);
748  newFile += nBinaryFiles;
749  mille = new Mille(newFile);
750  milleFileSize = 0;
751  nBinaryFiles++;
752 }
753 
757  layer = wireData->getLayer();
758  setIModILay(sector,module,layer); // calculate iMod and iLay
762  if(layerPart == 1) {
763  layPosShCurr += 20;
764  layNegShCurr += 20;
765  }
766 }
767 
769  setIModILay();
770  Int_t nParam = 0;
771 //!!! Int_t nMdcPar = mShFlag[sector][module][6] ? 7 : 6;
772  Int_t nMdcPar = cellThicknFree ? 7 : 6;
773  for(Int_t p=0;p<nMdcPar;p++) {
774  globLabel[p] = packLabel(sector,module,-1,p); //!!!
775  globDer[p] = calcGlobDer(p);
776  nParam++;
777  }
778 
779 // Double_t cs[6]={0.766044,0.939693,1.,1.,0.939693,0.766044};
780  if(shiftType == 1 && !useMdcShParOnly) {
781  if(layerPart == 0 || !doLayP2Align) { // Layer part 1
782  globLabel[nParam] = packLabel(sector,module,layer,1); // in mdc sys
783  globDer[nParam] = globDer[1]; // go to central wire number
784  // ? globDer[nParam] = globDer[1]*cs[layer];
785  nParam++;
786  // ? globLabel[nParam] = packLabel(sector,module,layer,2);
787  // ? globDer[nParam] = globDer[2];
788  // ? nParam++;
789  globLabel[nParam] = packLabel(sector,module,layer,5); // in mdc sys
790  globDer[nParam] = globDer[5]; // go to the wire orientation
791  nParam++;
792  } else if(layerPart == 1) { // Layer part 2 (it can be in mdc3 & 4 only)
793  globLabel[nParam] = packLabel(sector,module,layer,7); // in mdc sys
794  globDer[nParam] = globDer[1]; // go to central wire number
795  nParam++;
796  globLabel[nParam] = packLabel(sector,module,layer,8); // in mdc sys
797  globDer[nParam] = globDer[5]; // go to the wire orientation
798  nParam++;
799  }
800  }
801 
802 //printf("0Glob: Der %.3g %.3g %.3g %.3g %.3g %.3g %.3g %.3g\n",
803 //globDer[0],globDer[1],globDer[2],globDer[3],globDer[4],globDer[5],globDer[6],globDer[7]);
804  return nParam;
805 }
806 
808  setIModILay();
809 //!!! Int_t nMdcPar = mShFlag[sector][module][6] ? 7 : 6;
810  Int_t nMdcPar = cellThicknFree ? 7 : 6;
811  Int_t nParam = nMdcPar;
812  for(Int_t p=0;p<nMdcPar;p++) {
813  globLabel[p] = packLabel(sector,module,-1,p);
814  Double_t dTm;
815  Double_t corr = calcDerTofCorr(p,&dTm);
816  globDer[p] = (dTm + tofDerCorr[iMod][p])*derNorm[p];
817  if(useMdcShParOnly || shiftType!=1 || (p!=1 && p!=5) ) continue;
818  // For HMdcLayerGeomPar:
819 
820  Int_t pInd = p==1 ? nMdcPar : nMdcPar+1;
821  Int_t po = p;
822  if(layerPart==1 && doLayP2Align) {
823  po = p==1 ? 7 : 8;
824  pInd += 2;
825  }
826  globLabel[pInd] = packLabel(sector,module,layer,po);
827  globDer[pInd] = (dTm + corr)*derNorm[p];
828  nParam++;
829  }
830  return nParam;
831 }
832 
833 Double_t HMdcMille::calcDerTofCorr(Int_t p,Double_t* dDrTm) {
834  Double_t dDrT = calcDDriftTime(p);
835  if(dDrTm) *dDrTm = dDrT;
836  Double_t wtNorm = wireData->getWtNorm();
837  Double_t sumWtNorm = finalParam->getSumWtNorm(wireData->getModIndex());
838  if(sumWtNorm > 0.) return -dDrT*wtNorm/sumWtNorm;
839  return 0.;
840 }
841 
842 int HMdcMille::packLabel(Int_t s,Int_t m,Int_t l,Int_t parn) {
843  // if lay<0 - it is shift of mdc, otherwise shift of layer
844  if(l<0) l = -1;
846  if(l<0) mShFlag[s][m][parn] = kTRUE;
847  else lShFlag[s][m][l][parn] = kTRUE;
848  }
849  int label = (((shiftType*10 + s+1)*10 + m+1)*10 + l+1)*10 + parn+1;
850  if(beamTimeId>0 && l<0 && parn!=6) label += beamTimeId*100000;
851  return label;
852 }
853 
854 Bool_t HMdcMille::unpackLabel(Int_t label,Int_t& btId,Int_t& sys,
855  Int_t& sec,Int_t& mod,Int_t& lay,Int_t& parn) {
856  // shSys = 0 shift in respect sector system
857  // shSys = 1 shift in respect mdc system
858  // if lay==-1 - it is shift of mdc, otherwise shift of layer
859  // parn = 0 - shift along X
860  // parn = 1 - shift along Y
861  // parn = 2 - shift along Z
862  // parn = 3 - rotation around X
863  // parn = 4 - rotation around Y
864  // parn = 5 - rotation around Z
865  // parn = 6 - cell thickness
866  // parn = 7 - shift along Y for second layer part
867  // parn = 8 - rotation around Z for second layer part
868  btId = label/100000; // beam time Id
869  sys = (label%100000)/10000;
870  sec = (label%10000)/1000 - 1;
871  mod = (label%1000)/100 - 1;
872  lay = (label%100)/10 - 1;
873  parn = (label%10) - 1;
874  if(sec<0||sec>5 || mod<0||mod>3 || lay<-1||lay>5 || parn<0) return kFALSE;
875  return kTRUE;
876 }
877 
878 Double_t HMdcMille::calcGlobDer(Int_t parNum) {
879  Double_t dDrT = calcDDriftTime(parNum);
880  if(parNum==6) return dDrT*derLThNorm[layer];
881  else return dDrT*derNorm[parNum]; // derNorm=1/parSteps
882 }
883 
884 Double_t HMdcMille::calcDDriftTime(Int_t p) {
885  Int_t dSign1 = {0};
886  Int_t dSign2 = {0};
887  Double_t drTime1 = calcDriftTime(layPosShCurr[p],dSign1);
888  Double_t drTime2 = calcDriftTime(layNegShCurr[p],dSign2);
889  if(dSign1 == dSign2) return (drTime1-drTime2)*0.5;
890  Int_t dSignT = wireData->getDistanceSign();
891  if(dSignT==dSign1) return drTime1 - wireData->getDrTime();
892  return wireData->getDrTime() - drTime2; // dSignT==dSign2
893 }
894 
895 Double_t HMdcMille::calcDriftTime(const HGeomTransform& laySys,Int_t& distSign) {
896  HGeomVector p1l,p2l;
897  Int_t parSec = finalParam->getSec();
898  if(parSec==sector || parSec<0) {
899  p1l = p1;
900  p2l = p2;
901  } else {
903  pSCells->transLineToOtherSec(*finalParam,sector,p1l,p2l);
904  }
905  p1l = laySys.transTo(p1l);
906  p2l = laySys.transTo(p2l);
907  return wireData->calcDriftTimeForAlign(p1l,p2l,distSign);
908 }
909 
912  if(pedeInParFName.Length() > 0) createPedeInFile = kTRUE;
913  else Error("setPedeInFileName","File name is not specified!!!");
914 }
915 
917  FILE* file = fopen(pedeInParFName.Data(),"r");
918  if(file != NULL) {
919  printf("******CreatePedeInParamFile******\n");
920  printf("* File %s exist!\n",pedeInParFName.Data());
921  printf("* If you want to change this file\n");
922  printf("* do it by hand or delete file.\n");
923  printf("*********************************\n");
924  fclose(file);
925  return;
926  }
927  file = fopen(pedeInParFName.Data(),"w");
928  if(file == NULL) Error("creatPedeInParamFile","Can't open file %s",pedeInParFName.Data());
929  else {
930  TString sys;
931  if(shiftType == 0) sys = "All shifts done in sector sys.";
932  if(shiftType == 1) sys = "All shifts done in mdc sys.";
933  if(shiftType == 2) sys = "All shifts done in lab. sys.";
934  fprintf(file," Parameter ! %s\n",sys.Data());
935  for(Int_t s=0;s<6;s++) for(Int_t m=0;m<4;m++) for(Int_t l=-1;l<6;l++) {
936  if(l>=0 && useMdcShParOnly) continue;
937  for(Int_t p=0;p<nParMax;p++) {
938  if(l<0 && !mShFlag[s][m][p]) continue;
939  if(l>=0 && !lShFlag[s][m][l][p]) continue;
940  if( (l<0 && mFixFlag[s][m][p]) || (l>=0 && lFixFlag[s][m][l][p]) )
941  fprintf(file,"%11i 0.0 -1.0000 !",packLabel(s,m,l,p));
942  else fprintf(file,"%11i 0.0 0.0000 !",packLabel(s,m,l,p));
943 
944  if(l<0) fprintf(file," %is. %im.: mdc",s+1,m+1);
945  else {
946  fprintf(file," %is. %im. %il.: layer",s+1,m+1,l+1);
947  if(p==7 || p==8) fprintf(file," part II");
948  }
949  fprintf(file," %s\n",shitsInfo[p].Data()); // Print description
950  }
951  }
952  fclose(file);
953  }
954 }
955 
957  HRuntimeDb* rtdb = gHades->getRuntimeDb();
958  if(rtdb == 0) return kFALSE;
959  HParAsciiFileIo* inputFile = new HParAsciiFileIo;
960  if(!inputFile->open(const_cast<char*>(geomParFileName.Data()),"in")) return kTRUE; // No file!
961  HMdcGeomStruct* pMdcGeomStruct = (HMdcGeomStruct*)rtdb->getContainer("MdcGeomStruct");
962  if( !((HParSet*)pMdcGeomStruct)->init(inputFile) ) return kFALSE;
963  pMdcGeomStruct->setStatic();
964  pMdcGeomStruct->setInputVersion(1,1);
965 
966  HSpecGeomPar* pSpecGeomPar = (HSpecGeomPar*)rtdb->getContainer("SpecGeomPar");
967  HMdcGeomPar* pMdcGeomPar = (HMdcGeomPar*)rtdb->getContainer("MdcGeomPar");
968  HMdcLayerGeomPar* pMdcLGeomPar = (HMdcLayerGeomPar*)rtdb->getContainer("MdcLayerGeomPar");
969  HMdcLayerCorrPar* pMdcLCorrPar = (HMdcLayerCorrPar*)rtdb->getContainer("MdcLayerCorrPar");
970 
971  if( !((HParSet*)pSpecGeomPar)->init(inputFile) ||
972  !((HParSet*)pMdcGeomPar )->init(inputFile) ||
973  !((HParSet*)pMdcLGeomPar)->init(inputFile) ||
974  !((HParSet*)pMdcLCorrPar)->init(inputFile)) return kFALSE;
975 
976  pSpecGeomPar->setInputVersion(1,1);
977  pMdcGeomPar ->setInputVersion(1,1);
978  pMdcLGeomPar->setInputVersion(1,1);
979  pMdcLCorrPar->setInputVersion(1,1);
980 
981  pSpecGeomPar->setStatic();
982  pMdcGeomPar ->setStatic();
983  pMdcLGeomPar->setStatic();
984  pMdcLCorrPar->setStatic();
985 
986  inputFile->close();
987  delete inputFile;
988  isGeomFileExist = kTRUE;
989  return kTRUE;
990 }
991 
993  if(fileName == 0) printErrorAndExit("readPedeResFile","File name of Pede results not specified!");
994  // Open file:
995  FILE* file = fopen(fileName,"r");
996  if(file == NULL) {
997  if(iteration==0) {
998  for(Int_t s=0;s<6;s++) for(Int_t m=0;m<4;m++) for(Int_t p=0;p<9;p++) mShFlag[s][m][p] = kTRUE;
999  return;
1000  }
1001  printErrorAndExit("readPedeResFile","No pede result file!");
1002  }
1003 
1004  // Read file:
1005  TString buffer;
1006  if( !buffer.Gets(file) ) printErrorAndExit("readPedeResFile","Wrong format of file %s",fileName);
1007  Ssiz_t pos = buffer.Index("Parameter");
1008  if(pos<0 || pos>10) printErrorAndExit("readPedeResFile","Wrong format of file %s",fileName);
1009  while(buffer.Gets(file) && buffer.Length()>16 && addShifts(buffer));
1010 
1011  fclose(file);
1012 }
1013 
1015  if( !isGeomChanged ) return kTRUE;
1016  if( !isGeomFileExist ) return kFALSE;
1017  HMdcLayerGeomPar* pMdcLGeomPar = pGetCont->getMdcLayerGeomPar();
1018  HMdcLayerCorrPar* pMdcLCorrPar = pGetCont->getMdcLayerCorrPar();
1019  Bool_t reCalcLayerTransf = kFALSE;
1020  HGeomTransform shift;
1021  HGeomTransform trans;
1022  for(Int_t s=0;s<6;s++) for(Int_t m=0;m<4;m++) {
1023  if(pMdcDet->getModule(s,m) == 0) continue;
1024  HMdcLayerGeomParMod& fLayerParMod = (*pMdcLGeomPar)[s][m];
1025  if(!useMdcShParOnly) for(Int_t l=0;l<6;l++) {
1026  if(shiftSysLay[s][m][l] != 1 ) continue; //!!!!!!!!!!!
1027  // Tol'ko dlya Y seychas!!!!!!!!!!!!!!!!!!!!!!!!
1028 //V kakoy posledovatel'nosti????
1029  HMdcLayerGeomParLay& fLayerParLay = fLayerParMod[l];
1030  Double_t pitch = fLayerParLay.getPitch();
1031  Double_t centWrOld = fLayerParLay.getCentralWireNr();
1032  Double_t wireOrOld = fLayerParLay.getWireOrient();
1033  Double_t wireOrNew = wireOrOld + shiftsLay[s][m][l][5]; // Layer rotation
1034  Double_t wireOrR = wireOrNew*TMath::DegToRad();
1035  Double_t layShift = shiftsLay[s][m][l][1]*TMath::Cos(wireOrR) -
1036  shiftsLay[s][m][l][0]*TMath::Sin(wireOrR); // Layer shift
1037  Double_t centWrNew = centWrOld - layShift/pitch;
1038  fLayerParLay.setWireOrient(wireOrNew);
1039  fLayerParLay.setCentralWireNr(centWrNew);
1040 
1041  // ------- Layer second part: -------
1042  if(m<2 || !doLayP2Align) continue;
1043  Int_t firstWrPII;
1044  Float_t layShiftPIIold, wrOrCorrPIIold;
1045  if(!pMdcLCorrPar->getLayerCorrPar(s,m,l, firstWrPII,layShiftPIIold,wrOrCorrPIIold)) continue;
1046  if(firstWrPII>300) pMdcLCorrPar->addLayerShift(s,m,l,firstWrPII,0.,0.);
1047  else {
1048 // Double_t wireOrPIIold = wireOrOld - wrOrCorrPIIold;
1049 // Double_t wireOrPIInew = wireOrPIIold + shiftsLay[s][m][l][8]; // Layer PII rotation
1050 //
1051 // Double_t yShiftPIIold = layShiftPIIold/TMath::Cos(wireOrPIIold);
1052 // Double_t yShiftPIInew = yShiftPIIold + shiftsLay[s][m][l][7]; // Layer PII y-shift
1053 //
1054 // Float_t wrOrCorrPIInew = wireOrNew - wireOrPIInew;
1055 // Float_t layShiftPIInew = yShiftPIInew*TMath::Cos(wireOrPIInew*TMath::DegToRad()) - layShift;
1056  Float_t wrOrCorrPIInew = wrOrCorrPIIold + shiftsLay[s][m][l][5] - shiftsLay[s][m][l][8];
1057  Double_t cosWrOrPIInew = TMath::Cos((wireOrNew - wrOrCorrPIInew)*TMath::DegToRad());
1058  Float_t layShiftPIInew = layShiftPIIold + shiftsLay[s][m][l][7]*cosWrOrPIInew - layShift;
1059 
1060  pMdcLCorrPar->addLayerShift(s,m,l,firstWrPII,layShiftPIInew,wrOrCorrPIInew);
1061  }
1062  }
1063  if(shiftSysMdc[s][m]<0) continue;
1064  // Cell thickness:
1065  if(shiftsMdc[s][m][6] != 0.) {
1066  for(Int_t l=0;l<6;l++) {
1067  HMdcLayerGeomParLay& fLayerParLay = fLayerParMod[l];
1068  Float_t catDist = fLayerParLay.getCatDist();
1069  fLayerParLay.setCatDist(catDist + shiftsMdc[s][m][6]);
1070  }
1071  reCalcLayerTransf = kTRUE;
1072  }
1073 
1074  // Shift MDC:
1075  HMdcSizesCells::setTransform(shiftsMdc[s][m],shift); //first 6 par. are used
1076  HModGeomPar* fMod = pGetCont->getModGeomPar(s,m);
1077  if(shiftSysMdc[s][m] == 0) { // shifts in sector coor.system
1078  if( !pGetCont->getSecTransMod(trans,s,m) ) return kFALSE;
1079  trans.transFrom(shift);
1080  trans.transFrom(pGetCont->getLabTransSec(s));
1081  fMod->getLabTransform() = trans;
1082  } else if(shiftSysMdc[s][m] == 1) { // shifts in mdc coor.system
1083  if(!pGetCont->getLabTransMod(trans,s,m)) return kFALSE;
1084  shift.transFrom(trans);
1085  fMod->getLabTransform() = shift;
1086  } else if(shiftSysMdc[s][m] == 2) { // shifts in lab.coor.system
1087  if( !pGetCont->getLabTransMod(trans,s,m) ) return kFALSE;
1088  trans.transFrom(shift);
1089  fMod->getLabTransform() = trans;
1090  }
1091  }
1092  if(reCalcLayerTransf) pMdcLGeomPar->calcLayerTransformations();
1093  if(isGeomChanged) printf("\nMDC geometry has changed !!!\n\n");
1094  return kTRUE;
1095 }
1096 
1098  TString outGeomFileName(geomParFileName + ".last");
1099  HParAsciiFileIo* outputFile = new HParAsciiFileIo;
1100  if( !outputFile->open(const_cast<char*>(outGeomFileName.Data()),"out") )
1101  return kFALSE;
1102  HRuntimeDb* rtdb = gHades->getRuntimeDb();
1103  if(rtdb == 0) return kFALSE;
1104  if(rtdb->getContainer("MdcGeomStruct")->write(outputFile)<0) return kFALSE;
1105  if(rtdb->getContainer("SpecGeomPar")->write(outputFile)<0) return kFALSE;
1106  if(rtdb->getContainer("MdcGeomPar")->write(outputFile)<0) return kFALSE;
1107  if(rtdb->getContainer("MdcLayerGeomPar")->write(outputFile)<0) return kFALSE;
1108  if(rtdb->getContainer("MdcLayerCorrPar")->write(outputFile)<0) return kFALSE;
1109 // if(rtdb->getContainer("MdcCalParRaw")->write(outputFile)<0) return kFALSE;
1110  outputFile->close();
1111  delete outputFile;
1112  if(doCopyGeomFile) copyFile("cp",outGeomFileName);
1113  if(!isGeomFileExist && !isGeomChanged) {
1114  TString cpFile("cp "+outGeomFileName+" "+geomParFileName);
1115  gSystem->Exec(cpFile);
1116  }
1117  return kTRUE;
1118 }
1119 
1120 void HMdcMille::fixFullMod(Int_t s,Int_t m) {
1121  for(Int_t p=0;p<7;p++) fixModPar(s,m,p);
1122 }
1123 
1124 void HMdcMille::fixModPar(Int_t s,Int_t m,Int_t p) {
1125  if(s>=0&&s<6 && m>=0&&m<4 && p>=0&&p<7) {
1126  mFixFlag[s][m][p] = kTRUE;
1127  if(sigmaMdc[s][m][p] == 0.) sigmaMdc[s][m][p] = -1.;
1128  }
1129 }
1130 
1131 void HMdcMille::fixFullLay(Int_t s,Int_t m,Int_t l) {
1132  for(Int_t p=0;p<nParMax;p++) fixLayPar(s,m,l,p);
1133 }
1134 
1135 void HMdcMille::fixLayPar(Int_t s,Int_t m,Int_t l,Int_t p) {
1136  if(s<0 || s>5 || m<0 || m>3 || l<0 || l>5 || p<0 || p>8) return;
1137  lFixFlag[s][m][l][p] = kTRUE;
1138  if(sigmaLay[s][m][l][p] == 0.) sigmaLay[s][m][l][p] = -1.;
1139 }
1140 
1141 void HMdcMille::fixPar(Int_t s,Int_t m,Int_t l,Int_t p) {
1142  if(l<0) {
1143  if(p<0) fixFullMod(s,m);
1144  else fixModPar(s,m,p);
1145  } else {
1146  if(p<0) fixFullLay(s,m,l);
1147  else fixLayPar(s,m,l,p);
1148  }
1149 }
1150 
1152  // Open file:
1153  FILE* file = fopen(sumOfShiftsFName.Data(),"r");
1154  if(file == NULL) return;
1155 
1156  // Read file:
1157  TString buffer;
1158  if( !buffer.Gets(file) ) printErrorAndExit("readSumOfShiftsFile",
1159  "Wrong format of file %s",sumOfShiftsFName.Data()); // Stop!
1160  Ssiz_t pos = buffer.Index("SumOfShifts:");
1161  if(pos<0 || pos>3) printErrorAndExit("readSumOfShiftsFile",
1162  "Wrong format of file %s",sumOfShiftsFName.Data()); // Stop!
1163 
1164  Int_t sys;
1165  Int_t n = sscanf(buffer.Data(),"%*s %i %*s %*s %i",&iteration,&sys);
1166  if(n!=2 || iteration<0) printErrorAndExit("readSumOfShiftsFile",
1167  "Wrong format of file %s",sumOfShiftsFName.Data()); // Stop!
1168  iteration++;
1169  if( !buffer.Gets(file) ) printErrorAndExit("readSumOfShiftsFile",
1170  "Wrong format of file %s",sumOfShiftsFName.Data()); // Stop!
1171  pos = buffer.Index("Statistics:");
1172  if(pos<0 || pos>3) printErrorAndExit("readSumOfShiftsFile",
1173  "Wrong format of file %s",sumOfShiftsFName.Data()); // Stop!
1174 
1175  while(buffer.Gets(file) && buffer.Length()>16 && addShifts(buffer));
1176 
1177  fclose(file);
1178 }
1179 
1180 Bool_t HMdcMille::addShifts(TString& buffer) {
1181  TObjArray* arr = buffer.Tokenize(" ");
1182  Bool_t exitFlag = kFALSE;
1183  if(arr->GetEntries()>=3) {
1184  Int_t label = ((TObjString*)(arr->At(0)))->GetString().Atoi();
1185  Double_t shift = ((TObjString*)(arr->At(1)))->GetString().Atof();
1186  Double_t coef = ((TObjString*)(arr->At(2)))->GetString().Atof();
1187  Int_t btId,sys,sec,mod,lay,parn;
1188  if(!unpackLabel(label,btId,sys,sec,mod,lay,parn))
1189  printErrorAndExit("addShifts","Wrong format of input file!!! Stop!");
1190  exitFlag = kTRUE;
1191  if((btId == 0 || btId == beamTimeId) && pMdcDet->getModule(sec,mod)) {
1192  mdcMods[sec][mod] = 1;
1193  if(lay < 0) {
1194  shiftsMdc[sec][mod][parn] += shift;
1195  if(coef >= 0.) shiftSysMdc[sec][mod] = sys; // Otherwise = -1
1196  if(shift != 0.) isGeomChanged = kTRUE;
1197  mShFlag[sec][mod][parn] = kTRUE;
1198  sigmaMdc[sec][mod][parn] = coef;
1199  } else {
1200  shiftsLay[sec][mod][lay][parn] += shift;
1201  if(coef >= 0.) shiftSysLay[sec][mod][lay] = sys;
1202  if(shift != 0.) isGeomChanged = kTRUE;
1203  lShFlag[sec][mod][lay][parn] = kTRUE;
1204  sigmaLay[sec][mod][lay][parn] = coef;
1205  }
1206  }
1207  }
1208  arr->Delete();
1209  delete arr;
1210  return exitFlag;
1211 }
1212 
1213 void HMdcMille::printErrorAndExit(const char* func,const char* form,const char* str) {
1214 #warning the following warning because of string literal can be ignored!
1215  if(str == NULL) Error(func,form);
1216  else Error(func,form,str);
1217  exit(1);
1218 }
1219 
1221  FILE* file = fopen(sumOfShiftsFName.Data(),"w");
1222  if(file == NULL) printErrorAndExit("creatSumOfShiftsFile",
1223  "Can't open file %s",sumOfShiftsFName.Data());
1224  TString sys;
1225  if(shiftType == 0) sys = "all shifts done in sector sys.";
1226  if(shiftType == 1) sys = "all shifts done in mdc sys.";
1227  if(shiftType == 2) sys = "all shifts done in lab. sys.";
1228  Double_t mNWr = Double_t(nWiresTot)/nTracks;
1229  Double_t mNLy = Double_t(nLayersTot)/nTracks;
1230  fprintf(file,"SumOfShifts: %i iteration. System %i - %s\n",iteration,shiftType,sys.Data());
1231  fprintf(file,"Statistics: %i tracks; %8.4f <wires/tr.>; %8.4f <layers/tr.>\n",nTracks,mNWr,mNLy);
1232  for(Int_t s=0;s<6;s++) for(Int_t m=0;m<4;m++) for(Int_t l=-1;l<6;l++) {
1233  if(l>=0 && useMdcShParOnly) continue;
1234  for(Int_t p=0;p<nParMax;p++) {
1235  if(l<0 && !mShFlag[s][m][p]) continue;
1236  if(l>=0 && !lShFlag[s][m][l][p]) continue;
1237 
1238  fprintf(file,"%11i",packLabel(s,m,l,p)); // Print label
1239  if(l < 0) fprintf(file," %18.10e %9.4e ! %is. %im.: mdc",
1240  shiftsMdc[s][m][p],sigmaMdc[s][m][p],s+1,m+1);
1241  else {
1242  fprintf(file," %18.10e %9.4e ! %is. %im. %il.: layer",
1243  shiftsLay[s][m][l][p],sigmaLay[s][m][l][p],s+1,m+1,l+1);
1244  if(p==7 || p==8) fprintf(file," part II");
1245  }
1246  fprintf(file," %s\n",shitsInfo[p].Data()); // Print description
1247  }
1248  }
1249  fclose(file);
1251  if(doCopyResFile) copyFile("cp",pedeResFName,".res");
1252  if(doCopyHisFile) copyFile("mv",pedeResFName,".his");
1253  if(doCopyLogFile) copyFile("mv",pedeResFName,".log");
1254 }
1255 
1256 void HMdcMille::copyFile(const char* op,TString& file,const char* ext) {
1257  TString cpFile(op);
1258  if(ext) cpFile += " "+file+ext+" "+file+ext+stepIter;
1259  else cpFile += " "+file+" "+file+stepIter;
1260  gSystem->Exec(cpFile);
1261 }
1262 
1265  if(pedeTaskFileName.Length()<=0) printErrorAndExit("creatPedeTaskFile",
1266  "File name of pede task list is not specified!");
1267  crPedeTaskFile = kTRUE;
1268  // Set default method:
1269  method = "inversion ";
1270  numOfIter = 3;
1271  accuracy = 0.001;
1272  mthDescr = "Gauss matrix inversion";
1273  bandwidth = -1;
1274 }
1275 
1276 void HMdcMille::useInversionMethod(Int_t nit,Float_t acc) {
1277  method = "inversion ";
1278  numOfIter = nit;
1279  accuracy = acc;
1280  mthDescr = "Gauss matrix inversion";
1281  bandwidth = -1;
1282 }
1283 
1284 void HMdcMille::useBandCholeskyMethod(Int_t nit,Float_t acc,Int_t bwidth) {
1285  method = "bandcholesky ";
1286  numOfIter = nit;
1287  accuracy = acc;
1288  mthDescr = "Cholesky";
1289  bandwidth = bwidth;
1290 }
1291 
1292 void HMdcMille::useCholeskyMethod(Int_t nit,Float_t acc,Int_t bwidth) {
1293  method = "cholesky ";
1294  numOfIter = nit;
1295  accuracy = acc;
1296  mthDescr = "Cholesky";
1297  bandwidth = bwidth;
1298 }
1299 
1300 void HMdcMille::useSparseGMRESMethod(Int_t nit,Float_t acc) {
1301  method = "sparseGMRES ";
1302  numOfIter = nit;
1303  accuracy = acc;
1304  mthDescr = "minimal residual";
1305  bandwidth = -1;
1306 }
1307 
1308 void HMdcMille::useFullGMRESMethod(Int_t nit,Float_t acc) {
1309  method = "fullGMRES ";
1310  numOfIter = nit;
1311  accuracy = acc;
1312  mthDescr = "minimal residual";
1313  bandwidth = -1;
1314 }
1315 
1316 void HMdcMille::useDiagonalizationMethod(Int_t nit,Float_t acc) {
1317  method = "diagonalization";
1318  numOfIter = nit;
1319  accuracy = acc;
1320  mthDescr = "diagonalization";
1321  bandwidth = -1;
1322 }
1323 
1324 void HMdcMille::setParConstrainFile(const char* file) {
1325  parConstrainFile = file;
1326 }
1327 
1329  if(!crPedeTaskFile) return;
1330  if(pedeInParFName.Length() <= 0) printErrorAndExit("writePedeTaskFile",
1331  "File name of parameters list not specified!");
1332  FILE* file = fopen(pedeTaskFileName.Data(),"w");
1333  if(file == NULL) printErrorAndExit("writePedeTaskFile","Can't open file %s",pedeTaskFileName.Data());
1334  fprintf(file,"Cfiles !\n");
1335 
1336  if(parConstrainFile.Length() > 0)
1337  fprintf(file,"%s ! constraints text file\n",parConstrainFile.Data());
1338  else fprintf(file,"!parConstrain.txt ! constraints text file\n");
1339 
1340  fprintf(file,"%s ! parameter text file\n",pedeInParFName.Data());
1341 
1342  for(Int_t nf=0;nf<nBinaryFiles;nf++) fprintf(file,
1343  "%s%i ! binary data file\n",milleFileName.Data(),nf);
1344 
1345  if(bandwidth > 0) fprintf(file,"bandwidth %i ! width of band matrix\n",bandwidth);
1346  fprintf(file,"method %s %i %g ! %s\n",method.Data(),numOfIter,accuracy,mthDescr.Data());
1347  fclose(file);
1348 }
1349 
1351  for(Int_t s=0;s<6;s++) nWiresCut[s] = wc[s];
1352 }
1353 
1354 void HMdcMille::setMaxNumWiresCut(Int_t cs1,Int_t cs2,Int_t cs3,Int_t cs4,Int_t cs5,Int_t cs6) {
1355  nWiresCut[0] = cs1;
1356  nWiresCut[1] = cs2;
1357  nWiresCut[2] = cs3;
1358  nWiresCut[3] = cs4;
1359  nWiresCut[4] = cs5;
1360  nWiresCut[5] = cs6;
1361 }
1362 
1364  useMdcShParOnly = fl;
1365  if(useMdcShParOnly) {
1366  for(Int_t s=0;s<6;s++) for(Int_t m=0;m<4;m++) for(Int_t p=6;p<9;p++) mShFlag[s][m][p] = kFALSE;
1367  }
1368 }
static HMdcEvntListCells * getExObject(void)
HMdcList24GroupCells & getOutputListCells3(void)
Definition: hmdcwiredata.h:331
HGeomTransform * layPosShCurr
Definition: hmdcmille.h:81
void setMaxNumWiresCut(Int_t *wc)
Definition: hmdcmille.cc:1350
Bool_t init(void)
Definition: hmdcmille.cc:169
Int_t nMods[6]
Definition: hmdcmille.h:53
Double_t z2(void) const
Definition: hmdcgeomobj.h:255
Bool_t finalize(void)
Definition: hmdcmille.cc:370
void fixPar(Int_t s, Int_t m, Int_t l, Int_t p)
Definition: hmdcmille.cc:1141
Int_t getNumOfGoodWires(void) const
Bool_t doCopyLogFile
Definition: hmdcmille.h:35
Double_t sigmaMdc[6][4][9]
Definition: hmdcmille.h:104
HGeomTransform negShifts[ASIZE]
Definition: hmdcmille.h:61
UChar_t getSec(void) const
Definition: hmdcclus.h:135
void readSumOfShiftsFile(void)
Definition: hmdcmille.cc:1151
const HGeomTransform & getLabTransSec(Int_t sec, Bool_t init=kFALSE)
static void setAnotherFit(HMdc12Fit *fit)
Int_t getNLayers(void) const
Double_t y2(void) const
Definition: hmdcgeomobj.h:254
void fixLayPar(Int_t s, Int_t m, Int_t l, Int_t p)
Definition: hmdcmille.cc:1135
TString pedeResFName
Definition: hmdcmille.h:93
const HGeomTransform * getRotLayP2SysRSec(void) const
Int_t getGlobalDerWTof(void)
Definition: hmdcmille.cc:807
void transLineToOtherSec(const HMdcLineParam &ln, Int_t sec, HGeomVector &p1, HGeomVector &p2)
Double_t sigmaLay[6][4][6][9]
Definition: hmdcmille.h:107
Double_t getWtNorm(void) const
Definition: hmdcwiredata.h:155
static HMdcGetContainers * getObject()
static Bool_t fPrint(void)
void creatPedeInParamFile(void)
Definition: hmdcmille.cc:916
Bool_t doCopyHisFile
Definition: hmdcmille.h:36
HMdcTrackParam * getFinalParam(void)
Double_t thetaCut
Definition: hmdcmille.h:55
TString milleFileName
Definition: hmdcmille.h:41
HMdcWireData * getWireData(Int_t n)
Definition: hmdcwiredata.h:304
Bool_t getFitStatus(void) const
void setTrShiftInModSys(void)
Definition: hmdcmille.cc:270
Int_t sector
Definition: hmdcmille.h:73
TString mthDescr
Definition: hmdcmille.h:127
void setShiftTransformation(void)
Definition: hmdcmille.cc:247
UChar_t getIOSeg(void) const
Definition: hmdcclus.h:136
HCategory * fClusCat
Definition: hmdc12fit.h:54
Int_t bandwidth
Definition: hmdcmille.h:128
Double_t calcDriftTimeForAlign(const HGeomVector &p1, const HGeomVector &p2, Int_t &distSign)
Int_t shiftType
Definition: hmdcmille.h:23
HMdcDetector * getMdcDetector()
const HGeomTransform * getRotLayP1SysRMod(void) const
UInt_t milleFileSize
Definition: hmdcmille.h:43
void setDef(void)
Definition: hmdcmille.cc:74
Bool_t mShFlag[6][4][9]
Definition: hmdcmille.h:112
Double_t derLThNorm[6]
Definition: hmdcmille.h:64
Bool_t getPrintFlag(void) const
Int_t getLayerNParts(void) const
Int_t iMod
Definition: hmdcmille.h:77
int globLabel[ASIZE]
Definition: hmdcmille.h:90
Bool_t doCopyResFile
Definition: hmdcmille.h:34
Int_t getGlobalDer(void)
Definition: hmdcmille.cc:768
Int_t getSector3(void) const
Definition: hmdcwiredata.h:339
Bool_t doCopyGeomFile
Definition: hmdcmille.h:32
Bool_t makeShifts(void)
Definition: hmdcmille.cc:1014
void setPrintFlag(Bool_t prnt)
Int_t getNDrTimesMod(Int_t m) const
Bool_t lFixFlag[6][4][6][9]
Definition: hmdcmille.h:115
Bool_t printDebugFlag
Definition: hmdcmille.h:131
Int_t nWiresMinTr
Definition: hmdcmille.h:54
HRuntimeDb * getRuntimeDb(void)
Definition: hades.h:111
HMdcTrackFitter * fitter
Definition: hmdc12fit.h:51
const Int_t kSkipEvent
Definition: haddef.h:30
TString pedeInParFName
Definition: hmdcmille.h:111
Int_t iteration
Definition: hmdcmille.h:94
Double_t z1(void) const
Definition: hmdcgeomobj.h:252
Bool_t isGeomFileExist
Definition: hmdcmille.h:98
Mille * mille
Definition: hmdcmille.h:22
Int_t nParMax
Definition: hmdcmille.h:102
virtual Int_t getModule(Int_t sector, Int_t mod)
Definition: hdetector.cc:85
Int_t nLayers[6][4]
Definition: mdcEfficiency.C:37
Double_t getChi2(void) const
Bool_t lShFlag[6][4][6][9]
Definition: hmdcmille.h:113
const HGeomTransform & getLabTransMod(Int_t sec, Int_t mod, Bool_t init=kFALSE)
HMdcGetContainers * pGetCont
Definition: hmdcmille.h:39
void setMdcShiftParOnly(Bool_t fl=kTRUE)
Definition: hmdcmille.cc:1363
Bool_t fitSeg(HMdcClus *fClst, Int_t arrInd=0)
Definition: hmdc12fit.cc:677
HGeomTransform negLThSh[6]
Definition: hmdcmille.h:63
void calcShiftInSecSys(const HGeomTransform *laySec, const HGeomTransform &shift, HGeomTransform &laySh)
Definition: hmdcmille.cc:349
Int_t n
void setInputVersion(Int_t v=-1, Int_t i=0)
Definition: hparset.h:31
static HMdcSizesCells * getExObject(void)
static HMdcSizesCells * getObject(void)
void setILay(Int_t l)
Definition: hmdcmille.h:190
Bool_t writeParAsciiFile(void)
Definition: hmdcmille.cc:1097
Bool_t fitNSectors(HMdcClus *cl1, HMdcClus *cl2, HMdcClus *cl3, HMdcClus *cl4, HMdcClus *cl5, HMdcClus *cl6, HMdcClus *cl7, HMdcClus *cl8)
Definition: hmdc12fit.cc:797
void useDiagonalizationMethod(Int_t nit=5, Float_t acc=0.01)
Definition: hmdcmille.cc:1316
Double_t getThetaDeg(void) const
Definition: hmdcgeomobj.h:261
void readPedeResFile(const char *fileName)
Definition: hmdcmille.cc:992
Int_t getNDrTimes(void) const
float locDer[4]
Definition: hmdcmille.h:88
void setParConstrainFile(const char *file)
Definition: hmdcmille.cc:1324
HMdcWiresArr & getWiresArr(void)
Int_t getSector(void) const
Definition: hmdcwiredata.h:307
HModGeomPar * getModGeomPar(Int_t sec, Int_t mod, Bool_t init=kFALSE)
Bool_t open(const Text_t *fname, const Text_t *status="in")
Int_t getLayerPart(Int_t c) const
Int_t getEventCounter()
Definition: hades.h:96
int end()
Definition: Mille.cc:131
HMdcList24GroupCells & getOutputListCells2(void)
Definition: hmdcwiredata.h:330
Bool_t cellThicknFree
Definition: hmdcmille.h:44
static Bool_t unpackLabel(Int_t label, Int_t &btId, Int_t &sys, Int_t &sec, Int_t &mod, Int_t &lay, Int_t &parn)
Definition: hmdcmille.cc:854
HMdcWireData * wireData
Definition: hmdcmille.h:72
Bool_t addLayerShift(Int_t s, Int_t m, Int_t l, Int_t fstWr, Float_t sh, Float_t orCorr=0.)
HGeomTransform layPosSh[144][ASIZE]
Definition: hmdcmille.h:65
HMdcLayerCorrPar * getMdcLayerCorrPar(Bool_t init=kFALSE)
Double_t x2(void) const
Definition: hmdcgeomobj.h:253
virtual Int_t write()
Definition: hparset.cc:119
void fixModPar(Int_t s, Int_t m, Int_t p)
Definition: hmdcmille.cc:1124
Double_t x1(void) const
Definition: hmdcgeomobj.h:250
HMdcList24GroupCells & getOutputListCells(void)
Definition: hmdcwiredata.h:329
HMdcTrackFitInOut fitpar
Definition: hmdc12fit.h:49
void setIMod(Int_t s, Int_t m)
Definition: hmdcmille.h:189
Bool_t reinit(void)
Definition: hmdc12fit.cc:249
Int_t getNLayersMod(Int_t mod) const
Int_t beamTimeId
Definition: hmdcmille.h:38
HIterator * iterClus
Definition: hmdc12fit.h:55
Int_t layerPart
Definition: hmdcmille.h:79
void useCholeskyMethod(Int_t nit=5, Float_t acc=0.001, Int_t bwidth=6)
Definition: hmdcmille.cc:1292
HGeomVector p2
Definition: hmdcmille.h:85
Int_t module
Definition: hmdcmille.h:74
Bool_t init(void)
Definition: hmdc12fit.cc:183
Double_t func(Double_t *x, Double_t *par)
Bool_t isGeomChanged
Definition: hmdcmille.h:99
Int_t getTypeClFinder(void) const
Definition: hmdcclus.h:134
Bool_t fit2Sectors(HMdcClus *fClst1, HMdcClus *fClst2, HMdcClus *fClst3, HMdcClus *fClst4)
Definition: hmdc12fit.cc:758
Int_t getDistanceSign(void) const
Definition: hmdcwiredata.h:149
UInt_t nTracks
Definition: hmdcmille.h:29
TString pedeTaskFileName
Definition: hmdcmille.h:123
Int_t getSector4(void) const
Definition: hmdcwiredata.h:340
void useSparseGMRESMethod(Int_t nit=5, Float_t acc=0.01)
Definition: hmdcmille.cc:1300
TString fileName("be1834507002801.hld")
Bool_t getLayerCorrPar(Int_t s, Int_t m, Int_t l, Int_t &fstWr, Float_t &sh, Float_t &orCorr) const
Double_t y1(void) const
Definition: hmdcgeomobj.h:251
Int_t layer
Definition: hmdcmille.h:75
void creatSumOfShiftsFile(void)
Definition: hmdcmille.cc:1220
HMdcSizesCellsLayer * pSCLayer[144]
Definition: hmdcmille.h:59
void setXYZ(const Double_t xx, const Double_t yy, const Double_t zz)
Definition: hgeomvector.h:25
void setTrShiftInSecSys(void)
Definition: hmdcmille.cc:293
Float_t accuracy
Definition: hmdcmille.h:126
Double_t tofDerCorr[24][9]
Definition: hmdcmille.h:116
Definition: hparset.h:9
void setChi2PerNdfCut(Double_t cut=50.)
void setCatDist(Float_t c)
Bool_t getSecTransMod(HGeomTransform &trans, Int_t sec, Int_t mod, Bool_t init=kFALSE)
void setTrShiftInLabSys(void)
Definition: hmdcmille.cc:314
Int_t iSec
Definition: hmdcmille.h:76
Hades * gHades
Definition: hades.cc:1213
~HMdcMille(void)
Definition: hmdcmille.cc:165
ClassImp(HMdcMille) HMdcMille
Definition: hmdcmille.cc:42
Double_t getSumWtNorm(Int_t m) const
void writePedeTaskFile(void)
Definition: hmdcmille.cc:1328
static void deleteCont()
Bool_t getAnalytDeriv(Float_t *der, HMdcTrackParam *par=0)
HMdcList24GroupCells & getOutputListCells4(void)
Definition: hmdcwiredata.h:332
Int_t numOfIter
Definition: hmdcmille.h:125
Bool_t doLayP2Align
Definition: hmdcmille.h:80
HMdcSizesCellsMod * pSCMod[24]
Definition: hmdcmille.h:58
Bool_t mFixFlag[6][4][9]
Definition: hmdcmille.h:114
void getIndexRegChilds(Int_t &first, Int_t &last) const
Definition: hmdcclus.h:160
Bool_t fixFitTOffset
Definition: hmdcmille.h:26
Int_t getModule(void) const
Definition: hmdcwiredata.h:136
Int_t getModIndex(void) const
Definition: hmdcwiredata.h:137
Double_t getDrTime(void) const
Definition: hmdcwiredata.h:145
void mille(int NLC, const float *derLc, int NGL, const float *derGl, const int *label, float rMeas, float sigma)
Definition: Mille.cc:44
Int_t execute(void)
Definition: hmdcmille.cc:380
void fitAlgorithmForMille(void)
Definition: hmdcmille.cc:428
HGeomTransform layNegSh[144][ASIZE]
Definition: hmdcmille.h:66
void printErrorAndExit(const char *func, const char *form, const char *str=NULL)
Definition: hmdcmille.cc:1213
HParSet * getContainer(const Text_t *)
Definition: hruntimedb.cc:124
void setCentralWireNr(Float_t e)
HGeomVector transFrom(const HGeomVector &p) const
Int_t iLay
Definition: hmdcmille.h:78
Bool_t eventFilter(void)
Definition: hmdcmille.cc:397
const HGeomTransform * getRotLayP1SysRSec(void) const
void fixFullMod(Int_t s, Int_t m)
Definition: hmdcmille.cc:1120
UChar_t filteringFlag
Definition: hmdcmille.h:48
HGeomTransform posShifts[ASIZE]
Definition: hmdcmille.h:60
Double_t shiftsMdc[6][4][9]
Definition: hmdcmille.h:103
void creatPedeTaskFile(const char *fileName="pedeTask.txt")
Definition: hmdcmille.cc:1263
Int_t histInd[6]
Definition: hmdcmille.h:119
Char_t getHitStatus(void) const
Definition: hmdcwiredata.h:153
void fixFullLay(Int_t s, Int_t m, Int_t l)
Definition: hmdcmille.cc:1131
Bool_t useMdcShParOnly
Definition: hmdcmille.h:37
Char_t fname[]
Definition: drawAccCuts.C:3
Bool_t useSector[6]
Definition: hmdcmille.h:45
Bool_t crPedeTaskFile
Definition: hmdcmille.h:122
TString sumOfShiftsFName
Definition: hmdcmille.h:96
void calcShiftInLabSys(const HGeomTransform *secSys, const HGeomTransform *laySysMod, const HGeomTransform &shift, HGeomTransform &laySh)
Definition: hmdcmille.cc:358
Int_t mdcMods[6][4]
Definition: hmdcmille.h:100
Int_t getSector(void) const
Definition: hmdcwiredata.h:135
void openNewBinaryFile(void)
Definition: hmdcmille.cc:744
Double_t getPhiDeg(void) const
Definition: hmdcgeomobj.h:262
Bool_t doCopySumShFile
Definition: hmdcmille.h:33
void sendToMille(void)
Definition: hmdcmille.cc:650
HMdcTrackFitter * getFitter(Int_t ind)
Definition: hmdc12fit.h:37
TString stepIter
Definition: hmdcmille.h:95
Bool_t fitMixedClus(HMdcClus *fClst, Int_t arrInd=0)
Definition: hmdc12fit.cc:916
void setStatic(Bool_t flag=kTRUE)
Definition: hparset.h:38
void createHists(Int_t sec)
Definition: hmdcmille.cc:217
Double_t calcDDriftTime(Int_t p)
Definition: hmdcmille.cc:884
Bool_t reinit(void)
Definition: hmdcmille.cc:173
HMdcTrackParam * finalParam
Definition: hmdcmille.h:83
HMdcFittersArray fittersArr[2]
Definition: hmdc12fit.h:50
static void setTransform(Double_t *par, HGeomTransform &trans)
const HGeomTransform * getSecTrans(void) const
void fitAlgForMilleCosmic(void)
Definition: hmdcmille.cc:495
Double_t parSteps[ASIZE]
Definition: hmdcmille.h:27
Bool_t isInCell(void) const
Definition: hmdcwiredata.h:151
Int_t getSector2(void) const
Definition: hmdcwiredata.h:338
UInt_t nLayersTot
Definition: hmdcmille.h:31
int packLabel(Int_t s, Int_t m, Int_t l, Int_t parn)
Definition: hmdcmille.cc:842
Int_t getIndexParent(void) const
Definition: hmdcclus.h:159
Int_t getOutputNLayers(void) const
Int_t nWiresMin
Definition: hmdcmille.h:52
Int_t nBinaryFiles
Definition: hmdcmille.h:42
void copyFile(const char *op, TString &file, const char *ext=0)
Definition: hmdcmille.cc:1256
TString parConstrainFile
Definition: hmdcmille.h:129
Double_t getErrTdcTime(void) const
Definition: hmdcwiredata.h:144
HMdcLayerGeomPar * getMdcLayerGeomPar(Bool_t init=kFALSE)
Int_t getLayer(void) const
Definition: hmdcwiredata.h:138
Bool_t doHists
Definition: hmdcmille.h:132
Definition: Mille.h:26
HGeomTransform & getLabTransform()
Definition: hdetgeompar.h:21
void useFullGMRESMethod(Int_t nit=5, Float_t acc=0.01)
Definition: hmdcmille.cc:1308
HGeomTransform posLThSh[6]
Definition: hmdcmille.h:62
Double_t getDev(void) const
Definition: hmdcwiredata.h:146
HGeomVector p1
Definition: hmdcmille.h:84
Double_t calcGlobDer(Int_t parNum)
Definition: hmdcmille.cc:878
Int_t getActiveModule(void) const
Bool_t loadGeometryPar(void)
Definition: hmdcmille.cc:956
TString geomParFileName
Definition: hmdcmille.h:97
void reset(void)
Definition: hmdc12fit.h:39
HMdcDetector * pMdcDet
Definition: hmdcmille.h:40
float globDer[ASIZE]
Definition: hmdcmille.h:89
UInt_t nWiresTot
Definition: hmdcmille.h:30
Int_t clusFitAlg
Definition: hmdc12fit.h:82
const HGeomTransform * getRotLayP2SysRMod(void) const
void calcShiftInModSys(const HGeomTransform *mdcSys, const HGeomTransform *laySysMod, const HGeomTransform &shift, HGeomTransform &laySh)
Definition: hmdcmille.cc:338
Double_t shiftsLay[6][4][6][9]
Definition: hmdcmille.h:106
void useInversionMethod(Int_t nit=3, Float_t acc=0.001)
Definition: hmdcmille.cc:1276
TString shitsInfo[9]
Definition: hmdcmille.h:67
Int_t getSegIndex(void) const
Double_t calcDerTofCorr(Int_t p, Double_t *dDrTm=0)
Definition: hmdcmille.cc:833
UChar_t nLayerParts[144]
Definition: hmdcmille.h:68
void setWireOrient(Float_t d)
Bool_t addShifts(TString &buffer)
Definition: hmdcmille.cc:1180
Bool_t fitSec(HMdcClus *fClst1, HMdcClus *fClst2)
Definition: hmdc12fit.cc:717
Int_t nHists
Definition: hmdcmille.h:118
Int_t shiftSysMdc[6][4]
Definition: hmdcmille.h:101
TList * fHistograms
HGeomVector transTo(const HGeomVector &p) const
UInt_t nClustersCut
Definition: hmdcmille.h:46
Int_t nWiresCut[6]
Definition: hmdcmille.h:47
void useBandCholeskyMethod(Int_t nit=8, Float_t acc=0.01, Int_t bwidth=6)
Definition: hmdcmille.cc:1284
Int_t getSec(void) const
Definition: hmdcgeomobj.h:243
Bool_t createPedeInFile
Definition: hmdcmille.h:110
HGeomTransform * layNegShCurr
Definition: hmdcmille.h:82
Int_t getNDriftTimes(void) const
Definition: hmdcwiredata.h:309
Int_t shiftSysLay[6][4][6]
Definition: hmdcmille.h:105
Double_t calcDriftTime(const HGeomTransform &laySys, Int_t &distSign)
Definition: hmdcmille.cc:895
Int_t getCell(void) const
Definition: hmdcwiredata.h:139
void setPedeInFileName(const char *fname)
Definition: hmdcmille.cc:910
TString method
Definition: hmdcmille.h:124
Bool_t isMdcInAlign[6][4]
Definition: hmdcmille.h:69
void setIModILay(void)
Definition: hmdcmille.cc:754
HMdcTrkCand * fillTrkCandISeg(void)
Definition: hmdc12fit.cc:1282
Double_t derNorm[ASIZE]
Definition: hmdcmille.h:28