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  ObjId(0), Id(0),Nhit(0), Nparj(0), Npare(0),
18  AlignableChi2(0.), AlignableNdof(0)
19 {
20  treename = "T9";
21  treetxt = "HIP User Variables";
22 
23  for (int i=0;i<nparmax*(nparmax+1)/2;++i)
24  Jtvj[i] = 0.;
25  for (int i=0;i<nparmax;++i)
26  Jtve[i] = 0.;
27 }
28 
29 // ----------------------------------------------------------------------------
30 
32 {
33  tree->Branch("Id", &Id, "Id/i");
34  tree->Branch("ObjId", &ObjId, "ObjId/I");
35 
36  tree->Branch("Nhit", &Nhit, "Nhit/I");
37  tree->Branch("Nparj", &Nparj, "Nparj/I");
38  tree->Branch("Jtvj", &Jtvj, "Jtvj[Nparj]/D");
39  tree->Branch("Npare", &Npare, "Npare/I");
40  tree->Branch("Jtve", &Jtve, "Jtve[Npare]/D");
41  tree->Branch("AlignableChi2", &AlignableChi2, "AlignableChi2/D");
42  tree->Branch("AlignableNdof", &AlignableNdof, "AlignableNdof/i");
43 }
44 
45 // ----------------------------------------------------------------------------
46 
48 {
49  tree->SetBranchAddress("Id", &Id);
50  tree->SetBranchAddress("ObjId", &ObjId);
51 
52  tree->SetBranchAddress("Nhit", &Nhit);
53  tree->SetBranchAddress("Nparj", &Nparj);
54  tree->SetBranchAddress("Jtvj", &Jtvj);
55  tree->SetBranchAddress("Npare", &Npare);
56  tree->SetBranchAddress("Jtve", &Jtve);
57  tree->SetBranchAddress("AlignableChi2", &AlignableChi2);
58  tree->SetBranchAddress("AlignableNdof", &AlignableNdof);
59 
60 }
61 
62 // ----------------------------------------------------------------------------
63 // find tree entry based on detID and typeID
64 
65 int HIPUserVariablesIORoot::findEntry(unsigned int detId,int comp)
66 {
67  if (newopen) { // we're here for the first time
68  edm::LogInfo("Alignment") <<"[HIPUserVariablesIORoot::findEntry] fill map ...";
69  treemap.erase(treemap.begin(),treemap.end());
70  for (int ev = 0;ev<tree->GetEntries();ev++) {
71  tree->GetEntry(ev);
72  treemap[std::make_pair(Id,ObjId)]=ev;
73  }
74  newopen=false;
75  }
76 
77  // now we have filled the map
78  treemaptype::iterator imap = treemap.find(std::make_pair(detId,comp));
79  int result=-1;
80  if (imap != treemap.end()) result=(*imap).second;
81  return result;
82 
83 
84  //double noAliPar = tree->GetEntries();
85  //for (int ev = 0;ev<noAliPar;ev++) {
86  // tree->GetEntry(ev);
87  // if(Id==detId&&comp==ObjId) return (ev);
88  //}
89  //return(-1);
90 }
91 
92 // ----------------------------------------------------------------------------
93 
95 {
97 
98  if ((ap->userVariables())==0) {
99  edm::LogError("Alignment") <<"UserVariables not found!";
100  return -1;
101  }
102 
103  HIPUserVariables* uvar =
104  dynamic_cast<HIPUserVariables*>(ap->userVariables());
105 
106  AlgebraicSymMatrix jtvj = uvar->jtvj;
107  AlgebraicVector jtve = uvar->jtve;
108  int nhit=uvar->nhit;
109  int np=jtve.num_row();
110 
111  Nhit=nhit;
112  Npare=np;
113  Nparj=np*(np+1)/2;
114  int count=0;
115  for(int row=0;row<np;row++){
116  Jtve[row]=jtve[row];
117  for(int col=0;col<np;col++){
118  if(row-1<col){Jtvj[count]=jtvj[row][col];count++;}
119  }
120  }
121  Id = ali->id();
122  ObjId = ali->alignableObjectId();
123 
124  //Chi^2 of alignable
125  AlignableChi2= uvar->alichi2 ;
126  AlignableNdof= uvar->alindof ;
127 
128  tree->Fill();
129  return 0;
130 }
131 
132 // ----------------------------------------------------------------------------
133 
135  int& ierr)
136 {
137  ierr=0;
138  HIPUserVariables* uvar;
139 
140  int entry = findEntry(ali->id(), ali->alignableObjectId());
141  if(entry!=-1) {
142  tree->GetEntry(entry);
143 
144  int np=Npare;
145  AlgebraicVector jtve(np,0);
146  AlgebraicSymMatrix jtvj(np,0);
147  int count=0;
148  for(int row=0;row<np;row++) {
149  jtve[row]=Jtve[row];
150  for(int col=0; col < np;col++) {
151  if(row-1<col) {jtvj[row][col]=Jtvj[count];count++;}
152  }
153  }
154 
155  uvar = new HIPUserVariables(np);
156  uvar->jtvj=jtvj;
157  uvar->jtve=jtve;
158  uvar->nhit=Nhit;
159 
160  //Chi2n
161  uvar->alichi2=AlignableChi2;
162  uvar->alindof=AlignableNdof;
163 
164  return uvar;
165  }
166 
167  // ierr=-1;
168  return 0 ;
169 }
170 
171 //-----------------------------------------------------------------------------
172 
173 void
175  const char* filename, int iter, bool validCheck, int& ierr)
176 {
177  ierr=0;
178  int iret;
179  iret = open(filename,iter,true);
180  if (iret!=0) { ierr=-1; return;}
181  iret = write(alivec,validCheck);
182  if (iret!=0) { ierr=-2; return;}
183  iret = close();
184  if (iret!=0) { ierr=-3; return;}
185 }
186 
187 //-----------------------------------------------------------------------------
188 
189 std::vector<AlignmentUserVariables*>
191  const char* filename, int iter, int& ierr)
192 {
193  std::vector<AlignmentUserVariables*> result;
194  ierr=0;
195  int iret;
196  iret = open(filename,iter,false);
197  if (iret!=0) { ierr=-1; return result;}
198  result = read(alivec,iret);
199  if (iret!=0) { ierr=-2; return result;}
200  iret = close();
201  if (iret!=0) { ierr=-3; return result;}
202 
203  return result;
204 }
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 i
Definition: DBlmapReader.cc:9
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.
int np
Definition: AMPTWrapper.h:33
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
int col
Definition: cuy.py:1008