summaryrefslogtreecommitdiff
path: root/generator/g4gen.py
blob: 9650b68eb18354045ad8601413804b9b5137752e (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
from Geant4 import *
import g4py.ezgeom
import g4py.NISTmaterials
import g4py.ParticleGun
import pyublas
import numpy as np
import chroma.event as event

try:
    import G4chroma
except:
    # Try building the module
    import subprocess
    import sys, os
    module_dir = os.path.split(os.path.realpath(__file__))[0]
    print >>sys.stderr, 'Compiling G4chroma.so...'
    retcode = subprocess.call('g++ -o \'%s/G4chroma.so\' -shared \'%s/G4chroma.cc\' -fPIC `geant4-config --cflags --libs` `python-config --cflags --libs --ldflags` -lboost_python' % (module_dir, module_dir), shell=True)
    assert retcode == 0
    import G4chroma

class G4Generator(object):
    def __init__(self, material, seed=None):
        '''Create generator to produce photons inside the specified material.

           material: chroma.geometry.Material object with density, 
                     composition dict and refractive_index.

                     composition dictionary should be 
                        { element_symbol : fraction_by_weight, ... }.
           seed: Random number generator seed for HepRandom. If None,
                 generator is not seeded.
        '''
        if seed is not None:
            HepRandom.setTheSeed(seed)

        g4py.NISTmaterials.Construct()
        g4py.ezgeom.Construct()
        self.physics_list = G4chroma.ChromaPhysicsList()
        gRunManager.SetUserInitialization(self.physics_list)
        self.particle_gun = g4py.ParticleGun.Construct()
        
        self.world_material = self.create_g4material(material)
        g4py.ezgeom.SetWorldMaterial(self.world_material)
        g4py.ezgeom.ResizeWorld(100*m, 100*m, 100*m)

        self.world = g4py.ezgeom.G4EzVolume('world')
        self.world.CreateBoxVolume(self.world_material, 100*m, 100*m, 100*m)
        self.world.PlaceIt(G4ThreeVector(0,0,0))

        self.tracking_action = G4chroma.PhotonTrackingAction()
        gRunManager.SetUserAction(self.tracking_action)
        gRunManager.Initialize()
        
        #preinitialize the process by running a simple event
        self.generate_photons(event.Event(event_id=0, particle_name='e-', 
                                          gen_position=(0,0,0), 
                                          gen_direction=(1,0,0), 
                                          gen_total_energy=1.0))
        
    def create_g4material(self, material):
        g4material = G4Material('world_material', material.density * g / cm3,
                                len(material.composition))

        # Add elements
        for element_name, element_frac_by_weight in material.composition.items():
            g4material.AddElement(G4Element.GetElement(element_name, True),
                                  element_frac_by_weight)

        # Set index of refraction
        prop_table = G4MaterialPropertiesTable()
        # Reverse entries so they are in ascending energy order rather
        # than wavelength
        energy = list((2*pi*hbarc / (material.refractive_index[::-1,0] * nanometer)).astype(float))
        values = list(material.refractive_index[::-1, 1].astype(float))
        prop_table.AddProperty('RINDEX', energy, values)

        # Load properties
        g4material.SetMaterialPropertiesTable(prop_table)
        return g4material

    def _extract_photons_from_tracking_action(self):
        n = self.tracking_action.GetNumPhotons()        
        pos = np.zeros(shape=(n,3), dtype=np.float32)
        pos[:,0] = self.tracking_action.GetX()
        pos[:,1] = self.tracking_action.GetY()
        pos[:,2] = self.tracking_action.GetZ()

        dir = np.zeros(shape=(n,3), dtype=np.float32)
        dir[:,0] = self.tracking_action.GetDirX()
        dir[:,1] = self.tracking_action.GetDirY()
        dir[:,2] = self.tracking_action.GetDirZ()

        pol = np.zeros(shape=(n,3), dtype=np.float32)
        pol[:,0] = self.tracking_action.GetPolX()
        pol[:,1] = self.tracking_action.GetPolY()
        pol[:,2] = self.tracking_action.GetPolZ()
        
        wavelengths = self.tracking_action.GetWavelength().astype(np.float32)
        times = self.tracking_action.GetT0().astype(np.float32)

        return event.Photons(positions=pos, directions=dir, polarizations=pol, times=times, wavelengths=wavelengths)

    def generate_photons(self, ev):
        '''Use GEANT4 to generate photons produced by the given particle.
           
               ev: a generator.event.Event object with the particle
                   properties set.  If it contains subtracks, those
                   will be used to create the photon vertices rather
                   than the main particle.

               Returns an instance of event.Photons containing the
               generated photon vertices for the primary particle or
               all the subtracks, if present.
        '''
        photons = []
        if ev.subtracks:
            subtracks = ev.subtracks
        else:
            # Create temporary subtrack for single primary particle
            subtracks = [event.Subtrack(particle_name=ev.particle_name,
                                        position=ev.gen_position,
                                        direction=ev.gen_direction,
                                        start_time=0.0,
                                        total_energy=ev.gen_total_energy)]

        for subtrack in subtracks:
            self.particle_gun.SetParticleByName(subtrack.particle_name)
            self.particle_gun.SetParticleEnergy(subtrack.total_energy * MeV)
            self.particle_gun.SetParticlePosition(G4ThreeVector(*subtrack.position)*m)
            self.particle_gun.SetParticleMomentumDirection(G4ThreeVector(*subtrack.direction).unit())
            
            self.tracking_action.Clear()
            gRunManager.BeamOn(1)
            photons.append(self._extract_photons_from_tracking_action())

        # Merge all photon lists into one big list
        return event.concatenate_photons(photons)

if __name__ == '__main__':
    import time
    import optics
    gen = G4Generator(optics.water)

    start = time.time()
    n = 0
    for i in xrange(100):
        photons = gen.generate_photons(event.Event('mu-', (0,0,0), (1,0,0), 1.0))
        n += len(photons.times)
        print photons.positions[0].min(), photons.positions[0].max()
    stop = time.time()
    print stop - start, 'sec'
    print n / (stop-start), 'photons/sec'