Source

enzo-timing-woc / src / enzo / calc_rates.src

The trunk branch has multiple heads

Full commit
   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
#include "fortran.def"
c=======================================================================
c//////////////////////  SUBROUTINE CALC_RATES  \\\\\\\\\\\\\\\\\\\\\\\\\
c
      subroutine calc_rates(
     &             nratec, aye, temstart, temend, alpha0, f3, iradtype,
     &             utem, uxyz, uaye, urho, utim,
     &             ceHIa, ceHeIa, ceHeIIa, ciHIa, ciHeIa, 
     &             ciHeISa, ciHeIIa, reHIIa, reHeII1a, 
     &             reHeII2a, reHeIIIa, brema, compa,
     &             piHI, piHeI, piHeII,
     &             hyd01ka, h2k01a, vibha, rotha, rotla, 
     &             gpldl, gphdl, hdltea, hdlowa,
     &             k1a, k2a, k3a, k4a, k5a, k6a, k7a, k8a, k9a, k10a,
     &             k11a, k12a, k13a, k13dda, k14a, k15a, k16a, k17a,
     &             k18a, k19a, k20a, k21a, k22a,
     &             k24, k25, k26, k27, k28, k29, k30, k31,
     &             k50a, k51a, k52a, k53a, k54a, k55a, k56a, ioutput
     &                     )
c
c  COMPUTE MULTISPECIES RATE LOOKUP TABLE
c
c  written by: Yu Zhang, Peter Anninos and Tom Abel
c  date:       
c  modified1: January, 1996 by Greg Bryan; adapted to KRONOS
c  modified2: October, 1996 by GB; moved to AMR
c  modified3: February, 2000 by GB; added Toms new collisional rates
c
c  PURPOSE:
c    Construct tables for rate coefficients and some cooling functions 
c    versus the logarithm of temperature.
c
c  UNITS:
c    Most rate coefficients are computed in cgs units and then converted
c    into more convenient units.  This would usually be code units
c    (see units.src) but in many cases the natural code units are not
c    near unity (e.g. particles/box volume), so we include units and
c    constants from the equations themselves.  This is explained for
c    each of two types used: rate coefficients, cooling coefficients.
c    Note also that the densities in the rate equations are true, not
c    comoving densities and so we must include a factor of a^3.
c    NOTE: if the spectrum changes, this routine must be called again.
c
c  PARAMETERS:
c    nnu      = number of frequency bins
c    everg    = ergs per eV (Rydberg constant) in cgs
c    evhz     = everg / h
c    mh       = mass of hydrogen atom
c    temstart = table start temperature (K)
c    temend   = table end temperature (K)
c    tevk     = conversion constant [eV] = [K] / tevk
c    ireco    = Recombination cooling flag (1 = on)
c    ibrco    = Brehmmstrahlung cooling flag (1 = on)
c    icico    = Collisional ionization flag (1 = on)
c    iceco    = Collisional excitation flag (1 = on)
c    iphoto_mol = Molecular photo-dissociation flag (1 = on)
c
c  INPUTS
c    nratec   = number of temperature bins in rate equations
c
c    tbase1   = code time units
c    xbase1   = code length units
c    kunit    = conversion factor xbase1**3/tbase1
c
c    alpha0   = power law slope of radiation flux (?)
c    f3       = amplitude of external flux at z = 3
c    iradtype = radiation field type (0 - none,
c                8 - hard spectrum w/ absorption,
c                otherwise ignored)
c
c  RATES:
c
c  (New rates numbering as of Feb/2000, agrees with original
c   Abel etal 1997)
c    old   -> new              old   -> new
c    k1-10 -> k1-10            k16   -> k14
c    k11   -> k13              k17   -> k15
c    k12   -> k11              k18   -> k16
c    k13   -> removed          k19   -> k17
c    k14   -> k12              k20   -> k18
c    k15   -> removed          k21   -> k19
c
C ------- 1:    HI    + e   -> HII   + 2e
C ------- 2:    HII   + e   -> H     + p
C ------- 3:    HeI   + e   -> HeII  + 2e
C ------- 4:    HeII  + e   -> HeI   + p
C ------- 5:    HeII  + e   -> HeIII + 2e
C ------- 6:    HeIII + e   -> HeII  + p
C ------- 7:    HI    + e   -> HM    + p
C ------- 8:    HM    + HI  -> H2I*  + e
C ------- 9:    HI    + HII -> H2II  + p
C ------- 10:   H2II  + HI  -> H2I*  + HII
c
c
C --------------  old  ---------------------
C ------- 11:   H2I   + H   -> 3H
C ------- 12:   H2I   + HII -> H2II  + H
C ------- 13:   H2I   + e   -> HI    + HM
C ------- 14:   H2I   + e   -> 2HI   + e
C ------- 15:   H2I   + H2I -> H2I   + 2HI
C ------- 16:   HM    + e   -> HI    + 2e
C ------- 17:   HM    + HI  -> 2H    + e
C ------- 18:   HM    + HII -> 2HI
C ------- 19:   HM    + HII -> H2II  + e
C ------- 20:   H2II  + e   -> 2HI
C ------- 21:   H2II  + HM  -> HI    + H2I
C --------------  old  ---------------------
c
c --------------  new  ---------------------
C ---11--       H2I   + HII -> H2II  + H
C ---12--       H2I   + e   -> 2HI   + e
C ---13--       H2I   + H   -> 3H
C ---14--       HM    + e   -> HI    + 2e
C ---15--       HM    + HI  -> 2H    + e
C ---16--       HM    + HII -> 2HI
C ---17--       HM    + HII -> H2II  + e
C ---18--       H2II  + e   -> 2HI
C ---19--       H2II  + HM  -> HI    + H2I
c (20-21) - unused
c --------------  new  ---------------------
c
c
C ------- 22:   2H    + H   -> H2I   + H
C ------- 24:   HI    + p   -> HII   + e
C ------- 25:   HeII  + p   -> HeIII + e
C ------- 26:   HeI   + p   -> HeII  + e
C ------- 27:   HM    + p   -> HI    + e
C ------- 28:   H2II  + p   -> HI    + HII
C ------- 29:   H2I   + p   -> H2II  + e
C ------- 30:   H2II  + p   -> 2HII  + e
C ------- 31:   H2I   + p   -> 2HI
C
c ------- 50-56: deuterium rates (given below)
c
c-----------------------------------------------------------------------
c
c  This define indicates if the Tom Abels new (as of Feb/2000) collisional
c    rates (k1-k19) should be used.
c
#define USE_NEW_COLLISIONAL_RATES
c
      implicit NONE
