CMS 3D CMS Logo

AlignmentParametersIORoot.cc
Go to the documentation of this file.
1 // this class's header
3 
7 
9 
10 #include "TTree.h"
11 
12 // ----------------------------------------------------------------------------
13 // constructor
15  treename = "AlignmentParameters";
16  treetxt = "Alignment Parameters";
17 }
18 
19 // ----------------------------------------------------------------------------
21  tree->Branch("parSize", &theCovRang, "CovRang/I");
22  tree->Branch("Id", &theId, "Id/i");
23  tree->Branch("paramType", &theParamType, "paramType/I");
24  tree->Branch("Par", &thePar, "Par[CovRang]/D");
25  tree->Branch("covarSize", &theCovarRang, "CovarRang/I");
26  tree->Branch("Cov", &theCov, "Cov[CovarRang]/D");
27  tree->Branch("ObjId", &theObjId, "ObjId/I");
28  tree->Branch("HieraLevel", &theHieraLevel, "HieraLevel/I");
29 }
30 
31 // ----------------------------------------------------------------------------
33  tree->SetBranchAddress("parSize", &theCovRang);
34  tree->SetBranchAddress("covarSize", &theCovarRang);
35  tree->SetBranchAddress("Id", &theId);
36  tree->SetBranchAddress("Par", &thePar);
37  tree->SetBranchAddress("paramType", &theParamType);
38  tree->SetBranchAddress("Cov", &theCov);
39  tree->SetBranchAddress("ObjId", &theObjId);
40  tree->SetBranchAddress("HieraLevel", &theHieraLevel);
41 }
42 
43 // ----------------------------------------------------------------------------
45  const AlignmentParameters* ap = ali->alignmentParameters();
46  const AlgebraicVector& params = ap->parameters();
47  const AlgebraicSymMatrix& cov = ap->covariance();
48 
49  theCovRang = params.num_row();
50  theCovarRang = theCovRang * (theCovRang + 1) / 2;
51  int count = 0;
52  for (int row = 0; row < theCovRang; row++) {
53  thePar[row] = params[row];
54  for (int col = 0; col < theCovRang; col++) {
55  if (row - 1 < col) {
56  theCov[count] = cov[row][col];
57  count++;
58  }
59  }
60  }
61 
62  theId = ali->id();
63  theParamType = ap->type();
64  theObjId = ali->alignableObjectId();
66 
67  tree->Fill();
68  return 0;
69 }
70 
71 // ----------------------------------------------------------------------------
73  if (tree->GetEntryWithIndex(ali->id(), ali->alignableObjectId()) > 0) {
74  int covsize = theCovRang;
75  int count = 0;
76  AlgebraicVector par(covsize, 0);
77  AlgebraicSymMatrix cov(covsize, 0);
78  for (int row = 0; row < covsize; row++) {
79  par[row] = thePar[row];
80  for (int col = 0; col < covsize; col++) {
81  if (row - 1 < col) {
82  cov[row][col] = theCov[count];
83  count++;
84  }
85  }
86  }
87 
88  using namespace AlignmentParametersFactory;
90  AlignmentParameters* alipar1;
91  if (ali->alignmentParameters()) {
92  const std::vector<bool>& sel = ali->alignmentParameters()->selector();
93  alipar1 = createParameters(ali, parType, sel);
94  } else {
95  const std::vector<bool> sel(theCovRang, true);
96  alipar1 = createParameters(ali, parType, sel);
97  }
98  AlignmentParameters* alipar = alipar1->clone(par, cov);
99  alipar->setValid(true);
100  ierr = 0;
101  delete alipar1;
102  return alipar;
103  }
104 
105  ierr = -1;
106  return (nullptr);
107 }
108 
110  if (bWrite) {
111  int nIndices = tree->BuildIndex("Id", "ObjId");
112  edm::LogInfo("Alignment") << "@SUB=AlignmentParametersIORoot::setBranchAddresses"
113  << "number of indexed entries: " << nIndices;
114  }
115 
116  return closeRoot();
117 }
AlignmentParameters * alignmentParameters() const
Get the AlignmentParameters.
Definition: Alignable.h:58
const AlgebraicSymMatrix & covariance(void) const
Get parameter covariance matrix.
const AlgebraicVector & parameters(void) const
Get alignment parameters.
int writeOne(Alignable *ali) override
Write AlignmentParameters of one Alignable.
ParametersType
enums for all available AlignmentParameters
double theCov[nParMax *(nParMax+1)/2]
void setValid(bool v)
Set validity flag.
void createBranches(void) override
Create all branches and give names.
virtual StructureType alignableObjectId() const =0
Return the alignable type identifier.
AlignmentParameters * createParameters(Alignable *ali, ParametersType parType, const std::vector< bool > &sel)
Log< level::Info, false > LogInfo
CLHEP::HepVector AlgebraicVector
align::ID id() const
Return the ID of Alignable, i.e. DetId of &#39;first&#39; component GeomDet(Unit).
Definition: Alignable.h:180
virtual unsigned int hierarchyLevel() const
virtual int type() const =0
tell type (AlignmentParametersFactory::ParametersType - but no circular dependency) ...
int closeRoot(void)
close IO
virtual AlignmentParameters * clone(const AlgebraicVector &par, const AlgebraicSymMatrix &cov) const =0
Enforce clone methods in derived classes.
const std::vector< bool > & selector(void) const
Get alignment parameter selector vector.
void setBranchAddresses(void) override
Set branch adresses.
col
Definition: cuy.py:1009
CLHEP::HepSymMatrix AlgebraicSymMatrix
ParametersType parametersType(const std::string &typeString)
convert string to ParametersType - exception if not known
Definition: tree.py:1
int close(void) override
Close IO.
AlignmentParameters * readOne(Alignable *ali, int &ierr) override
Read AlignmentParameters of one Alignable.