CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CocoaDBMgr.cc
Go to the documentation of this file.
11 
21 
30 
34 
36 
38 
39 //----------------------------------------------------------------------
41 {
42  if(!instance) {
43  instance = new CocoaDBMgr;
44  }
45  return instance;
46 }
47 
48 //----------------------------------------------------------------------
50 {
51 }
52 
53 //-----------------------------------------------------------------------
55 {
57 
59  int nrcd;
60 
61  cond::Time_t appendTime = Fit::nEvent+1;
62  if(gomgr->GlobalOptions()["writeDBOptAlign"] > 0 ) {
63 
64  //----- Build OpticalAlignments
66 
67  //--- Dump OpticalAlignments
68  nrcd = optalign->opticalAlignments_.size();
69  if( !myDbService.isAvailable() ){
70  throw cms::Exception("CocoaDBMgr::DumpCocoaResults DB not available");
71  }
72  // try{
73  if ( myDbService->isNewTagRequest( "OpticalAlignmentsRcd" ) ) {
74  std::cout << " new OA to DB " << "begin " << myDbService->beginOfTime() << " current " << myDbService->currentTime() << " end " << myDbService->endOfTime() << std::endl;
75  myDbService->createNewIOV<OpticalAlignments>(optalign,
76  myDbService->beginOfTime(),
77  myDbService->endOfTime(),
78  // myDbService->endOfTime(),
79  "OpticalAlignmentsRcd");
80  } else {
81  std::cout << " old OA to DB " << " current " << myDbService->currentTime() << " end " << myDbService->endOfTime() << std::endl;
82  myDbService->appendSinceTime<OpticalAlignments>( optalign,
83  // myDbService->endOfTime(),
84  appendTime,
85  // myDbService->currentTime(),
86  "OpticalAlignmentsRcd");
87  }
88 
89 
90  /* }catch(const cond::Exception& er) {
91  std::cout<<er.what()<<std::endl;
92  }catch(const std::exception& er){
93  std::cout<<"caught std::exception "<<er.what()<<std::endl;
94  }catch(...){
95  std::cout<<"Funny error"<<std::endl;
96  } */
97 
98  if(ALIUtils::debug >= 2) std::cout << "OpticalAlignmentsRcd WRITTEN TO DB : "<< nrcd << std::endl;
99  }
100 
101  if( gomgr->GlobalOptions()["writeDBAlign"] > 0) {
102 
103  // Build DT alignments and errors
104  std::pair< Alignments*,AlignmentErrorsExtended*> dtali = BuildAlignments(1);
105  Alignments* dt_Alignments = dtali.first;
106  AlignmentErrorsExtended* dt_AlignmentErrors = dtali.second;
107 
108  // Dump DT alignments and errors
109  nrcd = dt_Alignments->m_align.size();
110  if ( myDbService->isNewTagRequest( "DTAlignmentRcd" ) ) {
111  myDbService->createNewIOV<Alignments>(&(*dt_Alignments),
112  myDbService->beginOfTime(),
113  myDbService->endOfTime(),
114  "DTAlignmentRcd");
115  } else {
116  myDbService->appendSinceTime<Alignments>( &(*dt_Alignments),
117  appendTime,
118  // myDbService->currentTime(),
119  "DTAlignmentRcd");
120  }
121  if(ALIUtils::debug >= 2) std::cout << "DTAlignmentRcd WRITTEN TO DB : "<< nrcd << std::endl;
122 
123  nrcd = dt_AlignmentErrors->m_alignError.size();
124  if ( myDbService->isNewTagRequest( "DTAlignmentErrorExtendedRcd" ) ) {
125  myDbService->createNewIOV<AlignmentErrorsExtended>(&(*dt_AlignmentErrors),
126  myDbService->beginOfTime(),
127  myDbService->endOfTime(),
128  "DTAlignmentErrorExtendedRcd");
129  } else {
130  myDbService->appendSinceTime<AlignmentErrorsExtended>( &(*dt_AlignmentErrors),
131  appendTime,
132  "DTAlignmentErrorExtendedRcd");
133  }
134  if(ALIUtils::debug >= 2) std::cout << "DTAlignmentErrorExtendedRcd WRITTEN TO DB : "<< nrcd << std::endl;
135 
136  // Build CSC alignments and errors
137  std::pair< Alignments*,AlignmentErrorsExtended*> cscali = BuildAlignments(0);
138  Alignments* csc_Alignments = cscali.first;
139  AlignmentErrorsExtended* csc_AlignmentErrors = cscali.second;
140 
141  // Dump CSC alignments and errors
142  nrcd = csc_Alignments->m_align.size();
143  if ( myDbService->isNewTagRequest( "CSCAlignmentRcd" ) ) {
144  myDbService->createNewIOV<Alignments>(&(*csc_Alignments),
145  myDbService->beginOfTime(),
146  myDbService->endOfTime(),
147  "CSCAlignmentRcd");
148  } else {
149  myDbService->appendSinceTime<Alignments>( &(*csc_Alignments),
150  appendTime,
151  "CSCAlignmentRcd");
152  }
153  if(ALIUtils::debug >= 2) std::cout << "CSCAlignmentRcd WRITTEN TO DB : "<< nrcd << std::endl;
154 
155  nrcd = csc_AlignmentErrors->m_alignError.size();
156  if ( myDbService->isNewTagRequest( "CSCAlignmentErrorExtendedRcd" ) ) {
157  myDbService->createNewIOV<AlignmentErrorsExtended>(&(*csc_AlignmentErrors),
158  myDbService->beginOfTime(),
159  myDbService->endOfTime(),
160  "CSCAlignmentErrorExtendedRcd");
161  } else {
162  myDbService->appendSinceTime<AlignmentErrorsExtended>( &(*csc_AlignmentErrors),
163  appendTime,
164  "CSCAlignmentErrorExtendedRcd");
165  }
166  if(ALIUtils::debug >= 2) std::cout << "CSCAlignmentErrorExtendedRcd WRITTEN TO DB : "<< nrcd << std::endl;
167 
168  //? gives unreadable error??? std::cout << "@@@@ OPTICALALIGNMENTS WRITTEN TO DB " << *optalign << std::endl;
169 
170  return TRUE;
171  }
172 
173  return TRUE;
174 }
175 
176 
177 //@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@
179 {
180  std::cout << " CocoaDBMgr::GetOptAlignInfoFromOptO " << opto->name() << std::endl;
182  data.ID_=opto->getCmsswID();
183  data.type_=opto->type();
184  data.name_=opto->name();
185 
186  //----- Centre in local coordinates
187  CLHEP::Hep3Vector centreLocal = opto->centreGlob() - opto->parent()->centreGlob();
188  CLHEP::HepRotation parentRmGlobInv = inverseOf( opto->parent()->rmGlob() );
189  centreLocal = parentRmGlobInv * centreLocal;
190 
191  const std::vector< Entry* > theCoordinateEntryVector = opto->CoordinateEntryList();
192  std::cout << " CocoaDBMgr::GetOptAlignInfoFromOptO starting coord " <<std::endl;
193 
194  data.x_.value_= centreLocal.x() / 100.; // in cm
195  std::cout << " matrix " << Fit::GetAtWAMatrix() << std::endl;
196  std::cout << " matrix " << Fit::GetAtWAMatrix()->Mat() << " " << theCoordinateEntryVector[0]->fitPos() << std::endl;
197  data.x_.error_= GetEntryError( theCoordinateEntryVector[0] ) / 100.; // in cm
198 
199  data.y_.value_= centreLocal.y() / 100.; // in cm
200  std::cout << " matrix " << Fit::GetAtWAMatrix()->Mat() << " " << theCoordinateEntryVector[1]->fitPos() << std::endl;
201  data.y_.error_= GetEntryError( theCoordinateEntryVector[1] ) / 100.; // in cm
202 
203  data.z_.value_= centreLocal.z() / 100.; // in cm
204  std::cout << " matrix " << Fit::GetAtWAMatrix()->Mat() << " " << theCoordinateEntryVector[2]->fitPos() << std::endl;
205  data.z_.error_= GetEntryError( theCoordinateEntryVector[2] ) / 100.; // in cm
206 
207  //----- angles in local coordinates
208  std::vector<double> anglocal = opto->getLocalRotationAngles( theCoordinateEntryVector );
209 
210  data.angx_.value_= anglocal[0] *180./M_PI; // in deg
211  std::cout << " matrix " << Fit::GetAtWAMatrix()->Mat() << theCoordinateEntryVector[3]->fitPos() << std::endl;
212  data.angx_.error_= GetEntryError( theCoordinateEntryVector[3] ) * 180./M_PI; // in deg;
213 
214  data.angy_.value_= anglocal[1] * 180./M_PI; // in deg
215  std::cout << " matrix " << Fit::GetAtWAMatrix()->Mat() << theCoordinateEntryVector[4]->fitPos() << std::endl;
216  data.angy_.error_= GetEntryError( theCoordinateEntryVector[4] ) * 180./M_PI; // in deg;;
217 
218  data.angz_.value_= anglocal[2] *180./M_PI; // in deg
219  std::cout << " matrix " << Fit::GetAtWAMatrix()->Mat() << theCoordinateEntryVector[5]->fitPos() << std::endl;
220  data.angz_.error_= GetEntryError( theCoordinateEntryVector[5] ) * 180./M_PI; // in deg;
221 
222 
223  const std::vector< Entry* > theExtraEntryVector = opto->ExtraEntryList(); std::cout << " CocoaDBMgr::GetOptAlignInfoFromOptO starting entry " << std::endl;
224 
225  std::vector< Entry* >::const_iterator ite;
226  for( ite = theExtraEntryVector.begin(); ite != theExtraEntryVector.end(); ite++ ) {
227  OpticalAlignParam extraEntry;
228  extraEntry.name_ = (*ite)->name();
229  extraEntry.dim_type_ = (*ite)->type();
230  extraEntry.value_ = (*ite)->value();
231  extraEntry.error_ = (*ite)->sigma();
232  extraEntry.quality_ = (*ite)->quality();
233  data.extraEntries_.push_back( extraEntry );
234  std::cout << " CocoaDBMgr::GetOptAlignInfoFromOptO done extra entry " << extraEntry.name_ << std::endl;
235 
236  }
237 
238  return data;
239 }
240 
241 
242 //-----------------------------------------------------------------------
243 double CocoaDBMgr::GetEntryError( const Entry* entry )
244 {
245  if( entry->quality() > 0 ) {
246  return sqrt(Fit::GetAtWAMatrix()->Mat()->me[entry->fitPos()][entry->fitPos()]);
247  } else { //entry not fitted, return original error
248  return entry->sigma();
249  }
250 }
251 
252 
253 //-----------------------------------------------------------------------
254 double CocoaDBMgr::GetEntryError( const Entry* entry1, const Entry* entry2 )
255 {
256  if( entry1 == entry2 ) return GetEntryError( entry1 );
257 
258  if( entry1->quality() > 0 && entry2->quality() > 0 ) {
259  return sqrt(Fit::GetAtWAMatrix()->Mat()->me[entry1->fitPos()][entry2->fitPos()]);
260  } else { //entries not fitted, correlation is 0
261  return 0.;
262  }
263 }
264 
265 
266 //-----------------------------------------------------------------------
268 {
269  OpticalAlignments* optalign= new OpticalAlignments;
270 
271  static std::vector< OpticalObject* > optolist = Model::OptOList();
272  static std::vector< OpticalObject* >::const_iterator ite;
273  for(ite = optolist.begin(); ite != optolist.end(); ite++ ){
274  if( (*ite)->type() == "system" ) continue;
276  optalign->opticalAlignments_.push_back(data);
277  if(ALIUtils::debug >= 5) {
278  std::cout << "@@@@ OPTALIGNINFO TO BE WRITTEN TO DB "
279  << data
280  << std::endl;
281  }
282  }
283  return optalign;
284 }
285 
286 
287 //-----------------------------------------------------------------------
288 std::pair< Alignments*,AlignmentErrorsExtended*> CocoaDBMgr::BuildAlignments(bool bDT)
289 {
290  Alignments* alignments = new Alignments;
291  AlignmentErrorsExtended* alignmentErrors = new AlignmentErrorsExtended;
292 
293  //read
294  static std::vector< OpticalObject* > optolist = Model::OptOList();
295  static std::vector< OpticalObject* >::const_iterator ite;
296  for(ite = optolist.begin(); ite != optolist.end(); ite++ ){
297  if( (*ite)->type() == "system" ) continue;
298  std::cout << "CocoaDBMgr::BuildAlignments getCmsswID " << (*ite) << std::endl;
299  std::cout << "CocoaDBMgr::BuildAlignments getCmsswID " << (*ite)->getCmsswID() << std::endl;
300  //check CMSSW ID
301  if( (*ite)->getCmsswID() > 0 ) { //put the numbers of DT or CSC objects
302  std::cout << " cal fill alignments " << std::endl;
303  alignments->m_align.push_back( *(GetAlignInfoFromOptO( *ite )));
304  std::cout << " fill alignments " << std::endl;
305  // AlignTransformErrorExtended* err =
306  //GetAlignInfoErrorFromOptO( *ite );
307  alignmentErrors->m_alignError.push_back(*(GetAlignInfoErrorFromOptO( *ite )));
308  std::cout << "CocoaDBMgr::BuildAlignments add alignmentError " << alignmentErrors->m_alignError.size() << std::endl;
309  }
310  }
311 
312  if(ALIUtils::debug >= 4) std::cout << "CocoaDBMgr::BuildAlignments end with n alignment " << alignments->m_align.size() << " n alignmentError " << alignmentErrors->m_alignError.size() << std::endl;
313  return std::pair< Alignments*,AlignmentErrorsExtended*>(alignments,alignmentErrors);
314 }
315 
316 
317 //-----------------------------------------------------------------------
319 {
320  if(ALIUtils::debug >= 3) std::cout << "@@@ CocoaDBMgr::GetAlignInfoFromOptO " << opto->name() << std::endl;
321 
324  align::ID cmsswID = opto->getCmsswID();
325 
326  std::cout << "@@@ CocoaDBMgr::GetAlignInfoFromOptO buildalign" << opto->name() << std::endl;
327  AlignTransform* align = new AlignTransform( trans, rot, cmsswID );
328 
329  std::cout << "@@@ CocoaDBMgr::GetAlignInfoFromOptO alig built " << opto->name() << std::endl;
330 
331  return align;
332  // return dd;
333 }
334 
335 //-----------------------------------------------------------------------
337 {
338  if(ALIUtils::debug >= 3) std::cout << "@@@ CocoaDBMgr::GetAlignInfoErrorFromOptO " << opto->name() << std::endl;
339 
340  align::ID cmsswID = opto->getCmsswID();
341 
342  GlobalError gerr(1.,
343  0.,
344  1.,
345  0.,
346  0.,
347  1.);
348  //double(dx*dx), 0., double(dy*dy), 0., 0., double(dz*dz) ) ;
349  CLHEP::HepSymMatrix errms = asHepMatrix(gerr.matrix());
350  AlignTransformErrorExtended* alignError = new AlignTransformErrorExtended( errms, cmsswID );
351  return alignError;
352 
353  CLHEP::HepMatrix errm(3,3);
354  const std::vector< Entry* > theCoordinateEntryVector = opto->CoordinateEntryList();
355 std::cout << "@@@ CocoaDBMgr::GetAlignInfoFromOptOfill errm " << opto->name() << std::endl;
356  errm(0,0) = GetEntryError( theCoordinateEntryVector[0] ) / 100.; // in cm
357  errm(1,1) = GetEntryError( theCoordinateEntryVector[1] ) / 100.; // in cm
358  errm(2,2) = GetEntryError( theCoordinateEntryVector[2] ) / 100.; // in cm
359  errm(0,1) = GetEntryError( theCoordinateEntryVector[0], theCoordinateEntryVector[1] ) / 100.; // in cm
360  errm(0,2) = GetEntryError( theCoordinateEntryVector[0], theCoordinateEntryVector[2] ) / 100.; // in cm
361  errm(1,2) = GetEntryError( theCoordinateEntryVector[1], theCoordinateEntryVector[2] ) / 100.; // in cm
362  // errm(1,0) = errm(0,1);
363  // errm(2,0) = errm(0,2);
364  // errm(2,1) = errm(1,2);
365 
366 std::cout << "@@@ CocoaDBMgr::GetAlignInfoFromOptO errm filled" << opto->name() << std::endl;
367 // CLHEP::HepSymMatrix errms(3);
368 // errms.assign(errm);
369 
370 std::cout << "@@@ CocoaDBMgr::GetAlignInfoFromOptO errms filled " << opto->name() << std::endl;
371 // AlignTransformErrorExtended* alignError = new AlignTransformErrorExtended( errms, cmsswID );
372 // AlignTransformErrorExtended* alignError = 0;
373 
374  std::cout << alignError << "@@@ CocoaDBMgr::GetAlignInfoFromOptO error built " << opto->name() << std::endl;
375  //t return alignError;
376  return (AlignTransformErrorExtended*)(0);
377 }
378 
379 
const ALIuint getCmsswID() const
CLHEP::HepMatrix asHepMatrix(const ROOT::Math::SMatrix< double, N1, N2, typename ROOT::Math::MatRepStd< double, N1, N2 > > &rm)
Definition: Migration.h:49
#define TRUE
Definition: scimark2.h:12
OpticalAlignParam x_
AlignTransformErrorExtended * GetAlignInfoErrorFromOptO(OpticalObject *opto)
Definition: CocoaDBMgr.cc:336
Definition: Entry.h:18
const std::vector< Entry * > & ExtraEntryList() const
Definition: OpticalObject.h:69
static ALIMatrix * GetAtWAMatrix()
Definition: Fit.h:146
uint32_t ID
Definition: Definitions.h:26
bool DumpCocoaResults()
Definition: CocoaDBMgr.cc:54
CLHEP::Hep3Vector Translation
const AlgebraicSymMatrix33 & matrix() const
OpticalAlignments * BuildOpticalAlignments()
Definition: CocoaDBMgr.cc:267
const std::vector< Entry * > & CoordinateEntryList() const
Definition: OpticalObject.h:65
static ALIint debug
Definition: ALIUtils.h:35
static GlobalOptionMgr * getInstance()
std::vector< AlignTransform > m_align
Definition: Alignments.h:19
const CLHEP::HepRotation & rmGlob() const
static CocoaDBMgr * instance
Definition: CocoaDBMgr.h:62
OpticalAlignParam angx_
std::vector< double > getLocalRotationAngles(const std::vector< Entry * > &entries) const
OpticalAlignParam y_
void appendSinceTime(T *payloadObj, cond::Time_t sinceTime, const std::string &recordName, bool withlogging=false)
unsigned long long Time_t
Definition: Time.h:16
bool isNewTagRequest(const std::string &recordName)
std::vector< OpticalAlignInfo > opticalAlignments_
T sqrt(T t)
Definition: SSEVec.h:48
const OpticalObject * parent() const
Definition: OpticalObject.h:62
bool isAvailable() const
Definition: Service.h:46
const MAT * Mat() const
static ALIuint nEvent
Definition: Fit.h:204
AlignTransform * GetAlignInfoFromOptO(OpticalObject *opto)
Definition: CocoaDBMgr.cc:318
double GetEntryError(const Entry *entry)
Definition: CocoaDBMgr.cc:243
#define M_PI
OpticalAlignInfo GetOptAlignInfoFromOptO(OpticalObject *opto)
Definition: CocoaDBMgr.cc:178
void createNewIOV(T *firstPayloadObj, cond::Time_t firstSinceTime, cond::Time_t firstTillTime, const std::string &recordName, bool withlogging=false)
OpticalAlignParam angz_
ALIint fitPos() const
Definition: Entry.h:60
static std::vector< OpticalObject * > & OptOList()
Definition: Model.h:71
OpticalAlignParam z_
std::vector< AlignTransformErrorExtended > m_alignError
static CocoaDBMgr * getInstance()
Definition: CocoaDBMgr.cc:40
std::vector< OpticalAlignParam > extraEntries_
const CLHEP::Hep3Vector & centreGlob() const
Definition: OpticalObject.h:85
OpticalAlignParam angy_
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
ALIint quality() const
Definition: Entry.h:59
tuple cout
Definition: gather_cfg.py:121
const ALIstring & name() const
Definition: OpticalObject.h:60
std::pair< Alignments *, AlignmentErrorsExtended * > BuildAlignments(bool bDT)
Definition: CocoaDBMgr.cc:288
std::string dim_type_
std::map< ALIstring, ALIdouble, std::less< ALIstring > > & GlobalOptions()
unsigned int ID_
const ALIstring & type() const
Definition: OpticalObject.h:61
ALIdouble sigma() const
Definition: Entry.h:57
CLHEP::HepRotation Rotation