00001
00002 #include "SimG4Core/Notification/interface/TrackInformation.h"
00003 #include "DetectorDescription/Core/interface/DDFilter.h"
00004 #include "DetectorDescription/Core/interface/DDFilteredView.h"
00005 #include "DetectorDescription/Core/interface/DDSolid.h"
00006 #include "DetectorDescription/Core/interface/DDSplit.h"
00007 #include "DetectorDescription/Core/interface/DDValue.h"
00008
00009 #include "G4LogicalVolumeStore.hh"
00010 #include "G4LogicalVolume.hh"
00011 #include "G4Step.hh"
00012 #include "G4Track.hh"
00013 #include "G4VProcess.hh"
00014 #include "G4Poisson.hh"
00015
00016
00017 #include "FWCore/ServiceRegistry/interface/Service.h"
00018 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00019 #include <TTree.h>
00020
00021
00022 #include "SimG4CMS/CherenkovAnalysis/interface/DreamSD.h"
00023 #include "SimG4CMS/CherenkovAnalysis/interface/PMTResponse.h"
00024
00025
00026 DreamSD::DreamSD(G4String name, const DDCompactView & cpv,
00027 SensitiveDetectorCatalog & clg,
00028 edm::ParameterSet const & p, const SimTrackManager* manager) :
00029 CaloSD(name, cpv, clg, p, manager) {
00030
00031 edm::ParameterSet m_EC = p.getParameter<edm::ParameterSet>("ECalSD");
00032 useBirk= m_EC.getParameter<bool>("UseBirkLaw");
00033 doCherenkov_ = m_EC.getParameter<bool>("doCherenkov");
00034 birk1 = m_EC.getParameter<double>("BirkC1")*(g/(MeV*cm2));
00035 birk2 = m_EC.getParameter<double>("BirkC2");
00036 birk3 = m_EC.getParameter<double>("BirkC3");
00037 slopeLY= m_EC.getParameter<double>("SlopeLightYield");
00038 readBothSide_ = m_EC.getUntrackedParameter<bool>("ReadBothSide", false);
00039
00040 edm::LogInfo("EcalSim") << "Constructing a DreamSD with name " << GetName() << "\n"
00041 << "DreamSD:: Use of Birks law is set to "
00042 << useBirk << " with three constants kB = "
00043 << birk1 << ", C1 = " << birk2 << ", C2 = "
00044 << birk3 << "\n"
00045 << " Slope for Light yield is set to "
00046 << slopeLY << "\n"
00047 << " Parameterization of Cherenkov is set to "
00048 << doCherenkov_ << " and readout both sides is "
00049 << readBothSide_;
00050
00051 initMap(name,cpv);
00052
00053
00054 edm::Service<TFileService> tfile;
00055
00056 if ( !tfile.isAvailable() )
00057 throw cms::Exception("BadConfig") << "TFileService unavailable: "
00058 << "please add it to config file";
00059
00060 ntuple_ = tfile->make<TTree>("tree","Cherenkov photons");
00061 if (doCherenkov_) {
00062 ntuple_->Branch("nphotons",&nphotons_,"nphotons/I");
00063 ntuple_->Branch("px",px_,"px[nphotons]/F");
00064 ntuple_->Branch("py",py_,"py[nphotons]/F");
00065 ntuple_->Branch("pz",pz_,"pz[nphotons]/F");
00066 ntuple_->Branch("x",x_,"x[nphotons]/F");
00067 ntuple_->Branch("y",y_,"y[nphotons]/F");
00068 ntuple_->Branch("z",z_,"z[nphotons]/F");
00069 }
00070
00071 }
00072
00073
00074 bool DreamSD::ProcessHits(G4Step * aStep, G4TouchableHistory *) {
00075
00076 if (aStep == NULL) {
00077 return true;
00078 } else {
00079 side = 1;
00080 if (getStepInfo(aStep)) {
00081 if (hitExists() == false && edepositEM+edepositHAD>0.)
00082 currentHit = createNewHit();
00083 if (readBothSide_) {
00084 side = -1;
00085 getStepInfo(aStep);
00086 if (hitExists() == false && edepositEM+edepositHAD>0.)
00087 currentHit = createNewHit();
00088 }
00089 }
00090 }
00091 return true;
00092 }
00093
00094
00095
00096 bool DreamSD::getStepInfo(G4Step* aStep) {
00097
00098 preStepPoint = aStep->GetPreStepPoint();
00099 theTrack = aStep->GetTrack();
00100 G4String nameVolume = preStepPoint->GetPhysicalVolume()->GetName();
00101
00102
00103 double weight = 1.;
00104 weight *= curve_LY(aStep, side);
00105 if (useBirk) weight *= getAttenuation(aStep, birk1, birk2, birk3);
00106 edepositEM = aStep->GetTotalEnergyDeposit() * weight;
00107 LogDebug("EcalSim") << "DreamSD:: " << nameVolume << " Side " << side
00108 <<" Light Collection Efficiency " << weight
00109 << " Weighted Energy Deposit " << edepositEM/MeV
00110 << " MeV";
00111
00112 if ( doCherenkov_ ) {
00113 edepositHAD = cherenkovDeposit_( aStep );
00114 } else {
00115 edepositHAD = 0;
00116 }
00117
00118 double time = (aStep->GetPostStepPoint()->GetGlobalTime())/nanosecond;
00119 unsigned int unitID= setDetUnitId(aStep);
00120 if (side < 0) unitID++;
00121 TrackInformation * trkInfo = (TrackInformation *)(theTrack->GetUserInformation());
00122 int primaryID;
00123
00124 if (trkInfo)
00125 primaryID = trkInfo->getIDonCaloSurface();
00126 else
00127 primaryID = 0;
00128
00129 if (primaryID == 0) {
00130 edm::LogWarning("EcalSim") << "CaloSD: Problem with primaryID **** set by "
00131 << "force to TkID **** "
00132 << theTrack->GetTrackID() << " in Volume "
00133 << preStepPoint->GetTouchable()->GetVolume(0)->GetName();
00134 primaryID = theTrack->GetTrackID();
00135 }
00136
00137 bool flag = (unitID > 0);
00138 G4TouchableHistory* touch =(G4TouchableHistory*)(theTrack->GetTouchable());
00139 if (flag) {
00140 currentID.setID(unitID, time, primaryID, 0);
00141
00142 LogDebug("EcalSim") << "CaloSD:: GetStepInfo for"
00143 << " PV " << touch->GetVolume(0)->GetName()
00144 << " PVid = " << touch->GetReplicaNumber(0)
00145 << " MVid = " << touch->GetReplicaNumber(1)
00146 << " Unit " << currentID.unitID()
00147 << " Edeposit = " << edepositEM << " " << edepositHAD;
00148 } else {
00149 LogDebug("EcalSim") << "CaloSD:: GetStepInfo for"
00150 << " PV " << touch->GetVolume(0)->GetName()
00151 << " PVid = " << touch->GetReplicaNumber(0)
00152 << " MVid = " << touch->GetReplicaNumber(1)
00153 << " Unit " << std::hex << unitID << std::dec
00154 << " Edeposit = " << edepositEM << " " << edepositHAD;
00155 }
00156 return flag;
00157
00158 }
00159
00160
00161
00162 void DreamSD::initRun() {
00163
00164
00165 DimensionMap::const_iterator ite = xtalLMap.begin();
00166 G4LogicalVolume* lv = (ite->first);
00167 G4Material* material = lv->GetMaterial();
00168 edm::LogInfo("EcalSim") << "DreamSD::initRun: Initializes for material "
00169 << material->GetName() << " in " << lv->GetName();
00170 materialPropertiesTable = material->GetMaterialPropertiesTable();
00171 if ( !materialPropertiesTable ) {
00172 if ( !setPbWO2MaterialProperties_( material ) ) {
00173 edm::LogWarning("EcalSim") << "Couldn't retrieve material properties table\n"
00174 << " Material = " << material->GetName();
00175 }
00176 materialPropertiesTable = material->GetMaterialPropertiesTable();
00177 }
00178 }
00179
00180
00181
00182 uint32_t DreamSD::setDetUnitId(G4Step * aStep) {
00183 const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
00184 uint32_t id = (touch->GetReplicaNumber(1))*10 + (touch->GetReplicaNumber(0));
00185 LogDebug("EcalSim") << "DreamSD:: ID " << id;
00186 return id;
00187 }
00188
00189
00190
00191 void DreamSD::initMap(G4String sd, const DDCompactView & cpv) {
00192
00193 G4String attribute = "ReadOutName";
00194 DDSpecificsFilter filter;
00195 DDValue ddv(attribute,sd,0);
00196 filter.setCriteria(ddv,DDSpecificsFilter::equals);
00197 DDFilteredView fv(cpv);
00198 fv.addFilter(filter);
00199 fv.firstChild();
00200
00201 const G4LogicalVolumeStore * lvs = G4LogicalVolumeStore::GetInstance();
00202 std::vector<G4LogicalVolume *>::const_iterator lvcite;
00203 bool dodet=true;
00204 while (dodet) {
00205 const DDSolid & sol = fv.logicalPart().solid();
00206 std::vector<double> paras(sol.parameters());
00207 G4String name = sol.name().name();
00208 G4LogicalVolume* lv=0;
00209 for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++)
00210 if ((*lvcite)->GetName() == name) {
00211 lv = (*lvcite);
00212 break;
00213 }
00214 LogDebug("EcalSim") << "DreamSD::initMap (for " << sd << "): Solid "
00215 << name << " Shape " << sol.shape() <<" Parameter 0 = "
00216 << paras[0] << " Logical Volume " << lv;
00217 double length = 0, width = 0;
00218
00219 std::sort( paras.begin(), paras.end() );
00220 length = 2.0*paras.back();
00221 width = 2.0*paras.front();
00222 xtalLMap.insert( std::pair<G4LogicalVolume*,Doubles>(lv,Doubles(length,width)) );
00223 dodet = fv.next();
00224 }
00225 LogDebug("EcalSim") << "DreamSD: Length Table for " << attribute << " = "
00226 << sd << ":";
00227 DimensionMap::const_iterator ite = xtalLMap.begin();
00228 int i=0;
00229 for (; ite != xtalLMap.end(); ite++, i++) {
00230 G4String name = "Unknown";
00231 if (ite->first != 0) name = (ite->first)->GetName();
00232 LogDebug("EcalSim") << " " << i << " " << ite->first << " " << name
00233 << " L = " << ite->second.first
00234 << " W = " << ite->second.second;
00235 }
00236 }
00237
00238
00239 double DreamSD::curve_LY(G4Step* aStep, int flag) {
00240
00241 G4StepPoint* stepPoint = aStep->GetPreStepPoint();
00242 G4LogicalVolume* lv = stepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
00243 G4String nameVolume= lv->GetName();
00244
00245 double weight = 1.;
00246 G4ThreeVector localPoint = setToLocal(stepPoint->GetPosition(),
00247 stepPoint->GetTouchable());
00248 double crlength = crystalLength(lv);
00249 double localz = localPoint.x();
00250 double dapd = 0.5 * crlength - flag*localz;
00251 if (dapd >= -0.1 || dapd <= crlength+0.1) {
00252 if (dapd <= 100.)
00253 weight = 1.0 + slopeLY - dapd * 0.01 * slopeLY;
00254 } else {
00255 edm::LogWarning("EcalSim") << "DreamSD: light coll curve : wrong distance "
00256 << "to APD " << dapd << " crlength = "
00257 << crlength << " crystal name = " << nameVolume
00258 << " z of localPoint = " << localz
00259 << " take weight = " << weight;
00260 }
00261 LogDebug("EcalSim") << "DreamSD, light coll curve : " << dapd
00262 << " crlength = " << crlength
00263 << " crystal name = " << nameVolume
00264 << " z of localPoint = " << localz
00265 << " take weight = " << weight;
00266 return weight;
00267 }
00268
00269
00270 const double DreamSD::crystalLength(G4LogicalVolume* lv) const {
00271
00272 double length= -1.;
00273 DimensionMap::const_iterator ite = xtalLMap.find(lv);
00274 if (ite != xtalLMap.end()) length = ite->second.first;
00275 return length;
00276
00277 }
00278
00279
00280 const double DreamSD::crystalWidth(G4LogicalVolume* lv) const {
00281
00282 double width= -1.;
00283 DimensionMap::const_iterator ite = xtalLMap.find(lv);
00284 if (ite != xtalLMap.end()) width = ite->second.second;
00285 return width;
00286
00287 }
00288
00289
00290
00291
00292
00293 double DreamSD::cherenkovDeposit_( G4Step* aStep ) {
00294
00295 double cherenkovEnergy = 0;
00296 if (!materialPropertiesTable) return cherenkovEnergy;
00297 G4Material* material = aStep->GetTrack()->GetMaterial();
00298
00299
00300 const G4MaterialPropertyVector* Rindex = materialPropertiesTable->GetProperty("RINDEX");
00301 if ( Rindex == NULL ) {
00302 edm::LogWarning("EcalSim") << "Couldn't retrieve refractive index";
00303 return cherenkovEnergy;
00304 }
00305
00306 LogDebug("EcalSim") << "Material properties: " << "\n"
00307 << " Pmin = " << Rindex->GetMinPhotonEnergy()
00308 << " Pmax = " << Rindex->GetMaxPhotonEnergy();
00309
00310
00311 G4StepPoint* pPreStepPoint = aStep->GetPreStepPoint();
00312 G4StepPoint* pPostStepPoint = aStep->GetPostStepPoint();
00313 G4ThreeVector x0 = pPreStepPoint->GetPosition();
00314 G4ThreeVector p0 = aStep->GetDeltaPosition().unit();
00315 const G4DynamicParticle* aParticle = aStep->GetTrack()->GetDynamicParticle();
00316 const double charge = aParticle->GetDefinition()->GetPDGCharge();
00317
00318 double beta = 0.5*( pPreStepPoint->GetBeta() + pPostStepPoint->GetBeta() );
00319 double BetaInverse = 1.0/beta;
00320
00321 LogDebug("EcalSim") << "Particle properties: " << "\n"
00322 << " charge = " << charge
00323 << " beta = " << beta;
00324
00325
00326 double meanNumberOfPhotons = getAverageNumberOfPhotons_( charge, beta, material, Rindex );
00327 if ( meanNumberOfPhotons <= 0.0 ) {
00328 LogDebug("EcalSim") << "Mean number of photons is zero: " << meanNumberOfPhotons
00329 << ", stopping here";
00330 return cherenkovEnergy;
00331 }
00332
00333
00334 meanNumberOfPhotons *= aStep->GetStepLength();
00335
00336
00337 int numPhotons = static_cast<int>( G4Poisson(meanNumberOfPhotons) );
00338
00339 if ( numPhotons <= 0 ) {
00340 LogDebug("EcalSim") << "Poission number of photons is zero: " << numPhotons
00341 << ", stopping here";
00342 return cherenkovEnergy;
00343 }
00344
00345
00346 double Pmin = Rindex->GetMinPhotonEnergy();
00347 double Pmax = Rindex->GetMaxPhotonEnergy();
00348 double dp = Pmax - Pmin;
00349 double maxCos = BetaInverse / Rindex->GetMaxProperty();
00350 double maxSin2 = (1.0 - maxCos) * (1.0 + maxCos);
00351
00352
00353 for ( int iPhoton = 0; iPhoton<numPhotons; ++iPhoton ) {
00354
00355
00356 double randomNumber;
00357 double sampledMomentum, sampledRI;
00358 double cosTheta, sin2Theta;
00359
00360
00361 do {
00362 randomNumber = G4UniformRand();
00363 sampledMomentum = Pmin + randomNumber * dp;
00364 sampledRI = Rindex->GetProperty(sampledMomentum);
00365 cosTheta = BetaInverse / sampledRI;
00366
00367 sin2Theta = (1.0 - cosTheta)*(1.0 + cosTheta);
00368 randomNumber = G4UniformRand();
00369
00370 } while (randomNumber*maxSin2 > sin2Theta);
00371
00372
00373
00374 randomNumber = G4UniformRand();
00375
00376 double phi = twopi*randomNumber;
00377 double sinPhi = sin(phi);
00378 double cosPhi = cos(phi);
00379
00380
00381
00382
00383 double sinTheta = sqrt(sin2Theta);
00384 double px = sinTheta*cosPhi;
00385 double py = sinTheta*sinPhi;
00386 double pz = cosTheta;
00387 G4ThreeVector photonDirection(px, py, pz);
00388
00389
00390 photonDirection.rotateUz(p0);
00391
00392
00393 randomNumber = G4UniformRand();
00394 G4ThreeVector photonPosition = x0 + randomNumber * aStep->GetDeltaPosition();
00395 G4ThreeVector photonMomentum = sampledMomentum*photonDirection;
00396
00397
00398 cherenkovEnergy += getPhotonEnergyDeposit_( photonMomentum, photonPosition, aStep );
00399
00400
00401 nphotons_ = numPhotons;
00402 px_[iPhoton] = photonMomentum.x();
00403 py_[iPhoton] = photonMomentum.y();
00404 pz_[iPhoton] = photonMomentum.z();
00405 x_[iPhoton] = photonPosition.x();
00406 y_[iPhoton] = photonPosition.y();
00407 z_[iPhoton] = photonPosition.z();
00408 }
00409
00410
00411 ntuple_->Fill();
00412
00413
00414 return cherenkovEnergy;
00415
00416 }
00417
00418
00419
00420
00421
00422 const double DreamSD::getAverageNumberOfPhotons_( const double charge,
00423 const double beta,
00424 const G4Material* aMaterial,
00425 const G4MaterialPropertyVector* Rindex ) const
00426 {
00427 const G4double rFact = 369.81/(eV * cm);
00428
00429 if( beta <= 0.0 ) return 0.0;
00430
00431 double BetaInverse = 1./beta;
00432
00433
00434
00435
00436
00437
00438 double Pmin = Rindex->GetMinPhotonEnergy();
00439 double Pmax = Rindex->GetMaxPhotonEnergy();
00440
00441
00442 double nMin = Rindex->GetMinProperty();
00443 double nMax = Rindex->GetMaxProperty();
00444
00445
00446 double CAImax = chAngleIntegrals_->GetMaxValue();
00447
00448 double dp = 0., ge = 0., CAImin = 0.;
00449
00450
00451 if ( nMax < BetaInverse) { }
00452
00453
00454 else if (nMin > BetaInverse) {
00455 dp = Pmax - Pmin;
00456 ge = CAImax;
00457 }
00458
00459
00460
00461
00462
00463 else {
00464 Pmin = Rindex->GetPhotonEnergy(BetaInverse);
00465 dp = Pmax - Pmin;
00466
00467
00468 bool isOutRange;
00469 double CAImin = chAngleIntegrals_->GetValue(Pmin, isOutRange);
00470 ge = CAImax - CAImin;
00471
00472 }
00473
00474
00475 double numPhotons = rFact * charge/eplus * charge/eplus *
00476 (dp - ge * BetaInverse*BetaInverse);
00477
00478 LogDebug("EcalSim") << "@SUB=getAverageNumberOfPhotons"
00479 << "CAImin = " << CAImin << "\n"
00480 << "CAImax = " << CAImax << "\n"
00481 << "dp = " << dp << ", ge = " << ge << "\n"
00482 << "numPhotons = " << numPhotons;
00483
00484
00485
00486 return numPhotons;
00487
00488 }
00489
00490
00491
00492
00493
00494 bool DreamSD::setPbWO2MaterialProperties_( G4Material* aMaterial ) {
00495
00496 std::string pbWO2Name("E_PbWO4");
00497 if ( pbWO2Name != aMaterial->GetName() ) {
00498 edm::LogWarning("EcalSim") << "This is not the right material: "
00499 << "expecting " << pbWO2Name
00500 << ", got " << aMaterial->GetName();
00501 return false;
00502 }
00503
00504 G4MaterialPropertiesTable* table = new G4MaterialPropertiesTable();
00505
00506
00507
00508 const int nEntries = 14;
00509 double PhotonEnergy[nEntries] = { 1.7712*eV, 1.8368*eV, 1.90745*eV, 1.98375*eV, 2.0664*eV,
00510 2.15625*eV, 2.25426*eV, 2.3616*eV, 2.47968*eV, 2.61019*eV,
00511 2.75521*eV, 2.91728*eV, 3.09961*eV, 3.30625*eV };
00512 double RefractiveIndex[nEntries] = { 2.17728, 2.18025, 2.18357, 2.18753, 2.19285,
00513 2.19813, 2.20441, 2.21337, 2.22328, 2.23619,
00514 2.25203, 2.27381, 2.30282, 2.34666 };
00515
00516 table->AddProperty( "RINDEX", PhotonEnergy, RefractiveIndex, nEntries );
00517 aMaterial->SetMaterialPropertiesTable(table);
00518
00519
00520
00521 chAngleIntegrals_ =
00522 std::auto_ptr<G4PhysicsOrderedFreeVector>( new G4PhysicsOrderedFreeVector() );
00523
00524 int index = 0;
00525 double currentRI = RefractiveIndex[index];
00526 double currentPM = PhotonEnergy[index];
00527 double currentCAI = 0.0;
00528 chAngleIntegrals_->InsertValues(currentPM, currentCAI);
00529 double prevPM = currentPM;
00530 double prevCAI = currentCAI;
00531 double prevRI = currentRI;
00532 while ( ++index < nEntries ) {
00533 currentRI = RefractiveIndex[index];
00534 currentPM = PhotonEnergy[index];
00535 currentCAI = 0.5*(1.0/(prevRI*prevRI) + 1.0/(currentRI*currentRI));
00536 currentCAI = prevCAI + (currentPM - prevPM) * currentCAI;
00537
00538 chAngleIntegrals_->InsertValues(currentPM, currentCAI);
00539
00540 prevPM = currentPM;
00541 prevCAI = currentCAI;
00542 prevRI = currentRI;
00543 }
00544
00545 LogDebug("EcalSim") << "Material properties set for " << aMaterial->GetName();
00546
00547 return true;
00548
00549 }
00550
00551
00552
00553
00554
00555
00556
00557 const double DreamSD::getPhotonEnergyDeposit_( const G4ThreeVector& p,
00558 const G4ThreeVector& x,
00559 const G4Step* aStep )
00560 {
00561
00562 double energy = 0;
00563
00564
00565
00566
00567
00568
00569
00570
00571
00572 G4StepPoint* stepPoint = aStep->GetPreStepPoint();
00573 G4LogicalVolume* lv = stepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
00574 G4String nameVolume= lv->GetName();
00575
00576 double crlength = crystalLength(lv);
00577 double crwidth = crystalWidth(lv);
00578 double dapd = 0.5 * crlength - x.x();
00579 double y = p.y()/p.x()*dapd;
00580
00581 LogDebug("EcalSim") << "Distance to APD: " << dapd
00582 << " - y at APD: " << y;
00583
00584
00585 if ( fabs(y)>crwidth*0.5 ) {
00586
00587 }
00588
00589
00590 double waveLength = p.mag()*1.239e8;
00591
00592
00593 energy = p.mag()*PMTResponse::getEfficiency(waveLength);
00594
00595 LogDebug("EcalSim") << "Wavelength: " << waveLength << " - Energy: " << energy;
00596
00597 return energy;
00598
00599 }