CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

NavVolume6Faces Class Reference

#include <NavVolume6Faces.h>

Inheritance diagram for NavVolume6Faces:
NavVolume MagVolume GloballyPositioned< float > MagneticField

List of all members.

Public Member Functions

virtual VolumeCrossReturnType crossToNextVolume (const TrajectoryStateOnSurface &currentState, const Propagator &prop) const
 Cross this volume and point at the next.
virtual const std::vector
< VolumeSide > & 
faces () const
 Access to volume faces.
bool inside (const GlobalPoint &gp, double tolerance) const
bool isIron () const
 Access to Iron/Air information:
 NavVolume6Faces (const MagVolume &magvol, const bool isIron=false)
 A NavVolume6Faces that corresponds exactly to a MagVolume.
 NavVolume6Faces (const PositionType &pos, const RotationType &rot, DDSolidShape shape, const std::vector< NavVolumeSide > &faces, const MagneticFieldProvider< float > *mfp)
virtual Container nextSurface (const NavVolume::LocalPoint &pos, const NavVolume::LocalVector &mom, double charge, PropagationDirection propDir=alongMomentum) const
 Give a sorted list of possible surfaces to propagate to.
virtual Container nextSurface (const NavVolume::LocalPoint &pos, const NavVolume::LocalVector &mom, double charge, PropagationDirection propDir, ConstReferenceCountingPointer< Surface > NotThisSurfaceP) const
 Same, giving lowest priority to the surface we are on now (=NotThisSurface)

Private Member Functions

void computeBounds (const std::vector< NavVolumeSide > &faces)
BoundscomputeBounds (int index, const std::vector< NavVolumeSide > &faces)
BoundscomputeBounds (int index, const std::vector< const Plane * > &bpc)

Private Attributes

bool isThisIron
std::vector< VolumeSidetheFaces
Container theNavSurfaces

Detailed Description

Definition at line 13 of file NavVolume6Faces.h.


Constructor & Destructor Documentation

NavVolume6Faces::NavVolume6Faces ( const PositionType pos,
const RotationType rot,
DDSolidShape  shape,
const std::vector< NavVolumeSide > &  faces,
const MagneticFieldProvider< float > *  mfp 
)

Definition at line 25 of file NavVolume6Faces.cc.

References computeBounds(), i, and theFaces.

                                                                            :
  NavVolume(pos,rot,shape,mfp) 
{
  for (std::vector<NavVolumeSide>::const_iterator i=faces.begin();
       i != faces.end(); i++) {
    theFaces.push_back( VolumeSide( const_cast<Surface*>(&(i->surface().surface())), 
                                    i->globalFace(), i->surfaceSide()));
    //  LogDebug("NavGeometry") << " or actually this is where we have side " << i->surfaceSide() << " and face " << i->globalFace() ;
  }



    computeBounds(faces);
}
NavVolume6Faces::NavVolume6Faces ( const MagVolume magvol,
const bool  isIron = false 
) [explicit]

A NavVolume6Faces that corresponds exactly to a MagVolume.

Definition at line 44 of file NavVolume6Faces.cc.

References computeBounds(), MagVolume::faces(), and i.

                                                                      :
  NavVolume( magvol.position(), magvol.rotation(), magvol.shapeType(), magvol.provider()),
  theFaces(magvol.faces()), isThisIron(isIron)
{
  std::vector<NavVolumeSide> navSides;
  std::vector<VolumeSide> magSides( magvol.faces());
  NavSurfaceBuilder navBuilder;

  for (std::vector<VolumeSide>::const_iterator i=magSides.begin();
       i != magSides.end(); i++) {
    NavSurface* navSurface = navBuilder.build( i->surface());
    navSides.push_back( NavVolumeSide( navSurface, i->globalFace(), i->surfaceSide()));
  }
  computeBounds(navSides);
}

Member Function Documentation

void NavVolume6Faces::computeBounds ( const std::vector< NavVolumeSide > &  faces) [private]

Definition at line 73 of file NavVolume6Faces.cc.

References NavSurface::addVolume(), NavSurface::bounds(), i, and theNavSurfaces.

