C File name: svx.f ( was previously inr.f) C ------xxx---------- 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 support lampshadeshade c (SISR) with pipes and cables. Added cooling tubes to the barrel. c Hubert: added a 'dummy silicon plane' in the location of South and North Mutr c stations 1. This allows simple cuts on surviving particles. c Switched on/off with new phnx.par variable station1_dummy. Aug 2005 c Hubert: added option to install only 1/4 of the North endcaps, switched with c variable quarter_section in the phenx.par file c Hubert: Bug fix: declared SIP(S,M,B) to be 'MANY', since they are overlapping c Hubert: Added switch (sili_endcap_type) to make 'flat' endcaps. c 0=no endcaps, 1=umbrella endcaps, 1=flat endcaps. 06 Jan 2006 (w/ Sasha) c 3=titlted endcap pixel ('ldrd') pixel planes, 4=flat ldrd planes c Sasha L: Added misalignment for the barrel (07 September 2006) c Hubert: Added global in/out coordinates to the output (18 sep 2006) c Sasha L: staggered geometry for barrel strips (layers 3 and 4) c Hubert: removed 'station1' option, replaced by 'stagger' option to give stations c small, staggered rotations. (13 Dec 2006). c Hubert: Fixed S support cone orientation (sisr#2), moved cooling tubes to central c up/down sectors. c Hubert: Removed cooling tubes from SISR, 1 of 2 SJCC, SIIC, SBX1,2,3. Enlarged c SISP, thinned SICB,M,S 3->2mm. Enlarged SISP's. SICT 3->4mm. Jan 2007 c Added support posts SISQ for SISP's, and cooling manifolds SJSP and SJSQ c Removed code for sili_endcap_type = 1 c Introduced material 'honeycomb-1' for support structures. c Hubert: Enlarged the master volume SIEN to hold the readout wheels at +-Z. Added c readout wheel volume SIWH, with pc boards, connectors, cables rohacell. c *===================================================================================== 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) /10*0.0/ ! z-pos. of the Cage corners Real sili_cg_rmx(10)/10*0.0/ ! Outer SIEN radii at the corners, cm Real sili_cg_xdisp /0.0/ ! x-displacement of SIEN in HALL, cm Real sili_cg_ydisp /0.0/ ! y-displacement of SIEN in HALL, cm Real sili_cg_zdisp /0.0/ ! z-displacement of SIEN in HALL, cm namelist /sili_cg_par/ sili_cg_npcon,sili_cg_z $ ,sili_cg_rmn,sili_cg_rmx,sili_cg_thck,sili_cg_inthck $ ,sili_cg_xdisp,sili_cg_ydisp,sili_cg_zdisp $ ,sili_cg_tempc c """"""""""""""""""""""""""""""""""""""""""""""""""""""""""""" c--- VTX Barrel parameters c """"""""""""""""""""""""""""""""""""""""""""""""""""""""""""" Integer sili_br_nlayers /8/ ! Number of barrel layers Real sili_br_snhalfx(20) /.696,19*1.71555/ ! Si sensor width, cm Real sili_br_snhalfy(20) /0.01,19*0.02/ ! Si sensor thickness, cm Real sili_br_snhalfz(20) /2.836,19*3.2291/ ! Si sensor length, cm Real sili_br_x0add(20) /20*0.01/ ! Passive material thickness c ! in the ladders added on top c ! of Si sensor, RadLength X0 Real sili_br_snzgap(20) /20*0./ ! Sensor z-gap (if >=0) or c ! sensor z-overlap (if <0), cm c Real sili_br_tilt(20) /7.8,7.0,7.3,7.0,16*0./ ! OLD Ladder tilts Real sili_br_tilt(20) /7.8,7.0,0.0,0.0,16*0./ ! Ladder tilts Integer sili_br_nsn(20) /4,19*5/ ! Number of Si sensors/ladder Real sili_br_r(20) /2.5,5.,10.,14.,0.3,15*0./ ! 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 c ! The fifth parameter is staggering c ! in strip layers 3 and 4. c ! Set it to 0 to get old geometry 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/ real sili_br_misalignment(8) /8*0./ ! misalignment in X-Y plane (cm), east arm first real sili_br_shiftstag logical lbarrel /.true./ ! allows the barrel to be switched on/off integer sili_wheel_yesno /1/ ! model readout wheels c ====================================================================== c The parameters are for support rings at the end of the ladders 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 Real sili_br_manz(20) /-17.,17.,18*0./ ! Z-positions cooling manifold Real sili_br_manr(20) /20*2.4/ ! mean radius cooling manifold Real sili_br_mandia(20) /20*10./ ! diameter cooling manifold real siwh_pcb_thick ! pcb half-thickness = 1/32" real siwh_pcb_z(10) ! pcb center positions real siwh_pcb_rmax ! pcb max radius real siwh_pcb_connz ! pcb connector half-z real siwh_pcb_suppz(10) ! pcb support/cooling half-z real siwh_roha1_thick ! inner rohacell cover half-thickness 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 & ,sili_br_manz,sili_br_manr ,sili_br_mandia $ ,sili_br_misalignment, & sili_wheel_yesno, siwh_pcb_thick, siwh_pcb_z, siwh_pcb_rmax, & siwh_pcb_connz, siwh_pcb_suppz, siwh_roha1_thick c--- VTX Endcap parameters ---------------------------------------------- integer sili_endcap_type /2/ ! no endcaps (0) cones (1) or flat disks (2) ! tilted squares (3) or flat squares (4) integer sili_endcap_layers /8/ ! 4 south, 4 north real stagger /0.0/ ! turn on endcap staggering in phi real stag_ang(12) /12*0.0/ ! small rotations of endcap planes integer quarter_sector /0/ ! build only 1/4 of the endcaps real sili_endcap_z(8) / -35.16, -29.16, -23.16, -17.87, & 17.87, 23.16, 29.16, 35.16/ real sili_srvc_pos, sili_srvc_rmn, ct_radius, ct_zlayer1, & cc_thick, sili_endcap_angle, sili_endcap_zshade(8), & sili_endcap_zflat(8), sili_endcap_zldrd(4), panthk real yldrd real zldrd(4) / 20.0, 26.0, 32.0, 38.0 / namelist /sili_endcap_par/ & sili_endcap_type, sili_endcap_layers, & sili_endcap_zshade, sili_endcap_zflat, sili_endcap_zldrd, & sili_srvc_pos, sili_srvc_rmn, & ct_radius, ct_zlayer1, cc_thick, panthk, & stagger, quarter_sector 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 siliSupport /'SISP'/ ! Barrel support rings 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 /21/ ! Number of hit components integer*4 inrNBITSH(21) /21*32/ ! Bits for packing the hits c Hit component names character*4 inrNMSH(21) /'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 & ,'XIGL','YIGL','ZIGL','XOGL','YOGL','ZOGL'/ ! global entry & exit c Default setting of offsets and gains REAL inrORIG(21) /3*1000.,3*0.,3*1000.,6*1000.,6*1000./ ! offsets REAL inrFACT(21) /3*100000.,1.E7,1.e12,1.0,3*100000. & ,6*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 integer sili_med_coolant /124/ ! C5F12 liquid coolant integer sili_med_honeycomb /125/ ! 1/4" honeycomb, .5mm c-c skin, Al core 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(1),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(30),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 plength, plength2, plen3, plen4, & sbx1_thk, strip3, sbx2_thk, strip6, strip7, gap1, & chiplen, strip_z3, strip_z5, strip_z6, & chipwid, chipthk, sithk, silen, silen2, silen3, silen4, silen5, & silen6, silen7, silen8, silen9, rinner, zsep, zdist, & zdist2, sangle2, d, d2, dd, halfz, halfz2, deg, rad, & pangle, sangle, stripz, stripix, plen6, stripz2, c--------1---------1---------1---------1---------1---------1---------1-- & l1, stripz3, stripoxb, stripoy, stripz4, stripx4, stripy4, & par(20), parb(20), pars(20), par_s1(4), par_s2(4), & par_s4(4), par_s5(4), delz, dely, shade_z, shade_zmid, z_shade, & plen8, plength3, d3, halfz3, stripz5, parc(20), & sili_sr_rmx(10), srvc_zmid, sili_sr_rmn, srvc_z, sili_sr_thck, & zangle_srvc, plength4, d4, pard(9), tlength(4), stripozb, & stripoxm, stripozm, ztoshade, ytoshade, rmax_shades(8), rmax, & xcool, ycool, rcool, a_cool(2), z_cool(2), sien_delr, ystrip, & gap, rcool0, rladd0, posx, posy, a_honey(2), z_honey(2), & siwh_rmin, siwh_rmax real zangle /68.0/ ! lampshade angle wrt z=axis ( = 90-22) integer ivol1, ii, i, j, k, l, ir1, ir2, irx(24), irsh(12), ir, & ipanel, ierr, ishade, irshade, icopy, panels, ipstart, ipend, & shade_start, nparb, npars, nstr, irsisi, irottop, & irotbot character*4 sil_name, panel_name real rmyres,rmyrndm,rmydum integer rmyrndmi logical ltest ! xxx data ltest /.false./ ! WARNING do NOT run events with ltest .true.!!!! ! for display only (volume overlap will kill hits) C ================================================================ cxx integer ipart, itrtyp real amass, charge, tlife, ub, nwb character*12 chpart cxx 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 write(*,'(a16,5a5)') 'SVX CVOLU_OPT = ', & CVOLU_OPT(1,3),CVOLU_OPT(2,3),CVOLU_OPT(3,3), & CVOLU_OPT(4,3),CVOLU_OPT(5,3) IF(CVOLU_OPT(1,3).EQ.'BARR') continue write(*,*) ' sili_endcap_type = ', sili_endcap_type if (sili_endcap_type.eq.1) then stop '*** svx.f: type 1 no longer supported ***' endif IF(CVOLU_OPT(1,3).EQ.'FULL'.OR. & CVOLU_OPT(1,3).EQ.'VOLS'.OR. & CVOLU_OPT(1,3).EQ.'BARR') 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 sienxxx 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 ! number of user words in GSMATE calls 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) v_m_name = 'HALL' v_i_name = siliEnvelope ! = SIEN 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, a cylinder with a hole for the beampipe. dim_sili(1) = 0.0 dim_sili(2) = 360.0 dim_sili(3) = 2 dim_sili(4) = sili_cg_z(1) + sili_cg_thck dim_sili(5) = sili_cg_rmn dim_sili(6) = sili_cg_rmx(3) - sili_cg_thck dim_sili(7) = sili_cg_z(6) - sili_cg_thck dim_sili(8) = sili_cg_rmn dim_sili(9) = sili_cg_rmx(3) - sili_cg_thck C************************************************************************** c Define material 'AIRCOLD' at 0C nmat = sili_med_coldair CALL GFMATE(nmat,sil_name,anumb,znumb, ! First check if this material & dens,radl,absl,ubuf,nwbuf) ! number has already been used: * write (6,*)' aircold:' if (anumb.ne. -1.0) goto 993 ! if so, abort. 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. ubuf(1) = sili_cg_tempc nwbuf = 1 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,9,ivolu) call gsatt(v_i_name,'SEEN',1) call gspos(v_i_name,nr,v_m_name,0.,0.,0.,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 ! 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. anumb = 18.14 ! NOTE that an unsuccessful znumb = 9.065 ! call to GFMATE can clobber dens = 1.68 ! the returned variables, so we radl = 25. ! set them all AFTER calling GFMATE absl = 56.7 ! Jan 2007 HvH ubuf(1) = 0. nwbuf = 1 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 CALL GFMATE(nmat,sil_name,anumb,znumb, ! First check if this material & dens,radl,absl,ubuf,nwbuf) ! number has already been used: * write (6,*)' carbon:',ubuf,nwbuf if (anumb.ne. -1.0) goto 993 ! if so, abort. anumb = 12.01 ! from call to gpmate(0) znumb = 6.00 ! dens = 1.78 ! from Hytec report (default=2.265) radl = 23.9 ! scaled up from density absl = 63.5 ! scaled up from density ubuf(1) = 0. nwbuf = 1 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) c Define material/media 'Freon-coolant' ! C5F12 (from Hytec report) nmat = sili_med_coolant! = 124 CALL GFMATE(nmat,sil_name,anumb,znumb, ! First check if this material & dens,radl,absl,ubuf,nwbuf) ! number has already been used: * write (6,*)' coolant:' if (anumb.ne. -1.0) goto 993 ! if so, abort. If not, define: a_cool(1) = 12.010 ! Carbon z_cool(1) = 6.000 ! wmat (1) = 5 ! a_cool(2) = 18.998 ! Fluor z_cool(2) = 9.000 ! wmat (2) = 12 ! dens = 1.597 ! from Hytec report CALL GSMIXT(nmat,'Freon-coolant$', a_cool, z_cool, dens,-2,wmat) nmed = sili_med_coolant ! use parameters from carbon-carbon CALL GSTMED(nmed,'Freon-coolant$',nmat,isvol,ifield,fieldm & ,tmaxfd,dmaxms,deemax,epsil,stmin,ubuf,nwbuf) * write (6,*)' honeycomb:' nmat = sili_med_honeycomb ! = 125 a_honey(1) = 12.010 ! Carbon z_honey(1) = 6.000 ! wmat (1) = 0.97 ! 97% by weigh a_honey(2) = 26.98 ! Aluminum z_honey(2) = 13.000 ! wmat (2) = 0.03 ! 3% by weight dens = 0.251 ! scaled from c-c CALL GSMIXT(nmat,'Honeycomb-1$', a_honey, z_honey, dens,2,wmat) nmed = sili_med_honeycomb ! use parameters from carbon-carbon CALL GSTMED(nmed,'Honeycomb-1$',nmat,isvol,ifield,fieldm & ,tmaxfd,dmaxms,deemax,epsil,stmin,ubuf,nwbuf) cxx call gpmate(0) ! print all materials cxx cxx call gptmed(0) ! print all media nmed = sili_med_honeycomb ! used to be carbon-carbon v_i_name = siliSupport ! SISP are support rings v_m_name = siliCage call gsvolu(v_i_name,'TUBE',nmed,dim_sili,0,ivolu) call gsatt(v_i_name,'SEEN',1) ! SISQ are posts that connect rings to outside frame call gsvolu('SISQ','TRD1',nmed,dim_sili,0,ivolu) call gsatt(v_i_name,'SEEN',1) nmed = sili_med_coolant call gsvolu('SJSP','TUBE',nmed,dim_sili,0,ivolu) ! cooling tube circle call gsatt(v_i_name,'SEEN',1) call gsvolu('SJSQ','TUBE',nmed,dim_sili,0,ivolu) ! cooling tube feeds call gsatt(v_i_name,'SEEN',1) npar = 3 irot = irot+2 ! for top SISQ posts call gsrotm(irot,90.,0.,180.,0.,90.,90.) irottop = irot irot = irot+1 ! for bottom SISQ posts call gsrotm(irot,90.,0.,0.,0.,90.,-90.) irotbot = irot Do nr = 1, sili_br_nspring ! loop over 8 rings dim_sili(1) = sili_br_sprin(nr) ! support rings dim_sili(2) = sili_br_sprout(nr) dim_sili(3) = sili_br_spthck(nr) if (lbarrel) then call gsposp(v_i_name,nr,v_m_name,0.,0., & sili_br_spz(nr), & irotnull,'ONLY',dim_sili,npar) ! endif dim_sili(1) = 0.14* sili_br_sprout(nr) ! ring support posts/wedges dim_sili(2) = 2.56 ! = sili_cg_rmx(1)*tan(8 degrees) dim_sili(3) = sili_br_spthck(nr) dim_sili(4) = (sili_cg_rmx(3) -0.5 -sili_br_sprout(nr))/2-0.12 if (nr.eq.1.or.nr.eq.5) then ! special for innermost rings dim_sili(2) = 0.14*sili_br_sprin(2) ! shared support dim_sili(4) = (sili_br_sprin(2) - sili_br_sprout(1))/2-0.05 endif posy = sili_br_sprout(nr) + dim_sili(4) if (lbarrel) then call gsposp('SISQ',nr,v_m_name,0.,posy, & sili_br_spz(nr), & irottop,'ONLY',dim_sili,4) ! call gsposp('SISQ',nr+8,v_m_name,0.,-posy, & sili_br_spz(nr), & irotbot,'ONLY',dim_sili,4) ! endif dim_sili(1) = sili_br_manr(nr) - 0.443*sili_br_mandia(nr) ! cooling manifolds - rings dim_sili(2) = sili_br_manr(nr) + 0.443*sili_br_mandia(nr) dim_sili(3) = 0.443*sili_br_mandia(nr) if (lbarrel) then call gsposp('SJSP',nr,v_m_name,0.,0., & sili_br_manz(nr), & irotnull,'ONLY',dim_sili,npar) ! endif dim_sili(1) = 0.0 ! cooling manifolds - supply pipes dim_sili(2) = 0.5*sili_br_mandia(nr) dim_sili(3) = (sili_cg_rmx(3) - sili_br_manr(nr) - & sili_cg_thck - 0.443*sili_br_mandia(nr))/2.0-0.1 posx = 1.0 posy = sili_br_manr(nr) +0.443*sili_br_mandia(nr)+ dim_sili(3) if (lbarrel) then call gsposp('SJSQ',nr ,v_m_name, posx, posy, & sili_br_manz(nr), & irottop,'ONLY',dim_sili,4) ! call gsposp('SJSQ',nr+8 ,v_m_name, posx,-posy, & sili_br_manz(nr), & irotbot,'ONLY',dim_sili,4) ! call gsposp('SJSQ',nr+16,v_m_name,-posx, posy, & sili_br_manz(nr), & irottop,'ONLY',dim_sili,4) ! call gsposp('SJSQ',nr+24,v_m_name,-posx,-posy, & sili_br_manz(nr), & irotbot,'ONLY',dim_sili,4) ! endif Enddo C (4b) cooling tube (copied from further down, since copies are used in the barrel) call gsvolu( 'SICT', 'TUBE', sili_med_coolant, par, 0, ivol1) 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 c nmed = sili_med_passive ! this is Aluminum nmed = sili_med_carbon ! changed ladder support to c-c 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 ) --------------------------------------- if(sili_br_r(5).gt.0.) then write(*,*) ' SVX BARREL uses Staggered Geometry: ', * sili_br_r(5) endif icopy = 0 ! for cooling tube SICT 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 (SI01, SI02, SI03, SI04) 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 c Cooling tube parameters: par(1) = 0.0 par(2) = ct_radius par(3) = ladd_halfz ! this depends on the layer Do nr = 1, sili_br_nsn(iLayer) ! loop over sensors/stave c Sensor ! = 4, 4, 5, 6 v_i_name = siliBrSensor dim_sili(1) = sen_y(mod(nr-1,2)+1) dim_sili(4) = sili_br_snhalfy(iLayer) c active area is 3.6 mm narrower than passive layer if((iLayer.eq.3).or.(iLayer.eq.4)) then dim_sili(3) = sili_br_snhalfx(iLayer)-0.18 endif 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 c active area is 3.6 mm narrower than passive layer if((iLayer.eq.3).or.(iLayer.eq.4)) then dim_sili(3) = sili_br_snhalfx(iLayer) endif 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) ! partial thickness yladd = rladd + xladd*COS(dphit) xladd = xladd*SIN(dphit) rladd = sqrt(xladd*xladd + yladd*yladd) rcool = rladd + ladd_halfthck + ct_radius ! tube on the outside of the ladder 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 ! Note: we are inside of a loop over 4 layers Do isec = 1, sili_br_nsec(iLayer) ! 2 sectors: East, West nladd = sili_br_nlad(isec,iLayer) ! ladders per sector 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 ! loop over ladders in this sector C C------------ staggered geometry ----------------------------------------- C rladd0 = rladd + sili_br_shiftstag rcool0 = rcool + sili_br_shiftstag if((iLayer.eq.3).or.(iLayer.eq.4)) then sili_br_shiftstag = sili_br_r(5) * (-1.0)**iladd rladd0 = rladd + sili_br_shiftstag rcool0 = rcool + sili_br_shiftstag endif C C------------------------------------------------------------------------- C c c------------------------------------------------------------------------- c add misalignment rmyres = sili_br_misalignment((isec-1)*4+iLayer) ! resolution in cm if(rmyres.gt.0.) then rmydum = 10000. call gpoiss(rmydum,rmyrndmi,1) ! gaussian with mean 10000 and sigma 100 rmyrndm = rmyrndmi ! convert to real rmyrndm = (rmyrndm-10000.)/100. ! gaussian with mean 0 and sigma 1 write(*,*) ' SVX BARREL IS MISALLIGNED : ', * isec, iLayer, iladd, rmyrndm rmyrndm = rmyrndm * rmyres ! random shift in cm, assuming rmyres alignment accuracy rmyrndm = rmyrndm / sili_br_r(iLayer) ! random shift in radians ... philadd = philadd + rmyrndm ! shifted angle for each ladder endif c c------------------------------------------------------------------------- c xladd = rladd0*cos(philadd) yladd = rladd0*sin(philadd) xcool = rcool0 * cos(philadd) ycool = rcool0 * sin(philadd) irot = irot + 1 CALL GSROTM(irot,90.,phirotladd,90. * ,phirotladd+90.,0.,0.) if (lbarrel) then call gspos(v_i_name,nr,v_m_name ! place the ladder * ,xladd,yladd,sili_br_z(iLayer),irot,'ONLY') endif philadd = philadd + dphir phirotladd = phirotladd + dphid nr = nr + 1 ! now place cooling tube: icopy = icopy + 1 if (lbarrel) then call gsposp('SICT', icopy, v_m_name, & xcool,ycool,sili_br_z(iLayer), & irotnull, 'ONLY', par, 3) endif 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========================= Readout wheels code ===================c c See http://p25ext.lanl.gov/~hubert/phenix/silicon/simulations/jan07/wheel.html c======================================================================c if (sili_wheel_yesno.eq.1) then siwh_rmin = sili_cg_rmx(3) siwh_rmax = sili_cg_rmx(1) ! note Helium bag HEB2 starts at 53.6 cm nmed = sili_med_coldair dim_sili( 1) = 0 ! Holder volume for the readout wheels. dim_sili( 2) = 360 ! this is a hollow cylinder fitting dim_sili( 3) = 6 ! just around the Silicon Enclosure SIEN ! It consists of 2 sections, with a 0-thickness dim_sili( 4) = sili_cg_z(1)! central section (to clear the Helium bag). dim_sili( 5) = siwh_rmin dim_sili( 6) = siwh_rmax dim_sili( 7) = sili_cg_z(2) ! note in helium_bag.f, I moved HEB1 in to 26.0cm dim_sili( 8) = siwh_rmin dim_sili( 9) = siwh_rmax dim_sili(10) = sili_cg_z(3) dim_sili(11) = siwh_rmin dim_sili(12) = siwh_rmin dim_sili(13) = sili_cg_z(4) dim_sili(14) = siwh_rmin dim_sili(15) = siwh_rmin dim_sili(16) = sili_cg_z(5) dim_sili(17) = siwh_rmin dim_sili(18) = siwh_rmax dim_sili(19) = sili_cg_z(6) dim_sili(20) = siwh_rmin dim_sili(21) = siwh_rmax call gsvolu('SIWH', 'PCON', nmed, dim_sili, 21, ivolu) call gsatt ('SIWH', 'SEEN', 1) call gspos('SIWH', 1, 'SIEN', 0., 0., 0. ,irotnull,'ONLY') * wheels at 27.42 31.89 34.34 36.78 40.0 dim_sili( 1) = 0 ! dim_sili( 2) = 360 ! dim_sili( 3) = 2 ! ! dim_sili( 4) = -siwh_pcb_thick ! dim_sili( 5) = siwh_rmin dim_sili( 6) = siwh_pcb_rmax dim_sili( 7) = +siwh_pcb_thick ! dim_sili( 8) = siwh_rmin dim_sili( 9) = siwh_pcb_rmax nmed = sili_med_carbon call gsvolu('SWRB', 'PCON', nmed, dim_sili, 9, ivolu) call gsatt ('SWRB', 'SEEN', 1) call gsatt ('SWRB', 'COLO', 2) do i=1,10 call gspos('SWRB', i, 'SIWH', 0., 0., siwh_pcb_z(i), & irotnull,'ONLY') enddo * connectors on the boards, modeled as solid plastic rings: dim_sili( 1) = 0 ! dim_sili( 2) = 360 ! dim_sili( 3) = 2 ! ! dim_sili( 4) = -siwh_pcb_connz ! dim_sili( 5) = siwh_pcb_rmax - 1.0 dim_sili( 6) = siwh_pcb_rmax dim_sili( 7) = +siwh_pcb_connz ! dim_sili( 8) = siwh_pcb_rmax - 1.0 dim_sili( 9) = siwh_pcb_rmax nmed = 802 ! = plastic call gsvolu('SWCN', 'PCON', nmed, dim_sili, 9, ivolu) call gsatt ('SWCN', 'SEEN', 1) call gsatt ('SWCN', 'COLO', 4) do i=1,10 if (i.le.4.or.i.eq.6) then call gspos('SWCN', i, 'SIWH', 0., 0., & siwh_pcb_z(i) + siwh_pcb_thick + siwh_pcb_connz, & irotnull,'ONLY') else call gspos('SWCN', i, 'SIWH', 0., 0., & siwh_pcb_z(i) - siwh_pcb_thick - siwh_pcb_connz, & irotnull,'ONLY') endif enddo * support/cooling planes, one behind each pc board: nmed = sili_med_honeycomb call gsvolu('SWSP', 'PCON', nmed, dim_sili, 0, ivolu) call gsatt ('SWSP', 'SEEN', 1) call gsatt ('SWSP', 'COLO', 3) dim_sili( 1) = 0 ! dim_sili( 2) = 360 ! dim_sili( 3) = 2 ! ! * dim_sili( 4) = -siwh_pcb_suppz ! defined in each gposp() dim_sili( 5) = siwh_rmin dim_sili( 6) = siwh_pcb_rmax * dim_sili( 7) = +siwh_pcb_suppz ! dim_sili( 8) = siwh_rmin dim_sili( 9) = siwh_pcb_rmax do i=1,10 dim_sili( 4) = -siwh_pcb_suppz(i) ! dim_sili( 7) = +siwh_pcb_suppz(i) ! if (i.le.4.or.i.eq.6) then call gsposp('SWSP', i, 'SIWH', 0., 0., & siwh_pcb_z(i) - siwh_pcb_thick - siwh_pcb_suppz(i), & irotnull,'ONLY',dim_sili,9) write (6,*)' minus: ',i, siwh_pcb_suppz(i) else call gsposp('SWSP', i, 'SIWH', 0., 0., & siwh_pcb_z(i) + siwh_pcb_thick + siwh_pcb_suppz(i), & irotnull,'ONLY',dim_sili,9) write (6,*)' plus: ',i, siwh_pcb_suppz(i) endif enddo * one rohacell thermal disk as an inner-z cover: nmed = sili_med_cg call gsvolu('SWRH', 'PCON', nmed, dim_sili, 0, ivolu) call gsatt ('SWRH', 'SEEN', 1) call gsatt ('SWRH', 'COLO', 7) dim_sili( 1) = 0 ! dim_sili( 2) = 360 ! dim_sili( 3) = 2 ! ! dim_sili( 4) = -siwh_roha1_thick ! defined in each gposp() dim_sili( 5) = siwh_rmin dim_sili( 6) = siwh_pcb_rmax dim_sili( 7) = +siwh_roha1_thick ! dim_sili( 8) = siwh_rmin dim_sili( 9) = siwh_pcb_rmax call gsposp('SWRH', 1, 'SIWH', 0., 0., & siwh_pcb_z(5) + siwh_pcb_thick + 2.*siwh_pcb_suppz(i) + 0.2 & + siwh_roha1_thick, irotnull,'ONLY',dim_sili,9) call gsposp('SWRH', 2, 'SIWH', 0., 0., & siwh_pcb_z(6) - siwh_pcb_thick - 2.*siwh_pcb_suppz(i) - 0.2 & - siwh_roha1_thick, irotnull,'ONLY',dim_sili,9) * cables coming off the perimeter: dim_sili( 1) = 0 ! dim_sili( 2) = 360 ! dim_sili( 3) = 2 ! ! dim_sili( 4) = -5.75 dim_sili( 5) = siwh_pcb_rmax dim_sili( 6) = siwh_rmax dim_sili( 7) = +5.75 dim_sili( 8) = siwh_pcb_rmax dim_sili( 9) = siwh_pcb_rmax +0.5 nmed = 802 ! = plastic call gsvolu('SWCB', 'PCON', nmed, dim_sili, 9, ivolu) call gsatt ('SWCB', 'SEEN', 1) call gsatt ('SWCB', 'COLO', 6) call gspos('SWCB', 1, 'SIWH', 0., 0., & -34.25, irotnull,'ONLY') irot = irot+1 call gsrotm(irot,90.,180., 90.,90., 180.,0.) call gspos('SWCB', 2, 'SIWH', 0., 0., & +34.25, irot,'ONLY') nmed = sili_med_coldair ! reset to be sure endif ! readout wheels modeled or not c======================================================================c c============================= Endcaps code ===================c c======================================================================c panels=24 ! total # of carbon panels per "cone" section * if (ltest) panels = 1 plength =14.40 ! overall length of each large carbon panel plength2=10.55 ! overall length of each medium carbon panel plength3= 6.80 ! overall length of each small carbon panel plength4=14.0 ! '' service panel * plen3=.94 ! length along x of small end on large carbon panel plen3=.92 ! shrunk by .2mm so there is no overlap in the 'flat' mode plen4=4.45 ! length along x of large end on large carbon panel * panthk=.3 ! now in namelist: total thickness of a carbon panel (along z) sbx1_thk=.10 ! thickness (y) of inner carbon strip strip3=.38 ! width (x) of one carbon strip (same for all 3) sbx2_thk=.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 medium 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=.030 ! total thickness of one chip (along y) Dec 2006 HvH strip_z5=5.0*chiplen ! length (z) of inner-radius carbon strips (big) strip_z6=6.0*chiplen ! length (z) of outer-radius carbon strips strip_z3=3.0*chiplen ! ... outer-radius ... (medium) sithk=.0300 ! total thickness of silicon panel, doe detectors if (sili_endcap_type>=3) sithk=0.0300 ! ldrd detectors 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 if ( sili_endcap_type .eq. 3.) then zangle = 68.0 ! angle each panel makes with z-axis elseif (sili_endcap_type .eq. 2 .or. sili_endcap_type .eq. 4) then zangle = 90.0 ! flat planes endif 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 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 pangle = (360/24.0) ! The total angle of each carbon panel sangle = 3.47433 ! angle of carbon strips stripz = -plength /2 + (gap1+strip_z5/2)*COS(sangle*deg) ! z-position of strips stripz2 = -plength2/2 + (gap1+strip_z5/2)*COS(sangle*deg) stripz3 = -plength3/2 + (gap1+strip_z5/2)*COS(sangle*deg) stripix = plen3/4+(gap1+strip_z5/2)*SIN(sangle*deg) ! x-position of inner SI, SBX1 plen6 = 2*plength2*TAN(2*sangle*deg)+plen3 plen8 = 2*plength3*TAN(2*sangle*deg)+plen3 stripz5 = -plength3/2 + (gap1+strip_z5/2)*COS(sangle*deg) l1 = rinner/SIN(zangle*deg) stripoxb = plen3/4 +(strip6+strip_z6/2)*SIN(sangle*deg) stripoxm = plen3/4 +(strip6+strip_z3/2)*SIN(sangle*deg) stripoy = (panthk+sbx2_thk)/2 stripozb = -plength /2 + (strip6+strip_z6/2)*COS(sangle*deg) stripozm = -plength2/2 + (strip6+strip_z3/2)*COS(sangle*deg) stripx4 = plen3/4 +(strip6+strip7/2)*SIN(sangle*deg) stripy4 = (panthk+sbx2_thk)/2 stripz4 = -plength2/2 + (strip6+strip7/2)*COS(sangle*deg) do i=1,8 if (sili_endcap_type .eq. 2) then ! z-positions for flat disks: sili_endcap_z(i) = sili_endcap_zflat(i) endif enddo c============= Now define the endcap volumes: =============================c C--- Some items common to DOE and LDRD are needed first: 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 * write (15,'(''small si trapezoid: '',4f6.3)') (par_s1(i),i=1,4) if (stagger.gt.0) then stag_ang( 5) = 0. ! staggering angle in +phi direction stag_ang( 6) = 0.9375 ! 1x 360 / 96 / 4 stag_ang( 7) = 1.8750 ! 2x stag_ang( 8) = 2.8125 ! 3x stagger >0 PATTERN: stag_ang( 9) = 0. ! 0x (1 2 3 4 - 1 2 3 4) stag_ang(10) = 0.9735 ! 1x stag_ang(11) = 1.8750 ! 2x stag_ang(12) = 2.8125 ! 3x elseif(stagger.lt.0) then stag_ang( 8) = 0. ! staggering angle in +phi direction stag_ang( 6) = 0.9375 ! 1x 360 / 96 / 4 stag_ang( 5) = 1.8750 ! 2x stag_ang( 7) = 2.8125 ! 3x stagger <0 PATTERN: stag_ang( 9) = 0. ! 0x (3 2 4 1 - 1 4 2 3) stag_ang(11) = 0.9735 ! 1x stag_ang(12) = 1.8750 ! 2x stag_ang(10) = 2.8125 ! 3x endif do ishade=5,12 ! make a rot matrix for each shade irot = irot + 1 ! such that all face in the same direction, irsh(ishade) = irot ! and phi=0 starts at y=0, and goes positive. CALL GSROTM( irsh(ishade), & 90., 90. +stag_ang(ishade)*abs(stagger), & 90., 0. +stag_ang(ishade)*abs(stagger), & 180., 0. ) cxx & 90., 90. -stag_ang(ishade)*abs(stagger), ! reversed cxx & 90., 0. -stag_ang(ishade)*abs(stagger), cxx & 180., 0. ) enddo if (sili_endcap_type.eq.3 .or. sili_endcap_type.ge.4) goto 666 C---- This is the section for flat and conical full endcaps, sili_endcap_type 1 and 2: 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+sbx1_thk+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_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_passive, PAR, 3, IVOL1) ! med is Al 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) = plen3/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+sbx1_thk+panthk/2 ! half length along y PAR(4) = plength2/2 ! half length along z CALL GSVOLU( 'SIPM', 'TRD1 ', sili_med_carbon, PAR, 4, IVOL1) PAR(1) = plen3/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) 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) = plen3/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+sbx1_thk+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) = plen3/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) c---- PCON mother volume to hold lampshades: --------------------------!--------- shade_z = 1.0 ! lampshade thickness (down from 1.0) if (sili_endcap_type .eq. 2) then ! flat silicon planes parb( 1) = 0 ! start angle parb( 2) = 360 ! full 360 cone parb( 3) = 2 ! 2 z's will be given parb( 4) = -shade_z/2.0 ! z1 parb( 5) = sili_cg_rmn parb( 6) = sili_cg_rmx(3) - 0.1 - 0.5 parb( 7) = -parb(4) ! z2 parb( 8) = parb(5) ! at inner radius parb( 9) = parb(6) ! at inner radius nparb = 9 endif ! Next for the small shade mother volume: rmax = sili_cg_rmx(3) - sili_cg_thck - (plength-plength3) +0.6-0.5 delz = (rmax - sili_cg_rmn)/tan(zangle*deg) dely = shade_z * tan(zangle*deg) shade_zmid = (shade_z + delz) / 2.0 if (sili_endcap_type .eq. 2) then ! flat plane, zangle=90 pars( 1) = 0 ! start angle pars( 2) = 360 ! full 360 cone pars( 3) = 2 ! 2 z's will be given pars(4) = -shade_z/2.0 ! z1 pars(5) = sili_cg_rmn ! outer radius pars(6) = 10.5 ! outer radius (number by hand) pars(7) = shade_z/2.0 ! z2 pars(8) = sili_cg_rmn pars(9) = pars(6) - ytoshade npars = 9 endif if (quarter_sector.eq.1) then shade_start = 9 else shade_start = 5 endif do ishade=shade_start,12 ! define 4, 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,npars,ivolu) else ! big lampshades CALL GSVOLU(sil_name,'PCON',sili_med_coldair,parb,nparb,ivolu) endif enddo C================= Position the detectors ======================c 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+0.5)/panels) , & 90+zangle , (90. -360*(I+0.5)/panels) , & zangle , (90 -360*(I+0.5)/panels)) enddo if (sili_endcap_type.eq.2) then ! place endcaps do ishade = shade_start,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 = irsh(ishade) ! rotation of this shade 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 else ! big panel parameters: panel_name = 'SIPB' ! big panel name dd = d ! center endif if (quarter_sector.eq.1) then ! North quarter test configuration ipstart = 17 ! ipend = 22 ! 0-5 if (ishade.lt.9) ipend = -1 ! only on the N side else ! full shades, N and S ipstart = 0 ! changed from 0 for test only xxx ipend = panels-1 ! 0-23 endif ! do i = ipstart, ipend ! place panels in lampshades, 'MANY' !!! CALL GSPOS (panel_name , i+1, sil_name , & dd*SIN(2*PI*(i+0.5)/panels), & dd*COS(2*PI*(i+0.5)/panels), 0.0 , irx(i+1), 'MANY') enddo ! end placing 24 panels in the shade enddo ! loop over all lampshades endif ! if normal endcaps c------------------------------------------------------!---------------!------------- * 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.) irot = irot+1 ! and another copy on the other side ir2 = irot CALL GSROTM( ir2, 90+sangle, 0. , 90. , 90. , -sangle, 180.) * Place the inner-radius silicon inside SIPB CALL GSPOSP( 'SISI' , 1, 'SIPB' , -stripix , & panthk/2+sbx1_thk+chipthk+sithk/2 , stripz ,ir1,'ONLY', & par_s1, 4) CALL GSPOSP ( 'SISI' , 2, 'SIPB' , stripix , & -(panthk/2)-sbx1_thk-chipthk-sithk/2 , stripz ,ir2,'ONLY', & par_s1, 4) ! Place the outer-radius silicon in SIPB CALL GSPOSP ( 'SISI' , 3, 'SIPB' , stripoxb , panthk/2+ & sithk/2+sbx2_thk +chipthk, stripozb, ir2,'ONLY', & par_s2,4) CALL GSPOSP ( 'SISI' , 4, 'SIPB' , -stripoxb, -panthk/2- & (chipthk+sithk)/2-sbx2_thk -chipthk/2, stripozb, 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 + strip_z6/2 + i*chiplen/2) * sin(sangle*deg), & (panthk/2) + sbx2_thk + chipthk/2, -plength/2 + & (strip6 + strip_z6/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 c--------1---------1---------1---------1---------1---------1---------1-- CALL GSPOS( 'SCHP', (i+19)/2, 'SIPB' , -plen3/4 - & (strip6 + strip_z6/2 + i*chiplen/2) * sin(sangle*deg), & -(panthk/2) - sbx2_thk - chipthk/2, -plength/2 + & (strip6 + strip_z6/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' , & -stripix - i*chiplen*sin(sangle*deg) , & (panthk/2) + sbx1_thk + 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' , & stripix + i*chiplen*sin(sangle*deg) , & -(panthk/2) - sbx1_thk - 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 inner(radius) silicon in medium panel volume SIPM: CALL GSPOSP ( 'SISI' , 1, 'SIPM' , -stripix , & panthk/2 + sbx1_thk + chipthk + sithk/2 , stripz2,ir1,'ONLY', & par_s1,4) CALL GSPOSP ( 'SISI' , 2, 'SIPM' , stripix , & -panthk/2 - sbx1_thk - chipthk - sithk/2 , stripz2,ir2,'ONLY', & par_s1,4) * Place outer radius silicon directly in medium panel volume SIPM: c--------1---------1---------1---------1---------1---------1---------1-- CALL GSPOSP ( 'SISI' , 3, 'SIPM' , stripx4 , & panthk/2 + sbx2_thk + chipthk + sithk/2, stripz4, ir2,'ONLY', & par_s4,4) CALL GSPOSP ( 'SISI' , 4, 'SIPM' , -stripx4 , & -panthk/2 - sbx2_thk - 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' , & -stripix - i*chiplen*sin(sangle*deg) , & (panthk/2) + sbx1_thk + 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' , & stripix + i*chiplen*sin(sangle*deg) , & -(panthk/2) - sbx1_thk - chipthk/2 , & stripz2 + i*chiplen*cos(sangle*deg), & ir2,'ONLY') enddo ! outer-radius chips 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) + sbx2_thk + 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) - sbx2_thk - 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 inner(radius) silicon in small panel volume SIPS: c--------1---------1---------1---------1---------1---------1---------1-- CALL GSPOSP ( 'SISI' , 1, 'SIPS' , -stripix , & panthk/2 + sbx1_thk + chipthk + sithk/2 , stripz5,ir1,'ONLY', & par_s1,4) CALL GSPOSP ( 'SISI' , 2, 'SIPS' , stripix , & -panthk/2 - sbx1_thk - 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' , & -stripix - i*chiplen*sin(sangle*deg) , & (panthk/2) + sbx1_thk + 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' , & stripix + i*chiplen*sin(sangle*deg) , & -(panthk/2) - sbx1_thk - 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, M or S do ishade = shade_start,12 write (namesw(4),'(''SI'',I2.2)') ishade namesw(5) = 'SIPB' if (ishade.eq.8 .or. ishade.eq. 9) namesw(5) = 'SIPS' * write (6,*)' gsdet 1:' 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 c---- Hide some of the volumes, and set colors ----------------------!--------------------------- 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 do ishade = shade_start,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( 'SCHP ', 'SEEN ', 1) ! readout chips ! 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 *------ end of big cone/flat endcaps ------------------------------------------------------ 666 continue ! here are volumes common to the full DOE and the small LDRD detector: C (4) This is one service panel. PAR(1) = 2.0*pi*sili_srvc_rmn/48.0 ! half length along x at -z PAR(2) = 2.0*pi*(sili_cg_rmx(3) - sili_cg_thck -0.5)/48.0 -0.1 ! half length along x at +z srvc_z = 0.7 ! width of cone base (down from 1.5 par(3) = 0.15 ! half-thickness (down from 0.5) zangle_srvc = 50.0 ! angle of service shade with z-axis plength4 = (sili_cg_rmx(3) - sili_cg_thck & - sili_srvc_rmn)/sin(zangle_srvc*deg) & - PAR(3)/cos(zangle_srvc*DEG) -0.1 -0.5 PAR(4) = plength4/2.0 ! half length along z d4 = (plength4/2.0 + 0.5*PAR(3)/cos(zangle_srvc*DEG)) & *SIN(zangle_srvc*DEG)+sili_srvc_rmn CALL GSVOLU( 'SIPT', 'TRD1 ', sili_med_coldair, PAR, 4, IVOL1) CALL GSVOLU( 'SIQT', 'TRD1 ', sili_med_coldair, PAR, 4, IVOL1) C (4a) Thin carbon panel inside service panel par(3) = cc_thick/2.0 ! 2mm thick Carbon-Carbon CALL GSVOLU( 'SICC', 'TRD1 ', sili_med_carbon, PAR, 4, IVOL1) C (4b) cooling tube (moved further up, since copies are used in the barrel) c call gsvolu( 'SICT', 'TUBE', sili_med_coolant, par, 0, ivol1) C (5) This is one small service panel. PAR(1) = 2.0*pi*sili_cg_rmn /48.0 ! half length along x at +z PAR(2) = 2.0*pi*sili_srvc_rmn/48.0 ! half length along x at -z par(3) = 0.15 ! half thickness = y doen from 0.5 PAR(4) = 0.5*(sili_srvc_rmn-sili_cg_rmn) ! half length along z CALL GSVOLU( 'SJPT', 'TRD1 ', sili_med_coldair, PAR, 4, IVOL1) CALL GSVOLU( 'SJQT', 'TRD1 ', sili_med_coldair, PAR, 4, IVOL1) C Thin carbon panel inside service panel par(3) = cc_thick/2.0 ! 2mm thick Carbon-Carbon CALL GSVOLU( 'SJCC', 'TRD1 ', sili_med_carbon, PAR, 4, IVOL1) if (ltest) then PAR(1) = 9.0 ! test box to show +x, +y and +z directions: PAR(2) = 9.0 PAR(3) = (sili_cg_z(6)-sili_cg_thck)/2. CALL GSVOLU( 'TEST', 'BOX ', sili_med_coldair , PAR, 3, IVOL1) call gspos ('TEST',1,'SICG', par(1), par(2), par(3), & 0, 'ONLY') ! test box showing +x,+y,+z endif ! place carbon sheets in support panel: call gspos('SICC',1,'SIPT', 0.0,-cc_thick/2,0.0, irotnull, 'ONLY') call gspos('SICC',2,'SIPT', 0.0, cc_thick/2,0.0, irotnull, 'ONLY') call gspos('SICC',3,'SIQT', 0.0,-cc_thick/2,0.0, irotnull, 'ONLY') call gspos('SICC',4,'SIQT', 0.0, cc_thick/2,0.0, irotnull, 'ONLY') call gspos('SJCC',2,'SJPT', 0.0, cc_thick/2,0.0, irotnull, 'ONLY') call gspos('SJCC',4,'SJQT', 0.0, cc_thick/2,0.0, irotnull, 'ONLY') tlength(1) = plength4 tlength(2) = plength4 tlength(3) = 9.5 tlength(4) = 4.2 par(1) = 0.0 ! cooling tubes for layers 1-4: par(2) = ct_radius icopy = 0 do i=1,4 par(3) = tlength(i)/2.0 icopy = icopy+1 if (i.ne.1) then icopy = icopy+1 endif enddo par(3) = 0.5*(sili_srvc_rmn-sili_cg_rmn) ! tube for layer-1 pixels c(x) this is the barrel support structure cone -----------------------!---------------------- rmax = sili_cg_rmx(3) - sili_cg_thck + 0.15 - 0.7 delz = (rmax - sili_srvc_rmn)/tan(zangle_srvc*deg) dely = srvc_z * tan(zangle_srvc*deg) srvc_zmid = (srvc_z + delz) / 2.0 parc( 1) = 0 ! start angle parc( 2) = 360 ! full 360 cone parc( 3) = 5 ! 5 z's will be givenc parc( 4) = -srvc_zmid ! z1 parc( 5) = sili_cg_rmn parc( 6) = sili_srvc_rmn parc( 7) = parc(4) + srvc_z ! z2 parc( 8) = sili_cg_rmn ! at inner radius parc( 9) = sili_srvc_rmn + dely parc(10) = parc(4) + srvc_z ! z3 parc(11) = sili_srvc_rmn ! at inner radius parc(12) = sili_srvc_rmn + dely parc(13) = parc(7) + delz - srvc_z ! z4 parc(14) = rmax -dely parc(15) = rmax ! outer radius parc(16) = parc(4) + delz + srvc_z ! z5 parc(17) = rmax ! at outer radius parc(18) = rmax CALL GSVOLU('SISR','PCON', sili_med_coldair,parc,18,ivolu) do i = 0, panels-1 ! place support panels in cone *xxx if (i.ne.0.and.i.ne.11.and.i.ne.12.and.i.ne.23) then swap dec 2006 hvh if (i.eq.0. or.i.eq.11. or.i.eq.12. or.i.eq.23) then irot = irot+1 ! will be used several times CALL GSROTM(irot, & 90. , ( -360*(I+0.5)/panels), & 90.+zangle_srvc, (90. -360*(I+0.5)/panels), & zangle_srvc, (90. - 360*(I+0.5)/panels) ) CALL GSPOS ('SIPT' , i+1, 'SISR' , ! volume & d4*SIN(2*PI*(i+0.5)/panels), & d4*COS(2*PI*(i+0.5)/panels), 0.00 , irot, 'ONLY') irot = irot+1 ! will be used several times CALL GSROTM(irot, 90.0, ( -360.*(I+0.5)/panels), & 0.0, 0.0, & 90.0, (90.-360.*(I+0.5)/panels) ) CALL GSPOS ('SJPT' , i+1, 'SISR' , ! volume & 0.5*(sili_srvc_rmn+sili_cg_rmn)*SIN(2*PI*(i+0.5)/panels), & 0.5*(sili_srvc_rmn+sili_cg_rmn)*COS(2*PI*(i+0.5)/panels), & -srvc_zmid+srvc_z/2.0 , irot, 'ONLY') else ! 4 panels without tubes irot = irot+1 ! will be used several times CALL GSROTM(irot, 90. , ( -360*(I+0.5)/panels), & 90+zangle_srvc, (90. -360*(I+0.5)/panels), & zangle_srvc, (90. -360*(I+0.5)/panels) ) CALL GSPOS ('SIQT' , i+1, 'SISR' , & d4*SIN(2*PI*(i+0.5)/panels), & d4*COS(2*PI*(i+0.5)/panels), 0.0 , irot, 'ONLY') irot = irot+1 ! will be used several times CALL GSROTM(irot, 90.0, ( -360.*I/panels), & 0.0, 0.0, & 90.0, (90.-360.*I/panels) ) CALL GSPOS ('SJQT' , i+1, 'SISR' , & 0.5*(sili_srvc_rmn+sili_cg_rmn)*SIN(2*PI*i/panels), & 0.5*(sili_srvc_rmn+sili_cg_rmn)*COS(2*PI*i/panels), & -srvc_zmid+srvc_z/2.0 , irot, 'ONLY') endif ! tubes / no tubes enddo ! end placing 24 panels in the shade if (lbarrel) then CALL GSPOS ('SISR', 1, 'SICG', ! Place the services - north cone & 0., 0., sili_srvc_pos , irotnull, 'ONLY') irot = irot+1 ! rotate 180 around y CALL GSROTM(irot, 90.,180.,90.,90.,180.,0.) CALL GSPOS ('SISR', 2, 'SICG', ! Place the services - south cone & 0., 0., -sili_srvc_pos , irot, 'ONLY') endif c--- end build the common support structure cones ----------------------!----------------- c--- build the LDRD detector -------------------------------------------!----------------- if(sili_endcap_type.ge.3 .and. sili_endcap_type.le.5) then ! cxxx write (6,*)' install ldrd planes' par_s5(1) = 3.68 ! silicon parameters par_s5(2) = 3.68 ! the funny order is to duplicate the DOE par_s5(3) = sithk/2.0 ! (=type 1,2) silicon orientation par_s5(4) = 0.32 ! gap = 0.56 ! si-to-si gap delz = 0.1 ! c-to-c z-spacing between sensitive silicon irot = irot+1 ! z up locally in the silicon irsisi = irot ! y down the beampipe call gsrotm(irsisi, 90., 180., 0., 0., 90., 90.) ! x sideways irot = irot+1 ! zangle = 68 degrees call gsrotm(irot, 90.0, 0.0, & 180.0 -zangle, 90.0, & 90.0-zangle, 90.0) do ishade=9,12 ! place only on the North side if (sili_endcap_type.eq.4 .and. ishade.eq.9) then PAR(1) = par_s5(1) ! width master volume box nstr = 6 PAR(2) = ((par_s5(4)*2+gap)*nstr -gap/2 + par_s5(4))/2 ! height = 7.24 for 2x6 PAR(3) = 0.2 ! thickness (arbitrary) yldrd = 3.5 + (delz + sithk/2.0)*cos(zangle*deg) & + par(2)*sin(zangle*deg) write (sil_name, '(''SI'',I2.2)') ishade CALL GSVOLU(sil_name, 'BOX ', sili_med_coldair,PAR,3,IVOL1) call gsatt (sil_name, 'COLO',2) call gspos (sil_name,1,'SICG', & 0.0, yldrd, sili_endcap_zldrd(ishade-8), irot, 'ONLY') do i=1,nstr ystrip = (i-1)*1.2 - par(2) + par_s5(4) call gsposp('SISI', i, sil_name, * 0.0, ystrip, -delz, irsisi,'ONLY', * par_s5, 4) enddo do i=1,nstr ystrip = (i-1)*1.2 - par(2) + 2*par_s5(4) + gap/2 call gsposp('SISI', i+nstr, sil_name, * 0.0, ystrip, +delz, irsisi,'ONLY', * par_s5, 4) enddo else ! shades 10, 11, 12 ------------------------ PAR(1) = par_s5(1) ! width master volume box nstr = 10 ! # strips per side PAR(2) = ((par_s5(4)*2+gap)*nstr -gap/2 + par_s5(4))/2 ! height = 7.22 for 2x12 PAR(3) = 0.2 ! thickness yldrd = 3.5 + (delz + sithk/2.0)*cos(zangle*deg) & + par(2)*sin(zangle*deg) write (sil_name, '(''SI'',I2.2)') ishade CALL GSVOLU(sil_name, 'BOX ', sili_med_coldair,PAR,3,IVOL1) call gsatt (sil_name, 'COLO',2) call gspos (sil_name,1,'SICG', & 0.0, yldrd, sili_endcap_zldrd(ishade-8), irot, 'ONLY') do i=1,nstr ystrip = (i-1)*1.2 - par(2) + par_s5(4) call gsposp('SISI', i, sil_name, * 0.0, ystrip, -delz, irsisi,'ONLY', * par_s5, 4) enddo do i=1,nstr ystrip = (i-1)*1.2 - par(2) + 2*par_s5(4) + gap/2 call gsposp('SISI', i+nstr, sil_name, * 0.0, ystrip, +delz, irsisi,'ONLY', * par_s5, 4) enddo endif ! small(9) or big(10,11,12) shade enddo *------ declare the sensitive volumes -------------------------------- !--------------------------- * namesw = HALL, SIEN, SICG, SIxx, SISI, where xx=09-12 do ishade = 9,12 write (namesw(4),'(''SI'',I2.2)') ishade write (6,*)' gsdet 4:' call gsdet (set_id, namesw(4), 5, namesw, nbitsv, idtype, & nwpa, nwsa, iset, idet) call gsdeth(set_id, namesw(4), nhh,inrNMSH,inrNBITSH, & inrORIG,inrFACT) enddo endif !-------- end build the LDRD detector ----------------!--------------------------- c---- Hide some of the volumes, and set colors ----------------------!--------------------------- 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 CALL GSATT( 'SISI', 'SEEN ', 1) ! silicon CALL GSATT( 'SISI', 'COLO', 6) ! silicon CALL GSATT( 'SIPT', 'COLO', 7) ! service panel CALL GSATT( 'SIQT', 'COLO', 7) ! service panel CALL GSATT( 'SJPT', 'COLO', 7) ! service panel CALL GSATT( 'SJQT', 'COLO', 7) ! service panel CALL GSATT( 'SICC', 'COLO', 4) ! service panel CALL GSATT( 'SJCC', 'COLO', 4) ! service panel CALL GSATT( 'SISI','WORK',1) ! Make volumes sensitive c---- end colors and other attributes ENDIF ! this was a check on volume character from line ~420 above 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+1) 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: --------- if (sili_endcap_type.eq.0) goto 667 ! 0=no endcaps write(15,*,err=994) 'SVX endcap parameters:' write(15,*,err=994) sili_endcap_layers write (15,'(''small Si trapezoid: '',4f7.4)') (par_s1(i),i=1,4) write (15,'('' big Si trapezoid: '',4f7.4)') (par_s2(i),i=1,4) write(15,*,err=994) & 'SVX Endcaps sensor rotation matrices and translation vectors' if (sili_endcap_type.eq.2) then * write (6,*)' shade_start A: ', shade_start do ishade = shade_start,12 ! loop over 4 or 8 lampshades write (sil_name, '(''SI'',I2.2)') ishade ! call uctoh(sil_name, LNAM(4), 4, 4) ! lnum(4) = 1 * write (6,*)' ishade, ipstart, ipend: ', ishade, ipstart, ipend do ipanel = ipstart+1, ipend+1 ! loop over 6 or 24 panels per shade lnum(5) = ipanel ! sil_name = 'SIPB' ! big, small panels 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 ! 6, 24 panels per lampshade enddo ! loop over 4, 8 lampshades elseif (sili_endcap_type.eq.3 .or. sili_endcap_type.eq.4) then do ishade = 9,12 ! loop over 4 planes write (sil_name, '(''SI'',I2.2)') ishade ! 1,2,3 = HALL, SIEN, SICG call uctoh(sil_name, LNAM(4), 4, 4) ! 4 = SI09 (10,11,12) lnum(4) = 1 nsensors = nstr*2 ! 2x10 or 12 if (ishade.eq.9) nsensors = 12 ! 2x6 do isen = 1, nsensors ! lnum(5) = isen ! copy number of the silicon write (15,'(5i3)', err=994) ishade, isen call uctoh('SISI', LNAM(5),4,4) nlevel=0 ! for next call to glvolu: call glvolu(5, 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 ! loop over 4, 8 lampshades endif nlevel = 0 667 continue ! end skip endcaps 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 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 c HvH: these parameters are obsolete do ipart=1,52 call GFPART(ipart,chpart,itrtyp,amass,charge,tlife,ub,nwb) write (51,*)'GFPART:',ipart,chpart,itrtyp,amass,charge,tlife enddo write (51,*) '...' do ipart=1049,1075 call GFPART(ipart,chpart,itrtyp,amass,charge,tlife,ub,nwb) write (51,*)'GFPART:',ipart,chpart,itrtyp,amass,charge,tlife 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