c
c  Arguments
c
      integer nratec, iradtype
      real    aye, temstart, temend, alpha0, f3,
     &        utem, uxyz, uaye, urho, utim
      real    ceHIa(nratec), ceHeIa(nratec), ceHeIIa(nratec),
     &        ciHIa(nratec), ciHeIa(nratec), ciHeISa(nratec), 
     &        ciHeIIa(nratec), reHIIa(nratec), reHeII1a(nratec), 
     &        reHeII2a(nratec), reHeIIIa(nratec), brema(nratec)
      real    hyd01ka(nratec), h2k01a(nratec), vibha(nratec), 
     &        rotha(nratec), rotla(nratec), gpldl(nratec),
     &        gphdl(nratec), hdltea(nratec), hdlowa(nratec)
      real    compa, piHI, piHeI, piHeII
      real    k1a (nratec), k2a (nratec), k3a (nratec), k4a (nratec), 
     &        k5a (nratec), k6a (nratec), k7a (nratec), k8a (nratec), 
     &        k9a (nratec), k10a(nratec), k11a(nratec), k12a(nratec), 
     &        k13a(nratec), k14a(nratec), k15a(nratec), k16a(nratec), 
     &        k17a(nratec), k18a(nratec), k19a(nratec), k20a(nratec), 
     &        k21a(nratec), k24, k25, k26, k27, k28, k29, k30, k31,
     &        k22a(nratec), k50a(nratec), k51a(nratec), k52a(nratec), 
     &        k53a(nratec), k54a(nratec), k55a(nratec), k56a(nratec)
      integer ioutput
      real    k13dda(nratec, 7)
c
c  Parameters
c
      integer nnu
      parameter (nnu = 400)
      integer ireco, ibrco, icico, iceco, iphoto_mol
      parameter (ireco = 1, ibrco = 1, icico = 1, iceco = 1, 
     &           iphoto_mol = 0)

      double precision everg, evhz, tevk, mh, pi, tiny_cgs, dhuge, kb
      parameter (everg = 1.60184d-12, evhz  = 2.41838d+14, 
     &           tevk = 1.1605d+4, mh = 1.67d-24, pi=3.14159d0,
     &           tiny_cgs = 1.0d-37, dhuge = 1.0d+30,
     &           kb = 1.380658d-16)

c
c         Set various cutoff values in eV.
c
      double precision e24,e25,e26,e27,e28a,e28b,e28c,e29a,e29b,e29c,
     &                 e30a,e30b
        PARAMETER(
     &    e24  = 13.6d0
     &   ,e25  = 54.4d0
     &   ,e26  = 24.6d0
     &   ,e27  = 0.755d0
     &   ,e28a = 2.65d0
     &   ,e28b = 11.27d0
     &   ,e28c = 21.0d0
     &   ,e29a = 15.42d0
     &   ,e29b = 16.5d0
     &   ,e29c = 17.7d0
     &   ,e30a = 30.0d0
     &   ,e30b = 70.0d0)

c  Locals
c
      integer i,j,n
      double precision nu, delnu, hnu, logttt, ttt, tev, logtev, 
     &                 xx, dum, tbase1, xbase1, kunit, coolunit,
     &                 dbase1, dlogtem, kunit_3bdy, x, y
      double precision sigma24(nnu), sigma25(nnu), sigma26(nnu), 
     &                 sigma27(nnu), sigma28(nnu), sigma29(nnu), 
     &                 sigma30(nnu), sigma31(nnu), fext(nnu)
      real tm
      double precision HDLR, HDLV, lt, t3
c
c\\\\\\\\\\\\\\\\\\\\\\\\\\\\\\\\\\/////////////////////////////////////
c=======================================================================
c
c  Get conversion units
c
c     t/x/dbase1 is the number (z dependant) that converts from the
c       dimensionless code units to physical units.  Also, in the
c       code aye = 1 at z=zinit, so to convert the usual a (=1 at z=0)
c       to a~ (written in the code as aye), we use a = a~*[a] 
c
      tbase1 = utim
      xbase1 = uxyz/(aye*uaye)    ! uxyz is [x]*a     
      dbase1 = urho*(aye*uaye)**3 ! urho is [dens]/a^3
