ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/CMSSW/Alignment/CommonAlignmentAlgorithm/src/AlignmentParametersIORoot.cc
Revision: 1.5
Committed: Thu Jul 12 15:08:28 2007 UTC (17 years, 9 months ago) by flucke
Content type: text/plain
Branch: MAIN
CVS Tags: V01-04-00-08, V01-04-00-07, V01-04-00-06
Changes since 1.4: +5 -8 lines
Log Message:
adjust to vanishing CompositeRigidBodyAlignmentParameters

File Contents

# User Rev Content
1 flucke 1.2 // this class's header
2     #include "Alignment/CommonAlignmentAlgorithm/interface/AlignmentParametersIORoot.h"
3    
4     #include <TTree.h>
5    
6 flucke 1.5 #include "Alignment/CommonAlignmentParametrization/interface/RigidBodyAlignmentParameters.h"
7 flucke 1.2
8 fronga 1.1 #include "Alignment/TrackerAlignment/interface/TrackerAlignableId.h"
9    
10     // ----------------------------------------------------------------------------
11     // constructor
12    
13     AlignmentParametersIORoot::AlignmentParametersIORoot()
14     {
15     treename = "AlignmentParameters";
16     treetxt = "Alignment Parameters";
17     }
18    
19     // ----------------------------------------------------------------------------
20    
21     void AlignmentParametersIORoot::createBranches(void)
22     {
23     tree->Branch("parSize", &theCovRang, "CovRang/I");
24     tree->Branch("Id", &theId, "Id/I");
25     tree->Branch("Par", &thePar, "Par[CovRang]/D");
26     tree->Branch("covarSize", &theCovarRang, "CovarRang/I");
27     tree->Branch("Cov", &theCov, "Cov[CovarRang]/D");
28     tree->Branch("ObjId", &theObjId, "ObjId/I");
29 flucke 1.3 tree->Branch("HieraLevel",&theHieraLevel,"HieraLevel/I");
30 fronga 1.1 }
31    
32     // ----------------------------------------------------------------------------
33    
34     void AlignmentParametersIORoot::setBranchAddresses(void)
35     {
36     tree->SetBranchAddress("parSize", &theCovRang);
37     tree->SetBranchAddress("covarSize", &theCovarRang);
38     tree->SetBranchAddress("Id", &theId);
39     tree->SetBranchAddress("Par", &thePar);
40     tree->SetBranchAddress("Cov", &theCov);
41     tree->SetBranchAddress("ObjId", &theObjId);
42 flucke 1.3 tree->SetBranchAddress("HieraLevel",&theHieraLevel);
43 fronga 1.1 }
44    
45     // ----------------------------------------------------------------------------
46     // find tree entry based on detID and typeID
47    
48     int AlignmentParametersIORoot::findEntry(unsigned int detId,int comp)
49     {
50     double noAliPar = tree->GetEntries();
51     for (int ev = 0;ev<noAliPar;ev++) {
52     tree->GetEntry(ev);
53     if ( theId==detId && theObjId==comp ) return (ev);
54     }
55     return -1;
56     }
57    
58     // ----------------------------------------------------------------------------
59     int AlignmentParametersIORoot::writeOne(Alignable* ali)
60     {
61 flucke 1.3 const AlignmentParameters* ap =ali->alignmentParameters();
62     const AlgebraicVector& params = ap->parameters();
63     const AlgebraicSymMatrix& cov = ap->covariance();
64 fronga 1.1
65     theCovRang = params.num_row();
66     theCovarRang = theCovRang*(theCovRang+1)/2;
67     int count=0;
68     for(int row=0;row<theCovRang;row++){
69     thePar[row]=params[row];
70     for(int col=0;col<theCovRang;col++){
71     if(row-1<col) { theCov[count] = cov[row][col]; count++; }
72     }
73     }
74    
75     TrackerAlignableId converter;
76     theId = converter.alignableId(ali);
77     theObjId = converter.alignableTypeId(ali);
78 flucke 1.3 theHieraLevel = ap->hierarchyLevel();
79 fronga 1.1
80     tree->Fill();
81     return 0;
82     }
83    
84     // ----------------------------------------------------------------------------
85    
86 flucke 1.2 AlignmentParameters* AlignmentParametersIORoot::readOne( Alignable* ali, int& ierr )
87 fronga 1.1 {
88    
89     TrackerAlignableId converter;
90     int obj = converter.alignableTypeId(ali);
91     unsigned int detInt = converter.alignableId(ali);
92    
93 flucke 1.5 AlignmentParameters* alipar = 0;
94 fronga 1.1 AlgebraicVector par(nParMax,0);
95     AlgebraicSymMatrix cov(nParMax,0);
96 flucke 1.5 const std::vector<bool> &sel = ali->alignmentParameters()->selector();
97 fronga 1.1
98     int entry = findEntry(detInt,obj);
99     if( entry != -1 )
100     {
101     tree->GetEntry(entry);
102     int covsize = theCovRang;
103     int count=0;
104     for(int row=0;row<covsize;row++)
105     {
106     par[row]=thePar[row];
107     for(int col=0; col < covsize;col++) {
108     if(row-1<col) {cov[row][col]=theCov[count];count++;}
109     }
110     }
111 flucke 1.5 // FIXME: In future should check which kind of parameters to construct...
112     alipar = new RigidBodyAlignmentParameters(ali,par,cov,sel);
113 fronga 1.1 alipar->setValid(true);
114     ierr=0;
115     return alipar;
116     }
117    
118     ierr=-1;
119     return(0);
120     }
121    
122