CMS 3D CMS Logo

MillePedeVariablesIORoot.cc
Go to the documentation of this file.
1 
8 // this class's header
10 
12 
16 
17 #include "TTree.h"
18 
19 // -------------------------------------------------------------------------------------------------
21  myId(0), myObjId(0), myNumPar(0),
22  myHitsX(0), myHitsY(0), myLabel(0), myName(""), myNamePtr(&myName)
23 {
24  treename = "MillePedeUser";
25  treetxt = "MillePede User Variables";
26  for (unsigned int i=0;i<kMaxNumPar;i++) {
27  myIsValid[i] = 0;
28  myDiffBefore[i] = 0.;
29  myGlobalCor[i] = 0.;
30  myPreSigma[i] = 0.;
31  myParameter[i] = 0.;
32  mySigma[i] = 0.;
33  }
34 }
35 
36 // -------------------------------------------------------------------------------------------------
38 (const align::Alignables& alivec, const char *filename, int iter, bool validCheck, int &ierr)
39 {
40  ierr = 0;
41 
42  int iret = this->open(filename, iter, true);
43  if (iret != 0) {
44  ierr = -1;
45  } else {
46  iret = this->write(alivec, validCheck);
47  tree->BuildIndex("Id", "ObjId");
48  if (iret != 0) {
49  ierr = -2;
50  } else {
51  iret = this->close();
52  if (iret != 0) {
53  ierr = -3;
54  }
55  }
56  }
57 
58  return;
59 }
60 
61 // -------------------------------------------------------------------------------------------------
62 std::vector<AlignmentUserVariables*> MillePedeVariablesIORoot::readMillePedeVariables
63 (const align::Alignables& alivec, const char *filename, int iter, int &ierr)
64 {
65  std::vector<AlignmentUserVariables*> result;
66  ierr = 0;
67  int iret = this->open(filename, iter, false);
68  if (iret != 0) {
69  ierr = -1;
70  } else {
71  result = this->read(alivec, iret);
72  if (iret != 0) {
73  ierr = -2;
74  } else {
75  iret = this->close();
76  if (iret != 0) {
77  ierr = -3;
78  }
79  }
80  }
81 
82  return result;
83 }
84 
85 // -------------------------------------------------------------------------------------------------
87 {
88  if (!ali || !ali->alignmentParameters()
89  || !dynamic_cast<MillePedeVariables*>(ali->alignmentParameters()->userVariables())) {
90  edm::LogError("Alignment") << "@SUB=MillePedeVariablesIORoot::writeOne"
91  << "No MillePedeVariables found!";
92  return -1;
93  }
94 
95  const MillePedeVariables *mpVar =
96  static_cast<MillePedeVariables*>(ali->alignmentParameters()->userVariables());
97  myNumPar = mpVar->size();
98  if (myNumPar >= kMaxNumPar) {
99  edm::LogError("Alignment") << "@SUB=MillePedeVariablesIORoot::writeOne"
100  << "Ignoring parameters " << static_cast<int>(kMaxNumPar) << " to " << myNumPar-1;
102  }
103 
104  for (unsigned int iPar = 0; iPar < myNumPar; ++iPar) {
105  myIsValid[iPar] = mpVar->isValid()[iPar];
106  myDiffBefore[iPar] = mpVar->diffBefore()[iPar];
107  myGlobalCor[iPar] = mpVar->globalCor()[iPar];
108  myPreSigma[iPar] = mpVar->preSigma()[iPar];
109  myParameter[iPar] = mpVar->parameter()[iPar];
110  mySigma[iPar] = mpVar->sigma()[iPar];
111  }
112  myHitsX = mpVar->hitsX();
113  myHitsY = mpVar->hitsY();
114  myLabel = mpVar->label();
115  myName = mpVar->name();
116 
117  myId = ali->id();
118  myObjId = ali->alignableObjectId();
119 
120  tree->Fill();
121 
122  return 0;
123 }
124 
125 // -------------------------------------------------------------------------------------------------
127 {
128  ierr = 0;
129 
130  if (tree->GetEntryWithIndex(ali->id(), ali->alignableObjectId()) < 0) {
131  edm::LogError("Alignment") << "@SUB=MillePedeVariablesIORoot::readOne"
132  << "No index for id/type = (" << ali->id() << "/"
133  << ali->alignableObjectId() << ") found!";
134  ierr = 1;
135  return nullptr;
136  }
137 
139  for (unsigned int iPar = 0; iPar < myNumPar; ++iPar) {
140  mpVar->isValid()[iPar] = myIsValid[iPar];
141  mpVar->diffBefore()[iPar] = myDiffBefore[iPar];
142  mpVar->globalCor()[iPar] = myGlobalCor[iPar];
143  mpVar->preSigma()[iPar] = myPreSigma[iPar];
144  mpVar->parameter()[iPar] = myParameter[iPar];
145  mpVar->sigma()[iPar] = mySigma[iPar];
146  }
147  mpVar->setHitsX(myHitsX);
148  mpVar->setHitsY(myHitsY);
149 
150  return mpVar;
151 }
152 
153 // -------------------------------------------------------------------------------------------------
155 {
156  tree->Branch("Id", &myId, "Id/i");
157  tree->Branch("ObjId", &myObjId, "ObjId/I");
158  tree->Branch("NumPar", &myNumPar, "NumPar/i");
159  tree->Branch("IsValid", myIsValid, "IsValid[NumPar]/b");
160  tree->Branch("DiffBefore", myDiffBefore,"DiffBefore[NumPar]/F");
161  tree->Branch("GlobalCor", myGlobalCor, "GlobalCor[NumPar]/F");
162  tree->Branch("PreSigma", myPreSigma, "PreSigma[NumPar]/F");
163  tree->Branch("Par", myParameter, "Par[NumPar]/F"); // name as in AlignmentParametersIORoot
164  tree->Branch("Sigma", mySigma, "Sigma[NumPar]/F");
165  tree->Branch("HitsX", &myHitsX, "HitsX/i");
166  tree->Branch("HitsY", &myHitsY, "HitsY/i");
167  tree->Branch("Label", &myLabel, "Label/i");
168  tree->Branch("Name", &myNamePtr);
169 }
170 
171 // -------------------------------------------------------------------------------------------------
173 {
174  tree->SetBranchAddress("Id", &myId);
175  tree->SetBranchAddress("ObjId", &myObjId);
176  tree->SetBranchAddress("NumPar", &myNumPar);
177  tree->SetBranchAddress("IsValid", myIsValid);
178  tree->SetBranchAddress("DiffBefore", myDiffBefore);
179  tree->SetBranchAddress("GlobalCor", myGlobalCor);
180  tree->SetBranchAddress("PreSigma", myPreSigma);
181  tree->SetBranchAddress("Par", myParameter);
182  tree->SetBranchAddress("Sigma", mySigma);
183  tree->SetBranchAddress("HitsX", &myHitsX);
184  tree->SetBranchAddress("HitsY", &myHitsY);
185  tree->SetBranchAddress("Label", &myLabel);
186  tree->SetBranchAddress("Name", &myNamePtr);
187 }
unsigned int hitsX() const
get number of hits for x-measurement
unsigned int label() const
get alignable label as used by pede
const std::vector< float > & globalCor() const
get global correlation array
align::ID id() const
Return the ID of Alignable, i.e. DetId of &#39;first&#39; component GeomDet(Unit).
Definition: Alignable.h:189
unsigned int size() const
number of parameters
void setBranchAddresses() override
set root branche addresses
const std::vector< float > & parameter() const
get array of parameters
void writeMillePedeVariables(const align::Alignables &alivec, const char *filename, int iter, bool validCheck, int &ierr)
AlignmentParameters * alignmentParameters() const
Get the AlignmentParameters.
Definition: Alignable.h:61
const std::string & name() const
get alignable name
AlignmentUserVariables * userVariables(void) const
Get pointer to user variables.
AlignmentUserVariables * readOne(Alignable *ali, int &ierr) override
virtual StructureType alignableObjectId() const =0
Return the alignable type identifier.
void setHitsY(unsigned int hitsY)
const std::vector< float > & preSigma() const
get array of presigmas (<= 0: means fixed)
unsigned int hitsY() const
get number of hits for y-measurement
const std::vector< float > & sigma() const
get array of sigmas
std::vector< AlignmentUserVariables * > readMillePedeVariables(const align::Alignables &alivec, const char *filename, int iter, int &ierr)
(Abstract) Base class for alignment algorithm user variables
const std::vector< float > & diffBefore() const
get array of differences to start value
const std::vector< bool > & isValid() const
get valid flag array
std::vector< Alignable * > Alignables
Definition: Utilities.h:32
int write(const align::Alignables &alivec, bool validCheck)
std::vector< AlignmentUserVariables * > read(const align::Alignables &alivec, int &ierr)
Definition: tree.py:1
int open(const char *filename, int iteration, bool writemode) override
void createBranches() override
create root branches
int writeOne(Alignable *ali) override
void setHitsX(unsigned int hitsX)