CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

cms::FastMuPropagator Class Reference

#include <FastMuPropagator.h>

Inheritance diagram for cms::FastMuPropagator:
Propagator

List of all members.

Public Member Functions

virtual FastMuPropagatorclone () const
 FastMuPropagator (const MagneticField *mf, FmpConst *fmp, PropagationDirection dir=alongMomentum)
 FastMuPropagator (const MagneticField *mf, PropagationDirection dir=alongMomentum)
virtual const MagneticFieldmagneticField () const
TrajectoryStateOnSurface propagate (const FreeTrajectoryState &fts, const Cylinder &bound) const
TrajectoryStateOnSurface propagate (const FreeTrajectoryState &fts, const Plane &) const
TrajectoryStateOnSurface propagate (const FreeTrajectoryState &fts, const Surface &surface) const
virtual std::pair
< TrajectoryStateOnSurface,
double > 
propagateWithPath (const FreeTrajectoryState &state, const Plane &bc) const
virtual std::pair
< TrajectoryStateOnSurface,
double > 
propagateWithPath (const FreeTrajectoryState &state, const Cylinder &bc) const
virtual ~FastMuPropagator ()

Private Member Functions

bool checkfts (const FreeTrajectoryState &fts) const

Private Attributes

const MagneticFieldfield
FmpConsttheFmpConst

Detailed Description

Definition at line 18 of file FastMuPropagator.h.


Constructor & Destructor Documentation

cms::FastMuPropagator::FastMuPropagator ( const MagneticField mf,
PropagationDirection  dir = alongMomentum 
) [inline]

Definition at line 20 of file FastMuPropagator.h.

References field, and theFmpConst.

Referenced by clone().

                                    {theFmpConst=new FmpConst(); field = mf;} 
cms::FastMuPropagator::FastMuPropagator ( const MagneticField mf,
FmpConst fmp,
PropagationDirection  dir = alongMomentum 
) [inline]

Definition at line 23 of file FastMuPropagator.h.

References field, and theFmpConst.

                                    {theFmpConst=fmp; field = mf;}
virtual cms::FastMuPropagator::~FastMuPropagator ( ) [inline, virtual]

Definition at line 28 of file FastMuPropagator.h.

                               {
 //   delete theFmpConst;
  }

Member Function Documentation

bool cms::FastMuPropagator::checkfts ( const FreeTrajectoryState fts) const [private]

Definition at line 283 of file FastMuPropagator.cc.

References abs, CastorDataFrameFilter_impl::check(), gather_cfg::cout, FreeTrajectoryState::parameters(), PV3DBase< T, PVType, FrameType >::perp(), GlobalTrajectoryParameters::position(), alignCSCRings::r, PV3DBase< T, PVType, FrameType >::z(), and z.

                                                                        {
          bool check=true;
          double z=fts.parameters().position().z();
          double r=fts.parameters().position().perp();
          float mubarrelrad=theFmpConst->mubarrelrad;
          float muforwardrad=theFmpConst->muforwardrad;
          float muoffset=theFmpConst->muoffset;
          mubarrelrad=350.; muforwardrad=400.;
#ifdef PROPAGATOR_DB
        cout<<"checkfts::r,z="<<r<<" "<<z<<" "<<check<<endl;
#endif    
          if(r<mubarrelrad-muoffset&&abs(z)<muforwardrad-muoffset){
          check=false;
#ifdef PROPAGATOR_DB
        cout<<"checkfts::false="<<check<<endl;
#endif    
          }   
          return check;
}
virtual FastMuPropagator* cms::FastMuPropagator::clone ( void  ) const [inline, virtual]

Implements Propagator.

Definition at line 38 of file FastMuPropagator.h.

References alongMomentum, dir, FastMuPropagator(), and field.

virtual const MagneticField* cms::FastMuPropagator::magneticField ( ) const [inline, virtual]

Implements Propagator.

Definition at line 60 of file FastMuPropagator.h.

References field.

{return field;}
TrajectoryStateOnSurface cms::FastMuPropagator::propagate ( const FreeTrajectoryState state,
const Surface sur 
) const [inline, virtual]

Propagate from a free state (e.g. position and momentum in in global cartesian coordinates) to a surface. Only use the generic method if the surface type (plane or cylinder) is not known at the calling point.

Reimplemented from Propagator.

Definition at line 44 of file FastMuPropagator.h.

References propagate().

                                                                   {
    return Propagator::propagate( fts, surface);
  }
TrajectoryStateOnSurface cms::FastMuPropagator::propagate ( const FreeTrajectoryState fts,
const Plane boundplane 
) const [virtual]

Implements Propagator.

Definition at line 173 of file FastMuPropagator.cc.

