File indexing completed on 2024-04-06 11:55:58
0001 #include "FWCore/ServiceRegistry/interface/Service.h"
0002 #include "CondFormats/DataRecord/interface/OpticalAlignmentsRcd.h"
0003 #include "CondFormats/OptAlignObjects/interface/OpticalAlignMeasurementInfo.h"
0004 #include <DD4hep/DD4hepUnits.h>
0005 #include "DetectorDescription/DDCMS/interface/DDCompactView.h"
0006 #include "DetectorDescription/DDCMS/interface/DDFilteredView.h"
0007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0008 #include "FWCore/Framework/interface/ESHandle.h"
0009 #include "FWCore/Framework/interface/ESTransientHandle.h"
0010 #include "FWCore/Framework/interface/Event.h"
0011 #include "FWCore/Framework/interface/EventSetup.h"
0012 #include "FWCore/Framework/interface/MakerMacros.h"
0013 #include "Geometry/Records/interface/IdealGeometryRecord.h"
0014
0015 #include "Alignment/CocoaUtilities/interface/ALIUtils.h"
0016 #include "Alignment/CocoaModel/interface/Model.h"
0017 #include "Alignment/CocoaFit/interface/Fit.h"
0018 #include "Alignment/CocoaModel/interface/Entry.h"
0019 #include "Alignment/CocoaUtilities/interface/ALIFileOut.h"
0020 #include "Alignment/CocoaModel/interface/CocoaDaqReaderRoot.h"
0021 #include "Alignment/CocoaModel/interface/OpticalObject.h"
0022 #include "Alignment/CocoaUtilities/interface/GlobalOptionMgr.h"
0023 #include "Alignment/CocoaFit/interface/CocoaDBMgr.h"
0024 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0025 #include "DetectorDescription/DDCMS/interface/DDSpecParRegistry.h"
0026 #include "CondFormats/OptAlignObjects/interface/OpticalAlignments.h"
0027 #include "CondFormats/OptAlignObjects/interface/OpticalAlignMeasurements.h"
0028
0029 class CocoaAnalyzer : public edm::one::EDAnalyzer<edm::one::SharedResources> {
0030 public:
0031 explicit CocoaAnalyzer(edm::ParameterSet const& p);
0032 explicit CocoaAnalyzer(int i) {}
0033 ~CocoaAnalyzer() override {}
0034
0035 void beginJob() override;
0036 void analyze(const edm::Event& e, const edm::EventSetup& c) override;
0037
0038 private:
0039 void readXMLFile(const edm::EventSetup& evts);
0040
0041 std::vector<OpticalAlignInfo> readCalibrationDB(const edm::EventSetup& evts);
0042 void correctAllOpticalAlignments(std::vector<OpticalAlignInfo>& allDBOpticalAlignments);
0043 void correctOpticalAlignmentParameter(OpticalAlignParam& myXMLParam, const OpticalAlignParam& myDBParam);
0044
0045 void runCocoa();
0046
0047 private:
0048 edm::ESGetToken<OpticalAlignments, OpticalAlignmentsRcd> optAliToken;
0049 edm::ESGetToken<cms::DDCompactView, IdealGeometryRecord> ddcvToken;
0050 OpticalAlignments oaList_;
0051 OpticalAlignMeasurements measList_;
0052 std::string theCocoaDaqRootFileName_;
0053 };
0054
0055 CocoaAnalyzer::CocoaAnalyzer(edm::ParameterSet const& pset) : optAliToken(esConsumes()), ddcvToken(esConsumes()) {
0056 theCocoaDaqRootFileName_ = pset.getParameter<std::string>("cocoaDaqRootFile");
0057 int maxEvents = pset.getParameter<int32_t>("maxEvents");
0058 GlobalOptionMgr::getInstance()->setDefaultGlobalOptions();
0059 GlobalOptionMgr::getInstance()->setGlobalOption("maxEvents", maxEvents);
0060 GlobalOptionMgr::getInstance()->setGlobalOption("writeDBAlign", 1);
0061 GlobalOptionMgr::getInstance()->setGlobalOption("writeDBOptAlign", 1);
0062 usesResource("CocoaAnalyzer");
0063 }
0064
0065 void CocoaAnalyzer::beginJob() {}
0066
0067 void CocoaAnalyzer::analyze(const edm::Event& evt, const edm::EventSetup& evts) {
0068 ALIUtils::setDebugVerbosity(5);
0069
0070
0071 readXMLFile(evts);
0072
0073
0074 std::vector<OpticalAlignInfo> oaListCalib = readCalibrationDB(evts);
0075 correctAllOpticalAlignments(oaListCalib);
0076
0077
0078 runCocoa();
0079 }
0080
0081
0082
0083
0084
0085
0086 void CocoaAnalyzer::readXMLFile(const edm::EventSetup& evts) {
0087 edm::ESTransientHandle<cms::DDCompactView> myCompactView = evts.getTransientHandle(ddcvToken);
0088 const cms::DDDetector* mySystem = myCompactView->detector();
0089
0090 if (mySystem) {
0091
0092 const dd4hep::Volume& worldVolume = mySystem->worldVolume();
0093
0094 if (ALIUtils::debug >= 3) {
0095 edm::LogInfo("Alignment") << "CocoaAnalyzer::ReadXML: world object = " << worldVolume.name();
0096 }
0097
0098 OpticalAlignInfo worldInfo;
0099 worldInfo.ID_ = 0;
0100 worldInfo.name_ = worldVolume.name();
0101 worldInfo.type_ = "system";
0102 worldInfo.parentName_ = "";
0103 worldInfo.x_.value_ = 0.;
0104 worldInfo.x_.error_ = 0.;
0105 worldInfo.x_.quality_ = 0;
0106 worldInfo.y_.value_ = 0.;
0107 worldInfo.y_.error_ = 0.;
0108 worldInfo.y_.quality_ = 0;
0109 worldInfo.z_.value_ = 0.;
0110 worldInfo.z_.error_ = 0.;
0111 worldInfo.z_.quality_ = 0;
0112 worldInfo.angx_.value_ = 0.;
0113 worldInfo.angx_.error_ = 0.;
0114 worldInfo.angx_.quality_ = 0;
0115 worldInfo.angy_.value_ = 0.;
0116 worldInfo.angy_.error_ = 0.;
0117 worldInfo.angy_.quality_ = 0;
0118 worldInfo.angz_.value_ = 0.;
0119 worldInfo.angz_.error_ = 0.;
0120 worldInfo.angz_.quality_ = 0;
0121 oaList_.opticalAlignments_.emplace_back(worldInfo);
0122
0123
0124
0125
0126
0127 const cms::DDSpecParRegistry& allSpecParSections = myCompactView->specpars();
0128
0129
0130
0131
0132 cms::DDFilteredView myFilteredView(mySystem, worldVolume);
0133
0134 cms::DDSpecParRefs cocoaParameterSpecParSections;
0135
0136 const std::string cocoaParameterAttribute = "COCOA";
0137 const std::string cocoaParameterValue = "COCOA";
0138
0139
0140 allSpecParSections.filter(cocoaParameterSpecParSections, cocoaParameterAttribute, cocoaParameterValue);
0141
0142
0143
0144 myFilteredView.mergedSpecifics(cocoaParameterSpecParSections);
0145
0146
0147 int nObjects = 0;
0148 bool doCOCOA = myFilteredView.firstChild();
0149
0150
0151 while (doCOCOA) {
0152 ++nObjects;
0153
0154 OpticalAlignInfo oaInfo;
0155 OpticalAlignParam oaParam;
0156 OpticalAlignMeasurementInfo oaMeas;
0157
0158
0159 const dd4hep::PlacedVolume& myPlacedVolume = myFilteredView.volume();
0160 const std::string& name = myPlacedVolume.name();
0161 const std::string& nodePath = myFilteredView.path();
0162 oaInfo.name_ = nodePath;
0163
0164
0165 oaInfo.parentName_ = nodePath.substr(0, nodePath.rfind('/', nodePath.length()));
0166
0167 if (ALIUtils::debug >= 4) {
0168 edm::LogInfo("Alignment") << " CocoaAnalyzer::ReadXML reading object " << name;
0169 edm::LogInfo("Alignment") << " @@ Name built= " << oaInfo.name_ << " short_name= " << name
0170 << " parent= " << oaInfo.parentName_;
0171 }
0172
0173
0174
0175
0176
0177 const dd4hep::Direction& transl = myPlacedVolume.position();
0178
0179 if (ALIUtils::debug >= 4) {
0180 edm::LogInfo("Alignment") << "Local translation in dd4hep units = " << transl;
0181 }
0182
0183
0184
0185 oaInfo.x_.name_ = "X";
0186 oaInfo.x_.dim_type_ = "centre";
0187 oaInfo.x_.value_ = transl.x() / dd4hep::m;
0188 oaInfo.x_.error_ = cms::getParameterValueFromSpecParSections<double>(allSpecParSections,
0189 nodePath,
0190 "centre_X_sigma",
0191 0) /
0192 dd4hep::m;
0193 oaInfo.x_.quality_ = static_cast<int>(
0194 cms::getParameterValueFromSpecParSections<double>(allSpecParSections, nodePath, "centre_X_quality", 0));
0195
0196 oaInfo.y_.name_ = "Y";
0197 oaInfo.y_.dim_type_ = "centre";
0198 oaInfo.y_.value_ = transl.y() / dd4hep::m;
0199 oaInfo.y_.error_ = cms::getParameterValueFromSpecParSections<double>(allSpecParSections,
0200 nodePath,
0201 "centre_Y_sigma",
0202 0) /
0203 dd4hep::m;
0204 oaInfo.y_.quality_ = static_cast<int>(
0205 cms::getParameterValueFromSpecParSections<double>(allSpecParSections, nodePath, "centre_Y_quality", 0));
0206
0207 oaInfo.z_.name_ = "Z";
0208 oaInfo.z_.dim_type_ = "centre";
0209 oaInfo.z_.value_ = transl.z() / dd4hep::m;
0210 oaInfo.z_.error_ = cms::getParameterValueFromSpecParSections<double>(allSpecParSections,
0211 nodePath,
0212 "centre_Z_sigma",
0213 0) /
0214 dd4hep::m;
0215 oaInfo.z_.quality_ = static_cast<int>(
0216 cms::getParameterValueFromSpecParSections<double>(allSpecParSections, nodePath, "centre_Z_quality", 0));
0217
0218
0219
0220
0221
0222
0223
0224
0225
0226
0227
0228 const TGeoHMatrix parentToChild = myPlacedVolume.matrix();
0229
0230 const TGeoHMatrix& childToParent = parentToChild.Inverse();
0231
0232
0233
0234
0235 const Double_t* rot = childToParent.GetRotationMatrix();
0236 const double xx = rot[0];
0237 const double xy = rot[1];
0238 const double xz = rot[2];
0239 const double yx = rot[3];
0240 const double yy = rot[4];
0241 const double yz = rot[5];
0242 const double zx = rot[6];
0243 const double zy = rot[7];
0244 const double zz = rot[8];
0245 if (ALIUtils::debug >= 4) {
0246 edm::LogInfo("Alignment") << "Local rotation = ";
0247 edm::LogInfo("Alignment") << xx << " " << xy << " " << xz;
0248 edm::LogInfo("Alignment") << yx << " " << yy << " " << yz;
0249 edm::LogInfo("Alignment") << zx << " " << zy << " " << zz;
0250 }
0251 const CLHEP::Hep3Vector colX(xx, yx, zx);
0252 const CLHEP::Hep3Vector colY(xy, yy, zy);
0253 const CLHEP::Hep3Vector colZ(xz, yz, zz);
0254 const CLHEP::HepRotation rotclhep(colX, colY, colZ);
0255 const std::vector<double>& angles = ALIUtils::getRotationAnglesFromMatrix(rotclhep, 0., 0., 0.);
0256
0257
0258
0259 oaInfo.angx_.name_ = "X";
0260 oaInfo.angx_.dim_type_ = "angles";
0261 oaInfo.angx_.value_ =
0262 cms::getParameterValueFromSpecParSections<double>(allSpecParSections, nodePath, "angles_X_value", 0);
0263 oaInfo.angx_.error_ =
0264 cms::getParameterValueFromSpecParSections<double>(allSpecParSections, nodePath, "angles_X_sigma", 0);
0265 oaInfo.angx_.quality_ = static_cast<int>(
0266 cms::getParameterValueFromSpecParSections<double>(allSpecParSections, nodePath, "angles_X_quality", 0));
0267
0268 oaInfo.angy_.name_ = "Y";
0269 oaInfo.angy_.dim_type_ = "angles";
0270 oaInfo.angy_.value_ =
0271 cms::getParameterValueFromSpecParSections<double>(allSpecParSections, nodePath, "angles_Y_value", 0);
0272 oaInfo.angy_.error_ =
0273 cms::getParameterValueFromSpecParSections<double>(allSpecParSections, nodePath, "angles_Y_sigma", 0);
0274 oaInfo.angy_.quality_ = static_cast<int>(
0275 cms::getParameterValueFromSpecParSections<double>(allSpecParSections, nodePath, "angles_Y_quality", 0));
0276
0277 oaInfo.angz_.name_ = "Z";
0278 oaInfo.angz_.dim_type_ = "angles";
0279 oaInfo.angz_.value_ =
0280 cms::getParameterValueFromSpecParSections<double>(allSpecParSections, nodePath, "angles_Z_value", 0);
0281 oaInfo.angz_.error_ =
0282 cms::getParameterValueFromSpecParSections<double>(allSpecParSections, nodePath, "angles_Z_sigma", 0);
0283 oaInfo.angz_.quality_ = static_cast<int>(
0284 cms::getParameterValueFromSpecParSections<double>(allSpecParSections, nodePath, "angles_Z_quality", 0));
0285
0286 oaInfo.type_ =
0287 cms::getParameterValueFromSpecParSections<std::string>(allSpecParSections, nodePath, "cocoa_type", 0);
0288
0289 oaInfo.ID_ = static_cast<int>(
0290 cms::getParameterValueFromSpecParSections<double>(allSpecParSections, nodePath, "cmssw_ID", 0));
0291
0292 if (ALIUtils::debug >= 4) {
0293 edm::LogInfo("Alignment") << "CocoaAnalyzer::ReadXML OBJECT " << oaInfo.name_ << " pos/angles read ";
0294 }
0295
0296
0297
0298 if (fabs(oaInfo.angx_.value_ - angles[0]) > 1.E-9 || fabs(oaInfo.angy_.value_ - angles[1]) > 1.E-9 ||
0299 fabs(oaInfo.angz_.value_ - angles[2]) > 1.E-9) {
0300 edm::LogError("Alignment") << " WRONG ANGLE IN OBJECT " << oaInfo.name_ << oaInfo.angx_.value_ << " =? "
0301 << angles[0] << oaInfo.angy_.value_ << " =? " << angles[1] << oaInfo.angz_.value_
0302 << " =? " << angles[2];
0303 }
0304
0305
0306
0307 const std::vector<std::string>& names =
0308 cms::getAllParameterValuesFromSpecParSections<std::string>(allSpecParSections, nodePath, "extra_entry");
0309 const std::vector<std::string>& dims =
0310 cms::getAllParameterValuesFromSpecParSections<std::string>(allSpecParSections, nodePath, "dimType");
0311 const std::vector<double>& values =
0312 cms::getAllParameterValuesFromSpecParSections<double>(allSpecParSections, nodePath, "value");
0313 const std::vector<double>& errors =
0314 cms::getAllParameterValuesFromSpecParSections<double>(allSpecParSections, nodePath, "sigma");
0315 const std::vector<double>& quality =
0316 cms::getAllParameterValuesFromSpecParSections<double>(allSpecParSections, nodePath, "quality");
0317
0318 if (ALIUtils::debug >= 4) {
0319 edm::LogInfo("Alignment") << " CocoaAnalyzer::ReadXML: Fill extra entries with read parameters ";
0320 }
0321
0322 if (names.size() == dims.size() && dims.size() == values.size() && values.size() == errors.size() &&
0323 errors.size() == quality.size()) {
0324 for (size_t i = 0; i < names.size(); ++i) {
0325 double dimFactor = 1.;
0326 const std::string& type = dims[i];
0327 if (type == "centre" || type == "length") {
0328 dimFactor =
0329 1. / dd4hep::m;
0330 } else if (type == "angles" || type == "angle" || type == "nodim") {
0331 dimFactor = 1.;
0332 }
0333 oaParam.value_ = values[i] * dimFactor;
0334 oaParam.error_ = errors[i] * dimFactor;
0335 oaParam.quality_ = static_cast<int>(quality[i]);
0336 oaParam.name_ = names[i];
0337 oaParam.dim_type_ = dims[i];
0338 oaInfo.extraEntries_.emplace_back(oaParam);
0339 oaParam.clear();
0340 }
0341
0342 oaList_.opticalAlignments_.emplace_back(oaInfo);
0343 } else {
0344 edm::LogInfo("Alignment") << "WARNING FOR NOW: sizes of extra parameters (names, dimType, value, quality) do"
0345 << " not match! Did not add " << nObjects << " item to OpticalAlignments.";
0346 }
0347
0348
0349 const std::vector<std::string>& measNames =
0350 cms::getAllParameterValuesFromSpecParSections<std::string>(allSpecParSections, nodePath, "meas_name");
0351 const std::vector<std::string>& measTypes =
0352 cms::getAllParameterValuesFromSpecParSections<std::string>(allSpecParSections, nodePath, "meas_type");
0353
0354 std::map<std::string, std::vector<std::string>> measObjectNames;
0355 std::map<std::string, std::vector<std::string>> measParamNames;
0356 std::map<std::string, std::vector<double>> measParamValues;
0357 std::map<std::string, std::vector<double>> measParamSigmas;
0358 std::map<std::string, std::vector<double>> measIsSimulatedValue;
0359 for (const auto& name : measNames) {
0360 measObjectNames[name] = cms::getAllParameterValuesFromSpecParSections<std::string>(
0361 allSpecParSections, nodePath, "meas_object_name_" + name);
0362 measParamNames[name] = cms::getAllParameterValuesFromSpecParSections<std::string>(
0363 allSpecParSections, nodePath, "meas_value_name_" + name);
0364 measParamValues[name] =
0365 cms::getAllParameterValuesFromSpecParSections<double>(allSpecParSections, nodePath, "meas_value_" + name);
0366 measParamSigmas[name] =
0367 cms::getAllParameterValuesFromSpecParSections<double>(allSpecParSections, nodePath, "meas_sigma_" + name);
0368 measIsSimulatedValue[name] = cms::getAllParameterValuesFromSpecParSections<double>(
0369 allSpecParSections, nodePath, "meas_is_simulated_value_" + name);
0370 }
0371
0372 if (ALIUtils::debug >= 4) {
0373 edm::LogInfo("Alignment") << " CocoaAnalyzer::ReadXML: Fill measurements with read parameters ";
0374 }
0375
0376 if (measNames.size() == measTypes.size()) {
0377 for (size_t i = 0; i < measNames.size(); ++i) {
0378 oaMeas.ID_ = i;
0379 oaMeas.name_ = measNames[i];
0380 oaMeas.type_ = measTypes[i];
0381 oaMeas.measObjectNames_ = measObjectNames[oaMeas.name_];
0382 if (measParamNames.size() == measParamValues.size() && measParamValues.size() == measParamSigmas.size()) {
0383 for (size_t i2 = 0; i2 < measParamNames[oaMeas.name_].size(); i2++) {
0384 oaParam.name_ = measParamNames[oaMeas.name_][i2];
0385 oaParam.value_ = measParamValues[oaMeas.name_][i2];
0386 oaParam.error_ = measParamSigmas[oaMeas.name_][i2];
0387 oaParam.quality_ = 2;
0388 if (oaMeas.type_ == "SENSOR2D" || oaMeas.type_ == "COPS" || oaMeas.type_ == "DISTANCEMETER" ||
0389 oaMeas.type_ == "DISTANCEMETER!DIM" || oaMeas.type_ == "DISTANCEMETER3DIM") {
0390 oaParam.dim_type_ = "length";
0391 } else if (oaMeas.type_ == "TILTMETER") {
0392 oaParam.dim_type_ = "angle";
0393 } else {
0394 edm::LogError("Alignment") << "CocoaAnalyzer::readXMLFile. Invalid measurement type: " << oaMeas.type_;
0395 }
0396
0397 oaMeas.values_.emplace_back(oaParam);
0398 oaMeas.isSimulatedValue_.emplace_back(measIsSimulatedValue[oaMeas.name_][i2]);
0399 if (ALIUtils::debug >= 5) {
0400 edm::LogInfo("Alignment") << oaMeas.name_ << " copying issimu "
0401 << oaMeas.isSimulatedValue_[oaMeas.isSimulatedValue_.size() - 1] << " = "
0402 << measIsSimulatedValue[oaMeas.name_][i2];
0403 }
0404 oaParam.clear();
0405 }
0406 } else {
0407 if (ALIUtils::debug >= 2) {
0408 edm::LogWarning("Alignment") << "WARNING FOR NOW: sizes of measurement parameters (name, value, sigma) do"
0409 << " not match! for measurement " << oaMeas.name_
0410 << " !Did not fill parameters for this measurement ";
0411 }
0412 }
0413 measList_.oaMeasurements_.emplace_back(oaMeas);
0414 if (ALIUtils::debug >= 5) {
0415 edm::LogInfo("Alignment") << "CocoaAnalyser: MEASUREMENT " << oaMeas.name_ << " extra entries read "
0416 << oaMeas;
0417 }
0418 oaMeas.clear();
0419 }
0420
0421 } else {
0422 if (ALIUtils::debug >= 2) {
0423 edm::LogWarning("Alignment") << "WARNING FOR NOW: sizes of measurements (names, types do"
0424 << " not match! Did not add " << nObjects << " item to XXXMeasurements";
0425 }
0426 }
0427
0428 oaInfo.clear();
0429 doCOCOA = myFilteredView.firstChild();
0430 }
0431
0432 if (ALIUtils::debug >= 3) {
0433 edm::LogInfo("Alignment") << "CocoaAnalyzer::ReadXML: Finished building " << nObjects + 1
0434 << " OpticalAlignInfo objects"
0435 << " and " << measList_.oaMeasurements_.size()
0436 << " OpticalAlignMeasurementInfo objects ";
0437 }
0438 if (ALIUtils::debug >= 5) {
0439 edm::LogInfo("Alignment") << " @@@@@@ OpticalAlignments " << oaList_;
0440 edm::LogInfo("Alignment") << " @@@@@@ OpticalMeasurements " << measList_;
0441 }
0442 }
0443 }
0444
0445
0446
0447
0448
0449 std::vector<OpticalAlignInfo> CocoaAnalyzer::readCalibrationDB(const edm::EventSetup& evts) {
0450 if (ALIUtils::debug >= 3) {
0451 edm::LogInfo("Alignment") << "$$$ CocoaAnalyzer::readCalibrationDB: ";
0452 }
0453
0454 using namespace edm::eventsetup;
0455 const OpticalAlignments* pObjs = &evts.getData(optAliToken);
0456 const std::vector<OpticalAlignInfo>& infoFromDB = pObjs->opticalAlignments_;
0457
0458 if (ALIUtils::debug >= 5) {
0459 edm::LogInfo("Alignment") << "CocoaAnalyzer::readCalibrationDB: Number of OpticalAlignInfo READ "
0460 << infoFromDB.size();
0461 for (const auto& myInfoFromDB : infoFromDB) {
0462 edm::LogInfo("Alignment") << "CocoaAnalyzer::readCalibrationDB: OpticalAlignInfo READ " << myInfoFromDB;
0463 }
0464 }
0465
0466 return infoFromDB;
0467 }
0468
0469
0470
0471
0472 void CocoaAnalyzer::correctAllOpticalAlignments(std::vector<OpticalAlignInfo>& allDBOpticalAlignments) {
0473 if (ALIUtils::debug >= 3) {
0474 edm::LogInfo("Alignment") << "$$$ CocoaAnalyzer::correctAllOpticalAlignments: ";
0475 }
0476
0477 for (const auto& myDBInfo : allDBOpticalAlignments) {
0478 if (ALIUtils::debug >= 5) {
0479 edm::LogInfo("Alignment") << "CocoaAnalyzer::findOpticalAlignInfoXML: Looking for OAI " << myDBInfo.name_;
0480 }
0481
0482 std::vector<OpticalAlignInfo>& allXMLOpticalAlignments = oaList_.opticalAlignments_;
0483 const auto& myXMLInfo = std::find_if(
0484 allXMLOpticalAlignments.begin(), allXMLOpticalAlignments.end(), [&](const auto& myXMLInfoCandidate) {
0485 return myXMLInfoCandidate.name_ == myDBInfo.name_;
0486 });
0487
0488 if (myXMLInfo != allXMLOpticalAlignments.end()) {
0489 if (ALIUtils::debug >= 4) {
0490 edm::LogInfo("Alignment") << "CocoaAnalyzer::findOpticalAlignInfoXML: OAI found " << myXMLInfo->name_;
0491 edm::LogInfo("Alignment")
0492 << "CocoaAnalyzer::correctAllOpticalAlignments: correcting data from XML with DB info.";
0493 }
0494 correctOpticalAlignmentParameter(myXMLInfo->x_, myDBInfo.x_);
0495 correctOpticalAlignmentParameter(myXMLInfo->y_, myDBInfo.y_);
0496 correctOpticalAlignmentParameter(myXMLInfo->z_, myDBInfo.z_);
0497 correctOpticalAlignmentParameter(myXMLInfo->angx_, myDBInfo.angx_);
0498 correctOpticalAlignmentParameter(myXMLInfo->angy_, myDBInfo.angy_);
0499 correctOpticalAlignmentParameter(myXMLInfo->angz_, myDBInfo.angz_);
0500
0501
0502 const std::vector<OpticalAlignParam>& allDBExtraEntries = myDBInfo.extraEntries_;
0503 std::vector<OpticalAlignParam>& allXMLExtraEntries = myXMLInfo->extraEntries_;
0504 for (const auto& myDBExtraEntry : allDBExtraEntries) {
0505 const auto& myXMLExtraEntry = std::find_if(
0506 allXMLExtraEntries.begin(), allXMLExtraEntries.end(), [&](const auto& myXMLExtraEntryCandidate) {
0507 return myXMLExtraEntryCandidate.name_ == myDBExtraEntry.name_;
0508 });
0509
0510 if (myXMLExtraEntry != allXMLExtraEntries.end()) {
0511 correctOpticalAlignmentParameter(*myXMLExtraEntry, myDBExtraEntry);
0512 } else {
0513 if (myDBExtraEntry.name_ != "None") {
0514 if (ALIUtils::debug >= 2) {
0515 edm::LogError("Alignment")
0516 << "CocoaAnalyzer::correctAllOpticalAlignments: extra entry read from DB is not present in XML "
0517 << myDBExtraEntry << " in object " << myDBInfo;
0518 }
0519 }
0520 }
0521 }
0522
0523 if (ALIUtils::debug >= 5) {
0524 edm::LogInfo("Alignment") << "CocoaAnalyzer::correctAllOpticalAlignments: corrected OpticalAlingInfo "
0525 << oaList_;
0526 }
0527 } else {
0528 if (ALIUtils::debug >= 2) {
0529 edm::LogError("Alignment") << "CocoaAnalyzer::correctAllOpticalAlignments: OpticalAlignInfo read from DB "
0530 << myDBInfo << " is not present in XML.";
0531 }
0532 }
0533 }
0534 }
0535
0536
0537
0538
0539 void CocoaAnalyzer::correctOpticalAlignmentParameter(OpticalAlignParam& myXMLParam,
0540 const OpticalAlignParam& myDBParam) {
0541 if (myDBParam.value_ != -9.999E9) {
0542 const std::string& type = myDBParam.dimType();
0543 double dimFactor = 1.;
0544
0545 if (type == "centre" || type == "length") {
0546 dimFactor = 1. / dd4hep::m;
0547 } else if (type == "angles" || type == "angle" || type == "nodim") {
0548 dimFactor = 1.;
0549 } else {
0550 edm::LogError("Alignment") << "Incorrect OpticalAlignParam type = " << type;
0551 }
0552
0553 const double correctedValue = myDBParam.value_ * dimFactor;
0554 if (ALIUtils::debug >= 4) {
0555 edm::LogInfo("Alignment") << "CocoaAnalyzer::correctOpticalAlignmentParameter old value= " << myXMLParam.value_
0556 << " new value= " << correctedValue;
0557 }
0558 myXMLParam.value_ = correctedValue;
0559 }
0560 }
0561
0562
0563
0564
0565 void CocoaAnalyzer::runCocoa() {
0566 if (ALIUtils::debug >= 3) {
0567 edm::LogInfo("Alignment") << "$$$ CocoaAnalyzer::runCocoa: ";
0568 }
0569
0570
0571 Model& model = Model::getInstance();
0572 model.BuildSystemDescriptionFromOA(oaList_);
0573
0574 if (ALIUtils::debug >= 3) {
0575 edm::LogInfo("Alignment") << "$$ CocoaAnalyzer::runCocoa: geometry built ";
0576 }
0577
0578
0579 model.BuildMeasurementsFromOA(measList_);
0580
0581 if (ALIUtils::debug >= 3) {
0582 edm::LogInfo("Alignment") << "$$ CocoaAnalyzer::runCocoa: measurements built ";
0583 }
0584
0585
0586 Fit::getInstance();
0587 Fit::startFit();
0588
0589 if (ALIUtils::debug >= 0)
0590 edm::LogInfo("Alignment") << "............ program ended OK";
0591 if (ALIUtils::report >= 1) {
0592 ALIFileOut& fileout = ALIFileOut::getInstance(Model::ReportFName());
0593 fileout << "............ program ended OK";
0594 }
0595 }
0596
0597 DEFINE_FWK_MODULE(CocoaAnalyzer);