c
c  1) Set the dimensions of the (non-radiative) rate coefficients.  
c    Note that we have included the units that convert density to 
c    number density, so the rate equations should look like 
c    (in dimensionless units, hence the primes):
c
c       d(d0~)/dt~ = k~ * d1~ * d2~ / a~^3
c
c    where k~ is the dimenionless rate coefficients and d0-2~ are three
c     dimensionless densities (i.e. d = [dens]*d~) and a~ is the 
c     dimensionless expansion coefficient (see above).
c
c    rate eqn        : delta(n0)  = k  * n1        * n2        * dt     / a^3
c    rate eqn units  : [dens]/mh  = k  * [dens]/mh * [dens]/mh * [time] / [a]^3
c    rate eqn dimless: delta(n0~) = k~ * n1~       * n2~       * dt~    / a~^3
c    so: k = [k] * k~  where [k] = ( [a]^3 * mh ) / ( [dens] * [time] )  (~)
c    reminder: the number densities here are normalized with [dens] which
c              is not a constant (it has a factor a^3), so the number
c              densities must be converted from comoving to proper.
c
      kunit   = (uaye**3 * mh) / (dbase1 * tbase1)
      kunit_3bdy  = kunit * (uaye**3 * mh) / dbase1
c
c  2) Set the dimension of the cooling coefficients (including constants)
c     (this equation has a rho because e is the specific energy, not
c      energy/unit volume).
c       delta(e)  = L     * n1        * n2        * dt     / dens   / a^3
c       [e]       = L     * [dens]/mh * [dens]/mh * [time] / [dens] / [a]^3
c       delta(e~) = L~    * n1~       * n2~       * dt~    / dens~  / a~^3 [~]
c     so L = [L] * L~ where [L] = [e] * mh**2 * [a]^3 / ([dens] * [time]) [~]
c       but [e] = ([a]*[x])**2 / [time]**2 and ([a] = 1 / (1 + zri) )
c      [L] = ([a]**5 * [x]**2 * mh**2) / ([dens] * [time]**3)
c
      coolunit = (uaye**5 * xbase1**2 * mh**2) / (tbase1**3 * dbase1)
c
c    Note: some of the coffiecients have only one power of n.  These
c          do not have the /a^3 factor, also they have units
c          [L1] = ([a]**2 * [x]**2 * mh) / [time]**3
c               = [L] * [dens] * [a]**3 / mh
c          This is done through the dom variable in cool.src
c         (some have three powers of n and they are different by the
c          reciprocal of the above factor multiplying [L]).
c
c  3) the units for the radiative rate coefficients is just 1/[time]
c
c  Compute log spacing in temperature
c
      ttt    = temstart
      logttt = log(ttt)
      dlogtem= (log(temend) - log(temstart))/real(nratec-1)
c
c  Initialize constants to tiny
c
      do i = 1, nratec
c
        k1a (i) = tiny
        k2a (i) = tiny
        k3a (i) = tiny
        k4a (i) = tiny
        k5a (i) = tiny
        k6a (i) = tiny
        k7a (i) = tiny
        k8a (i) = tiny
        k9a (i) = tiny
        k10a(i) = tiny
        k11a(i) = tiny
        k12a(i) = tiny
        k13a(i) = tiny
        k14a(i) = tiny
        k15a(i) = tiny
        k16a(i) = tiny
        k17a(i) = tiny
        k18a(i) = tiny
        k19a(i) = tiny
        k20a(i) = tiny
        k21a(i) = tiny
        k22a(i) = tiny
        k50a(i) = tiny
        k51a(i) = tiny
        k52a(i) = tiny
        k53a(i) = tiny
        k54a(i) = tiny
        k55a(i) = tiny
        k56a(i) = tiny
c
        do j = 1, 7
           k13dda(i,j) = tiny
        enddo
c
        ceHIa(i)    = tiny
        ceHeIa(i)   = tiny
        ceHeIIa(i)  = tiny
        ciHIa(i)    = tiny
        ciHeIa(i)   = tiny
        ciHeISa(i)  = tiny
        ciHeIIa(i)  = tiny
        reHIIa(i)   = tiny
        reHeII1a(i) = tiny
        reHeII2a(i) = tiny
        reHeIIIa(i) = tiny
        brema(i)    = tiny
        compa       = tiny
c
        hyd01ka(i)  = tiny
        h2k01a(i)   = tiny
        vibha(i)    = tiny
        rotha(i)    = tiny
        rotla(i)    = tiny
        hdltea(i)   = tiny
        hdlowa(i)   = tiny
c
      enddo
c
c  Fill in tables over the range temstart to temend
c
c -------------------------------------------------
c  1) rate coefficients (excluding external radiation field)
c
      do i = 1, nratec
c
c       Compute temperature of this bin (in eV)
c
        logttt = log(temstart) + real(i-1)*dlogtem
        ttt = exp(logttt)
        tev = ttt/tevk
        logtev = log(tev)
c
#ifdef USE_NEW_COLLISIONAL_RATES
c
c       Call Tom Abels routine (from his web page, Feb/2000)
c
        call coll_rates(ttt, k1a(i), k2a(i), k3a(i), k4a(i), k5a(i),
     &                  k6a(i), k7a(i), k8a(i), k9a(i), k10a(i),
     &                  k11a(i), k12a(i), k13a(i), k14a(i), k15a(i), 
     &                  k16a(i), k17a(i), k18a(i), k19a(i), kunit)
c
#else /* USE_NEW_COLLISIONAL_RATES */
c
c       Compute rates in cgs: cm^3/s (note: tev, logtev is double)
c        (the k1-19 rates below are an older set)
c 
        if (tev .gt. 0.8) then
          k1a(i) = exp(-32.71396786375 
     &           + 13.53655609057*logtev
     &           - 5.739328757388*logtev**2 
     &           + 1.563154982022*logtev**3
     &           - 0.2877056004391*logtev**4
     &           + 0.03482559773736999*logtev**5
     &           - 0.00263197617559*logtev**6
     &           + 0.0001119543953861*logtev**7
     &           - 2.039149852002d-6*logtev**8) / kunit
