File: fluorescence_mod.f90

package info (click to toggle)
mocassin 2.02.73.2-1
  • links: PTS, VCS
  • area: main
  • in suites: bookworm, bullseye, forky, sid, trixie
  • size: 42,116 kB
  • sloc: f90: 18,400; makefile: 75
file content (1871 lines) | stat: -rw-r--r-- 81,395 bytes parent folder | download | duplicates (3)
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
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
1502
1503
1504
1505
1506
1507
1508
1509
1510
1511
1512
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
1523
1524
1525
1526
1527
1528
1529
1530
1531
1532
1533
1534
1535
1536
1537
1538
1539
1540
1541
1542
1543
1544
1545
1546
1547
1548
1549
1550
1551
1552
1553
1554
1555
1556
1557
1558
1559
1560
1561
1562
1563
1564
1565
1566
1567
1568
1569
1570
1571
1572
1573
1574
1575
1576
1577
1578
1579
1580
1581
1582
1583
1584
1585
1586
1587
1588
1589
1590
1591
1592
1593
1594
1595
1596
1597
1598
1599
1600
1601
1602
1603
1604
1605
1606
1607
1608
1609
1610
1611
1612
1613
1614
1615
1616
1617
1618
1619
1620
1621
1622
1623
1624
1625
1626
1627
1628
1629
1630
1631
1632
1633
1634
1635
1636
1637
1638
1639
1640
1641
1642
1643
1644
1645
1646
1647
1648
1649
1650
1651
1652
1653
1654
1655
1656
1657
1658
1659
1660
1661
1662
1663
1664
1665
1666
1667
1668
1669
1670
1671
1672
1673
1674
1675
1676
1677
1678
1679
1680
1681
1682
1683
1684
1685
1686
1687
1688
1689
1690
1691
1692
1693
1694
1695
1696
1697
1698
1699
1700
1701
1702
1703
1704
1705
1706
1707
1708
1709
1710
1711
1712
1713
1714
1715
1716
1717
1718
1719
1720
1721
1722
1723
1724
1725
1726
1727
1728
1729
1730
1731
1732
1733
1734
1735
1736
1737
1738
1739
1740
1741
1742
1743
1744
1745
1746
1747
1748
1749
1750
1751
1752
1753
1754
1755
1756
1757
1758
1759
1760
1761
1762
1763
1764
1765
1766
1767
1768
1769
1770
1771
1772
1773
1774
1775
1776
1777
1778
1779
1780
1781
1782
1783
1784
1785
1786
1787
1788
1789
1790
1791
1792
1793
1794
1795
1796
1797
1798
1799
1800
1801
1802
1803
1804
1805
1806
1807
1808
1809
1810
1811
1812
1813
1814
1815
1816
1817
1818
1819
1820
1821
1822
1823
1824
1825
1826
1827
1828
1829
1830
1831
1832
1833
1834
1835
1836
1837
1838
1839
1840
1841
1842
1843
1844
1845
1846
1847
1848
1849
1850
1851
1852
1853
1854
1855
1856
1857
1858
1859
1860
1861
1862
1863
1864
1865
1866
1867
1868
1869
1870
1871
! Copyright (C) 2007 Barbara Ercolano
!
! Version 3.00
module fluorescence_mod

    use common_mod
    use constants_mod
    use continuum_mod
    use grid_mod
    use interpolation_mod
    use pathIntegration_mod
    use vector_mod


    type(vector), parameter :: origin=vector(0.,0.,0.)  ! origin of the cartesian grid axes

    integer     , parameter :: safeLim = 10000          ! safety limit for the loops

    contains

      recursive subroutine fluorescencePacketRun(grid, fType, rvec, xP, yP, zP, gP)
        implicit none


        real :: number
        real, save :: ionPhot = 0.

        type(vector), intent(inout)                  :: rVec        !

        integer, dimension(2), intent(inout)              :: xP, yP, &
             & zP                                            ! cartesian axes indeces
        ! 1= mother; 2=sub

        integer, intent(inout)           :: gP               ! grid index
        integer                          :: igpr             ! grid pointer 1= mother 2=sub
        integer                          :: difSourceL(3)    ! cell indeces
        integer                          :: err              ! allocation error status
        integer                          :: i, j             ! counters
        integer                          :: idirP, idirT     ! direction cosines

        type(grid_type), dimension(:), intent(inout) :: grid        ! the grid(s)

        type(photon_packet)              :: fluoPacket         ! the energu packet

        character(len=4), intent(inout)                   :: fType ! what resonance line?


        if (fType == 'diffuse') then
           ! the energy packet has beenid absorbed - reemit a new packet from this position
           call energyPacketRun(ftype, rVec, enPacket%xP(1:2), enPacket%yP(1:2), &
                & enPacket%zP(1:2), gP)
           return
        end if

        if (gP==1) then
           igpr = 1
        else if (gP>1) then
           igpr = 2
        else
           print*,  "! fluorescencePacketRun: insane grid index"
           stop
        end if

        difSourceL(1) = xP(igpr)
        difSourceL(2) = yP(igpr)
        difSourceL(3) = zP(igpr)

        fluoPacket = newFluorescencePacket(chType=fType, gP=gP, difSource=difSourceL, rvec=rvec)

        ! compute the next segment of trajectory
        call fluoPathSegment(fluoPacket)

        return


      contains


        ! this function initializes a photon packet
        function initFluorescencePacket(nuP,  position, lgLine, lgStellar, xP, yP, zP, gP)
          implicit none


          real                     :: random            ! random number
          type(vector), intent(in) :: position          ! the position at which the photon

          integer, intent(in)      :: nuP               ! the frequency of the photon packet
          integer, intent(in),dimension(2) :: xP, yP, &
               & zP                                     ! indeces of position on the x, y and z axes
          integer, intent(in)      :: gP                ! grid index
          integer                  :: igpi              ! grid pointer 1=mother, 2=sub
          integer                  :: i, irepeat        ! counter

          type(photon_packet)      :: initFluorescencePacket  ! the photon packet

          logical, intent(in)      :: lgLine, lgStellar ! line, stellar packet?



          initFluorescencePacket%position = position

          initFluorescencePacket%iG  = gP

          if (gP==1) then
             igpi=1
          else if (gp>1) then
             igpi=2
          else
             print*, "! initFluorescencePacket: insane gridp pointer"
             stop
          end if

          initFluorescencePacket%nuP      = nuP

          initFluorescencePacket%lgStellar = lgStellar

          ! check if photon packen is line or continuum photon
          if ( lgLine ) then
             ! line photon
             initFluorescencePacket%nu       = 0.
             initFluorescencePacket%lgLine   = .true.
          else
             ! continuum photon
             initFluorescencePacket%nu       = nuArray(nuP)
             initFluorescencePacket%lgLine   = .false.
          end if

          initFluorescencePacket%xP  = xP
          initFluorescencePacket%yP  = yP
          initFluorescencePacket%zP  = zP


          ! cater for plane parallel ionization case
          if (initFluorescencePacket%lgStellar .and. lgPlaneIonization) then

             ! get position

             ! x-direction
             call random_number(random)
             random = 1. - random
             initFluorescencePacket%position%x = &
                  & -(grid(gP)%xAxis(2)-grid(gP)%xAxis(1))/2. + random*( &
                  & (grid(gP)%xAxis(2)-grid(gP)%xAxis(1))/2.+&
                  & (grid(gP)%xAxis(grid(gP)%nx)-grid(gP)%xAxis(grid(gP)%nx-1))/2.+&
                  & grid(gP)%xAxis(grid(gP)%nx))
             if (initFluorescencePacket%position%x<grid(gP)%xAxis(1)) &
                  & initFluorescencePacket%position%x=grid(gP)%xAxis(1)
             if (initFluorescencePacket%position%x>grid(gP)%xAxis(grid(gP)%nx)) &
                  initFluorescencePacket%position%x=grid(gP)%xAxis(grid(gP)%nx)

             call locate(grid(gP)%xAxis, initFluorescencePacket%position%x, initFluorescencePacket%xP(igpi))
             if (initFluorescencePacket%xP(igpi) < grid(gP)%nx) then
                if (initFluorescencePacket%position%x >= (grid(gP)%xAxis(initFluorescencePacket%xP(igpi))+&
                     & grid(gP)%xAxis(initFluorescencePacket%xP(igpi)+1))/2.) &
                     & initFluorescencePacket%xP(igpi) = initFluorescencePacket%xP(igpi)+1
             end if

             ! y-direction
             initFluorescencePacket%position%y = 0.
             initFluorescencePacket%yP(igpi) = 1

             ! z-direction
             call random_number(random)
             random = 1. - random
             initFluorescencePacket%position%z = -(grid(gP)%zAxis(2)-grid(gP)%zAxis(1))/2. + random*( &
                  & (grid(gP)%zAxis(2)-grid(gP)%zAxis(1))/2.+&
                  & (grid(gP)%zAxis(grid(gP)%nz)-grid(gP)%zAxis(grid(gP)%nz-1))/2.+&
                  & grid(gP)%zAxis(grid(gP)%nz))
             if (initFluorescencePacket%position%z<grid(gP)%zAxis(1)) &
                  & initFluorescencePacket%position%z=grid(gP)%zAxis(1)
             if (initFluorescencePacket%position%z>grid(gP)%zAxis(grid(gP)%nz)) &
                  & initFluorescencePacket%position%z=grid(gP)%zAxis(grid(gP)%nz)

             call locate(grid(gP)%zAxis, initFluorescencePacket%position%z, initFluorescencePacket%zP(igpi))
             if (initFluorescencePacket%zP(igpi) < grid(gP)%nz) then
                if (initFluorescencePacket%position%z >= (grid(gP)%xAxis(initFluorescencePacket%zP(igpi))+&
                     & grid(gP)%zAxis(initFluorescencePacket%zP(igpi)+1))/2.) initFluorescencePacket%zP(igpi) =&
                     & initFluorescencePacket%zP(igpi)+1
             end if

             if (initFluorescencePacket%xP(igpi)<1) initFluorescencePacket%xP(igpi)=1
             if (initFluorescencePacket%zP(igpi)<1) initFluorescencePacket%zP(igpi)=1

             ! direction is parallel to y-axis direction
             initFluorescencePacket%direction%x = 0.
             initFluorescencePacket%direction%y = 1.
             initFluorescencePacket%direction%z = 0.

             if (initFluorescencePacket%xP(igpi) >  grid(gP)%xAxis(grid(gP)%nx) .or. &
                  & initFluorescencePacket%zP(igpi) >  grid(gP)%zAxis(grid(gP)%nz)) then
                print*, "! initFluorescencePacket: insanity in planeIonisation init"
                print*, igpi, initFluorescencePacket%xP(igpi),  grid(gP)%xAxis(grid(gP)%nx), &
                     & initFluorescencePacket%zP(igpi), grid(gP)%zAxis(grid(gP)%nz),  &
                     &random, initFluorescencePacket%position%z

                stop
             end if

             planeIonDistribution(initFluorescencePacket%xP(igpi),initFluorescencePacket%zP(igpi)) = &
                  & planeIonDistribution(initFluorescencePacket%xP(igpi),initFluorescencePacket%zP(igpi)) + 1

          else

             do irepeat = 1, 1000000
                ! get a random direction
                initFluorescencePacket%direction = randomUnitVector()
                if (initFluorescencePacket%direction%x/=0. .and. &
                     & initFluorescencePacket%direction%y/=0. .and. &
                     & initFluorescencePacket%direction%z/=0.) exit
             end do
          end if

          initFluorescencePacket%origin(1) = gP
          initFluorescencePacket%origin(2) = grid(gP)%active(initFluorescencePacket%xP(igpi),&
               & initFluorescencePacket%yP(igpi), initFluorescencePacket%zP(igpi))

        end function initFluorescencePacket


        ! this function creates a new fluorescence packet
        function newFluorescencePacket(chType, gP, difSource, rvec)

            type(photon_packet)                :: newFluorescencePacket! the photon packet to be created
            type(vector), intent(in)           :: rVec           !
            real                               :: random         ! random number

            type(vector)                       :: positionLoc    ! the position of the photon
                                                                 ! packet
            integer                            :: nuP            ! the frequency index of the photon packet
            integer, dimension(2)              :: orX,orY,orZ    ! dummy
            integer, intent(in)                :: difSource(3)   ! grid and cell indeces
            integer, intent(inout)             :: gP
            integer                            :: igpn           ! grid pointe 1=motehr, 2=sub

            character(len=4), intent(in)       :: chType         ! what line?




            if (gP==1) then
               igpn = 1
            else if (gp>1) then
               igpn = 2
            else
               print*,  "! newFluorescencePacket: insane grid pointer"
               stop
            end if

            positionLoc%x = rVec%x
            positionLoc%y = rVec%y
            positionLoc%z = rVec%z

            ! initialize the new fluorescence packet
            orX(igPn) = difSource(1)
            orY(igPn) = difSource(2)
            orZ(igPn) = difSource(3)

            ! initialize the new fluorescence packet
            if (grid(gP)%active(orX(igpn), orY(igpn), orZ(igpn)) < 0.) then
               print*, "! newFluorescencePacket: new packet cannot be emitted from re-mapped cell -3-"
               print*, "chType, nuP,  .false., .true., orX,orY,orZ, gp"
               print*, chType, nuP,  .false., .true., orX,orY,orZ, gp
               stop
            end if

            select case (chType)
            ! if the photon is stellar
            case ("FeKa")

               newFluorescencePacket = initFluorescencePacket(FeKaColdP, positionLoc, .false., .false., orX,&
                    & orY, orZ, gP)

            case ("FeL1")

               newFluorescencePacket = initFluorescencePacket(FeL1P, positionLoc, .false., .false., orX,&
                    & orY, orZ, gP)
            case ("FeL2")

               newFluorescencePacket = initFluorescencePacket(FeL2P, positionLoc, .false., .false., orX,&
                    & orY, orZ, gP)

            case ("CKa")

               newFluorescencePacket = initFluorescencePacket(CKaP, positionLoc, .false., .false., orX,&
                    & orY, orZ, gP)

            case ("NKa")

               newFluorescencePacket = initFluorescencePacket(NKaP, positionLoc, .false., .false., orX,&
                    & orY, orZ, gP)

            case ("OKa")

               newFluorescencePacket = initFluorescencePacket(OKaP, positionLoc, .false., .false., orX,&
                    & orY, orZ, gP)

            case ("NeKa")

               newFluorescencePacket = initFluorescencePacket(NeKaP, positionLoc, .false., .false., orX,&
                    & orY, orZ, gP)

            case ("MgKa")

               newFluorescencePacket = initFluorescencePacket(MgKaP, positionLoc, .false., .false., orX,&
                    & orY, orZ, gP)

            case ("AlKa")

               newFluorescencePacket = initFluorescencePacket(AlKaP, positionLoc, .false., .false., orX,&
                    & orY, orZ, gP)

            case ("SiKa")

               newFluorescencePacket = initFluorescencePacket(SiKaP, positionLoc, .false., .false., orX,&
                    & orY, orZ, gP)

            case ("SKa")

               newFluorescencePacket = initFluorescencePacket(SKaP, positionLoc, .false., .false., orX,&
                    & orY, orZ, gP)

            case ("ArKa")

               newFluorescencePacket = initFluorescencePacket(ArKaP, positionLoc, .false., .false., orX,&
                    & orY, orZ, gP)

            case ("CaKa")

               newFluorescencePacket = initFluorescencePacket(CaKaP, positionLoc, .false., .false., orX,&
                    & orY, orZ, gP)

            case ("NiKa")

               newFluorescencePacket = initFluorescencePacket(NiKaP, positionLoc, .false., .false., orX,&
                    & orY, orZ, gP)

            ! if the photon packet type is wrong or missing
            case default

                print*, "! newFluorescencePacket: unsupported fluorescent transition.", chType
                stop

             end select

           end function newFluorescencePacket

           subroutine fluoPathSegment(enPacket)
             implicit none


             ! local variables
             type(vector)                    :: vHat     ! direction vector
             type(vector)                    :: rVec     ! position vector

             real                            :: absTau   ! optical depth
             real                            :: dlLoc    ! local displacement
             real                            :: dx, dy, dz
             real                            :: dSx, dSy, dSz
                                                      ! distances from x,y and z wall
             real                            :: dS       ! distance from nearest wall
             real                            :: dV       ! lume of this cell
             real                            :: passProb ! prob of passing the next segment
             real                            :: probSca  ! prob that the packet scatters
             real                            :: radius   ! radius
             real                            :: random   ! random number
             real                            :: tauCell  ! local tau

             integer                         :: idirT,idirP ! direction cosine counters
             integer                         :: i, j, nS ! counter
             integer                         :: xP,yP,zP ! cartesian axes indeces
             integer                         :: gP       ! grid index
             integer                         :: igpp     ! grid index 1=mother 2=sub
             integer                         :: safeLimit =1000 ! safe limit for the loop

             type(photon_packet), intent(inout) :: enPacket ! the energy packet

             character(len=7)                :: packetType ! what line?

             logical                         :: lgScattered ! is the packet scattering with dust?
             logical                         :: lgReturn

             if (enPacket%iG == 1) then
                igpp = 1
             else if (enPacket%iG>1) then
                igpp = 2
             else
                print*, "! fluoPathSegment: insane grid index"
                stop
             end if

             ! check that the input position is not outside the grid
             if ( (enPacket%iG <= 0).or.(enPacket%iG > nGrids) ) then
                print*, "! fluoPathSegment: starting position not in any defined gridhses",&
                     & enPacket
                stop
             else if ( (enPacket%xP(igpp) <= 0).or.&
                  &(enPacket%xP(igpp) > grid(enPacket%iG)%nx) ) then
                print*, "! fluoPathSegment: starting position in x is outside the grid",&
                     & enPacket
                stop
             else if ( (enPacket%yP(igpp) <= 0).or. &
                  & (enPacket%yP(igpp) > grid(enPacket%iG)%ny) ) then
                print*, "! fluoPathSegment: starting position in y is outside the grid",&
                     & enPacket
                stop
             else if ( (enPacket%zP(igpp) <= 0).or.&
                  & (enPacket%zP(igpp) > grid(enPacket%iG)%nz) ) then
                print*, "! fluoPathSegment: starting position in z is outside the grid",&
                     & enPacket
                stop
             end if

             ! define vHat and rVec
             rVec = enPacket%position
             vHat = enPacket%direction

             ! initialize xP, yP,zP
             xP = enPacket%xP(igpp)
             yP = enPacket%yP(igpp)
             zP = enPacket%zP(igpp)
             gP = enPacket%iG

             ! initialise distance from walls
             dSx = 0.
             dSy = 0.
             dSz = 0.

             if (lg1D) then
                radius = 1.e10*sqrt((rVec%x/1.e10)*(rVec%x/1.e10) + &
                     &                               (rVec%y/1.e10)*(rVec%y/1.e10) + &
                     &                               (rVec%z/1.e10)*(rVec%z/1.e10))
                call locate(grid(1)%xAxis, radius, xP)
                if (nGrids > 1 .or. gP >1) then
                   print*, " ! fluorescencePacketRun: multiple grids are not allowed in a 1D simulation"
                   stop
                end if
             end if

             ! initialize optical depth
             absTau = 0.

             ! get a random number
             call random_number(random)

             ! calculate the probability
             passProb = -log(1.-random)

             ! speed up photons that my be trapped
             if (lgPlaneIonization) then
                safeLimit=5000
             else
