CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_8_patch3/src/Alignment/LaserAlignment/plugins/TkLasBeamFitter.cc

Go to the documentation of this file.
00001 
00016 // system include files
00017 #include <memory>
00018 
00019 // framework include files
00020 #include "FWCore/Framework/interface/Frameworkfwd.h"
00021 #include "FWCore/Framework/interface/EDProducer.h"
00022 #include "FWCore/Framework/interface/ESHandle.h"
00023 #include "FWCore/Framework/interface/Event.h"
00024 #include "FWCore/Framework/interface/EventSetup.h"
00025 #include "FWCore/Framework/interface/MakerMacros.h"
00026 #include "FWCore/Framework/interface/Run.h"
00027 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00028 #include "FWCore/Utilities/interface/InputTag.h"
00029 #include "FWCore/ServiceRegistry/interface/Service.h"
00030 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00031 
00032 #include <Geometry/CommonDetUnit/interface/GeomDetUnit.h>
00033 #include <Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h>
00034 #include <Geometry/Records/interface/TrackerDigiGeometryRecord.h>
00035 #include <MagneticField/Engine/interface/MagneticField.h>
00036 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00037 
00038 // data formats
00039 // for edm::InRun
00040 #include "DataFormats/Provenance/interface/BranchType.h"
00041 
00042 // laser data formats
00043 #include "DataFormats/Alignment/interface/TkLasBeam.h"
00044 #include "DataFormats/Alignment/interface/TkFittedLasBeam.h"
00045 
00046 // further includes
00047 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00048 #include "Alignment/LaserAlignment/interface/TsosVectorCollection.h"
00049 #include "TrackingTools/TrajectoryState/interface/FreeTrajectoryState.h"
00050 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00051 
00052 #include <iostream>
00053 #include "TMinuit.h"
00054 #include "TGraphErrors.h"
00055 #include "TF1.h"
00056 #include "TH1.h"
00057 #include "TH2.h"
00058 
00059 using namespace edm;
00060 using namespace std;
00061 
00062 
00063 //
00064 // class declaration
00065 //
00066 
00067 class TkLasBeamFitter : public edm::EDProducer {
00068 public:
00069   explicit TkLasBeamFitter(const edm::ParameterSet &config);
00070   ~TkLasBeamFitter();
00071   
00072   //virtual void beginJob(const edm::EventSetup& /*access deprecated*/) {}
00073   virtual void produce(edm::Event &event, const edm::EventSetup &setup);
00074   // virtual void beginRun(edm::Run &run, const edm::EventSetup &setup);
00075   virtual void endRun(edm::Run &run, const edm::EventSetup &setup);
00076   //virtual void endJob() {}
00077 
00078 private:
00081   void getLasBeams(TkFittedLasBeam &beam,vector<TrajectoryStateOnSurface> &tsosLas);
00082   void getLasHits(TkFittedLasBeam &beam, const SiStripLaserRecHit2D &hit, 
00083                   vector<const GeomDetUnit*> &gd, vector<GlobalPoint> &globHit,
00084                   unsigned int &hitsAtTecPlus); 
00085 //   void fillVectors(TkFittedLasBeam &beam);
00086 
00087   // need static functions to be used in fitter(..);
00088   // all parameters used therein have to be static, as well (see below)
00089   static double tecPlusFunction(double *x, double *par);
00090   static double tecMinusFunction(double *x, double *par);
00091   static double atFunction(double *x, double *par);
00092 
00093   void fitter(TkFittedLasBeam &beam, AlgebraicSymMatrix &covMatrix, 
00094               unsigned int &hitsAtTecPlus, unsigned int &nFitParams,
00095               std::vector<double> &hitPhi, std::vector<double> &hitPhiError, std::vector<double> &hitZprimeError,
00096               double &zMin, double &zMax, double &bsAngleParam,
00097               double &offset, double &offsetError, double &slope, double &slopeError,
00098               double &phiAtMinusParam, double &phiAtPlusParam,
00099               double &atThetaSplitParam);
00100 
00101   void trackPhi(TkFittedLasBeam &beam, unsigned int &hit,
00102                 double &trackPhi, double &trackPhiRef,
00103                 double &offset, double &slope, double &bsAngleParam, 
00104                 double &phiAtMinusParam, double &phiAtPlusParam,
00105                 double &atThetaSplitParam, std::vector<GlobalPoint> &globHit);
00106 
00107   void globalTrackPoint(TkFittedLasBeam &beam, 
00108                         unsigned int &hit, unsigned int &hitsAtTecPlus,
00109                         double &trackPhi, double &trackPhiRef, 
00110                         std::vector<GlobalPoint> &globHit, std::vector<GlobalPoint> &globPtrack,
00111                         GlobalPoint &globPref, std::vector<double> &hitPhiError);
00112 
00113   void buildTrajectory(TkFittedLasBeam &beam, unsigned int &hit,
00114                        vector<const GeomDetUnit*> &gd, std::vector<GlobalPoint> &globPtrack,
00115                        vector<TrajectoryStateOnSurface> &tsosLas, GlobalPoint &globPref);
00116 
00117   bool fitBeam(TkFittedLasBeam &beam, AlgebraicSymMatrix &covMatrix, 
00118                unsigned int &hitsAtTecPlus, unsigned int &nFitParams,
00119                double &offset, double &slope, vector<GlobalPoint> &globPtrack, 
00120                double &bsAngleParam, double &chi2);
00121 
00122   // ----------member data ---------------------------
00123   const edm::InputTag src_;
00124   bool fitBeamSplitters_;
00125   unsigned int nAtParameters_;
00126 
00127   edm::Service<TFileService> fs;
00128 
00129   // static parameters used in static parametrization functions
00130   static vector<double> gHitZprime;
00131   static vector<double> gBarrelModuleRadius;
00132   static vector<double> gBarrelModuleOffset;
00133   static float gTIBparam;
00134   static float gTOBparam;
00135   static double gBeamR;
00136   static double gBeamZ0;
00137   static double gBeamSplitterZprime;
00138   static unsigned int gHitsAtTecMinus;
00139   static double gBSparam;
00140   static bool gFitBeamSplitters;
00141   static bool gIsInnerBarrel;
00142 
00143   // histograms
00144   TH1F *h_bsAngle, *h_hitX, *h_hitXTecPlus, *h_hitXTecMinus,
00145     *h_hitXAt, *h_chi2, *h_chi2ndof, *h_pull, *h_res, 
00146     *h_resTecPlus, *h_resTecMinus, *h_resAt;
00147   TH2F *h_bsAngleVsBeam, *h_hitXvsZTecPlus, *h_hitXvsZTecMinus,
00148     *h_hitXvsZAt, *h_resVsZTecPlus, *h_resVsZTecMinus, *h_resVsZAt,
00149     *h_resVsHitTecPlus, *h_resVsHitTecMinus, *h_resVsHitAt;
00150 }; 
00151 
00152 //
00153 // constants, enums and typedefs
00154 //
00155 
00156 
00157 //
00158 // static data member definitions
00159 //
00160 
00161 // static parameters used in parametrization functions
00162 vector<double> TkLasBeamFitter::gHitZprime;
00163 vector<double> TkLasBeamFitter::gBarrelModuleRadius;
00164 vector<double> TkLasBeamFitter::gBarrelModuleOffset;
00165 float TkLasBeamFitter::gTIBparam = 0.097614; // = abs(r_offset/r_module) (nominal!)
00166 float TkLasBeamFitter::gTOBparam = 0.034949; // = abs(r_offset/r_module) (nominal!)
00167 double TkLasBeamFitter::gBeamR = 0.0;
00168 double TkLasBeamFitter::gBeamZ0 = 0.0; 
00169 double TkLasBeamFitter::gBeamSplitterZprime = 0.0;
00170 unsigned int TkLasBeamFitter::gHitsAtTecMinus = 0;
00171 double TkLasBeamFitter::gBSparam = 0.0;
00172 bool TkLasBeamFitter::gFitBeamSplitters = 0;
00173 bool TkLasBeamFitter::gIsInnerBarrel = 0;
00174 
00175 // handles
00176 Handle<TkLasBeamCollection> laserBeams;
00177 ESHandle<MagneticField> fieldHandle;
00178 ESHandle<TrackerGeometry> geometry;
00179 
00180 //
00181 // constructors and destructor
00182 //
00183 TkLasBeamFitter::TkLasBeamFitter(const edm::ParameterSet &iConfig) :
00184   src_(iConfig.getParameter<edm::InputTag>("src")),
00185   fitBeamSplitters_(iConfig.getParameter<bool>("fitBeamSplitters")),
00186   nAtParameters_(iConfig.getParameter<unsigned int>("numberOfFittedAtParameters")),
00187   h_bsAngle(0), h_hitX(0), h_hitXTecPlus(0), h_hitXTecMinus(0),
00188   h_hitXAt(0), h_chi2(0), h_chi2ndof(0), h_pull(0), h_res(0), 
00189   h_resTecPlus(0), h_resTecMinus(0), h_resAt(0),
00190   h_bsAngleVsBeam(0), h_hitXvsZTecPlus(0), h_hitXvsZTecMinus(0),
00191   h_hitXvsZAt(0), h_resVsZTecPlus(0), h_resVsZTecMinus(0), h_resVsZAt(0),
00192   h_resVsHitTecPlus(0), h_resVsHitTecMinus(0), h_resVsHitAt(0)
00193 {
00194   // declare the products to produce
00195   this->produces<TkFittedLasBeamCollection, edm::InRun>();
00196   this->produces<TsosVectorCollection, edm::InRun>();
00197   
00198   //now do what ever other initialization is needed
00199 }
00200 
00201 //---------------------------------------------------------------------------------------
00202 TkLasBeamFitter::~TkLasBeamFitter()
00203 {
00204  
00205    // do anything here that needs to be done at desctruction time
00206    // (e.g. close files, deallocate resources etc.)
00207 
00208 }
00209 
00210 
00211 //
00212 // member functions
00213 //
00214 
00215 //---------------------------------------------------------------------------------------
00216 // ------------ method called to produce the data  ------------
00217 void TkLasBeamFitter::produce(edm::Event &iEvent, const edm::EventSetup &setup)
00218 {
00219   // Nothing per event!
00220 }
00221 
00222 //---------------------------------------------------------------------------------------
00223 // ------------ method called at end of each run  ---------------------------------------
00224 void TkLasBeamFitter::endRun(edm::Run &run, const edm::EventSetup &setup)
00225 {
00226 // }
00227 // // FIXME!
00228 // // Indeed, that should be in endRun(..) - as soon as AlignmentProducer can call
00229 // // the algorithm's endRun correctly!
00230 //
00231 //
00232 // void TkLasBeamFitter::beginRun(edm::Run &run, const edm::EventSetup &setup)
00233 // {
00234 
00235   // book histograms
00236   h_hitX = fs->make<TH1F>("hitX","local x of LAS hits;local x [cm];N",100,-0.5,0.5);
00237   h_hitXTecPlus = fs->make<TH1F>("hitXTecPlus","local x of LAS hits in TECplus;local x [cm];N",100,-0.5,0.5);
00238   h_hitXTecMinus = fs->make<TH1F>("hitXTecMinus","local x of LAS hits in TECminus;local x [cm];N",100,-0.5,0.5);
00239   h_hitXAt = fs->make<TH1F>("hitXAt","local x of LAS hits in ATs;local x [cm];N",100,-2.5,2.5);
00240   h_hitXvsZTecPlus = fs->make<TH2F>("hitXvsZTecPlus","local x vs z in TECplus;z [cm];local x [cm]",80,120,280,100,-0.5,0.5);
00241   h_hitXvsZTecMinus = fs->make<TH2F>("hitXvsZTecMinus","local x vs z in TECMinus;z [cm];local x [cm]",80,-280,-120,100,-0.5,0.5);
00242   h_hitXvsZAt = fs->make<TH2F>("hitXvsZAt","local x vs z in ATs;z [cm];local x [cm]",200,-200,200,100,-0.5,0.5);
00243   h_chi2 = fs->make<TH1F>("chi2","#chi^{2};#chi^{2};N",100,0,2000);
00244   h_chi2ndof = fs->make<TH1F>("chi2ndof","#chi^{2} per degree of freedom;#chi^{2}/N_{dof};N",100,0,300);
00245   h_pull = fs->make<TH1F>("pull","pulls of #phi residuals;pull;N",50,-10,10);
00246   h_res = fs->make<TH1F>("res","#phi residuals;#phi_{track} - #phi_{hit} [rad];N",60,-0.0015,0.0015);
00247   h_resTecPlus = fs->make<TH1F>("resTecPlus","#phi residuals in TECplus;#phi_{track} - #phi_{hit} [rad];N",30,-0.0015,0.0015);
00248   h_resTecMinus = fs->make<TH1F>("resTecMinus","#phi residuals in TECminus;#phi_{track} - #phi_{hit} [rad];N",60,-0.0015,0.0015);
00249   h_resAt = fs->make<TH1F>("resAt","#phi residuals in ATs;#phi_{track} - #phi_{hit} [rad];N",30,-0.0015,0.0015);
00250   h_resVsZTecPlus = fs->make<TH2F>("resVsZTecPlus","phi residuals vs. z in TECplus;z [cm];#phi_{track} - #phi_{hit} [rad]",
00251                          80,120,280,100,-0.0015,0.0015);
00252   h_resVsZTecMinus = fs->make<TH2F>("resVsZTecMinus","phi residuals vs. z in TECminus;z [cm];#phi_{track} - #phi_{hit} [rad]",
00253                           80,-280,-120,100,-0.0015,0.0015);
00254   h_resVsZAt = fs->make<TH2F>("resVsZAt","#phi residuals vs. z in ATs;N;#phi_{track} - #phi_{hit} [rad]",
00255                     200,-200,200,100,-0.0015,0.0015);
00256   h_resVsHitTecPlus = fs->make<TH2F>("resVsHitTecPlus","#phi residuals vs. hits in TECplus;hit no.;#phi_{track} - #phi_{hit} [rad]",
00257                            144,0,144,100,-0.0015,0.0015);
00258   h_resVsHitTecMinus = fs->make<TH2F>("resVsHitTecMinus","#phi residuals vs. hits in TECminus;hit no.;#phi_{track} - #phi_{hit} [rad]",
00259                             144,0,144,100,-0.0015,0.0015);
00260   h_resVsHitAt = fs->make<TH2F>("resVsHitAt","#phi residuals vs. hits in ATs;hit no.;#phi_{track} - #phi_{hit} [rad]",
00261                       176,0,176,100,-0.0015,0.0015);
00262   h_bsAngle = fs->make<TH1F>("bsAngle","fitted beam splitter angle;BS angle [rad];N",40,-0.004,0.004);
00263   h_bsAngleVsBeam = fs->make<TH2F>("bsAngleVsBeam","fitted beam splitter angle per beam;Beam no.;BS angle [rad]",
00264                          40,0,300,100,-0.004,0.004);
00265 
00266   // Create output collections - they are parallel.
00267   // (edm::Ref etc. and thus edm::AssociationVector are not supported for edm::Run...)
00268   std::auto_ptr<TkFittedLasBeamCollection> fittedBeams(new TkFittedLasBeamCollection);
00269   // One std::vector<TSOS> for each TkFittedLasBeam:
00270   std::auto_ptr<TsosVectorCollection> tsosesVec(new TsosVectorCollection);
00271 
00272   // get TkLasBeams, Tracker geometry, magnetic field
00273   run.getByLabel( "LaserAlignment", "tkLaserBeams", laserBeams );
00274   setup.get<TrackerDigiGeometryRecord>().get(geometry);
00275   setup.get<IdealMagneticFieldRecord>().get(fieldHandle);
00276 
00277   // hack for fixed BSparams (ugly!)
00278 //   double bsParams[34] = {-0.000266,-0.000956,-0.001205,-0.000018,-0.000759,0.002554,
00279 //                       0.000465,0.000975,0.001006,0.002027,-0.001263,-0.000763,
00280 //                       -0.001702,0.000906,-0.002120,0.001594,0.000661,-0.000457,
00281 //                       -0.000447,0.000347,-0.002266,-0.000446,0.000659,0.000018,
00282 //                       -0.001630,-0.000324,
00283 //                       // ATs
00284 //                       -999.,-0.001709,-0.002091,-999.,
00285 //                       -0.001640,-999.,-0.002444,-0.002345};
00286 
00287   double bsParams[40] = {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,
00288                          0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
00289 
00290   // beam counter
00291   unsigned int beamNo(0);
00292   // fit BS? If false, values from bsParams are taken
00293   gFitBeamSplitters = fitBeamSplitters_;
00294   if(fitBeamSplitters_) cout << "Fitting BS!" << endl;
00295   else cout << "BS fixed, not fitted!" << endl;
00296 
00297   // loop over LAS beams
00298   for(TkLasBeamCollection::const_iterator iBeam = laserBeams->begin(), iEnd = laserBeams->end();
00299        iBeam != iEnd; ++iBeam){
00300 
00301     TkFittedLasBeam beam(*iBeam);
00302     vector<TrajectoryStateOnSurface> tsosLas;
00303 
00304     // set BS param for fit
00305     gBSparam = bsParams[beamNo];
00306 
00307     // call main function; all other functions are called inside getLasBeams(..)
00308     this->getLasBeams(beam, tsosLas);
00309     
00310     // fill output products
00311     fittedBeams->push_back(beam);
00312     tsosesVec->push_back(tsosLas);
00313 
00314 //     if(!this->fitBeam(fittedBeams->back(), tsosesVec->back())){
00315 //       edm::LogError("BadFit") 
00316 //       << "Problems fitting TkLasBeam, id " << fittedBeams->back().getBeamId() << ".";
00317 //       fittedBeams->pop_back(); // remove last entry added just before
00318 //       tsosesVec->pop_back();   // dito
00319 //     }
00320 
00321     beamNo++;
00322   }
00323   
00324   // finally, put fitted beams and TSOS vectors into run
00325   run.put(fittedBeams);
00326   run.put(tsosesVec);
00327 }
00328 
00329 // methods for las data processing
00330 
00331 // -------------- loop over beams, call functions ----------------------------
00332 void TkLasBeamFitter::getLasBeams(TkFittedLasBeam &beam, vector<TrajectoryStateOnSurface> &tsosLas)
00333 {
00334   cout << "---------------------------------------" << endl;
00335   cout << "beam id: " << beam.getBeamId() // << " isTec: " << (beam.isTecInternal() ? "Y" : "N") 
00336        << " isTec+: " << (beam.isTecInternal(1) ? "Y" : "N") << " isTec-: " << (beam.isTecInternal(-1) ? "Y" : "N")
00337        << " isAt: " << (beam.isAlignmentTube() ? "Y" : "N") << " isR6: " << (beam.isRing6() ? "Y" : "N")
00338        << endl;
00339  
00340   // reset static variables 
00341   gHitsAtTecMinus = 0;
00342   gHitZprime.clear();
00343   gBarrelModuleRadius.clear();
00344   gBarrelModuleOffset.clear();
00345 
00346   // set right beam radius
00347   gBeamR = beam.isRing6() ? 84.0 : 56.4;
00348   
00349   vector<const GeomDetUnit*> gd;
00350   vector<GlobalPoint> globHit;
00351   unsigned int hitsAtTecPlus(0);
00352   double sumZ(0.);
00353 
00354   // loop over hits
00355   for( TkLasBeam::const_iterator iHit = beam.begin(); iHit < beam.end(); ++iHit ){
00356     // iHit is a SiStripLaserRecHit2D
00357 
00358     const SiStripLaserRecHit2D hit(*iHit);
00359 
00360     this->getLasHits(beam, hit, gd, globHit, hitsAtTecPlus);
00361     sumZ += globHit.back().z();
00362 
00363     // fill histos
00364     h_hitX->Fill(hit.localPosition().x());
00365     // TECplus
00366     if(beam.isTecInternal(1)){
00367       h_hitXTecPlus->Fill(hit.localPosition().x());
00368       h_hitXvsZTecPlus->Fill(globHit.back().z(),hit.localPosition().x());
00369     }
00370     // TECminus
00371     else if(beam.isTecInternal(-1)){
00372       h_hitXTecMinus->Fill(hit.localPosition().x());
00373       h_hitXvsZTecMinus->Fill(globHit.back().z(),hit.localPosition().x());
00374     }
00375     // ATs
00376     else{
00377       h_hitXAt->Fill(hit.localPosition().x());
00378       h_hitXvsZAt->Fill(globHit.back().z(),hit.localPosition().x());
00379     }  
00380   }
00381 
00382   gBeamZ0 = sumZ / globHit.size();
00383   double zMin(0.), zMax(0.);
00384   // TECplus
00385   if(beam.isTecInternal(1)){
00386     gBeamSplitterZprime = 205.75 - gBeamZ0;
00387     zMin = 120.0 - gBeamZ0;
00388     zMax = 280.0 - gBeamZ0;
00389   }
00390   // TECminus
00391   else if(beam.isTecInternal(-1)){
00392     gBeamSplitterZprime = -205.75 - gBeamZ0;
00393     zMin = -280.0 - gBeamZ0;
00394     zMax = -120.0 - gBeamZ0;
00395   }
00396   // AT
00397   else{
00398     gBeamSplitterZprime = 112.3 - gBeamZ0;
00399     zMin = -200.0 - gBeamZ0;
00400     zMax = 200.0 - gBeamZ0;
00401   }
00402 
00403   // fill vectors for fitted quantities
00404   vector<double> hitPhi, hitPhiError, hitZprimeError;
00405 
00406   for(unsigned int hit = 0; hit < globHit.size(); ++hit){
00407     hitPhi.push_back(static_cast<double>(globHit[hit].phi()));
00408     // localPositionError[hit] or assume 0.003, 0.006
00409     hitPhiError.push_back( 0.003 / globHit[hit].perp());
00410     // no errors on z, fill with zeros
00411     hitZprimeError.push_back(0.0);
00412     // barrel-specific values
00413     if(beam.isAlignmentTube() && abs(globHit[hit].z()) < 112.3){
00414       gBarrelModuleRadius.push_back(globHit[hit].perp());
00415       gBarrelModuleOffset.push_back(gBarrelModuleRadius.back() - gBeamR);
00416       // TIB/TOB flag
00417       if(gBarrelModuleOffset.back() < 0.0){
00418         gIsInnerBarrel = 1;
00419       }
00420       else{
00421         gIsInnerBarrel = 0;
00422       }
00423       gHitZprime.push_back(globHit[hit].z() - gBeamZ0 - abs(gBarrelModuleOffset.back())); 
00424     }
00425     // non-barrel z'
00426     else{
00427       gHitZprime.push_back(globHit[hit].z() - gBeamZ0);
00428     }
00429   }
00430 
00431   // number of fit parameters, 3 for TECs (always!); 3, 5, or 6 for ATs
00432   unsigned int tecParams(3), atParams(0);
00433   if(nAtParameters_ == 3) atParams = 3;
00434   else if(nAtParameters_ == 5) atParams = 5;
00435   else atParams = 6;                            // <-- default value, recommended
00436   unsigned int nFitParams(0);
00437   if(!fitBeamSplitters_ || 
00438      (hitsAtTecPlus == 0 && beam.isAlignmentTube() ) ){
00439     tecParams = tecParams - 1;
00440     atParams = atParams - 1;
00441   }
00442   if(beam.isTecInternal()){
00443     nFitParams = tecParams;
00444   }
00445   else{
00446     nFitParams = atParams;
00447   }
00448   
00449   // fit parameter definitions
00450   double offset(0.), offsetError(0.), slope(0.), slopeError(0.),
00451     bsAngleParam(0.), phiAtMinusParam(0.), phiAtPlusParam(0.),
00452     atThetaSplitParam(0.);
00453   AlgebraicSymMatrix covMatrix;
00454   if(!fitBeamSplitters_ || (beam.isAlignmentTube() && hitsAtTecPlus == 0)){
00455     covMatrix = AlgebraicSymMatrix(nFitParams, 1);
00456   }
00457   else{
00458     covMatrix = AlgebraicSymMatrix(nFitParams - 1, 1);
00459   }
00460 
00461   this->fitter(beam, covMatrix, 
00462                hitsAtTecPlus, nFitParams, 
00463                hitPhi, hitPhiError, hitZprimeError, 
00464                zMin, zMax, bsAngleParam,
00465                offset, offsetError, slope, slopeError,
00466                phiAtMinusParam, phiAtPlusParam,
00467                atThetaSplitParam); 
00468   
00469   vector<GlobalPoint> globPtrack;
00470   GlobalPoint globPref;
00471   double chi2(0.);
00472 
00473   for(unsigned int hit = 0; hit < gHitZprime.size(); ++hit){
00474     
00475     // additional phi value (trackPhiRef) for trajectory calculation
00476     double trackPhi(0.), trackPhiRef(0.); 
00477     
00478     this->trackPhi(beam, hit, trackPhi, trackPhiRef,
00479                    offset, slope, bsAngleParam, 
00480                    phiAtMinusParam, phiAtPlusParam,
00481                    atThetaSplitParam, globHit);
00482 
00483     cout << "track phi = " << trackPhi 
00484          << ", hit phi = " << hitPhi[hit] 
00485          << ", zPrime = " << gHitZprime[hit] 
00486          << " r = " << globHit[hit].perp() << endl;
00487 
00488     this->globalTrackPoint(beam, hit, hitsAtTecPlus, 
00489                            trackPhi, trackPhiRef, 
00490                            globHit, globPtrack, globPref, 
00491                            hitPhiError);
00492 
00493     // calculate residuals = pred - hit (in global phi)
00494     const double phiResidual = globPtrack[hit].phi() - globHit[hit].phi();
00495     // pull calculation (FIX!!!)
00496     const double phiResidualPull = phiResidual / hitPhiError[hit];
00497     //       sqrt(hitPhiError[hit]*hitPhiError[hit] + 
00498     //     (offsetError*offsetError + globPtrack[hit].z()*globPtrack[hit].z() * slopeError*slopeError));
00499 
00500     // calculate chi2
00501     chi2 += phiResidual*phiResidual / (hitPhiError[hit]*hitPhiError[hit]);
00502     
00503     // fill histos
00504     h_res->Fill(phiResidual);
00505     // TECplus
00506     if(beam.isTecInternal(1)){
00507       h_pull->Fill(phiResidualPull);
00508       h_resTecPlus->Fill(phiResidual);
00509       h_resVsZTecPlus->Fill(globPtrack[hit].z(), phiResidual);
00510       // Ring 6
00511       if(beam.isRing6()){
00512         h_resVsHitTecPlus->Fill(hit+(beam.getBeamId()-1)/10*9+72, phiResidual);
00513       }
00514       // Ring 4
00515       else{
00516         h_resVsHitTecPlus->Fill(hit+beam.getBeamId()/10*9, phiResidual);
00517       }
00518     }
00519     // TECminus
00520     else if(beam.isTecInternal(-1)){
00521       h_pull->Fill(phiResidualPull);
00522       h_resTecMinus->Fill(phiResidual);
00523       h_resVsZTecMinus->Fill(globPtrack[hit].z(), phiResidual);
00524       // Ring 6
00525       if(beam.isRing6()){
00526         h_resVsHitTecMinus->Fill(hit+(beam.getBeamId()-101)/10*9+72, phiResidual);
00527       }
00528       // Ring 4
00529       else{
00530         h_resVsHitTecMinus->Fill(hit+(beam.getBeamId()-100)/10*9, phiResidual);
00531       }
00532     }
00533     // ATs
00534     else{
00535       h_pull->Fill(phiResidualPull);
00536       h_resAt->Fill(phiResidual);
00537       h_resVsZAt->Fill(globPtrack[hit].z(), phiResidual);
00538       h_resVsHitAt->Fill(hit+(beam.getBeamId()-200)/10*22, phiResidual);
00539     }
00540 
00541     this->buildTrajectory(beam, hit, gd, globPtrack, tsosLas, globPref);
00542   }
00543 
00544   cout << "chi^2 = " << chi2 << ", chi^2/ndof = " << chi2/(gHitZprime.size() - nFitParams) << endl;
00545   this->fitBeam(beam, covMatrix, hitsAtTecPlus, nFitParams,
00546                 offset, slope, globPtrack, bsAngleParam, chi2);
00547   
00548   cout << "bsAngleParam = " << bsAngleParam << endl;
00549 
00550   // fill histos
00551   // include slope, offset, covariance plots here
00552   h_chi2->Fill(chi2);
00553   h_chi2ndof->Fill(chi2/(gHitZprime.size() - nFitParams));
00554   if(bsAngleParam != 0.0){
00555     h_bsAngle->Fill(2.0*atan(0.5*bsAngleParam));
00556     h_bsAngleVsBeam->Fill(beam.getBeamId(), 2.0*atan(0.5*bsAngleParam));
00557   }
00558 }
00559 
00560 
00561 // --------- get hits, convert to global coordinates ---------------------------
00562 void TkLasBeamFitter::getLasHits(TkFittedLasBeam &beam, const SiStripLaserRecHit2D &hit, 
00563                                  vector<const GeomDetUnit*> &gd, vector<GlobalPoint> &globHit,
00564                                  unsigned int &hitsAtTecPlus)
00565 { 
00566   // get global position of LAS hits
00567   gd.push_back(geometry->idToDetUnit(hit.getDetId()));
00568   GlobalPoint globPtemp(gd.back()->toGlobal(hit.localPosition()));
00569 
00570   // testing: globPtemp should be right
00571   globHit.push_back(globPtemp);
00572 
00573   if(beam.isAlignmentTube()){
00574     if(abs(globPtemp.z()) > 112.3){
00575       if(globPtemp.z() < 112.3) gHitsAtTecMinus++ ;
00576       else hitsAtTecPlus++ ;
00577     }
00578   }
00579 }
00580 
00581 
00582 // ------------ parametrization functions for las beam fits ------------
00583 double TkLasBeamFitter::tecPlusFunction(double *x, double *par)
00584 {
00585   double z = x[0]; // 'primed'? -> yes!!!
00586 
00587   if(z < gBeamSplitterZprime){
00588     return par[0] + par[1] * z;
00589   } 
00590   else{
00591     if(gFitBeamSplitters){
00592       // par[2] = 2*tan(BeamSplitterAngle/2.0)
00593       return par[0] + par[1] * z - par[2] * (z - gBeamSplitterZprime)/gBeamR; 
00594     }
00595     else{
00596       return par[0] + par[1] * z - gBSparam * (z - gBeamSplitterZprime)/gBeamR; 
00597     }
00598   }
00599 }
00600 
00601 
00602 double TkLasBeamFitter::tecMinusFunction(double *x, double *par)
00603 {
00604   double z = x[0]; // 'primed'? -> yes!!!
00605 
00606   if(z > gBeamSplitterZprime){
00607     return par[0] + par[1] * z;
00608   } 
00609   else{
00610     if(gFitBeamSplitters){
00611       // par[2] = 2*tan(BeamSplitterAngle/2.0)
00612       return par[0] + par[1] * z + par[2] * (z - gBeamSplitterZprime)/gBeamR; 
00613     }
00614     else{
00615       return par[0] + par[1] * z + gBSparam * (z - gBeamSplitterZprime)/gBeamR; 
00616     }
00617   }
00618 }
00619  
00620 double TkLasBeamFitter::atFunction(double *x, double *par)
00621 {
00622   double z = x[0]; // 'primed'? -> yes!!!
00623   // TECminus
00624   if(z < -gBeamSplitterZprime - 2.0*gBeamZ0){
00625     return par[0] + par[1] * z;
00626   }
00627   // BarrelMinus
00628   else if(-gBeamSplitterZprime - 2.0*gBeamZ0 < z && z < -gBeamZ0){
00629     // z value includes module offset from main beam axis
00630     // TOB
00631     if(!gIsInnerBarrel){
00632       return par[0] + par[1] * z + gTOBparam * (par[2] + par[4]);
00633     }
00634     // TIB
00635     else{
00636       return par[0] + par[1] * z - gTIBparam * (par[2] - par[4]);
00637     }
00638   }
00639   // BarrelPlus
00640   else if(-gBeamZ0 < z && z < gBeamSplitterZprime){
00641     // z value includes module offset from main beam axis
00642     // TOB
00643     if(!gIsInnerBarrel){
00644       return par[0] + par[1] * z + gTOBparam * (par[3] - par[4]);
00645     }
00646     // TIB
00647     else{
00648       return par[0] + par[1] * z - gTIBparam * (par[3] + par[4]);
00649     }
00650   }
00651   // TECplus
00652   else{
00653     if(gFitBeamSplitters){
00654       // par[2] = 2*tan(BeamSplitterAngle/2.0)
00655       return par[0] + par[1] * z - par[5] * (z - gBeamSplitterZprime)/gBeamR; // BS par: 5, 4, or 2
00656     }
00657     else{
00658       return par[0] + par[1] * z - gBSparam * (z - gBeamSplitterZprime)/gBeamR;
00659     }
00660   }
00661 }
00662 
00663 
00664 // ------------ perform fit of beams ------------------------------------
00665 void TkLasBeamFitter::fitter(TkFittedLasBeam &beam, AlgebraicSymMatrix &covMatrix, 
00666                              unsigned int &hitsAtTecPlus, unsigned int &nFitParams,
00667                              vector<double> &hitPhi, vector<double> &hitPhiError, vector<double> &hitZprimeError,
00668                              double &zMin, double &zMax, double &bsAngleParam, 
00669                              double &offset, double &offsetError, double &slope, double &slopeError,
00670                              double &phiAtMinusParam, double &phiAtPlusParam, double &atThetaSplitParam)
00671 {
00672   TGraphErrors *lasData = new TGraphErrors(gHitZprime.size(), 
00673                                            &(gHitZprime[0]), &(hitPhi[0]), 
00674                                            &(hitZprimeError[0]), &(hitPhiError[0]));
00675   
00676   // do fit (R = entire range)
00677   if(beam.isTecInternal(1)){
00678     TF1 tecPlus("tecPlus", tecPlusFunction, zMin, zMax, nFitParams );
00679     tecPlus.SetParameter( 1, 0 ); // slope
00680     tecPlus.SetParameter( nFitParams - 1, 0 ); // BS 
00681     lasData->Fit("tecPlus", "R"); // "R", "RV" or "RQ"
00682   }
00683   else if(beam.isTecInternal(-1)){
00684     TF1 tecMinus("tecMinus", tecMinusFunction, zMin, zMax, nFitParams );
00685     tecMinus.SetParameter( 1, 0 ); // slope
00686     tecMinus.SetParameter( nFitParams - 1, 0 ); // BS 
00687     lasData->Fit("tecMinus", "R");
00688   }
00689   else{
00690     TF1 at("at", atFunction, zMin, zMax, nFitParams );
00691     at.SetParameter( 1, 0 ); // slope
00692     at.SetParameter( nFitParams - 1, 0 ); // BS 
00693     lasData->Fit("at","R");
00694   }
00695   
00696   // get values and errors for offset and slope
00697   gMinuit->GetParameter(0, offset, offsetError);
00698   gMinuit->GetParameter(1, slope, slopeError);
00699 
00700   // additional AT parameters
00701   // define param errors that are not used later
00702   double bsAngleParamError(0.), phiAtMinusParamError(0.),
00703     phiAtPlusParamError(0.), atThetaSplitParamError(0.);
00704 
00705   if(beam.isAlignmentTube()){
00706     gMinuit->GetParameter(2, phiAtMinusParam, phiAtMinusParamError);
00707     gMinuit->GetParameter(3, phiAtPlusParam, phiAtPlusParamError);
00708     gMinuit->GetParameter(4, atThetaSplitParam, atThetaSplitParamError);
00709   }
00710   // get Beam Splitter parameters
00711   if(fitBeamSplitters_){
00712     if(beam.isAlignmentTube() && hitsAtTecPlus == 0){
00713       bsAngleParam = gBSparam;
00714     }
00715     else{
00716       gMinuit->GetParameter( nFitParams - 1 , bsAngleParam, bsAngleParamError);
00717     }
00718   }
00719   else{
00720     bsAngleParam = gBSparam;
00721   }
00722 
00723   // fill covariance matrix
00724   vector<double> vec( covMatrix.num_col() * covMatrix.num_col() );
00725   gMinuit->mnemat( &vec[0], covMatrix.num_col() );
00726   for(int col = 0; col < covMatrix.num_col(); col++){
00727     for(int row = 0; row < covMatrix.num_col(); row++){
00728       covMatrix[col][row] = vec[row + covMatrix.num_col()*col];
00729     }
00730   }
00731   // compute correlation between parameters
00732 //   double corr01 = covMatrix[1][0]/(offsetError*slopeError);
00733 
00734   delete lasData;
00735 }
00736 
00737 
00738 // -------------- calculate track phi value ----------------------------------
00739 void TkLasBeamFitter::trackPhi(TkFittedLasBeam &beam, unsigned int &hit, 
00740                                double &trackPhi, double &trackPhiRef,
00741                                double &offset, double &slope, double &bsAngleParam,
00742                                double &phiAtMinusParam, double &phiAtPlusParam,
00743                                double &atThetaSplitParam, vector<GlobalPoint> &globHit)
00744 {
00745   // TECplus
00746   if(beam.isTecInternal(1)){
00747     if(gHitZprime[hit] < gBeamSplitterZprime){
00748       trackPhi = offset + slope * gHitZprime[hit];
00749       trackPhiRef = offset + slope * (gHitZprime[hit] + 1.0);
00750     }
00751     else{
00752       trackPhi = offset + slope * gHitZprime[hit] 
00753         - bsAngleParam * (gHitZprime[hit] - gBeamSplitterZprime)/gBeamR;
00754       trackPhiRef = offset + slope * (gHitZprime[hit] + 1.0)
00755         - bsAngleParam * ((gHitZprime[hit] + 1.0) - gBeamSplitterZprime)/gBeamR;
00756     }
00757   }
00758   // TECminus
00759   else if(beam.isTecInternal(-1)){
00760     if(gHitZprime[hit] > gBeamSplitterZprime){
00761       trackPhi = offset + slope * gHitZprime[hit];
00762       trackPhiRef = offset + slope * (gHitZprime[hit] + 1.0);
00763     }
00764     else{
00765       trackPhi = offset + slope * gHitZprime[hit] 
00766         + bsAngleParam * (gHitZprime[hit] - gBeamSplitterZprime)/gBeamR;
00767       trackPhiRef = offset + slope * (gHitZprime[hit] + 1.0)
00768         + bsAngleParam * ((gHitZprime[hit] + 1.0) - gBeamSplitterZprime)/gBeamR;
00769     }
00770   }
00771   // ATs
00772   else{
00773     // TECminus
00774     if(gHitZprime[hit] < -gBeamSplitterZprime - 2.0*gBeamZ0){
00775       trackPhi = offset + slope * gHitZprime[hit];
00776       trackPhiRef = offset + slope * (gHitZprime[hit] + 1.0);
00777     }
00778     // BarrelMinus
00779     else if(-gBeamSplitterZprime - 2.0*gBeamZ0 < gHitZprime[hit] && gHitZprime[hit] < -gBeamZ0){
00780       if(!gIsInnerBarrel){
00781         trackPhi = offset + slope * gHitZprime[hit] + gTOBparam * (phiAtMinusParam + atThetaSplitParam); 
00782       }
00783       else{
00784         trackPhi = offset + slope * gHitZprime[hit] - gTIBparam * (phiAtMinusParam - atThetaSplitParam);
00785       }
00786       trackPhiRef = offset + slope * (gHitZprime[hit] + abs(gBarrelModuleOffset[hit-gHitsAtTecMinus]));
00787     }
00788     // BarrelPlus
00789     else if(-gBeamZ0 < gHitZprime[hit] && gHitZprime[hit] < gBeamSplitterZprime){
00790       if(!gIsInnerBarrel){
00791         trackPhi = offset + slope * gHitZprime[hit] + gTOBparam * (phiAtPlusParam - atThetaSplitParam);
00792       }
00793       else{
00794         trackPhi = offset + slope * gHitZprime[hit] - gTIBparam * (phiAtPlusParam + atThetaSplitParam);
00795       }
00796       trackPhiRef = offset + slope * (gHitZprime[hit] + abs(gBarrelModuleOffset[hit-gHitsAtTecMinus]));
00797     }
00798     // TECplus
00799     else{
00800       trackPhi = offset + slope * gHitZprime[hit] 
00801         - bsAngleParam * (gHitZprime[hit] - gBeamSplitterZprime)/gBeamR;
00802       trackPhiRef = offset + slope * (gHitZprime[hit] + 1.0)
00803         - bsAngleParam * ((gHitZprime[hit] + 1.0) - gBeamSplitterZprime)/gBeamR;
00804     }
00805   }
00806 }
00807 
00808 
00809 // -------------- calculate global track points, hit residuals, chi2 ----------------------------------
00810 void TkLasBeamFitter::globalTrackPoint(TkFittedLasBeam &beam, 
00811                                        unsigned int &hit, unsigned int &hitsAtTecPlus,
00812                                        double &trackPhi, double &trackPhiRef, 
00813                                        vector<GlobalPoint> &globHit, vector<GlobalPoint> &globPtrack,
00814                                        GlobalPoint &globPref, vector<double> &hitPhiError)
00815 {
00816   // TECs
00817   if(beam.isTecInternal(0)){
00818     globPtrack.push_back(GlobalPoint(GlobalPoint::Cylindrical(gBeamR, trackPhi, globHit[hit].z())));
00819     globPref = GlobalPoint(GlobalPoint::Cylindrical(gBeamR, trackPhiRef, globHit[hit].z() + 1.0));
00820   }
00821   // ATs
00822   else{
00823     // TECminus
00824     if(hit < gHitsAtTecMinus){ // gHitZprime[hit] < -gBeamSplitterZprime - 2.0*gBeamZ0
00825       globPtrack.push_back(GlobalPoint(GlobalPoint::Cylindrical(gBeamR, trackPhi, globHit[hit].z())));
00826       globPref = GlobalPoint(GlobalPoint::Cylindrical(gBeamR, trackPhiRef, globHit[hit].z() + 1.0));
00827     }
00828     // TECplus
00829     else if(hit > gHitZprime.size() - hitsAtTecPlus - 1){ // gHitZprime[hit] > gBeamSplitterZprime
00830       globPtrack.push_back(GlobalPoint(GlobalPoint::Cylindrical(gBeamR, trackPhi, globHit[hit].z())));
00831       globPref = GlobalPoint(GlobalPoint::Cylindrical(gBeamR, trackPhiRef, globHit[hit].z() + 1.0));
00832     }
00833     // Barrel
00834     else{
00835       globPtrack.push_back(GlobalPoint(GlobalPoint::Cylindrical(globHit[hit].perp(), trackPhi, globHit[hit].z())));
00836       globPref = GlobalPoint(GlobalPoint::Cylindrical(gBeamR, trackPhiRef, globHit[hit].z()));
00837     }
00838   }
00839 }
00840 
00841 
00842 // ----------- create TrajectoryStateOnSurface for each track hit ----------------------------------------------
00843 void TkLasBeamFitter::buildTrajectory(TkFittedLasBeam &beam, unsigned int &hit, 
00844                                       vector<const GeomDetUnit*> &gd, vector<GlobalPoint> &globPtrack,
00845                                       vector<TrajectoryStateOnSurface> &tsosLas, GlobalPoint &globPref)
00846 {
00847   const MagneticField* magneticField = fieldHandle.product();
00848   GlobalVector trajectoryState;
00849   
00850   // TECplus
00851   if(beam.isTecInternal(1)){
00852     trajectoryState = GlobalVector(globPref-globPtrack[hit]);
00853   }
00854   // TECminus
00855   else if(beam.isTecInternal(-1)){
00856     trajectoryState = GlobalVector(globPtrack[hit]-globPref);
00857   }
00858   // ATs
00859   else{
00860     // TECminus
00861     if(gHitZprime[hit] < -gBeamSplitterZprime - 2.0*gBeamZ0){
00862       trajectoryState = GlobalVector(globPtrack[hit]-globPref);
00863     }
00864     // TECplus
00865     else if(gHitZprime[hit] > gBeamSplitterZprime){
00866       trajectoryState = GlobalVector(globPref-globPtrack[hit]);
00867     }
00868     // Barrel
00869     else{
00870       trajectoryState = GlobalVector(globPtrack[hit]-globPref);
00871     }
00872   }     
00873 //   cout << "trajectory: " << trajectoryState << endl;                 
00874   const FreeTrajectoryState ftsLas = FreeTrajectoryState(globPtrack[hit],trajectoryState,0,magneticField);
00875   tsosLas.push_back(TrajectoryStateOnSurface(ftsLas,gd[hit]->surface(),
00876                                              SurfaceSideDefinition::beforeSurface));
00877 }
00878 
00879 
00880 //---------------------- set beam parameters for fittedBeams ---------------------------------
00881 bool TkLasBeamFitter::fitBeam(TkFittedLasBeam &beam, AlgebraicSymMatrix &covMatrix, 
00882                               unsigned int &hitsAtTecPlus, unsigned int &nFitParams,
00883                               double &offset, double &slope, vector<GlobalPoint> &globPtrack,
00884                               double &bsAngleParam, double &chi2)
00885 {
00886   // set beam parameters for beam output
00887   unsigned int paramType(0);
00888   if(!fitBeamSplitters_) paramType = 1;
00889   if(beam.isAlignmentTube() && hitsAtTecPlus == 0) paramType = 0;
00890 //   const unsigned int nPedeParams = nFitParams + paramType;
00891 
00892   // test without BS params
00893   const unsigned int nPedeParams(nFitParams);
00894 //   cout << "number of Pede parameters: " << nPedeParams << endl;
00895 
00896   std::vector<TkFittedLasBeam::Scalar> params(nPedeParams);
00897   params[0] = offset;
00898   params[1] = slope;
00899   // no BS parameter for AT beams without TECplus hits
00900 //   if(beam.isTecInternal() || hitsAtTecPlus > 0) params[2] = bsAngleParam;
00901 
00902   AlgebraicMatrix derivatives(gHitZprime.size(), nPedeParams);
00903   // fill derivatives matrix with local track derivatives
00904   for(unsigned int hit = 0; hit < gHitZprime.size(); ++hit){
00905     
00906     // d(delta phi)/d(offset) is identical for every hit
00907     derivatives[hit][0] = 1.0;
00908 
00909     // d(delta phi)/d(slope) and d(delta phi)/d(bsAngleParam) depend on parametrizations
00910     // TECplus
00911     if(beam.isTecInternal(1)){
00912       derivatives[hit][1] = globPtrack[hit].z();
00913 //       if(gHitZprime[hit] < gBeamSplitterZprime){
00914 //      derivatives[hit][2] = 0.0;
00915 //       }
00916 //       else{
00917 //      derivatives[hit][2] = - (globPtrack[hit].z() - gBeamSplitterZprime) / gBeamR;
00918 //       }
00919     }
00920     // TECminus
00921     else if(beam.isTecInternal(-1)){
00922       derivatives[hit][1] = globPtrack[hit].z();
00923 //       if(gHitZprime[hit] > gBeamSplitterZprime){
00924 //      derivatives[hit][2] = 0.0;
00925 //       }
00926 //       else{
00927 //      derivatives[hit][2] = (globPtrack[hit].z() - gBeamSplitterZprime) / gBeamR;
00928 //       }
00929     }
00930     // ATs
00931     else{
00932       // TECminus
00933       if(gHitZprime[hit] < -gBeamSplitterZprime - 2.0*gBeamZ0){
00934         derivatives[hit][1] = globPtrack[hit].z();
00935 //      if(hitsAtTecPlus > 0){
00936 //        derivatives[hit][2] = 0.0;
00937 //      }
00938       }
00939       // TECplus
00940       else if(gHitZprime[hit] > gBeamSplitterZprime){
00941         derivatives[hit][1] = globPtrack[hit].z();
00942 //      if(hitsAtTecPlus > 0){
00943 //        derivatives[hit][2] = - (globPtrack[hit].z() - gBeamSplitterZprime) / gBeamR;
00944 //      }
00945       }
00946       // Barrel
00947       else{
00948         derivatives[hit][1] = globPtrack[hit].z() - gBarrelModuleOffset[hit-gHitsAtTecMinus];
00949 //      if(hitsAtTecPlus > 0){
00950 //        derivatives[hit][2] = 0.0;
00951 //      }
00952       }
00953     }
00954   }
00955 
00956   unsigned int firstFixedParam(covMatrix.num_col()); // FIXME! --> no, is fine!!!
00957 //   unsigned int firstFixedParam = nPedeParams - 1;
00958 //   if(beam.isAlignmentTube() && hitsAtTecPlus == 0) firstFixedParam = nPedeParams;
00959 //   cout << "first fixed parameter: " << firstFixedParam << endl;
00960   // set fit results
00961   beam.setParameters(paramType, params, covMatrix, derivatives, firstFixedParam, chi2);
00962 
00963   return true; // return false in case of problems
00964 }
00965 
00966 
00967 //---------------------------------------------------------------------------------------
00968 //define this as a plug-in
00969 DEFINE_FWK_MODULE(TkLasBeamFitter);