Referenced by NavVolume6Faces().

{
    bool allPlanes = true;
    // bool allPlanes = false; // for TESTS ONLY!!!
    std::vector<const Plane*> planes;
    for (std::vector<NavVolumeSide>::const_iterator iface=faces.begin(); iface!=faces.end(); iface++) {
        const Plane* plane = dynamic_cast<const Plane*>(&(iface->surface()));
        if (plane != 0) {
            planes.push_back(plane);
        }
        else allPlanes = false;
    }
    
    for (unsigned int i=0; i<faces.size(); i++) {

      // FIXME: who owns the new NavSurface? memory leak???

        NavSurface& navSurf = faces[i].mutableSurface();
        Bounds* myBounds = 0;
        if (allPlanes) {        
            myBounds = computeBounds( i, planes);
        }
        else {
            myBounds = computeBounds( i, faces);
        }
        navSurf.addVolume( this, myBounds, faces[i].surfaceSide());
        delete myBounds; // since navSurf now owns a copy

// this is tricky: we want to avoid multiple copies of the Bounds; the NavSurface owns
// a copy of Bounds for each touching Volume (instantiated in the call to addVolume).
// We would like to keep a pointer to the same Bounds in the NavVolume, so we have to ASK
// the NavSurface for the Bounds* of the Bounds we just gave it!
        //LogDebug("NavGeometry") << "Adding a Volume Side with center " << navSurf.surface().position() << " side "<< faces[i].surfaceSide() << " and face " << faces[i].globalFace();
        theNavSurfaces.push_back( SurfaceAndBounds(&navSurf, navSurf.bounds(this), faces[i].surfaceSide(), faces[i].globalFace()));
    }
}
Bounds * NavVolume6Faces::computeBounds ( int  index,
const std::vector< const Plane * > &  bpc 
) [private]

Definition at line 110 of file NavVolume6Faces.cc.

References gather_cfg::cout, ThreePlaneCrossing::crossing(), ExpressReco_HICollisions_FallBack::e, i, j, Plane::localZ(), Plane::normalVector(), SurfaceOrientation::onSurface, GloballyPositioned< T >::position(), Plane::side(), and GloballyPositioned< T >::toLocal().

{
  const Plane* plane( bpc[index]);

  // find the 4 intersecting planes
  int startIndex = 2*(1+index/2); // 2, 4, 6
  std::vector<const Plane*> crossed; crossed.reserve(4);
  for (int j = startIndex; j <  startIndex+4; j++) {
    crossed.push_back(bpc[j%6]);
  }

  // compute intersection corners of the plane triplets
  std::vector<GlobalPoint> corners; corners.reserve(4);
  ThreePlaneCrossing crossing;
  for ( int i=0; i<2; i++) {
    for ( int j=2; j<4; j++) {
      GlobalPoint corner( crossing.crossing( *plane, *crossed[i], *crossed[j]));
      corners.push_back(corner);

#ifdef DEBUG
      std::cout << "Crossing of planes is " << corner << std::endl;
      std::cout << "NormalVectors of the planes are " << plane->normalVector()
           << " " << crossed[i]->normalVector() << " " << crossed[j]->normalVector() << std::endl;
      std::cout << "Positions of planes are " << plane->position()
           << " " << crossed[i]->position() << " " << crossed[j]->position() << std::endl;
      if (plane->side( corner, 1.e-5) == SurfaceOrientation::onSurface &&
          crossed[i]->side( corner, 1.e-5) == SurfaceOrientation::onSurface &&
          crossed[j]->side( corner, 1.e-5) == SurfaceOrientation::onSurface) {
          std::cout << "Crossing is really on all three surfaces" << std::endl;
      }
      else {
          std::cout << "CROSSING IS NOT ON SURFACES!!!" << std::endl;
          std::cout << plane->localZ(corner) << std::endl;
          std::cout << crossed[i]->localZ(corner) << std::endl;
          std::cout << crossed[j]->localZ(corner) << std::endl;
       }
#endif

    }
  }

  // put corners in cyclic sequence (2 and 3 swapped)
  return new FourPointPlaneBounds( plane->toLocal( corners[0]), plane->toLocal( corners[1]),
                                   plane->toLocal( corners[3]), plane->toLocal( corners[2]));
}
Bounds * NavVolume6Faces::computeBounds ( int  index,
const std::vector< NavVolumeSide > &  faces 
) [private]

