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