CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/RecoEcal/EgammaClusterProducers/src/UnifiedSCCollectionProducer.cc

Go to the documentation of this file.
00001 // C/C++ headers
00002 #include <iostream>
00003 #include <vector>
00004 #include <memory>
00005 
00006 // Framework
00007 #include "FWCore/Framework/interface/Event.h"
00008 #include "FWCore/Framework/interface/EventSetup.h"
00009 #include "DataFormats/Common/interface/Handle.h"
00010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00011 #include "FWCore/Framework/interface/ESHandle.h"
00012 // Reconstruction Classes
00013 #include "DataFormats/EcalRecHit/interface/EcalRecHit.h"
00014 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
00015 #include "DataFormats/EcalDetId/interface/EBDetId.h"
00016 #include "DataFormats/EgammaReco/interface/BasicCluster.h"
00017 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
00018 #include "DataFormats/EgammaReco/interface/SuperClusterFwd.h"
00019 #include "DataFormats/EgammaReco/interface/BasicClusterFwd.h"
00020 #include "CondFormats/DataRecord/interface/EcalChannelStatusRcd.h"
00021 #include "CondFormats/EcalObjects/interface/EcalChannelStatus.h"
00022 
00023 // Class header file
00024 #include "RecoEcal/EgammaClusterProducers/interface/UnifiedSCCollectionProducer.h"
00025 #include "RecoEcal/EgammaCoreTools/interface/PositionCalc.h"
00026 #include "DataFormats/CaloRecHit/interface/CaloCluster.h"
00027 
00028 
00029 /*
00030 UnifiedSCCollectionProducer:
00031 ^^^^^^^^^^^^^^^^^^^^^^^^^^^^
00032 
00033 Takes  as  input the cleaned  and the  uncleaned collection of SC
00034 and  produces two collections of SC: one with the clean SC, but flagged
00035 such that with the algoID value one can identify the SC that are also
00036 in the unclean collection and a collection with the unclean only SC.
00037 This collection has the algoID enumeration of the SC altered
00038 such that:
00039 flags = 0   (cleanedOnly)     cluster is only in the cleaned collection
00040 flags = 100 (common)          cluster is common in both collections
00041 flags = 200 (uncleanedOnly)   cluster is only in the uncleaned collection
00042 
00043 In that way the user can get hold of objects from the
00044 -  cleaned   collection only if they choose flags <  200
00045 -  uncleaned collection only if they choose flags >= 100
00046 
00047 18 Aug 2010
00048 Nikolaos Rompotis and Chris Seez  - Imperial College London
00049 many thanks to David Wardrope, Shahram Rahatlou and Federico Ferri
00050 */
00051 
00052 
00053 UnifiedSCCollectionProducer::UnifiedSCCollectionProducer(const edm::ParameterSet& ps)
00054 {
00055 
00056         // get the parameters
00057         // the cleaned collection:
00058         cleanBcCollection_ = ps.getParameter<edm::InputTag>("cleanBcCollection");
00059         cleanScCollection_ = ps.getParameter<edm::InputTag>("cleanScCollection");
00060         // the uncleaned collection
00061         uncleanBcCollection_ = ps.getParameter<edm::InputTag>("uncleanBcCollection");
00062         uncleanScCollection_ = ps.getParameter<edm::InputTag>("uncleanScCollection");
00063         // the names of the products to be produced:
00064         //
00065         // the clean collection: this is as it was before, but labeled
00066         bcCollection_ = ps.getParameter<std::string>("bcCollection");
00067         scCollection_ = ps.getParameter<std::string>("scCollection");
00068         // the unclean only collection: SC unique to the unclean collection
00069         bcCollectionUncleanOnly_ = ps.getParameter<std::string>("bcCollectionUncleanOnly");
00070         scCollectionUncleanOnly_ = ps.getParameter<std::string>("scCollectionUncleanOnly");
00071         // the products:
00072         produces< reco::BasicClusterCollection >(bcCollection_);
00073         produces< reco::SuperClusterCollection >(scCollection_);
00074         produces< reco::BasicClusterCollection >(bcCollectionUncleanOnly_);
00075         produces< reco::SuperClusterCollection >(scCollectionUncleanOnly_);
00076 
00077 }
00078 
00079 
00080 UnifiedSCCollectionProducer::~UnifiedSCCollectionProducer() {;}
00081 
00082 
00083 void UnifiedSCCollectionProducer::produce(edm::Event& evt, 
00084                                           const edm::EventSetup& es)
00085 {
00086 
00087   edm::LogInfo("UnifiedSC")<< ">>>>> Entering UnifiedSCCollectionProducer <<<<<";
00088   // get the input collections
00089   // __________________________________________________________________________
00090   //
00091   // cluster collections:
00092   edm::Handle<reco::BasicClusterCollection> pCleanBC;
00093   edm::Handle<reco::SuperClusterCollection> pCleanSC;
00094   //
00095   edm::Handle<reco::BasicClusterCollection> pUncleanBC;
00096   edm::Handle<reco::SuperClusterCollection> pUncleanSC;
00097   
00098   evt.getByLabel(cleanScCollection_, pCleanSC);
00099   if (!(pCleanSC.isValid())) 
00100     {
00101       
00102       edm::LogError("MissingInput") << "could not find clean super clusters";
00103       return;
00104     }
00105   
00106   evt.getByLabel(cleanBcCollection_, pCleanBC);
00107   if (!(pCleanBC.isValid())) 
00108     {
00109       
00110       edm::LogError("MissingInput") << "could not find " << cleanBcCollection_;
00111       return;
00112     }
00113 
00114   evt.getByLabel(uncleanBcCollection_, pUncleanBC);
00115   if (!(pUncleanBC.isValid())) 
00116     {
00117       
00118       edm::LogError("MissingInput") << "could not find " <<  uncleanBcCollection_;
00119       return;
00120     }
00121 
00122 
00123 
00124         evt.getByLabel(uncleanScCollection_, pUncleanSC);
00125         if (!(pUncleanSC.isValid())) 
00126         {
00127 
00128           edm::LogError("MissingInput")<< "could not handle unclean super clusters!" ;
00129           return;
00130         }
00131 
00132         // the collections to be produced ___________________________________________
00133         reco::BasicClusterCollection basicClusters;
00134         reco::SuperClusterCollection superClusters;
00135         //
00136         reco::BasicClusterCollection basicClustersUncleanOnly;
00137         reco::SuperClusterCollection superClustersUncleanOnly;
00138         //
00139         // run over the uncleaned SC and check how many of them are matched to 
00140         // the cleaned ones
00141         // if you find a matched one, then keep the info that it is matched 
00142         //    along with which clean SC was matched + its basic clusters
00143         // if you find an unmatched one, keep the info and store its basic clusters
00144         //
00145         // 
00146         int uncleanSize = pUncleanSC->size();
00147         int cleanSize =   pCleanSC->size();
00148         
00149         LogTrace("UnifiedSC") << "Size of Clean Collection: " << cleanSize 
00150                         << ", uncleanSize: " << uncleanSize;
00151 
00152 
00153         // keep the indices
00154         std::vector<int> inUncleanOnlyInd;      // counting the unclean
00155         std::vector<int> inCleanInd;            // counting the unclean
00156         std::vector<int> inCleanOnlyInd;        // counting the clean
00157         std::vector<DetId> scUncleanSeedDetId;  // counting the unclean
00158         std::vector<DetId> scCleanSeedDetId;    // counting the clean
00159         // ontains the index of the SC that owns that BS
00160         // first basic cluster index, second: 0 for unclean and 1 for clean
00161         std::vector< std::pair<int, int> > basicClusterOwner; 
00162         std::vector< std::pair<int, int> > basicClusterOwnerUncleanOnly; 
00163         // if this basic cluster is a seed it is 1
00164         std::vector<int> uncleanBasicClusterIsSeed;
00165 
00166         // loop over unclean SC _____________________________________________________
00167         for (int isc =0; isc< uncleanSize; ++isc) {
00168                 reco::SuperClusterRef unscRef( pUncleanSC, isc);    
00169                 const std::vector< std::pair<DetId, float> > & uhits = unscRef->hitsAndFractions();
00170                 int uhitsSize = uhits.size();
00171                 bool foundTheSame = false;
00172                 for (int jsc=0; jsc < cleanSize; ++jsc) { // loop over the cleaned SC
00173                         reco::SuperClusterRef cscRef( pCleanSC, jsc );
00174                         const std::vector<std::pair<DetId,float> > & chits = cscRef->hitsAndFractions();
00175                         int chitsSize =  chits.size();
00176                         foundTheSame = false;
00177                         if (unscRef->seed()->seed()==cscRef->seed()->seed() && chitsSize == uhitsSize) { 
00178                                 // if the clusters are exactly the same then because the clustering
00179                                 // algorithm works in a deterministic way, the order of the rechits
00180                                 // will be the same
00181                                 foundTheSame = true;
00182                                 for (int i=0; i< chitsSize; ++i) {
00183                                         if (uhits[i].first != chits[i].first ) { 
00184                                           foundTheSame=false;    
00185                                           break;
00186                                         }
00187                                 }
00188                                 
00189                         }
00190                         if (foundTheSame) { // ok you have found it:
00191                                 // this supercluster belongs to both collections
00192                                 inUncleanOnlyInd.push_back(0);
00193                                 inCleanInd.push_back(jsc); // keeps the index of the clean SC
00194                                 scUncleanSeedDetId.push_back(unscRef->seed()->seed());
00195                                 //
00196                                 // keep its basic clusters:
00197                                 for (reco::CaloCluster_iterator bciter = unscRef->clustersBegin(); bciter != unscRef->clustersEnd(); ++bciter) {
00198                                         // the basic clusters
00199                                         basicClusters.push_back(**bciter);
00200                                         // index of the unclean SC
00201                                         basicClusterOwner.push_back( std::make_pair(isc,0) ); 
00202                                 } 
00203                                 break; // break the loop over unclean sc
00204                         }
00205                 }
00206                 if (not foundTheSame) { // this SC is only in the unclean collection
00207                         // mark it as unique in the uncleaned
00208                         inUncleanOnlyInd.push_back(1);
00209                         scUncleanSeedDetId.push_back(unscRef->seed()->seed());
00210                         // keep all its basic clusters
00211                         for (reco::CaloCluster_iterator bciter = unscRef->clustersBegin(); bciter != unscRef->clustersEnd(); ++bciter) {
00212                                 // the basic clusters
00213                                 basicClustersUncleanOnly.push_back(**bciter);
00214                                 basicClusterOwnerUncleanOnly.push_back( std::make_pair(isc,0) );
00215                         }
00216                 }
00217         } // loop over the unclean SC _______________________________________________
00218         //
00219         int inCleanSize = inCleanInd.size();
00220         //
00221         // loop over the clean SC, check that are not in common with the unclean
00222         // ones and then store their SC as before ___________________________________
00223         for (int jsc =0; jsc< cleanSize; ++jsc) {
00224                 // check whether this index is already in the common collection
00225                 bool takenAlready = false;
00226                 for (int j=0; j< inCleanSize; ++j) {
00227                         if (jsc == inCleanInd[j]) { takenAlready = true ;break;}
00228                 }
00229                 if (takenAlready) {
00230                         inCleanOnlyInd.push_back(0);
00231                         scCleanSeedDetId.push_back(DetId(0));
00232                         continue;
00233                 }
00234                 inCleanOnlyInd.push_back(1);
00235                 reco::SuperClusterRef cscRef( pCleanSC, jsc );
00236                 scCleanSeedDetId.push_back(cscRef->seed()->seed());
00237                 for (reco::CaloCluster_iterator bciter = cscRef->clustersBegin(); bciter != cscRef->clustersEnd(); ++bciter) {
00238                         // the basic clusters
00239                         basicClusters.push_back(**bciter);
00240                         basicClusterOwner.push_back( std::make_pair(jsc,1) );
00241                 }
00242         } // end loop over clean SC _________________________________________________
00243         //
00244         //
00245 
00246         // Final check: in the endcap BC may exist that are not associated to SC,
00247         // we need to recover them as well (e.g. multi5x5 algo)
00248         // This is should be optimized (SA, 20110621)
00249 
00250 
00251         // loop on original clean BC collection and see if the BC is missing from the new one
00252         for (reco::BasicClusterCollection::const_iterator bc = pCleanBC->begin();
00253              bc!=pCleanBC->end();
00254              ++bc){
00255           
00256           bool foundTheSame = false;
00257           for (reco::BasicClusterCollection::const_iterator cleanonly_bc = basicClusters.begin();
00258                cleanonly_bc!=basicClusters.end();
00259                ++cleanonly_bc){     
00260             
00261             const std::vector<std::pair<DetId,float> > & chits = bc->hitsAndFractions();
00262             int chitsSize =  chits.size();
00263             
00264             const std::vector<std::pair<DetId,float> > & uhits = cleanonly_bc->hitsAndFractions();
00265             int uhitsSize =  uhits.size();
00266             
00267             
00268             if (cleanonly_bc->seed()==bc->seed() && chitsSize == uhitsSize) {
00269                 foundTheSame = true;
00270                 for (int i=0; i< chitsSize; ++i) {
00271                   if (uhits[i].first != chits[i].first ) { 
00272                     foundTheSame=false;    
00273                     break;
00274                   }
00275                 }
00276             }
00277                     
00278           } // loop on new clean BC collection
00279           
00280           // clean basic cluster is not associated to SC and does not belong to the 
00281           // new collection, add it
00282           if (!foundTheSame){         
00283             basicClusters.push_back(*bc);
00284             LogTrace("UnifiedSC")<<"found BC to add that was not associated to any SC";
00285           }
00286           
00287         } // loop on original clean BC collection
00288 
00289 
00290         // at this point we have the basic cluster collection ready
00291         // Up to index   basicClusterOwner.size() we have the BC owned by a SC
00292         // The remaining are BCs not owned by a SC
00293 
00294         int bcSize = (int) basicClusterOwner.size();
00295         int bcSizeUncleanOnly = (int) basicClustersUncleanOnly.size();
00296 
00297         LogTrace("UnifiedSC") << "Found cleaned SC: " << cleanSize 
00298                               <<  " uncleaned SC: "   << uncleanSize ;
00299         //
00300         // export the clusters to the event from the clean clusters
00301         std::auto_ptr< reco::BasicClusterCollection> 
00302                 basicClusters_p(new reco::BasicClusterCollection);
00303         basicClusters_p->assign(basicClusters.begin(), basicClusters.end());
00304         edm::OrphanHandle<reco::BasicClusterCollection> bccHandle =  
00305                 evt.put(basicClusters_p, bcCollection_);
00306         if (!(bccHandle.isValid())) {
00307               
00308           edm::LogWarning("MissingInput")<< "could not handle the new BasicClusters!";
00309           return;
00310         }
00311         reco::BasicClusterCollection basicClustersProd = *bccHandle;
00312 
00313         LogTrace("UnifiedSC")<< "Got the BasicClusters from the event again" ;
00314         //
00315         // export the clusters to the event: from the unclean only clusters
00316         std::auto_ptr< reco::BasicClusterCollection> 
00317                 basicClustersUncleanOnly_p(new reco::BasicClusterCollection);
00318         basicClustersUncleanOnly_p->assign(basicClustersUncleanOnly.begin(), 
00319                                            basicClustersUncleanOnly.end());
00320         edm::OrphanHandle<reco::BasicClusterCollection> bccHandleUncleanOnly =  
00321                 evt.put(basicClustersUncleanOnly_p, bcCollectionUncleanOnly_);
00322         if (!(bccHandleUncleanOnly.isValid())) {
00323 
00324           edm::LogWarning("MissingInput")<< "could not handle the new BasicClusters (Unclean Only)!" ;
00325           return;
00326         }
00327         reco::BasicClusterCollection basicClustersUncleanOnlyProd = *bccHandleUncleanOnly;
00328         LogTrace("UnifiedSC")<< "Got the BasicClusters from the event again  (Unclean Only)" ;
00329         //
00330 
00331         // now we can build the SC collection
00332         //
00333         // start again from the unclean collection
00334         // all the unclean SC will become members of the new collection
00335         // with different algoIDs ___________________________________________________
00336         for (int isc=0; isc< uncleanSize; ++isc) {
00337                 reco::CaloClusterPtrVector clusterPtrVector;
00338                 // the seed is the basic cluster with the highest energy
00339                 reco::CaloClusterPtr seed; 
00340                 if (inUncleanOnlyInd[isc] == 1) { // unclean SC Unique in Unclean
00341                         for (int jbc=0; jbc< bcSizeUncleanOnly; ++jbc) {
00342                                 std::pair<int, int> theBcOwner = basicClusterOwnerUncleanOnly[jbc];
00343                                 if (theBcOwner.first == isc && theBcOwner.second == 0) {
00344                                         reco::CaloClusterPtr currentClu=reco::CaloClusterPtr(bccHandleUncleanOnly,jbc);
00345                                         clusterPtrVector.push_back(currentClu);
00346                                         if (scUncleanSeedDetId[isc] == currentClu->seed()) {
00347                                                 seed = currentClu;
00348                                         }
00349                                 }
00350                         }
00351 
00352                 }
00353                 else { // unclean SC common in clean and unclean
00354                         for (int jbc=0; jbc< bcSize; ++jbc) {
00355                                 std::pair<int, int> theBcOwner = basicClusterOwner[jbc];
00356                                 if (theBcOwner.first == isc && theBcOwner.second == 0) {
00357                                         reco::CaloClusterPtr currentClu=reco::CaloClusterPtr(bccHandle,jbc);
00358                                         clusterPtrVector.push_back(currentClu);
00359                                         if (scUncleanSeedDetId[isc] == currentClu->seed()) {
00360                                                 seed = currentClu;
00361                                         }
00362                                 }
00363                         }
00364                 }
00365                 //std::cout << "before getting the uncl" << std::endl;
00366                 reco::SuperClusterRef unscRef( pUncleanSC, isc ); 
00367                 reco::SuperCluster newSC(unscRef->energy(), unscRef->position(), 
00368                                          seed, clusterPtrVector );
00369                 // now set the algoID for this SC again
00370                 if (inUncleanOnlyInd[isc] == 1) {
00371                         // set up the quality to unclean only .............
00372                         newSC.setFlags(reco::CaloCluster::uncleanOnly);
00373                         superClustersUncleanOnly.push_back(newSC);
00374                 }
00375                 else {
00376                         // set up the quality to common  .............
00377                         newSC.setFlags(reco::CaloCluster::common);
00378                         superClusters.push_back(newSC);
00379                 }
00380                 // now you can store your SC
00381 
00382         } // end loop over unclean SC _______________________________________________
00383         //  flags numbering scheme
00384         //  flags =   0 = cleanedOnly     cluster is only in the cleaned collection
00385         //  flags = 100 = common          cluster is common in both collections
00386         //  flags = 200 = uncleanedOnly   cluster is only in the uncleaned collection
00387 
00388         // now loop over the clean SC and do the same but now you have to avoid the
00389         // the duplicated ones ______________________________________________________
00390         for (int jsc=0; jsc< cleanSize; ++jsc) {
00391                 //std::cout << "working in cl #" << jsc << std::endl;
00392                 // check that the SC is not in the unclean collection
00393                 if (inCleanOnlyInd[jsc] == 0) continue;
00394                 reco::CaloClusterPtrVector clusterPtrVector;
00395                 // the seed is the basic cluster with the highest energy
00396                 reco::CaloClusterPtr seed; 
00397                 for (int jbc=0; jbc< bcSize; ++jbc) {
00398                         std::pair<int, int> theBcOwner = basicClusterOwner[jbc];
00399                         if (theBcOwner.first == jsc && theBcOwner.second == 1) {
00400                                 reco::CaloClusterPtr currentClu=reco::CaloClusterPtr(bccHandle,jbc);
00401                                 clusterPtrVector.push_back(currentClu);
00402                                 if (scCleanSeedDetId[jsc] == currentClu->seed()) {
00403                                         seed = currentClu;
00404                                 }
00405                         }
00406                 }
00407                 reco::SuperClusterRef cscRef( pCleanSC, jsc ); 
00408                 reco::SuperCluster newSC(cscRef->energy(), cscRef->position(),
00409                                          seed, clusterPtrVector );
00410                 newSC.setFlags(reco::CaloCluster::cleanOnly);
00411 
00412                 // add it to the collection:
00413                 superClusters.push_back(newSC);
00414 
00415         } // end loop over clean SC _________________________________________________
00416 
00417 
00418 
00419 
00420 
00421 
00422 
00423         LogTrace("UnifiedSC")<< "New SC collection was created";
00424 
00425         std::auto_ptr< reco::SuperClusterCollection> 
00426                 superClusters_p(new reco::SuperClusterCollection);
00427         superClusters_p->assign(superClusters.begin(), superClusters.end());
00428 
00429         evt.put(superClusters_p, scCollection_);
00430         
00431         LogTrace("UnifiedSC") << "Clusters (Basic/Super) added to the Event! :-)";
00432 
00433         std::auto_ptr< reco::SuperClusterCollection> 
00434                 superClustersUncleanOnly_p(new reco::SuperClusterCollection);
00435         superClustersUncleanOnly_p->assign(superClustersUncleanOnly.begin(), 
00436                                            superClustersUncleanOnly.end());
00437 
00438         evt.put(superClustersUncleanOnly_p, scCollectionUncleanOnly_);
00439 
00440         // ----- debugging ----------
00441         // print the new collection SC quantities
00442 
00443         // print out the clean collection SC
00444         LogTrace("UnifiedSC") << "Clean Collection SC ";
00445         for (int i=0; i < cleanSize; ++i) {
00446           reco::SuperClusterRef cscRef( pCleanSC, i );
00447           LogTrace("UnifiedSC") << " >>> clean    #" << i << "; Energy: " << cscRef->energy()
00448                                 << " eta: " << cscRef->eta() 
00449                                 << " sc seed detid: " << cscRef->seed()->seed().rawId()
00450                                 << std::endl;
00451                 }
00452                 // the unclean SC
00453                 LogTrace("UnifiedSC") << "Unclean Collection SC ";
00454                 for (int i=0; i < uncleanSize; ++i) {
00455                         reco::SuperClusterRef uscRef( pUncleanSC, i );
00456                         LogTrace("UnifiedSC") << " >>> unclean  #" << i << "; Energy: " << uscRef->energy()
00457                                 << " eta: " << uscRef->eta() 
00458                                 << " sc seed detid: " << uscRef->seed()->seed().rawId();
00459                 }
00460                 // the new collection
00461                 LogTrace("UnifiedSC")<< "The new SC clean collection with size "<< superClusters.size()  << std::endl;
00462 
00463                 int new_unclean = 0, new_clean=0;
00464                 for (int i=0; i < (int) superClusters.size(); ++i) {
00465                         const reco::SuperCluster & nsc = superClusters[i];
00466                         LogTrace("UnifiedSC") << "SC was got" << std::endl
00467                                 << " ---> energy: " << nsc.energy() << std::endl
00468                                 << " ---> eta: " << nsc.eta() << std::endl
00469                                 << " ---> inClean: " << nsc.isInClean() << std::endl
00470                                 << " ---> id: " << nsc.seed()->seed().rawId() << std::endl
00471                                 << " >>> newSC    #" << i << "; Energy: " << nsc.energy()
00472                                 << " eta: " << nsc.eta()  << " isClean=" 
00473                                 << nsc.isInClean() << " isUnclean=" << nsc.isInUnclean()
00474                                 << " sc seed detid: " << nsc.seed()->seed().rawId();
00475 
00476                         if (nsc.isInUnclean()) ++new_unclean;
00477                         if (nsc.isInClean()) ++new_clean;
00478                 }
00479                 LogTrace("UnifiedSC")<< "The new SC unclean only collection with size "<< superClustersUncleanOnly.size();
00480                 for (int i=0; i < (int) superClustersUncleanOnly.size(); ++i) {
00481                         const reco::SuperCluster nsc = superClustersUncleanOnly[i];
00482                         LogTrace ("UnifiedSC") << " >>> newSC    #" << i << "; Energy: " << nsc.energy()
00483                                 << " eta: " << nsc.eta()  << " isClean=" 
00484                                 << nsc.isInClean() << " isUnclean=" << nsc.isInUnclean()
00485                                 << " sc seed detid: " << nsc.seed()->seed().rawId();
00486                         if (nsc.isInUnclean()) ++new_unclean;
00487                         if (nsc.isInClean()) ++new_clean;      
00488                 }
00489                 if ( (new_unclean != uncleanSize) || (new_clean != cleanSize) ) {
00490                         LogTrace("UnifiedSC") << ">>>>!!!!!! MISMATCH: new unclean/ old unclean= " 
00491                                 << new_unclean << " / " << uncleanSize 
00492                                 << ", new clean/ old clean" << new_clean << " / " << cleanSize;
00493                 }
00494         
00495 }