COSMOS v7.655  COSMOSv7655
(AirShowerMC)
cmyEfield.f
Go to the documentation of this file.
1  subroutine cmyefield(aTrack, Efout)
2 ! This is a template for giving electric field which is
3 ! complex and cannot be expressed by the default method.
4 ! This should be copied to your application area and
5 ! modified as you like. For compilaton, you may simply
6 ! change a line in chook.mk: i.e
7 ! objs = chook.o
8 ! to
9 ! objs = chook.o cmyEfield.o
10 
11 ! Before making your application, the library must be
12 ! updated; For that
13 ! Define MYEFIELD in Cosmos/cosmos/Zcondc.h
14 ! remove all in Cosmos/lib/$ARCH/
15 ! make veryclean
16 ! make
17 ! For actual run, you must specify HowEfield=2 in the
18 ! namelist ($HPARAM part).
19 !
20  implicit none
21 #include "Zglobalc.h"
22 #include "Ztrack.h"
23 #include "Zobs.h"
24 #include "Zobsv.h"
25 #include "Zincidentv.h"
26  type(track):: aTrack ! input. Track info in E-xyz.
27  ! position and "where" info (or others ) are used.
28  real(8),intent(out):: Efout(3) ! Elecric field vector in
29  ! the E-xyz system. unit is V/m.
30 ! This following information may be used to fix Efout
31 !
32 ! particle position.
33 ! aTrack.pos.height: vertical height in m from the sealevel
34 ! aTrack.pos.depth: depth in kg/m2.
35 ! aTrack.pos.xyz.r(1:3): xyz pos in E-xyz system.
36 ! to convert to other system,
37 ! see cdet2xyzFE.f and cxyz2det.f in Cosmos/Tracking.
38 ! Use routines cxyz2***. *** depends on the requirment.
39 !
40 ! aTrack.where: the observation level to which the ptcl
41 ! is directed.
42 ! The next observation point
43 ! ObsSites(aTrack.where).pos.height ! height in m
44 ! ObsSites(aTrack.where).pos.depth ! height in m
45 ! Deepest observation point
46 ! ObsSites(NoOfSites).pos.height
47 ! ObsSites(NoOfSites).pos.depth
48 ! distance r to the shower axis (perpendicular to the axis). r in m.
49 ! real(8)::r
50 ! call cgetPcoreDist(aTrack, r)
51 ! distance r to the shower axis (horizontal to the axis). r in m
52 ! call cgetHcoreDist(aTrack, r)
53 ! time information
54 ! aTrack.t/c*Tonsec in ns.
55 ! (c and Tonsec is deifined in Zglobalc.h; 0 is first col.
56 ! point)
57 ! to convert Eelectric field (Ef) in detector or primary system
58 ! to the Exyz system.
59 ! real(8):: Ef(3)
60 ! call cdet2xyzD(Ef, Efout) ! for detector to Exyz
61 ! call cprim2xyzD(Ef,Efout) ! for primary to Exyz
62  real(8)::h
63  real(8)::x,y
64  real(8):: Efz, Ef(3)
65  type(coord):: xyz
66 
67  h = atrack.pos.height
68 
69  if(h > 7e3 ) then
70  efz = 0.
71  else
72 ! convert x,y. origin is at the deepest det. pos
73  call cxyz2det( obssites(noofsites).pos.xyz,
74  * atrack.pos.xyz, xyz)
75  x = xyz.r(1)
76  y = xyz.r(2)
77 ! E eists only in some area
78  if( x > 0. .and. y > 0.) then
79  efz= 360d3*sin(3.1415/2.*h/2d3)*exp(-h/5.d3)
80  else
81  efz = 0.
82  endif
83  endif
84  ef(:) =(/0.d0, 0.d0, efz/)
85  call cdet2xyzd(ef, efout)
86  end subroutine cmyefield
subroutine cmyefield(aTrack, Efout)
Definition: cmyEfield.f:2
*Zfirst p fm *Zfirst p Zfirst p Zfirst p *Zfirst p *Zfirst pos xyz r
Definition: ZavoidUnionMap.h:1
block data cblkIncident data *Za1ry *HeightOfInj d3
Definition: cblkIncident.h:7
Definition: Ztrack.h:44
block data cblkEvhnp ! currently usable models data RegMdls ad *special data *Cekaon d0
Definition: cblkEvhnp.h:5
*Zfirst p fm *Zfirst p Zfirst p Zfirst p *Zfirst p *Zfirst pos xyz Zfirst pos xyz *Zfirst pos *Zfirst pos Zfirst pos height
Definition: ZavoidUnionMap.h:1
subroutine cdet2xyzd(a, b)
Definition: cxyz2det.f:101
block data cblkElemag data *AnihiE e3
Definition: cblkElemag.h:7
Definition: Zcoord.h:43
subroutine cxyz2det(det, a, b)
Definition: cxyz2det.f:7