/*
<img src="picts/AliPMDv1Class.gif">
*/
//End_Html
#include <Riostream.h>
#include <TGeoManager.h>
#include <TGeoGlobalMagField.h>
#include <TVirtualMC.h>
#include "AliConst.h"
#include "AliLog.h"
#include "AliMC.h"
#include "AliMagF.h"
#include "AliPMDv2008.h"
#include "AliRun.h"
const Int_t AliPMDv2008::fgkNcolUM1 = 48;
const Int_t AliPMDv2008::fgkNcolUM2 = 96;
const Int_t AliPMDv2008::fgkNrowUM1 = 96;
const Int_t AliPMDv2008::fgkNrowUM2 = 48;
const Float_t AliPMDv2008::fgkCellRadius = 0.25;
const Float_t AliPMDv2008::fgkCellWall = 0.02;
const Float_t AliPMDv2008::fgkCellDepth = 0.50;
const Float_t AliPMDv2008::fgkThBase = 0.2;
const Float_t AliPMDv2008::fgkThBKP = 0.1;
const Float_t AliPMDv2008::fgkThAir = 1.03;
const Float_t AliPMDv2008::fgkThPCB = 0.16;
const Float_t AliPMDv2008::fgkThLead = 1.5;
const Float_t AliPMDv2008::fgkThSteel = 0.5;
const Float_t AliPMDv2008::fgkGap = 0.025;
const Float_t AliPMDv2008::fgkZdist = 361.5;
const Float_t AliPMDv2008::fgkSqroot3 = 1.7320508;
const Float_t AliPMDv2008::fgkSqroot3by2 = 0.8660254;
const Float_t AliPMDv2008::fgkSSBoundary = 0.3;
const Float_t AliPMDv2008::fgkThSS = 1.03;
const Float_t AliPMDv2008::fgkThG10 = 1.03;
ClassImp(AliPMDv2008)
AliPMDv2008::AliPMDv2008():
fSMthick(0.),
fDthick(0.),
fSMLengthax(0.),
fSMLengthay(0.),
fSMLengthbx(0.),
fSMLengthby(0.),
fMedSens(0)
{
for (Int_t i = 0; i < 3; i++)
{
fDboxmm1[i] = 0.;
fDboxmm12[i] = 0.;
fDboxmm2[i] = 0.;
fDboxmm22[i] = 0.;
}
}
AliPMDv2008::AliPMDv2008(const char *name, const char *title):
AliPMD(name,title),
fSMthick(0.),
fDthick(0.),
fSMLengthax(0.),
fSMLengthay(0.),
fSMLengthbx(0.),
fSMLengthby(0.),
fMedSens(0)
{
for (Int_t i = 0; i < 3; i++)
{
fDboxmm1[i] = 0.;
fDboxmm12[i] = 0.;
fDboxmm2[i] = 0.;
fDboxmm22[i] = 0.;
}
}
void AliPMDv2008::CreateGeometry()
{
GetParameters();
CreateSupermodule();
CreatePMD();
}
void AliPMDv2008::CreateSupermodule()
{
Int_t i,j;
Int_t number;
Int_t ihrotm,irotdm;
Float_t xb, yb, zb;
Int_t *idtmed = fIdtmed->GetArray()-599;
AliMatrix(ihrotm, 90., 30., 90., 120., 0., 0.);
AliMatrix(irotdm, 90., 180., 90., 270., 180., 0.);
Float_t hexd2[10] = {0.,360.,6,2,-0.25,0.,0.23,0.25,0.,0.23};
hexd2[4] = -fgkCellDepth/2.;
hexd2[7] = fgkCellDepth/2.;
hexd2[6] = fgkCellRadius - fgkCellWall;
hexd2[9] = fgkCellRadius - fgkCellWall;
TVirtualMC::GetMC()->Gsvolu("ECAR", "PGON", idtmed[604], hexd2,10);
Float_t hexd1[10] = {0.,360.,6,2,-0.25,0.,0.25,0.25,0.,0.25};
hexd1[4] = -fgkCellDepth/2.;
hexd1[7] = fgkCellDepth/2.;
hexd1[6] = fgkCellRadius;
hexd1[9] = fgkCellRadius;
TVirtualMC::GetMC()->Gsvolu("ECCU", "PGON", idtmed[614], hexd1,10);
TVirtualMC::GetMC()->Gspos("ECAR", 1, "ECCU", 0., 0., 0., 0, "ONLY");
Float_t dbox1[3];
dbox1[0] = fgkCellRadius/fgkSqroot3by2;
dbox1[1] = fgkNrowUM1*fgkCellRadius;
dbox1[2] = fgkCellDepth/2.;
TVirtualMC::GetMC()->Gsvolu("EST1","BOX", idtmed[698], dbox1, 3);
Float_t dbox2[3];
dbox2[1] = fgkNrowUM2*fgkCellRadius;
dbox2[0] = dbox1[0];
dbox2[2] = dbox1[2];
TVirtualMC::GetMC()->Gsvolu("EST2","BOX", idtmed[698], dbox2, 3);
xb = 0.;
zb = 0.;
yb = (dbox1[1]) - fgkCellRadius;
for (i = 1; i <= fgkNrowUM1; ++i)
{
number = i;
TVirtualMC::GetMC()->Gspos("ECCU", number, "EST1", xb,yb,zb, 0, "ONLY");
yb -= (fgkCellRadius*2.);
}
xb = 0.;
zb = 0.;
yb = (dbox2[1]) - fgkCellRadius;
for (i = 1; i <= fgkNrowUM2; ++i)
{
number = i;
TVirtualMC::GetMC()->Gspos("ECCU", number, "EST2", xb,yb,zb, 0, "ONLY");
yb -= (fgkCellRadius*2.);
}
Float_t dbox3[3];
dbox3[0] = (dbox1[0]*fgkNcolUM1)-(fgkCellRadius*fgkSqroot3*(fgkNcolUM1-1)/6.);
dbox3[1] = dbox1[1]+fgkCellRadius/2.;
dbox3[2] = fgkCellDepth/2.;
TVirtualMC::GetMC()->Gsvolu("EHC1","BOX", idtmed[698], dbox3, 3);
xb = dbox3[0]-dbox1[0];
for (j = 1; j <= fgkNcolUM1; ++j)
{
if(j%2 == 0)
{
yb = -fgkCellRadius/2.0;
}
else
{
yb = fgkCellRadius/2.0;
}
number = j;
TVirtualMC::GetMC()->Gspos("EST1",number, "EHC1", xb, yb , 0. , 0, "MANY");
xb = (dbox3[0]-dbox1[0])-j*fgkCellRadius*fgkSqroot3;
}
dbox3[0] = (dbox1[0]*fgkNcolUM1)-(fgkCellRadius*fgkSqroot3*(fgkNcolUM1-1)/6.);
dbox3[1] = dbox1[1]+fgkCellRadius/2.;
dbox3[2] = fgkCellDepth/2.;
Float_t dbox4[3];
dbox4[0] =(dbox2[0]*fgkNcolUM2)-(fgkCellRadius*fgkSqroot3*(fgkNcolUM2-1)/6.);
dbox4[1] = dbox2[1] + fgkCellRadius/2.;
dbox4[2] = dbox3[2];
TVirtualMC::GetMC()->Gsvolu("EHC2","BOX", idtmed[698], dbox4, 3);
xb = dbox4[0]-dbox2[0];
for (j = 1; j <= fgkNcolUM2; ++j)
{
if(j%2 == 0)
{
yb = -fgkCellRadius/2.0;
}
else
{
yb = +fgkCellRadius/2.0;
}
number = j;
TVirtualMC::GetMC()->Gspos("EST2",number, "EHC2", xb, yb , 0. ,0, "MANY");
xb = (dbox4[0]-dbox2[0])-j*fgkCellRadius*fgkSqroot3;
}
Float_t dboxPcbA[3];
dboxPcbA[0] = dbox3[0];
dboxPcbA[1] = dbox3[1];
dboxPcbA[2] = fgkThPCB/2.;
TVirtualMC::GetMC()->Gsvolu("EPCA","BOX", idtmed[607], dboxPcbA, 3);
Float_t dboxBPlaneA[3];
dboxBPlaneA[0] = dbox3[0];
dboxBPlaneA[1] = dbox3[1];
dboxBPlaneA[2] = fgkThBKP/2.;
TVirtualMC::GetMC()->Gsvolu("EBKA","BOX", idtmed[607], dboxBPlaneA, 3);
Float_t dboxAir3A[3];
dboxAir3A[0] = dbox3[0]+(2.0*fgkGap);
dboxAir3A[1] = dbox3[1]+(2.0*fgkGap);
dboxAir3A[2] = fgkThAir/2.;
TVirtualMC::GetMC()->Gsvolu("ECGA","BOX", idtmed[698], dboxAir3A, 3);
TVirtualMC::GetMC()->Gsvolu("ECVA","BOX", idtmed[698], dboxAir3A, 3);
Float_t dboxGGA[3];
dboxGGA[0] = dboxAir3A[0]+(2.0*fgkGap);
dboxGGA[1] = dboxAir3A[1]+(2.0*fgkGap);
dboxGGA[2] = fgkThG10/2.;
TVirtualMC::GetMC()->Gsvolu("EDGA","BOX", idtmed[607], dboxGGA, 3);
TVirtualMC::GetMC()->Gsvolu("EDVA","BOX", idtmed[607], dboxGGA, 3);
Float_t dboxSS1[3];
dboxSS1[0] = dboxGGA[0]+fgkSSBoundary;
dboxSS1[1] = dboxGGA[1]+fgkSSBoundary;
dboxSS1[2] = fgkThSS/2.;
TVirtualMC::GetMC()->Gsvolu("ESSA","BOX", idtmed[618], dboxSS1, 3);
TVirtualMC::GetMC()->Gsvolu("ESVA","BOX", idtmed[618], dboxSS1, 3);
TVirtualMC::GetMC()->Gspos("EDVA", 1, "ESVA", 0., 0., 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("ECVA", 1, "EDVA", 0., 0., 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EDGA", 1, "ESSA", 0., 0., 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("ECGA", 1, "EDGA", 0., 0., 0., 0, "ONLY");
Float_t zbpcb = -dboxAir3A[2] + (2.0*fgkGap) + fgkThPCB/2.;
TVirtualMC::GetMC()->Gspos("EPCA", 1, "ECVA", 0., 0., zbpcb, 0, "ONLY");
Float_t zhc = zbpcb + fgkThPCB/2. + fgkCellDepth/2.;
TVirtualMC::GetMC()->Gspos("EHC1", 1, "ECVA", 0., 0., zhc, 0, "ONLY");
Float_t ztpcb = zhc + fgkCellDepth/2 + fgkThPCB/2.;
TVirtualMC::GetMC()->Gspos("EPCA", 2, "ECVA", 0., 0., ztpcb, 0, "ONLY");
zbpcb = -dboxAir3A[2] + fgkThPCB + fgkThPCB/2.;
TVirtualMC::GetMC()->Gspos("EPCA", 1, "ECGA", 0., 0., zbpcb, 0, "ONLY");
zhc = zbpcb + fgkThPCB/2. + fgkCellDepth/2.;
TVirtualMC::GetMC()->Gspos("EHC1", 1, "ECGA", 0., 0., zhc, 0, "ONLY");
ztpcb = zhc + fgkCellDepth/2 + fgkThPCB/2.;
TVirtualMC::GetMC()->Gspos("EPCA", 2, "ECGA", 0., 0., ztpcb, 0, "ONLY");
Float_t dboxUM1[3];
dboxUM1[0] = dboxSS1[0];
dboxUM1[1] = dboxSS1[1];
dboxUM1[2] = fgkThSS/2. +0.15;
TVirtualMC::GetMC()->Gsvolu("EUM1","BOX", idtmed[698], dboxUM1, 3);
TVirtualMC::GetMC()->Gsvolu("EUV1","BOX", idtmed[698], dboxUM1, 3);
Float_t dboxBaseA[3];
dboxBaseA[0] = dboxSS1[0];
dboxBaseA[1] = dboxSS1[1];
dboxBaseA[2] = fgkThBase/2.;
TVirtualMC::GetMC()->Gsvolu("EBPA","BOX", idtmed[607], dboxBaseA, 3);
Float_t zbaseplate = -dboxUM1[2] + fgkThBase/2.;
TVirtualMC::GetMC()->Gspos("EBPA", 1, "EUV1", 0., 0., zbaseplate, 0, "ONLY");
Float_t zss = zbaseplate + fgkThBase/2. + fgkThSS/2.;
TVirtualMC::GetMC()->Gspos("ESVA", 1, "EUV1", 0., 0., zss, 0, "ONLY");
Float_t zbkp = zss + fgkThSS/2. + fgkThBKP/2.;
TVirtualMC::GetMC()->Gspos("EBKA", 1, "EUV1", 0., 0., zbkp, 0, "ONLY");
zbkp = -dboxUM1[2] + fgkThBKP/2.;
TVirtualMC::GetMC()->Gspos("EBKA", 1, "EUM1", 0., 0., zbkp, 0, "ONLY");
zss = zbkp + fgkThBKP/2. + fgkThSS/2.;
TVirtualMC::GetMC()->Gspos("ESSA", 1, "EUM1", 0., 0., zss, 0, "ONLY");
zbaseplate = zss + fgkThSS/2 + fgkThBase/2.;
TVirtualMC::GetMC()->Gspos("EBPA", 1, "EUM1", 0., 0., zbaseplate, 0, "ONLY");
Float_t dboxPcbB[3];
dboxPcbB[0] = dbox4[0];
dboxPcbB[1] = dbox4[1];
dboxPcbB[2] = fgkThPCB/2.;
TVirtualMC::GetMC()->Gsvolu("EPCB","BOX", idtmed[607], dboxPcbB, 3);
Float_t dboxBPlaneB[3];
dboxBPlaneB[0] = dbox4[0];
dboxBPlaneB[1] = dbox4[1];
dboxBPlaneB[2] = fgkThBKP/2.;
TVirtualMC::GetMC()->Gsvolu("EBKB","BOX", idtmed[607], dboxBPlaneB, 3);
Float_t dboxAir3B[3];
dboxAir3B[0] = dbox4[0]+(2.0*fgkGap);
dboxAir3B[1] = dbox4[1]+(2.0*fgkGap);
dboxAir3B[2] = fgkThAir/2.;
TVirtualMC::GetMC()->Gsvolu("ECGB","BOX", idtmed[698], dboxAir3B, 3);
TVirtualMC::GetMC()->Gsvolu("ECVB","BOX", idtmed[698], dboxAir3B, 3);
Float_t dboxGGB[3];
dboxGGB[0] = dboxAir3B[0]+(2.0*fgkGap);
dboxGGB[1] = dboxAir3B[1]+(2.0*fgkGap);
dboxGGB[2] = fgkThG10/2.;
TVirtualMC::GetMC()->Gsvolu("EDGB","BOX", idtmed[607], dboxGGB, 3);
TVirtualMC::GetMC()->Gsvolu("EDVB","BOX", idtmed[607], dboxGGB, 3);
Float_t dboxSS2[3];
dboxSS2[0] = dboxGGB[0] + fgkSSBoundary;
dboxSS2[1] = dboxGGB[1] + fgkSSBoundary;
dboxSS2[2] = fgkThSS/2.;
TVirtualMC::GetMC()->Gsvolu("ESSB","BOX", idtmed[618], dboxSS2, 3);
TVirtualMC::GetMC()->Gsvolu("ESVB","BOX", idtmed[618], dboxSS2, 3);
TVirtualMC::GetMC()->Gspos("EDGB", 1, "ESSB", 0., 0., 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("ECGB", 1, "EDGB", 0., 0., 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EDVB", 1, "ESVB", 0., 0., 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("ECVB", 1, "EDVB", 0., 0., 0., 0, "ONLY");
Float_t zbpcb2 = -dboxAir3B[2] + (2.0*fgkGap) + fgkThPCB/2.;
TVirtualMC::GetMC()->Gspos("EPCB", 1, "ECVB", 0., 0., zbpcb2, 0, "ONLY");
Float_t zhc2 = zbpcb2 + fgkThPCB/2. + fgkCellDepth/2.;
TVirtualMC::GetMC()->Gspos("EHC2", 1, "ECVB", 0., 0., zhc2, 0, "ONLY");
Float_t ztpcb2 = zhc2 + fgkCellDepth/2 + fgkThPCB/2.;
TVirtualMC::GetMC()->Gspos("EPCB", 2, "ECVB", 0., 0., ztpcb2, 0, "ONLY");
zbpcb2 = -dboxAir3B[2] + fgkThPCB + fgkThPCB/2.;
TVirtualMC::GetMC()->Gspos("EPCB", 1, "ECGB", 0., 0., zbpcb2, 0, "ONLY");
zhc2 = zbpcb2 + fgkThPCB/2. + fgkCellDepth/2.;
TVirtualMC::GetMC()->Gspos("EHC2", 1, "ECGB", 0., 0., zhc2, 0, "ONLY");
ztpcb2 = zhc2 + fgkCellDepth/2 + fgkThPCB/2.;
TVirtualMC::GetMC()->Gspos("EPCB", 2, "ECGB", 0., 0., ztpcb2, 0, "ONLY");
Float_t dboxUM2[3];
dboxUM2[0] = dboxSS2[0];
dboxUM2[1] = dboxSS2[1];
dboxUM2[2] = fgkThSS/2. +0.15;
TVirtualMC::GetMC()->Gsvolu("EUM2","BOX", idtmed[698], dboxUM2, 3);
TVirtualMC::GetMC()->Gsvolu("EUV2","BOX", idtmed[698], dboxUM2, 3);
Float_t dboxBaseB[3];
dboxBaseB[0] = dboxSS2[0];
dboxBaseB[1] = dboxSS2[1];
dboxBaseB[2] = fgkThBase/2.;
TVirtualMC::GetMC()->Gsvolu("EBPB","BOX", idtmed[607], dboxBaseB, 3);
Float_t zbaseplate2 = -dboxUM2[2] + fgkThBase/2.;
TVirtualMC::GetMC()->Gspos("EBPB", 1, "EUV2", 0., 0., zbaseplate2, 0, "ONLY");
Float_t zss2 = zbaseplate2 + fgkThBase/2. + fgkThSS/2.;
TVirtualMC::GetMC()->Gspos("ESVB", 1, "EUV2", 0., 0., zss2, 0, "ONLY");
Float_t zbkp2 = zss2 + fgkThSS/2. + fgkThBKP/2.;
TVirtualMC::GetMC()->Gspos("EBKB", 1, "EUV2", 0., 0., zbkp2, 0, "ONLY");
zbkp2 = -dboxUM2[2] + fgkThBKP/2.;
TVirtualMC::GetMC()->Gspos("EBKB", 1, "EUM2", 0., 0., zbkp2, 0, "ONLY");
zss2 = zbkp2 + fgkThBKP/2. + fgkThSS/2.;
TVirtualMC::GetMC()->Gspos("ESSB", 1, "EUM2", 0., 0., zss2, 0, "ONLY");
zbaseplate2 = zss2 + fgkThSS/2 + fgkThBase/2.;
TVirtualMC::GetMC()->Gspos("EBPB", 1, "EUM2", 0., 0., zbaseplate2, 0, "ONLY");
Float_t dboxPba[3];
dboxPba[0] = dboxUM1[0];
dboxPba[1] = dboxUM1[1];
dboxPba[2] = fgkThLead/2.;
TVirtualMC::GetMC()->Gsvolu("EPB1","BOX", idtmed[600], dboxPba, 3);
Float_t dboxPbb[3];
dboxPbb[0] = dboxUM2[0];
dboxPbb[1] = dboxUM2[1];
dboxPbb[2] = fgkThLead/2.;
TVirtualMC::GetMC()->Gsvolu("EPB2","BOX", idtmed[600], dboxPbb, 3);
Float_t dboxSM1[3];
dboxSM1[0] = 3.0*dboxUM1[0] + (2.0*0.075);
dboxSM1[1] = 2.0*dboxUM1[1] + 0.05;
dboxSM1[2] = dboxUM1[2];
TVirtualMC::GetMC()->Gsvolu("ESMA","BOX", idtmed[698], dboxSM1, 3);
TVirtualMC::GetMC()->Gsvolu("EMVA","BOX", idtmed[698], dboxSM1, 3);
Float_t xa1,xa2,xa3,ya1,ya2;
xa1 = dboxSM1[0] - dboxUM1[0];
xa2 = xa1 - dboxUM1[0] - 0.15 - dboxUM1[0];
xa3 = xa2 - dboxUM1[0] - 0.15 - dboxUM1[0];
ya1 = dboxSM1[1] - dboxUM1[1];
ya2 = ya1 - dboxUM1[1] - 0.1 - dboxUM1[1];
TVirtualMC::GetMC()->Gspos("EUM1", 2, "ESMA", xa2, ya1, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EUM1", 3, "ESMA", xa3, ya1, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EUM1", 4, "ESMA", xa1, ya2, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EUM1", 5, "ESMA", xa2, ya2, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EUM1", 6, "ESMA", xa3, ya2, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EUV1", 1, "EMVA", xa1, ya1, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EUV1", 2, "EMVA", xa2, ya1, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EUV1", 3, "EMVA", xa3, ya1, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EUV1", 4, "EMVA", xa1, ya2, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EUV1", 5, "EMVA", xa2, ya2, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EUV1", 6, "EMVA", xa3, ya2, 0., 0, "ONLY");
Float_t dboxSM2[3];
dboxSM2[0] = 2.0*dboxUM2[0] + 0.075;
dboxSM2[1] = 3.0*dboxUM2[1] + (2.0*0.05);
dboxSM2[2] = dboxUM2[2];
TVirtualMC::GetMC()->Gsvolu("ESMB","BOX", idtmed[698], dboxSM2, 3);
TVirtualMC::GetMC()->Gsvolu("EMVB","BOX", idtmed[698], dboxSM2, 3);
Float_t xb1,xb2,yb1,yb2,yb3;
xb1 = dboxSM2[0] - dboxUM2[0];
xb2 = xb1 - dboxUM2[0] - 0.15 - dboxUM2[0];
yb1 = dboxSM2[1] - dboxUM2[1];
yb2 = yb1 - dboxUM2[1] - 0.1 - dboxUM2[1];
yb3 = yb2 - dboxUM2[1] - 0.1 - dboxUM2[1];
TVirtualMC::GetMC()->Gspos("EUM2", 3, "ESMB", xb1, yb2, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EUM2", 4, "ESMB", xb2, yb2, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EUM2", 5, "ESMB", xb1, yb3, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EUM2", 6, "ESMB", xb2, yb3, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EUV2", 1, "EMVB", xb1, yb1, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EUV2", 2, "EMVB", xb2, yb1, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EUV2", 3, "EMVB", xb1, yb2, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EUV2", 4, "EMVB", xb2, yb2, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EUV2", 5, "EMVB", xb1, yb3, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EUV2", 6, "EMVB", xb2, yb3, 0., 0, "ONLY");
Float_t dboxSMPb1[3];
dboxSMPb1[0] = 3.0*dboxUM1[0] + (2.0*0.075);
dboxSMPb1[1] = 2.0*dboxUM1[1] + 0.05;
dboxSMPb1[2] = fgkThLead/2.;
TVirtualMC::GetMC()->Gsvolu("ESPA","BOX", idtmed[698], dboxSMPb1, 3);
Float_t xpa1,xpa2,xpa3,ypa1,ypa2;
xpa1 = -dboxSMPb1[0] + dboxUM1[0];
xpa2 = xpa1 + dboxUM1[0] + 0.15 + dboxUM1[0];
xpa3 = xpa2 + dboxUM1[0] + 0.15 + dboxUM1[0];
ypa1 = dboxSMPb1[1] - dboxUM1[1];
ypa2 = ypa1 - dboxUM1[1] - 0.1 - dboxUM1[1];
TVirtualMC::GetMC()->Gspos("EPB1", 1, "ESPA", xpa1, ypa1, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EPB1", 2, "ESPA", xpa2, ypa1, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EPB1", 3, "ESPA", xpa3, ypa1, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EPB1", 4, "ESPA", xpa1, ypa2, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EPB1", 5, "ESPA", xpa2, ypa2, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EPB1", 6, "ESPA", xpa3, ypa2, 0., 0, "ONLY");
Float_t dboxSMPb2[3];
dboxSMPb2[0] = 2.0*dboxUM2[0] + 0.075;
dboxSMPb2[1] = 3.0*dboxUM2[1] + (2.0*0.05);
dboxSMPb2[2] = fgkThLead/2.;
TVirtualMC::GetMC()->Gsvolu("ESPB","BOX", idtmed[698], dboxSMPb2, 3);
Float_t xpb1,xpb2,ypb1,ypb2,ypb3;
xpb1 = -dboxSMPb2[0] + dboxUM2[0];
xpb2 = xpb1 + dboxUM2[0] + 0.15 + dboxUM2[0];
ypb1 = dboxSMPb2[1] - dboxUM2[1];
ypb2 = ypb1 - dboxUM2[1] - 0.1 - dboxUM2[1];
ypb3 = ypb2 - dboxUM2[1] - 0.1 - dboxUM2[1];
TVirtualMC::GetMC()->Gspos("EPB2", 1, "ESPB", xpb1, ypb1, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EPB2", 2, "ESPB", xpb2, ypb1, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EPB2", 3, "ESPB", xpb1, ypb2, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EPB2", 4, "ESPB", xpb2, ypb2, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EPB2", 5, "ESPB", xpb1, ypb3, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EPB2", 6, "ESPB", xpb2, ypb3, 0., 0, "ONLY");
Float_t dboxFEE[3];
dboxFEE[0] = 0.05;
dboxFEE[1] = 3.50;
dboxFEE[2] = 1.20;
TVirtualMC::GetMC()->Gsvolu("EFEE","BOX", idtmed[607], dboxFEE, 3);
Float_t dboxFEEBPlaneA[3];
dboxFEEBPlaneA[0] = dboxBPlaneA[0];
dboxFEEBPlaneA[1] = dboxBPlaneA[1];
dboxFEEBPlaneA[2] = 1.2;
TVirtualMC::GetMC()->Gsvolu("EFBA","BOX", idtmed[698], dboxFEEBPlaneA, 3);
Float_t dboxFEEBPlaneB[3];
dboxFEEBPlaneB[0] = dboxBPlaneB[0];
dboxFEEBPlaneB[1] = dboxBPlaneB[1];
dboxFEEBPlaneB[2] = 1.2;
TVirtualMC::GetMC()->Gsvolu("EFBB","BOX", idtmed[698], dboxFEEBPlaneB, 3);
Float_t xFee;
Float_t yFee;
Float_t zFee = 0.0;
Float_t xA = 0.25;
Float_t yA = 4.00;
Float_t xSepa = 1.70;
Float_t ySepa = 8.00;
number = 1;
yFee = dboxFEEBPlaneA[1] - yA;
for (i = 1; i <= 6; ++i)
{
xFee = -dboxFEEBPlaneA[0] + xA;
for (j = 1; j <= 12; ++j)
{
TVirtualMC::GetMC()->Gspos("EFEE", number, "EFBA", xFee,yFee,zFee, 0, "ONLY");
xFee += xSepa;
number += 1;
}
yFee -= ySepa;
}
number = 1;
yFee = dboxFEEBPlaneB[1] - yA;
for (i = 1; i <= 3; ++i)
{
xFee = -dboxFEEBPlaneB[0] + xA;
for (j = 1; j <= 24; ++j)
{
TVirtualMC::GetMC()->Gspos("EFEE", number, "EFBB", xFee,yFee,zFee, 0, "ONLY");
xFee += xSepa;
number += 1;
}
yFee -= ySepa;
}
Float_t dboxEFSA[3];
dboxEFSA[0] = 3.0*dboxFEEBPlaneA[0] + 0.92;
dboxEFSA[1] = 2.0*dboxFEEBPlaneA[1] + (0.95/2.0);
dboxEFSA[2] = dboxFEEBPlaneA[2];
TVirtualMC::GetMC()->Gsvolu("EFSA","BOX", idtmed[698],dboxEFSA, 3);
Float_t dboxEFSB[3];
dboxEFSB[0] = 2.0*dboxFEEBPlaneB[0] + (0.938/2.0);
dboxEFSB[1] = 3.0*dboxFEEBPlaneB[1] + 1.05;
dboxEFSB[2] = dboxFEEBPlaneB[2];
TVirtualMC::GetMC()->Gsvolu("EFSB","BOX", idtmed[698],dboxEFSB, 3);
Float_t xfs1,xfs2,xfs3,yfs1,yfs2,yfs3;
xfs1 = -dboxEFSA[0] + dboxFEEBPlaneA[0];
xfs2 = xfs1 + dboxFEEBPlaneA[0] + 0.92 + dboxFEEBPlaneA[0];
xfs3 = xfs2 + dboxFEEBPlaneA[0] + 0.92 + dboxFEEBPlaneA[0];
yfs1 = dboxEFSA[1] - dboxFEEBPlaneA[1];
yfs2 = yfs1 - dboxFEEBPlaneA[1] - 0.95 - dboxFEEBPlaneA[1];
TVirtualMC::GetMC()->Gspos("EFBA", 2, "EFSA", xfs2, yfs1, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EFBA", 3, "EFSA", xfs3, yfs1, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EFBA", 4, "EFSA", xfs1, yfs2, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EFBA", 5, "EFSA", xfs2, yfs2, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EFBA", 6, "EFSA", xfs3, yfs2, 0., 0, "ONLY");
xfs1 = -dboxEFSB[0] + dboxFEEBPlaneB[0];
xfs2 = xfs1 + dboxFEEBPlaneB[0] + 0.938 + dboxFEEBPlaneB[0];
yfs1 = dboxEFSB[1] - dboxFEEBPlaneB[1];
yfs2 = yfs1 - dboxFEEBPlaneB[1] - 1.05 - dboxFEEBPlaneB[1];
yfs3 = yfs2 - dboxFEEBPlaneB[1] - 1.05 - dboxFEEBPlaneB[1];
TVirtualMC::GetMC()->Gspos("EFBB", 3, "EFSB", xfs1, yfs2, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EFBB", 4, "EFSB", xfs2, yfs2, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EFBB", 5, "EFSB", xfs1, yfs3, 0., 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EFBB", 6, "EFSB", xfs2, yfs3, 0., 0, "ONLY");
}
void AliPMDv2008::CreatePMD()
{
Float_t zp;
Int_t jhrot12,jhrot13, irotdm;
Int_t *idtmed = fIdtmed->GetArray()-599;
Float_t dboxFea[3];
dboxFea[0] = fSMLengthax;
dboxFea[1] = fSMLengthay;
dboxFea[2] = fgkThSteel/2.;
TVirtualMC::GetMC()->Gsvolu("EFEA","BOX", idtmed[618], dboxFea, 3);
Float_t dboxFeb[3];
dboxFeb[0] = fSMLengthbx;
dboxFeb[1] = fSMLengthby;
dboxFeb[2] = fgkThSteel/2.;
TVirtualMC::GetMC()->Gsvolu("EFEB","BOX", idtmed[618], dboxFeb, 3);
AliMatrix(irotdm, 90., 0., 90., 90., 180., 0.);
AliMatrix(jhrot12, 90., 180., 90., 270., 0., 0.);
AliMatrix(jhrot13, 90., 240., 90., 330., 0., 0.);
Float_t gaspmd[3];
gaspmd[0] = fSMLengthax;
gaspmd[1] = fSMLengthay;
gaspmd[2] = fSMthick;
TVirtualMC::GetMC()->Gsvolu("EPM1", "BOX", idtmed[698], gaspmd, 3);
TVirtualMC::GetMC()->Gsvolu("EPM2", "BOX", idtmed[698], gaspmd, 3);
Float_t zpsa,zpba,zfea,zcva,zfee;
zfee=-gaspmd[2] + 1.2;
zcva = zfee + 1.2 + fDthick;
zfea = zcva + fDthick + fgkThSteel/2.;
TVirtualMC::GetMC()->Gspos("EFEA", 1, "EPM1", 0., 0., zfea, 0, "ONLY");
zpba=zfea+fgkThSteel/2.+ fgkThLead/2.;
TVirtualMC::GetMC()->Gspos("ESPA", 1, "EPM1", 0., 0., zpba, 0, "ONLY");
zpsa = zpba + fgkThLead/2. + fDthick;
TVirtualMC::GetMC()->Gspos("ESMA", 1, "EPM1", 0., 0., zpsa, 0, "ONLY");
zfee=zpsa + fDthick + 1.2;
TVirtualMC::GetMC()->Gspos("EFSA", 3, "EPM1", 0., 0., zfee, 0, "ONLY");
gaspmd[0] = fSMLengthbx;
gaspmd[1] = fSMLengthby;
gaspmd[2] = fSMthick;
TVirtualMC::GetMC()->Gsvolu("EPM3", "BOX", idtmed[698], gaspmd, 3);
TVirtualMC::GetMC()->Gsvolu("EPM4", "BOX", idtmed[698], gaspmd, 3);
Float_t zpsb,zpbb,zfeb,zcvb;
zfee=-gaspmd[2] + 1.2;
zcvb= zfee + 1.2 + fDthick;
zfeb= zcvb + fDthick + fgkThSteel/2.;
TVirtualMC::GetMC()->Gspos("EFEB", 4, "EPM4", 0., 0., zfeb, 0, "ONLY");
zpbb= zfeb + fgkThSteel/2.+ fgkThLead/2.;
TVirtualMC::GetMC()->Gspos("ESPB", 4, "EPM4", 0., 0., zpbb, 0, "ONLY");
zpsb = zpbb + fgkThLead/2.+ fDthick;
TVirtualMC::GetMC()->Gspos("ESMB", 4, "EPM4", 0., 0., zpsb, jhrot12, "ONLY");
zfee=zpsb + fDthick + 1.2;
TVirtualMC::GetMC()->Gspos("EFSB", 8, "EPM4", 0., 0., zfee, jhrot12, "ONLY");
zp = fgkZdist;
Float_t xfinal,yfinal;
Float_t xsmb,ysmb;
Float_t xsma,ysma;
xfinal = fSMLengthax + 0.48/2 + fSMLengthbx;
yfinal = fSMLengthay + 0.20/2 + fSMLengthby;
xsma = xfinal - fSMLengthax;
ysma = yfinal - fSMLengthay;
xsmb = -xfinal + fSMLengthbx;
ysmb = yfinal - fSMLengthby;
TVirtualMC::GetMC()->Gspos("EPM1", 1, "ALIC", xsma,ysma,zp, 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EPM2", 1, "ALIC", -xsma,-ysma,zp, 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EPM3", 1, "ALIC", xsmb,ysmb,zp, 0, "ONLY");
TVirtualMC::GetMC()->Gspos("EPM4", 1, "ALIC", -xsmb,-ysmb,zp, 0, "ONLY");
}
void AliPMDv2008::CreateMaterials()
{
Int_t isxfld = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->Integ();
Float_t sxmgmx = ((AliMagF*)TGeoGlobalMagField::Instance()->GetField())->Max();
AliMaterial(1, "Pb $", 207.19, 82., 11.35, .56, 18.5);
Float_t dAr = 0.001782;
Float_t x0Ar = 19.55 / dAr;
AliMaterial(2, "Argon$", 39.95, 18., dAr, x0Ar, 6.5e4);
Float_t aCO2[2] = { 12.,16. };
Float_t zCO2[2] = { 6.,8. };
Float_t wCO2[2] = { 1.,2. };
Float_t dCO2 = 0.001977;
AliMixture(3, "CO2 $", aCO2, zCO2, dCO2, -2, wCO2);
AliMaterial(4, "Al $", 26.98, 13., 2.7, 8.9, 18.5);
Float_t aArCO2[3] = {39.948,12.0107,15.9994};
Float_t zArCO2[3] = {18.,6.,8.};
Float_t wArCO2[3] = {0.7,0.08,0.22};
Float_t dArCO2 = dAr * 0.7 + dCO2 * 0.3;
AliMixture(5, "ArCO2$", aArCO2, zArCO2, dArCO2, 3, wArCO2);
AliMaterial(6, "Fe $", 55.85, 26., 7.87, 1.76, 18.5);
Float_t aG10[4]={1.,12.011,15.9994,28.086};
Float_t zG10[4]={1.,6.,8.,14.};
Float_t wG10[4]={0.15201,0.10641,0.49444,0.24714};
AliMixture(8,"G10",aG10,zG10,1.7,4,wG10);
AliMaterial(15, "Cu $", 63.54, 29., 8.96, 1.43, 15.);
Float_t aSteel[4] = { 55.847,51.9961,58.6934,28.0855 };
Float_t zSteel[4] = { 26.,24.,28.,14. };
Float_t wSteel[4] = { .715,.18,.1,.005 };
Float_t dSteel = 7.88;
AliMixture(19, "STAINLESS STEEL$", aSteel, zSteel, dSteel, 4, wSteel);
Float_t aAir[4]={12.0107,14.0067,15.9994,39.948};
Float_t zAir[4]={6.,7.,8.,18.};
Float_t wAir[4]={0.000124,0.755267,0.231781,0.012827};
Float_t dAir1 = 1.20479E-10;
Float_t dAir = 1.20479E-3;
AliMixture(98, "Vacum$", aAir, zAir, dAir1, 4, wAir);
AliMixture(99, "Air $", aAir, zAir, dAir , 4, wAir);
AliMedium(1, "Pb conv.$", 1, 0, 0, isxfld, sxmgmx, 1., .1, .01, .1);
AliMedium(4, "Al $", 4, 0, 0, isxfld, sxmgmx, .1, .1, .01, .1);
AliMedium(5, "ArCO2 $", 5, 1, 0, isxfld, sxmgmx, .1, .1, .10, .1);
AliMedium(6, "Fe $", 6, 0, 0, isxfld, sxmgmx, .1, .1, .01, .1);
AliMedium(8, "G10plate$", 8, 0, 0, isxfld, sxmgmx, 1., .1, .01, .1);
AliMedium(15, "Cu $", 15, 0, 0, isxfld, sxmgmx, .1, .1, .01, .1);
AliMedium(19, "S steel$", 19, 0, 0, isxfld, sxmgmx, 1., .1, .01, .1);
AliMedium(98, "Vacuum $", 98, 0, 0, isxfld, sxmgmx, 1., .1, .10, 10);
AliMedium(99, "Air gaps$", 99, 0, 0, isxfld, sxmgmx, 1., .1, .10, .1);
AliDebug(1,"Outside create materials");
}
void AliPMDv2008::Init()
{
AliDebug(2,"Inside Init");
AliDebug(2,"PMD simulation package (v1) initialised");
AliDebug(2,"parameters of pmd");
AliDebug(2,Form("%10.2f %10.2f %10.2f %10.2f\n",
fgkCellRadius,fgkCellWall,fgkCellDepth,fgkZdist));
Int_t *idtmed = fIdtmed->GetArray()-599;
fMedSens=idtmed[605-1];
gGeoManager->SetVolumeAttribute("ECAR", "SEEN", 0);
gGeoManager->SetVolumeAttribute("ECCU", "SEEN", 0);
gGeoManager->SetVolumeAttribute("ECCU", "COLO", 4);
gGeoManager->SetVolumeAttribute("EST1", "SEEN", 0);
gGeoManager->SetVolumeAttribute("EST2", "SEEN", 0);
gGeoManager->SetVolumeAttribute("EHC1", "SEEN", 0);
gGeoManager->SetVolumeAttribute("EHC2", "SEEN", 0);
gGeoManager->SetVolumeAttribute("EPCA", "SEEN", 0);
gGeoManager->SetVolumeAttribute("EBKA", "SEEN", 0);
gGeoManager->SetVolumeAttribute("ECGA", "SEEN", 0);
gGeoManager->SetVolumeAttribute("ECVA", "SEEN", 0);
gGeoManager->SetVolumeAttribute("EDGA", "SEEN", 0);
gGeoManager->SetVolumeAttribute("EDVA", "SEEN", 0);
gGeoManager->SetVolumeAttribute("ESSA", "SEEN", 0);
gGeoManager->SetVolumeAttribute("ESVA", "SEEN", 0);
gGeoManager->SetVolumeAttribute("EUM1", "SEEN", 0);
gGeoManager->SetVolumeAttribute("EUV1", "SEEN", 0);
gGeoManager->SetVolumeAttribute("EBPA", "SEEN", 0);
gGeoManager->SetVolumeAttribute("EPCB", "SEEN", 0);
gGeoManager->SetVolumeAttribute("EBKB", "SEEN", 0);
gGeoManager->SetVolumeAttribute("ECGB", "SEEN", 0);
gGeoManager->SetVolumeAttribute("ECVB", "SEEN", 0);
gGeoManager->SetVolumeAttribute("EDGB", "SEEN", 0);
gGeoManager->SetVolumeAttribute("EDVB", "SEEN", 0);
gGeoManager->SetVolumeAttribute("ESSB", "SEEN", 0);
gGeoManager->SetVolumeAttribute("ESVB", "SEEN", 0);
gGeoManager->SetVolumeAttribute("EUM2", "SEEN", 0);
gGeoManager->SetVolumeAttribute("EUV2", "SEEN", 0);
gGeoManager->SetVolumeAttribute("EBPB", "SEEN", 0);
gGeoManager->SetVolumeAttribute("EPB1", "SEEN", 0);
gGeoManager->SetVolumeAttribute("EPB2", "SEEN", 0);
gGeoManager->SetVolumeAttribute("ESMA", "SEEN", 0);
gGeoManager->SetVolumeAttribute("EMVA", "SEEN", 0);
gGeoManager->SetVolumeAttribute("ESMB", "SEEN", 0);
gGeoManager->SetVolumeAttribute("EMVB", "SEEN", 0);
gGeoManager->SetVolumeAttribute("ESPA", "SEEN", 0);
gGeoManager->SetVolumeAttribute("ESPB", "SEEN", 0);
gGeoManager->SetVolumeAttribute("EFEE", "SEEN", 0);
gGeoManager->SetVolumeAttribute("EFEE", "COLO", 4);
gGeoManager->SetVolumeAttribute("EFBA", "SEEN", 0);
gGeoManager->SetVolumeAttribute("EFBB", "SEEN", 0);
gGeoManager->SetVolumeAttribute("EFSA", "SEEN", 0);
gGeoManager->SetVolumeAttribute("EFSB", "SEEN", 0);
gGeoManager->SetVolumeAttribute("EFEA", "SEEN", 0);
gGeoManager->SetVolumeAttribute("EFEB", "SEEN", 0);
gGeoManager->SetVolumeAttribute("EPM1", "SEEN", 1);
gGeoManager->SetVolumeAttribute("EPM2", "SEEN", 1);
gGeoManager->SetVolumeAttribute("EPM3", "SEEN", 1);
gGeoManager->SetVolumeAttribute("EPM4", "SEEN", 1);
}
void AliPMDv2008::StepManager()
{
Int_t copy;
Float_t hits[5], destep;
Float_t center[3] = {0,0,0};
Int_t vol[6];
if(TVirtualMC::GetMC()->CurrentMedium() == fMedSens && (destep = TVirtualMC::GetMC()->Edep())) {
TVirtualMC::GetMC()->CurrentVolID(copy);
vol[0] = copy;
TVirtualMC::GetMC()->CurrentVolOffID(1,copy);
vol[1] = copy;
TVirtualMC::GetMC()->CurrentVolOffID(2,copy);
vol[2] = copy;
TVirtualMC::GetMC()->CurrentVolOffID(3,copy);
vol[3] = copy;
TVirtualMC::GetMC()->CurrentVolOffID(4,copy);
vol[4] = copy;
TVirtualMC::GetMC()->CurrentVolOffID(5,copy);
vol[5] = copy;
TVirtualMC::GetMC()->Gdtom(center,hits,1);
hits[3] = destep*1e9;
hits[4] = TVirtualMC::GetMC()->TrackTime();
AddHit(gAlice->GetMCApp()->GetCurrentTrackNumber(), vol, hits);
}
}
void AliPMDv2008::GetParameters()
{
fSMLengthax = 32.7434;
fSMLengthbx = 42.5886;
fSMLengthay = 49.1;
fSMLengthby = 37.675;
fDthick = fgkThSS/2. +0.15;
fSMthick = 2.0*(fgkThSS/2. +0.15)
+fgkThSteel/2.+fgkThLead/2.0 + 2.4;
}
void AliPMDv2008::AddAlignableVolumes() const
{
SetSectorAlignable();
}
void AliPMDv2008::SetSectorAlignable() const
{
TString vpsector = "ALIC_1/EPM";
TString vpappend = "_1";
TString snsector="PMD/Sector";
TString volpath, symname;
for(Int_t cnt=1; cnt<=4; cnt++){
volpath = vpsector;
volpath += cnt;
volpath += vpappend;
symname = snsector;
symname += cnt;
if(!gGeoManager->SetAlignableEntry(symname.Data(),volpath.Data()))
{
AliFatal("Unable to set alignable entry!");
}
}
}