00001
00011 #include "Alignment/MillePedeAlignmentAlgorithm/interface/MillePedeAlignmentAlgorithm.h"
00012
00013
00014 #include "FWCore/Framework/interface/ESHandle.h"
00015 #include "FWCore/Framework/interface/EventSetup.h"
00016 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00017
00018 #include "TrackingTools/PatternTools/interface/Trajectory.h"
00019
00020
00021
00022 #include "Alignment/MillePedeAlignmentAlgorithm/interface/MillePedeMonitor.h"
00023 #include "Alignment/MillePedeAlignmentAlgorithm/interface/MillePedeVariables.h"
00024 #include "Alignment/MillePedeAlignmentAlgorithm/interface/MillePedeVariablesIORoot.h"
00025 #include "Mille.h"
00026 #include "PedeSteerer.h"
00027 #include "PedeReader.h"
00028 #include "Alignment/MillePedeAlignmentAlgorithm/interface/PedeLabelerBase.h"
00029 #include "Alignment/MillePedeAlignmentAlgorithm/interface/PedeLabelerPluginFactory.h"
00030
00031 #include "Alignment/ReferenceTrajectories/interface/TrajectoryFactoryBase.h"
00032 #include "Alignment/ReferenceTrajectories/interface/TrajectoryFactoryPlugin.h"
00033
00034 #include "Alignment/CommonAlignmentAlgorithm/interface/AlignmentParameterStore.h"
00035 #include "Alignment/CommonAlignmentAlgorithm/interface/AlignmentIORoot.h"
00036 #include "Alignment/CommonAlignmentAlgorithm/interface/IntegratedCalibrationBase.h"
00037
00038 #include "Alignment/CommonAlignment/interface/AlignmentParameters.h"
00039 #include "Alignment/CommonAlignment/interface/AlignableNavigator.h"
00040 #include "Alignment/CommonAlignment/interface/AlignableDetOrUnitPtr.h"
00041
00042
00043 #include "Alignment/TrackerAlignment/interface/AlignableTracker.h"
00044 #include "Alignment/MuonAlignment/interface/AlignableMuon.h"
00045 #include "Alignment/CommonAlignment/interface/AlignableExtras.h"
00046
00047 #include "DataFormats/CLHEP/interface/AlgebraicObjects.h"
00048 #include "DataFormats/TrackReco/interface/Track.h"
00049 #include "DataFormats/SiStripDetId/interface/SiStripDetId.h"
00050 #include "DataFormats/Alignment/interface/TkFittedLasBeam.h"
00051 #include "Alignment/LaserAlignment/interface/TsosVectorCollection.h"
00052
00053 #include <Geometry/CommonDetUnit/interface/GeomDetUnit.h>
00054 #include <Geometry/CommonDetUnit/interface/GeomDetType.h>
00055 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
00056 #include "Geometry/Records/interface/IdealGeometryRecord.h"
00057
00058 #include "DataFormats/TrackerRecHit2D/interface/ProjectedSiStripRecHit2D.h"
00059
00060 #include <fstream>
00061 #include <sstream>
00062 #include <algorithm>
00063
00064 #include <TMath.h>
00065 #include <TMatrixDSymEigen.h>
00066 typedef TransientTrackingRecHit::ConstRecHitContainer ConstRecHitContainer;
00067 typedef TransientTrackingRecHit::ConstRecHitPointer ConstRecHitPointer;
00068 typedef TrajectoryFactoryBase::ReferenceTrajectoryCollection RefTrajColl;
00069
00070
00071 #include <iostream>
00072 #include "Alignment/SurveyAnalysis/interface/SurveyPxbImage.h"
00073 #include "Alignment/SurveyAnalysis/interface/SurveyPxbImageLocalFit.h"
00074 #include "Alignment/SurveyAnalysis/interface/SurveyPxbImageReader.h"
00075 #include "Alignment/SurveyAnalysis/interface/SurveyPxbDicer.h"
00076 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
00077 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00078
00079 #include "Alignment/CommonAlignmentParametrization/interface/AlignmentParametersFactory.h"
00080
00081
00082
00083 MillePedeAlignmentAlgorithm::MillePedeAlignmentAlgorithm(const edm::ParameterSet &cfg) :
00084 AlignmentAlgorithmBase(cfg),
00085 theConfig(cfg), theMode(this->decodeMode(theConfig.getUntrackedParameter<std::string>("mode"))),
00086 theDir(theConfig.getUntrackedParameter<std::string>("fileDir")),
00087 theAlignmentParameterStore(0), theAlignables(), theAlignableNavigator(0),
00088 theMonitor(0), theMille(0), thePedeLabels(0), thePedeSteer(0),
00089 theTrajectoryFactory(0),
00090 theMinNumHits(cfg.getParameter<unsigned int>("minNumHits")),
00091 theMaximalCor2D(cfg.getParameter<double>("max2Dcorrelation")),
00092 theLastWrittenIov(0)
00093 {
00094 if (!theDir.empty() && theDir.find_last_of('/') != theDir.size()-1) theDir += '/';
00095 edm::LogInfo("Alignment") << "@SUB=MillePedeAlignmentAlgorithm" << "Start in mode '"
00096 << theConfig.getUntrackedParameter<std::string>("mode")
00097 << "' with output directory '" << theDir << "'.";
00098 if (this->isMode(myMilleBit)) {
00099 theMille = new Mille((theDir + theConfig.getParameter<std::string>("binaryFile")).c_str());
00100 }
00101 }
00102
00103
00104
00105 MillePedeAlignmentAlgorithm::~MillePedeAlignmentAlgorithm()
00106 {
00107 delete theAlignableNavigator;
00108 theAlignableNavigator = 0;
00109 delete theMille;
00110 theMille = 0;
00111 delete theMonitor;
00112 theMonitor = 0;
00113 delete thePedeSteer;
00114 thePedeSteer = 0;
00115 delete thePedeLabels;
00116 thePedeLabels = 0;
00117 delete theTrajectoryFactory;
00118 theTrajectoryFactory = 0;
00119 }
00120
00121
00122
00123 void MillePedeAlignmentAlgorithm::initialize(const edm::EventSetup &setup,
00124 AlignableTracker *tracker, AlignableMuon *muon, AlignableExtras *extras,
00125 AlignmentParameterStore *store)
00126 {
00127 if (muon) {
00128 edm::LogWarning("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::initialize"
00129 << "Running with AlignabeMuon not yet tested.";
00130 }
00131
00132
00133 edm::ESHandle<TrackerTopology> tTopoHandle;
00134 setup.get<IdealGeometryRecord>().get(tTopoHandle);
00135 const TrackerTopology* const tTopo = tTopoHandle.product();
00136
00137 theAlignableNavigator = new AlignableNavigator(extras, tracker, muon);
00138 theAlignmentParameterStore = store;
00139 theAlignables = theAlignmentParameterStore->alignables();
00140
00141 edm::ParameterSet pedeLabelerCfg(theConfig.getParameter<edm::ParameterSet>("pedeLabeler"));
00142 edm::VParameterSet RunRangeSelectionVPSet(theConfig.getUntrackedParameter<edm::VParameterSet>("RunRangeSelection"));
00143 pedeLabelerCfg.addUntrackedParameter<edm::VParameterSet>("RunRangeSelection",
00144 RunRangeSelectionVPSet);
00145
00146 std::string labelerPlugin = "PedeLabeler";
00147 if (RunRangeSelectionVPSet.size()>0) {
00148 labelerPlugin = "RunRangeDependentPedeLabeler";
00149 if (pedeLabelerCfg.exists("plugin")) {
00150 std::string labelerPluginCfg = pedeLabelerCfg.getParameter<std::string>("plugin");
00151 if ((labelerPluginCfg!="PedeLabeler" && labelerPluginCfg!="RunRangeDependentPedeLabeler") ||
00152 pedeLabelerCfg.getUntrackedParameter<edm::VParameterSet>("parameterInstances").size()>0) {
00153 throw cms::Exception("BadConfig")
00154 << "MillePedeAlignmentAlgorithm::initialize"
00155 << "both RunRangeSelection and generic labeler specified in config file. "
00156 << "Please get rid of either one of them.\n";
00157 }
00158 }
00159 } else {
00160 if (pedeLabelerCfg.exists("plugin")) {
00161 labelerPlugin = pedeLabelerCfg.getParameter<std::string>("plugin");
00162 }
00163 }
00164
00165 if (!pedeLabelerCfg.exists("plugin")) {
00166 pedeLabelerCfg.addUntrackedParameter<std::string>("plugin", labelerPlugin);
00167 }
00168
00169 edm::LogInfo("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::initialize"
00170 << "Using plugin '" << labelerPlugin << "' to generate labels.";
00171
00172 thePedeLabels = PedeLabelerPluginFactory::get()->create(labelerPlugin,
00173 PedeLabelerBase::TopLevelAlignables(tracker, muon, extras),
00174 pedeLabelerCfg);
00175
00176
00177 edm::ParameterSet pedeSteerCfg(theConfig.getParameter<edm::ParameterSet>("pedeSteerer"));
00178 thePedeSteer = new PedeSteerer(tracker, muon, extras,
00179 theAlignmentParameterStore, thePedeLabels,
00180 pedeSteerCfg, theDir, !this->isMode(myPedeSteerBit));
00181
00182
00183
00184 const std::vector<edm::ParameterSet> mprespset
00185 (theConfig.getParameter<std::vector<edm::ParameterSet> >("pedeReaderInputs"));
00186 if (!mprespset.empty()) {
00187 edm::LogInfo("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::initialize"
00188 << "Apply " << mprespset.end() - mprespset.begin()
00189 << " previous MillePede constants from 'pedeReaderInputs'.";
00190 }
00191
00192
00193
00194 RunRange runrange(cond::timeTypeSpecs[cond::runnumber].beginValue,
00195 cond::timeTypeSpecs[cond::runnumber].endValue);
00196 for (std::vector<edm::ParameterSet>::const_iterator iSet = mprespset.begin(), iE = mprespset.end();
00197 iSet != iE; ++iSet) {
00198
00199
00200 if (!this->readFromPede((*iSet), false, runrange)) {
00201 throw cms::Exception("BadConfig")
00202 << "MillePedeAlignmentAlgorithm::initialize: Problems reading input constants of "
00203 << "pedeReaderInputs entry " << iSet - mprespset.begin() << '.';
00204 }
00205 theAlignmentParameterStore->applyParameters();
00206
00207 theAlignmentParameterStore->resetParameters();
00208 }
00209
00210
00211 thePedeSteer->buildSubSteer(tracker, muon, extras);
00212
00213
00214 this->buildUserVariables(theAlignables);
00215
00216 if (this->isMode(myMilleBit)) {
00217 if (!theConfig.getParameter<std::vector<std::string> >("mergeBinaryFiles").empty() ||
00218 !theConfig.getParameter<std::vector<std::string> >("mergeTreeFiles").empty()) {
00219 throw cms::Exception("BadConfig")
00220 << "'vstring mergeTreeFiles' and 'vstring mergeBinaryFiles' must be empty for "
00221 << "modes running mille.";
00222 }
00223 const std::string moniFile(theConfig.getUntrackedParameter<std::string>("monitorFile"));
00224 if (moniFile.size()) theMonitor = new MillePedeMonitor(tTopo, (theDir + moniFile).c_str());
00225
00226
00227 const edm::ParameterSet fctCfg(theConfig.getParameter<edm::ParameterSet>("TrajectoryFactory"));
00228 const std::string fctName(fctCfg.getParameter<std::string>("TrajectoryFactoryName"));
00229 theTrajectoryFactory = TrajectoryFactoryPlugin::get()->create(fctName, fctCfg);
00230 }
00231
00232 if (this->isMode(myPedeSteerBit)) {
00233
00234 const edm::ParameterSet pxbSurveyCfg(theConfig.getParameter<edm::ParameterSet>("surveyPixelBarrel"));
00235 theDoSurveyPixelBarrel = pxbSurveyCfg.getParameter<bool>("doSurvey");
00236 if (theDoSurveyPixelBarrel) this->addPxbSurvey(pxbSurveyCfg);
00237 }
00238 }
00239
00240
00241 bool MillePedeAlignmentAlgorithm::addCalibrations(const std::vector<IntegratedCalibrationBase*> &iCals)
00242 {
00243 theCalibrations.insert(theCalibrations.end(), iCals.begin(), iCals.end());
00244 thePedeLabels->addCalibrations(iCals);
00245 return true;
00246 }
00247
00248
00249 bool MillePedeAlignmentAlgorithm::setParametersForRunRange(const RunRange &runrange)
00250 {
00251 if (this->isMode(myPedeReadBit)) {
00252
00253 theAlignmentParameterStore->restoreCachedTransformations();
00254
00255
00256 theAlignmentParameterStore->resetParameters();
00257
00258 this->buildUserVariables(theAlignables);
00259
00260 if (!this->readFromPede(theConfig.getParameter<edm::ParameterSet>("pedeReader"), true, runrange)) {
00261 edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::setParametersForRunRange"
00262 << "Problems reading pede result, but applying!";
00263 }
00264 theAlignmentParameterStore->applyParameters();
00265
00266 this->doIO(++theLastWrittenIov);
00267 }
00268
00269 return true;
00270 }
00271
00272
00273
00274 void MillePedeAlignmentAlgorithm::terminate(const edm::EventSetup& iSetup)
00275 {
00276 delete theMille;
00277 theMille = 0;
00278
00279 std::vector<std::string> files;
00280 if (this->isMode(myMilleBit) || !theConfig.getParameter<std::string>("binaryFile").empty()) {
00281 files.push_back(theDir + theConfig.getParameter<std::string>("binaryFile"));
00282 } else {
00283 const std::vector<std::string> plainFiles
00284 (theConfig.getParameter<std::vector<std::string> >("mergeBinaryFiles"));
00285 for (std::vector<std::string>::const_iterator i = plainFiles.begin(), iEnd = plainFiles.end();
00286 i != iEnd; ++i) {
00287 files.push_back(theDir + *i);
00288 }
00289 }
00290
00291
00292 theAlignmentParameterStore->cacheTransformations();
00293
00294 const std::string masterSteer(thePedeSteer->buildMasterSteer(files));
00295 if (this->isMode(myPedeRunBit)) {
00296 thePedeSteer->runPede(masterSteer);
00297 }
00298
00299
00300
00301 this->doIO(0);
00302 theLastWrittenIov = 0;
00303 }
00304
00305
00306
00307 void MillePedeAlignmentAlgorithm::run(const edm::EventSetup &setup, const EventInfo &eventInfo)
00308 {
00309 if (!this->isMode(myMilleBit)) return;
00310 const ConstTrajTrackPairCollection &tracks = eventInfo.trajTrackPairs_;
00311
00312 if (theMonitor) {
00313 for (ConstTrajTrackPairCollection::const_iterator iTrajTrack = tracks.begin();
00314 iTrajTrack != tracks.end(); ++iTrajTrack) {
00315 theMonitor->fillTrack((*iTrajTrack).second);
00316 }
00317 }
00318
00319 const RefTrajColl trajectories(theTrajectoryFactory->trajectories(setup, tracks, eventInfo.beamSpot_));
00320
00321
00322 unsigned int refTrajCount = 0;
00323 for (RefTrajColl::const_iterator iRefTraj = trajectories.begin(), iRefTrajE = trajectories.end();
00324 iRefTraj != iRefTrajE; ++iRefTraj, ++refTrajCount) {
00325
00326 RefTrajColl::value_type refTrajPtr = *iRefTraj;
00327 if (theMonitor) theMonitor->fillRefTrajectory(refTrajPtr);
00328
00329 const std::pair<unsigned int, unsigned int> nHitXy
00330 = this->addReferenceTrajectory(setup, eventInfo, refTrajPtr);
00331
00332 if (theMonitor && (nHitXy.first || nHitXy.second)) {
00333
00334
00335 const reco::Track *trackPtr =
00336 (trajectories.size() == tracks.size() ? tracks[refTrajCount].second : 0);
00337 theMonitor->fillUsedTrack(trackPtr, nHitXy.first, nHitXy.second);
00338 }
00339
00340 }
00341 }
00342
00343
00344
00345 std::pair<unsigned int, unsigned int>
00346 MillePedeAlignmentAlgorithm::addReferenceTrajectory(const edm::EventSetup &setup,
00347 const EventInfo &eventInfo,
00348 const RefTrajColl::value_type &refTrajPtr)
00349 {
00350 std::pair<unsigned int, unsigned int> hitResultXy(0,0);
00351 if (refTrajPtr->isValid()) {
00352
00353
00354 std::vector<AlignmentParameters*> parVec(refTrajPtr->recHits().size());
00355
00356 std::vector<bool> validHitVecY(refTrajPtr->recHits().size(), false);
00357
00358 for (unsigned int iHit = 0; iHit < refTrajPtr->recHits().size(); ++iHit) {
00359 const int flagXY = this->addMeasurementData(setup, eventInfo, refTrajPtr, iHit, parVec[iHit]);
00360
00361 if (flagXY < 0) {
00362 hitResultXy.first = 0;
00363 break;
00364 } else {
00365 if (flagXY >= 1) ++hitResultXy.first;
00366 validHitVecY[iHit] = (flagXY >= 2);
00367 }
00368 }
00369
00370
00371 for (unsigned int iVirtualMeas = 0; iVirtualMeas < refTrajPtr->numberOfVirtualMeas(); ++iVirtualMeas) {
00372 this->addVirtualMeas(refTrajPtr, iVirtualMeas);
00373 }
00374
00375
00376 if (hitResultXy.first == 0 || hitResultXy.first < theMinNumHits) {
00377 theMille->kill();
00378 hitResultXy.first = hitResultXy.second = 0;
00379 } else {
00380 theMille->end();
00381
00382
00383 hitResultXy.second = this->addHitCount(parVec, validHitVecY);
00384 }
00385 }
00386
00387 return hitResultXy;
00388 }
00389
00390
00391 unsigned int
00392 MillePedeAlignmentAlgorithm::addHitCount(const std::vector<AlignmentParameters*> &parVec,
00393 const std::vector<bool> &validHitVecY) const
00394 {
00395
00396 unsigned int nHitY = 0;
00397 for (unsigned int iHit = 0; iHit < validHitVecY.size(); ++iHit) {
00398 Alignable *ali = (parVec[iHit] ? parVec[iHit]->alignable() : 0);
00399
00400
00401
00402 while (ali) {
00403 AlignmentParameters *pars = ali->alignmentParameters();
00404 if (pars) {
00405
00406 MillePedeVariables *mpVar = static_cast<MillePedeVariables*>(pars->userVariables());
00407
00408 mpVar->increaseHitsX();
00409 if (validHitVecY[iHit]) {
00410 mpVar->increaseHitsY();
00411 if (pars == parVec[iHit]) ++nHitY;
00412 }
00413 }
00414 ali = ali->mother();
00415 }
00416 }
00417
00418 return nHitY;
00419 }
00420
00421
00422 void MillePedeAlignmentAlgorithm::endRun(const EventInfo &eventInfo, const EndRunInfo &runInfo,
00423 const edm::EventSetup &setup)
00424 {
00425 if(runInfo.tkLasBeams_ && runInfo.tkLasBeamTsoses_){
00426
00427 this->addLaserData(eventInfo, *(runInfo.tkLasBeams_), *(runInfo.tkLasBeamTsoses_));
00428 }
00429 }
00430
00431
00432 int MillePedeAlignmentAlgorithm::addMeasurementData(const edm::EventSetup &setup,
00433 const EventInfo &eventInfo,
00434 const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr,
00435 unsigned int iHit,
00436 AlignmentParameters *¶ms)
00437 {
00438 params = 0;
00439 theFloatBufferX.clear();
00440 theFloatBufferY.clear();
00441 theIntBuffer.clear();
00442
00443 const TrajectoryStateOnSurface &tsos = refTrajPtr->trajectoryStates()[iHit];
00444 const ConstRecHitPointer &recHitPtr = refTrajPtr->recHits()[iHit];
00445
00446 if (!recHitPtr->isValid()) return 0;
00447
00448
00449
00450 this->globalDerivativesCalibration(recHitPtr, tsos, setup, eventInfo,
00451 theFloatBufferX, theFloatBufferY, theIntBuffer);
00452
00453
00454 AlignableDetOrUnitPtr alidet(theAlignableNavigator->alignableFromDetId(recHitPtr->geographicalId()));
00455
00456 if (!this->globalDerivativesHierarchy(eventInfo,
00457 tsos, alidet, alidet, theFloatBufferX,
00458 theFloatBufferY, theIntBuffer, params)) {
00459 return -1;
00460 } else if (theFloatBufferX.empty()) {
00461 return 0;
00462 } else {
00463 return this->callMille(refTrajPtr, iHit, theIntBuffer, theFloatBufferX, theFloatBufferY);
00464 }
00465 }
00466
00467
00468 bool MillePedeAlignmentAlgorithm
00469 ::globalDerivativesHierarchy(const EventInfo &eventInfo,
00470 const TrajectoryStateOnSurface &tsos,
00471 Alignable *ali, const AlignableDetOrUnitPtr &alidet,
00472 std::vector<float> &globalDerivativesX,
00473 std::vector<float> &globalDerivativesY,
00474 std::vector<int> &globalLabels,
00475 AlignmentParameters *&lowestParams) const
00476 {
00477
00478 if (!ali) return true;
00479
00480 if (false && theMonitor && alidet != ali) theMonitor->fillFrameToFrame(alidet, ali);
00481
00482 AlignmentParameters *params = ali->alignmentParameters();
00483
00484 if (params) {
00485 if (!lowestParams) lowestParams = params;
00486
00487 bool hasSplitParameters = thePedeLabels->hasSplitParameters(ali);
00488 const unsigned int alignableLabel = thePedeLabels->alignableLabel(ali);
00489
00490 if (0 == alignableLabel) {
00491 edm::LogWarning("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::globalDerivativesHierarchy"
00492 << "Label not found, skip Alignable.";
00493 return false;
00494 }
00495
00496 const std::vector<bool> &selPars = params->selector();
00497 const AlgebraicMatrix derivs(params->derivatives(tsos, alidet));
00498
00499
00500 for (unsigned int iSel = 0; iSel < selPars.size(); ++iSel) {
00501 if (selPars[iSel]) {
00502 globalDerivativesX.push_back(derivs[iSel][kLocalX]
00503 /thePedeSteer->cmsToPedeFactor(iSel));
00504 if (hasSplitParameters==true) {
00505 globalLabels.push_back(thePedeLabels->parameterLabel(ali, iSel, eventInfo, tsos));
00506 } else {
00507 globalLabels.push_back(thePedeLabels->parameterLabel(alignableLabel, iSel));
00508 }
00509 globalDerivativesY.push_back(derivs[iSel][kLocalY]
00510 /thePedeSteer->cmsToPedeFactor(iSel));
00511 }
00512 }
00513
00514 if (thePedeSteer->isNoHiera(ali)) return true;
00515 }
00516
00517 return this->globalDerivativesHierarchy(eventInfo,
00518 tsos, ali->mother(), alidet,
00519 globalDerivativesX, globalDerivativesY,
00520 globalLabels, lowestParams);
00521 }
00522
00523
00524
00525 void MillePedeAlignmentAlgorithm::
00526 globalDerivativesCalibration(const TransientTrackingRecHit::ConstRecHitPointer &recHit,
00527 const TrajectoryStateOnSurface &tsos,
00528 const edm::EventSetup &setup, const EventInfo &eventInfo,
00529 std::vector<float> &globalDerivativesX,
00530 std::vector<float> &globalDerivativesY,
00531 std::vector<int> &globalLabels) const
00532 {
00533 std::vector<IntegratedCalibrationBase::ValuesIndexPair> derivs;
00534 for (auto iCalib = theCalibrations.begin(); iCalib != theCalibrations.end(); ++iCalib) {
00535
00536 (*iCalib)->derivatives(derivs, *recHit, tsos, setup, eventInfo);
00537 for (auto iValuesInd = derivs.begin(); iValuesInd != derivs.end(); ++iValuesInd) {
00538
00539 globalLabels.push_back(thePedeLabels->calibrationLabel(*iCalib, iValuesInd->second));
00540 globalDerivativesX.push_back(iValuesInd->first.first);
00541 globalDerivativesY.push_back(iValuesInd->first.second);
00542 }
00543 }
00544 }
00545
00546
00547
00548
00549
00550
00551
00552
00553
00554
00555
00556
00557
00558
00559
00560
00561
00562
00563
00564
00565
00566
00567
00568
00569
00570
00571
00572
00573
00574
00575
00576
00577
00578
00579
00580
00581 bool MillePedeAlignmentAlgorithm::is2D(const ConstRecHitPointer &recHit) const
00582 {
00583
00584
00585 if (recHit->dimension() < 2) {
00586 return false;
00587 } else if (recHit->detUnit()) {
00588 return recHit->detUnit()->type().isTrackerPixel();
00589 } else {
00590 if (dynamic_cast<const ProjectedSiStripRecHit2D*>(recHit->hit())) {
00591
00592 return false;
00593 } else {
00594 return true;
00595 }
00596 }
00597 }
00598
00599
00600 bool MillePedeAlignmentAlgorithm::readFromPede(const edm::ParameterSet &mprespset, bool setUserVars,
00601 const RunRange &runrange)
00602 {
00603 bool allEmpty = this->areEmptyParams(theAlignables);
00604
00605 PedeReader reader(mprespset, *thePedeSteer, *thePedeLabels, runrange);
00606 std::vector<Alignable*> alis;
00607 bool okRead = reader.read(alis, setUserVars);
00608 bool numMatch = true;
00609
00610 std::stringstream out;
00611 out << "Read " << alis.size() << " alignables";
00612 if (alis.size() != theAlignables.size()) {
00613 out << " while " << theAlignables.size() << " in store";
00614 numMatch = false;
00615 }
00616 if (!okRead) out << ", but problems in reading";
00617 if (!allEmpty) out << ", possibly overwriting previous settings";
00618 out << ".";
00619
00620 if (okRead && allEmpty) {
00621 if (numMatch) {
00622 edm::LogInfo("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::readFromPede" << out.str();
00623 } else if (alis.size()) {
00624 edm::LogWarning("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::readFromPede" << out.str();
00625 } else {
00626 edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::readFromPede" << out.str();
00627 return false;
00628 }
00629 return true;
00630 }
00631
00632 edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::readFromPede" << out.str();
00633 return false;
00634 }
00635
00636
00637 bool MillePedeAlignmentAlgorithm::areEmptyParams(const std::vector<Alignable*> &alignables) const
00638 {
00639
00640 for (std::vector<Alignable*>::const_iterator iAli = alignables.begin();
00641 iAli != alignables.end(); ++iAli) {
00642 const AlignmentParameters *params = (*iAli)->alignmentParameters();
00643 if (params) {
00644 const AlgebraicVector &parVec(params->parameters());
00645 const AlgebraicMatrix &parCov(params->covariance());
00646 for (int i = 0; i < parVec.num_row(); ++i) {
00647 if (parVec[i] != 0.) return false;
00648 for (int j = i; j < parCov.num_col(); ++j) {
00649 if (parCov[i][j] != 0.) return false;
00650 }
00651 }
00652 }
00653 }
00654
00655 return true;
00656 }
00657
00658
00659 unsigned int MillePedeAlignmentAlgorithm::doIO(int loop) const
00660 {
00661 unsigned int result = 0;
00662
00663 const std::string outFilePlain(theConfig.getParameter<std::string>("treeFile"));
00664 if (outFilePlain.empty()) {
00665 edm::LogInfo("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::doIO"
00666 << "treeFile parameter empty => skip writing for 'loop' " << loop;
00667 return result;
00668 }
00669
00670 const std::string outFile(theDir + outFilePlain);
00671
00672 AlignmentIORoot aliIO;
00673 int ioerr = 0;
00674 if (loop == 0) {
00675 aliIO.writeAlignableOriginalPositions(theAlignables, outFile.c_str(), loop, false, ioerr);
00676 if (ioerr) {
00677 edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::doIO"
00678 << "Problem " << ioerr << " in writeAlignableOriginalPositions";
00679 ++result;
00680 }
00681 } else if (loop == 1) {
00682
00683 const std::vector<std::string> inFiles
00684 (theConfig.getParameter<std::vector<std::string> >("mergeTreeFiles"));
00685 const std::vector<std::string> binFiles
00686 (theConfig.getParameter<std::vector<std::string> >("mergeBinaryFiles"));
00687 if (inFiles.size() != binFiles.size()) {
00688 edm::LogWarning("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::doIO"
00689 << "'vstring mergeTreeFiles' and 'vstring mergeBinaryFiles' "
00690 << "differ in size";
00691 }
00692 this->addHitStatistics(0, outFile, inFiles);
00693 }
00694 MillePedeVariablesIORoot millePedeIO;
00695 millePedeIO.writeMillePedeVariables(theAlignables, outFile.c_str(), loop, false, ioerr);
00696 if (ioerr) {
00697 edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::doIO"
00698 << "Problem " << ioerr << " writing MillePedeVariables";
00699 ++result;
00700 }
00701
00702 aliIO.writeOrigRigidBodyAlignmentParameters(theAlignables, outFile.c_str(), loop, false, ioerr);
00703 if (ioerr) {
00704 edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::doIO" << "Problem " << ioerr
00705 << " in writeOrigRigidBodyAlignmentParameters, " << loop;
00706 ++result;
00707 }
00708 aliIO.writeAlignableAbsolutePositions(theAlignables, outFile.c_str(), loop, false, ioerr);
00709 if (ioerr) {
00710 edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::doIO" << "Problem " << ioerr
00711 << " in writeAlignableAbsolutePositions, " << loop;
00712 ++result;
00713 }
00714
00715 return result;
00716 }
00717
00718
00719 void MillePedeAlignmentAlgorithm::buildUserVariables(const std::vector<Alignable*> &alis) const
00720 {
00721 for (std::vector<Alignable*>::const_iterator iAli = alis.begin(); iAli != alis.end(); ++iAli) {
00722 AlignmentParameters *params = (*iAli)->alignmentParameters();
00723 if (!params) {
00724 throw cms::Exception("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::buildUserVariables"
00725 << "No parameters for alignable";
00726 }
00727 MillePedeVariables *userVars = dynamic_cast<MillePedeVariables*>(params->userVariables());
00728 if (userVars) {
00729 for (unsigned int iPar = 0; iPar < userVars->size(); ++iPar) {
00730
00731
00732
00733 userVars->setAllDefault(iPar);
00734
00735 }
00736 } else {
00737 userVars = new MillePedeVariables(params->size(), thePedeLabels->alignableLabel(*iAli));
00738 params->setUserVariables(userVars);
00739 }
00740 }
00741 }
00742
00743
00744 unsigned int MillePedeAlignmentAlgorithm::decodeMode(const std::string &mode) const
00745 {
00746 if (mode == "full") {
00747 return myMilleBit + myPedeSteerBit + myPedeRunBit + myPedeReadBit;
00748 } else if (mode == "mille") {
00749 return myMilleBit;
00750 } else if (mode == "pede") {
00751 return myPedeSteerBit + myPedeRunBit + myPedeReadBit;
00752 } else if (mode == "pedeSteer") {
00753 return myPedeSteerBit;
00754 } else if (mode == "pedeRun") {
00755 return myPedeSteerBit + myPedeRunBit + myPedeReadBit;
00756 } else if (mode == "pedeRead") {
00757 return myPedeReadBit;
00758 }
00759
00760 throw cms::Exception("BadConfig")
00761 << "Unknown mode '" << mode
00762 << "', use 'full', 'mille', 'pede', 'pedeRun', 'pedeSteer' or 'pedeRead'.";
00763
00764 return 0;
00765 }
00766
00767
00768 bool MillePedeAlignmentAlgorithm::addHitStatistics(int fromIov, const std::string &outFile,
00769 const std::vector<std::string> &inFiles) const
00770 {
00771 bool allOk = true;
00772 int ierr = 0;
00773 MillePedeVariablesIORoot millePedeIO;
00774 for (std::vector<std::string>::const_iterator iFile = inFiles.begin();
00775 iFile != inFiles.end(); ++iFile) {
00776 const std::string inFile(theDir + *iFile);
00777 const std::vector<AlignmentUserVariables*> mpVars =
00778 millePedeIO.readMillePedeVariables(theAlignables, inFile.c_str(), fromIov, ierr);
00779 if (ierr || !this->addHits(theAlignables, mpVars)) {
00780 edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::addHitStatistics"
00781 << "Error " << ierr << " reading from " << inFile
00782 << ", tree " << fromIov << ", or problems in addHits";
00783 allOk = false;
00784 }
00785 for (std::vector<AlignmentUserVariables*>::const_iterator i = mpVars.begin();
00786 i != mpVars.end(); ++i) {
00787 delete *i;
00788 }
00789 }
00790
00791 return allOk;
00792 }
00793
00794
00795 bool MillePedeAlignmentAlgorithm::addHits(const std::vector<Alignable*> &alis,
00796 const std::vector<AlignmentUserVariables*> &mpVars) const
00797 {
00798 bool allOk = (mpVars.size() == alis.size());
00799 std::vector<AlignmentUserVariables*>::const_iterator iUser = mpVars.begin();
00800 for (std::vector<Alignable*>::const_iterator iAli = alis.begin();
00801 iAli != alis.end() && iUser != mpVars.end(); ++iAli, ++iUser) {
00802 MillePedeVariables *mpVarNew = dynamic_cast<MillePedeVariables*>(*iUser);
00803 AlignmentParameters *ps = (*iAli)->alignmentParameters();
00804 MillePedeVariables *mpVarOld = (ps ? dynamic_cast<MillePedeVariables*>(ps->userVariables()) : 0);
00805 if (!mpVarNew || !mpVarOld || mpVarOld->size() != mpVarNew->size()) {
00806 allOk = false;
00807 continue;
00808 }
00809
00810 mpVarOld->increaseHitsX(mpVarNew->hitsX());
00811 mpVarOld->increaseHitsY(mpVarNew->hitsY());
00812 }
00813
00814 return allOk;
00815 }
00816
00817
00818 void MillePedeAlignmentAlgorithm::makeGlobDerivMatrix(const std::vector<float> &globalDerivativesx,
00819 const std::vector<float> &globalDerivativesy,
00820 TMatrixF &aGlobalDerivativesM)
00821 {
00822
00823 for (unsigned int i = 0; i < globalDerivativesx.size(); ++i) {
00824 aGlobalDerivativesM(0,i) = globalDerivativesx[i];
00825 aGlobalDerivativesM(1,i) = globalDerivativesy[i];
00826 }
00827 }
00828
00829
00830 void MillePedeAlignmentAlgorithm::diagonalize
00831 (TMatrixDSym &aHitCovarianceM, TMatrixF &aLocalDerivativesM, TMatrixF &aHitResidualsM,
00832 TMatrixF &aGlobalDerivativesM) const
00833 {
00834 TMatrixDSymEigen myDiag(aHitCovarianceM);
00835 TMatrixD aTranfoToDiagonalSystem = myDiag.GetEigenVectors();
00836 TMatrixD aTranfoToDiagonalSystemInv = myDiag.GetEigenVectors( );
00837 TMatrixF aTranfoToDiagonalSystemInvF = myDiag.GetEigenVectors( );
00838 TMatrixD aMatrix = aTranfoToDiagonalSystemInv.Invert() * aHitCovarianceM * aTranfoToDiagonalSystem;
00839
00840
00841
00842
00843
00844 aHitCovarianceM = TMatrixDSym(2, aMatrix.GetMatrixArray());
00845 aTranfoToDiagonalSystemInvF.Invert();
00846
00847 aLocalDerivativesM = aTranfoToDiagonalSystemInvF * aLocalDerivativesM;
00848
00849
00850 aHitResidualsM = aTranfoToDiagonalSystemInvF * aHitResidualsM;
00851 aGlobalDerivativesM = aTranfoToDiagonalSystemInvF * aGlobalDerivativesM;
00852 }
00853
00854
00855 void MillePedeAlignmentAlgorithm
00856 ::addRefTrackVirtualMeas1D(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr,
00857 unsigned int iVirtualMeas, TMatrixDSym &aHitCovarianceM,
00858 TMatrixF &aHitResidualsM, TMatrixF &aLocalDerivativesM)
00859 {
00860
00861
00862 const unsigned int xIndex = iVirtualMeas + refTrajPtr->numberOfHitMeas();
00863
00864
00865
00866 aHitCovarianceM(0,0)=refTrajPtr->measurementErrors()[xIndex][xIndex];
00867
00868
00869 aHitResidualsM(0,0)= refTrajPtr->measurements()[xIndex];
00870
00871
00872 const AlgebraicMatrix &locDerivMatrix = refTrajPtr->derivatives();
00873
00874
00875
00876 for (int i = 0; i < locDerivMatrix.num_col(); ++i) {
00877 aLocalDerivativesM(0,i) = locDerivMatrix[xIndex][i];
00878 }
00879 }
00880
00881
00882 void MillePedeAlignmentAlgorithm
00883 ::addRefTrackData2D(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr,
00884 unsigned int iTrajHit, TMatrixDSym &aHitCovarianceM,
00885 TMatrixF &aHitResidualsM, TMatrixF &aLocalDerivativesM)
00886 {
00887
00888
00889 const unsigned int xIndex = iTrajHit*2;
00890 const unsigned int yIndex = iTrajHit*2+1;
00891
00892
00893
00894 aHitCovarianceM(0,0)=refTrajPtr->measurementErrors()[xIndex][xIndex];
00895 aHitCovarianceM(0,1)=refTrajPtr->measurementErrors()[xIndex][yIndex];
00896 aHitCovarianceM(1,0)=refTrajPtr->measurementErrors()[yIndex][xIndex];
00897 aHitCovarianceM(1,1)=refTrajPtr->measurementErrors()[yIndex][yIndex];
00898
00899
00900 aHitResidualsM(0,0)= refTrajPtr->measurements()[xIndex] - refTrajPtr->trajectoryPositions()[xIndex];
00901 aHitResidualsM(1,0)= refTrajPtr->measurements()[yIndex] - refTrajPtr->trajectoryPositions()[yIndex];
00902
00903
00904 const AlgebraicMatrix &locDerivMatrix = refTrajPtr->derivatives();
00905
00906
00907
00908 for (int i = 0; i < locDerivMatrix.num_col(); ++i) {
00909 aLocalDerivativesM(0,i) = locDerivMatrix[xIndex][i];
00910 aLocalDerivativesM(1,i) = locDerivMatrix[yIndex][i];
00911 }
00912 }
00913
00914
00915 int MillePedeAlignmentAlgorithm
00916 ::callMille(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr,
00917 unsigned int iTrajHit, const std::vector<int> &globalLabels,
00918 const std::vector<float> &globalDerivativesX,
00919 const std::vector<float> &globalDerivativesY)
00920 {
00921 const ConstRecHitPointer aRecHit(refTrajPtr->recHits()[iTrajHit]);
00922
00923 if((aRecHit)->dimension() == 1) {
00924 return this->callMille1D(refTrajPtr, iTrajHit, globalLabels, globalDerivativesX);
00925 } else {
00926 return this->callMille2D(refTrajPtr, iTrajHit, globalLabels,
00927 globalDerivativesX, globalDerivativesY);
00928 }
00929 }
00930
00931
00932
00933 int MillePedeAlignmentAlgorithm
00934 ::callMille1D(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr,
00935 unsigned int iTrajHit, const std::vector<int> &globalLabels,
00936 const std::vector<float> &globalDerivativesX)
00937 {
00938 const ConstRecHitPointer aRecHit(refTrajPtr->recHits()[iTrajHit]);
00939 const unsigned int xIndex = iTrajHit*2;
00940
00941
00942 const AlgebraicMatrix &locDerivMatrix = refTrajPtr->derivatives();
00943 const int nLocal = locDerivMatrix.num_col();
00944 std::vector<float> localDerivatives(nLocal);
00945 for (unsigned int i = 0; i < localDerivatives.size(); ++i) {
00946 localDerivatives[i] = locDerivMatrix[xIndex][i];
00947 }
00948
00949
00950 float residX = refTrajPtr->measurements()[xIndex] - refTrajPtr->trajectoryPositions()[xIndex];
00951 float hitErrX = TMath::Sqrt(refTrajPtr->measurementErrors()[xIndex][xIndex]);
00952
00953
00954 const int nGlobal = globalDerivativesX.size();
00955
00956
00957
00958 theMille->mille(nLocal, &(localDerivatives[0]), nGlobal, &(globalDerivativesX[0]),
00959 &(globalLabels[0]), residX, hitErrX);
00960
00961 if (theMonitor) {
00962 theMonitor->fillDerivatives(aRecHit, &(localDerivatives[0]), nLocal,
00963 &(globalDerivativesX[0]), nGlobal, &(globalLabels[0]));
00964 theMonitor->fillResiduals(aRecHit, refTrajPtr->trajectoryStates()[iTrajHit],
00965 iTrajHit, residX, hitErrX, false);
00966 }
00967
00968 return 1;
00969 }
00970
00971
00972 int MillePedeAlignmentAlgorithm
00973 ::callMille2D(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr,
00974 unsigned int iTrajHit, const std::vector<int> &globalLabels,
00975 const std::vector<float> &globalDerivativesx,
00976 const std::vector<float> &globalDerivativesy)
00977 {
00978 const ConstRecHitPointer aRecHit(refTrajPtr->recHits()[iTrajHit]);
00979
00980 if((aRecHit)->dimension() != 2) {
00981 edm::LogError("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::callMille2D"
00982 << "You try to call method for 2D hits for a "
00983 << (aRecHit)->dimension()
00984 << "D Hit. Hit gets ignored!";
00985 return -1;
00986 }
00987
00988 TMatrixDSym aHitCovarianceM(2);
00989 TMatrixF aHitResidualsM(2,1);
00990 TMatrixF aLocalDerivativesM(2, refTrajPtr->derivatives().num_col());
00991
00992 this->addRefTrackData2D(refTrajPtr, iTrajHit, aHitCovarianceM,aHitResidualsM,aLocalDerivativesM);
00993 TMatrixF aGlobalDerivativesM(2,globalDerivativesx.size());
00994 this->makeGlobDerivMatrix(globalDerivativesx, globalDerivativesy, aGlobalDerivativesM);
00995
00996
00997
00998
00999 const double corr = aHitCovarianceM(0,1) / sqrt(aHitCovarianceM(0,0) * aHitCovarianceM(1,1));
01000 if (theMonitor) theMonitor->fillCorrelations2D(corr, aRecHit);
01001 bool diag = false;
01002 switch(aRecHit->geographicalId().subdetId()) {
01003 case SiStripDetId::TID:
01004 case SiStripDetId::TEC:
01005 if (aRecHit->geographicalId().det() == DetId::Tracker && TMath::Abs(corr) > theMaximalCor2D) {
01006 this->diagonalize(aHitCovarianceM, aLocalDerivativesM, aHitResidualsM, aGlobalDerivativesM);
01007 diag = true;
01008 }
01009 break;
01010 default:;
01011 }
01012
01013 float newResidX = aHitResidualsM(0,0);
01014 float newResidY = aHitResidualsM(1,0);
01015 float newHitErrX = TMath::Sqrt(aHitCovarianceM(0,0));
01016 float newHitErrY = TMath::Sqrt(aHitCovarianceM(1,1));
01017 float *newLocalDerivsX = aLocalDerivativesM[0].GetPtr();
01018 float *newLocalDerivsY = aLocalDerivativesM[1].GetPtr();
01019 float *newGlobDerivsX = aGlobalDerivativesM[0].GetPtr();
01020 float *newGlobDerivsY = aGlobalDerivativesM[1].GetPtr();
01021 const int nLocal = aLocalDerivativesM.GetNcols();
01022 const int nGlobal = aGlobalDerivativesM.GetNcols();
01023
01024 if (diag && (newHitErrX > newHitErrY)) {
01025
01026 std::swap(newResidX, newResidY);
01027 std::swap(newHitErrX, newHitErrY);
01028 std::swap(newLocalDerivsX, newLocalDerivsY);
01029 std::swap(newGlobDerivsX, newGlobDerivsY);
01030 }
01031
01032
01033
01034 theMille->mille(nLocal, newLocalDerivsX, nGlobal, newGlobDerivsX,
01035 &(globalLabels[0]), newResidX, newHitErrX);
01036
01037 if (theMonitor) {
01038 theMonitor->fillDerivatives(aRecHit, newLocalDerivsX, nLocal, newGlobDerivsX, nGlobal,
01039 &(globalLabels[0]));
01040 theMonitor->fillResiduals(aRecHit, refTrajPtr->trajectoryStates()[iTrajHit],
01041 iTrajHit, newResidX, newHitErrX, false);
01042 }
01043 const bool isReal2DHit = this->is2D(aRecHit);
01044 if (isReal2DHit) {
01045 theMille->mille(nLocal, newLocalDerivsY, nGlobal, newGlobDerivsY,
01046 &(globalLabels[0]), newResidY, newHitErrY);
01047 if (theMonitor) {
01048 theMonitor->fillDerivatives(aRecHit, newLocalDerivsY, nLocal, newGlobDerivsY, nGlobal,
01049 &(globalLabels[0]));
01050 theMonitor->fillResiduals(aRecHit, refTrajPtr->trajectoryStates()[iTrajHit],
01051 iTrajHit, newResidY, newHitErrY, true);
01052 }
01053 }
01054
01055 return (isReal2DHit ? 2 : 1);
01056 }
01057
01058
01059 void MillePedeAlignmentAlgorithm
01060 ::addVirtualMeas(const ReferenceTrajectoryBase::ReferenceTrajectoryPtr &refTrajPtr, unsigned int iVirtualMeas)
01061 {
01062 TMatrixDSym aHitCovarianceM(1);
01063 TMatrixF aHitResidualsM(1,1);
01064 TMatrixF aLocalDerivativesM(1, refTrajPtr->derivatives().num_col());
01065
01066 this->addRefTrackVirtualMeas1D(refTrajPtr, iVirtualMeas, aHitCovarianceM, aHitResidualsM, aLocalDerivativesM);
01067
01068
01069 TMatrixF aGlobalDerivativesM(1,1);
01070 aGlobalDerivativesM(0,0) = 0;
01071
01072 float newResidX = aHitResidualsM(0,0);
01073 float newHitErrX = TMath::Sqrt(aHitCovarianceM(0,0));
01074 float *newLocalDerivsX = aLocalDerivativesM[0].GetPtr();
01075 float *newGlobDerivsX = aGlobalDerivativesM[0].GetPtr();
01076 const int nLocal = aLocalDerivativesM.GetNcols();
01077 const int nGlobal = 0;
01078
01079 theMille->mille(nLocal, newLocalDerivsX, nGlobal, newGlobDerivsX,
01080 &nGlobal, newResidX, newHitErrX);
01081 }
01082
01083
01084 void MillePedeAlignmentAlgorithm::addLaserData(const EventInfo &eventInfo,
01085 const TkFittedLasBeamCollection &lasBeams,
01086 const TsosVectorCollection &lasBeamTsoses)
01087 {
01088 TsosVectorCollection::const_iterator iTsoses = lasBeamTsoses.begin();
01089 for(TkFittedLasBeamCollection::const_iterator iBeam = lasBeams.begin(), iEnd = lasBeams.end();
01090 iBeam != iEnd; ++iBeam, ++iTsoses){
01091
01092 edm::LogInfo("Alignment") << "@SUB=MillePedeAlignmentAlgorithm::addLaserData"
01093 << "Beam " << iBeam->getBeamId() << " with "
01094 << iBeam->parameters().size() << " parameters and "
01095 << iBeam->getData().size() << " hits.\n There are "
01096 << iTsoses->size() << " TSOSes.";
01097
01098 this->addLasBeam(eventInfo, *iBeam, *iTsoses);
01099 }
01100 }
01101
01102
01103 void MillePedeAlignmentAlgorithm::addLasBeam(const EventInfo &eventInfo,
01104 const TkFittedLasBeam &lasBeam,
01105 const std::vector<TrajectoryStateOnSurface> &tsoses)
01106 {
01107 AlignmentParameters *dummyPtr = 0;
01108 std::vector<float> lasLocalDerivsX;
01109 const unsigned int beamLabel = thePedeLabels->lasBeamLabel(lasBeam.getBeamId());
01110
01111 for (unsigned int iHit = 0; iHit < tsoses.size(); ++iHit) {
01112 if (!tsoses[iHit].isValid()) continue;
01113
01114 theFloatBufferX.clear();
01115 theFloatBufferY.clear();
01116 theIntBuffer.clear();
01117 lasLocalDerivsX.clear();
01118
01119 const SiStripLaserRecHit2D &hit = lasBeam.getData()[iHit];
01120 AlignableDetOrUnitPtr lasAli(theAlignableNavigator->alignableFromDetId(hit.getDetId()));
01121 this->globalDerivativesHierarchy(eventInfo,
01122 tsoses[iHit], lasAli, lasAli,
01123 theFloatBufferX, theFloatBufferY, theIntBuffer, dummyPtr);
01124
01125 for (unsigned int nFitParams = 0;
01126 nFitParams < static_cast<unsigned int>(lasBeam.parameters().size());
01127 ++nFitParams) {
01128 const float derivative = lasBeam.derivatives()[iHit][nFitParams];
01129 if (nFitParams < lasBeam.firstFixedParameter()) {
01130 lasLocalDerivsX.push_back(derivative);
01131 } else {
01132 const unsigned int numPar = nFitParams - lasBeam.firstFixedParameter();
01133 theIntBuffer.push_back(thePedeLabels->parameterLabel(beamLabel, numPar));
01134 theFloatBufferX.push_back(derivative);
01135 }
01136 }
01137
01138 const float residual = hit.localPosition().x() - tsoses[iHit].localPosition().x();
01139
01140 const float error = 0.003;
01141
01142 theMille->mille(lasLocalDerivsX.size(), &(lasLocalDerivsX[0]), theFloatBufferX.size(),
01143 &(theFloatBufferX[0]), &(theIntBuffer[0]), residual, error);
01144 }
01145
01146 theMille->end();
01147 }
01148
01149 void MillePedeAlignmentAlgorithm::addPxbSurvey(const edm::ParameterSet &pxbSurveyCfg)
01150 {
01151
01152 const bool doOutputOnStdout(pxbSurveyCfg.getParameter<bool>("doOutputOnStdout"));
01153 if (doOutputOnStdout) std::cout << "# Output from addPxbSurvey follows below because doOutputOnStdout is set to True" << std::endl;
01154
01155
01156 SurveyPxbDicer dicer(pxbSurveyCfg.getParameter<std::vector<edm::ParameterSet> >("toySurveyParameters"), pxbSurveyCfg.getParameter<unsigned int>("toySurveySeed"));
01157 std::ofstream outfile(pxbSurveyCfg.getUntrackedParameter<std::string>("toySurveyFile").c_str());
01158
01159
01160 std::vector<SurveyPxbImageLocalFit> measurements;
01161 std::string filename(pxbSurveyCfg.getParameter<edm::FileInPath>("infile").fullPath());
01162 SurveyPxbImageReader<SurveyPxbImageLocalFit> reader(filename, measurements, 800);
01163
01164
01165 for(std::vector<SurveyPxbImageLocalFit>::size_type i=0; i!=measurements.size(); i++)
01166 {
01167 if (doOutputOnStdout) std::cout << "Module " << i << ": ";
01168
01169
01170 AlignableDetOrUnitPtr mod1(theAlignableNavigator->alignableFromDetId(measurements[i].getIdFirst()));
01171 AlignableDetOrUnitPtr mod2(theAlignableNavigator->alignableFromDetId(measurements[i].getIdSecond()));
01172 const AlignableSurface& surf1 = mod1->surface();
01173 const AlignableSurface& surf2 = mod2->surface();
01174
01175
01176 const LocalPoint fidpoint0(-0.91,+3.30);
01177 const LocalPoint fidpoint1(+0.91,+3.30);
01178 const LocalPoint fidpoint2(+0.91,-3.30);
01179 const LocalPoint fidpoint3(-0.91,-3.30);
01180
01181
01182
01183
01184 const GlobalPoint surf2point0(surf2.toGlobal(fidpoint0));
01185 const GlobalPoint surf2point1(surf2.toGlobal(fidpoint1));
01186 const LocalPoint fidpoint0inSurf1frame(surf1.toLocal(surf2point0));
01187 const LocalPoint fidpoint1inSurf1frame(surf1.toLocal(surf2point1));
01188
01189
01190 SurveyPxbImageLocalFit::fidpoint_t fidpointvec;
01191 fidpointvec.push_back(fidpoint0inSurf1frame);
01192 fidpointvec.push_back(fidpoint1inSurf1frame);
01193 fidpointvec.push_back(fidpoint2);
01194 fidpointvec.push_back(fidpoint3);
01195
01196
01197 if (pxbSurveyCfg.getParameter<bool>("doToySurvey"))
01198 {
01199 dicer.doDice(fidpointvec,measurements[i].getIdPair(), outfile);
01200 }
01201
01202
01203 measurements[i].doFit(fidpointvec, thePedeLabels->alignableLabel(mod1), thePedeLabels->alignableLabel(mod2));
01204 SurveyPxbImageLocalFit::localpars_t a;
01205 a = measurements[i].getLocalParameters();
01206 const SurveyPxbImageLocalFit::value_t chi2 = measurements[i].getChi2();
01207
01208
01209 if (doOutputOnStdout)
01210 {
01211 std::cout << "a: " << a[0] << ", " << a[1] << ", " << a[2] << ", " << a[3]
01212 << " S= " << sqrt(a[2]*a[2]+a[3]*a[3])
01213 << " phi= " << atan(a[3]/a[2])
01214 << " chi2= " << chi2 << std::endl;
01215 }
01216 if (theMonitor)
01217 {
01218 theMonitor->fillPxbSurveyHistsChi2(chi2);
01219 theMonitor->fillPxbSurveyHistsLocalPars(a[0],a[1],sqrt(a[2]*a[2]+a[3]*a[3]),atan(a[3]/a[2]));
01220 }
01221
01222
01223 for(SurveyPxbImageLocalFit::count_t j=0; j!=SurveyPxbImageLocalFit::nMsrmts; j++)
01224 {
01225 theMille->mille((int)measurements[i].getLocalDerivsSize(),
01226 measurements[i].getLocalDerivsPtr(j),
01227 (int)measurements[i].getGlobalDerivsSize(),
01228 measurements[i].getGlobalDerivsPtr(j),
01229 measurements[i].getGlobalDerivsLabelPtr(j),
01230 measurements[i].getResiduum(j),
01231 measurements[i].getSigma(j));
01232 }
01233 theMille->end();
01234 }
01235 outfile.close();
01236 }
01237