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