9 #include "G4LogicalVolumeStore.hh"
10 #include "G4LogicalVolume.hh"
13 #include "G4VProcess.hh"
14 #include "G4Poisson.hh"
25 #include "G4SystemOfUnits.hh"
26 #include "G4PhysicalConstants.hh"
32 CaloSD(name, cpv, clg, p, manager) {
43 edm::LogInfo(
"EcalSim") <<
"Constructing a DreamSD with name " << GetName() <<
"\n"
44 <<
"DreamSD:: Use of Birks law is set to "
45 << useBirk <<
" with three constants kB = "
46 << birk1 <<
", C1 = " << birk2 <<
", C2 = "
48 <<
" Slope for Light yield is set to "
50 <<
" Parameterization of Cherenkov is set to "
51 << doCherenkov_ <<
" and readout both sides is "
61 <<
"please add it to config file";
63 ntuple_ = tfile->
make<TTree>(
"tree",
"Cherenkov photons");
103 G4String nameVolume =
preStepPoint->GetPhysicalVolume()->GetName();
110 LogDebug(
"EcalSim") <<
"DreamSD:: " << nameVolume <<
" Side " <<
side
111 <<
" Light Collection Efficiency " << weight
112 <<
" Weighted Energy Deposit " <<
edepositEM/MeV
121 double time = (aStep->GetPostStepPoint()->GetGlobalTime())/nanosecond;
123 if (
side < 0) unitID++;
132 if (primaryID == 0) {
133 edm::LogWarning(
"EcalSim") <<
"CaloSD: Problem with primaryID **** set by "
134 <<
"force to TkID **** "
135 <<
theTrack->GetTrackID() <<
" in Volume "
136 <<
preStepPoint->GetTouchable()->GetVolume(0)->GetName();
140 bool flag = (unitID > 0);
141 G4TouchableHistory* touch =(G4TouchableHistory*)(
theTrack->GetTouchable());
145 LogDebug(
"EcalSim") <<
"CaloSD:: GetStepInfo for"
146 <<
" PV " << touch->GetVolume(0)->GetName()
147 <<
" PVid = " << touch->GetReplicaNumber(0)
148 <<
" MVid = " << touch->GetReplicaNumber(1)
152 LogDebug(
"EcalSim") <<
"CaloSD:: GetStepInfo for"
153 <<
" PV " << touch->GetVolume(0)->GetName()
154 <<
" PVid = " << touch->GetReplicaNumber(0)
155 <<
" MVid = " << touch->GetReplicaNumber(1)
156 <<
" Unit " << std::hex << unitID << std::dec
168 DimensionMap::const_iterator ite =
xtalLMap.begin();
169 G4LogicalVolume* lv = (ite->first);
170 G4Material* material = lv->GetMaterial();
171 edm::LogInfo(
"EcalSim") <<
"DreamSD::initRun: Initializes for material "
172 << material->GetName() <<
" in " << lv->GetName();
176 edm::LogWarning(
"EcalSim") <<
"Couldn't retrieve material properties table\n"
177 <<
" Material = " << material->GetName();
186 const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
187 uint32_t
id = (touch->GetReplicaNumber(1))*10 + (touch->GetReplicaNumber(0));
188 LogDebug(
"EcalSim") <<
"DreamSD:: ID " << id;
196 G4String attribute =
"ReadOutName";
204 const G4LogicalVolumeStore * lvs = G4LogicalVolumeStore::GetInstance();
205 std::vector<G4LogicalVolume *>::const_iterator lvcite;
211 G4LogicalVolume* lv=0;
212 for (lvcite = lvs->begin(); lvcite != lvs->end(); lvcite++)
213 if ((*lvcite)->GetName() ==
name) {
217 LogDebug(
"EcalSim") <<
"DreamSD::initMap (for " << sd <<
"): Solid "
218 << name <<
" Shape " << sol.
shape() <<
" Parameter 0 = "
219 << paras[0] <<
" Logical Volume " << lv;
220 double length = 0,
width = 0;
223 length = 2.0*paras.back();
224 width = 2.0*paras.front();
228 LogDebug(
"EcalSim") <<
"DreamSD: Length Table for " << attribute <<
" = "
230 DimensionMap::const_iterator ite =
xtalLMap.begin();
232 for (; ite !=
xtalLMap.end(); ite++, i++) {
233 G4String
name =
"Unknown";
234 if (ite->first != 0) name = (ite->first)->GetName();
235 LogDebug(
"EcalSim") <<
" " << i <<
" " << ite->first <<
" " << name
236 <<
" L = " << ite->second.first
237 <<
" W = " << ite->second.second;
244 G4StepPoint* stepPoint = aStep->GetPreStepPoint();
245 G4LogicalVolume* lv = stepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
246 G4String nameVolume= lv->GetName();
249 G4ThreeVector localPoint =
setToLocal(stepPoint->GetPosition(),
250 stepPoint->GetTouchable());
252 double localz = localPoint.x();
253 double dapd = 0.5 * crlength - flag*localz;
254 if (dapd >= -0.1 || dapd <= crlength+0.1) {
258 edm::LogWarning(
"EcalSim") <<
"DreamSD: light coll curve : wrong distance "
259 <<
"to APD " << dapd <<
" crlength = "
260 << crlength <<
" crystal name = " << nameVolume
261 <<
" z of localPoint = " << localz
262 <<
" take weight = " <<
weight;
264 LogDebug(
"EcalSim") <<
"DreamSD, light coll curve : " << dapd
265 <<
" crlength = " << crlength
266 <<
" crystal name = " << nameVolume
267 <<
" z of localPoint = " << localz
268 <<
" take weight = " <<
weight;
276 DimensionMap::const_iterator ite =
xtalLMap.find(lv);
277 if (ite !=
xtalLMap.end()) length = ite->second.first;
286 DimensionMap::const_iterator ite =
xtalLMap.find(lv);
287 if (ite !=
xtalLMap.end()) width = ite->second.second;
298 double cherenkovEnergy = 0;
300 G4Material* material = aStep->GetTrack()->GetMaterial();
304 if ( Rindex ==
NULL ) {
306 return cherenkovEnergy;
311 int Rlength = Rindex->GetVectorLength() - 1;
312 double Pmin = Rindex->Energy(0);
313 double Pmax = Rindex->Energy(Rlength);
314 LogDebug(
"EcalSim") <<
"Material properties: " <<
"\n"
315 <<
" Pmin = " << Pmin
316 <<
" Pmax = " << Pmax;
319 G4StepPoint* pPreStepPoint = aStep->GetPreStepPoint();
320 G4StepPoint* pPostStepPoint = aStep->GetPostStepPoint();
321 G4ThreeVector x0 = pPreStepPoint->GetPosition();
322 G4ThreeVector p0 = aStep->GetDeltaPosition().unit();
323 const G4DynamicParticle* aParticle = aStep->GetTrack()->GetDynamicParticle();
324 const double charge = aParticle->GetDefinition()->GetPDGCharge();
326 double beta = 0.5*( pPreStepPoint->GetBeta() + pPostStepPoint->GetBeta() );
327 double BetaInverse = 1.0/
beta;
329 LogDebug(
"EcalSim") <<
"Particle properties: " <<
"\n"
330 <<
" charge = " << charge
331 <<
" beta = " <<
beta;
335 if ( meanNumberOfPhotons <= 0.0 ) {
336 LogDebug(
"EcalSim") <<
"Mean number of photons is zero: " << meanNumberOfPhotons
337 <<
", stopping here";
338 return cherenkovEnergy;
342 meanNumberOfPhotons *= aStep->GetStepLength();
345 int numPhotons =
static_cast<int>( G4Poisson(meanNumberOfPhotons) );
347 if ( numPhotons <= 0 ) {
348 LogDebug(
"EcalSim") <<
"Poission number of photons is zero: " << numPhotons
349 <<
", stopping here";
350 return cherenkovEnergy;
354 double dp = Pmax - Pmin;
355 double maxCos = BetaInverse / (*Rindex)[Rlength];
356 double maxSin2 = (1.0 - maxCos) * (1.0 + maxCos);
359 for (
int iPhoton = 0; iPhoton<numPhotons; ++iPhoton ) {
363 double sampledMomentum, sampledRI;
364 double cosTheta, sin2Theta;
368 randomNumber = G4UniformRand();
369 sampledMomentum = Pmin + randomNumber *
dp;
370 sampledRI = Rindex->Value(sampledMomentum);
371 cosTheta = BetaInverse / sampledRI;
373 sin2Theta = (1.0 - cosTheta)*(1.0 + cosTheta);
374 randomNumber = G4UniformRand();
376 }
while (randomNumber*maxSin2 > sin2Theta);
380 randomNumber = G4UniformRand();
382 double phi = twopi*randomNumber;
383 double sinPhi =
sin(phi);
384 double cosPhi =
cos(phi);
389 double sinTheta =
sqrt(sin2Theta);
390 double px = sinTheta*cosPhi;
391 double py = sinTheta*sinPhi;
392 double pz = cosTheta;
393 G4ThreeVector photonDirection(px, py, pz);
396 photonDirection.rotateUz(p0);
399 randomNumber = G4UniformRand();
400 G4ThreeVector photonPosition = x0 + randomNumber * aStep->GetDeltaPosition();
401 G4ThreeVector photonMomentum = sampledMomentum*photonDirection;
408 px_[iPhoton] = photonMomentum.x();
409 py_[iPhoton] = photonMomentum.y();
410 pz_[iPhoton] = photonMomentum.z();
411 x_[iPhoton] = photonPosition.x();
412 y_[iPhoton] = photonPosition.y();
413 z_[iPhoton] = photonPosition.z();
420 return cherenkovEnergy;
430 const G4Material* aMaterial,
431 G4MaterialPropertyVector* Rindex )
433 const G4double rFact = 369.81/(eV * cm);
435 if( beta <= 0.0 )
return 0.0;
437 double BetaInverse = 1./
beta;
444 int Rlength = Rindex->GetVectorLength() - 1;
445 double Pmin = Rindex->Energy(0);
446 double Pmax = Rindex->Energy(Rlength);
449 double nMin = (*Rindex)[0];
450 double nMax = (*Rindex)[Rlength];
455 double dp = 0., ge = 0., CAImin = 0.;
458 if ( nMax < BetaInverse) { }
461 else if (nMin > BetaInverse) {
471 Pmin = Rindex->Value(BetaInverse);
476 ge = CAImax - CAImin;
481 double numPhotons = rFact * charge/eplus * charge/eplus *
482 (dp - ge * BetaInverse*BetaInverse);
484 LogDebug(
"EcalSim") <<
"@SUB=getAverageNumberOfPhotons"
485 <<
"CAImin = " << CAImin <<
"\n"
486 <<
"CAImax = " << CAImax <<
"\n"
487 <<
"dp = " << dp <<
", ge = " << ge <<
"\n"
488 <<
"numPhotons = " << numPhotons;
503 if ( pbWO2Name != aMaterial->GetName() ) {
505 <<
"expecting " << pbWO2Name
506 <<
", got " << aMaterial->GetName();
510 G4MaterialPropertiesTable*
table =
new G4MaterialPropertiesTable();
514 const int nEntries = 14;
515 double PhotonEnergy[nEntries] = { 1.7712*eV, 1.8368*eV, 1.90745*eV, 1.98375*eV, 2.0664*eV,
516 2.15625*eV, 2.25426*eV, 2.3616*eV, 2.47968*eV, 2.61019*eV,
517 2.75521*eV, 2.91728*eV, 3.09961*eV, 3.30625*eV };
518 double RefractiveIndex[nEntries] = { 2.17728, 2.18025, 2.18357, 2.18753, 2.19285,
519 2.19813, 2.20441, 2.21337, 2.22328, 2.23619,
520 2.25203, 2.27381, 2.30282, 2.34666 };
522 table->AddProperty(
"RINDEX", PhotonEnergy, RefractiveIndex, nEntries );
523 aMaterial->SetMaterialPropertiesTable(table);
528 std::auto_ptr<G4PhysicsOrderedFreeVector>(
new G4PhysicsOrderedFreeVector() );
531 double currentRI = RefractiveIndex[
index];
532 double currentPM = PhotonEnergy[
index];
533 double currentCAI = 0.0;
535 double prevPM = currentPM;
536 double prevCAI = currentCAI;
537 double prevRI = currentRI;
538 while ( ++index < nEntries ) {
539 currentRI = RefractiveIndex[
index];
540 currentPM = PhotonEnergy[
index];
541 currentCAI = 0.5*(1.0/(prevRI*prevRI) + 1.0/(currentRI*currentRI));
542 currentCAI = prevCAI + (currentPM - prevPM) * currentCAI;
547 prevCAI = currentCAI;
551 LogDebug(
"EcalSim") <<
"Material properties set for " << aMaterial->GetName();
564 const G4ThreeVector&
x,
565 const G4Step* aStep )
578 G4StepPoint* stepPoint = aStep->GetPreStepPoint();
579 G4LogicalVolume* lv = stepPoint->GetTouchable()->GetVolume(0)->GetLogicalVolume();
580 G4String nameVolume= lv->GetName();
584 double dapd = 0.5 * crlength - x.x();
585 double y = p.y()/p.x()*dapd;
587 LogDebug(
"EcalSim") <<
"Distance to APD: " << dapd
588 <<
" - y at APD: " <<
y;
591 if ( fabs(y)>crwidth*0.5 ) {
596 double waveLength = p.mag()*1.239e8;
601 LogDebug(
"EcalSim") <<
"Wavelength: " << waveLength <<
" - Energy: " <<
energy;
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
const std::vector< double > & parameters(void) const
Give the parameters of the solid.
const DDLogicalPart & logicalPart() const
The logical-part of the current node in the filtered-view.
void addFilter(const DDFilter &, log_op op=AND)
virtual G4bool getStepInfo(G4Step *aStep)
double getPhotonEnergyDeposit_(const G4ParticleMomentum &p, const G4ThreeVector &x, const G4Step *aStep)
Returns energy deposit for a given photon.
std::pair< double, double > Doubles
virtual bool ProcessHits(G4Step *step, G4TouchableHistory *tHistory)
Sin< T >::type sin(const T &t)
double getAverageNumberOfPhotons_(const double charge, const double beta, const G4Material *aMaterial, G4MaterialPropertyVector *rIndex)
Returns average number of photons created by track.
bool setPbWO2MaterialProperties_(G4Material *aMaterial)
Sets material properties at run-time...
T * make(const Args &...args) const
make new ROOT object
std::auto_ptr< G4PhysicsOrderedFreeVector > chAngleIntegrals_
Table of Cherenkov angle integrals vs photon momentum.
const DDSolid & solid(void) const
Returns a reference object of the solid being the shape of this LogicalPart.
virtual uint32_t setDetUnitId(G4Step *)
type of data representation of DDCompactView
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
A DDSolid represents the shape of a part.
bool next()
set current node to the next node in the filtered tree
Cos< T >::type cos(const T &t)
DreamSD(G4String, const DDCompactView &, SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
DDSolidShape shape(void) const
The type of the solid.
double getAttenuation(G4Step *aStep, double birk1, double birk2, double birk3)
G4MaterialPropertiesTable * materialPropertiesTable
const double crystalWidth(G4LogicalVolume *) const
void setID(uint32_t unitID, double timeSlice, int trackID, uint16_t depth=0)
G4StepPoint * preStepPoint
double cherenkovDeposit_(G4Step *aStep)
Returns the total energy due to Cherenkov radiation.
bool firstChild()
set the current node to the first child ...
void setCriteria(const DDValue &nameVal, comp_op, log_op l=AND, bool asString=true, bool merged=true)
const double crystalLength(G4LogicalVolume *) const
void initMap(G4String, const DDCompactView &)
const std::string & name() const
Returns the name.
double curve_LY(G4Step *, int)
static const double getEfficiency(const double &waveLengthNm)
Return efficiency for given photon wavelength (in nm)
CaloG4Hit * createNewHit()
The DDGenericFilter is a runtime-parametrized Filter looking on DDSpecifcs.
G4ThreeVector setToLocal(const G4ThreeVector &, const G4VTouchable *)