C File name: svx.f ( was previously inr.f) C ---------------- C C Original author: Shaheen Tonse (LLNL) C Creation date: March 18, 1993 C C Purpose: Set up the Silicon Vertex Detector (SVX) C C Revision History: c """""""""""""""" C C.F. Maguire June 28, 2001 Convert to use for PHENIX upgrades group c c V. L. Rykov 26-Aug-2003: c 1) The vertex detector encloser ("envelope", cooling bag) c is introduced. c 2) For the barrel VTX, "plain" Si cylinders are replaced with c something more realistic (ladders and sensors). c 3) Set of hit components extended with the "wish list", but, c for the time being, only 8 the first components are stored. c 4) The gains for position resolution set for resolution 0.1 mkm. c V. L. Rykov 03-Sep-2003: c Using the full (extended) set for hit components. c V. L. Rykov 19-Feb-2004: c inrORIG & inrFACT adjusted for geting rid of the bad diagnostic c due to endcap, which still is not splitted into sensors. c V. L. Rykov 29-Mar-2004: c Fixed the bug in the passive layer definition. c V. L. Rykov 16-Apr-2004: c SVX cage and barrel parameters along with the sensor rotation c matrices and translation vestors are written into the c svxPISA.par file. c V. L. Rykov 21-Apr-2004: c Fixed some strange volume interference of SVX & MUI first c reported by Hua Pei on 04/19/2004. c c Hubert van Hecke, Aug/Sep 2004: replaced endcap geometry with detailed lampshade c description from Gerd Kunde, Michael Malik and Jan Boissevain c Moved the material/medium numbers into the MVD range c Hubert van Hecke, Jan 2005: shrunk middle 4 endcaps, added pipes and cables c c The tree structure for the barrel is: c c HALL c | c SIEN (VTX Envelope/cage) c | c SICG (VTX cage inner surface) c | c -------------------------------------------------------------- c | | | | | ....| c SI01(1) SI01(2)... SI02(1) SI02(2) ... (B-ladders) (Endcaps) c | | | c ------ ------ ------------------ c | |... | |... | | c SISN(1) SISN(2) ... (B-sensors) * *==================== The tree structure for the endcaps is: ============================ * * HALL * | * SIEN (VTX Envelope/Cage) * | * SICG (VTX cage inner surface) * | * ++-----++--------------+----+--+--------------------+--+-------------++++ * || || | | | | |||| * SI05,06, 11,12 SI07,SI10 SI08,SI09 (barrel) * | | | * SIPB(24) SIPM(24) SIPS(24) * | | | * +-----+--------+--------+-------+ | | * | | | | | | | * SICB SISI(4)* SCHP(22) SBX1(2) SBX2(2) | | * | | * big shades^ +-----+--------+-------+-------+ | * | | | | | | * medium shades>> SICC SISI(4)* SCHP(16)SBX1(2) SBX2(2) | * | * +-----+--------+-------+-------+ * | | | | | * small shades>> SICC SISI(4)* SCHP(16)SBX1(2) SBX2(2) * * SI05 - SI12 are the 8 lampshades, 4 big, 2 medium and 2 small * each lampshade containing 24 panels, SIPB, SIPM or SIPS * Each panel consists of carbon support, silicon(4), readout chips(22, 16 or 12) and * standoff strips(2x2) * *===================================================================================== SUBROUTINE SVX(FULL,NH) C Implicit none C C--- Formal Argument Declarations C ---------------------------- c... Input (?): character*4 full ! set before call in gugeom C... Output: number of components of a hit integer*4 nh ! set before call in gugeom C C--- External Functions C ------------------ C C ================================================================ C--- Global Declarations C ================================================================ include 'guphnx.inc' include 'gclist.inc' include 'gconst.inc' include 'gcflag.inc' include 'gugeom.inc' include 'gcvolu.inc' C C need to access zebra to write parameters to FZOUT file C include 'fstore.inc' include 'sublink.inc' include 'fpdlink.inc' C C ================================================================ C--- Local declarations of the input data from phnx.par file C ================================================================ c--- VTX Envelope/Cage parameters: Volumes SIEN(outer)/SICG(inner) c """"""""""""""""""""""""""""""""""""""""""""""""""""""""""""" Real sili_cg_rmn /2.2/ ! Inner cage radius, cm Real sili_cg_thck /0.5/ ! Cage wall thickness, cm Real sili_cg_inthck /0.2/ ! Thickness of the beam pipe ins., cm Real sili_cg_tempc /0.0/ ! Temperature inside the Cage, deg. C Integer sili_cg_npcon /6/ ! Number of corners for SIEN's PCON Real sili_cg_z(10) /-40.0,-25.6,-13.0 ! z-pos. of the Cage corners: * ,13.0,25.6,40.0,4*0./ ! MUST: z(1)=0) or c ! sensor z-overlap (if <0), cm Real sili_br_tilt(20) /10.0,6.0 ! Ladder tilt in layers, deg. * ,5.5,5.0 * ,16*0./ Integer sili_br_nsn(20) /4,19*5/ ! Number of Si sensors/ladder Real sili_br_r(20) /2.5,6.,8.,10.,16*20./ ! Radial positions ... c ! (read carefully!!!) c ! of the center line, c ! x=y=0 (loca Si sensor c ! coordinates), of the Si c ! sensors in the layers Real sili_br_z(20) /20*0./ ! Z-pos. of the ladder centers in SICG Real sili_br_dphi(20) /29.40,29.40 ! Azim. spacing of ladders, deg. * ,22.87,18.71 * ,16*0./ Integer sili_br_nsec(20) /20*2/ ! Number of cont. phi-sections/layer Real sili_br_phic(5,20) /0.,180.,3*0. ! Sect. center azimuth pos. * ,0.,180.,3*0. * ,0.,180.,3*0. * ,0.,180.,3*0. * ,80*0./ Integer sili_br_nlad(5,20) /5*5 ! Number of ladders/phi-section * ,5*5 * ,5*7 * ,5*9 * ,80*1/ c ====================================================================== c The parameters between === are for the temporary fake support emulation c with "plain" discs introduced for accomodating the request of c some people. You are welcome to spread this brand new "artificial c butter" over your piece of bread. Just, please, do not eat it! - VR Integer sili_br_nspring /2/ ! Number of rings Real sili_br_spz(20) /-17.,17.,18*0./ ! Z-positions, cm Real sili_br_sprin(20) /20*2.4/ ! Inner radii, cm Real sili_br_sprout(20) /20*10./ ! Outer radii, cm Real sili_br_spthck(20) /20*0.25/ ! Half-hickness, cm c ====================================================================== namelist /sili_br_par/ sili_br_nlayers,sili_br_r,sili_br_z $ ,sili_br_nsn,sili_br_snhalfx,sili_br_snhalfy $ ,sili_br_snhalfz,sili_br_x0add,sili_br_snzgap $ ,sili_br_nsec,sili_br_nlad,sili_br_phic,sili_br_dphi $ ,sili_br_tilt $ ,sili_br_nspring,sili_br_spz,sili_br_sprin $ ,sili_br_sprout,sili_br_spthck c--- VTX Endcap parameters ---------------------------------------------- integer sili_endcap_type /1/ ! cones (1) or flat disks (2) integer sili_endcap_layers /8/ ! 4 south, 4 north real sili_endcap_z(8) / -35.16, -29.16, -23.16, -17.87, & 17.87, 23.16, 29.16, 35.16/ namelist /sili_endcap_par/ & sili_endcap_type, sili_endcap_layers, sili_endcap_z C ================================================================ C--- Local definitions C ================================================================ c Input filename character*50 par_file /'phnx.par'/ ! Input filename character*50 svxpar_file /'svxPISA.par'/ ! Output file for parameters character*4 set_id /'SVX '/ ! Detector/hit set ID c """"""""""""""""""""""""""""""""""""""""""""""""""""""""""""" c VTX cage volume names character*4 siliEnvelope/'SIEN'/ ! Envelope (outer cage surface) character*4 siliCage /'SICG'/ ! Cage inner surface character*4 siliInCage /'SIIC'/ ! Thermal insulation around beam-pipe character*4 siliSupport /'SISP'/ ! Fake barrel support c VTX layer names (for barrel, names of ladders in a layer) character*4 siliNames(12) /'SI01', 'SI02', 'SI03', 'SI04', & 'SI05', 'SI06', 'SI07', 'SI08', & 'SI09', 'SI10', 'SI11', 'SI12'/ character*4 siliBrSensor /'SISN'/ ! Barrel sensor name character*4 siliBrPassive /'SIPV'/ ! Barrel passive material name c """"""""""""""""""""""""""""""""""""""""""""""""""""""""""""" c The following are used in GSDET: namesv(4) is to be replaced with c 'SI01', 'SI02', ... at running time c Integer nwpa /500/ ! Init. size of HITS banks Integer nwsa /500/ ! Init. size of DIGI banks Integer idtype /2001/ ! User def. detector type Integer nbrv /5/ ! Num. of br. vol. desc. Integer necv /6/ ! Num. of endcap vol. desc. Integer nv /6/ ! max(nbrv,necv) Integer nbitsv(6) /6*8/ ! Bits to pack vol. copy # Character*4 namesv(6) /'HALL','SIEN','SICG' * ,'SIxx','SISN','SDUM'/ ! Volume names, barrel Character*4 namesw(6) /'HALL','SIEN','SICG' * ,'SISx','SIPx','SISI'/ ! Volume names, endcaps c """"""""""""""""""""""""""""""""""""""""""""""""""""""""""""" c The following are used in GSDETH c """""""""""""""""""""""""""""""" c Hit parameters will be global position(3), energy loss, time of flight, c particle type and c entry momentum(3), local in & out positions(6), c Integer nhh /15/ ! Number of hit components integer*4 inrNBITSH(15) /15*32/ ! Bits for packing the hits c Hit component names character*4 inrNMSH(15) /'POSX','POSY','POSZ' ! Global positions & ,'DELE','TOFL' ! Energy loss & TOF & ,'P_ID','MOMX', 'MOMY', 'MOMZ' ! Particle ID & Entry mom. & ,'XILC','YILC','ZILC','XOLC','YOLC','ZOLC'/ ! Local entry & exit c Default setting of offsets and gains REAL inrORIG(15) /3*1000.,3*0.,3*1000.,6*1000./ ! offsets REAL inrFACT(15) /3*100000.,1.E7,1.e12,1.0,3*100000. & ,6*100000./ ! gains c c The above gains give c - 0.1 keV energy deposition resolution c - 0.0001 mm position resolution c - 0.01 MeV/c momentum resolution c """"""""""""""""""""""""""""""""""""""""""""""""""""""""""""" c Rotation matrices c """"""""""""""""""""""""""""""""""""""""""""""""""""""""""""" Integer irot_cage ! Cage (SIEN) rotation matrix # C ================================================================ C--- Local definitions of materials/media C ================================================================ Integer sili_med_silicon /10/ ! Sensitive silicon nmed (nmat=50,Si) Integer sili_med_passive /26/ ! Ladder passive nmed (nmat=09,Al) Integer sili_med_cg /120/ ! Cage (SIEN-SICG) nmat/nmed (Rohacell) Integer sili_med_coldair /121/ ! Gas inside the SICG (cold air) Integer sili_med_gfrp /122/ ! GFRP for the fake support integer sili_med_carbon /123/ ! carbon-carbon composite c Material for the cage, volumes SIEN-SICG c """"""""""""""""""""""""""""""""""""""""""""""""""""""""""""" C Copied from pisa200/src/ver/vermatdef.f, Rev 1.2 C Rohacell(H11-C8-N1-O12) Mixture Parameters for the vertex detector: REAL AROHA(4)/ 1.008 , 12.01 , 14.008 , 16. / ! A for Rohacell REAL ZROHA(4)/ 1. , 6. , 7. , 8. / ! Z for Rohacell REAL WROHA(4)/ 11. , 8. , 1. , 2. / ! Weights for Rohacell c """"""""""""""""""""""""""""""""""""""""""""""""""""""""""""" C ================================================================ c--- Local work variables (counters, etc) C ================================================================ CHARACTER*20 NAMTMED INTEGER NMAT, ISVOL, IFIELD, NWBUF, LNAM(6) integer LNUM(6) /1,1,1,1,1,1/ REAL FIELDM,TMAXFD,DMAXMS,DEEMAX,EPSIL,STMIN & ,WMAT(3),UBUF(10),ANUMB,ZNUMB,DENS,RADL,ABSL CHARACTER*80 CHFORM character*4 v_m_name,v_i_name, si_name Integer iLayer,iPoint,icnt,nr,npar,nmed,ivolu & ,iset,idet,iod Integer nladd, iladd, isec, isen, nsensors Real dim_sili(20),philadd,dphid,dphir,dphit,sgn & ,phi1,phi2,dthck,phirotladd,rladd,xladd,yladd & ,ladd_halfx,ladd_halfthck,ladd_halfz,pasv_radl,pasv_halfy & ,sen_y(2) real panels, plength, plength2, plen3, plen4, panthk, strip, & strip2, strip3, strip4, strip5, strip6, strip7, gap1, chiplen, & chipwid, chipthk, sithk, silen, silen2, silen3, silen4, silen5, & silen6, silen7, silen8, silen9, rinner, zangle, zsep, zdist, & zdist2, sangle2, panthick, d, d2, dd, halfz, halfz2, deg, rad, & pangle, sangle, stripz, stripx, plen5, plen6, stripz2, stripx2, & l1, stripz3, stripx3, stripy3, stripz4, stripx4, stripy4, & par(20), parb(20), parm(20), pars(20), par_s1(4), par_s2(4), & par_s4(4), delz, dely, shade_z, shade_zmid, z_shade, plen7, & plen8, plength3, d3, halfz3, stripx5, stripz5 integer ivol1, ii, i, j, k, l, ir1, ir2, irx(24), irsh, ir, & ipanel, ierr, ishade, irshade character*4 sil_name, panel_name C ================================================================ C ================================================================ C C Executable code C =============== C c>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> c BEGIN c>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> c c Read the geometery file segments c open(unit=15,file=par_file,status='OLD',err=999) read(15,nml=sili_cg_par,err=998) read(15,nml=sili_br_par,err=997) read(15,nml=sili_endcap_par,err=996) close(unit=15) C C only book volumes if input parameters are OK C IF(CVOLU_OPT(1,3).EQ.'FULL'.OR.CVOLU_OPT(1,3).EQ.'VOLS')THEN NH = nhh ! Number of hit components to NH output parameter nv = max(nbrv, necv) ! Number of volume descriptors C... Create outer surface of the cage, SIEN dim_sili(1) = 0. dim_sili(2) = 360. dim_sili(3) = sili_cg_npcon npar = 3 Do icnt = 1, sili_cg_npcon npar = npar + 1 dim_sili(npar) = sili_cg_z(icnt) npar = npar + 1 dim_sili(npar) = sili_cg_rmn npar = npar + 1 dim_sili(npar) = sili_cg_rmx(icnt) Enddo nwbuf = 1 C************************************************************************** C Define mixture ROHACELL for the SIEN C Copied from pisa200/src/ver/vermatdef.f, Rev 1.2 C NMAT = sili_med_cg ! Rohacell CALL GSMIXT (NMAT,' ROHACELL$',AROHA,ZROHA,0.075,-4,WROHA) C C Tracking media # sili_cg_wall - Rohacell C NMED = sili_med_cg ! Rohacell ISVOL = 0 ! Not sensitive IFIELD = 1 ! Magnetic field FIELDM = 10.0 ! max field TMAXFD = 45.0 ! maximum angle due to field (one step) in deg DMAXMS = 0.2 ! max disp. due to mulsct. in one step (cm) DEEMAX = 0.1 ! max fractional energy loss in one step EPSIL = .001 ! tracking precision (cm) STMIN = 0.5 ! min step due to e loss or mulsct. (cm) UBUF(1) = 0. ! tracking stop switch CALL GSTMED(NMED,'Rohacell$',NMAT,ISVOL,IFIELD 1 ,FIELDM,TMAXFD,DMAXMS,DEEMAX,EPSIL,STMIN,UBUF,NWBUF) C************************************************************************** v_m_name = 'HALL' v_i_name = siliEnvelope nr = 1 call gsvolu(v_i_name,'PCON',nmed,dim_sili,npar,ivolu) call gsatt(v_i_name,'SEEN',1) call gspos(v_i_name,nr,v_m_name * ,sili_cg_xdisp,sili_cg_ydisp,sili_cg_zdisp * ,irotnull,'ONLY') C... Create inner volume of the cage SICG npar = 3 sgn = 1. If(sili_cg_z(1).GT.sili_cg_z(sili_cg_npcon)) sgn = -1. phi1 = 1.570796327 Do icnt = 1, sili_cg_npcon npar = npar + 3 If(icnt .LT. sili_cg_npcon) THEN phi2 = ATAN2(sili_cg_rmx(icnt+1)-sili_cg_rmx(icnt), * sili_cg_z(icnt+1)-sili_cg_z(icnt)) Else phi2 = -1.570796327 Endif dphir = 0.5*(phi1 + phi2) dthck = sili_cg_thck/COS(0.5*(phi1-phi2)) dim_sili(npar-2) = sili_cg_z(icnt) + dthck*SIN(dphir) dim_sili(npar-1) = sili_cg_rmn ! Inner radius dim_sili(npar) = sili_cg_rmx(icnt) - dthck*COS(dphir) If(dim_sili(npar) .LT. sili_cg_rmn) THEN If(icnt .EQ. 1) THEN dphir = phi2 Elseif(icnt .EQ. sili_cg_npcon) THEN dphir = phi1 Endif dim_sili(npar) = sili_cg_rmn dim_sili(npar-2) = (sili_cg_rmx(icnt)-sili_cg_rmn)/ * TAN(dphir) dim_sili(npar-2) = sili_cg_z(icnt) + * sgn*(sili_cg_thck/SIN(dphir) - dim_sili(npar-2)) Endif phi1 = phi2 Enddo ! loop over sili_cg_npcon C************************************************************************** c Define material 'AIRCOLD' at 0C ubuf(1) = sili_cg_tempc nmat = sili_med_coldair anumb = 14.61 znumb = 7.3 dens = 1.293e-3*273./(sili_cg_tempc+273.) radl = 0.283e5*(sili_cg_tempc+273.)/273 absl = 0.696e5*(sili_cg_tempc+273.)/273. CALL GFMATE(nmat,sil_name,anumb,znumb, ! First check if this material & dens,radl,absl,ubuf,nwbuf) ! number has already been used: if (anumb.ne. -1.0) goto 993 ! if so, abort. CALL GSMATE(nmat,'AIRCOLD',anumb,znumb, 1 dens,radl,absl,ubuf,nwbuf) nmed = sili_med_coldair ! Air at temperature sili_cg_tempc C isvol = 0 ! Not sensitive ifield = 1 ! magnetic field fieldm = 10. ! max field, kGs tmaxfd = 0.2 ! max angle due to field (one step) in degrees dmaxms = 0.1 ! max disp. due to mulsct. in one step, cm deemax = 0.01 ! max fractional energy loss in one step epsil = 0.001 ! tracking precision, cm stmin = 0.1 ! min step due to e-loss or multsct., cm ubuf(1) = 0. ! tracking stop switch CALL GSTMED(nmed,'AIRCOLD$',nmat,isvol,ifield,fieldm * ,tmaxfd,dmaxms,deemax,epsil,stmin,ubuf,nwbuf) C************************************************************************** v_i_name = siliCage v_m_name = siliEnvelope nr = 1 call gsvolu(v_i_name,'PCON',nmed,dim_sili,npar,ivolu) call gsatt(v_i_name,'SEEN',1) call gspos(v_i_name,nr,v_m_name,0.,0.,0.,irotnull,'ONLY') c Create and position thermal insulation of the beam-pipe made of Rohacell v_i_name = siliInCage v_m_name = siliEnvelope nmed = sili_med_cg ! Rohacell nr = 1 npar = 3 dim_sili(1) = sili_cg_rmn - sili_cg_inthck dim_sili(2) = sili_cg_rmn dim_sili(3) = 0.5*(sili_cg_z(sili_cg_npcon)-sili_cg_z(1)) dim_sili(4) = 0.5*(sili_cg_z(sili_cg_npcon)+sili_cg_z(1)) call gsvolu(v_i_name,'TUBE',nmed,dim_sili,npar,ivolu) call gsatt(v_i_name,'SEEN',1) call gspos(v_i_name,nr,v_m_name,0.,0.,dim_sili(4) * ,irotnull,'ONLY') c ====================================================================== c The parameters between === are for the temporary fake support emulation c with "plain" discs introduced for accomodating the request of c some people. You are welcome to spread this brand new "artificial c butter" over your piece of bread. Just, please, do not eat it! - VR c Define material/media 'GFRP' (slightly modified G10 from ../itr/pc1gem.f) nmat = sili_med_gfrp anumb = 18.14 znumb = 9.065 dens = 1.68 radl = 25. absl = 56.7 CALL GFMATE(nmat,sil_name,anumb,znumb, ! First check if this material & dens,radl,absl,ubuf,nwbuf) ! number has already been used: if (anumb.ne. -1.0) goto 993 ! if so, abort. CALL GSMATE(nmat,'GFRP',anumb,znumb,dens 1 ,radl,absl,ubuf,nwbuf) nmed = sili_med_gfrp isvol = 0 ! Not sensitive ifield = 1 ! magnetic field fieldm = 10. ! max field, kGs tmaxfd = 0.2 ! max angle due to field (one step) in degrees dmaxms = 0.1 ! max disp. due to mulsct. in one step, cm deemax = 0.01 ! max fractional energy loss in one step epsil = 0.001 ! tracking precision, cm stmin = 0.1 ! min step due to e-loss or multsct., cm ubuf(1) = 0. ! tracking stop switch CALL GSTMED(nmed,'GFRP$',nmat,isvol,ifield,fieldm * ,tmaxfd,dmaxms,deemax,epsil,stmin,ubuf,nwbuf) c Define material/media 'Carbon-carbon', for endcap composite panels nmat = sili_med_carbon ! = 123 anumb = 12.010 ! from the call to gpmate(0) znumb = 6.000 ! dens = 1.78 ! from Hytec report (default=2.265) radl = 23.9 ! scaled up from density absl = 63.5 ! scaled up from density CALL GFMATE(nmat,sil_name,anumb,znumb, ! First check if this material & dens,radl,absl,ubuf,nwbuf) ! number has already been used: if (anumb.ne. -1.0) goto 993 ! if so, abort. CALL GSMATE(nmat,'Carbon-carbon',anumb,znumb,dens & ,radl,absl,ubuf,nwbuf) nmed = sili_med_carbon isvol = 0 ! Not sensitive ifield = 1 ! magnetic field fieldm = 10. ! max field, kGs tmaxfd = 0.2 ! max angle due to field (one step) in degrees dmaxms = 0.1 ! max disp. due to mulsct. in one step, cm deemax = 0.01 ! max fractional energy loss in one step epsil = 0.001 ! tracking precision, cm stmin = 0.1 ! min step due to e-loss or multsct., cm ubuf(1) = 0. ! tracking stop switch CALL GSTMED(nmed,'Carbon-carbon$',nmat,isvol,ifield,fieldm * ,tmaxfd,dmaxms,deemax,epsil,stmin,ubuf,nwbuf) v_i_name = siliSupport v_m_name = siliCage call gsvolu(v_i_name,'TUBE',nmed,dim_sili,0,ivolu) call gsatt(v_i_name,'SEEN',1) npar = 3 Do nr = 1, sili_br_nspring dim_sili(1) = sili_br_sprin(nr) dim_sili(2) = sili_br_sprout(nr) dim_sili(3) = sili_br_spthck(nr) call gsposp(v_i_name,nr,v_m_name,0.,0.,sili_br_spz(nr) * ,irotnull,'ONLY',dim_sili,npar) Enddo c ====================================================================== c c... Build barrel VTX c Get passive material radiation length nmed = sili_med_passive CALL GFTMED(nmed,NAMTMED,nmat,isvol,ifield,fieldm * ,tmaxfd,dmaxms,deemax,epsil,stmin,ubuf,nwbuf) CALL GFMATE(nmat,NAMTMED,anumb,znumb,dens,pasv_radl,absl * ,ubuf,nwbuf) c Create a layer of passive material behind the sensor nmed = sili_med_passive npar = 0 call gsvolu(siliBrPassive,'BOX ',nmed,dim_sili,npar,ivolu) call GSATT(siliBrPassive,'SEEN',1) call GSATT(siliBrPassive,'COLO',3) c Create Si sensors SISN nmed = sili_med_silicon npar = 0 call gsvolu(siliBrSensor,'BOX ',nmed,dim_sili,npar,ivolu) call GSATT(siliBrSensor,'SEEN',1) call GSATT(siliBrSensor,'COLO',2) CALL GSATT(siliBrSensor,'WORK',1) ! Make volume sensitive c--- Cycle over the layers ( 4 barrel layers ) --------------------------------------- Do iLayer = 1, sili_br_nlayers c... Define ladder volume, SInn v_m_name = siliNames(iLayer) c Define ladder dimensions pasv_halfy = 0.5*sili_br_x0add(iLayer)*pasv_radl ! Al Thick/2 ladd_halfx = sili_br_snhalfx(iLayer) ! Width/2 dim_sili(1) = ladd_halfx ladd_halfthck = pasv_halfy + sili_br_snhalfy(iLayer) ! Thick/2 If(sili_br_snzgap(iLayer) .LT. 0.) Then ladd_halfthck = 2.*ladd_halfthck Endif dim_sili(2) = ladd_halfthck ladd_halfz = sili_br_nsn(iLayer)*sili_br_snhalfz(iLayer) * + 0.5*sili_br_snzgap(iLayer)*(sili_br_nsn(iLayer)-1) ! Length/2 dim_sili(3) = ladd_halfz c Create a ladder volume SInn made of cold air npar = 3 nmed = sili_med_coldair call gsvolu(v_m_name,'BOX ',nmed,dim_sili,npar,ivolu) call GSATT(v_m_name,'SEEN',0) call GSATT(v_m_name,'COLO',4) c Position sensors SISN and passive layers in the ladder SInn dim_sili(2) = sili_br_snhalfz(iLayer) - ladd_halfz npar = 3 dim_sili(3) = sili_br_snhalfx(iLayer) dim_sili(5) = sili_br_snhalfz(iLayer) sen_y(1) = sili_br_snhalfy(iLayer) - ladd_halfthck If(sili_br_snzgap(iLayer) .GE. 0.) Then sen_y(2) = sen_y(1) Else sen_y(2) = sen_y(1) + * 2.*(sili_br_snhalfy(iLayer) + pasv_halfy) Endif Do nr = 1, sili_br_nsn(iLayer) c Sensor v_i_name = siliBrSensor dim_sili(1) = sen_y(mod(nr-1,2)+1) dim_sili(4) = sili_br_snhalfy(iLayer) call gsposp(v_i_name,nr,v_m_name * ,0.,dim_sili(1),dim_sili(2),irotnull,'ONLY' * ,dim_sili(3),npar) c Passive layer v_i_name = siliBrPassive dim_sili(1) = dim_sili(1) + dim_sili(4) + pasv_halfy dim_sili(4) = pasv_halfy call gsposp(v_i_name,nr,v_m_name * ,0.,dim_sili(1),dim_sili(2),irotnull,'ONLY' * ,dim_sili(3),npar) dim_sili(2) = dim_sili(2) + 2.*sili_br_snhalfz(iLayer) * + sili_br_snzgap(iLayer) Enddo c--- Place ladders SInn in the cage SICG rladd = sili_br_r(iLayer) dphit = DEGRAD*sili_br_tilt(iLayer) xladd = ladd_halfthck - sili_br_snhalfy(iLayer) yladd = rladd + xladd*COS(dphit) xladd = xladd*SIN(dphit) rladd = sqrt(xladd*xladd + yladd*yladd) dphit = RADDEG*ATAN(xladd/rladd) dphid = sili_br_dphi(iLayer) dphir = DEGRAD*dphid v_i_name = v_m_name v_m_name = siliCage nr = 1 Do isec = 1, sili_br_nsec(iLayer) nladd = sili_br_nlad(isec,iLayer) philadd = sili_br_phic(isec,iLayer) + dphit * - (0.5*dphid)*(nladd-1) phirotladd = -90.+sili_br_tilt(iLayer) + philadd philadd = DEGRAD*philadd DO iladd = 1, nladd xladd = rladd*cos(philadd) yladd = rladd*sin(philadd) irot = irot + 1 CALL GSROTM(irot,90.,phirotladd,90. * ,phirotladd+90.,0.,0.) call gspos(v_i_name,nr,v_m_name * ,xladd,yladd,sili_br_z(iLayer),irot,'ONLY') philadd = philadd + dphir phirotladd = phirotladd + dphid nr = nr + 1 Enddo Enddo c put SInn in set 'SVX' c Note for the barrel v_i_name = SI01,2,3,4 namesv(4) = v_i_name call gsdet(set_id,v_i_name,nv,namesv,nbitsv,idtype,nwpa & ,nwsa,iset,idet) call gsdeth(set_id,v_i_name,nhh,inrNMSH,inrNBITSH & ,inrORIG,inrFACT) Enddo ! loop over 4 barrel layers c======================================================================c c============================= Endcaps code ===================c c======================================================================c panels=24. ! total # of carbon panels per "cone" section plength =14.4 ! overall length of each large carbon panel plength2=10.55 ! overall length of each medium carbon panel plength3= 6.8 ! overall length of each small carbon panel plen3=.94 ! length along x of small end on large carbon panel plen4=4.45 ! length along x of large end on large carbon panel panthk=.3 ! total thickness of a carbon panel (along y) strip=6.5 ! length (z) of inner(radius) carbon strips (= 5 x 1.3) strip2=.10 ! thickness (y) of inner carbon strip strip3=.38 ! width (x) of one carbon strip (same for all 3) strip4=7.8 ! length (z) of outer(radius) carbon strips (= 6 x 1.3) strip5=.025 ! thickness (y) of outer carbon strip (changed from 0.05) strip6=6.41 ! distance of outer carbon strips from sm end of panel strip7=3.9 ! length (z) of outer radius strip on small panel gap1=.11 ! tiny gap between carbon strip and edge of panel chiplen=1.3 ! total length of one chip (along z) chipwid=.38 ! total length of one chip (along x) chipthk=.025 ! total thickness of one chip (along y) sithk=.03 ! total thickness of silicon panel silen=6.6 ! overall length of INNER silicon panel silen2=1.6 ! length along x of large end of INNER Si panel silen3=.54 ! length along x of small end of INNER Si panel silen4=7.9 ! overall length of the larger, OUTER silicon panel silen5=2.7 ! length along x of large end of larger OUTER Si panel silen6=1.6 ! length along x of small end of larger OUTER Si panel silen7=4.0 ! overall length of the smaller, OUTER silicon panel silen8=2.16 ! length along x of large end of smaller OUTER Si panel silen9=1.6 ! length along x of small end of smaller OUTER Si panel rinner=3.5 ! inner "radius" of each "cone" dist from z-axis zangle=68. ! angle each panel makes with z-axis zsep=6. ! distance between "cone" sections zdist=19.86+18 ! z dist of outermost edge on outermost "cone" ! Note: zdist depends on zsep. zdist2=19.86 ! z dist of outermost edge on small "cones" sangle2=0.0 ! just a place holder for now C ----These values (below) are calculated from those entered above----- C C d = distance from z-axis to the midpoint of each large carbon panel C d2 = dist from z-axis to midpoint of each medium carbon panel C d3 = dist from z-axis to midpoint of each small carbon panel C halfz = the half length along z for each large carbon panel C halfz2 = half length along z for each medium carbon panel c 3 small panthick = panthk/2 ! half-z thickness of each carbon panel deg = 2*PI/360 ! To convert degrees into radians d = (plength /2)*SIN(zangle*DEG)+rinner d2 = (plength2/2)*SIN(zangle*DEG)+rinner d3 = (plength3/2)*SIN(zangle*DEG)+rinner halfz = (plength/2)*COS(zangle*DEG) halfz2 = (plength2/2)*COS(zangle*DEG) halfz3 = (plength3/2)*COS(zangle*DEG) rad = 360/(2*PI) ! To convert radians into degrees pangle = (360/panels) ! The total angle of each carbon panel sangle = 3.47433 ! angle of carbon strips stripz = -plength/2 + (gap1+strip/2)*COS(sangle*deg) !z-position of strips stripx = plen3/4+(gap1+strip/2)*SIN(sangle*deg) !x-position of strips plen5 = plen3 ! length along x of small end on medium carbon panel plen6 = 2*plength2*TAN(2*sangle*deg)+plen3 plen7 = plen3 ! length along x of small end on small carbon panel plen8 = 2*plength3*TAN(2*sangle*deg)+plen3 stripz2 = -plength2/2 + (gap1+strip/2)*COS(sangle*deg) stripx2 = plen5/4+(gap1+strip/2)*SIN(sangle*deg) stripz5 = -plength3/2 + (gap1+strip/2)*COS(sangle*deg) stripx5 = plen7/4+(gap1+strip/2)*SIN(sangle*deg) l1 = rinner/SIN(zangle*deg) stripx3 = plen3/4 +(strip6+strip4/2)*SIN(sangle*deg) stripy3 = (panthk+strip5)/2 stripz3 = -plength/2 + (strip6+strip4/2)*COS(sangle*deg) stripx4 = plen3/4 +(strip6+strip7/2)*SIN(sangle*deg) stripy4 = (panthk+strip5)/2 stripz4 = -plength2/2 + (strip6+strip7/2)*COS(sangle*deg) c============= Now define the endcap volumes: =============================c C (1) This is one large mother panel. The extras in x (.15 and .3 are so C that the overhanging silicon does not protrude from this mother volume. PAR(1) = plen3/2 + .15 ! half length along x at -z PAR(2) = plen4/2 + .3 ! half length along x at +z PAR(3) = sithk+chipthk+strip2+panthk/2 ! half thickness (y) PAR(4) = plength/2 ! half length along z CALL GSVOLU( 'SIPB', 'TRD1 ', sili_med_coldair, PAR, 4, IVOL1) PAR(1) = plen3/2 ! This is one large carbon panel inside SIPB PAR(2) = plen4/2 PAR(3) = panthk/2 PAR(4) = plength/2 CALL GSVOLU( 'SICB', 'TRD1 ', sili_med_carbon, PAR, 4, IVOL1) PAR_s1(1) = silen3/2 ! This is one inner silicon panel PAR_s1(2) = silen2/2 ! found in SIPB, SIPM and SIPS PAR_s1(3) = sithk/2 PAR_s1(4) = silen/2 CALL GSVOLU( 'SISI', 'TRD1 ', sili_med_silicon, PAR_s1, 0, IVOL1) ! 0 parameters PAR_s2(1) = silen6/2 ! This is one outer silicon panel in SIPB PAR_s2(2) = silen5/2 PAR_s2(3) = sithk/2 PAR_s2(4) = silen4/2 ! no GSVOLU call, done with GSPOSP PAR_s4(1) = silen9/2 ! This is one outer silicon panel in SIPM PAR_s4(2) = silen8/2 PAR_s4(3) = sithk/2 PAR_s4(4) = silen7/2 ! no GSVOLU call, done with GSPOSP PAR(1) = chipwid/2 ! This is one computer chip inside SIPB,M,S PAR(2) = chipthk/2 PAR(3) = chiplen/2 CALL GSVOLU( 'SCHP', 'BOX ', sili_med_silicon, PAR, 3, IVOL1) C (2) This is one medium mother panel. The extras in x (.15 and .3) are so C that the overhanging silicon does not protrude from this mother volume. PAR(1) = plen5/2 + 0.15 ! half length along x at -z PAR(2) = plen6/2 + 0.30 ! half length along x at +z PAR(3) = sithk+chipthk+strip2+panthk/2 ! half length along y PAR(4) = plength2/2 ! half length along z CALL GSVOLU( 'SIPM', 'TRD1 ', sili_med_coldair, PAR, 4, IVOL1) PAR(1) = plen5/2 ! This is one small carbon panel inside SIPM PAR(2) = plen6/2 PAR(3) = panthk/2 PAR(4) = plength2/2 CALL GSVOLU( 'SICM', 'TRD1 ', sili_med_carbon, PAR, 4, IVOL1) PAR(1) = strip3/2 ! This is the inner(radius)thin carbon PAR(2) = strip2/2 ! strip on the carbon panels PAR(3) = strip/2 CALL GSVOLU( 'SBX1', 'BOX ', sili_med_carbon, PAR, 3, IVOL1) PAR(1) = strip3/2 ! This is the outer(radius)thin carbon PAR(2) = strip5/2 ! strip on the carbon panels PAR(3) = strip4/2 CALL GSVOLU( 'SBX2', 'BOX ', sili_med_carbon, PAR, 3, IVOL1) C (3) This is one small mother panel. The extras in x (.15 and .3) are so C that the overhanging silicon does not protrude from this mother volume. PAR(1) = plen7/2 + 0.15 ! half length along x at -z PAR(2) = plen8/2 + 0.30 ! half length along x at +z PAR(3) = sithk+chipthk+strip2+panthk/2 ! half length along y PAR(4) = plength3/2 ! half length along z CALL GSVOLU( 'SIPS', 'TRD1 ', sili_med_coldair, PAR, 4, IVOL1) PAR(1) = plen7/2 ! This is one small carbon panel inside SIPS PAR(2) = plen8/2 PAR(3) = panthk/2 PAR(4) = plength3/2 CALL GSVOLU( 'SICS', 'TRD1 ', sili_med_carbon, PAR, 4, IVOL1) cx PAR(1) = 9.0 ! test box to show +x, +y and +z directions: cx PAR(2) = 9.0 cx PAR(3) = (sili_cg_z(6)-sili_cg_thck)/2. cx CALL GSVOLU( 'TEST', 'BOX ', , PAR, 3, IVOL1) c---- PCON mother volume to hold lampshades: ------------------------------------ shade_z = 1.4 ! lampshade thickness delz = (sili_cg_rmx(3) - sili_cg_thck & - sili_cg_rmn)/tan(zangle*deg) dely = shade_z * tan(zangle*deg) shade_zmid = (shade_z + delz) / 2.0 parb( 1) = 0 ! start angle parb( 2) = 360 ! full 360 cone parb( 3) = 4 ! 4 z's will be given parb( 4) = -shade_zmid ! z1 parb( 5) = sili_cg_rmn parb( 6) = sili_cg_rmn parb( 7) = parb(4) + shade_z ! z2 parb( 8) = sili_cg_rmn ! at inner radius parb( 9) = sili_cg_rmn + dely parb(10) = parb(7) + delz - shade_z ! z3 parb(11) = sili_cg_rmx(3) - sili_cg_thck - dely parb(12) = sili_cg_rmx(3) - sili_cg_thck ! outer radius parb(13) = parb(4) + delz + shade_z ! z4 parb(14) = sili_cg_rmx(3) - sili_cg_thck ! at outer radius parb(15) = sili_cg_rmx(3) - sili_cg_thck ! no GSVOLU call, done with GSPOSP ! Next for the medium shade mother volume: delz = (sili_cg_rmx(3) - sili_cg_thck - (plength-plength2) & - sili_cg_rmn)/tan(zangle*deg) dely = shade_z * tan(zangle*deg) shade_zmid = (shade_z + delz) / 2.0 parm( 1) = 0 ! start angle parm( 2) = 360 ! full 360 cone parm( 3) = 4 ! 4 z's will be given parm( 4) = -shade_zmid ! z1 parm( 5) = sili_cg_rmn parm( 6) = sili_cg_rmn parm( 7) = parm(4) + shade_z ! z2 parm( 8) = sili_cg_rmn ! at inner radius parm( 9) = sili_cg_rmn + dely parm(10) = parm(7) + delz - shade_z ! z3 parm(11) = sili_cg_rmx(3) - sili_cg_thck -dely -(plength-plength2) parm(12) = sili_cg_rmx(3) - sili_cg_thck -(plength-plength2) ! outer radius parm(13) = parm(4) + delz + shade_z ! z4 parm(14) = sili_cg_rmx(3) - sili_cg_thck -(plength-plength2) ! at outer radius parm(15) = sili_cg_rmx(3) - sili_cg_thck - (plength-plength2) ! Next for the small shade mother volume: delz = (sili_cg_rmx(3) - sili_cg_thck - (plength-plength3) & - sili_cg_rmn)/tan(zangle*deg) dely = shade_z * tan(zangle*deg) shade_zmid = (shade_z + delz) / 2.0 pars( 1) = 0 ! start angle pars( 2) = 360 ! full 360 cone pars( 3) = 4 ! 4 z's will be given pars( 4) = -shade_zmid ! z1 pars( 5) = sili_cg_rmn pars( 6) = sili_cg_rmn pars( 7) = pars(4) + shade_z ! z2 pars( 8) = sili_cg_rmn ! at inner radius pars( 9) = sili_cg_rmn + dely pars(10) = pars(7) + delz - shade_z ! z3 pars(11) = sili_cg_rmx(3) - sili_cg_thck -dely -(plength-plength3) pars(12) = sili_cg_rmx(3) - sili_cg_thck -(plength-plength3) ! outer radius pars(13) = pars(4) + delz + shade_z ! z4 pars(14) = sili_cg_rmx(3) - sili_cg_thck -(plength-plength3) ! at outer radius pars(15) = sili_cg_rmx(3) - sili_cg_thck - (plength-plength3) do ishade=5,12 ! define 8 copies SI05 - SI12 write (sil_name, '(''SI'',I2.2)') ishade if (ishade.eq.8 .or. ishade.eq.9) then ! small lampshades CALL GSVOLU(sil_name,'PCON', sili_med_coldair,pars,15,ivolu) elseif (ishade.eq.7 .or. ishade.eq.10) then ! medium lampshades CALL GSVOLU(sil_name,'PCON', sili_med_coldair,parm,15,ivolu) else ! big lampshades CALL GSVOLU(sil_name,'PCON', sili_med_coldair,parb,15,ivolu) endif enddo c------------------------------------------------------------------------------- c Hide some of the volumes c -2 = only that volume seen, no descendants c -1 = none visible c 0 = volume is not visible but descendants are c 1 = both volume and descendants are visible cx call gsatt( 'TEST', 'SEEN',1) ! test box showing +x+y+z do ishade=5,12 ! SI05, SI06 ... SI12 write (sil_name, '(''SI'',I2.2)') ishade call gsatt(sil_name, 'SEEN', 1) ! call gsatt(sil_name, 'COLO', 7) ! enddo CALL GSATT( 'SIPB', 'SEEN ', 0) ! big panels in big lampshade CALL GSATT( 'SIPM', 'SEEN ', 0) ! medium '' medium '' CALL GSATT( 'SIPS', 'SEEN ', 0) ! small '' small '' CALL GSATT( 'SICB', 'SEEN ', 1) ! Carbon panels big CALL GSATT( 'SICM', 'SEEN ', 1) ! medium CALL GSATT( 'SICS', 'SEEN ', 1) ! small CALL GSATT( 'SISI', 'SEEN ', 1) ! silicon CALL GSATT( 'SCHP ', 'SEEN ', 1) ! readout chips CALL GSATT( 'SBX1', 'SEEN ', 1) ! Carbon strips CALL GSATT( 'SBX2', 'SEEN ', 1) ! standoffs ! Add color to individual pieces CALL GSATT( 'SICB', 'COLO', 1) ! Carbon 1=black CALL GSATT( 'SICM', 'COLO', 1) ! CALL GSATT( 'SICS', 'COLO', 1) ! 2=red 5=yellow 8=white CALL GSATT( 'SBX1', 'COLO', 2) ! 3=green 6=magenta 17-166: CALL GSATT( 'SBX2', 'COLO', 2) ! 4=blue 7=lightblue many shades CALL GSATT( 'SISI', 'COLO', 6) ! silicon CALL GSATT( 'SCHP', 'COLO', 7) ! readout chips CALL GSATT( 'SISI','WORK',1) ! Make volumes sensitive C================= Position the detectors ======================c irot = irot + 1 ! for rotation of shades about y-axis irsh = irot ! x -> -x and z -> -z CALL GSROTM( irsh, 90., 180., 90., 90., 180., 0. ) do i = 0, panels-1 ! make 24 matrices, these irot = irot+1 ! will be used several times irx(i+1) = irot CALL GSROTM( irx(i+1), 90., (-360*I/panels) , 90+zangle , (90. + -360*I/panels) , zangle , (90-360*I/panels)) enddo do ishade=5,12 ! Place shades: loop north to south write (sil_name, '(''SI'',I2.2)') ishade ! SI05, SI06 ... SI12 z_shade = sili_endcap_z(ishade-4) ! z from par file irshade = 0 ! no rotation if (ishade.gt.8) irshade = irsh ! flip shades on North side CALL GSPOS (sil_name, 1, 'SICG', ! Place this shade & 0., 0., z_shade , irshade, 'ONLY') ! z = -zdist+(ishade-5)*zsep+halfz ! Now fill them with panels: if (ishade.eq.8 .or. ishade.eq.9) then ! small panel parameters: panel_name = 'SIPS' ! name dd = d3 ! center distance elseif (ishade.eq.7 .or. ishade.eq.10) then ! medium panel parameters: panel_name = 'SIPM' ! name dd = d2 ! center distance else ! big panel parameters: panel_name = 'SIPB' ! big panel name dd = d ! center endif do i = 0, panels-1 ! place panels in lampshades CALL GSPOS (panel_name , i+1, sil_name , & dd*SIN(2*PI*i/panels), & dd*COS(2*PI*i/panels), 0.0 , irx(i+1), 'ONLY') enddo ! end placing 24 panels in the shade enddo ! loop over all lampshades *----------------------------------------------------------------------------------* cx call gspos ('TEST',1,'SICG',9.0, 9.0, cx & (sili_cg_z(6)+sili_cg_thck)/2.0, 0, 'ONLY') ! test box showing +x,+y,+z * Place the inner(radius) carbon strips on the big carbon panels irot = irot+1 ir1 = irot CALL GSROTM( ir1, 90-sangle, 0. , 90. , 90. , sangle, 180.) CALL GSPOS( 'SBX1' , 1, 'SIPB' , -stripx , + (panthk/2)+strip2/2 , stripz ,ir1,'ONLY') irot = irot+1 ! and another copy on the other side ir2 = irot CALL GSROTM( ir2, 90+sangle, 0. , 90. , 90. , -sangle, 180.) CALL GSPOS( 'SBX1' , 2, 'SIPB' , stripx , + -(panthk/2)-strip2/2 , stripz ,ir2,'ONLY') * Place the outer(radius) carbon strips on the outer carbon surface CALL GSPOS( 'SBX2' , 1, 'SIPB' , -stripx3 , -stripy3 , + stripz3 ,ir1,'ONLY') CALL GSPOS( 'SBX2' , 2, 'SIPB' , stripx3 , stripy3 , + stripz3 ,ir2,'ONLY') * Place the inner-radius silicon inside SIPB CALL GSPOSP( 'SISI' , 1, 'SIPB' , -stripx , & panthk/2+strip2+chipthk+sithk/2 , stripz ,ir1,'ONLY', & par_s1, 4) CALL GSPOSP ( 'SISI' , 2, 'SIPB' , stripx , & -(panthk/2)-strip2-chipthk-sithk/2 , stripz ,ir2,'ONLY', & par_s1, 4) ! Place the outer-radius silicon in SIPB CALL GSPOSP ( 'SISI' , 3, 'SIPB' , stripx3 , (panthk/2)+ & (chipthk+sithk)/2+strip5 +chipthk/2, stripz3, ir2,'ONLY', & par_s2,4) CALL GSPOSP ( 'SISI' , 4, 'SIPB' , -stripx3, -panthk/2- & (chipthk+sithk)/2-strip5 -chipthk/2, stripz3, ir1,'ONLY', & par_s2,4) ! Place the readout chips in SIPB, do i=-5,5,2 ! 6 under one outer si chip: copies 1-6 CALL GSPOS( 'SCHP', (i+7)/2, 'SIPB' , & plen3/4 + (strip6 + strip4/2 + i*chiplen/2) * sin(sangle*deg), & (panthk/2) + strip5 + chipthk/2, & -plength/2 + (strip6 + strip4/2 + i*chiplen/2)*cos(sangle*deg), & ir2,'ONLY') enddo ! and 6 under the other outer si chip: do i=-5,5,2 ! copies 7-12 CALL GSPOS( 'SCHP', (i+19)/2, 'SIPB' , & -plen3/4 - (strip6 + strip4/2 + i*chiplen/2) * sin(sangle*deg), & -(panthk/2) - strip5 - chipthk/2, & -plength/2 + (strip6 + strip4/2 + i*chiplen/2)*cos(sangle*deg), & ir1,'ONLY') enddo ! chips under the inner-radius silicon do i=-2,2 ! copies 13-17 CALL GSPOS( 'SCHP' , 15+i, 'SIPB' , & -stripx - i*chiplen*sin(sangle*deg) , & (panthk/2) + strip2 + chipthk/2 , & stripz + i*chiplen*cos(sangle*deg), & ir1,'ONLY') enddo do i=-2,2 ! copies 18-22, inner-radius, other side CALL GSPOS( 'SCHP' , 20+i, 'SIPB' , & stripx + i*chiplen*sin(sangle*deg) , & -(panthk/2) - strip2 - chipthk/2 , & stripz + i*chiplen*cos(sangle*deg), & ir2,'ONLY') enddo C Place the Large carbon panels in the mother panel volume CALL GSPOS( 'SICB' , 1, 'SIPB', 0. , 0. , 0. ,0, 'ONLY') c---- now fill the medium panel mother volume: --------------------------------- * Place carbon panel in small lampshade panel mother volume CALL GSPOS( 'SICM' , 1, 'SIPM', 0. , 0. , 0. ,0,'ONLY') * Place the small carbon mount strips on the carbon panels (inner) CALL GSPOS( 'SBX1' , 3, 'SIPM' , -stripx2 , + panthk/2+strip2/2 , stripz2 ,ir1,'ONLY') * Place the small carbon strips on the carbon panels (outer) CALL GSPOS( 'SBX1' , 4, 'SIPM' , stripx2 , + -panthk/2-strip2/2 , stripz2 ,ir2,'ONLY') * Place the inner(radius) silicon in medium panel volume SIPM: CALL GSPOSP ( 'SISI' , 1, 'SIPM' , -stripx2 , & panthk/2 + strip2 + chipthk + sithk/2 , stripz2 ,ir1,'ONLY', & par_s1,4) CALL GSPOSP ( 'SISI' , 2, 'SIPM' , stripx2 , & -panthk/2 - strip2 - chipthk - sithk/2 , stripz2 ,ir2,'ONLY', & par_s1,4) * Place outer radius silicon directly in medium panel volume SIPM: CALL GSPOSP ( 'SISI' , 3, 'SIPM' , stripx4 , & panthk/2 + chipthk + sithk/2, stripz4, ir2,'ONLY', & par_s4,4) CALL GSPOSP ( 'SISI' , 4, 'SIPM' , -stripx4 , & -panthk/2 - chipthk - sithk/2, stripz4, ir1,'ONLY', & par_s4,4) ! Readout chips under the inner-radius silicon do i=-2,2 ! copies 23-27 CALL GSPOS( 'SCHP' , 25+i, 'SIPM' , & -stripx - i*chiplen*sin(sangle*deg) , & (panthk/2) + strip2 + chipthk/2 , & stripz2 + i*chiplen*cos(sangle*deg), & ir1,'ONLY') enddo do i=-2,2 ! inner, on the other (x) side, copies 28-32 CALL GSPOS( 'SCHP' , 30+i, 'SIPM' , & stripx + i*chiplen*sin(sangle*deg) , & -(panthk/2) - strip2 - chipthk/2 , & stripz2 + i*chiplen*cos(sangle*deg), & ir2,'ONLY') enddo ! outer-radius silicon in small panel SIPM do i=-1,1 ! copies 33-35 CALL GSPOS( 'SCHP' , 34+i, 'SIPM' , & stripx4 + i*chiplen*sin(sangle*deg) , & (panthk/2) + chipthk/2 , & stripz4 + i*chiplen*cos(sangle*deg), & ir2,'ONLY') enddo do i=-1,1 ! same, other (x) side, copies 36-38 CALL GSPOS( 'SCHP' , 37+i, 'SIPM' , & -stripx4 - i*chiplen*sin(sangle*deg) , & -(panthk/2) - chipthk/2, & stripz4 + i*chiplen*cos(sangle*deg), & ir1,'ONLY') enddo c---- now fill the small panel mother volume: --------------------------------- * Place carbon panel in small lampshade panel mother volume CALL GSPOS( 'SICS' , 1, 'SIPS', 0. , 0. , 0. ,0,'ONLY') * Place the small carbon mount strips on the carbon panels (inner) CALL GSPOS( 'SBX1' , 3, 'SIPS' , -stripx5 , + panthk/2+strip2/2 , stripz3 ,ir1,'ONLY') CALL GSPOS( 'SBX1' , 4, 'SIPS' , stripx5 , + -panthk/2-strip2/2 , stripz3 ,ir2,'ONLY') * Place the inner(radius) silicon in small panel volume SIPS: CALL GSPOSP ( 'SISI' , 1, 'SIPS' , -stripx5 , & panthk/2 + strip2 + chipthk + sithk/2 , stripz5 ,ir1,'ONLY', & par_s1,4) CALL GSPOSP ( 'SISI' , 2, 'SIPS' , stripx5 , & -panthk/2 - strip2 - chipthk - sithk/2 , stripz5 ,ir2,'ONLY', & par_s1,4) ! Readout chips under the inner-radius silicon do i=-2,2 ! copies 39-43 CALL GSPOS( 'SCHP' , 41+i, 'SIPS' , & -stripx - i*chiplen*sin(sangle*deg) , & (panthk/2) + strip2 + chipthk/2 , & stripz5 + i*chiplen*cos(sangle*deg), & ir1,'ONLY') enddo do i=-2,2 ! inner, on the other (x) side, copies 44-48 CALL GSPOS( 'SCHP' , 46+i, 'SIPS' , & stripx + i*chiplen*sin(sangle*deg) , & -(panthk/2) - strip2 - chipthk/2 , & stripz5 + i*chiplen*cos(sangle*deg), & ir2,'ONLY') enddo c================= declare the sensitive volumes ====================c * namesw = HALL, SIEN, SICG, SIxx, SIPy, SISI, where xx=05-12, y = B or S do ishade=5,12 write (namesw(4),'(''SI'',I2.2)') ishade namesw(5) = 'SIPB' if (ishade.ge.8 .and. ishade.le.9) namesw(5) = 'SIPM' call gsdet (set_id, namesw(4), 6, namesw, nbitsv, idtype, & nwpa, nwsa, iset, idet) call gsdeth(set_id, namesw(4), nhh,inrNMSH,inrNBITSH, & inrORIG,inrFACT) enddo ENDIF ! check on volume character from line 367 CVOLU_OPT... c--- Write parameters into the svxPISA.par file (temporary solution) c """"""""""""""""""""""""""""""""""""""""""""""""""""""""""""""" open(unit=15,file=svxpar_file,status='UNKNOWN',err=995) c SVX cage parameters write(15,*,err=994) 'SVX cage parameters:' write(15,*,err=994) sili_cg_rmn write(15,*,err=994) sili_cg_thck, sili_cg_inthck write(15,*,err=994) sili_cg_tempc write(15,*,err=994) sili_cg_npcon write(15,*,err=994) (sili_cg_z(iLayer), iLayer=1,sili_cg_npcon) write(15,*,err=994) (sili_cg_rmx(iLayer), iLayer=1,sili_cg_npcon) write(15,*,err=994) sili_cg_xdisp, sili_cg_ydisp, sili_cg_zdisp c SVX barrel parameters write(15,*,err=994) 'SVX barrel parameters:' write(15,*,err=994) sili_br_nlayers write(15,*,err=994) * (sili_br_r(iLayer), iLayer=1,sili_br_nlayers) write(15,*,err=994) * (sili_br_z(iLayer), iLayer=1,sili_br_nlayers) write(15,*,err=994) * (sili_br_nsn(iLayer), iLayer=1,sili_br_nlayers) write(15,*,err=994) * (sili_br_snhalfx(iLayer), iLayer=1,sili_br_nlayers) write(15,*,err=994) * (sili_br_snhalfy(iLayer), iLayer=1,sili_br_nlayers) write(15,*,err=994) * (sili_br_snhalfz(iLayer), iLayer=1,sili_br_nlayers) write(15,*,err=994) * (sili_br_x0add(iLayer), iLayer=1,sili_br_nlayers) write(15,*,err=994) * (sili_br_snzgap(iLayer), iLayer=1,sili_br_nlayers) write(15,*,err=994) * (sili_br_dphi(iLayer), iLayer=1,sili_br_nlayers) write(15,*,err=994) * (sili_br_tilt(iLayer), iLayer=1,sili_br_nlayers) do iLayer = 1, sili_br_nlayers write(15,*,err=994) sili_br_nsec(iLayer) write(15,*,err=994) * (sili_br_nlad(isec,iLayer),isec=1,sili_br_nsec(iLayer)) write(15,*,err=994) * (sili_br_phic(isec,iLayer),isec=1,sili_br_nsec(iLayer)) enddo c BARREL: Write sensor rotation matrices and translation vectors write(15,*,err=994) * 'SVX barrel sensor rotation matrices and translation vectors' call uctoh(namesv(1), LNAM(1), 4, 4) LNUM(1) = 1 call uctoh(namesv(2), LNAM(2), 4, 4) LNUM(2) = 1 call uctoh(namesv(3), LNAM(3), 4, 4) LNUM(3) = 1 do iLayer = 1, sili_br_nlayers ! loop over barrels 1-4 call uctoh(siliNames(iLayer), LNAM(4), 4, 4) nladd = 0 do isec = 1, sili_br_nsec(iLayer) nladd = nladd + sili_br_nlad(isec,iLayer) enddo do iladd = 1, nladd ! loop over ALL ladders LNUM(4) = iladd call uctoh(namesv(5), LNAM(5), 4, 4) do isen = 1, sili_br_nsn(iLayer) ! 4 Si dets par ladder LNUM(5) = isen write(15,*,err=994) 0, iLayer, iladd, isen NLEVEL = 0 ! nlevel is in /GCVOLU/ call glvolu(5, LNAM, LNUM, isec) ! fills /GCVOLU/ write (66,'(''5'',5a5,7i3)') namesv(1),namesv(2), & namesv(3),siliNames(iLayer),namesv(5), LNUM, isec write(15,*,err=994) (gtran(isec,NLEVEL), isec=1,3) ! in GCVOLU: xyz pos write(15,*,err=994) (grmat(isec,NLEVEL), isec=1,3) ! in GCVOLU: rotation write(15,*,err=994) (grmat(isec,NLEVEL), isec=4,6) ! matrix write(15,*,err=994) (grmat(isec,NLEVEL), isec=7,9) ! '' enddo enddo enddo NLEVEL = 0 c--------------------------------------------------! same for the endcaps: --------- write(15,*,err=994) 'SVX endcap parameters:' write(15,*,err=994) sili_endcap_layers write(15,*,err=994) & 'SVX Endcaps sensor rotation matrices and translation vectors' do ishade = 5,12 ! loop over 8 lampshades write (sil_name, '(''SI'',I2.2)') ishade ! call uctoh(sil_name, LNAM(4), 4, 4) ! lnum(4) = 1 do ipanel = 1,24 ! loop over 24 panels per shade lnum(5) = ipanel ! sil_name = 'SIPB' ! big and small panels if (ishade.eq.7 .or. ishade.eq.10) sil_name = 'SIPM' if (ishade.eq.8 .or. ishade.eq. 9) sil_name = 'SIPS' call uctoh(sil_name, LNAM(5), 4, 4) ! nsensors = 4 ! if (sil_name.eq.'SIPS') nsensors = 2 ! No outer silicon do isen = 1, nsensors ! front, back, and inner, outer radius lnum(6) = isen ! of the carbon panel write (15,'(5i3)', err=994) ishade, ipanel, isen call uctoh('SISI', LNAM(6),4,4) nlevel=0 ! for next call to glvolu: call glvolu(6, LNAM, LNUM, ierr) ! fills /GCVOLU/ write (66,'(''6'',6a5,7i3)') (lnam(i),i=1,6),LNUM,ierr write(15,*,err=994) (gtran(i, nlevel), i = 1,3) ! in GCVOLU: xyz pos write(15,*,err=994) (grmat(i, nlevel), i = 1,3) ! in GCVOLU: rotation write(15,*,err=994) (grmat(i, nlevel), i = 4,6) ! matrix write(15,*,err=994) (grmat(i, nlevel), i = 7,9) ! '' enddo ! 4 si detectors per panel enddo ! 24 panels per lampshade enddo ! loop over 8 lampshades nlevel = 0 c------------------------! end glvolu calls for the endcaps --------------- close(unit=15) ! file svxPISA.par c--- c--- Fill 'PARA' zebra-bank c--- CHFORM = '5I -F' ! 4 integer counts, then use all float call mzform('PARA',CHFORM,iod) ! book characteristic c c write the parameters to a zebra bank. later they will go to output file c--- Counting number of parameters npar = 1 ! Number of hit components c Contribution from SIEN/SICG/SIIC npar = npar + 8 + 2*sili_cg_npcon c Contribution from Barrel npar = npar + 2 + 11*sili_br_nlayers Do iLayer = 1, sili_br_nlayers npar = npar + 2*sili_br_nsec(iLayer) Enddo c Contribution from Endcap npar = npar + 2 + 11*sili_endcap_layers call mzbook(ixdiv_fr, lFD_PARA, lFD_PARA, 1, & 'PARA', 0, 0, npar, iod, 0) c C fill the bank c c Two first integers: numbers of layers in Barrel & Endcap iPoint = 1 iqf(lfd_para + iPoint) = sili_br_nlayers iPoint = iPoint + 1 iqf(lfd_para + iPoint) = sili_endcap_layers c Number of hit components iPoint = iPoint + 1 iqf(lfd_para + iPoint) = nhh c Number of barrel & endcap volume descriptors iPoint = iPoint + 1 iqf(lfd_para + iPoint) = nbrv ! = 5 iPoint = iPoint + 1 iqf(lfd_para + iPoint) = necv ! = 6 c c Envelope/Cage parameters c iPoint = iPoint + 1 qf(lfd_para + iPoint) = sili_cg_rmn iPoint = iPoint + 1 qf(lfd_para + iPoint) = sili_cg_thck iPoint = iPoint + 1 qf(lfd_para + iPoint) = sili_cg_inthck iPoint = iPoint + 1 qf(lfd_para + iPoint) = sili_cg_tempc iPoint = iPoint + 1 qf(lfd_para + iPoint) = FLOAT(sili_cg_npcon) Do icnt = 1, sili_cg_npcon iPoint = iPoint + 1 qf(lfd_para + iPoint) = sili_cg_z(icnt) iPoint = iPoint + 1 qf(lfd_para + iPoint) = sili_cg_rmx(icnt) Enddo iPoint = iPoint + 1 qf(lfd_para + iPoint) = sili_cg_xdisp iPoint = iPoint + 1 qf(lfd_para + iPoint) = sili_cg_ydisp iPoint = iPoint + 1 qf(lfd_para + iPoint) = sili_cg_zdisp C C Barrel parameters C Do iLayer = 1, sili_br_nlayers iPoint = iPoint + 1 qf(lfd_para + iPoint) = sili_br_snhalfx(iLayer) iPoint = iPoint + 1 qf(lfd_para + iPoint) = sili_br_snhalfy(iLayer) iPoint = iPoint + 1 qf(lfd_para + iPoint) = sili_br_snhalfz(iLayer) iPoint = iPoint + 1 qf(lfd_para + iPoint) = sili_br_x0add(iLayer) iPoint = iPoint + 1 qf(lfd_para + iPoint) = sili_br_snzgap(iLayer) iPoint = iPoint + 1 qf(lfd_para + iPoint) = sili_br_tilt(iLayer) iPoint = iPoint + 1 qf(lfd_para + iPoint) = FLOAT(sili_br_nsn(iLayer)) iPoint = iPoint + 1 qf(lfd_para + iPoint) = sili_br_r(iLayer) iPoint = iPoint + 1 qf(lfd_para + iPoint) = sili_br_z(iLayer) iPoint = iPoint + 1 qf(lfd_para + iPoint) = sili_br_dphi(iLayer) iPoint = iPoint + 1 qf(lfd_para + iPoint) = FLOAT(sili_br_nsec(iLayer)) Do isec = 1, sili_br_nsec(iLayer) iPoint = iPoint + 1 qf(lfd_para + iPoint) = sili_br_phic(isec,iLayer) iPoint = iPoint + 1 qf(lfd_para + iPoint) = FLOAT(sili_br_nlad(isec,iLayer)) Enddo Enddo C C Endcap parameters cxxx HvH: these parameters are obsolete, but I leave them as a placeholder for now. * cxxx do iLayer = 1,sili_endcap_layers cxx iPoint = iPoint + 1 cx qf(lfd_para + iPoint) = sili_phi1_side(iLayer) c iPoint = iPoint + 1 c qf(lfd_para + iPoint) = sili_dph_side(iLayer) c iPoint = iPoint + 1 c qf(lfd_para + iPoint) = sili_z1_side(iLayer) c iPoint = iPoint + 1 c qf(lfd_para + iPoint) = sili_rmin1_side(iLayer) c iPoint = iPoint + 1 c qf(lfd_para + iPoint) = sili_rmax1_side(iLayer) c iPoint = iPoint + 1 c qf(lfd_para + iPoint) = sili_z2_side(iLayer) c iPoint = iPoint + 1 c qf(lfd_para + iPoint) = sili_rmin2_side(iLayer) c iPoint = iPoint + 1 c qf(lfd_para + iPoint) = sili_rmax2_side(iLayer) c iPoint = iPoint + 1 c qf(lfd_para + iPoint) = sili_npdv_side(iLayer) c iPoint = iPoint + 1 c qf(lfd_para + iPoint) = sili_nz_side(iLayer) cx iPoint = iPoint + 1 cxx qf(lfd_para + iPoint) = sili_zCenter_side(iLayer) cxxx enddo return ! from subroutine svx 993 stop ' Material number already used.' 994 stop ' svxPISA.par write error.' 995 stop ' Unable to open svxPISA.par.' 996 stop ' Read error in sili_endcap_par segment of phnx.par.' 997 stop ' Read error in sili_br_par segment of phnx.par.' 998 stop ' Read error in sili_cg_par segment of phnx.par.' 999 stop ' Unable to open phnx.par.' end ! end of subroutine svx c=============================================================================c