Definition at line 36 of file BladeShapeBuilderFromDet.cc.

References BoundingBox::corners(), funct::cos(), i, LogDebug, max(), min, phi, Geom::pi(), pi, csvReporter::r, funct::sin(), Surface::toGlobal(), GloballyPositioned< T >::toLocal(), z, and zPos.

{
  Surface::PositionType tmpPos = dets.front()->surface().position();


  float rmin(plane.toLocal(tmpPos).perp());
  float rmax(plane.toLocal(tmpPos).perp());
  float zmin(plane.toLocal(tmpPos).z());
  float zmax(plane.toLocal(tmpPos).z());
  float phimin(plane.toLocal(tmpPos).phi());
  float phimax(plane.toLocal(tmpPos).phi());

  for(vector<const GeomDet*>::const_iterator it=dets.begin(); it!=dets.end(); it++){
    vector<GlobalPoint> corners = BoundingBox().corners( (*it)->specificSurface() );

    for(vector<GlobalPoint>::const_iterator i=corners.begin();
        i!=corners.end(); i++) {

      float r   = plane.toLocal(*i).perp();
      float z   = plane.toLocal(*i).z();
      float phi = plane.toLocal(*i).phi();
      rmin = min( rmin, r);
      rmax = max( rmax, r);
      zmin = min( zmin, z);
      zmax = max( zmax, z);
      if ( PhiLess()( phi, phimin)) phimin = phi;
      if ( PhiLess()( phimax, phi)) phimax = phi;
    }
    // in addition to the corners we have to check the middle of the 
    // det +/- length/2, since the min (max) radius for typical fw
    // dets is reached there
        
    float rdet = (*it)->position().perp();
    float height  = (*it)->surface().bounds().width();
    rmin = min( rmin, rdet-height/2.F);
    rmax = max( rmax, rdet+height/2.F);  
    

  }

  if (!PhiLess()(phimin, phimax)) 
    edm::LogError("TkDetLayers") << " BladeShapeBuilderFromDet : " 
                                 << "Something went wrong with Phi Sorting !" ;
  float zPos = (zmax+zmin)/2.;
  float phiWin = phimax - phimin;
  float phiPos = (phimax+phimin)/2.;
  float rmed = (rmin+rmax)/2.;
  if ( phiWin < 0. ) {
    if ( (phimin < Geom::pi() / 2.) || (phimax > -Geom::pi()/2.) ){
      edm::LogError("TkDetLayers") << " something strange going on, please check " ;
    }
    //edm::LogInfo(TkDetLayers) << " Wedge at pi: phi " << phimin << " " << phimax << " " << phiWin 
    //   << " " << 2.*Geom::pi()+phiWin << " " ;
    phiWin += 2.*Geom::pi();
    phiPos += Geom::pi(); 
  }
  
  LocalVector localPos( rmed*cos(phiPos), rmed*sin(phiPos), zPos);

  LogDebug("TkDetLayers") << "localPos in computeBounds: " << localPos << "\n"
                          << "rmin:   " << rmin << "\n"
                          << "rmax:   " << rmax << "\n"
                          << "zmin:   " << zmin << "\n"
                          << "zmax:   " << zmax << "\n"
                          << "phiWin: " << phiWin ;

  return make_pair(DiskSectorBounds(rmin,rmax,zmin,zmax,phiWin),
                   plane.toGlobal(localPos) );

}
VolumeCrossReturnType NavVolume6Faces::crossToNextVolume ( const TrajectoryStateOnSurface currentState,
const Propagator prop 
) const [virtual]

Cross this volume and point at the next.

Implements NavVolume.

Definition at line 250 of file NavVolume6Faces.cc.

References alongMomentum, except, TrajectoryStateOnSurface::globalMomentum(), TrajectoryStateOnSurface::globalPosition(), LogDebug, nextSurface(), oppositeSide(), evf::utils::state, TrajectoryStateOnSurface::surface(), and GloballyPositioned< float >::toLocal().

