Go to the documentation of this file.00001
00002 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00003 #include "FWCore/Framework/interface/Event.h"
00004 #include "FWCore/Framework/interface/MakerMacros.h"
00005 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
00006 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
00007
00008 #include "FastSimulation/CaloRecHitsProducer/interface/CaloRecHitCopy.h"
00009
00010 #include <iostream>
00011
00012 CaloRecHitCopy::CaloRecHitCopy(edm::ParameterSet const & p)
00013
00014 {
00015
00016 theInputRecHitCollectionTypes = p.getParameter<std::vector<unsigned> >("InputRecHitCollectionTypes");
00017 theInputRecHitCollections = p.getParameter<std::vector<edm::InputTag> >("InputRecHitCollections");
00018 theOutputRecHitCollections = p.getParameter<std::vector<std::string> >("OutputRecHitCollections");
00019
00020 theOutputRecHitInstances.resize(theInputRecHitCollectionTypes.size());
00021
00022 for ( unsigned input=0; input<theInputRecHitCollectionTypes.size(); ++input ) {
00023
00024 theOutputRecHitInstances[input] =
00025 theOutputRecHitCollections[input] == "none" ?
00026 false : true;
00027
00028 switch ( theInputRecHitCollectionTypes[input] ) {
00029
00030 case 1:
00031 {
00032
00033 if ( !theOutputRecHitInstances[input] )
00034 produces<ESRecHitCollection>();
00035 else
00036 produces<ESRecHitCollection>(theOutputRecHitCollections[input]);
00037 }
00038 break;
00039
00040 case 2:
00041 {
00042
00043 if ( !theOutputRecHitInstances[input] )
00044 produces<EBRecHitCollection>();
00045 else
00046 produces<EBRecHitCollection>(theOutputRecHitCollections[input]);
00047 }
00048 break;
00049
00050 case 3:
00051 {
00052
00053 if ( !theOutputRecHitInstances[input] )
00054 produces<EERecHitCollection>();
00055 else
00056 produces<EERecHitCollection>(theOutputRecHitCollections[input]);
00057 }
00058 break;
00059
00060 case 4:
00061 {
00062
00063 if ( !theOutputRecHitInstances[input] )
00064 produces<HBHERecHitCollection>();
00065 else
00066 produces<HBHERecHitCollection>(theOutputRecHitCollections[input]);
00067 }
00068 break;
00069
00070 case 5:
00071 {
00072
00073 if ( !theOutputRecHitInstances[input] )
00074 produces<HORecHitCollection>();
00075 else
00076 produces<HORecHitCollection>(theOutputRecHitCollections[input]);
00077 }
00078 break;
00079
00080 case 6:
00081 {
00082
00083 if ( !theOutputRecHitInstances[input] )
00084 produces<HFRecHitCollection>();
00085 else
00086 produces<HFRecHitCollection>(theOutputRecHitCollections[input]);
00087 }
00088 break;
00089
00090 default:
00091
00092 break;
00093
00094 }
00095
00096 }
00097
00098
00099 }
00100
00101 CaloRecHitCopy::~CaloRecHitCopy() { }
00102
00103 void
00104 CaloRecHitCopy::endRun() { }
00105
00106 void
00107 CaloRecHitCopy::produce(edm::Event & iEvent, const edm::EventSetup & es)
00108 {
00109
00110
00111 for ( unsigned input=0; input<theInputRecHitCollectionTypes.size(); ++input ) {
00112
00113 switch ( theInputRecHitCollectionTypes[input] ) {
00114
00115 case 1:
00116 {
00117
00118 std::auto_ptr< ESRecHitCollection > copiedESRecHitCollection( new ESRecHitCollection );
00119 edm::Handle<ESRecHitCollection> ESRecHits;
00120 iEvent.getByLabel(theInputRecHitCollections[input],ESRecHits);
00121 ESRecHitCollection::const_iterator itES = ESRecHits->begin();
00122 ESRecHitCollection::const_iterator lastES = ESRecHits->end();
00123
00124 copiedESRecHitCollection->reserve(ESRecHits->size());
00125 for ( ; itES!=lastES; ++itES++ ) {
00126 EcalRecHit aHit(*itES);
00127 copiedESRecHitCollection->push_back(aHit);
00128 }
00129 if ( !theOutputRecHitInstances[input] )
00130 iEvent.put(copiedESRecHitCollection);
00131 else
00132 iEvent.put(copiedESRecHitCollection,theOutputRecHitCollections[input]);
00133 }
00134 break;
00135
00136 case 2:
00137 {
00138
00139 std::auto_ptr< EBRecHitCollection > copiedEBRecHitCollection( new EBRecHitCollection );
00140 edm::Handle<EBRecHitCollection> EBRecHits;
00141 iEvent.getByLabel(theInputRecHitCollections[input],EBRecHits);
00142 EBRecHitCollection::const_iterator itEB = EBRecHits->begin();
00143 EBRecHitCollection::const_iterator lastEB = EBRecHits->end();
00144
00145 copiedEBRecHitCollection->reserve(EBRecHits->size());
00146
00147 for ( ; itEB!=lastEB; ++itEB++ ) {
00148 EcalRecHit aHit(*itEB);
00149 copiedEBRecHitCollection->push_back(aHit);
00150 }
00151 if ( !theOutputRecHitInstances[input] )
00152 iEvent.put(copiedEBRecHitCollection);
00153 else
00154 iEvent.put(copiedEBRecHitCollection,theOutputRecHitCollections[input]);
00155 }
00156 break;
00157
00158 case 3:
00159 {
00160
00161 std::auto_ptr< EERecHitCollection > copiedEERecHitCollection( new EERecHitCollection );
00162 edm::Handle<EERecHitCollection> EERecHits;
00163 iEvent.getByLabel(theInputRecHitCollections[input],EERecHits);
00164 EERecHitCollection::const_iterator itEE = EERecHits->begin();
00165 EERecHitCollection::const_iterator lastEE = EERecHits->end();
00166
00167 copiedEERecHitCollection->reserve(EERecHits->size());
00168
00169 for ( ; itEE!=lastEE; ++itEE++ ) {
00170 EcalRecHit aHit(*itEE);
00171 copiedEERecHitCollection->push_back(aHit);
00172 }
00173 if ( !theOutputRecHitInstances[input] )
00174 iEvent.put(copiedEERecHitCollection);
00175 else
00176 iEvent.put(copiedEERecHitCollection,theOutputRecHitCollections[input]);
00177 }
00178 break;
00179
00180 case 4:
00181 {
00182
00183 std::auto_ptr< HBHERecHitCollection > copiedHBHERecHitCollection( new HBHERecHitCollection );
00184 edm::Handle<HBHERecHitCollection> HBHERecHits;
00185 iEvent.getByLabel(theInputRecHitCollections[input],HBHERecHits);
00186 HBHERecHitCollection::const_iterator itHBHE = HBHERecHits->begin();
00187 HBHERecHitCollection::const_iterator lastHBHE = HBHERecHits->end();
00188
00189 copiedHBHERecHitCollection->reserve(HBHERecHits->size());
00190
00191 for ( ; itHBHE!=lastHBHE; ++itHBHE++ ) {
00192 HBHERecHit aHit(*itHBHE);
00193 copiedHBHERecHitCollection->push_back(aHit);
00194 }
00195 if ( !theOutputRecHitInstances[input] )
00196 iEvent.put(copiedHBHERecHitCollection);
00197 else
00198 iEvent.put(copiedHBHERecHitCollection,theOutputRecHitCollections[input]);
00199 }
00200 break;
00201
00202 case 5:
00203 {
00204
00205 std::auto_ptr< HORecHitCollection > copiedHORecHitCollection( new HORecHitCollection );
00206 edm::Handle<HORecHitCollection> HORecHits;
00207 iEvent.getByLabel(theInputRecHitCollections[input],HORecHits);
00208 HORecHitCollection::const_iterator itHO = HORecHits->begin();
00209 HORecHitCollection::const_iterator lastHO = HORecHits->end();
00210
00211 copiedHORecHitCollection->reserve(HORecHits->size());
00212
00213 for ( ; itHO!=lastHO; ++itHO++ ) {
00214 HORecHit aHit(*itHO);
00215 copiedHORecHitCollection->push_back(aHit);
00216 }
00217 if ( !theOutputRecHitInstances[input] )
00218 iEvent.put(copiedHORecHitCollection);
00219 else
00220 iEvent.put(copiedHORecHitCollection,theOutputRecHitCollections[input]);
00221 }
00222 break;
00223
00224 case 6:
00225 {
00226
00227 std::auto_ptr< HFRecHitCollection > copiedHFRecHitCollection( new HFRecHitCollection );
00228 edm::Handle<HFRecHitCollection> HFRecHits;
00229 iEvent.getByLabel(theInputRecHitCollections[input],HFRecHits);
00230 HFRecHitCollection::const_iterator itHF = HFRecHits->begin();
00231 HFRecHitCollection::const_iterator lastHF = HFRecHits->end();
00232
00233 copiedHFRecHitCollection->reserve(HFRecHits->size());
00234
00235 for ( ; itHF!=lastHF; ++itHF++ ) {
00236 HFRecHit aHit(*itHF);
00237 copiedHFRecHitCollection->push_back(aHit);
00238 }
00239 if ( !theOutputRecHitInstances[input] )
00240 iEvent.put(copiedHFRecHitCollection);
00241 else
00242 iEvent.put(copiedHFRecHitCollection,theOutputRecHitCollections[input]);
00243 }
00244 break;
00245
00246 default:
00247
00248 break;
00249
00250 }
00251
00252 }
00253
00254 }
00255
00256 DEFINE_FWK_MODULE(CaloRecHitCopy);