References abs, GlobalTrajectoryParameters::charge(), DeDxDiscriminatorTools::charge(), funct::cos(), gather_cfg::cout, FreeTrajectoryState::curvilinearError(), m, CurvilinearTrajectoryError::matrix(), GlobalTrajectoryParameters::momentum(), FreeTrajectoryState::parameters(), PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::phi(), phi, GloballyPositioned< T >::position(), GlobalTrajectoryParameters::position(), ptmin, alignCSCRings::r, funct::sin(), launcher::step, funct::tan(), theta(), PV3DBase< T, PVType, FrameType >::theta(), PV3DBase< T, PVType, FrameType >::z(), and z.

{
  TrajectoryStateOnSurface badtsos;

#ifdef PROPAGATOR_DB 
  cout<<"FastMuPropagator::propagate::Start propagation in forward::zvert "<<theFmpConst->zvert<<endl;
#endif    
  if(!checkfts(fts)) return badtsos;
   
  // Extract information from muchambers

  int imin0, imax0;
  double ptmax,ptmin,pt,phi,theta,z,r,pl,plmin,plmax;
  double dfcalc,phnext;
  TrackCharge charge;
  float ptboun=theFmpConst->ptboun;
  float step=theFmpConst->step;

  
  GlobalVector moment=fts.parameters().momentum();
  GlobalPoint posit=fts.parameters().position();
  pt=moment.perp();

// Correct theta
          double zfts=fts.parameters().position().z()-theFmpConst->zvert;
          double rfts=fts.parameters().position().perp();
  theta = atan2(rfts,zfts);

#ifdef PROPAGATOR_DB
  cout<<"FastMuPropagator::propagate::Start propagation in forward::theta old, new "<<moment.theta()<<" "<<theta<<endl;
#endif

  phi=posit.phi();
  ptmax=pt+fts.curvilinearError().matrix()(1,1);
  ptmin=pt-fts.curvilinearError().matrix()(1,1);
  pl=pt/tan(theta);
  plmin=ptmin/tan(theta);
  plmax=ptmax/tan(theta);
  imax0=(int)((ptmax-ptboun)/step)+1;
  imin0=(int)((ptmin-ptboun)/step)+1;
#ifdef PROPAGATOR_DB     
  cout<<"FastMuPropagator::Look ptboun,step,imin0,imax0="<<ptboun<<" "<<step<<
    " "<<imin0<<" "<<imax0<<endl;
  cout<<"FastMuPropagator::Parameters (pt,theta,phi,ptmax,ptmin)="<<pt<<" "
      <<theta<<" "<<phi<<" "<<ptmax<<" "<<ptmin<<endl;   
#endif   
  if(imax0 > 1){
    if(imin0<1) imin0=1;
  
    z=boundplane.position().z();
    r=(z-theFmpConst->zvert)*tan(theta);
    charge=fts.parameters().charge();

    // pl evaluation  
  
    dfcalc=theFmpConst->forwparam[0]+
                   charge*theFmpConst->forwparam[1]/abs(pl);  
    phnext=dfcalc+phi;
    if(phnext>=twopi) phnext=phnext-twopi;
    if(phnext <0.) phnext=twopi+phnext;
#ifdef PROPAGATOR_DB    
    cout<<"FastMuPropagator::Phi in Muon Chamb="<<phi<<endl;
    cout<<"FastMuPropagator::Phnext determined="<<phnext<<" dfcalc="<<dfcalc<<endl;
#endif    

    int iwin;
    float awin,bwin,phbound;
    AlgebraicSymMatrix55 m;
      
    if(pt>40.) {
      phbound=theFmpConst->phiwinf[7];
    }else{ // r < bound
      iwin=(int)(pt/theFmpConst->ptstep);
      awin=(theFmpConst->phiwinf[iwin+1]-theFmpConst->phiwinf[iwin])
             /(theFmpConst->ptwmax[iwin]-theFmpConst->ptwmin[iwin]);
      bwin=theFmpConst->phiwinf[iwin]-awin*theFmpConst->ptwmin[iwin];
      phbound=awin*pt+bwin;
    }
#ifdef PROPAGATOR_DB
    cout<<"Forward::Size of window in phi="<<phbound<<endl;
#endif    
    m(0,0)=abs(plmax-plmin)/6.; m(1,1)=phbound/theFmpConst->sigf;
          m(2,2)=theFmpConst->zwin/theFmpConst->sigz;
    m(3,3)=phbound/(2.*theFmpConst->sigf);m(4,4)=theFmpConst->zwin/(2.*theFmpConst->sigz); 
    
    GlobalPoint aX(r*cos(phnext),r*sin(phnext),z);
    
    // Recalculate momentum with some trick...    
    double dfcalcmom=charge*theFmpConst->partrack*abs(z)/abs(pl);
    GlobalVector aP(pt*cos(phnext-dfcalcmom),pt*sin(phnext-dfcalcmom),pl); 
#ifdef PROPAGATOR_DB
    cout<<"dfcalcmom="<<charge<<" "<<theFmpConst->partrack<<" "<<z<<" "<<pl<<" "<<phnext<<endl;
    cout<<"phnextmom="<<phnext-dfcalcmom<<endl;

    cout<<"Plane aP="<<aP<<endl;
    cout<<"Before propagation="<<pt*cos(phi)<< " "<<pt*sin(phi)<<" "<<pl<<endl;
#endif       
       
    GlobalTrajectoryParameters gtp(aX,aP,charge,field);
    CurvilinearTrajectoryError cte(m);
   
    TrajectoryStateOnSurface tsos(gtp, CurvilinearTrajectoryError(m), boundplane);
    return tsos;                                                   
  }
    
  return badtsos;
}
TrajectoryStateOnSurface cms::FastMuPropagator::propagate ( const FreeTrajectoryState fts,
const Cylinder bound 
) const [virtual]