c
          k3a(i) = exp(-44.09864886561001
     &           + 23.91596563469*logtev
     &           - 10.75323019821*logtev**2
     &           + 3.058038757198*logtev**3
     &           - 0.5685118909884001*logtev**4
     &           + 0.06795391233790001*logtev**5
     &           - 0.005009056101857001*logtev**6
     &           + 0.0002067236157507*logtev**7
     &           - 3.649161410833d-6*logtev**8) / kunit
c
          k4a(i) = (1.54d-9*(1.+0.3/exp(8.099328789667/tev))
     &           / (exp(40.49664394833662/tev)*tev**1.5)
     &           + 3.92d-13/tev**0.6353) / kunit
c
          k5a(i) = exp(-68.71040990212001
     &           + 43.93347632635*logtev
     &           - 18.48066993568*logtev**2
     &           + 4.701626486759002*logtev**3
     &           - 0.7692466334492*logtev**4
     &           + 0.08113042097303*logtev**5
     &           - 0.005324020628287001*logtev**6
     &           + 0.0001975705312221*logtev**7
     &           - 3.165581065665d-6*logtev**8) / kunit
        else
          k1a(i) = tiny
          k3a(i) = tiny
          k4a(i) = 3.92d-13/tev**0.6353 / kunit
          k5a(i) = tiny
        endif
c
        if ( ttt .gt. 5500.0 ) then
          k2a(i) = exp(-28.61303380689232
     &           - 0.7241125657826851*logtev
     &           - 0.02026044731984691*logtev**2
     &           - 0.002380861877349834*logtev**3
     &           - 0.0003212605213188796*logtev**4
     &           - 0.00001421502914054107*logtev**5
     &           + 4.989108920299513d-6*logtev**6
     &           + 5.755614137575758d-7*logtev**7
     &           - 1.856767039775261d-8*logtev**8
     &           - 3.071135243196595d-9*logtev**9) / kunit
        else
          k2a(i) = k4a(i)
        endif

        k6a(i) = 3.36d-10/sqrt(ttt)/(ttt/1.d3)**0.2/(1+(ttt/1.d6)**0.7)
     &           / kunit
c
c       Molecular hydrogen and associated rates
c
        k7a(i) = 6.77d-15*tev**0.8779 / kunit
c
        if (tev .gt. 0.1) then
          k8a(i) = exp(-20.06913897587003
     &           + 0.2289800603272916*logtev
     &           + 0.03599837721023835*logtev**2
     &           - 0.004555120027032095*logtev**3
     &           - 0.0003105115447124016*logtev**4
     &           + 0.0001073294010367247*logtev**5
     &           - 8.36671960467864d-6*logtev**6
     &           + 2.238306228891639d-7*logtev**7) / kunit
        else
          k8a(i) = 1.43d-9 / kunit
        endif
c
        k9a(i) = 1.85d-23*ttt**1.8 / kunit
        if(ttt .gt. 6.7d3) 
     &     k9a(i) = 5.81d-16*(ttt/56200)**(-0.6657*log10(ttt/56200)) 
     &              / kunit
c
        k10a(i) = 6.0d-10 / kunit
c
        if (tev .gt. 0.3) then
          k13a(i) = 1.0670825d-10*tev**2.012/
     &              (exp(4.463/tev)*(1+0.2472*tev)**3.512) / kunit
c
          k11a(i) = exp(-24.24914687731536
     &            + 3.400824447095291*logtev
     &            - 3.898003964650152*logtev**2
     &            + 2.045587822403071*logtev**3
     &            - 0.5416182856220388*logtev**4
     &            + 0.0841077503763412*logtev**5
     &            - 0.007879026154483455*logtev**6
     &            + 0.0004138398421504563*logtev**7
     &            - 9.36345888928611d-6*logtev**8) / kunit

          k12a(i) = 4.38d-10*exp(-102000.0/ttt)*ttt**0.35 / kunit
        else
          k13a(i) = tiny
          k11a(i) = tiny
          k12a(i) = tiny
        endif
c
        if (tev .gt. 0.04) then
          k14a(i) = exp(-18.01849334273
     &            + 2.360852208681*logtev
     &            - 0.2827443061704*logtev**2
     &            + 0.01623316639567*logtev**3
     &            - 0.03365012031362999*logtev**4
     &            + 0.01178329782711*logtev**5
     &            - 0.001656194699504*logtev**6
     &            + 0.0001068275202678*logtev**7
     &            - 2.631285809207d-6*logtev**8) / kunit
        else
          k14a(i) =  tiny
        endif
c
        if (tev .gt. 0.1) then
          k15a(i) = exp(-20.37260896533324
     &            + 1.139449335841631*logtev
     &            - 0.1421013521554148*logtev**2
     &            + 0.00846445538663*logtev**3
     &            - 0.0014327641212992*logtev**4
     &            + 0.0002012250284791*logtev**5
     &            + 0.0000866396324309*logtev**6
     &            - 0.00002585009680264*logtev**7
     &            + 2.4555011970392d-6*logtev**8
     &            - 8.06838246118d-8*logtev**9) / kunit
        else
          k15a(i) = 2.56d-9*tev**1.78186 / kunit
        endif
c
        k16a(i) = 6.5d-9/sqrt(tev) / kunit
