CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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), myHitsX(0), myHitsY(0), myLabel(0), myName(""), myNamePtr(&myName) {
22  treename = "MillePedeUser";
23  treetxt = "MillePede User Variables";
24  for (unsigned int i = 0; i < kMaxNumPar; i++) {
25  myIsValid[i] = 0;
26  myDiffBefore[i] = 0.;
27  myGlobalCor[i] = 0.;
28  myPreSigma[i] = 0.;
29  myParameter[i] = 0.;
30  mySigma[i] = 0.;
31  }
32 }
33 
34 // -------------------------------------------------------------------------------------------------
36  const align::Alignables &alivec, const char *filename, int iter, bool validCheck, int &ierr) {
37  ierr = 0;
38 
39  int iret = this->open(filename, iter, true);
40  if (iret != 0) {
41  ierr = -1;
42  } else {
43  iret = this->write(alivec, validCheck);
44  tree->BuildIndex("Id", "ObjId");
45  if (iret != 0) {
46  ierr = -2;
47  } else {
48  iret = this->close();
49  if (iret != 0) {
50  ierr = -3;
51  }
52  }
53  }
54 
55  return;
56 }
57 
58 // -------------------------------------------------------------------------------------------------
59 std::vector<AlignmentUserVariables *> MillePedeVariablesIORoot::readMillePedeVariables(const align::Alignables &alivec,
60  const char *filename,
61  int iter,
62  int &ierr) {
63  std::vector<AlignmentUserVariables *> result;
64  ierr = 0;
65  int iret = this->open(filename, iter, false);
66  if (iret != 0) {
67  ierr = -1;
68  } else {
69  result = this->read(alivec, iret);
70  if (iret != 0) {
71  ierr = -2;
72  } else {
73  iret = this->close();
74  if (iret != 0) {
75  ierr = -3;
76  }
77  }
78  }
79 
80  return result;
81 }
82 
83 // -------------------------------------------------------------------------------------------------
85  if (!ali || !ali->alignmentParameters() ||
86  !dynamic_cast<MillePedeVariables *>(ali->alignmentParameters()->userVariables())) {
87  edm::LogError("Alignment") << "@SUB=MillePedeVariablesIORoot::writeOne"
88  << "No MillePedeVariables found!";
89  return -1;
90  }
91 
92  const MillePedeVariables *mpVar = static_cast<MillePedeVariables *>(ali->alignmentParameters()->userVariables());
93  myNumPar = mpVar->size();
94  if (myNumPar >= kMaxNumPar) {
95  edm::LogError("Alignment") << "@SUB=MillePedeVariablesIORoot::writeOne"
96  << "Ignoring parameters " << static_cast<int>(kMaxNumPar) << " to " << myNumPar - 1;
98  }
99 
100  for (unsigned int iPar = 0; iPar < myNumPar; ++iPar) {
101  myIsValid[iPar] = mpVar->isValid()[iPar];
102  myDiffBefore[iPar] = mpVar->diffBefore()[iPar];
103  myGlobalCor[iPar] = mpVar->globalCor()[iPar];
104  myPreSigma[iPar] = mpVar->preSigma()[iPar];
105  myParameter[iPar] = mpVar->parameter()[iPar];
106  mySigma[iPar] = mpVar->sigma()[iPar];
107  }
108  myHitsX = mpVar->hitsX();
109  myHitsY = mpVar->hitsY();
110  myLabel = mpVar->label();
111  myName = mpVar->name();
112 
113  myId = ali->id();
114  myObjId = ali->alignableObjectId();
115 
116  tree->Fill();
117 
118  return 0;
119 }
120 
121 // -------------------------------------------------------------------------------------------------
123  ierr = 0;
124 
125  if (tree->GetEntryWithIndex(ali->id(), ali->alignableObjectId()) < 0) {
126  edm::LogError("Alignment") << "@SUB=MillePedeVariablesIORoot::readOne"
127  << "No index for id/type = (" << ali->id() << "/" << ali->alignableObjectId()
128  << ") found!";
129  ierr = 1;
130  return nullptr;
131  }
132 
134  for (unsigned int iPar = 0; iPar < myNumPar; ++iPar) {
135  mpVar->isValid()[iPar] = myIsValid[iPar];
136  mpVar->diffBefore()[iPar] = myDiffBefore[iPar];
137  mpVar->globalCor()[iPar] = myGlobalCor[iPar];
138  mpVar->preSigma()[iPar] = myPreSigma[iPar];
139  mpVar->parameter()[iPar] = myParameter[iPar];
140  mpVar->sigma()[iPar] = mySigma[iPar];
141  }
142  mpVar->setHitsX(myHitsX);
143  mpVar->setHitsY(myHitsY);
144 
145  return mpVar;
146 }
147 
148 // -------------------------------------------------------------------------------------------------
150  tree->Branch("Id", &myId, "Id/i");
151  tree->Branch("ObjId", &myObjId, "ObjId/I");
152  tree->Branch("NumPar", &myNumPar, "NumPar/i");
153  tree->Branch("IsValid", myIsValid, "IsValid[NumPar]/b");
154  tree->Branch("DiffBefore", myDiffBefore, "DiffBefore[NumPar]/F");
155  tree->Branch("GlobalCor", myGlobalCor, "GlobalCor[NumPar]/F");
156  tree->Branch("PreSigma", myPreSigma, "PreSigma[NumPar]/F");
157  tree->Branch("Par", myParameter, "Par[NumPar]/F"); // name as in AlignmentParametersIORoot
158  tree->Branch("Sigma", mySigma, "Sigma[NumPar]/F");
159  tree->Branch("HitsX", &myHitsX, "HitsX/i");
160  tree->Branch("HitsY", &myHitsY, "HitsY/i");
161  tree->Branch("Label", &myLabel, "Label/i");
162  tree->Branch("Name", &myNamePtr);
163 }
164 
165 // -------------------------------------------------------------------------------------------------
167  tree->SetBranchAddress("Id", &myId);
168  tree->SetBranchAddress("ObjId", &myObjId);
169  tree->SetBranchAddress("NumPar", &myNumPar);
170  tree->SetBranchAddress("IsValid", myIsValid);
171  tree->SetBranchAddress("DiffBefore", myDiffBefore);
172  tree->SetBranchAddress("GlobalCor", myGlobalCor);
173  tree->SetBranchAddress("PreSigma", myPreSigma);
174  tree->SetBranchAddress("Par", myParameter);
175  tree->SetBranchAddress("Sigma", mySigma);
176  tree->SetBranchAddress("HitsX", &myHitsX);
177  tree->SetBranchAddress("HitsY", &myHitsY);
178  tree->SetBranchAddress("Label", &myLabel);
179  tree->SetBranchAddress("Name", &myNamePtr);
180 }
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:180
unsigned int size() const
number of parameters
void setBranchAddresses() override
set root branche addresses
const std::vector< float > & parameter() const
get array of parameters
Log< level::Error, false > LogError
void writeMillePedeVariables(const align::Alignables &alivec, const char *filename, int iter, bool validCheck, int &ierr)
AlignmentParameters * alignmentParameters() const
Get the AlignmentParameters.
Definition: Alignable.h:58
const std::string & name() const
get alignable name
tuple result
Definition: mps_fire.py:311
AlignmentUserVariables * userVariables(void) const
Get pointer to user variables.
AlignmentUserVariables * readOne(Alignable *ali, int &ierr) override
void setHitsY(unsigned int hitsY)
const std::vector< float > & preSigma() const
get array of presigmas (&lt;= 0: means fixed)
unsigned int hitsY() const
get number of hits for y-measurement
const std::vector< float > & sigma() const
get array of sigmas
virtual StructureType alignableObjectId() const =0
Return the alignable type identifier.
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:31
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
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)