CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/RecoLocalMuon/CSCSegment/src/CSCSegAlgoShowering.cc

Go to the documentation of this file.
00001 
00010 #include "RecoLocalMuon/CSCSegment/src/CSCSegAlgoShowering.h"
00011 
00012 #include "Geometry/CSCGeometry/interface/CSCLayer.h"
00013 #include <Geometry/CSCGeometry/interface/CSCChamber.h>
00014 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00015 
00016 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00017 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00018 
00019 #include <algorithm>
00020 #include <cmath>
00021 #include <iostream>
00022 #include <string>
00023 
00024 
00025 /* Constructor
00026  *
00027  */
00028 CSCSegAlgoShowering::CSCSegAlgoShowering(const edm::ParameterSet& ps) {
00029 //  debug                  = ps.getUntrackedParameter<bool>("CSCSegmentDebug");
00030   dRPhiFineMax           = ps.getParameter<double>("dRPhiFineMax");
00031   dPhiFineMax            = ps.getParameter<double>("dPhiFineMax");
00032   tanThetaMax            = ps.getParameter<double>("tanThetaMax");
00033   tanPhiMax              = ps.getParameter<double>("tanPhiMax");        
00034   maxRatioResidual       = ps.getParameter<double>("maxRatioResidualPrune");
00035 //  maxDR                  = ps.getParameter<double>("maxDR");
00036   maxDTheta              = ps.getParameter<double>("maxDTheta");
00037   maxDPhi                = ps.getParameter<double>("maxDPhi");
00038 }
00039 
00040 
00041 /* Destructor:
00042  *
00043  */
00044 CSCSegAlgoShowering::~CSCSegAlgoShowering(){
00045 
00046 }
00047 
00048 
00049 /* showerSeg
00050  *
00051  */
00052 CSCSegment CSCSegAlgoShowering::showerSeg( const CSCChamber* aChamber, ChamberHitContainer rechits ) {
00053 
00054   theChamber = aChamber;
00055   // Initialize parameters
00056   std::vector<float> x, y, gz;
00057   std::vector<int> n;
00058   
00059   
00060   for (int i = 0; i < 6; ++i) {
00061     x.push_back(0.);
00062     y.push_back(0.);
00063     gz.push_back(0.);
00064     n.push_back(0);
00065   }
00066 
00067   // Loop over hits to find center-of-mass position in each layer
00068   for (ChamberHitContainer::const_iterator it = rechits.begin(); it != rechits.end(); it++ ) {
00069     const CSCRecHit2D& hit = (**it);
00070     const CSCDetId id = hit.cscDetId();
00071     int l_id = id.layer();
00072     const CSCLayer* layer  = theChamber->layer(hit.cscDetId().layer());
00073     GlobalPoint gp         = layer->toGlobal(hit.localPosition());
00074     LocalPoint  lp         = theChamber->toLocal(gp);
00075 
00076     n[l_id -1]++;
00077     x[l_id -1] += lp.x();
00078     y[l_id -1] += lp.y();
00079     gz[l_id -1] += gp.z();
00080   }
00081 
00082 
00083   // Determine center of mass for each layer and average center of mass for chamber
00084   float avgChamberX = 0.;
00085   float avgChamberY = 0.;
00086   int n_lay = 0;
00087 
00088   for (unsigned i = 0; i < 6; ++i) {
00089     if (n[i] < 1 ) continue;
00090  
00091     x[i] = x[i]/n[i];
00092     y[i] = y[i]/n[i];
00093     avgChamberX += x[i];
00094     avgChamberY += y[i];
00095     n_lay++;
00096 
00097   }
00098 
00099   if ( n_lay > 0) {
00100     avgChamberX = avgChamberX / n_lay;
00101     avgChamberY = avgChamberY / n_lay;
00102   }
00103 
00104   // Create a FakeSegment origin that points back to the IP
00105   // Keep all this in global coordinates until last minute to avoid screwing up +/- signs !
00106 
00107   LocalPoint   lpCOM(avgChamberX, avgChamberY, 0.);
00108   GlobalPoint  gpCOM = theChamber->toGlobal(lpCOM);
00109   GlobalVector gvCOM(gpCOM.x(), gpCOM.y(), gpCOM.z());
00110 
00111   float Gdxdz = gpCOM.x()/gpCOM.z();
00112   float Gdydz = gpCOM.y()/gpCOM.z();
00113 
00114   // Figure out the intersection of this vector with each layer of the chamber
00115   // by projecting the vector
00116   std::vector<LocalPoint> layerPoints;
00117 
00118   for (unsigned i = 1; i < 7; i++) {
00119     // Get the layer z coordinates in global frame
00120     const CSCLayer* layer = theChamber->layer(i);
00121     LocalPoint temp(0., 0., 0.);
00122     GlobalPoint gp = layer->toGlobal(temp);
00123     float layer_Z = gp.z();
00124 
00125     // Then compute interesection of vector with that plane
00126     float layer_X = Gdxdz * layer_Z;
00127     float layer_Y = Gdydz * layer_Z;
00128     GlobalPoint Gintersect(layer_X, layer_Y, layer_Z);
00129     LocalPoint  Lintersect = theChamber->toLocal(Gintersect);
00130 
00131     float layerX = Lintersect.x();
00132     float layerY = Lintersect.y();
00133     float layerZ = Lintersect.z();
00134     LocalPoint layerPoint(layerX, layerY, layerZ);
00135     layerPoints.push_back(layerPoint);
00136   }
00137 
00138 
00139   std::vector<float> r_closest;
00140   std::vector<int> id;
00141   for (unsigned i = 0; i < 6; ++i ) {
00142     id.push_back(-1);
00143     r_closest.push_back(9999.);
00144   }
00145 
00146   int idx = 0;
00147 
00148   // Loop over all hits and find hit closest to com for that layer.
00149   for (ChamberHitContainer::const_iterator it = rechits.begin(); it != rechits.end(); it++ ) {    
00150     const CSCRecHit2D& hit = (**it);
00151     int layId = hit.cscDetId().layer();
00152     const CSCLayer* layer  = theChamber->layer(hit.cscDetId().layer());
00153     GlobalPoint gp         = layer->toGlobal(hit.localPosition());
00154     LocalPoint  lp         = theChamber->toLocal(gp);
00155 
00156     float d_x = lp.x() - layerPoints[layId-1].x();
00157     float d_y = lp.y() - layerPoints[layId-1].y();
00158 
00159     LocalPoint diff(d_x, d_y, 0.);
00160     
00161     if ( fabs(diff.mag() ) < r_closest[layId-1] ) {
00162        r_closest[layId-1] =  fabs(diff.mag());
00163        id[layId-1] = idx;
00164     }
00165     idx++;
00166   }
00167 
00168   // Now fill vector of rechits closest to center of mass:
00169   protoSegment.clear();
00170   idx = 0;
00171 
00172   // Loop over all hits and find hit closest to com for that layer.
00173   for (ChamberHitContainer::const_iterator it = rechits.begin(); it != rechits.end(); it++ ) {    
00174     const CSCRecHit2D& hit = (**it);
00175     int layId = hit.cscDetId().layer();
00176 
00177     if ( idx == id[layId-1] )protoSegment.push_back(*it);
00178 
00179     idx++;    
00180   }
00181 
00182   // Reorder hits in protosegment
00183   if ( gz[0] > 0. ) {
00184     if ( gz[0] > gz[5] ) { 
00185       reverse( protoSegment.begin(), protoSegment.end() );
00186     }    
00187   }
00188   else if ( gz[0] < 0. ) {
00189     if ( gz[0] < gz[5] ) {
00190       reverse( protoSegment.begin(), protoSegment.end() );
00191     }    
00192   }
00193 
00194 
00195   // Compute the segment properties
00196   updateParameters();
00197 
00198   // Clean up protosegment if there is one very bad hit on segment
00199   if (protoSegment.size() > 4) pruneFromResidual();
00200 
00201   // Look for better hits near segment  
00202   for (ChamberHitContainer::const_iterator it = rechits.begin(); it != rechits.end(); it++ ) {
00203 
00204     const CSCRecHit2D* h = *it;
00205     int layer = h->cscDetId().layer();
00206 
00207     if ( isHitNearSegment( h ) ) compareProtoSegment(h, layer);
00208   }
00209 
00210   // Prune worse hit if necessary
00211   if ( protoSegment.size() > 5 ) pruneFromResidual();
00212 
00213   // Update the parameters
00214   updateParameters();
00215 
00216   // Local direction
00217   double dz   = 1./sqrt(1. + protoSlope_u*protoSlope_u + protoSlope_v*protoSlope_v);
00218   double dx   = dz*protoSlope_u;
00219   double dy   = dz*protoSlope_v;
00220   LocalVector localDir(dx,dy,dz);
00221         
00222   // localDir may need sign flip to ensure it points outward from IP  
00223   double globalZpos    = ( theChamber->toGlobal( protoIntercept ) ).z();
00224   double globalZdir    = ( theChamber->toGlobal( localDir ) ).z();
00225   double directionSign = globalZpos * globalZdir;
00226   LocalVector protoDirection = (directionSign * localDir).unit();
00227 
00228   // Check to make sure the fitting didn't fuck things up
00229   GlobalVector protoGlobalDir = theChamber->toGlobal( protoDirection );  
00230   double protoTheta = protoGlobalDir.theta();
00231   double protoPhi = protoGlobalDir.phi();
00232   double simTheta = gvCOM.theta();
00233   double simPhi = gvCOM.phi();
00234   
00235   float dTheta = fabs(protoTheta - simTheta);
00236   float dPhi   = fabs(protoPhi - simPhi);
00237   //  float dR = sqrt(dEta*dEta + dPhi*dPhi);
00238   
00239   // Error matrix
00240   AlgebraicSymMatrix protoErrors = calculateError();     
00241   // but reorder components to match what's required by TrackingRecHit interface
00242   // i.e. slopes first, then positions
00243   flipErrors( protoErrors );
00244   
00245   //Flag the segment with chi12 of -1 of the segment isn't pointing toward origin      
00246   if (dTheta > maxDTheta || dPhi > maxDPhi) {
00247   CSCSegment temp(protoSegment, protoIntercept, protoDirection, protoErrors, -1); 
00248   return temp;
00249   }
00250   else {
00251   CSCSegment temp(protoSegment, protoIntercept, protoDirection, protoErrors, protoChi2);
00252   return temp;
00253   }
00254 
00255 } 
00256 
00257 
00258 
00259 
00260 /* isHitNearSegment
00261  *
00262  * Compare rechit with expected position from proto_segment
00263  */
00264 bool CSCSegAlgoShowering::isHitNearSegment( const CSCRecHit2D* hit) const {
00265 
00266   const CSCLayer* layer = theChamber->layer(hit->cscDetId().layer());
00267 
00268   // hit phi position in global coordinates
00269   GlobalPoint Hgp = layer->toGlobal(hit->localPosition());
00270   double Hphi = Hgp.phi();                                
00271   if (Hphi < 0.) Hphi += 2.*M_PI;
00272   LocalPoint Hlp = theChamber->toLocal(Hgp);
00273   double z = Hlp.z();  
00274 
00275   double LocalX = protoIntercept.x() + protoSlope_u * z;
00276   double LocalY = protoIntercept.y() + protoSlope_v * z;
00277   LocalPoint Slp(LocalX, LocalY, z);
00278   GlobalPoint Sgp = theChamber->toGlobal(Slp); 
00279   double Sphi = Sgp.phi();
00280   if (Sphi < 0.) Sphi += 2.*M_PI;
00281   double R = sqrt(Sgp.x()*Sgp.x() + Sgp.y()*Sgp.y());
00282   
00283   double deltaPhi = Sphi - Hphi;
00284   if (deltaPhi >  2.*M_PI) deltaPhi -= 2.*M_PI;
00285   if (deltaPhi < -2.*M_PI) deltaPhi += 2.*M_PI;
00286   if (deltaPhi < 0.) deltaPhi = -deltaPhi; 
00287 
00288   double RdeltaPhi = R * deltaPhi;
00289 
00290   if (RdeltaPhi < dRPhiFineMax && deltaPhi < dPhiFineMax ) return true;
00291 
00292   return false;
00293 }
00294 
00295 
00296 /* Method addHit
00297  *
00298  * Test if can add hit to proto segment. If so, try to add it.
00299  *
00300  */
00301 bool CSCSegAlgoShowering::addHit(const CSCRecHit2D* aHit, int layer) {
00302   
00303   // Return true if hit was added successfully and then parameters are updated.
00304   // Return false if there is already a hit on the same layer, or insert failed.
00305   
00306   bool ok = true;
00307   
00308   // Test that we are not trying to add the same hit again
00309   for ( ChamberHitContainer::const_iterator it = protoSegment.begin(); it != protoSegment.end(); it++ ) 
00310     if ( aHit == (*it)  ) return false;
00311   
00312   protoSegment.push_back(aHit);
00313 
00314   return ok;
00315 }    
00316 
00317 
00318 /* Method updateParameters
00319  *      
00320  * Perform a simple Least Square Fit on proto segment to determine slope and intercept
00321  *
00322  */   
00323 void CSCSegAlgoShowering::updateParameters() {
00324 
00325   // Compute slope from Least Square Fit    
00326   CLHEP::HepMatrix M(4,4,0);
00327   CLHEP::HepVector B(4,0);
00328 
00329   ChamberHitContainer::const_iterator ih;
00330   
00331   for (ih = protoSegment.begin(); ih != protoSegment.end(); ++ih) {
00332     
00333     const CSCRecHit2D& hit = (**ih);
00334     const CSCLayer* layer  = theChamber->layer(hit.cscDetId().layer());
00335     GlobalPoint gp         = layer->toGlobal(hit.localPosition());
00336     LocalPoint  lp         = theChamber->toLocal(gp); 
00337     
00338     double u = lp.x();
00339     double v = lp.y();
00340     double z = lp.z();
00341     
00342     // ptc: Covariance matrix of local errors 
00343     CLHEP::HepMatrix IC(2,2);
00344     IC(1,1) = hit.localPositionError().xx();
00345     IC(1,2) = hit.localPositionError().xy();
00346     IC(2,2) = hit.localPositionError().yy();
00347     IC(2,1) = IC(1,2); // since Cov is symmetric
00348     
00349     // ptc: Invert covariance matrix (and trap if it fails!)
00350     int ierr = 0;
00351     IC.invert(ierr); // inverts in place
00352     if (ierr != 0) {
00353       LogDebug("CSC") << "CSCSegment::fitSlopes: failed to invert covariance matrix=\n" << IC << "\n";      
00354     }
00355     
00356     M(1,1) += IC(1,1);
00357     M(1,2) += IC(1,2);
00358     M(1,3) += IC(1,1) * z;
00359     M(1,4) += IC(1,2) * z;
00360     B(1)   += u * IC(1,1) + v * IC(1,2);
00361     
00362     M(2,1) += IC(2,1);
00363     M(2,2) += IC(2,2);
00364     M(2,3) += IC(2,1) * z;
00365     M(2,4) += IC(2,2) * z;
00366     B(2)   += u * IC(2,1) + v * IC(2,2);
00367     
00368     M(3,1) += IC(1,1) * z;
00369     M(3,2) += IC(1,2) * z;
00370     M(3,3) += IC(1,1) * z * z;
00371     M(3,4) += IC(1,2) * z * z;
00372     B(3)   += ( u * IC(1,1) + v * IC(1,2) ) * z;
00373     
00374     M(4,1) += IC(2,1) * z;
00375     M(4,2) += IC(2,2) * z;
00376     M(4,3) += IC(2,1) * z * z;
00377     M(4,4) += IC(2,2) * z * z;
00378     B(4)   += ( u * IC(2,1) + v * IC(2,2) ) * z;
00379   }
00380   
00381   CLHEP::HepVector p = solve(M, B);
00382   
00383   // Update member variables 
00384   // Note that origin has local z = 0
00385 
00386   protoIntercept = LocalPoint(p(1), p(2), 0.);
00387   protoSlope_u = p(3);
00388   protoSlope_v = p(4);
00389 
00390   // Determine Chi^2 for the proto wire segment
00391   
00392   double chsq = 0.;
00393   
00394   for (ih = protoSegment.begin(); ih != protoSegment.end(); ++ih) {
00395     
00396     const CSCRecHit2D& hit = (**ih);
00397     const CSCLayer* layer  = theChamber->layer(hit.cscDetId().layer());
00398     GlobalPoint gp         = layer->toGlobal(hit.localPosition());
00399     LocalPoint lp          = theChamber->toLocal(gp);
00400     
00401     double u = lp.x();
00402     double v = lp.y();
00403     double z = lp.z();
00404     
00405     double du = protoIntercept.x() + protoSlope_u * z - u;
00406     double dv = protoIntercept.y() + protoSlope_v * z - v;
00407     
00408     CLHEP::HepMatrix IC(2,2);
00409     IC(1,1) = hit.localPositionError().xx();
00410     IC(1,2) = hit.localPositionError().xy();
00411     IC(2,2) = hit.localPositionError().yy();
00412     IC(2,1) = IC(1,2);
00413     
00414     // Invert covariance matrix
00415     int ierr = 0;
00416     IC.invert(ierr);
00417     if (ierr != 0) {
00418       LogDebug("CSC") << "CSCSegment::fillChiSquared: failed to invert covariance matrix=\n" << IC << "\n";      
00419     }
00420     chsq += du*du*IC(1,1) + 2.*du*dv*IC(1,2) + dv*dv*IC(2,2);
00421   }
00422   protoChi2 = chsq;
00423 }
00424 
00425 
00426 
00427 /* Method compareProtoSegment
00428  *      
00429  * For hit coming from the same layer of an existing hit within the proto segment
00430  * test if achieve better chi^2 by using this hit than the other
00431  *
00432  */ 
00433 void CSCSegAlgoShowering::compareProtoSegment(const CSCRecHit2D* h, int layer) {
00434   
00435   // Store old segment first
00436   double old_protoChi2                  = protoChi2;
00437   LocalPoint old_protoIntercept         = protoIntercept;
00438   float old_protoSlope_u                = protoSlope_u;
00439   float old_protoSlope_v                = protoSlope_v;
00440   LocalVector old_protoDirection        = protoDirection;
00441   ChamberHitContainer old_protoSegment  = protoSegment;
00442  
00443 
00444   // Try adding the hit to existing segment, and remove old one existing in same layer
00445   ChamberHitContainer::iterator it;
00446   for ( it = protoSegment.begin(); it != protoSegment.end(); ) {
00447     if ( (*it)->cscDetId().layer() == layer ) {
00448       it = protoSegment.erase(it);
00449     } else {
00450       ++it;
00451     }
00452   }
00453   bool ok = addHit(h, layer);
00454 
00455   if (ok) updateParameters();
00456   
00457   if ( (protoChi2 > old_protoChi2) || ( !ok ) ) {
00458     protoChi2       = old_protoChi2;
00459     protoIntercept  = old_protoIntercept;
00460     protoSlope_u    = old_protoSlope_u;
00461     protoSlope_v    = old_protoSlope_v;
00462     protoDirection  = old_protoDirection;
00463     protoSegment    = old_protoSegment;
00464   }
00465 }
00466 
00467 AlgebraicSymMatrix CSCSegAlgoShowering::weightMatrix() const {
00468 
00469   std::vector<const CSCRecHit2D*>::const_iterator it;
00470   int nhits = protoSegment.size();
00471   AlgebraicSymMatrix matrix(2*nhits, 0);
00472   int row = 0;
00473 
00474   for (it = protoSegment.begin(); it != protoSegment.end(); ++it) {
00475 
00476     const CSCRecHit2D& hit = (**it);
00477     ++row;
00478     matrix(row, row)   = hit.localPositionError().xx();
00479     matrix(row, row+1) = hit.localPositionError().xy();
00480     ++row;
00481     matrix(row, row-1) = hit.localPositionError().xy();
00482     matrix(row, row)   = hit.localPositionError().yy();
00483   }
00484   int ierr;
00485   matrix.invert(ierr);
00486   return matrix;
00487 }
00488 
00489 
00490 CLHEP::HepMatrix CSCSegAlgoShowering::derivativeMatrix() const {
00491 
00492   ChamberHitContainer::const_iterator it;
00493   int nhits = protoSegment.size();
00494   CLHEP::HepMatrix matrix(2*nhits, 4);
00495   int row = 0;
00496 
00497   for(it = protoSegment.begin(); it != protoSegment.end(); ++it) {
00498 
00499     const CSCRecHit2D& hit = (**it);
00500     const CSCLayer* layer = theChamber->layer(hit.cscDetId().layer());
00501     GlobalPoint gp = layer->toGlobal(hit.localPosition());
00502     LocalPoint lp = theChamber->toLocal(gp);
00503     float z = lp.z();
00504     ++row;
00505     matrix(row, 1) = 1.;
00506     matrix(row, 3) = z;
00507     ++row;
00508     matrix(row, 2) = 1.;
00509     matrix(row, 4) = z;
00510   }
00511   return matrix;
00512 }
00513 
00514 /* calculateError
00515  *
00516  */
00517 AlgebraicSymMatrix CSCSegAlgoShowering::calculateError() const {
00518 
00519   AlgebraicSymMatrix weights = weightMatrix();
00520   AlgebraicMatrix A = derivativeMatrix();
00521 
00522   // (AT W A)^-1
00523   // from https://www.phys.ufl.edu/~avery/fitting.html, part I
00524   int ierr;
00525   AlgebraicSymMatrix result = weights.similarityT(A);
00526   result.invert(ierr);
00527 
00528   // blithely assuming the inverting never fails...
00529   return result;
00530 }
00531 
00532 void CSCSegAlgoShowering::flipErrors( AlgebraicSymMatrix& a ) const {
00533 
00534   // The CSCSegment needs the error matrix re-arranged
00535 
00536   AlgebraicSymMatrix hold( a );
00537 
00538   // errors on slopes into upper left
00539   a(1,1) = hold(3,3);
00540   a(1,2) = hold(3,4);
00541   a(2,1) = hold(4,3);
00542   a(2,2) = hold(4,4);
00543 
00544   // errors on positions into lower right
00545   a(3,3) = hold(1,1);
00546   a(3,4) = hold(1,2);
00547   a(4,3) = hold(2,1);
00548   a(4,4) = hold(2,2);
00549 
00550   // off-diagonal elements remain unchanged
00551 
00552 }
00553 
00554 // Try to clean up segments by quickly looking at residuals
00555 void CSCSegAlgoShowering::pruneFromResidual(){
00556 
00557   // Only prune if have at least 5 hits 
00558   if ( protoSegment.size() < 5 ) return ;
00559 
00560 
00561   // Now Study residuals
00562       
00563   float maxResidual = 0.;
00564   float sumResidual = 0.;
00565   int nHits = 0;
00566   int badIndex = -1;
00567   int j = 0;
00568 
00569 
00570   ChamberHitContainer::const_iterator ih;
00571 
00572   for ( ih = protoSegment.begin(); ih != protoSegment.end(); ++ih ) {
00573     const CSCRecHit2D& hit = (**ih);
00574     const CSCLayer* layer  = theChamber->layer(hit.cscDetId().layer());
00575     GlobalPoint gp         = layer->toGlobal(hit.localPosition());
00576     LocalPoint lp          = theChamber->toLocal(gp);
00577 
00578     double u = lp.x();
00579     double v = lp.y();
00580     double z = lp.z();
00581 
00582     double du = protoIntercept.x() + protoSlope_u * z - u;
00583     double dv = protoIntercept.y() + protoSlope_v * z - v;
00584 
00585     float residual = sqrt(du*du + dv*dv);
00586 
00587     sumResidual += residual;
00588     nHits++;
00589     if ( residual > maxResidual ) {
00590       maxResidual = residual;
00591       badIndex = j;
00592     }
00593     j++;
00594   }
00595 
00596   float corrAvgResidual = (sumResidual - maxResidual)/(nHits -1);
00597 
00598   // Keep all hits 
00599   if ( maxResidual/corrAvgResidual < maxRatioResidual ) return;
00600 
00601 
00602   // Drop worse hit and recompute segment properties + fill
00603 
00604   ChamberHitContainer newProtoSegment;
00605 
00606   j = 0;
00607   for ( ih = protoSegment.begin(); ih != protoSegment.end(); ++ih ) {
00608     if ( j != badIndex ) newProtoSegment.push_back(*ih);
00609     j++;
00610   }
00611   
00612   protoSegment.clear();
00613 
00614   for ( ih = newProtoSegment.begin(); ih != newProtoSegment.end(); ++ih ) {
00615     protoSegment.push_back(*ih);
00616   }
00617 
00618   // Update segment parameters
00619   updateParameters();
00620 
00621 }
00622 
00623 
00624 
00625