c
        k17a(i) = 1.0d-8*ttt**(-0.4) / kunit
        if(ttt.gt.1.0d4)
     &     k17a(i)=4.0d-4*ttt**(-1.4)*exp(-15100.0/ttt) / kunit
c
        k18a(i) = 5.56396d-8/tev**0.6035 / kunit
        k19a(i) = 4.64d-8/sqrt(tev) / kunit
c
#endif /* USE_NEW_COLLISIONAL_RATES */
c
c       Compute the density-dependant collision H2 dissociation rates.
c          (this givens a 7-variable function that is used to compute
c           the log10 of the rate -- normalize by dividing by kunit).
c
        call colh2diss(ttt, k13dda(i,1), k13dda(i,2), k13dda(i,3), 
     &                 k13dda(i,4), k13dda(i,5), k13dda(i,6), 
     &                 k13dda(i,7))
        k13dda(i,1) = k13dda(i,1) - log10(kunit)
c
c       ------ 3-body H2 rate ----
c       The first bit is my fit to A.E. Orel 1987, J.Chem.Phys., 87, 
c       314. I then match it to the T^-1 of Palla etal (1983) 
c       Which is then 4 times smaller than the Palla rate ! Thats
c       molecule uncertainties ! :-)
c
        if (ttt .le. 300.0) then
           k22a(i) = 1.3d-32 * (ttt/300.0)**(-0.38) / kunit_3bdy
        else
           k22a(i) = 1.3d-32 * (ttt/300.0)**(-1.0) / kunit_3bdy
        endif
c
c       ------ Deuterium rates -----
c
c       50) H+  +  D  -> H  +  D+
c       51) H   +  D+ -> H+ +  D
c       52) H2  +  D+ -> HD +  H+
c       53) HD  +  H+ -> H2 +  D+
c       54) H2  +  D  -> HD +  H
c       55) HD  +  H  -> H2 +  D
c       56) D   +  H- -> HD +  e-
c      [57) D-  +  H  -> HD +  e-]  included by multiply 56 by 2
c
        k50a(i) = 1.0d-9 *exp(-4.1d1/ttt)  / kunit
        k51a(i) = 1.0d-9                   / kunit
        k52a(i) = 2.1d-9                   / kunit
        k53a(i) = 1.0d-9 *exp(-4.57d2/ttt) / kunit
        k54a(i) = 7.5d-11*exp(-3.82d3/ttt) / kunit
        k55a(i) = 7.5d-11*exp(-4.24d3/ttt) / kunit
        k56a(i) = 1.5d-9 *(ttt/300.0)**(-0.1) / kunit
c
      enddo
c
c  Write out cgs rate coefficients
c
c      open(10, file='k1-6.out', status='unknown')
c      do i = 1, nratec
c         write(10,1010) exp(log(temstart) + real(i-1)*dlogtem),
c     &               k1a(i),k2a(i),k3a(i),k4a(i),k5a(i),k6a(i)
c      enddo
c 1010 format(7(1pe11.4))
c      close(10)
c
c -------------------------------------------------
c  2) Cooling/heating rates (excluding external radiation field)
c
      do i = 1, nratec
c
        logttt = log(temstart) + real(i-1)*dlogtem
        ttt = exp(logttt)
c
c    a) Collisional excitations (Black 1981; Cen 1992)
c
      if ( iceco .eq. 1 ) then
        ceHIa(i)    = 7.5d-19*exp(-min(log(dhuge),118348/ttt))
     &              /(1.+sqrt(ttt/1.0d5)) / coolunit
        ceHeIa(i)   = 9.1d-27*exp(-min(log(dhuge),13179/ttt))
     &              *ttt**(-0.1687)/(1.+sqrt(ttt/1.0d5)) / coolunit
        ceHeIIa(i)  = 5.54d-17*exp(-min(log(dhuge),473638/ttt))
     &              *ttt**(-0.397)/(1.+sqrt(ttt/1.0d5)) / coolunit
      else
        ceHIa(i)    = tiny
        ceHeIa(i)   = tiny
        ceHeIIa(i)  = tiny
      endif
c
c    b) Collisional ionizations (Cen 1992 or Abel 1996)
c
      if ( icico .eq. 1 ) then
c
c        ciHIa(i)    = 1.27d-21*sqrt(ttt)/(1.+sqrt(ttt/1.0d5))
c     &              *exp(-min(log(dhuge),157809.1/ttt)) / coolunit
c        ciHeIa(i)   = 9.38d-22*sqrt(ttt)/(1.+sqrt(ttt/1.0d5))
c     &              *exp(-min(log(dhuge),285335.4/ttt)) / coolunit
c        ciHeIIa(i)  = 4.95d-22*sqrt(ttt)/(1.+sqrt(ttt/1.0d5))
c     &              *exp(-min(log(dhuge),631515.0/ttt)) / coolunit
c
        ciHeISa(i)  = 5.01d-27*(ttt)**(-0.1687)/(1.+sqrt(ttt/1.0d5))
     &              *exp(-min(log(dhuge),55338/ttt)) / coolunit
c
c       Collisional ionizations (polynomial fits from Tom Abel)
c
        ciHIa(i)    = 2.18d-11*k1a(i)*kunit / coolunit
        ciHeIa(i)   = 3.94d-11*k3a(i)*kunit / coolunit
        ciHeIIa(i)  = 8.72d-11*k5a(i)*kunit / coolunit
      else
        ciHIa(i)    = tiny
        ciHeIa(i)   = tiny
        ciHeIIa(i)  = tiny
        ciHeISa(i)  = tiny
      endif