{
  typedef TrajectoryStateOnSurface TSOS;
  typedef std::pair<TSOS,double> TSOSwithPath;

  NavVolume::Container nsc = nextSurface( toLocal( startingState.globalPosition()), 
                                          toLocal( startingState.globalMomentum()), -1,
                                          alongMomentum, &(startingState.surface()));
  int itry = 0;
  VolumeCrossReturnType VolumeCrossResult( 0, startingState, 0.0);

  for (NavVolume::Container::const_iterator isur = nsc.begin(); isur!=nsc.end(); isur++) {

    LogDebug("NavGeometry") <<  "crossToNextVolume: trying Surface no. " << itry ;
    TSOSwithPath state;
    
    try {
      state = isur->surface().propagateWithPath( prop, startingState);
    }
    catch (MagVolumeOutsideValidity& except) {
      LogDebug("NavGeometry") << "Ohoh... failed to stay inside magnetic field !! skip this surface " ;
        ++itry;
      continue;
    }
    
    if (!state.first.isValid()) {
      ++itry;
      continue;
    }

    LogDebug("NavGeometry") <<  "crossToNextVolume: reached Valid State at Surface no. " << itry ;
    LogDebug("NavGeometry") << " --> local position of Valid state is " <<  state.first.localPosition()  ;
    LogDebug("NavGeometry") << " --> global position of Valid state is " <<  state.first.globalPosition()  ;

    if (isur->bounds().inside(state.first.localPosition())) {
      //LogDebug("NavGeometry") << "crossToNextVolume: Surface containing destination point found at try " << itry ;
      // Found the exit surface !! Get pointer to next volume and save exit state:
      //VolumeCrossResult.first = isur->surface().nextVolume(state.localPosition(),oppositeSide(isur->side()));
      //VolumeCrossResult.second = state;
      //      exitSurface = &( isur->surface().surface() );
      //if(VolumeCrossResult.path() < 0.01) {
      //        LogDebug("NavGeometry") << " Stuck at  " << state.first.globalPosition() ;
      //}
      return VolumeCrossReturnType ( isur->surface().nextVolume(state.first.localPosition(), oppositeSide(isur->side())),
                            state.first, state.second );
      
      break;
    }
    else {
      LogDebug("NavGeometry") << "crossToNextVolume: BUT not inside the Bounds !! " ;
      ++itry;
    }
  }

  return VolumeCrossResult;
}
virtual const std::vector<VolumeSide>& NavVolume6Faces::faces ( ) const [inline, virtual]

Access to volume faces.

Implements MagVolume.

Definition at line 39 of file NavVolume6Faces.h.

References theFaces.

{ return theFaces; }
bool NavVolume6Faces::inside ( const GlobalPoint gp,
double  tolerance 
) const [virtual]

Implements MagVolume.

Definition at line 61 of file NavVolume6Faces.cc.

References i, SurfaceOrientation::onSurface, and theFaces.

{
  // check if the point is on the correct side of all delimiting surfaces
  for (std::vector<VolumeSide>::const_iterator i=theFaces.begin(); i!=theFaces.end(); i++) {
    Surface::Side side = i->surface().side( gp, tolerance);
    if ( side != i->surfaceSide() && side != SurfaceOrientation::onSurface) return false;
  }
  return true;
}
bool NavVolume6Faces::isIron ( ) const [inline, virtual]

Access to Iron/Air information:

Implements NavVolume.

Definition at line 42 of file NavVolume6Faces.h.

References isThisIron.

{ return isThisIron; }
NavVolume::Container NavVolume6Faces::nextSurface ( const NavVolume::LocalPoint pos,
const NavVolume::LocalVector mom,
double  charge,
PropagationDirection  propDir = alongMomentum 
) const [virtual]

Give a sorted list of possible surfaces to propagate to.

Implements NavVolume.

Definition at line 174 of file NavVolume6Faces.cc.

References alongMomentum, epsilon, i, query::result, theNavSurfaces, and GloballyPositioned< float >::toGlobal().

