CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
AlignmentCorrelationsIORoot.cc
Go to the documentation of this file.
1 #include "TTree.h"
2 
6 
7 // this class's header
9 
10 // ----------------------------------------------------------------------------
11 // constructor
12 
14 {
15  treename = "AlignmentCorrelations";
16  treetxt = "Correlations";
17 }
18 
19 // ----------------------------------------------------------------------------
20 
22 {
23  tree->Branch("Ali1Id", &Ali1Id, "Ali1Id/i");
24  tree->Branch("Ali2Id", &Ali2Id, "Ali2Id/i");
25  tree->Branch("Ali1ObjId", &Ali1ObjId, "Ali1ObjId/I");
26  tree->Branch("Ali2ObjId", &Ali2ObjId, "Ali2ObjId/I");
27  tree->Branch("corSize", &corSize, "corSize/I");
28  tree->Branch("CorMatrix", &CorMatrix, "CorMatrix[corSize]/D");
29 }
30 
31 // ----------------------------------------------------------------------------
32 
34 {
35  tree->SetBranchAddress("corSize", &corSize);
36  tree->SetBranchAddress("Ali1Id", &Ali1Id);
37  tree->SetBranchAddress("Ali2Id", &Ali2Id);
38  tree->SetBranchAddress("Ali1ObjId", &Ali1ObjId);
39  tree->SetBranchAddress("Ali2ObjId", &Ali2ObjId);
40  tree->SetBranchAddress("CorMatrix", &CorMatrix);
41 }
42 
43 // ----------------------------------------------------------------------------
44 
46 {
47  int icount=0;
48 
49  for(align::Correlations::const_iterator it=cor.begin();
50  it!=cor.end();++it) {
51  AlgebraicMatrix mat=(*it).second;
52  std::pair<Alignable*,Alignable*> Pair = (*it).first;
53  Alignable* ali1 = Pair.first;
54  Alignable* ali2 = Pair.second;
55  if( (ali1->alignmentParameters()->isValid()
56  && ali2->alignmentParameters()->isValid()) || !(validCheck)) {
57  Ali1ObjId = ali1->alignableObjectId();
58  Ali2ObjId = ali2->alignableObjectId();
59  Ali1Id = ali1->id();
60  Ali2Id = ali2->id();
61  int maxColumn = mat.num_row();
62  corSize = maxColumn*maxColumn;
63  for(int row = 0;row<maxColumn;row++)
64  for(int col = 0;col<maxColumn;col++)
65  CorMatrix[row+col*maxColumn] =mat[row][col];
66  tree->Fill();
67  icount++;
68  }
69  }
70  edm::LogInfo("AlignmentCorrelationsIORoot") << "Writing correlations: all,written: "
71  << cor.size() << "," << icount;
72  return 0;
73 }
74 
75 // ----------------------------------------------------------------------------
76 // read correlations for those alignables in vector given as argument
77 
80 {
81  align::Correlations theMap;
82 
83  // create ID map for all Alignables in alivec
84  align::Alignables::const_iterator it1;
85  std::map< std::pair<unsigned int,int>, Alignable* > idAlis;
86  for( it1=alivec.begin();it1!=alivec.end();++it1 )
87  idAlis[std::make_pair((*it1)->id(),(*it1)->alignableObjectId())] = (*it1);
88 
89  std::map<std::pair<unsigned int,int>,Alignable*>::const_iterator aliSearch1;
90  std::map<std::pair<unsigned int,int>,Alignable*>::const_iterator aliSearch2;
91  int nfound=0;
92  double maxEntry = tree->GetEntries();
93  for( int entry = 0;entry<maxEntry;entry++ )
94  {
95  tree->GetEntry(entry);
96  aliSearch1 = idAlis.find(std::make_pair(Ali1Id,Ali1ObjId));
97  aliSearch2 = idAlis.find(std::make_pair(Ali2Id,Ali2ObjId));
98  if (aliSearch1!=idAlis.end()
99  && aliSearch2!=idAlis.end())
100  {
101  // Alignables for this pair found
102  nfound++;
103  Alignable* myAli1 = (*aliSearch1).second;
104  Alignable* myAli2 = (*aliSearch2).second;
105  // FIXME: instead of nParMax in the next few lines one should probably
106  // use something like sqrt(corSize) - but take care of rounding!
107  // I have no time to test... :-( GF
109  for(int row = 0;row<nParMax;row++)
110  for(int col = 0;col<nParMax;col++)
111  mat[row][col] = CorMatrix[row+col*nParMax];
112  theMap[ std::make_pair(myAli1,myAli2) ] = mat;
113  }
114  }
115 
116  edm::LogInfo("AlignmentCorrelationsIORoot") << "Read correlations: all,read: "
117  << alivec.size() << "," << nfound;
118 
119  ierr=0;
120  return theMap;
121 }
align::ID Ali1Id
correlation tree
align::ID id() const
Return the ID of Alignable, i.e. DetId of &#39;first&#39; component GeomDet(Unit).
Definition: Alignable.h:185
static const int nParMax
align::Correlations read(const align::Alignables &alivec, int &ierr)
read correlations
std::map< std::pair< Alignable *, Alignable * >, AlgebraicMatrix > Correlations
Definition: Utilities.h:31
void createBranches(void)
create root branches
AlignmentParameters * alignmentParameters() const
Get the AlignmentParameters.
Definition: Alignable.h:57
CLHEP::HepMatrix AlgebraicMatrix
virtual StructureType alignableObjectId() const =0
Return the alignable type identifier.
std::vector< Alignable * > Alignables
Definition: Utilities.h:28
bool isValid(void) const
Get validity flag.
int write(const align::Correlations &cor, bool validCheck)
write correlations
void setBranchAddresses(void)
set root branches
int col
Definition: cuy.py:1008