!             safeLimit=500000
                safeLimit=1000
             end if

             do i = 1, safeLimit

                do j = 1, safeLimit

                   if (xP > grid(gP)%nx .or. xP < 1 .or. &
                        & yP > grid(gP)%ny .or. yP < 1 .or. &
                        & zP > grid(gP)%nz .or. zP < 1 ) then
                      print*, "! fluoPathSegment: insanity [gp,xp,yp,zp,j,i]", &
                           & gp, xp, yp, zp, j, i
                      stop
                   end if

                   if (grid(gP)%active(xP,yP,zP)<0) then

                      ! packet is entering a subgrid
                      enPacket%xP(1) = xP
                      enPacket%yP(1) = yP
                      enPacket%zP(1) = zP

                      gP = abs(grid(gP)%active(xP,yP,zP))

                      ! where is the packet in the sub-grid?

                      call locate(grid(gP)%xAxis, rVec%x, xP)
                      if (xP==0) xP = xP+1
                      if (xP< grid(gP)%nx) then
                         if (rVec%x >  (grid(gP)%xAxis(xP+1)+grid(gP)%xAxis(xP))/2.) &
                              & xP = xP + 1
                      end if

                      call locate(grid(gP)%yAxis, rVec%y, yP)
                      if (yP==0) yP=yP+1
                      if (yP< grid(gP)%ny) then
                         if (rVec%y >  (grid(gP)%yAxis(yP+1)+grid(gP)%yAxis(yP))/2.) &
                              & yP = yP + 1
                      end if

                      call locate(grid(gP)%zAxis, rVec%z, zP)
                      if (zP==0) zP=zP+1
                      if (zP< grid(gP)%nz) then
                         if (rVec%z >  (grid(gP)%zAxis(zP+1)+grid(gP)%zAxis(zP))/2.) &
                              & zP = zP + 1
                      end if

                   end if

                   enPacket%iG = gP
                   igpp = 2

                   ! find distances from all walls

                   if (lgSymmetricXYZ) then
                      if ( rVec%x <= grid(1)%xAxis(1) ) then
                         if (vHat%x<0.) vHat%x = -vHat%x
                         rVec%x = grid(1)%xAxis(1)
                      end if
                      if ( rVec%y <= grid(1)%yAxis(1) ) then
                         if (vHat%y<0.) vHat%y = -vHat%y
                         rVec%y = grid(1)%yAxis(1)
                      end if
                      if ( rVec%z <= grid(1)%zAxis(1) ) then
                         if (vHat%z<0.) vHat%z = -vHat%z
                         rVec%z = grid(1)%zAxis(1)
                      end if
                   end if

                   if (vHat%x>0.) then
                      if (xP<grid(gP)%nx) then

                         dSx = ( (grid(gP)%xAxis(xP+1)+grid(gP)%xAxis(xP))/2.-rVec%x)/vHat%x

                         if (abs(dSx)<1.e-10) then
                            rVec%x=(grid(gP)%xAxis(xP+1)+grid(gP)%xAxis(xP))/2.
                            xP = xP+1
                         end if
                      else
                         dSx = ( grid(gP)%xAxis(grid(gP)%nx)-rVec%x)/vHat%x
                         if (abs(dSx)<1.e-10) then
                            rVec%x=grid(gP)%xAxis(grid(gP)%nx)
                            if (.not.lgPlaneIonization .and. gP==1) return
                         end if
                      end if
                   else if (vHat%x<0.) then
                      if (xP>1) then
                         dSx = ( (grid(gP)%xAxis(xP)+grid(gP)%xAxis(xP-1))/2.-rVec%x)/vHat%x
                         if (abs(dSx)<1.e-10) then
                            rVec%x=(grid(gP)%xAxis(xP)+grid(gP)%xAxis(xP-1))/2.
                            xP = xP-1
                         end if
                      else
                         dSx = (grid(gP)%xAxis(1)-rVec%x)/vHat%x
                         if (abs(dSx)<1.e-10) then
                            rVec%x=grid(gP)%xAxis(1)
                         end if
                      end if
                   else if (vHat%x==0.) then
                      dSx = grid(gP)%xAxis(grid(gP)%nx)
                   end if

                   if (.not.lg1D) then
                      if (vHat%y>0.) then
                         if (yP<grid(gP)%ny) then
                            dSy = ( (grid(gP)%yAxis(yP+1)+grid(gP)%yAxis(yP))/2.-rVec%y)/vHat%y
                            if (abs(dSy)<1.e-10) then
                               rVec%y=(grid(gP)%yAxis(yP+1)+grid(gP)%yAxis(yP))/2.
                               yP = yP+1
                            end if
                         else
                            dSy = (  grid(gP)%yAxis(grid(gP)%ny)-rVec%y)/vHat%y
                            if (abs(dSy)<1.e-10) then
                               rVec%y=grid(gP)%yAxis(grid(gP)%ny)
                               if(gP==1) return
                            end if
                         end if
                      else if (vHat%y<0.) then
                         if (yP>1) then
                            dSy = ( (grid(gP)%yAxis(yP)+grid(gP)%yAxis(yP-1))/2.-rVec%y)/vHat%y
                            if (abs(dSy)<1.e-10) then
                               rVec%y=(grid(gP)%yAxis(yP)+grid(gP)%yAxis(yP-1))/2.
                               yP = yP-1
                            end if
                         else
                            dSy = ( grid(gP)%yAxis(1)-rVec%y)/vHat%y
                            if (abs(dSy)<1.e-10) then
                               rVec%y=grid(gP)%yAxis(1)
                            end if
                         end if
                      else if (vHat%y==0.) then
                         dSy = grid(gP)%yAxis(grid(gP)%ny)
                      end if

                      if (vHat%z>0.) then
                         if (zP<grid(gP)%nz) then
                            dSz = ( (grid(gP)%zAxis(zP+1)+grid(gP)%zAxis(zP))/2.-rVec%z)/vHat%z
                            if (abs(dSz)<1.e-10) then
                               rVec%z=(grid(gP)%zAxis(zP+1)+grid(gP)%zAxis(zP))/2.
                               zP = zP+1
                            end if
                         else
                            dSz = ( grid(gP)%zAxis(grid(gP)%nz)-rVec%z)/vHat%z
                            if (abs(dSz)<1.e-10) then
                               rVec%z=grid(gP)%zAxis(grid(gP)%nz)
                               if (.not.lgPlaneIonization .and. gP==1) return
                            end if
                         end if
                      else if (vHat%z<0.) then
                         if (zP>1) then
                            dSz = ( (grid(gP)%zAxis(zP)+grid(gP)%zAxis(zP-1))/2.-rVec%z)/vHat%z
                            if (abs(dSz)<1.e-10) then
                               rVec%z=(grid(gP)%zAxis(zP)+grid(gP)%zAxis(zP-1))/2.
                               zP = zP-1
                            end if
                         else
                            dSz = ( grid(gP)%zAxis(1)-rVec%z)/vHat%z
                            if (abs(dSz)<1.e-10) then
                               rVec%z=grid(gP)%zAxis(1)
                            end if
                         end if
                      else if (vHat%z==0.) then
                         dSz = grid(gP)%zAxis(grid(gP)%nz)
                      end if

                      if (xP > grid(gP)%nx .or. xP < 1 .or. &
                           & yP > grid(gP)%ny .or. yP < 1 .or. &
                           & zP > grid(gP)%nz .or. zP < 1 ) then
                         print*, "! fluoPathSegment: insanity -2- [gp,xp,yp,zp]", &
                              & gp, xp, yp, zp
                         stop
                      end if

                   end if

                   if (grid(gP)%active(xP,yP,zP)>=0) exit
                end do

                ! cater for cells on cell wall
                if ( abs(dSx)<1.e-10 ) dSx = grid(gP)%xAxis(grid(gP)%nx)
                if ( abs(dSy)<1.e-10 ) dSy = grid(gP)%yAxis(grid(gP)%ny)
                if ( abs(dSz)<1.e-10 ) dSz = grid(gP)%zAxis(grid(gP)%nz)

                ! find the nearest wall
                dSx = abs(dSx)
                dSy = abs(dSy)
                dSz = abs(dSz)

                if (dSx<=0.) then
                   print*, '! fluoPathSegment: [warning] dSx <= 0.',dSx
                   print*, 'grid(gP)%xAxis ', grid(gP)%xAxis
                   print*, 'gP,xP,grid(gP)%xAxis(xP), rVec%x, vHat%x'
                   print*, gP,xP,grid(gP)%xAxis(xP), rVec%x, vHat%x
                   dS = min(dSy, dSz)
                else if (dSy<=0.) then
                   print*, '! fluoPathSegment: [warning] dSy <= 0.', dSy
                   print*, 'grid(gP)%yAxis ', grid(gP)%yAxis
                   print*, 'gP,yP,grid(gP)%yAxis(yP), rVec%y, vHat%y'
                   print*, gP,yP,grid(gP)%yAxis(yP), rVec%y, vHat%y
                   dS = min(dSx, dSz)
                else if (dSz<=0.) then
                   print*, '! fluoPathSegment: [warning] dSz <= 0.', dSz
                   print*, 'grid(gP)%zAxis ', grid(gP)%zAxis
                   print*, 'gP,zP,grid(gP)%zAxis(zP), rVec%z, vHat%z'
                   print*, gP,zP,grid(gP)%zAxis(zP), rVec%z, vHat%z
                   dS = min(dSx, dSy)
                else
                   dS = min(dSx,dSy)
                   dS = min(dS, dSz)
                end if

                ! this should now never ever happen
                if (dS <= 0.) then
                   print*, 'fluoPathSegment: dS <= 0', dSx, dSy, dSz
                   print*, gP, rVec
                   stop
                end if

                ! calculate the optical depth to the next cell wall
                tauCell = dS*grid(gP)%opacity(grid(gP)%active(xP,yP,zP), enPacket%nuP)

                if (lg1D) then
                   if (nGrids>1) then
                      print*, '! getVolumeLoc: 1D option and multiple grids options are not compatible'
                      stop
                   end if

                   if (xP == 1) then

                      dV = 4.*Pi* ( (grid(gP)%xAxis(xP+1)/1.e15)**3)/3.


                   else if ( xP==grid(gP)%nx) then

                      dV = Pi* ( (3.*(grid(gP)%xAxis(xP)/1.e15)-(grid(gP)%xAxis(xP-1)/1.e15))**3 - &
                           & ((grid(gP)%xAxis(xP)/1.e15)+(grid(gP)%xAxis(xP-1)/1.e15))**3 ) / 6.

                   else

                      dV = Pi* ( ((grid(gP)%xAxis(xP+1)/1.e15)+(grid(gP)%xAxis(xP)/1.e15))**3 - &
                           & ((grid(gP)%xAxis(xP-1)/1.e15)+(grid(gP)%xAxis(xP)/1.e15))**3 ) / 6.

                   end if

                   dV = dV/8.

                else

                   if ( (xP>1) .and. (xP<grid(gP)%nx) ) then

                      dx = abs(grid(gP)%xAxis(xP+1)-grid(gP)%xAxis(xP-1))/2.
                   else if ( xP==1 ) then
                      if (lgSymmetricXYZ) then
                         dx = abs(grid(gP)%xAxis(xP+1)-grid(gP)%xAxis(xP))/2.
                      else
                         dx = abs(grid(gP)%xAxis(xP+1)-grid(gP)%xAxis(xP))
                      end if
                   else if ( xP==grid(gP)%nx ) then
                      dx = abs(grid(gP)%xAxis(xP)  -grid(gP)%xAxis(xP-1))
                   end if

                   if ( (yP>1) .and. (yP<grid(gP)%ny) ) then
                      dy = abs(grid(gP)%yAxis(yP+1)-grid(gP)%yAxis(yP-1))/2.
                   else if ( yP==1 ) then
                      if (lgSymmetricXYZ) then
                         dy = abs(grid(gP)%yAxis(yP+1)-grid(gP)%yAxis(yP))/2.
                      else
                         dy = abs(grid(gP)%yAxis(yP+1)-grid(gP)%yAxis(yP))
                      end if
                   else if ( yP==grid(gP)%ny ) then
                      dy = abs(grid(gP)%yAxis(yP)  -grid(gP)%yAxis(yP-1))
                   end if

                   if ( (zP>1) .and. (zP<grid(gP)%nz) ) then
                      dz = abs(grid(gP)%zAxis(zP+1)-grid(gP)%zAxis(zP-1))/2.
                   else if ( zP==1 ) then
                      if (lgSymmetricXYZ) then
                         dz = abs(grid(gP)%zAxis(zP+1)-grid(gP)%zAxis(zP))/2.
                      else
                         dz = abs(grid(gP)%zAxis(zP+1)-grid(gP)%zAxis(zP))
                      end if
                   else if ( zP==grid(gP)%nz ) then
                      dz = abs(grid(gP)%zAxis(zP)-grid(gP)%zAxis(zP-1))
                   end if

                   dx = dx/1.e15
                   dy = dy/1.e15
                   dz = dz/1.e15


                   ! calculate the volume
                   dV = dx*dy*dz

                end if

                ! check if the packet interacts within this cell
                if ((absTau+tauCell > passProb) .and. (grid(gP)%active(xP,yP,zP)>0)) then

                   ! packet interacts

                   ! calculate where within this cell the packet interacts
                   dlLoc = (passProb-absTau)/grid(gP)%opacity(grid(gP)%active(xP,yP,zP), enPacket%nuP)

                   ! update packet's position
                   rVec = rVec + dlLoc*vHat

                   if (lgSymmetricXYZ .and. gP==1) then
                      if ( rVec%x <= grid(gP)%xAxis(1) ) then
                         if (vHat%x<0.) vHat%x = -vHat%x
                         rVec%x = grid(gP)%xAxis(1)
                      end if
                      if ( rVec%y <= grid(gP)%yAxis(1) ) then
                         if (vHat%y<0.) vHat%y = -vHat%y
                         rVec%y = grid(gP)%yAxis(1)
                      end if
                      if ( rVec%z <= grid(gP)%zAxis(1) ) then
                         if (vHat%z<0.) vHat%z = -vHat%z
                         rVec%z = grid(gP)%zAxis(1)
                      end if
                   end if

                   ! add contribution of the packet to the radiation field
                   if (enPacket%lgStellar) then
                      grid(gP)%Jste(grid(gP)%active(xP,yP,zP),enPacket%nuP) = &
                           grid(gP)%Jste(grid(gP)%active(xP,yP,zP),enPacket%nuP) + dlLoc*deltaEUsed / dV
                   else ! if the energy packet is diffuse
                      if (lgDebug) then
                         grid(gP)%Jdif(grid(gP)%active(xP,yP,zP),enPacket%nuP) = &
                              & grid(gP)%Jdif(grid(gP)%active(xP,yP,zP),enPacket%nuP) + dlLoc*deltaEUsed / dV
                      else
                         grid(gP)%Jste(grid(gP)%active(xP,yP,zP),enPacket%nuP) = &
                              & grid(gP)%Jste(grid(gP)%active(xP,yP,zP),enPacket%nuP) + dlLoc*deltaEUsed / dV
                      end if
                   end if

                   ! check if the position within the cell is still within the outer radius
                   if ( sqrt( (rvec%x/1.e10)**2 + (rvec%y/1.e10)**2 + (rvec%z/1.e10)**2)*1.e10 >= R_out &
                        & .and. R_out > 0.) then

                      ! the packet escapes without further interaction
                      idirT = int(acos(enPacket%direction%z)/dTheta)+1
                      if (idirT>totangleBinsTheta) then
                         idirT=totangleBinsTheta
                      end if
                      if (idirT<1 .or. idirT>totAngleBinsTheta) then
                         print*, '! fluorescencePacketRun: error in theta direction cosine assignment',&
                              &  idirT, enPacket, dTheta, totAngleBinsTheta
                         stop
                      end if

                      if (enPacket%direction%x<1.e-35) then
                         idirP = 0
                      else
                         idirP = int(atan(enPacket%direction%y/enPacket%direction%x)/dPhi)
                      end if
                      if (idirP<0) idirP=totAngleBinsPhi+idirP
                      idirP=idirP+1
                      if (idirP>totangleBinsPhi) then
                         idirP=totangleBinsPhi
                      end if

                      if (idirP<1 .or. idirP>totAngleBinsPhi) then
                         print*, '! fluorescencePacketRun: error in Phi direction cosine assignment',&
                              &  idirP, enPacket, dPhi, totAngleBinsPhi
                         stop
                      end if

                      if (nAngleBins>0) then
                         if (viewPointPtheta(idirT) == viewPointPphi(idirP).or. &
                              & (viewPointTheta(viewPointPphi(idirP))==viewPointTheta(viewPointPtheta(idirT))) .or. &
                              & (viewPointPhi(viewPointPtheta(idirT))==viewPointPhi(viewPointPphi(idirP))) ) then
                            grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), enPacket%nuP,viewPointPtheta(idirT)) = &
                                 &grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                                 & enPacket%nuP,viewPointPtheta(idirT)) +  deltaEUsed
                            if (viewPointPtheta(idirT)/=0) grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                                 & enPacket%nuP,0) = &
                                 & grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                                 & enPacket%nuP,0) +  deltaEUsed
                         else
                            grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                                 & enPacket%nuP,0) = &
                                 & grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                                 & enPacket%nuP,0) +  deltaEUsed
                         end if
                      else

                         grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                              & enPacket%nuP,0) = &
                              & grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                              & enPacket%nuP,0) +  deltaEUsed

                      end if

                      return
                   end if


                   ! check if the packet is absorbed or scattered
                   ! NOTE : we are ignoring dust
                   probSca = KNsigmaT(enPacket%nuP)*&
                        &grid(gP)%Ne(grid(gP)%active(xp,yp,zp))/&
                        & (grid(gp)%opacity(grid(gp)%active(xp,yp,zp),&
                        &enPacket%nuP))

                   call random_number(random)

                   random = 1.-random

                   if (random > probSca) then
                      lgScattered = .false.
                   else if (random <= probSca) then
                      lgScattered = .true.
                   else
                      print*, '! fluoPathSegment: insanity occurred and scattering/absorption &
                           & decision stage.'
                      stop
                   end if

                   if (.not. lgScattered) then

                      absInt = absInt + 1.
                      exit

                   else

                      scaInt = scaInt + 1.
                      ! packet is compton scattered
                      ! calculate new direction
                      ! for now assume scattering is isotropic, when phase
                      ! function is introduced the following must be changed

                      enPacket%xP(igpp) = xP
                      enPacket%yP(igpp) = yP
                      enPacket%zP(igpp) = zP

                      if (grid(gP)%active(enPacket%xp(igpp), enPacket%yp(igpp), enPacket%zp(igpp)) < 0.) then
                         print*, "! fluoPathSegment: new packet cannot be emitted from re-mapped cell -1-"
                         print*, "nuP, .false., .true., xp,yp,zp, gp"
                         print*, nuP, .false., .true.,  xp,yp,zp, gp
                         stop
                      end if

                      call comptonScatter(enPacket)

                      vHat%x = enPacket%direction%x
                      vHat%y = enPacket%direction%y
                      vHat%z = enPacket%direction%z

                      ! initialize optical depth
                      absTau = 0.

                      ! get a random number
                      call random_number(random)

                      ! calculate the probability
                      passProb = -log(1.-random)

                   end if

                else

                   ! the packet is not absorbed within this cell
                   ! add contribution of the packet to the radiation field

                   if (enPacket%lgStellar) then
                      grid(gP)%Jste(grid(gP)%active(xP,yP,zP),enPacket%nuP) = &
                           grid(gP)%Jste(grid(gP)%active(xP,yP,zP),enPacket%nuP) + dS*deltaEUsed / dV
                   else ! if the energy packet is diffuse
                      if (lgDebug) then
                         grid(gP)%Jdif(grid(gP)%active(xP,yP,zP),enPacket%nuP) = &
                              & grid(gP)%Jdif(grid(gP)%active(xP,yP,zP),enPacket%nuP) + dS*deltaEUsed / dV
                      else
                         grid(gP)%Jste(grid(gP)%active(xP,yP,zP),enPacket%nuP) = &
                              & grid(gP)%Jste(grid(gP)%active(xP,yP,zP),enPacket%nuP) + dS*deltaEUsed / dV
                      end if
                   end if

                   ! update absTau
                   absTau = absTau+tauCell

                   ! update packet's position
                   rVec = rVec+dS*vHat

                   ! keep track of where you are on mother grid
                   if (gP>1) then

                      if (enPacket%xP(1) <= 0 .or. &
                           & enPacket%yP(1) <= 0 .or. &
                           & enPacket%zP(1) <= 0 ) then

                         ! locate where we are at on the mother grid
                         call locate(grid(grid(gP)%motherP)%xAxis, rvec%x, enPacket%xP(1))
                         if ( enPacket%xP(1)< grid(grid(gP)%motherP)%nx) then
                            if ( rvec%x > ( grid(grid(gP)%motherP)%xAxis(enPacket%xP(1)) + &
                                 & grid(grid(gP)%motherP)%xAxis(enPacket%xP(1)+1)/2.) ) &
                                 & enPacket%xP(1) = enPacket%xP(1)+1
                         end if

                         call locate( grid(grid(gP)%motherP)%yAxis, rvec%y, enPacket%yP(1))
                         if ( enPacket%yP(1)< grid(grid(gP)%motherP)%ny) then
                            if ( rvec%y > ( grid(grid(gP)%motherP)%yAxis(enPacket%yP(1)) + &
                                 & grid(grid(gP)%motherP)%yAxis(enPacket%yP(1)+1)/2.) ) &
                                 & enPacket%yP(1) = enPacket%yP(1)+1
                         end if

                         call locate(grid(grid(gP)%motherP)%zAxis, rvec%z, enPacket%zP(1))
                         if ( enPacket%zP(1)< grid(grid(gP)%motherP)%nz) then
                            if ( rvec%z > ( grid(grid(gP)%motherP)%zAxis(enPacket%zP(1)) + &
                                 & grid(grid(gP)%motherP)%zAxis(enPacket%zP(1)+1)/2.) ) &
                                 & enPacket%zP(1) = enPacket%zP(1)+1
                         end if


                      else

                         if (vHat%x>0.) then
                            if ( enPacket%xP(1) < grid(grid(gP)%motherP)%nx ) then
                               if ( rVec%x > (grid(grid(gP)%motherP)%xAxis(enPacket%xP(1))+&
                                    & grid(grid(gP)%motherP)%xAxis(enPacket%xP(1)+1))/2. ) then
                                  enPacket%xP(1) = enPacket%xP(1)+1
                               end if
                            else
                               if ( rVec%x > grid(grid(gP)%motherP)%xAxis(enPacket%xP(1))) then