Referenced by crossToNextVolume().

{
    typedef std::map<double,SurfaceAndBounds>   SortedContainer;

    GlobalPoint  gpos( toGlobal(pos));
    GlobalVector gmom( toGlobal(mom));
    GlobalVector gdir = (propDir == alongMomentum ? gmom : -gmom);

    SortedContainer sortedSurfaces;
    Container       verycloseSurfaces; // reachable surface with dist < epsilon !!
    Container       unreachableSurfaces;

    double epsilon = 0.01; // should epsilon be hard-coded or a variable in NavVolume?

    for (Container::const_iterator i=theNavSurfaces.begin(); i!=theNavSurfaces.end(); i++) {
        std::pair<bool,double> dist = i->surface().distanceAlongLine( gpos, gdir);
        if (dist.first) {
          if (dist.second > epsilon) sortedSurfaces[dist.second] = *i;
          else verycloseSurfaces.push_back(*i);
        } 
        else unreachableSurfaces.push_back(*i);
    }
    NavVolume::Container result;
    for (SortedContainer::const_iterator i=sortedSurfaces.begin(); i!=sortedSurfaces.end(); ++i) {
        result.push_back(i->second);
    }
    result.insert( result.end(), unreachableSurfaces.begin(), unreachableSurfaces.end());
    result.insert( result.end(), verycloseSurfaces.begin(), verycloseSurfaces.end());
    return result;
}
NavVolume::Container NavVolume6Faces::nextSurface ( const NavVolume::LocalPoint pos,
const NavVolume::LocalVector mom,
double  charge,
PropagationDirection  propDir,
ConstReferenceCountingPointer< Surface NotThisSurfaceP 
) const [virtual]

Same, giving lowest priority to the surface we are on now (=NotThisSurface)

Implements NavVolume.

Definition at line 208 of file NavVolume6Faces.cc.

References alongMomentum, epsilon, i, query::result, theNavSurfaces, and GloballyPositioned< float >::toGlobal().

{
    typedef std::map<double,SurfaceAndBounds>   SortedContainer;

    GlobalPoint  gpos( toGlobal(pos));
    GlobalVector gmom( toGlobal(mom));
    GlobalVector gdir = (propDir == alongMomentum ? gmom : -gmom);

    SortedContainer sortedSurfaces;
    Container       verycloseSurfaces; // reachable surface with dist < epsilon (if 6-surface check fails)
    Container       unreachableSurfaces;

    double epsilon = 0.01; // should epsilon be hard-coded or a variable in NavVolume?
    bool surfaceMatched = false;

    for (Container::const_iterator i=theNavSurfaces.begin(); i!=theNavSurfaces.end(); i++) {
      if (&(i->surface().surface()) == NotThisSurfaceP) surfaceMatched = true;
    }
    
    for (Container::const_iterator i=theNavSurfaces.begin(); i!=theNavSurfaces.end(); i++) {
        std::pair<bool,double> dist = i->surface().distanceAlongLine( gpos, gdir);
        if (dist.first) { 
          if ( &(i->surface().surface()) == NotThisSurfaceP || (!surfaceMatched && dist.second < epsilon)) 
            verycloseSurfaces.push_back(*i);
          else sortedSurfaces[dist.second] = *i;
        }
        else unreachableSurfaces.push_back(*i);
    }

    NavVolume::Container result;
    for (SortedContainer::const_iterator i=sortedSurfaces.begin(); i!=sortedSurfaces.end(); ++i) {
        result.push_back(i->second);
    }
    result.insert( result.end(), unreachableSurfaces.begin(), unreachableSurfaces.end());
    result.insert( result.end(), verycloseSurfaces.begin(), verycloseSurfaces.end());
    return result;
}

Member Data Documentation

Definition at line 51 of file NavVolume6Faces.h.

Referenced by isIron().

std::vector<VolumeSide> NavVolume6Faces::theFaces [private]

Definition at line 49 of file NavVolume6Faces.h.

Referenced by faces(), inside(), and NavVolume6Faces().

Definition at line 50 of file NavVolume6Faces.h.

Referenced by computeBounds(), and nextSurface().