CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch2/src/Alignment/CocoaApplication/src/CocoaAnalyzer.cc

Go to the documentation of this file.
00001 #include "FWCore/ServiceRegistry/interface/Service.h"
00002 #include "CondCore/DBCommon/interface/Exception.h"
00003 #include "CondCore/DBOutputService/interface/PoolDBOutputService.h"
00004 #include "FWCore/ParameterSet/interface/ParameterSet.h" 
00005 #include "FWCore/Framework/interface/ESHandle.h"
00006 #include "FWCore/Framework/interface/ESTransientHandle.h"
00007 
00008 #include "../interface/CocoaAnalyzer.h"
00009 #include "CondFormats/OptAlignObjects/interface/OpticalAlignMeasurementInfo.h" 
00010 #include "CondFormats/DataRecord/interface/OpticalAlignmentsRcd.h" 
00011 #include "DetectorDescription/Core/interface/DDFilteredView.h"
00012 #include "DetectorDescription/Core/interface/DDCompactView.h"
00013 #include "Geometry/Records/interface/IdealGeometryRecord.h"
00014 #include "FWCore/Framework/interface/Event.h"
00015 #include "FWCore/Framework/interface/EventSetup.h"
00016 
00017 #include "Alignment/CocoaUtilities/interface/ALIUtils.h"
00018 #include "Alignment/CocoaModel/interface/Model.h"
00019 #include "Alignment/CocoaFit/interface/Fit.h"
00020 #include "Alignment/CocoaModel/interface/Entry.h"
00021 #include "Alignment/CocoaUtilities/interface/ALIFileOut.h"
00022 #include "Alignment/CocoaModel/interface/CocoaDaqReaderRoot.h"
00023 #include "Alignment/CocoaModel/interface/OpticalObject.h"
00024 #include "Alignment/CocoaUtilities/interface/GlobalOptionMgr.h"
00025 #include "Alignment/CocoaFit/interface/CocoaDBMgr.h"
00026 
00027 //----------------------------------------------------------------------
00028 CocoaAnalyzer::CocoaAnalyzer(edm::ParameterSet const& pset) 
00029 {
00030   theCocoaDaqRootFileName = pset.getParameter< std::string >("cocoaDaqRootFile");
00031 
00032   int maxEvents = pset.getParameter< int32_t >("maxEvents");
00033   GlobalOptionMgr::getInstance()->setDefaultGlobalOptions();
00034   GlobalOptionMgr::getInstance()->setGlobalOption("maxEvents",maxEvents);
00035   GlobalOptionMgr::getInstance()->setGlobalOption("writeDBAlign",1);
00036   GlobalOptionMgr::getInstance()->setGlobalOption("writeDBOptAlign",1);
00037 }
00038 
00039 //----------------------------------------------------------------------
00040 void CocoaAnalyzer::beginJob()
00041 {
00042 }
00043 
00044 
00045 //------------------------------------------------------------------------
00046 void CocoaAnalyzer::RunCocoa()
00047 {
00048   if(ALIUtils::debug >= 3) {
00049     std::cout << std::endl << "$$$ CocoaAnalyzer::RunCocoa: " << std::endl;
00050   }
00051   //-ALIFileIn fin;
00052   //-  GlobalOptionMgr::getInstance()->setGlobalOption("debug_verbose",5, fin );
00053 
00054   //---------- Build the Model out of the system description text file
00055   Model& model = Model::getInstance();
00056 
00057   model.BuildSystemDescriptionFromOA( oaList_ );
00058 
00059   if(ALIUtils::debug >= 3) {
00060     std::cout << "$$ CocoaAnalyzer::RunCocoa: geometry built " << std::endl;
00061   }
00062 
00063   model.BuildMeasurementsFromOA( measList_ );
00064 
00065   if(ALIUtils::debug >= 3) {
00066     std::cout << "$$ CocoaAnalyzer::RunCocoa: measurements built " << std::endl;
00067   }
00068 
00069   Fit::getInstance();
00070 
00071   Fit::startFit();
00072 
00073   if(ALIUtils::debug >= 0) std::cout << "............ program ended OK" << std::endl;
00074   if( ALIUtils::report >=1 ) {
00075     ALIFileOut& fileout = ALIFileOut::getInstance( Model::ReportFName() );
00076     fileout << "............ program ended OK" << std::endl;
00077   }
00078 
00079 } // end of ::beginJob
00080 
00081 
00082 //-----------------------------------------------------------------------
00083 void CocoaAnalyzer::ReadXMLFile( const edm::EventSetup& evts ) 
00084 {
00085 
00086   // STEP ONE:  Initial COCOA objects will be built from a DDL geometry
00087   // description.  
00088   
00089   edm::ESTransientHandle<DDCompactView> cpv;
00090   evts.get<IdealGeometryRecord>().get(cpv);
00091 
00092   if(ALIUtils::debug >= 3) {
00093     std::cout << std::endl << "$$$ CocoaAnalyzer::ReadXML: root object= " << cpv->root() << std::endl;
00094   }
00095   
00096   //Build OpticalAlignInfo "system"
00097   const DDLogicalPart lv = cpv->root();
00098   
00099   OpticalAlignInfo oaInfo;
00100   oaInfo.ID_ = 0;
00101   //--- substract file name to object name
00102   oaInfo.name_ = lv.name().name();
00103   oaInfo.parentName_ = "";
00104   oaInfo.x_.quality_  = 0;    
00105   oaInfo.x_.value_ = 0.;
00106   oaInfo.x_.error_ = 0.;
00107   oaInfo.x_.quality_  = 0;
00108   oaInfo.y_.value_ = 0.;
00109   oaInfo.y_.error_ = 0.;
00110   oaInfo.y_.quality_  = 0;
00111   oaInfo.z_.value_ = 0.;
00112   oaInfo.z_.error_ = 0.;
00113   oaInfo.z_.quality_  = 0;
00114   oaInfo.angx_.value_ = 0.;
00115   oaInfo.angx_.error_ = 0.;
00116   oaInfo.angx_.quality_  = 0;
00117   oaInfo.angy_.value_ = 0.;
00118   oaInfo.angy_.error_ = 0.;
00119   oaInfo.angy_.quality_  = 0;
00120   oaInfo.angz_.value_ = 0.;
00121   oaInfo.angz_.error_ = 0.;
00122   oaInfo.angz_.quality_  = 0;
00123     
00124   oaInfo.type_ = "system";
00125 
00126   oaList_.opticalAlignments_.push_back(oaInfo);
00127   oaInfo.clear();
00128 
00129   // Example of traversing the whole optical alignment geometry.
00130   // At each node we get specpars as variables and use them in 
00131   // constructing COCOA objects. 
00132   //  It stores these objects in a private data member, opt
00133   std::string attribute = "COCOA"; 
00134   std::string value     = "COCOA";
00135   DDValue val(attribute, value, 0.0);
00136   
00137           // get all parts labelled with COCOA using a SpecPar
00138   DDSpecificsFilter filter;
00139   filter.setCriteria(val,  // name & value of a variable 
00140                      DDSpecificsFilter::matches,
00141                      DDSpecificsFilter::AND, 
00142                      true, // compare strings otherwise doubles
00143                      true  // use merged-specifics or simple-specifics
00144                      );
00145   DDFilteredView fv(*cpv);
00146   fv.addFilter(filter);
00147   bool doCOCOA = fv.firstChild();
00148   
00149   // Loop on parts
00150   int nObjects=0;
00151   OpticalAlignParam oaParam;
00152   OpticalAlignMeasurementInfo oaMeas;
00153 
00154   while ( doCOCOA ){
00155     ++nObjects;
00156     //    oaInfo.ID_ = nObjects;
00157     const DDsvalues_type params(fv.mergedSpecifics());
00158     
00159     const DDLogicalPart lv = fv.logicalPart();
00160     if(ALIUtils::debug >= 4) {
00161       std::cout << " CocoaAnalyzer::ReadXML reading object " << lv.name() << std::endl;
00162     }
00163 
00164     std::vector<DDExpandedNode> history = fv.geoHistory();
00165     oaInfo.parentName_ = "";
00166     size_t ii;
00167     for(ii = 0; ii < history.size()-1;ii++ ) {
00168       if( ii != 0 ) oaInfo.parentName_ += "/";
00169       std::string name = history[ii].logicalPart().name().name();
00170       oaInfo.parentName_ += name;
00171  //    oaInfo.parentName_ = (fv.geoHistory()[fv.geoHistory().size()-2]).logicalPart().name();
00172 //    icol = oaInfo.parentName_.find(":");
00173  //   oaInfo.parentName_ = oaInfo.parentName_.substr(icol+1,oaInfo.parentName_.length());
00174     }
00175 
00176     //--- build object name (= parent name + object name)
00177     std::string name = history[ii].logicalPart().name().name();
00178     //--- substract file name to object name
00179     oaInfo.name_ = oaInfo.parentName_ + "/" + name;
00180     if(ALIUtils::debug >= 5) {
00181       std::cout << " @@ Name built= " << oaInfo.name_ << " short_name= " << name << " parent= " << oaInfo.parentName_ << std::endl; 
00182     }
00183     //----- Read centre and angles
00184     oaInfo.x_.quality_  = int (myFetchDbl(params, "centre_X_quality", 0));
00185     DDTranslation transl = (fv.translation());
00186     DDRotationMatrix rot = (fv.rotation());
00187     DDExpandedNode parent = fv.geoHistory()[ fv.geoHistory().size()-2 ];
00188     DDTranslation parentTransl = parent.absTranslation();
00189     DDRotationMatrix parentRot = parent.absRotation();
00190     transl = parentRot.Inverse()*(transl - parentTransl );
00191     rot = parentRot.Inverse()*rot;
00192     rot = rot.Inverse(); //DDL uses opposite convention than COCOA
00193     /*    if(ALIUtils::debug >= 4) {
00194       ALIUtils::dumprm( rot, "local rotation ");
00195       ALIUtils::dump3v( transl, "local translation");
00196       } */
00197 
00198     oaInfo.x_.name_ = "X";
00199     oaInfo.x_.dim_type_ = "centre";
00200     oaInfo.x_.value_ = transl.x()*0.001; // CLHEP units are mm, COCOA are m 
00201     oaInfo.x_.error_ = myFetchDbl(params, "centre_X_sigma", 0)*0.001; // CLHEP units are mm, COCOA are m 
00202     oaInfo.x_.quality_  = int (myFetchDbl(params, "centre_X_quality", 0));
00203     
00204     oaInfo.y_.name_ = "Y";
00205     oaInfo.y_.dim_type_ = "centre";
00206     oaInfo.y_.value_ = transl.y()*0.001; // CLHEP units are mm, COCOA are m 
00207     oaInfo.y_.error_ = myFetchDbl(params, "centre_Y_sigma", 0)*0.001; // CLHEP units are mm, COCOA are m 
00208     oaInfo.y_.quality_  = int (myFetchDbl(params, "centre_Y_quality", 0));
00209 
00210     oaInfo.z_.name_ = "Z";
00211     oaInfo.z_.dim_type_ = "centre";
00212     oaInfo.z_.value_ = transl.z()*0.001; // CLHEP units are mm, COCOA are m 
00213     oaInfo.z_.error_ = myFetchDbl(params, "centre_Z_sigma", 0)*0.001; // CLHEP units are mm, COCOA are m 
00214     oaInfo.z_.quality_  = int (myFetchDbl(params, "centre_Z_quality", 0));
00215 
00216   //---- DDD convention is to use the inverse matrix, COCOA is the direct one!!!
00217     //---- convert it to CLHEP::Matrix
00218     double xx,xy,xz,yx,yy,yz,zx,zy,zz;
00219     rot.GetComponents (xx, xy, xz,
00220                  yx, yy, yz,
00221                  zx, zy, zz);
00222     CLHEP::Hep3Vector colX(xx,xy,xz);
00223     CLHEP::Hep3Vector colY(yx,yy,yz);
00224     CLHEP::Hep3Vector colZ(zx,zy,zz);
00225     CLHEP::HepRotation rotclhep( colX, colY, colZ );
00226     std::vector<double> angles = ALIUtils::getRotationAnglesFromMatrix( rotclhep,0., 0., 0. );
00227 
00228     oaInfo.angx_.name_ = "X";
00229     oaInfo.angx_.dim_type_ = "angles";
00230     //-    oaInfo.angx_.value_ = angles[0];
00231     oaInfo.angx_.value_ = myFetchDbl(params, "angles_X_value", 0);
00232     oaInfo.angx_.error_ = myFetchDbl(params, "angles_X_sigma", 0);
00233     oaInfo.angx_.quality_  = int (myFetchDbl(params, "angles_X_quality", 0));
00234 
00235     oaInfo.angy_.name_ = "Y";
00236     oaInfo.angy_.dim_type_ = "angles";
00237     //-    oaInfo.angy_.value_ = angles[1];
00238     oaInfo.angy_.value_ = myFetchDbl(params, "angles_Y_value", 0);
00239     oaInfo.angy_.error_ = myFetchDbl(params, "angles_Y_sigma", 0);
00240     oaInfo.angy_.quality_  = int (myFetchDbl(params, "angles_Y_quality", 0));
00241 
00242     oaInfo.angz_.name_ = "Z";
00243     oaInfo.angz_.dim_type_ = "angles";
00244     //    oaInfo.angz_.value_ = angles[2];
00245     oaInfo.angz_.value_ = myFetchDbl(params, "angles_Z_value", 0);
00246     oaInfo.angz_.error_ = myFetchDbl(params, "angles_Z_sigma", 0);
00247     oaInfo.angz_.quality_  = int (myFetchDbl(params, "angles_Z_quality", 0));
00248 
00249     oaInfo.type_ = myFetchString(params, "cocoa_type", 0);
00250 
00251     oaInfo.ID_ = int(myFetchDbl(params, "cmssw_ID", 0));
00252 
00253     if(ALIUtils::debug >= 4) {
00254       std::cout << "CocoaAnalyzer::ReadXML OBJECT " << oaInfo.name_ << " pos/angles read " << std::endl;
00255     }
00256 
00257     if( fabs( oaInfo.angx_.value_ - angles[0] ) > 1.E-9 || 
00258         fabs( oaInfo.angy_.value_ - angles[1] ) > 1.E-9 || 
00259         fabs( oaInfo.angz_.value_ - angles[2] ) > 1.E-9 ) {
00260       std::cerr << " WRONG ANGLE IN OBJECT " << oaInfo.name_<< 
00261         oaInfo.angx_.value_ << " =? " << angles[0] <<
00262         oaInfo.angy_.value_ << " =? " << angles[1] <<
00263         oaInfo.angz_.value_ << " =? " << angles[2] << std::endl;
00264     }
00265 
00266     //----- Read extra entries and measurements
00267     const std::vector<const DDsvalues_type *> params2(fv.specifics());
00268     std::vector<const DDsvalues_type *>::const_iterator spit = params2.begin();
00269     std::vector<const DDsvalues_type *>::const_iterator endspit = params2.end();
00270     //--- extra entries variables
00271     std::vector<std::string> names, dims;
00272     std::vector<double> values, errors, quality;
00273     //--- measurements variables
00274     std::vector<std::string> measNames;
00275     std::vector<std::string> measTypes;
00276     std::map<std::string, std::vector<std::string> > measObjectNames;
00277     std::map<std::string, std::vector<std::string> > measParamNames;
00278     std::map<std::string, std::vector<double> > measParamValues;
00279     std::map<std::string, std::vector<double> > measParamSigmas;
00280     std::map<std::string, std::vector<double> > measIsSimulatedValue;
00281 
00282     for ( ; spit != endspit; ++spit ) {
00283       DDsvalues_type::const_iterator sit = (**spit).begin();
00284       DDsvalues_type::const_iterator endsit = (**spit).end();
00285       for ( ; sit != endsit; ++sit ) {
00286         if (sit->second.name() == "extra_entry") {
00287           names = sit->second.strings();
00288         } else if (sit->second.name() == "dimType") {
00289           dims = sit->second.strings();
00290         } else if (sit->second.name() == "value") {
00291           values = sit->second.doubles();
00292         } else if (sit->second.name() == "sigma") {
00293           errors = sit->second.doubles();
00294         } else if (sit->second.name() == "quality") {
00295           quality = sit->second.doubles();
00296 
00297         } else if (sit->second.name() == "meas_name") {
00298           //-     std::cout << " meas_name found " << std::endl;
00299           measNames = sit->second.strings();
00300         } else if (sit->second.name() == "meas_type") {
00301           //- std::cout << " meas_type found " << std::endl;
00302           measTypes = sit->second.strings();
00303         }
00304         
00305       }
00306     }
00307 
00308     //---- loop again to look for the measurement object names, that have the meas name in the SpecPar title 
00309     //    <Parameter name="meas_object_name_SENSOR2D:OCMS/sens2" value="OCMS/laser1"  eval="false" /> 
00310     //   <Parameter name="meas_object_name_SENSOR2D:OCMS/sens2" value="OCMS/sens2"  eval="false" /> 
00311 
00312     std::vector<std::string>::iterator vsite;
00313     for ( spit = params2.begin(); spit != params2.end(); ++spit ) {
00314       //-  std::cout << "loop vector DDsvalues " << std::endl;
00315       DDsvalues_type::const_iterator sit = (**spit).begin();
00316       DDsvalues_type::const_iterator endsit = (**spit).end();
00317       for ( ; sit != endsit; ++sit ) {
00318         for( vsite = measNames.begin(); vsite != measNames.end(); vsite++ ){
00319           //- std::cout << "looping measObjectNames " << *vsite << std::endl;
00320           if (sit->second.name() == "meas_object_name_"+(*vsite)) {
00321             measObjectNames[*vsite] = sit->second.strings();
00322           }else if (sit->second.name() == "meas_value_name_"+(*vsite)) {
00323             measParamNames[*vsite] = sit->second.strings();
00324           }else if (sit->second.name() == "meas_value_"+(*vsite)) {
00325             measParamValues[*vsite] = sit->second.doubles();
00326           }else if (sit->second.name() == "meas_sigma_"+(*vsite)) {
00327             measParamSigmas[*vsite] = sit->second.doubles();
00328           }else if (sit->second.name() == "meas_is_simulated_value_"+(*vsite)) {
00329             measIsSimulatedValue[*vsite] = sit->second.doubles(); // this is not in OptAlignParam info
00330             if(ALIUtils::debug >= 5) {
00331               std::cout << *vsite << " setting issimu " << measIsSimulatedValue[*vsite][0] << std::endl;
00332             }
00333           }
00334           if(ALIUtils::debug >= 5) {
00335             std::cout << "CocoaAnalyser: looped measObjectNames " << "meas_object_name_"+(*vsite) << " n obj " << measObjectNames[*vsite].size() << std::endl;
00336           }
00337 
00338         }
00339         
00340       }
00341     }
00342     
00343     if(ALIUtils::debug >= 4) {
00344       std::cout << " CocoaAnalyzer::ReadXML:  Fill extra entries with read parameters " << std::endl;
00345     }
00346     //--- Fill extra entries with read parameters    
00347     if ( names.size() == dims.size() && dims.size() == values.size() 
00348          && values.size() == errors.size() && errors.size() == quality.size() ) {
00349       for ( size_t ind = 0; ind < names.size(); ++ind ) {
00350         double dimFactor = 1.;
00351         std::string type = oaParam.dimType();
00352         if( type == "centre" || type == "length" ) {
00353           dimFactor = 0.001; // in XML it is in mm 
00354         }else if ( type == "angles" || type == "angle" || type == "nodim" ){
00355           dimFactor = 1.;
00356         } 
00357         oaParam.value_ = values[ind]*dimFactor;
00358         oaParam.error_ = errors[ind]*dimFactor;
00359         oaParam.quality_ = int (quality[ind]);
00360         oaParam.name_ = names[ind];
00361         oaParam.dim_type_ = dims[ind];
00362         oaInfo.extraEntries_.push_back (oaParam);
00363         oaParam.clear();
00364       }
00365 
00366       //t      std::cout << names.size() << " OBJECT " << oaInfo.name_ << " extra entries read " << oaInfo << std::endl;
00367 
00368       oaList_.opticalAlignments_.push_back(oaInfo);
00369     } else {
00370       std::cout << "WARNING FOR NOW: sizes of extra parameters (names, dimType, value, quality) do"
00371                 << " not match!  Did not add " << nObjects << " item to OpticalAlignments." 
00372                 << std::endl;
00373     }
00374 
00375     if(ALIUtils::debug >= 4) {
00376       std::cout << " CocoaAnalyzer::ReadXML:  Fill measurements with read parameters " << std::endl;
00377     }
00378     //--- Fill measurements with read parameters    
00379     if ( measNames.size() == measTypes.size() ) {
00380       for ( size_t ind = 0; ind < measNames.size(); ++ind ) {
00381         oaMeas.ID_ = ind;
00382         oaMeas.name_ = measNames[ind];
00383         oaMeas.type_ = measTypes[ind];
00384         oaMeas.measObjectNames_ = measObjectNames[oaMeas.name_];
00385         if( measParamNames.size() == measParamValues.size() && measParamValues.size() == measParamSigmas.size() ) { 
00386           for( size_t ind2 = 0; ind2 < measParamNames[oaMeas.name_].size(); ind2++ ){
00387             oaParam.name_ = measParamNames[oaMeas.name_][ind2];
00388             oaParam.value_ = measParamValues[oaMeas.name_][ind2];
00389             oaParam.error_ = measParamSigmas[oaMeas.name_][ind2];
00390             if( oaMeas.type_ == "SENSOR2D" || oaMeas.type_ == "COPS" || oaMeas.type_ == "DISTANCEMETER" || oaMeas.type_ == "DISTANCEMETER!DIM" || oaMeas.type_ == "DISTANCEMETER3DIM" ) {
00391               oaParam.dim_type_ = "length";
00392             } else if( oaMeas.type_ == "TILTMETER" ) {
00393               oaParam.dim_type_ = "angle";
00394             } else {
00395               std::cerr << "CocoaAnalyzer::ReadXMLFile. Invalid measurement type: " <<  oaMeas.type_ << std::endl;
00396               std::exception();
00397             }
00398             
00399             oaMeas.values_.push_back( oaParam );
00400             oaMeas.isSimulatedValue_.push_back( measIsSimulatedValue[oaMeas.name_][ind2] );
00401             if(ALIUtils::debug >= 5) {
00402               std::cout << oaMeas.name_ << " copying issimu " << oaMeas.isSimulatedValue_[oaMeas.isSimulatedValue_.size()-1]  << " = " << measIsSimulatedValue[oaMeas.name_][ind2] << std::endl;
00403             //-           std::cout << ind2 << " adding meas value " << oaParam << std::endl;
00404             }
00405             oaParam.clear();    
00406           }
00407         } else {
00408           if(ALIUtils::debug >= 2) {
00409             std::cout << "WARNING FOR NOW: sizes of measurement parameters (name, value, sigma) do"
00410                       << " not match! for measurement " << oaMeas.name_ << " !Did not fill parameters for this measurement " << std::endl;
00411           }
00412         }
00413         measList_.oaMeasurements_.push_back (oaMeas);
00414         if(ALIUtils::debug >= 5) {
00415           std::cout << "CocoaAnalyser: MEASUREMENT " << oaMeas.name_ << " extra entries read " << oaMeas << std::endl;
00416         }
00417         oaMeas.clear();
00418       }
00419       
00420     } else {
00421       if(ALIUtils::debug >= 2) {
00422         std::cout << "WARNING FOR NOW: sizes of measurements (names, types do"
00423                   << " not match!  Did not add " << nObjects << " item to XXXMeasurements" 
00424                   << std::endl;
00425       }
00426     }
00427 
00428 //       std::cout << "sizes are values=" << values.size();
00429 //       std::cout << "  sigma(errors)=" << errors.size();
00430 //       std::cout << "  quality=" << quality.size();
00431 //       std::cout << "  names=" << names.size();
00432 //       std::cout << "  dimType=" << dims.size() << std::endl;
00433     oaInfo.clear();
00434     doCOCOA = fv.next(); // go to next part
00435   } // while (doCOCOA)
00436   if(ALIUtils::debug >= 3) {
00437     std::cout << "CocoaAnalyzer::ReadXML: Finished building " << nObjects+1 << " OpticalAlignInfo objects" << " and " <<  measList_.oaMeasurements_.size() << " OpticalAlignMeasurementInfo objects " << std::endl;
00438   }
00439   if(ALIUtils::debug >= 5) {
00440     std::cout << " @@@@@@ OpticalAlignments " << oaList_ << std::endl;
00441     std::cout << " @@@@@@ OpticalMeasurements " << measList_ << std::endl;
00442   }
00443 
00444 }
00445 
00446 //------------------------------------------------------------------------
00447 std::vector<OpticalAlignInfo> CocoaAnalyzer::ReadCalibrationDB( const edm::EventSetup& evts ) 
00448 {
00449   if(ALIUtils::debug >= 3) {
00450     std::cout<< std::endl <<"$$$ CocoaAnalyzer::ReadCalibrationDB: " << std::endl;
00451   }
00452   
00453   using namespace edm::eventsetup;
00454   edm::ESHandle<OpticalAlignments> pObjs;
00455   evts.get<OpticalAlignmentsRcd>().get(pObjs);
00456   const OpticalAlignments* dbObj = pObjs.product();
00457   
00458   if(ALIUtils::debug >= 5) {
00459     std::vector<OpticalAlignInfo>::const_iterator it;
00460     for( it=dbObj->opticalAlignments_.begin();it!=dbObj->opticalAlignments_.end(); ++it ){
00461       std::cout<<"CocoaAnalyzer::ReadCalibrationDB:  OpticalAlignInfo READ "<< *it << std::endl;
00462     }
00463   }
00464   
00465   if(ALIUtils::debug >= 4) {
00466     std::cout<<"CocoaAnalyzer::ReadCalibrationDB:  Number of OpticalAlignInfo READ "<< dbObj->opticalAlignments_.size() << std::endl;
00467   }
00468 
00469   return dbObj->opticalAlignments_;
00470 }
00471 
00472 
00473 //------------------------------------------------------------------------
00474 void CocoaAnalyzer::CorrectOptAlignments( std::vector<OpticalAlignInfo>& oaListCalib )
00475 { 
00476   if(ALIUtils::debug >= 3) {
00477     std::cout<< std::endl<< "$$$ CocoaAnalyzer::CorrectOptAlignments: " << std::endl;
00478   }
00479 
00480   std::vector<OpticalAlignInfo>::const_iterator it;
00481   for( it=oaListCalib.begin();it!=oaListCalib.end(); ++it ){
00482     OpticalAlignInfo oaInfoDB = *it;
00483     OpticalAlignInfo* oaInfoXML = FindOpticalAlignInfoXML( oaInfoDB );
00484     std::cerr << "error " << (*it).name_ << std::endl;
00485     if( oaInfoXML == 0 ) {
00486       if(ALIUtils::debug >= 2) {
00487         std::cerr << "@@@@@ WARNING CocoaAnalyzer::CorrectOptAlignments:  OpticalAlignInfo read from DB is not present in XML "<< *it << std::endl;
00488       }
00489     } else {
00490       //------ Correct info       
00491       if(ALIUtils::debug >= 5) {
00492         std::cout << "CocoaAnalyzer::CorrectOptAlignments: correcting data from DB info " << std::endl;
00493       }
00494       CorrectOaParam( &oaInfoXML->x_, oaInfoDB.x_ );
00495       CorrectOaParam( &oaInfoXML->y_, oaInfoDB.y_ );
00496       CorrectOaParam( &oaInfoXML->z_, oaInfoDB.z_ );
00497       CorrectOaParam( &oaInfoXML->angx_, oaInfoDB.angx_ );
00498       CorrectOaParam( &oaInfoXML->angy_, oaInfoDB.angy_ );
00499       CorrectOaParam( &oaInfoXML->angz_, oaInfoDB.angz_ );
00500       std::vector<OpticalAlignParam>::iterator itoap1, itoap2;
00501       std::vector<OpticalAlignParam> extraEntDB = oaInfoDB.extraEntries_;
00502       std::vector<OpticalAlignParam>* extraEntXML = &(oaInfoXML->extraEntries_);
00503       for( itoap1 = extraEntDB.begin(); itoap1 != extraEntDB.end(); itoap1++ ){
00504         bool pFound = false;
00505         //----- Look for the extra parameter in XML oaInfo that has the same name
00506         std::string oaName = (*itoap1).name_.substr( 1, (*itoap1).name_.size()-2 );
00507         for( itoap2 = extraEntXML->begin(); itoap2 != extraEntXML->end(); itoap2++ ){
00508           if( oaName == (*itoap2).name_ ) {
00509             CorrectOaParam( &(*itoap2), *itoap1 ); 
00510             pFound = true;
00511             break;
00512           }
00513         }
00514         if( !pFound && oaName != "None" ) {
00515           if(ALIUtils::debug >= 2) {
00516             std::cerr << "@@@@@ WARNING CocoaAnalyzer::CorrectOptAlignments:  extra entry read from DB is not present in XML "<< *itoap1 << " in object " << *it << std::endl;
00517           }
00518         }
00519 
00520       }
00521       if(ALIUtils::debug >= 5) {
00522         std::cout << "CocoaAnalyzer::CorrectOptAlignments: corrected OpticalAlingInfo " << oaList_ << std::endl;
00523       }
00524 
00525     }
00526   }
00527 
00528 }
00529 
00530 
00531 //------------------------------------------------------------------------
00532 OpticalAlignInfo* CocoaAnalyzer::FindOpticalAlignInfoXML( OpticalAlignInfo oaInfo )
00533 {
00534   OpticalAlignInfo* oaInfoXML = 0;
00535   std::vector<OpticalAlignInfo>::iterator it;
00536   for( it=oaList_.opticalAlignments_.begin();it!=oaList_.opticalAlignments_.end(); ++it ){
00537     std::string oaName = oaInfo.name_.substr( 1, oaInfo.name_.size()-2 );
00538 
00539     if(ALIUtils::debug >= 5) {
00540       std::cout << "CocoaAnalyzer::FindOpticalAlignInfoXML:  looking for OAI " << (*it).name_ << " =? " << oaName << std::endl;
00541     }
00542     if( (*it).name_ == oaName ) {
00543       oaInfoXML = &(*it);
00544       if(ALIUtils::debug >= 4) {
00545         std::cout << "CocoaAnalyzer::FindOpticalAlignInfoXML:  OAI found " << oaInfoXML->name_ << std::endl;
00546       }
00547       break;
00548     }
00549   }
00550 
00551   return oaInfoXML;
00552 }
00553 
00554 
00555 //------------------------------------------------------------------------
00556 bool CocoaAnalyzer::CorrectOaParam( OpticalAlignParam* oaParamXML, OpticalAlignParam oaParamDB )
00557 {
00558   if(ALIUtils::debug >= 4) {
00559     std::cout << "CocoaAnalyzer::CorrectOaParam  old value= " << oaParamXML->value_  << " new value= " << oaParamDB.value_ << std::endl;
00560   }
00561   if( oaParamDB.value_ == -9.999E9 ) return false;
00562  
00563   double dimFactor = 1.; 
00564   //loop for an Entry with equal type to entries to know which is the 
00565   std::string type = oaParamDB.dimType();
00566   if( type == "centre" || type == "length" ) {
00567     dimFactor = 0.01; // in DB it is in cm 
00568   }else if ( type == "angles" || type == "angle" || type == "nodim" ){
00569     dimFactor = 1.;
00570   }else {
00571     std::cerr << "!!! COCOA programming error: inform responsible: incorrect OpticalAlignParam type = " << type << std::endl;
00572     std::exception();
00573   }
00574 
00575   oaParamXML->value_ = oaParamDB.value_*dimFactor;
00576 
00577   return true;
00578 
00579 }
00580 
00581 
00582 //-#include "Alignment/CocoaUtilities/interface/GlobalOptionMgr.h"
00583 //-#include "Alignment/CocoaUtilities/interface/ALIFileIn.h"
00584 
00585 //-----------------------------------------------------------------------
00586 void CocoaAnalyzer::analyze(const edm::Event& evt, const edm::EventSetup& evts)
00587 {
00588   ALIUtils::setDebugVerbosity( 5 );
00589 
00590   ReadXMLFile( evts );
00591 
00592   std::vector<OpticalAlignInfo> oaListCalib = ReadCalibrationDB( evts );
00593 
00594   CorrectOptAlignments( oaListCalib );
00595 
00596   new CocoaDaqReaderRoot( theCocoaDaqRootFileName );
00597  
00598   /*-
00599   int nEvents = daqReader->GetNEvents();
00600   for( int ii = 0; ii < nEvents; ii++) {
00601     if( ! daqReader->ReadEvent( ii ) ) break;
00602   }
00603   */ 
00604 
00605   RunCocoa();
00606 
00607   //  std::cout << "!!!! NOT  DumpCocoaResults() " << std::endl;  
00608 
00609   //  CocoaDBMgr::getInstance()->DumpCocoaResults();  
00610 
00611   return;
00612 
00613 
00614   //  using namespace edm::eventsetup;
00615   //  std::cout <<" I AM IN RUN NUMBER "<<evt.id().run() <<std::endl;
00616   //  std::cout <<" ---EVENT NUMBER "<<evt.id().event() <<std::endl;
00617  
00618   // just a quick dump of the private OpticalAlignments object
00619   //  std::cout << oaList_ << std::endl;
00620  
00621   // STEP 2:
00622   // Get calibrated OpticalAlignments.  In this case I'm using
00623   // some sqlite database examples that are generated using
00624   // testOptAlignWriter.cc
00625   // from CondFormats/OptAlignObjects/test/
00626 
00627 //   edm::ESHandle<OpticalAlignments> oaESHandle;
00628 //   context.get<OpticalAlignmentsRcd>().get(oaESHandle);
00629 
00630 //   // This assumes they ALL come in together.  This may not be
00631 //   // the "real" case.  One could envision different objects coming
00632 //   // in and we would get by label each one (type).
00633 
00634 //   std::cout << "========== eventSetup data changes with IOV =========" << std::endl;
00635 //   std::cout << *oaESHandle << std::endl;
00636 //   //============== COCOA WORK!
00637 //   //  calibrated values should be used to "correct" the ones read in during beginJob
00638 //   //==============
00639   
00640 //   // 
00641 //   // to see how to iterate over the OpticalAlignments, please
00642 //   // refer to the << operator of OpticalAlignments, OpticalAlignInfo
00643 //   // and OpticalAlignParam.
00644 //   //     const OpticalAlignments* myoa=oa.product();
00645   
00646 //   // STEP 3:
00647 //   // This retrieves the Measurements
00648 //   // for each event, a new set of measurements is available.
00649 //  edm::Handle<OpticalAlignMeasurements> measHandle;
00650 //  evt.getByLabel("OptAlignGeneratedSource", measHandle); 
00651     
00652 
00653 //  std::cout << "========== event data product changes with every event =========" << std::endl;
00654 //  std::cout << *measHandle << std::endl;
00655 
00656   //============== COCOA WORK!
00657   //  Each set of optical alignment measurements can be used
00658   //  in whatever type of analysis COCOA does. 
00659   //==============
00660 
00661 } //end of ::analyze()
00662 
00663 // STEP 4:  one could use ::endJob() to write out the OpticalAlignments
00664 // generated by the analysis. Example code of writing is in
00665 // CondFormats/Alignment/test/testOptAlignWriter.cc
00666 
00667 
00668 //-----------------------------------------------------------------------
00669 double CocoaAnalyzer::myFetchDbl(const DDsvalues_type& dvst, 
00670                                       const std::string& spName,
00671                                       const size_t& vecInd ) {
00672   DDValue val(spName, 0.0);
00673   if (DDfetch(&dvst,val)) {
00674     if ( val.doubles().size() > vecInd ) {
00675       //          std::cout << "about to return: " << val.doubles()[vecInd] << std::endl;
00676       return val.doubles()[vecInd];
00677     } else {
00678       std::cout << "WARNING: OUT OF BOUNDS RETURNING 0 for index " << vecInd << " of SpecPar " << spName << std::endl;
00679     }
00680   }
00681   return 0.0;
00682 }
00683 
00684 //-----------------------------------------------------------------------
00685 std::string CocoaAnalyzer::myFetchString(const DDsvalues_type& dvst, 
00686                                       const std::string& spName,
00687                                       const size_t& vecInd ) {
00688   DDValue val(spName, 0.0);
00689   if (DDfetch(&dvst,val)) {
00690     if ( val.strings().size() > vecInd ) {
00691       //          std::cout << "about to return: " << val.doubles()[vecInd] << std::endl;
00692       return val.strings()[vecInd];
00693     } else {
00694       std::cout << "WARNING: OUT OF BOUNDS RETURNING 0 for index " << vecInd << " of SpecPar " << spName << std::endl;
00695     }
00696   }
00697   return "";
00698 }
00699