CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HIPUserVariablesIORoot.cc
Go to the documentation of this file.
1 #include "TTree.h"
2 
4 
8 
9 
10 // this class's header
12 
13 // ----------------------------------------------------------------------------
14 // constructor
15 
17 {
18  treename = "T9";
19  treetxt = "HIP User Variables";
20 }
21 
22 // ----------------------------------------------------------------------------
23 
25 {
26  tree->Branch("Id", &Id, "Id/i");
27  tree->Branch("ObjId", &ObjId, "ObjId/I");
28 
29  tree->Branch("Nhit", &Nhit, "Nhit/I");
30  tree->Branch("Nparj", &Nparj, "Nparj/I");
31  tree->Branch("Jtvj", &Jtvj, "Jtvj[Nparj]/D");
32  tree->Branch("Npare", &Npare, "Npare/I");
33  tree->Branch("Jtve", &Jtve, "Jtve[Npare]/D");
34  tree->Branch("AlignableChi2", &AlignableChi2, "AlignableChi2/D");
35  tree->Branch("AlignableNdof", &AlignableNdof, "AlignableNdof/i");
36 }
37 
38 // ----------------------------------------------------------------------------
39 
41 {
42  tree->SetBranchAddress("Id", &Id);
43  tree->SetBranchAddress("ObjId", &ObjId);
44 
45  tree->SetBranchAddress("Nhit", &Nhit);
46  tree->SetBranchAddress("Nparj", &Nparj);
47  tree->SetBranchAddress("Jtvj", &Jtvj);
48  tree->SetBranchAddress("Npare", &Npare);
49  tree->SetBranchAddress("Jtve", &Jtve);
50  tree->SetBranchAddress("AlignableChi2", &AlignableChi2);
51  tree->SetBranchAddress("AlignableNdof", &AlignableNdof);
52 
53 }
54 
55 // ----------------------------------------------------------------------------
56 // find tree entry based on detID and typeID
57 
58 int HIPUserVariablesIORoot::findEntry(unsigned int detId,int comp)
59 {
60  if (newopen) { // we're here for the first time
61  edm::LogInfo("Alignment") <<"[HIPUserVariablesIORoot::findEntry] fill map ...";
62  treemap.erase(treemap.begin(),treemap.end());
63  for (int ev = 0;ev<tree->GetEntries();ev++) {
64  tree->GetEntry(ev);
65  treemap[std::make_pair(Id,ObjId)]=ev;
66  }
67  newopen=false;
68  }
69 
70  // now we have filled the map
71  treemaptype::iterator imap = treemap.find(std::make_pair(detId,comp));
72  int result=-1;
73  if (imap != treemap.end()) result=(*imap).second;
74  return result;
75 
76 
77  //double noAliPar = tree->GetEntries();
78  //for (int ev = 0;ev<noAliPar;ev++) {
79  // tree->GetEntry(ev);
80  // if(Id==detId&&comp==ObjId) return (ev);
81  //}
82  //return(-1);
83 }
84 
85 // ----------------------------------------------------------------------------
86 
88 {
90 
91  if ((ap->userVariables())==0) {
92  edm::LogError("Alignment") <<"UserVariables not found!";
93  return -1;
94  }
95 
96  HIPUserVariables* uvar =
97  dynamic_cast<HIPUserVariables*>(ap->userVariables());
98 
99  AlgebraicSymMatrix jtvj = uvar->jtvj;
100  AlgebraicVector jtve = uvar->jtve;
101  int nhit=uvar->nhit;
102  int np=jtve.num_row();
103 
104  Nhit=nhit;
105  Npare=np;
106  Nparj=np*(np+1)/2;
107  int count=0;
108  for(int row=0;row<np;row++){
109  Jtve[row]=jtve[row];
110  for(int col=0;col<np;col++){
111  if(row-1<col){Jtvj[count]=jtvj[row][col];count++;}
112  }
113  }
114  Id = ali->id();
115  ObjId = ali->alignableObjectId();
116 
117  //Chi^2 of alignable
118  AlignableChi2= uvar->alichi2 ;
119  AlignableNdof= uvar->alindof ;
120 
121  tree->Fill();
122  return 0;
123 }
124 
125 // ----------------------------------------------------------------------------
126 
128  int& ierr)
129 {
130  ierr=0;
131  HIPUserVariables* uvar;
132 
133  int entry = findEntry(ali->id(), ali->alignableObjectId());
134  if(entry!=-1) {
135  tree->GetEntry(entry);
136 
137  int np=Npare;
138  AlgebraicVector jtve(np,0);
139  AlgebraicSymMatrix jtvj(np,0);
140  int count=0;
141  for(int row=0;row<np;row++) {
142  jtve[row]=Jtve[row];
143  for(int col=0; col < np;col++) {
144  if(row-1<col) {jtvj[row][col]=Jtvj[count];count++;}
145  }
146  }
147 
148  uvar = new HIPUserVariables(np);
149  uvar->jtvj=jtvj;
150  uvar->jtve=jtve;
151  uvar->nhit=Nhit;
152 
153  //Chi2n
154  uvar->alichi2=AlignableChi2;
155  uvar->alindof=AlignableNdof;
156 
157  return uvar;
158  }
159 
160  // ierr=-1;
161  return 0 ;
162 }
163 
164 //-----------------------------------------------------------------------------
165 
166 void
168  const char* filename, int iter, bool validCheck, int& ierr)
169 {
170  ierr=0;
171  int iret;
172  iret = open(filename,iter,true);
173  if (iret!=0) { ierr=-1; return;}
174  iret = write(alivec,validCheck);
175  if (iret!=0) { ierr=-2; return;}
176  iret = close();
177  if (iret!=0) { ierr=-3; return;}
178 }
179 
180 //-----------------------------------------------------------------------------
181 
182 std::vector<AlignmentUserVariables*>
184  const char* filename, int iter, int& ierr)
185 {
186  std::vector<AlignmentUserVariables*> result;
187  ierr=0;
188  int iret;
189  iret = open(filename,iter,false);
190  if (iret!=0) { ierr=-1; return result;}
191  result = read(alivec,iret);
192  if (iret!=0) { ierr=-2; return result;}
193  iret = close();
194  if (iret!=0) { ierr=-3; return result;}
195 
196  return result;
197 }
int writeOne(Alignable *ali)
align::ID id() const
Return the ID of Alignable, i.e. DetId of &#39;first&#39; component GeomDet(Unit).
Definition: Alignable.h:180
int findEntry(unsigned int detId, int comp)
int open(const char *filename, int iteration, bool writemode)
std::vector< Alignable * > Alignables
AlignmentUserVariables * readOne(Alignable *ali, int &ierr)
void setBranchAddresses(void)
set root branches
void createBranches(void)
create root branches
AlgebraicVector jtve
AlignmentParameters * alignmentParameters() const
Get the AlignmentParameters.
Definition: Alignable.h:57
void writeHIPUserVariables(const Alignables &alivec, const char *filename, int iter, bool validCheck, int &ierr)
AlgebraicSymMatrix jtvj
AlignmentUserVariables * userVariables(void) const
Get pointer to user variables.
std::vector< AlignmentUserVariables * > readHIPUserVariables(const Alignables &alivec, const char *filename, int iter, int &ierr)
std::pair< std::string, MonitorElement * > entry
Definition: ME_MAP.h:8
double Jtvj[nparmax *(nparmax+1)/2]
tuple result
Definition: query.py:137
virtual StructureType alignableObjectId() const =0
Return the alignable type identifier.
(Abstract) Base class for alignment algorithm user variables
CLHEP::HepVector AlgebraicVector
int write(const align::Alignables &alivec, bool validCheck)
std::vector< AlignmentUserVariables * > read(const align::Alignables &alivec, int &ierr)
tuple filename
Definition: lut2db_cfg.py:20
CLHEP::HepSymMatrix AlgebraicSymMatrix