Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:55:59

0001 #include "CondFormats/Alignment/interface/Alignments.h"
0002 #include "CondFormats/Alignment/interface/AlignmentErrorsExtended.h"
0003 #include "CondFormats/Alignment/interface/AlignTransform.h"
0004 #include "CondFormats/Alignment/interface/AlignTransformErrorExtended.h"
0005 #include "DataFormats/GeometryCommonDetAlgo/interface/GlobalError.h"
0006 #include <DD4hep/DD4hepUnits.h>
0007 #include "FWCore/ServiceRegistry/interface/Service.h"
0008 #include "CondCore/DBOutputService/interface/PoolDBOutputService.h"
0009 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0010 #include "FWCore/Framework/interface/ESHandle.h"
0011 
0012 #include "Alignment/CocoaFit/interface/CocoaDBMgr.h"
0013 #include "CondFormats/OptAlignObjects/interface/OpticalAlignMeasurementInfo.h"
0014 #include "CondFormats/DataRecord/interface/OpticalAlignmentsRcd.h"
0015 #include "CondFormats/AlignmentRecord/interface/DTAlignmentRcd.h"
0016 #include "Geometry/Records/interface/IdealGeometryRecord.h"
0017 #include "FWCore/Framework/interface/Event.h"
0018 #include "FWCore/Framework/interface/EventSetup.h"
0019 
0020 #include "Alignment/CocoaUtilities/interface/ALIUtils.h"
0021 #include "Alignment/CocoaModel/interface/Model.h"
0022 #include "Alignment/CocoaFit/interface/Fit.h"
0023 #include "Alignment/CocoaModel/interface/Entry.h"
0024 #include "Alignment/CocoaUtilities/interface/ALIFileOut.h"
0025 #include "Alignment/CocoaModel/interface/CocoaDaqReaderRoot.h"
0026 #include "Alignment/CocoaModel/interface/OpticalObject.h"
0027 #include "Alignment/CocoaUtilities/interface/GlobalOptionMgr.h"
0028 
0029 #include "CondFormats/OptAlignObjects/interface/OpticalAlignments.h"
0030 #include "CondFormats/OptAlignObjects/interface/OpticalAlignInfo.h"
0031 #include "CondFormats/OptAlignObjects/interface/OpticalAlignMeasurements.h"
0032 
0033 #include "CondCore/CondDB/interface/Serialization.h"
0034 
0035 CocoaDBMgr* CocoaDBMgr::instance = nullptr;
0036 
0037 //----------------------------------------------------------------------
0038 CocoaDBMgr* CocoaDBMgr::getInstance() {
0039   if (!instance) {
0040     instance = new CocoaDBMgr;
0041   }
0042   return instance;
0043 }
0044 
0045 //----------------------------------------------------------------------
0046 CocoaDBMgr::CocoaDBMgr() {}
0047 
0048 //-----------------------------------------------------------------------
0049 bool CocoaDBMgr::DumpCocoaResults() {
0050   edm::Service<cond::service::PoolDBOutputService> myDbService;
0051 
0052   GlobalOptionMgr* gomgr = GlobalOptionMgr::getInstance();
0053   int nrcd;
0054 
0055   cond::Time_t appendTime = Fit::nEvent + 1;
0056   if (gomgr->GlobalOptions()["writeDBOptAlign"] > 0) {
0057     //----- Build OpticalAlignments
0058     std::unique_ptr<OpticalAlignments> optalign = BuildOpticalAlignments();
0059 
0060     //--- Dump OpticalAlignments
0061     nrcd = optalign->opticalAlignments_.size();
0062     if (!myDbService.isAvailable()) {
0063       throw cms::Exception("CocoaDBMgr::DumpCocoaResults DB not available");
0064     }
0065     //    try{
0066     if (myDbService->isNewTagRequest("OpticalAlignmentsRcd")) {
0067       std::cout << " new OA to DB "
0068                 << "begin " << myDbService->beginOfTime() << " current " << myDbService->currentTime() << " end "
0069                 << myDbService->endOfTime() << std::endl;
0070       myDbService->createOneIOV<OpticalAlignments>(*optalign, myDbService->beginOfTime(), "OpticalAlignmentsRcd");
0071     } else {
0072       std::cout << " old OA to DB "
0073                 << " current " << myDbService->currentTime() << " end " << myDbService->endOfTime() << std::endl;
0074       myDbService->appendOneIOV<OpticalAlignments>(*optalign, appendTime, "OpticalAlignmentsRcd");
0075     }
0076 
0077     /*    }catch(const cond::Exception& er) {
0078       std::cout<<er.what()<<std::endl;
0079       }catch(const std::exception& er){
0080       std::cout<<"caught std::exception "<<er.what()<<std::endl;
0081       }catch(...){
0082       std::cout<<"Funny error"<<std::endl;
0083       } */
0084 
0085     if (ALIUtils::debug >= 2)
0086       std::cout << "OpticalAlignmentsRcd WRITTEN TO DB : " << nrcd << std::endl;
0087   }
0088 
0089   if (gomgr->GlobalOptions()["writeDBAlign"] > 0) {
0090     // Build DT alignments and errors
0091     const auto& dtali = BuildAlignments(true);
0092     auto& dt_Alignments = dtali.first;
0093     auto& dt_AlignmentErrors = dtali.second;
0094 
0095     // Dump DT alignments and errors
0096     nrcd = dt_Alignments->m_align.size();
0097     if (myDbService->isNewTagRequest("DTAlignmentRcd")) {
0098       myDbService->createOneIOV<Alignments>(*dt_Alignments, myDbService->beginOfTime(), "DTAlignmentRcd");
0099     } else {
0100       myDbService->appendOneIOV<Alignments>(*dt_Alignments, appendTime, "DTAlignmentRcd");
0101     }
0102     if (ALIUtils::debug >= 2)
0103       std::cout << "DTAlignmentRcd WRITTEN TO DB : " << nrcd << std::endl;
0104 
0105     nrcd = dt_AlignmentErrors->m_alignError.size();
0106     if (myDbService->isNewTagRequest("DTAlignmentErrorExtendedRcd")) {
0107       myDbService->createOneIOV<AlignmentErrorsExtended>(
0108           *dt_AlignmentErrors, myDbService->beginOfTime(), "DTAlignmentErrorExtendedRcd");
0109     } else {
0110       myDbService->appendOneIOV<AlignmentErrorsExtended>(
0111           *dt_AlignmentErrors, appendTime, "DTAlignmentErrorExtendedRcd");
0112     }
0113     if (ALIUtils::debug >= 2)
0114       std::cout << "DTAlignmentErrorExtendedRcd WRITTEN TO DB : " << nrcd << std::endl;
0115 
0116     // Build CSC alignments and errors
0117     const auto& cscali = BuildAlignments(false);
0118     auto& csc_Alignments = cscali.first;
0119     auto& csc_AlignmentErrors = cscali.second;
0120 
0121     // Dump CSC alignments and errors
0122     nrcd = csc_Alignments->m_align.size();
0123     if (myDbService->isNewTagRequest("CSCAlignmentRcd")) {
0124       myDbService->createOneIOV<Alignments>(*csc_Alignments, myDbService->beginOfTime(), "CSCAlignmentRcd");
0125     } else {
0126       myDbService->appendOneIOV<Alignments>(*csc_Alignments, appendTime, "CSCAlignmentRcd");
0127     }
0128     if (ALIUtils::debug >= 2)
0129       std::cout << "CSCAlignmentRcd WRITTEN TO DB : " << nrcd << std::endl;
0130 
0131     nrcd = csc_AlignmentErrors->m_alignError.size();
0132     if (myDbService->isNewTagRequest("CSCAlignmentErrorExtendedRcd")) {
0133       myDbService->createOneIOV<AlignmentErrorsExtended>(
0134           *csc_AlignmentErrors, myDbService->beginOfTime(), "CSCAlignmentErrorExtendedRcd");
0135     } else {
0136       myDbService->appendOneIOV<AlignmentErrorsExtended>(
0137           *csc_AlignmentErrors, appendTime, "CSCAlignmentErrorExtendedRcd");
0138     }
0139     if (ALIUtils::debug >= 2)
0140       std::cout << "CSCAlignmentErrorExtendedRcd WRITTEN TO DB : " << nrcd << std::endl;
0141 
0142     //? gives unreadable error???  std::cout << "@@@@ OPTICALALIGNMENTS WRITTEN TO DB " << *optalign << std::endl;
0143 
0144     return TRUE;
0145   }
0146 
0147   return TRUE;
0148 }
0149 
0150 //@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@
0151 OpticalAlignInfo CocoaDBMgr::GetOptAlignInfoFromOptO(OpticalObject* opto) {
0152   LogDebug("Alignment") << " CocoaDBMgr::GetOptAlignInfoFromOptO " << opto->name();
0153   OpticalAlignInfo data;
0154   data.ID_ = opto->getCmsswID();
0155   data.type_ = opto->type();
0156   data.name_ = opto->name();
0157 
0158   //----- Centre in local coordinates
0159   CLHEP::Hep3Vector centreLocal = opto->centreGlob() - opto->parent()->centreGlob();
0160   CLHEP::HepRotation parentRmGlobInv = inverseOf(opto->parent()->rmGlob());
0161   centreLocal = parentRmGlobInv * centreLocal;
0162 
0163   const std::vector<Entry*>& theCoordinateEntryVector = opto->CoordinateEntryList();
0164   LogDebug("Alignment") << " CocoaDBMgr::GetOptAlignInfoFromOptO starting coord ";
0165   if (theCoordinateEntryVector.size() == 6) {
0166     const Entry* const translationX = theCoordinateEntryVector.at(0);
0167     OpticalAlignParam translationXDataForDB;
0168     translationXDataForDB.name_ = translationX->name();
0169     translationXDataForDB.dim_type_ = translationX->type();
0170     translationXDataForDB.value_ = centreLocal.x() * dd4hep::m;              // m in COCOA, dd4hep unit in DB
0171     translationXDataForDB.error_ = GetEntryError(translationX) * dd4hep::m;  // m in COCOA, dd4hep unit in DB
0172     translationXDataForDB.quality_ = translationX->quality();
0173     data.x_ = translationXDataForDB;
0174 
0175     const Entry* const translationY = theCoordinateEntryVector.at(1);
0176     OpticalAlignParam translationYDataForDB;
0177     translationYDataForDB.name_ = translationY->name();
0178     translationYDataForDB.dim_type_ = translationY->type();
0179     translationYDataForDB.value_ = centreLocal.y() * dd4hep::m;              // m in COCOA, dd4hep unit in DB
0180     translationYDataForDB.error_ = GetEntryError(translationY) * dd4hep::m;  // m in COCOA, dd4hep unit in DB
0181     translationYDataForDB.quality_ = translationY->quality();
0182     data.y_ = translationYDataForDB;
0183 
0184     const Entry* const translationZ = theCoordinateEntryVector.at(2);
0185     OpticalAlignParam translationZDataForDB;
0186     translationZDataForDB.name_ = translationZ->name();
0187     translationZDataForDB.dim_type_ = translationZ->type();
0188     translationZDataForDB.value_ = centreLocal.z() * dd4hep::m;              // m in COCOA, dd4hep unit in DB
0189     translationZDataForDB.error_ = GetEntryError(translationZ) * dd4hep::m;  // m in COCOA, dd4hep unit in DB
0190     translationZDataForDB.quality_ = translationZ->quality();
0191     data.z_ = translationZDataForDB;
0192 
0193     //----- angles in local coordinates
0194     std::vector<double> anglocal = opto->getLocalRotationAngles(theCoordinateEntryVector);
0195     if (anglocal.size() == 3) {
0196       const Entry* const rotationX = theCoordinateEntryVector.at(3);
0197       OpticalAlignParam rotationXDataForDB;
0198       rotationXDataForDB.name_ = rotationX->name();
0199       rotationXDataForDB.dim_type_ = rotationX->type();
0200       rotationXDataForDB.value_ = anglocal.at(0);
0201       rotationXDataForDB.error_ = GetEntryError(rotationX);
0202       rotationXDataForDB.quality_ = rotationX->quality();
0203       data.angx_ = rotationXDataForDB;
0204 
0205       const Entry* const rotationY = theCoordinateEntryVector.at(4);
0206       OpticalAlignParam rotationYDataForDB;
0207       rotationYDataForDB.name_ = rotationY->name();
0208       rotationYDataForDB.dim_type_ = rotationY->type();
0209       rotationYDataForDB.value_ = anglocal.at(1);
0210       rotationYDataForDB.error_ = GetEntryError(rotationY);
0211       rotationYDataForDB.quality_ = rotationY->quality();
0212       data.angy_ = rotationYDataForDB;
0213 
0214       const Entry* const rotationZ = theCoordinateEntryVector.at(5);
0215       OpticalAlignParam rotationZDataForDB;
0216       rotationZDataForDB.name_ = rotationZ->name();
0217       rotationZDataForDB.dim_type_ = rotationZ->type();
0218       rotationZDataForDB.value_ = anglocal.at(2);
0219       rotationZDataForDB.error_ = GetEntryError(rotationZ);
0220       rotationZDataForDB.quality_ = rotationZ->quality();
0221       data.angz_ = rotationZDataForDB;
0222     }
0223   }
0224 
0225   std::cout << " CocoaDBMgr::GetOptAlignInfoFromOptO starting entry " << std::endl;
0226   for (const auto& myDBExtraEntry : opto->ExtraEntryList()) {
0227     OpticalAlignParam extraEntry;
0228     extraEntry.name_ = myDBExtraEntry->name();
0229     extraEntry.dim_type_ = myDBExtraEntry->type();
0230     extraEntry.value_ = myDBExtraEntry->value();
0231     extraEntry.error_ = myDBExtraEntry->sigma();
0232     if (extraEntry.dim_type_ == "centre" || extraEntry.dim_type_ == "length") {
0233       extraEntry.value_ *= dd4hep::m;  // m in COCOA, dd4hep unit in DB
0234       extraEntry.error_ *= dd4hep::m;  // m in COCOA, dd4hep unit in DB
0235     }
0236     extraEntry.quality_ = myDBExtraEntry->quality();
0237     data.extraEntries_.emplace_back(extraEntry);
0238     std::cout << " CocoaDBMgr::GetOptAlignInfoFromOptO done extra entry " << extraEntry.name_ << std::endl;
0239   }
0240 
0241   return data;
0242 }
0243 
0244 //-----------------------------------------------------------------------
0245 double CocoaDBMgr::GetEntryError(const Entry* entry) {
0246   if (entry->quality() > 0) {
0247     return sqrt(Fit::GetAtWAMatrix()->Mat()->me[entry->fitPos()][entry->fitPos()]);
0248   } else {  //entry not fitted, return original error
0249     return entry->sigma();
0250   }
0251 }
0252 
0253 //-----------------------------------------------------------------------
0254 double CocoaDBMgr::GetEntryError(const Entry* entry1, const Entry* entry2) {
0255   if (entry1 == entry2)
0256     return GetEntryError(entry1);
0257 
0258   if (entry1->quality() > 0 && entry2->quality() > 0) {
0259     return sqrt(Fit::GetAtWAMatrix()->Mat()->me[entry1->fitPos()][entry2->fitPos()]);
0260   } else {  //entries not fitted, correlation is 0
0261     return 0.;
0262   }
0263 }
0264 
0265 //-----------------------------------------------------------------------
0266 std::unique_ptr<OpticalAlignments> CocoaDBMgr::BuildOpticalAlignments() {
0267   std::unique_ptr<OpticalAlignments> optalign = std::make_unique<OpticalAlignments>();
0268 
0269   static std::vector<OpticalObject*> optolist = Model::OptOList();
0270   static std::vector<OpticalObject*>::const_iterator ite;
0271   for (ite = optolist.begin(); ite != optolist.end(); ++ite) {
0272     if ((*ite)->type() == "system")
0273       continue;
0274     OpticalAlignInfo data = GetOptAlignInfoFromOptO(*ite);
0275     optalign->opticalAlignments_.push_back(data);
0276     if (ALIUtils::debug >= 5) {
0277       std::cout << "@@@@ OPTALIGNINFO TO BE WRITTEN TO DB " << data << std::endl;
0278     }
0279   }
0280   return optalign;
0281 }
0282 
0283 //-----------------------------------------------------------------------
0284 std::pair<std::unique_ptr<Alignments>, std::unique_ptr<AlignmentErrorsExtended> > CocoaDBMgr::BuildAlignments(bool bDT) {
0285   std::unique_ptr<Alignments> alignments = std::make_unique<Alignments>();
0286   std::unique_ptr<AlignmentErrorsExtended> alignmentErrors = std::make_unique<AlignmentErrorsExtended>();
0287 
0288   //read
0289   static std::vector<OpticalObject*> optolist = Model::OptOList();
0290   static std::vector<OpticalObject*>::const_iterator ite;
0291   for (ite = optolist.begin(); ite != optolist.end(); ++ite) {
0292     if ((*ite)->type() == "system")
0293       continue;
0294     std::cout << "CocoaDBMgr::BuildAlignments getCmsswID " << (*ite) << std::endl;
0295     std::cout << "CocoaDBMgr::BuildAlignments getCmsswID " << (*ite)->getCmsswID() << std::endl;
0296     //check CMSSW ID
0297     if ((*ite)->getCmsswID() > 0) {  //put the numbers of DT or CSC objects
0298       std::cout << " cal fill alignments " << std::endl;
0299       alignments->m_align.push_back(*(GetAlignInfoFromOptO(*ite)));
0300       std::cout << " fill alignments " << std::endl;
0301       //      AlignTransformErrorExtended* err =
0302       //GetAlignInfoErrorFromOptO( *ite );
0303       alignmentErrors->m_alignError.push_back(*(GetAlignInfoErrorFromOptO(*ite)));
0304       std::cout << "CocoaDBMgr::BuildAlignments add alignmentError " << alignmentErrors->m_alignError.size()
0305                 << std::endl;
0306     }
0307   }
0308 
0309   if (ALIUtils::debug >= 4)
0310     std::cout << "CocoaDBMgr::BuildAlignments end with n alignment " << alignments->m_align.size()
0311               << " n alignmentError " << alignmentErrors->m_alignError.size() << std::endl;
0312 
0313   return std::make_pair(std::move(alignments), std::move(alignmentErrors));
0314 }
0315 
0316 //-----------------------------------------------------------------------
0317 AlignTransform* CocoaDBMgr::GetAlignInfoFromOptO(OpticalObject* opto) {
0318   if (ALIUtils::debug >= 3)
0319     std::cout << "@@@ CocoaDBMgr::GetAlignInfoFromOptO " << opto->name() << std::endl;
0320 
0321   const AlignTransform::Translation& trans = opto->centreGlob();
0322   const AlignTransform::Rotation& rot = opto->rmGlob();
0323   align::ID cmsswID = opto->getCmsswID();
0324 
0325   std::cout << "@@@ CocoaDBMgr::GetAlignInfoFromOptO buildalign" << opto->name() << std::endl;
0326   AlignTransform* align = new AlignTransform(trans, rot, cmsswID);
0327 
0328   std::cout << "@@@ CocoaDBMgr::GetAlignInfoFromOptO alig built " << opto->name() << std::endl;
0329 
0330   return align;
0331   //  return dd;
0332 }
0333 
0334 //-----------------------------------------------------------------------
0335 AlignTransformErrorExtended* CocoaDBMgr::GetAlignInfoErrorFromOptO(OpticalObject* opto) {
0336   if (ALIUtils::debug >= 3)
0337     std::cout << "@@@ CocoaDBMgr::GetAlignInfoErrorFromOptO " << opto->name() << std::endl;
0338 
0339   align::ID cmsswID = opto->getCmsswID();
0340 
0341   GlobalError gerr(1., 0., 1., 0., 0., 1.);
0342   //double(dx*dx),  0., double(dy*dy),     0., 0., double(dz*dz) ) ;
0343   CLHEP::HepSymMatrix errms = asHepMatrix(gerr.matrix());
0344   AlignTransformErrorExtended* alignError = new AlignTransformErrorExtended(errms, cmsswID);
0345   return alignError;
0346 
0347   CLHEP::HepMatrix errm(3, 3);
0348   const std::vector<Entry*>& theCoordinateEntryVector = opto->CoordinateEntryList();
0349   std::cout << "@@@ CocoaDBMgr::GetAlignInfoFromOptOfill errm " << opto->name() << std::endl;
0350   errm(0, 0) = GetEntryError(theCoordinateEntryVector[0]) * dd4hep::m;  // m in COCOA, dd4hep unit in DB
0351   errm(1, 1) = GetEntryError(theCoordinateEntryVector[1]) * dd4hep::m;  // m in COCOA, dd4hep unit in DB
0352   errm(2, 2) = GetEntryError(theCoordinateEntryVector[2]) * dd4hep::m;  // m in COCOA, dd4hep unit in DB
0353   errm(0, 1) = GetEntryError(theCoordinateEntryVector[0], theCoordinateEntryVector[1]) *
0354                dd4hep::m;  // m in COCOA, dd4hep unit in DB
0355   errm(0, 2) = GetEntryError(theCoordinateEntryVector[0], theCoordinateEntryVector[2]) *
0356                dd4hep::m;  // m in COCOA, dd4hep unit in DB
0357   errm(1, 2) = GetEntryError(theCoordinateEntryVector[1], theCoordinateEntryVector[2]) *
0358                dd4hep::m;  // m in COCOA, dd4hep unit in DB
0359   //   errm(1,0) = errm(0,1);
0360   // errm(2,0) = errm(0,2);
0361   // errm(2,1) = errm(1,2);
0362 
0363   std::cout << "@@@ CocoaDBMgr::GetAlignInfoFromOptO errm filled" << opto->name() << std::endl;
0364   //  CLHEP::HepSymMatrix errms(3);
0365   //  errms.assign(errm);
0366 
0367   std::cout << "@@@ CocoaDBMgr::GetAlignInfoFromOptO errms filled " << opto->name() << std::endl;
0368   //  AlignTransformErrorExtended* alignError = new AlignTransformErrorExtended( errms, cmsswID );
0369   //  AlignTransformErrorExtended* alignError = 0;
0370 
0371   std::cout << alignError << "@@@ CocoaDBMgr::GetAlignInfoFromOptO error built " << opto->name() << std::endl;
0372   //t  return alignError;
0373   return (AlignTransformErrorExtended*)nullptr;
0374 }