c
c    c) Recombinations (Cen 1992)
c
      if ( ireco .eq. 1 ) then
        reHIIa(i)   = 8.70d-27*sqrt(ttt)*(ttt/1000.0)**(-0.2)
     &              / (1.0 + (ttt/1.0d6)**(0.7)) / coolunit
        reHeII1a(i) = 1.55d-26*ttt**0.3647 / coolunit
        reHeII2a(i) = 1.24d-13*ttt**(-1.5)
     &              *exp(-min(log(dhuge),470000/ttt))
     &              *(1.+0.3*exp(-min(log(dhuge),94000/ttt))) 
     &              / coolunit
        reHeIIIa(i) = 3.48d-26*sqrt(ttt)*(ttt/1000.0)**(-0.2)
     &              / (1.0 + (ttt/1.0d6)**(0.7)) / coolunit
      else
        reHIIa(i)   = tiny
        reHeII1a(i) = tiny
        reHeII2a(i) = tiny
        reHeIIIa(i) = tiny
      endif
c
c    d) Bremsstrahlung (Black 1981)(Spitzer & Hart 1979)
c
      if ( ibrco .eq. 1 ) then
        brema(i)    = 1.43d-27*sqrt(ttt)
     &              *(1.1+0.34*exp(-(5.5-log10(ttt))**2/3.0)) / coolunit
      else
        brema(i)    = tiny
      endif
c
c         Bremsstrahlung (Shapiro & Kang 1987)(Spitzer 1978)
c
c        if(ttt .lt. 1.0*3.2d5) gaunt = 0.79464 + 0.1243*log10(ttt/1.)
c        if(ttt .ge. 1.0*3.2d5) gaunt = 2.13164 - 0.1240*log10(ttt/1.)
c        brem1a(i) = 1.426d-27*sqrt(ttt)*gaunt / coolunit
c        if(ttt .lt. 4.0*3.2d5) gaunt = 0.79464 + 0.1243*log10(ttt/4.)
c        if(ttt .ge. 4.0*3.2d5) gaunt = 2.13164 - 0.1240*log10(ttt/4.)
c        brem2a(i) = 1.426d-27*sqrt(ttt)*gaunt / coolunit
c
c    e) Molecular hydrogen cooling
c
c       (Which one is used is controlled by a flag in cool1d_mulit.src)
c
c    e part1) - Lepp & Shull rates
c
        xx = log10(ttt/1.0d4)
        vibha  (i) = 1.1d-18*exp(-min(log(dhuge),6744.0/ttt)) 
     &               / coolunit
c
        if( ttt .gt. 1635) then
           dum = 1.0d-12*sqrt(ttt)*exp(-1000.0/ttt)
        else
           dum = 1.4d-13*exp((ttt/125.) - (ttt/577.)**2)
        endif
        hyd01ka(i) = dum*exp(-min(log(dhuge),
     &                            8.152d-13/(1.38d-16*ttt) )) 
     &               / coolunit
c
        dum = 8.152d-13*(4.2/(1.38d-16*(ttt+1190))+1./(1.38d-16*ttt))
        h2k01a(i) = 1.45d-12*sqrt(ttt)*exp(-min(log(dhuge),dum)) 
     &              / coolunit
c
        if(ttt .gt. 4031)then
          rotla(i) = 1.38d-22*exp(-9243.0/ttt) / coolunit
        else
          rotla(i) = 10.0**(-22.9 - 0.553*xx - 1.148*xx*xx) / coolunit
        endif
c
        if(ttt .gt. 1087)then
          rotha(i) = 3.9d-19*exp(-6118.0/ttt) / coolunit
        else
          rotha(i) = 10.0**(-19.24 + 0.474*xx - 1.247*xx*xx) / coolunit
        endif
c
c     e part2) - Galli and Palla (1999) rates as fit by Tom Abel
c
        tm = max(ttt, 13.0d0)       ! no cooling below 13 Kelvin ...
        tm = min(tm, 1.e5)      ! fixes numerics
        lt = log10(tm)
c     low density limit from Galli and Palla
        gpldl(i)  = 10.**(-103.0+97.59*lt-48.05*lt**2+10.80*lt*lt*lt
     &       -0.9032*lt*lt*lt*lt) / coolunit
c     high density limit from HM79 (typo corrected Aug 30/2007)
        t3 = tm/1000.
        HDLR = ((9.5e-22*t3**3.76)/(1.+0.12*t3**2.1)*
     &            exp(-(0.13/t3)**3)+3.e-24*exp(-0.51/t3))
        HDLV = (6.7e-19*exp(-5.86/t3) + 1.6e-18*exp(-11.7/t3))
        gphdl(i)  = (HDLR + HDLV) / coolunit
c
c     f) HD cooling
c
c        HD Cooling Function (ergs cm3 /s)
c
c Fit to Lepp and Shull 1984 LTE (ergs/s) -> hdlte (ergs cm3/s)
c
      hdltea(i) = -35.6998d0 + 15.35716d0*dlog10(ttt) -
     &            5.58513d0   * (dlog10(ttt))**2 +
     &            0.8561149d0 * (dlog10(ttt))**3 - 
     &            1.75538d-2  * (dlog10(ttt))**4
      hdltea(i) = (10.0**min(hdltea(i), 15.0)) / coolunit
