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