Implements Propagator.

Definition at line 22 of file FastMuPropagator.cc.

References abs, GlobalTrajectoryParameters::charge(), DeDxDiscriminatorTools::charge(), funct::cos(), gather_cfg::cout, FreeTrajectoryState::curvilinearError(), m, CurvilinearTrajectoryError::matrix(), GlobalTrajectoryParameters::momentum(), FreeTrajectoryState::parameters(), PV3DBase< T, PVType, FrameType >::perp(), PV3DBase< T, PVType, FrameType >::phi(), phi, pi, GlobalTrajectoryParameters::position(), ptmin, Cylinder::radius(), funct::sin(), launcher::step, funct::tan(), theta(), PV3DBase< T, PVType, FrameType >::theta(), PV3DBase< T, PVType, FrameType >::z(), and z.

Referenced by propagate().

{
  TrajectoryStateOnSurface badtsos;

  if(!checkfts(fts)) return badtsos;

#ifdef PROPAGATOR_DB 
  cout<<"FastMuPropagator::propagate::Start propagation in barrel::zvert "<<theFmpConst->zvert<<endl;
#endif    

  // Extract information from muchambers

  int charge;

  int imin0, imax0;
  double ptmax,ptmin,pt,phi,theta,z,pl;
  double dfcalc,phnext,bound;
  float ptboun=theFmpConst->ptboun;
  float step=theFmpConst->step;
  
  GlobalVector moment=fts.parameters().momentum();
  GlobalPoint posit=fts.parameters().position();
  pt=moment.perp();

// Correct theta
          double zfts=fts.parameters().position().z()-theFmpConst->zvert;
          double rfts=fts.parameters().position().perp();
  theta = atan2(rfts,zfts);

#ifdef PROPAGATOR_DB
  cout<<"FastMuPropagator::propagate::Start propagation in barrel::theta old, new "<<moment.theta()<<" "<<theta<<endl;
#endif

  phi=posit.phi();
  ptmax=pt+fts.curvilinearError().matrix()(1,1);
  ptmin=pt-fts.curvilinearError().matrix()(1,1);
  bound=boundcyl.radius();
  charge=fts.parameters().charge();

  imax0=(int)((ptmax-ptboun)/step)+1;
  imin0=(int)((ptmin-ptboun)/step)+1;
#ifdef PROPAGATOR_DB     
  cout<<"FastMuPropagator::Look ptboun,step,imin0,imax0="<<ptboun<<" "<<step<<
    " "<<imin0<<" "<<imax0<<endl;
  cout<<"FastMuPropagator::Parameters (pt,theta,phi,ptmax,ptmin,bound)="<<pt<<" "
      <<theta<<" "<<phi<<" "<<ptmax<<" "<<ptmin<<" "<<bound<<endl;       
#endif   
  if(imax0 > 1){
    if(imin0<1) imin0=1;
    // ========================== phi ========================================    
    // new parametrisation (because of second mu-station)
    //

#ifdef PROPAGATOR_DB
    cout<<"FastMuPropagator::New parameters="<<theFmpConst->newparam[0]<<endl;
    cout<<"FastMuPropagator::New parameters for pt>40="<<
           theFmpConst->newparamgt40[0]<<endl;
#endif

    if(pt<40.) {
      dfcalc=charge/(theFmpConst->newparam[0]+theFmpConst->newparam[1]*pt
             +theFmpConst->newparam[2]*pt*pt);       
    }else{
      dfcalc=charge/(theFmpConst->newparamgt40[0]+
      theFmpConst->newparamgt40[1]*pt+theFmpConst->newparamgt40[2]*pt*pt);
    }
//
// !!!! temporarily decision till new parametrization appears.
//    
    if (abs(dfcalc) > 0.6 ) dfcalc = charge*0.6;
    
    phnext=dfcalc+phi;
    if(phnext>=twopi) phnext=phnext-twopi;
    if(phnext <0.) phnext=twopi+phnext;
#ifdef PROPAGATOR_DB    
    cout<<"FastMuPropagator::Phi in Muon Chamb="<<phi<<endl;
    cout<<"FastMuPropagator::Phnext determined="<<phnext<<" dfcalc="<<dfcalc<<endl;
#endif    
    // ==========================  Z  ========================================    
    if(abs(theta-pi/2.)>0.00001){
      z=bound/tan(theta)+theFmpConst->zvert;
      if(fabs(z)>140) {
// Temporary implementation =====      
//         if(z>0.) z = z-45.;
//       if(z<0.) z = z+45.;
// ==============================        
      }
    }else{
      z=0.;
    }
#ifdef PROPAGATOR_DB    
    cout<<"FastMuPropgator::Z coordinate="<<z<<"bound="<<bound<<"theta="<<theta<<" "<<
      abs(theta-pi/2.)<<endl;
#endif
    // ====================== fill global point/vector =======================    
    GlobalPoint aX(bound*cos(phnext),bound*sin(phnext),z);
    if(abs(theta-pi/2.)<0.00001){
      pl=0.;
    }else{
      pl=pt/tan(theta);
    }
    // Recalculate momentum with some trick...
    double dfcalcmom=charge*theFmpConst->partrack*bound/pt;    
    GlobalVector aP(pt*cos(phnext-dfcalcmom),pt*sin(phnext-dfcalcmom),pl);
     
#ifdef PROPAGATOR_DB
    cout<<"phnextmom="<<phnext-dfcalcmom<<endl;
    cout<<"Cylinder aP="<<aP<<endl;
    cout<<"Before propagation="<<pt*cos(phi)<< " "<<pt*sin(phi)<<" "<<pl<<endl;
#endif       
    // =======================================================================    
    
    GlobalTrajectoryParameters gtp(aX,aP,charge,field);
    AlgebraicSymMatrix55 m;
    int iwin;
    float awin,bwin,phbound;
    
    if(pt>40.) {
      iwin=8;
      phbound=theFmpConst->phiwinb[7];
    }else{
      iwin=(int)(pt/theFmpConst->ptstep);
      awin=(theFmpConst->phiwinb[iwin+1]-theFmpConst->phiwinb[iwin])
               /(theFmpConst->ptwmax[iwin]-theFmpConst->ptwmin[iwin]);
      bwin=theFmpConst->phiwinb[iwin]-awin*theFmpConst->ptwmin[iwin];
      phbound=awin*pt+bwin;

    }
#ifdef PROPAGATOR_DB
    cout<<"Size of window in phi="<<phbound<<" "<<iwin<<endl;
    cout<<theFmpConst->phiwinb[iwin+1]<<" "<<theFmpConst->phiwinb[iwin]
            <<" "<<theFmpConst->ptwmin[iwin]<<
      " "<<theFmpConst->ptwmax[iwin]<<" awin="<<awin<<" bwin="<<bwin<<endl;
#endif    
    m(0,0)=(ptmax-ptmin)/6.; m(1,1)=phbound/theFmpConst->sigf;
         m(2,2)=theFmpConst->zwin/theFmpConst->sigz;
    m(3,3)=phbound/(2.*theFmpConst->sigf);m(4,4)=theFmpConst->zwin/(2.*theFmpConst->sigz); 
    
    CurvilinearTrajectoryError cte(m);

    FreeTrajectoryState fts(gtp,cte);
    TrajectoryStateOnSurface tsos(gtp, CurvilinearTrajectoryError(m), boundcyl);
    return tsos;                           
    
  }
  return badtsos;
}
virtual std::pair< TrajectoryStateOnSurface, double> cms::FastMuPropagator::propagateWithPath ( const FreeTrajectoryState state,
const Cylinder bc 
) const [inline, virtual]

Implements Propagator.

Definition at line 55 of file FastMuPropagator.h.

                                                                                {
  std::pair<TrajectoryStateOnSurface,double> tp;
    return tp;
  }
virtual std::pair< TrajectoryStateOnSurface, double> cms::FastMuPropagator::propagateWithPath ( const FreeTrajectoryState state,
const Plane bc 
) const [inline, virtual]

Implements Propagator.

Definition at line 50 of file FastMuPropagator.h.

                                                                              {
  std::pair<TrajectoryStateOnSurface,double> tp;
    return tp;
  }

Member Data Documentation

Definition at line 65 of file FastMuPropagator.h.

Referenced by clone(), FastMuPropagator(), and magneticField().

Definition at line 64 of file FastMuPropagator.h.

Referenced by FastMuPropagator().