c      hdlte=10**lhdlte/den
c
c Galli and Palla 1998 low density limit (erg cm3 /s)
c  uses HD-He collisional data, so reduced by 1.27
c
      hdlowa(i) = ((3.0d0 * (4.4d-12 + 3.6d-13*ttt**0.77) *
     &            dexp(-128.d0/ttt) * 128.d0 +
     &            (5.0d0/3.0d0) * (4.1d-12+2.1e-13*ttt**0.92) *
     &            dexp(-255.d0/ttt) * 255.d0) * kb/1.27d0 ) / coolunit
c      hdcool=hdlte/(1.0d0+hdlte/hdlow)
c
      enddo
c
c     g) Compton cooling
c
      compa   = 5.65d-36 / coolunit
c
c -------------------------------------------------
c  3) External radiative processes
c
c     Initialize to tiny
c
      k24 = tiny
      k25 = tiny
      k26 = tiny
      k27 = tiny
      k28 = tiny
      k29 = tiny
      k30 = tiny
      k31 = tiny
c      
      piHI   = tiny
      piHeI  = tiny
      piHeII = tiny
c
c     Loop over all frequency bins, compute cross-sections
c        nu is in eV (range is 0.74 eV to 7.2 keV if nnu=400)
c
      do n = 1, nnu
        nu = 10**((n-14)*0.01)                    ! 0.74 ev -- 7.24d9 ev
c
c       24) HI photo-ionization cross-section
c
        if (nu .gt. e24) then
           dum = sqrt(nu/e24-1)
           sigma24(n) = 6.3d-18 * (e24/nu)**4
     &               *  exp(4.0-4.0*atan(dum)/dum)
     &               / (1-exp(-2.0*pi/dum))
        else
           sigma24(n) = 0.0
        endif
c
c       25) HeII photo-ionization cross-section
c
        if (nu .gt. e25) then
           dum = sqrt(nu/e25-1)
           sigma25(n) = 1.58d-18 * (e25/nu)**4
     &               * exp(4.0-4.0*atan(dum)/dum)
     &               / (1-exp(-2.0*pi/dum))
        else
           sigma25(n) = 0.0
        endif
c
c       26) HeI photo-ionization cross-section
c
        if (nu .gt. e26) then
c           sigma26(n) = 7.42d-18*(1.66*(nu/e26)**(-2.05)
c     &                          - 0.66*(nu/e26)**(-3.05))
c
c            Now From Verland, et al (1996)
c
           x = nu/13.61 - 0.4434
           y = sqrt(x**2 + 2.136**2)
           sigma26(n) = 9.492d-16*((x-1)**2 + 2.039**2) *
     &                  y**(0.5*3.188-5.5) *
     &                  (1.0 + sqrt(y/1.469))**(-3.188)        
	else
           sigma26(n) = 0.0
        endif
c
c       27) HM    + p   -> HI    + e
c
        if (nu .gt. e27) then
          sigma27(n) = 2.11d-16*(nu-e27)**1.5/nu**3
        else
          sigma27(n) = 0.0
        endif
c
c       28) H2II  + p   -> HI    + HII
c
        if (nu .gt. e28a .and. nu .le. e28b) then
          sigma28(n) = 10**(-40.97+6.03*nu-0.504*nu**2+1.387d-2*nu**3)
        elseif (nu .gt. e28b .and. nu .lt. e28c) then
          sigma28(n) = 10**(-30.26+2.79*nu-0.184*nu**2+3.535d-3*nu**3)
        else
          sigma28(n) = 0.0
        endif
c
c       29) H2I   + p   -> H2II  + e
c
        if (nu .gt. e29a .and. nu .le. e29b) then
          sigma29(n) = 6.2d-18*nu - 9.4d-17
        elseif (nu .gt. e29b .and. nu .le. e29c) then
          sigma29(n) = 1.4d-18*nu - 1.48d-17
        elseif (nu .gt. e29c) then
          sigma29(n) = 2.5d-14*nu**(-2.71)
        else
          sigma29(n) = 0.0
        endif
c
c       30) H2II  + p   -> 2HII  + e
c
        if (nu .ge. e30a .and. nu .lt. e30b) then
          sigma30(n) = 10**(-16.926-4.528d-2*nu+2.238d-4*nu**2
     &                      +4.245d-7*nu**3)
        else
          sigma30(n) = 0.0
        endif
c
c       31) H2I   + p   -> 2HI
c
c        sigma31(n) = 0.0

C	 BUG FIX: (n) subscripts below were (i)
C	 reported by Patrick McDonald, modified in revision 2310

        if (nu .ge. e28b .and. nu .lt. e24) THEN
           sigma31(n) = 3.71e-18
        else
           sigma31(n) = 0.0
        endif
c
c       Set the radiation field (only used in iradtype = 5 or 8)

        if (iradtype .eq. 5) then
c
c          5) A featureless power-low spectrum (quasar-like)
c         (f3 is the amplitude of flux at the HI ionization edge)
c
           fext(n) = f3*(nu/e24)**(-alpha0)
c
        elseif (iradtype .eq. 8) then
c
c          8) power law with absorbing neutral medium with
c              10^22 cm^-2 column density
c         (f3 is the amplitude of flux at the HI ionization edge,
c          the min is to prevent under-runs).
c
           fext(n) = f3*(nu/12.86)**(-1.0)*exp(-1.0*
     &       min(1.0d22*(sigma24(n) + 0.08*sigma26(n)), 100.0d0))
        else
c
c          otherwise) no radiation field
c
           fext(n) = 0.0
        endif
      enddo
c
c     Integrate over the frequency spectrum
c

!     call open_mpi_error_file( 'F89', 89, 'unknown' )

      do n = 1, nnu