!                            print*, '! fluoPathSegment: insanity occurred at mother grid transfer (x axis +)', &
!                                 & rVec%x, gP, grid(gP)%motherP
!                            stop
                               end if
                            end if
                         else
                            if ( enPacket%xP(1) > 1 ) then
                               if ( rVec%x <= (grid(grid(gP)%motherP)%xAxis(enPacket%xP(1)-1)+&
                                    & grid(grid(gP)%motherP)%xAxis(enPacket%xP(1)))/2. ) then
                                  enPacket%xP(1) = enPacket%xP(1)-1
                               end if
                            else
                               if (rVec%x < grid(grid(gP)%motherP)%xAxis(1)) then
!                            print*, '! fluoPathSegment: insanity occurred at mother grid transfer (x axis-)',&
!                                 & rVec%x, gP, grid(gP)%motherP
!                            stop
                               end if
                            end if
                         end if
                         if (vHat%y>0.) then
                            if (  enPacket%yP(1) < grid(grid(gP)%motherP)%ny ) then
                               if ( rVec%y > (grid(grid(gP)%motherP)%yAxis( enPacket%yP(1))+&
                                    & grid(grid(gP)%motherP)%yAxis(enPacket%yP(1)+1))/2. ) then
                                  enPacket%yP(1) =  enPacket%yP(1)+1
                               end if
                            else
                               if ( rVec%y > grid(grid(gP)%motherP)%yAxis( enPacket%yP(1))) then
