ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/CMSSW/Alignment/CommonAlignmentAlgorithm/src/AlignmentParametersIORoot.cc
Revision: 1.3
Committed: Fri Mar 16 16:35:04 2007 UTC (18 years, 1 month ago) by flucke
Content type: text/plain
Branch: MAIN
CVS Tags: CMSSW_1_4_10, CMSSW_1_4_9, CMSSW_1_4_8, CMSSW_1_4_7, CMSSW_1_5_4, CMSSW_1_5_3, CMSSW_1_4_6, CMSSW_1_5_2, CMSSW_1_6_0_pre1, CMSSW_1_5_1, CMSSW_1_4_5, CMSSW_1_5_0, V01-04-00-05, CMSSW_1_4_3g483, CMSSW_1_4_4, CMSSW_1_5_0_pre6, CMSSW_1_4_3, V01-04-00-03, CMSSW_1_5_0_pre5, CMSSW_1_4_2, CMSSW_1_4_1, CMSSW_1_5_0_pre4, CMSSW_1_5_0_pre3, CMSSW_1_4_0_DAQ1, V01-07-01, V01-04-00-02, CMSSW_1_4_0, CMSSW_1_5_0_pre2, CMSSW_1_4_0_pre7, CMSSW_1_4_0_pre6, CMSSW_1_4_0_pre5, V01-07-00, V01-04-00-01, CMSSW_1_5_0_pre1, CMSSW_1_4_0_pre4, V01-06-00, V01-05-00, V01-04-01, CMSSW_1_4_0_pre3, V01-04-00, V01-03-01, V01-03-00
Branch point for: V01-04-00-04
Changes since 1.2: +6 -5 lines
Log Message:
write hierarchy level of parameters

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