c          delnu = (10**((n-14)*0.01) - 10**((n-15)*0.01))*evhz
c          hnu   = 0.5*(10**((n-14)*0.01) + 10**((n-15)*0.01))*everg

          delnu = (10**((n-13.5)*0.01) - 10**((n-14.5)*0.01))*evhz
          hnu   = 10**((n-14)*0.01)*everg
c
c         Add to radiative rate coefficients
c
c            k24-k26) HI, HeII, HeI photo-ionization
c
          k24 = k24 + fext(n)*sigma24(n) * delnu/hnu * tbase1
          k25 = k25 + fext(n)*sigma25(n) * delnu/hnu * tbase1
          k26 = k26 + fext(n)*sigma26(n) * delnu/hnu * tbase1
c
c            k27-k31) HM, H2II, H2I, H2II, H2I photo-dissociation
c
          if (iphoto_mol .eq. 1) then
             k27 = k27 + fext(n)*sigma27(n) * delnu/hnu * tbase1
             k28 = k28 + fext(n)*sigma28(n) * delnu/hnu * tbase1
             k29 = k29 + fext(n)*sigma29(n) * delnu/hnu * tbase1
             k30 = k30 + fext(n)*sigma30(n) * delnu/hnu * tbase1
             k31 = k31 + fext(n)*sigma31(n) * delnu/hnu * tbase1
          endif
c
c            Add to HI,HeI,HeII photoionization heating coefficients
c
          dum = max(hnu - e24*everg, 0.0d0)
          piHI = piHI + fext(n)*sigma24(n)*delnu/hnu*dum/coolunit

!         write(89,1089) hnu, fext(n),sigma24(n),k24, dum, piHI

 1089     format(1p,10(e14.4,1x))
c
          dum = hnu - e26*everg
          piHeI = piHeI + fext(n)*sigma26(n)*delnu/hnu*dum/coolunit
c
          dum = hnu - e25*everg
          piHeII = piHeII + fext(n)*sigma25(n)*delnu/hnu*dum/coolunit
      enddo

!     call close_mpi_error_file( 89 )

c
c     The secondary electrons produced by steep spectral sources
c       are accounted for here using data from Shull & Steenberg.
c       Note that this is only good for ionization fractions < 1e-3 or so
c          (now down in the solvers to cover all ionization fractions)
c
#ifdef UNUSED
      if (iradtype .eq. 8) then
c
c        Add two terms which account for the amount of energy that
c          goes into secondary ionizations of HI and HeI (ignore HeII)
c          (We use piHI and piHeI which is the available energy and
c           then divide by the ionization energy of HI).  We also need
c           to convert from heating units (erg/s) to rate units (/s)
c          (The 0.08 accounts for the ratio n(HI)/n(HeI) where we are
c           assuming full neutrality)
c
         k24 = k24 + 0.39*(piHI + 0.08*piHeI)/(e24*everg)
     &         * coolunit * tbase1
         k26 = k26 + 0.055*(piHI/0.08 + piHeI)/(e26*everg)
     &         * coolunit * tbase1
c
c        Now correct the photo-heating rates to reflect that much of
c          the energy does not go into heating. (ignore HeII)
c
         piHI  = piHI *0.13
         piHeI = piHeI*0.13
c         
      endif
#endif /* UNUSED */

#define OUTPUT_COOLING_RATE_DATA
#ifdef OUTPUT_COOLING_RATE_DATA
c
c     Write out cooling rate data
c

      if (ioutput.eq.1) then
        open(10, file='cool_rates.out', status='unknown')

!      call open_mpi_error_file( 'cool_rates.out', 10, 'unknown' )

        do i=1, nratec
           logttt = log(temstart) + real(i-1)*dlogtem
           ttt = exp(logttt)
           write(10,1000) ttt, ceHIa(i), ceHeIa(i), ceHeIIa(i),
     &                    ciHIa(i), ciHeIa(i), ciHeISa(i),
     &                    ciHeIIa(i), reHIIa(i), reHeII1a(i),
     &                    reHeII2a(i), reHeIIIa(i), brema(i),
     &                    compa       
 1000      format(1p,30(e10.3,1x))
        enddo
        write(10,*) 'piHI=', piHI*coolunit
        write(10,*) 'piHeI=', piHeI*coolunit
        write(10,*) 'piHeII=', piHeII*coolunit

!      call close_mpi_error_file( 10 )

        close(10)

c
c     Write out species rate data
c

        open(10, file='rates.out', status='unknown')

!      call open_mpi_error_file( 'rates.out', 10, 'unknown' )

        do i=1, nratec
           logttt = log(temstart) + real(i-1)*dlogtem
           ttt = exp(logttt)
           write(10,1000) ttt, k1a(i), k2a(i), k3a(i), k4a(i), k5a(i),
     &                    k6a(i), k7a(i), k8a(i), k9a(i), k10a(i),
     &                    k11a(i), k12a(i), k13a(i), k14a(i), k15a(i), 
     &                    k16a(i), k17a(i), k18a(i), k19a(i), k22a(i)
        enddo
        write(10,*) 'k24=', k24/tbase1
        write(10,*) 'k25=', k25/tbase1
        write(10,*) 'k26=', k26/tbase1
        write(10,*) 'k27=', k27/tbase1
        write(10,*) 'k28=', k28/tbase1
        write(10,*) 'k29=', k29/tbase1
        write(10,*) 'k30=', k30/tbase1
        write(10,*) 'k31=', k31/tbase1

!      call close_mpi_error_file( 10 )

       close(10)

      endif

#endif /* OUTPUT_COOLING_RATE_DATA */

c
      return
      end