!                            print*, '! fluoPathSegment: insanity occurred at mother grid transfer (y axis +)',&
!                                 & rVec%y, gP, grid(gP)%motherP
!                            stop
                               end if
                            end if
                         else
                            if (  enPacket%yP(1) > 1 ) then
                               if ( rVec%y <= (grid(grid(gP)%motherP)%yAxis( enPacket%yP(1)-1)+&
                                    & grid(grid(gP)%motherP)%yAxis(enPacket%yP(1)))/2. ) then
                                  enPacket%yP(1) =  enPacket%yP(1)-1
                               end if
                            else
                               if (rVec%y < grid(grid(gP)%motherP)%yAxis(1)) then
!                            print*, '! fluoPathSegment: insanity occurred at mother grid transfer (y axis -)', &
!                                 & rVec%y, gP, grid(gP)%motherP
!                            stop
                               end if
                            end if
                         end if
                         if (vHat%z>0.) then
                            if (  enPacket%zP(1) < grid(grid(gP)%motherP)%nz ) then
                               if ( rVec%z > (grid(grid(gP)%motherP)%zAxis( enPacket%zP(1))+&
                                    & grid(grid(gP)%motherP)%zAxis(enPacket%zP(1)+1))/2. ) then
                                  enPacket%zP(1) =  enPacket%zP(1)+1
                               end if
                            else
                               if ( rVec%z > grid(grid(gP)%motherP)%zAxis( enPacket%zP(1))) then
!                            print*, '! fluoPathSegment: insanity occurred at mother grid transfer (z axis +)', &
!                                 & rVec%z, gP, grid(gP)%motherP
!                            stop
                               end if
                            end if
                         else
                            if (  enPacket%zP(1) > 1 ) then
                               if ( rVec%z <= (grid(grid(gP)%motherP)%zAxis( enPacket%zP(1)-1)+&
                                    & grid(grid(gP)%motherP)%zAxis(enPacket%zP(1)))/2. ) then
                                  enPacket%zP(1) =  enPacket%zP(1)-1
                               end if
                            else
                               if (rVec%z < grid(grid(gP)%motherP)%zAxis(1)) then
!                        print*, '! fluoPathSegment: insanity occurred at mother grid transfer (z axis -)', &
!                             & rVec%z, gP, grid(gP)%motherP
!                        stop
                               end if
                            end if
                         end if

                      end if
                   end if

                   if (.not.lg1D) then
                      if ( (dS == dSx) .and. (vHat%x > 0.)  ) then
                         xP = xP+1
                      else if ( (dS == dSx) .and. (vHat%x < 0.) ) then
                         xP = xP-1
                      else if ( (dS == dSy) .and. (vHat%y > 0.) ) then
                         yP = yP+1
                      else if ( (dS == dSy) .and. (vHat%y < 0.) ) then
                         yP = yP-1
                      else if ( (dS == dSz) .and. (vHat%z > 0.) ) then
                         zP = zP+1
                      else if ( (dS == dSz) .and. (vHat%z < 0.) ) then
                         zP = zP-1
                      else
                         print*, '! fluoPathSegment: insanity occurred in dS assignment &
                              & [dS,dSx,dSy,dSz,vHat]', dS,dSx,dSy,dSz,vHat
                      end if
                   else
                      radius = 1.e10*sqrt((rVec%x/1.e10)*(rVec%x/1.e10) + &
                           & (rVec%y/1.e10)*(rVec%y/1.e10) + &
                           & (rVec%z/1.e10)*(rVec%z/1.e10))
                      call locate(grid(gP)%xAxis, radius , xP)

                   end if

                   ! be 6/6/06
                   if(.not.lgPlaneIonization.and..not.lgSymmetricXYZ) then
                      lgReturn=.false.

                      if ( rVec%y <= grid(gP)%yAxis(1)-grid(gP)%geoCorrY .or. yP<1) then

                         ! the energy packet escapes this grid
                         if (gP==1) then
                            yP=1
                            lgReturn=.true.
                         else if (gP>1) then
                            xP = enPacket%xP(grid(gP)%motherP)
                            yP = enPacket%yP(grid(gP)%motherP)
                            zP = enPacket%zP(grid(gP)%motherP)
                            gP = grid(gP)%motherP
                         else
                            print*, '! fluoPathSegment: insanity occurred - invalid gP', gP
                            stop
                         end if

                      end if

                      if (rVec%y > grid(gP)%yAxis(grid(gP)%ny)+grid(gP)%geoCorrY .or. yP>grid(gP)%ny) then

                         if (gP==1) then
                            ! the energy packet escapes
                            yP = grid(gP)%ny
                            lgReturn=.true.
                         else if (gP>1) then
                            xP = enPacket%xP(grid(gP)%motherP)
                            yP =  enPacket%yP(grid(gP)%motherP)
                            zP =  enPacket%zP(grid(gP)%motherP)
                            gP = grid(gP)%motherP
                         else
                            print*, '! fluoPathSegment: insanity occurred - invalid gP', gP
                            stop
                         end if

                      end if

                      if ( (rVec%x <= grid(gP)%xAxis(1)-grid(gP)%geoCorrX .or. xP<1) .and. gP==1) then
                         xP=1
                         lgReturn=.true.
                      end if


                      if ( (rVec%x <= grid(gP)%xAxis(1)-grid(gP)%geoCorrX .or. xP<1) &
                           & .and. gP>1) then

                         xP = enPacket%xP(grid(gP)%motherP)
                         yP =  enPacket%yP(grid(gP)%motherP)
                         zP =  enPacket%zP(grid(gP)%motherP)
                         gP = grid(gP)%motherP

                      end if


                      if ( (rVec%x >=  grid(gP)%xAxis(grid(gP)%nx)+grid(gP)%geoCorrX &
                           & .or. xP>grid(gP)%nx) .and. gP==1 )then
                         xP = grid(gP)%nx
                         lgReturn=.true.
                      end if

                      if ( (rVec%x >=  grid(gP)%xAxis(grid(gP)%nx)+grid(gP)%geoCorrX&
                           & .or. xP>grid(gP)%nx) .and.  gP>1) then

                         xP = enPacket%xP(grid(gP)%motherP)
                         yP =  enPacket%yP(grid(gP)%motherP)
                         zP =  enPacket%zP(grid(gP)%motherP)
                         gP = grid(gP)%motherP

                      end if

                      if ( (rVec%z <= grid(gP)%zAxis(1)-grid(gP)%geoCorrZ .or.zP<1) &
                           & .and. gP==1) then
                         zP=1
                         lgReturn=.true.

                     end if

                     if ( (rVec%z <= grid(gP)%zAxis(1)-grid(gP)%geoCorrZ &
                          & .or.zP<1) .and. gP>1) then

                        xP = enPacket%xP(grid(gP)%motherP)
                        yP =  enPacket%yP(grid(gP)%motherP)
                        zP =  enPacket%zP(grid(gP)%motherP)
                        gP = grid(gP)%motherP

                     end if

                     if ( (rVec%z >=  grid(gP)%zAxis(grid(gP)%nz)+grid(gP)%geoCorrZ &
                          & .or. zP>grid(gP)%nz) &
                          & .and. gP==1) then

                        zP = grid(gP)%nz
                        lgReturn=.true.

                     end if

                     if ((rVec%z >=  grid(gP)%zAxis(grid(gP)%nz)+grid(gP)%geoCorrZ &
                          & .or. zP>grid(gP)%nz) .and. gP>1) then

                        xP = enPacket%xP(grid(gP)%motherP)
                        yP =  enPacket%yP(grid(gP)%motherP)
                        zP =  enPacket%zP(grid(gP)%motherP)
                        gP = grid(gP)%motherP

                     end if

                     if (lgReturn) then

                        ! the packet escapes without further interaction
                        idirT = int(acos(enPacket%direction%z)/dTheta)+1
                        if (idirT>totangleBinsTheta) then
                           idirT=totangleBinsTheta
                        end if
                        if (idirT<1 .or. idirT>totAngleBinsTheta) then
                           print*, '! fluorescencePacketRun: error in theta direction cosine assignment',&
                                &  idirT, enPacket, dTheta, totAngleBinsTheta
                           stop
                        end if


                        if (enPacket%direction%x<1.e-35) then
                           idirP = 0
                        else
                           idirP = int(atan(enPacket%direction%y/enPacket%direction%x)/dPhi)
                        end if
                        if (idirP<0) idirP=totAngleBinsPhi+idirP
                        idirP=idirP+1

                        if (idirP>totangleBinsPhi) then
                           idirP=totangleBinsPhi
                        end if
                        if (idirP<1 .or. idirP>totAngleBinsPhi) then
                           print*, '! fluorescencePacketRun: error in phi direction cosine assignment',&
                                &  idirP, enPacket, dPhi, totAngleBinsPhi
                           stop
                        end if


                        if (nAngleBins>0) then
                           if (viewPointPtheta(idirT) == viewPointPphi(idirP).or. &
                                & (viewPointTheta(viewPointPphi(idirP))==viewPointTheta(viewPointPtheta(idirT))) .or. &
                                & (viewPointPhi(viewPointPtheta(idirT))==viewPointPhi(viewPointPphi(idirP))) ) then
                              grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), enPacket%nuP,&
                                   & viewPointPtheta(idirT)) =grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2),&
                                   & enPacket%nuP,viewPointPtheta(idirT)) +  deltaEUsed

                              if (viewPointPtheta(idirT)/=0) grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                                   &enPacket%nuP,0) = &
                                   & grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                                   & enPacket%nuP,0) +  deltaEUsed

                           else
                              grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                                   & enPacket%nuP,0) = &
                                   & grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                                   & enPacket%nuP,0) +  deltaEUsed

                           end if

                        else

                           grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                                enPacket%nuP,0) = &
                                & grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                                & enPacket%nuP,0) +  deltaEUsed
                        end if

                        return
                     end if

                  end if

                  ! end be 6/6/06


                  if(lgPlaneIonization) then
                     lgReturn=.false.

                     if ( rVec%y <= grid(gP)%yAxis(1)-grid(gP)%geoCorrY .or. yP<1) then

                        ! the energy packet escapes this grid
                        if (gP==1) then
                           yP=1
                           lgReturn=.true.
                        else if (gP>1) then
                           xP = enPacket%xP(1)
                           yP = enPacket%yP(1)
                           zP = enPacket%zP(1)
                           gP = 1
                           igpp = 1
                        else
                           print*, '! fluoPathSegment: insanity occurred - invalid gP', gP
                           stop
                        end if

                     end if

                     if (rVec%y > grid(gP)%yAxis(grid(gP)%ny)+grid(gP)%geoCorrY .or. yP>grid(gP)%ny) then

                        if (gP==1) then
                           ! the energy packet escapes
                           yP = grid(gP)%ny
                           lgReturn=.true.
                        else if (gP>1) then
                           xP = enPacket%xP(1)
                           yP = enPacket%yP(1)
                           zP = enPacket%zP(1)
                           gP = 1
                           igpp = 1
                        else
                           print*, '! fluoPathSegment: insanity occurred - invalid gP', gP
                           stop
                        end if

                     end if

                     if ( (rVec%x <= grid(1)%xAxis(1) .or. xP<1) ) then
                        xP=1
                        rVec%x = grid(gP)%xAxis(1)
                        vHat%x = -vHat%x

                     end if

                     if ( (rVec%x <= grid(gP)%xAxis(1)-grid(gP)%geoCorrX .or. xP<1) &
                          & .and. gP>1) then

                        xP = enPacket%xP(1)
                        yP = enPacket%yP(1)
                        zP = enPacket%zP(1)
                        gP = 1
                        igpp = 1

                     end if

                     if ( (rVec%x >=  grid(1)%xAxis(grid(gP)%nx) &
                          & .or. xP>grid(gP)%nx)  )then
                        xP = grid(gP)%nx
                        rVec%x = grid(gP)%xAxis(grid(gP)%nx)
                        vHat%x = -vHat%x

                     end if

                     if ( (rVec%x >=  grid(gP)%xAxis(grid(gP)%nx)+grid(gP)%geoCorrX&
                          & .or. xP>grid(gP)%nx) .and.  gP>1) then

                        xP = enPacket%xP(1)
                        yP = enPacket%yP(1)
                        zP = enPacket%zP(1)
                        gP = 1
                        igpp = 1
                     end if

                     if ( (rVec%z <= grid(1)%zAxis(1) .or.zP<1) ) then
                        zP=1
                        rVec%z = grid(gP)%yAxis(1)
                        vHat%z = -vHat%z

                     end if

                     if ( (rVec%z <= grid(gP)%zAxis(1)-grid(gP)%geoCorrZ &
                          & .or.zP<1) .and. gP>1) then

                        xP = enPacket%xP(1)
                        yP = enPacket%yP(1)
                        zP = enPacket%zP(1)
                        gP = 1
                        igpp = 1

                     end if

                     if ( (rVec%z >=  grid(1)%zAxis(grid(gP)%nz) .or. zP>grid(gP)%nz) &
                          & ) then

                        zP = grid(gP)%nz
                        rVec%z = grid(gP)%zAxis(grid(gP)%nz)
                        vHat%z = -vHat%z

                     end if

                     if ((rVec%z >=  grid(gP)%zAxis(grid(gP)%nz)+grid(gP)%geoCorrZ &
                          & .or. zP>grid(gP)%nz) .and. gP>1) then

                        xP = enPacket%xP(1)
                        yP = enPacket%yP(1)
                        zP = enPacket%zP(1)
                        gP = 1
                        igpp = 1

                     end if


                     if (lgReturn) then

                        ! the packet escapes without further interaction

                        idirT = int(acos(enPacket%direction%z)/dTheta)+1
                        if (idirT>totangleBinsTheta) then
                           idirT=totangleBinsTheta
                        end if
                        if (idirT<1 .or. idirT>totAngleBinsTheta) then
                           print*, '! energyPacketRun: error in theta direction cosine assignment',&
                                &  idirT, enPacket, dTheta, totAngleBinsTheta
                           stop
                        end if


                        if (enPacket%direction%x<1.e-35) then
                           idirP = 0
                        else
                           idirP = int(atan(enPacket%direction%y/enPacket%direction%x)/dPhi)
                        end if
                        if (idirP<0) idirP=totAngleBinsPhi+idirP
                        idirP=idirP+1

                        if (idirP>totangleBinsPhi) then
                           idirP=totangleBinsPhi
                        end if
                        if (idirP<1 .or. idirP>totAngleBinsPhi) then
                           print*, '! energyPacketRun: error in phi direction cosine assignment',&
                                &  idirP, enPacket, dPhi, totAngleBinsPhi
                           stop
                        end if


                        if (nAngleBins>0) then
                           if (viewPointPtheta(idirT) == viewPointPphi(idirP).or. &
                                & (viewPointTheta(viewPointPphi(idirP))==viewPointTheta(viewPointPtheta(idirT))) .or. &
                                & (viewPointPhi(viewPointPtheta(idirT))==viewPointPhi(viewPointPphi(idirP))) ) then
                              grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), enPacket%nuP,viewPointPtheta(idirT)) = &
                                   &grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                                   & enPacket%nuP,viewPointPtheta(idirT)) +  deltaEUsed
                              if (viewPointPtheta(idirT)/=0) grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                                   enPacket%nuP,0) = &
                                   & grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                                   & enPacket%nuP,0) +  deltaEUsed

                           else
                              grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                                enPacket%nuP,0) = &
                                & grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                                & enPacket%nuP,0) +  deltaEUsed

                           end if

                        else

                           grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                                enPacket%nuP,0) = &
                                & grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                                & enPacket%nuP,0) +  deltaEUsed

                        end if

                        return
                     end if

                  end if

                  ! check if the path is still within the simulation region
                  radius = 1.e10*sqrt((rVec%x/1.e10)*(rVec%x/1.e10) + &
                       &                     (rVec%y/1.e10)*(rVec%y/1.e10) + &
                       &                     (rVec%z/1.e10)*(rVec%z/1.e10))

                  if (.not.lgPlaneIonization) then

                     if ( (.not.lgSymmetricXYZ .and. (rVec%x<=grid(1)%xAxis(1)-grid(1)%geoCorrX .or.&
                          & rVec%y<=grid(1)%yAxis(1)-grid(1)%geoCorrY .or. rVec%z<=grid(1)%zAxis(1)-&
                          & grid(1)%geoCorrZ)) .or.&
                          & (rVec%x >= grid(gP)%xAxis(grid(gP)%nx)+grid(gP)%geoCorrX) .or.&
                          &(rVec%y >= grid(gP)%yAxis(grid(gP)%ny)+grid(gP)%geoCorrY) .or.&
                          &(rVec%z >= grid(gP)%zAxis(grid(gP)%nz)+grid(gP)%geoCorrZ) .or. &
                          & xP>grid(gP)%nx .or. yP>grid(gP)%ny .or. zP>grid(gP)%nz ) then

                        if (gP==1) then

                           if (enPacket%xP(1) > grid(1)%nx) xP = grid(1)%nx
                           if (enPacket%yP(1) > grid(1)%ny) yP = grid(1)%ny
                           if (enPacket%zP(1) > grid(1)%nz) zP = grid(1)%nz
                           if (enPacket%xP(1) < 1) xP = 1
                           if (enPacket%yP(1) < 1) yP = 1
                           if (enPacket%zP(1) < 1) zP = 1


                           ! the energy packet escapes

                           ! the packet escapes without further interaction
                           idirT = int(acos(enPacket%direction%z)/dTheta)+1
                           if (idirT>totangleBinsTheta) then
                              idirT=totangleBinsTheta
                           end if
                           if (idirT<1 .or. idirT>totAngleBinsTheta) then
                              print*, '! energyPacketRun: error in theta direction cosine assignment',&
                                   &  idirT, enPacket, dTheta, totAngleBinsTheta
                              stop
                           end if

                           if (enPacket%direction%x<1.e-35) then
                              idirP = 0
                           else
                              idirP = int(atan(enPacket%direction%y/enPacket%direction%x)/dPhi)
                           end if
                           if (idirP<0) idirP=totAngleBinsPhi+idirP
                           idirP=idirP+1

                           if (idirP>totangleBinsPhi) then
                              idirP=totangleBinsPhi
                           end if

                           if (idirP<1 .or. idirP>totAngleBinsPhi) then
                              print*, '! energyPacketRun: error in phi direction cosine assignment',&
                                   &  idirP, enPacket, dPhi, totAngleBinsPhi
                              stop
                           end if

                           if (nAngleBins>0) then
                              if (viewPointPtheta(idirT) == viewPointPphi(idirP).or. &
                                   & (viewPointTheta(viewPointPphi(idirP))==viewPointTheta(viewPointPtheta(idirT))) .or. &
                                   & (viewPointPhi(viewPointPtheta(idirT))==viewPointPhi(viewPointPphi(idirP))) ) then
                                 grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), enPacket%nuP,&
                                      & viewPointPtheta(idirT)) = &
                                      &grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                                      & enPacket%nuP,viewPointPtheta(idirT)) +  deltaEUsed
                                 if (viewPointPtheta(idirT)/=0) grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                                      enPacket%nuP,0) = &
                                      & grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                                      & enPacket%nuP,0) +  deltaEUsed
                              else
                                 grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                                      enPacket%nuP,0) = &
                                      & grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                                      & enPacket%nuP,0) +  deltaEUsed
                              end if
                           else

                              grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                                   enPacket%nuP,0) = &
                                   & grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                                   & enPacket%nuP,0) +  deltaEUsed

                           end if

                           !b2.005
                           return

                        else if (gP>1) then

                           xP = enPacket%xP(1)
                           yP = enPacket%yP(1)
                           zP = enPacket%zP(1)
                           gP = 1
                           igpp = 1


                           if (gP/=1) then
                              print*, '! fluoPathSegment: nested multigrids still not implemented'
                              stop
                           end if

                           if ( (radius >= R_out .and. R_out >= 0.) .or.&
                                & (rVec%x >= grid(1)%xAxis(grid(1)%nx)+grid(1)%geoCorrX) .or.&
                                &(rVec%y >= grid(1)%yAxis(grid(1)%ny)+grid(1)%geoCorrY) .or.&
                                &(rVec%z >= grid(1)%zAxis(grid(1)%nz)+grid(1)%geoCorrZ) .or. &
                                & (.not.lgSymmetricXYZ .and.  (rVec%x<=grid(1)%xAxis(1)-grid(1)%geoCorrX .or.&
                                & rVec%y<=grid(1)%yAxis(1)-grid(1)%geoCorrY .or. rVec%z<=grid(1)%zAxis(1)-&
                                & grid(1)%geoCorrZ))) then

                              if (xP > grid(gP)%nx) xP = grid(gP)%nx
                              if (yP > grid(gP)%ny) yP = grid(gP)%ny
                              if (zP > grid(gP)%nz) zP = grid(gP)%nz

                              if (xP < 1) xP = 1
                              if (yP < 1) yP = 1
                              if (zP < 1) zP = 1

                              ! the energy packet escapes

                              ! the packet escapes without further interaction

                              idirT = int(acos(enPacket%direction%z)/dTheta)+1
                              if (idirT>totangleBinsTheta) then
                                 idirT=totangleBinsTheta
                              end if
                              if (idirT<1 .or. idirT>totAngleBinsTheta) then
                                 print*, '! energyPacketRun: error in theta direction cosine assignment',&
                                      &  idirT, enPacket, dTheta, totAngleBinsTheta
                                 stop
                              end if

                              if (enPacket%direction%x<1.e-35) then
                                 idirP = 0
                              else
                                 idirP = int(atan(enPacket%direction%y/enPacket%direction%x)/dPhi)
                              end if
                              if (idirP<0) idirP=totAngleBinsPhi+idirP
                              idirP=idirP+1

                              if (idirP>totangleBinsPhi) then
                                 idirP=totangleBinsPhi
                              end if

                              if (idirP<1 .or. idirP>totAngleBinsPhi) then
                                 print*, '! energyPacketRun: error in phi direction cosine assignment',&
                                      &  idirP, enPacket, dPhi, totAngleBinsPhi
                                 stop
                              end if

                              if (nAngleBins>0) then
                                 if (viewPointPtheta(idirT) == viewPointPphi(idirP).or. &
                                      &(viewPointTheta(viewPointPphi(idirP))==viewPointTheta(viewPointPtheta(idirT))).or. &
                                      & (viewPointPhi(viewPointPtheta(idirT))==viewPointPhi(viewPointPphi(idirP))) ) then
                                    grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), enPacket%nuP,&
                                         & viewPointPtheta(idirT)) = &
                                         &grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                                         & enPacket%nuP,viewPointPtheta(idirT)) +  deltaEUsed
                                    if (viewPointPtheta(idirT)/=0) grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                                         enPacket%nuP,0) = &
                                         & grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                                         & enPacket%nuP,0) +  deltaEUsed
                                 else
                                    grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                                         enPacket%nuP,0) = &
                                         & grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                                         & enPacket%nuP,0) +  deltaEUsed
                                 end if
                              else

                                 grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                                      enPacket%nuP,0) = &
                                      & grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                                      & enPacket%nuP,0) +  deltaEUsed

                              end if

                              !b2.005
                              return

                           end if
                        else
                           print*, '! fluoPathSegment: insanity occurred - invalid gP - ', gP
                           stop
                        end if

                     end if

                     if (lgSymmetricXYZ ) then
                         if (lgPlaneIonization) then
                            print*, '! fluoPathSegment: lgSymmetric and lgPlaneionization flags both raised'
                            stop
                         end if

                         if ( rVec%x <= grid(1)%xAxis(1) .or. (gP==1 .and. xP<1)) then
                            if (vHat%x<0.) vHat%x = -vHat%x
                            enPacket%xP(1) = 1
                            xP = 1
                            rVec%x = grid(gP)%xAxis(1)
                         end if
                         if ( rVec%y <= grid(1)%yAxis(1) .or. (gP==1 .and. yP<1)) then
                            if (vHat%y<0.) vHat%y = -vHat%y
                            enPacket%yP(1)=1
                            yP = 1
                            rVec%y = grid(gP)%yAxis(1)
                         end if
                         if ( rVec%z <= grid(1)%zAxis(1) .or. (gP==1 .and. zP<1)) then
                            if (vHat%z<0.) vHat%z = -vHat%z
                            enPacket%zP(1) = 1
                            zP=1
                            rVec%z = grid(1)%zAxis(1)
                         end if

                      end if

                   end if

                   if (gP>1) then
                      if ( ( (rVec%x <= grid(gP)%xAxis(1) &
                           &.or. xP<1) .and. vHat%x <=0.) .or. &
                           & ( (rVec%y <= grid(gP)%yAxis(1) &
                           & .or. yP<1) .and. vHat%y <=0.) .or. &
                           & ( (rVec%z <= grid(gP)%zAxis(1) &
                           &  .or. zP<1) .and. vHat%z <=0.) .or. &
                           & ( (rVec%x >= grid(gP)%xAxis(grid(gP)%nx) &
                           &.or. xP>grid(gP)%nx) .and. vHat%x >=0.) .or. &
                           & ( (rVec%y >= grid(gP)%yAxis(grid(gP)%ny) &
                           & .or. yP>grid(gP)%ny) .and. vHat%y >=0.) .or. &
                           & ( (rVec%z >= grid(gP)%zAxis(grid(gP)%nz) &
                           &  .or. zP>grid(gP)%nz) .and. vHat%z >=0.) ) then

                         ! go back to mother grid
                         xP = enPacket%xP(1)
                         yP = enPacket%yP(1)
                         zP = enPacket%zP(1)
                         gP = 1
                         igpp = 1

                      end if

                   end if

               end if

               if (.not. lgPlaneIonization .and. gP==1 .and. (xP > grid(gP)%nx  &
                    & .or. yP > grid(gP)%ny .or. zP > grid(gP)%nz) ) then

           ! the energy packet escapes

           ! the packet escapes without further interaction

           idirT = int(acos(enPacket%direction%z)/dTheta)+1
           if (idirT>totangleBinsTheta) then
              idirT=totangleBinsTheta
           end if
           if (idirT<1 .or. idirT>totAngleBinsTheta) then
              print*, '! energyPacketRun: error in theta direction cosine assignment',&
                   &  idirT, enPacket, dTheta, totAngleBinsTheta
              stop
           end if

           if (enPacket%direction%x<1.e-35) then
              idirP = 0
           else
              idirP = int(atan(enPacket%direction%y/enPacket%direction%x)/dPhi)
           end if
           if (idirP<0) idirP=totAngleBinsPhi+idirP
           idirP=idirP+1

           if (idirP>totangleBinsPhi) then
              idirP=totangleBinsPhi
           end if

           if (idirP<1 .or. idirP>totAngleBinsPhi) then
              print*, '! energyPacketRun: error in phi direction cosine assignment',&
                   &  idirP, enPacket, dPhi, totAngleBinsPhi
              stop
           end if

           if (nAngleBins>0) then
              if (viewPointPtheta(idirT) == viewPointPphi(idirP).or. &
                   & (viewPointTheta(viewPointPphi(idirP))==viewPointTheta(viewPointPtheta(idirT))) .or. &
                   & (viewPointPhi(viewPointPtheta(idirT))==viewPointPhi(viewPointPphi(idirP))) ) then
                 grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), enPacket%nuP,&
                      & viewPointPtheta(idirT)) = &
                      &grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                   & enPacket%nuP,viewPointPtheta(idirT)) +  deltaEUsed
                 if (viewPointPtheta(idirT)/=0) grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                      enPacket%nuP,0) = &
                      & grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                      & enPacket%nuP,0) +  deltaEUsed
              else
                 grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                      enPacket%nuP,0) = &
                      & grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                      & enPacket%nuP,0) +  deltaEUsed
              end if
           else

              grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                   enPacket%nuP,0) = &
                   & grid(enPacket%origin(1))%escapedPackets(enPacket%origin(2), &
                   & enPacket%nuP,0) +  deltaEUsed

           end if

           !b2.005
            return


         end if

      end do ! safelimit loop

      if (i>= safeLimit) then
         if (.not.lgPlaneIonization) then
            print*, '! fluoPathSegment: [warning] packet trajectory has exceeded&
                 &  maximum number of events', safeLimit, gP, xP,yP,zP, grid(gP)%active(xP,yP,zP), &
              & rvec, vhat, enPacket
         end if
         return

      end if

      if (gP==1) then
         igpp = 1
      else if (gP>1) then
         igpp = 2
      else
         print*, "! fluoPathSegment: insane grid index "
         stop
      end if

      enPacket%xP(igpp) = xP
      enPacket%yP(igpp) = yP
      enPacket%zP(igpp) = zP

      packetType = 'diffuse'

      ! the energy packet has beenid absorbed - reemit a new packet from this position
      call fluorescencePacketRun(grid, packetType, rvec, enPacket%xP(1:2), enPacket%yP(1:2), &
           & enPacket%zP(1:2), gP)

      return

 end subroutine fluoPathSegment

 ! determine energy and direction for compton redistribution of fluorescent lines
 subroutine comptonScatter(fPacket)
   implicit none


   real                              :: newNu,newTheta     ! new energy and theta
   real                              :: random             ! random number
   real                              :: u,v,w,t            ! direction units
   type(photon_packet),intent(inout) :: fPacket            ! the fluorescent packet

   integer                           :: newNuP,newThetaP   ! pointer to new energy and theta

   ! calculate new direction
   call getNu2(KNsigmaArray(fPacket%nuP,1:180),newThetaP)
   newTheta = real(newThetaP)*Pi/180.

   ! calculate u,v,w (Harries & Howarth, 1997)
   call random_number(random)
   w = 2.*random - 1.
   t = sqrt(1-w*w)
   u = t*cos(newTheta)
   v = t*sin(newTheta)

   fPacket%direction%x = u
   fPacket%direction%y = v
   fPacket%direction%z = w

   ! calculate new energy
   newNu = fPacket%nu*PcompArray(fPacket%nuP,newThetaP)

   ! find this in the nuArray
   call locate(nuArray,newNu,newNuP)

   fPacket%Nu  = newNu
   fPacket%NuP = newNuP

 end subroutine comptonScatter

 ! this subroutine determines the index in the PDFs
 ! same as getNu2
 subroutine getNu2(probDen, nuP)

   real, dimension(:), intent(in) :: probDen    ! probability density function

   real                           :: random     ! random number

   integer, intent(out)           :: nuP        ! frequency index of the new

   integer                        :: isearch,i  !

   ! get a random number
   call random_number(random)

   do i = 1, 10000
      if (random==0 .or. random==1.) then
         call random_number(random)
      else
         exit
      end if
   end do
   if (i>=10000) then
      print*, '! getNu2: problem with random number generator', random, i
      stop
   end if

   ! see what frequency random corresponds to
   nuP=1
   do isearch = 1, nbins
      if (random>=probDen(isearch)) then
         nuP=isearch
      else
         exit
      end if
   end do

   if (nuP>=nbins) then
      print*, 'random: ', random
      print*, 'probDen: ', probDen
   end if

 end subroutine getNu2

end subroutine fluorescencePacketRun


end module fluorescence_mod