File: equib_routines.f90

package info (click to toggle)
neat 2.3.2-2
  • links: PTS, VCS
  • area: main
  • in suites: bookworm, bullseye, sid, trixie
  • size: 4,108 kB
  • sloc: f90: 7,385; python: 211; makefile: 78; sh: 64
file content (1100 lines) | stat: -rw-r--r-- 34,631 bytes parent folder | download
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
! For input atomic data file formats see Readme file attached at the end
! **********************************************************************
!          Program EQUIB  (FORTRAN 77)
!
!    Programming history:
!
!    1981 May 3    IDH    Version 1
!    1981 May 5    IDH    Minibug fixed!
!    1981 May 7    IDH    Now takes collision rates or strengths
!    1981 Aug 3    SA     Interpolates collision strengths
!    1981 Aug 7    SA     Input method changed
!    1984 Nov 19   RESC   SA files entombed in scratch disk. Logical
!                         filenames given to SA's data files.
!    1995 Aug      DPR    Changed input file format. Increased matrices.
!    1996 Feb      XWL    Tidy up. SUBROUTINES SPLMAT, HGEN, CFY and CFD
!                         modified such that matrix sizes (i.e. maximum
!                         of Te and maximum no of levels) can now be cha
!                         by modifying the parameters NDIM1, NDIM2 and N
!                         in the Main program. EASY!
!                         Now takes collision rates as well.
!                         All variables are declared explicitly
!                         Generate two extra files (ionpop.lis and ionra
!                         of plain stream format for plotting
!    1996 June     CJP    Changed input data format for cases IBIG=1,2.
!                         Fixed readin bug for IBIG=2 case.
!                         Now reads reformatted upsilons (easier to see
!                         and the 0 0 0 data end is excluded for these c
!                         The A values have a different format for IBIG=
!    2009 April    RW     Converted to F90, inputs from cmd line, version
!                         written purely to do diagnostics.
!    2012 February RW     Tidying up, combining into simpler routines for NEAT
! ***** N.B!!  NO TRAPS FOR BAD DATA!!  TAKE CARE!! ****C
!
      module mod_equib
      use mod_globals
      contains

      subroutine get_diagnostic(ion,levu,levl,inratio,diagtype,fixedq,result,ndim2,ndim1,atomicdata,iion,tlower,tupper)
      use mod_atomicdata
      use mod_helium
      implicit none

      integer :: NDIM1, NDIM2, NDIM1T3, MAXND
                                                      !Maximum no of Te & levels
                                             !NDIM1T3 should be at least 3*NDIM1
                                                   !Maximum no. of Ne increments
      parameter (MAXND=100)
      integer :: G(NDIM2),                                                 &
     &  ITRANA(2,NDIM2),ITRANB(2,NDIM2),ITRANC(2,NDIM2),LOOP
      type(atomic_data),dimension(:),intent(in) :: atomicdata
      integer :: iion,nion
      real(kind=dp) :: N(NDIM2)
      real(kind=dp) :: & !TDRAT(2,MAXND)
     & TNIJ(NDIM2,NDIM2), FINTIJ(NDIM2,NDIM2),                          &
     & WAVA(NDIM2), WAVB(NDIM2), WAVC(NDIM2), CS(NDIM2,NDIM2),          &
     & QEFF(NDIM2,NDIM2), QQ(NDIM1),                                    &
     & QOM(NDIM1,NDIM2,NDIM2), A(NDIM2,NDIM2), E(NDIM2), T(NDIM1),      &
     & ROOTT(NDIM1), X(NDIM2,NDIM2), Y(NDIM2),                          &
     & X2(NDIM2,NDIM2), XKEEP(NDIM2,NDIM2), Y2(NDIM2), YKEEP(NDIM2),    &
     & HMH(NDIM1,NDIM1), D(NDIM1)
      character(len=20) :: LABEL(NDIM2)
      character(len=10) :: ION
      integer :: I, I1, I2, J, KK, LL, JT, JJD,                            &
     & NLEV, NTEMP, IBIG, IRATS,                                        &
     & NLEV1, INT, IND, IOPT, IT, IM1, JM1, IP1,                        &
     & IAPR, IBPR, ICPR, IKT, IA, IB, IC, IA1, IA2, IB1, IB2, IC1, IC2
      real(kind=dp) :: TEMPI, TINC, DENSI, DINC, DENS, DLOGD, TEMP, TLOGT,        &
     & TEMP2, DD, DELTEK, EXPE, VALUE, SUMN, TTT, TTP, AHB, EJI, WAV,   &
     & RLINT, FINT, SUMA, SUMB, FRAT, DEE

      real(kind=dp) :: fixedq
      real(kind=dp) :: inratio,result
      character(len=20) :: levu,levl
      character(len=1) :: diagtype
      real(kind=dp), DIMENSION(:,:), ALLOCATABLE :: RESULTS
      real(kind=dp) :: valtest(3)
      integer :: test
      real(kind=dp) :: tlower,tupper ! lower and upper limits for temperature calculations

!debugging
#ifdef CO
        print *,"subroutine: get_diagnostic. ion=",ion,levu,levl,inratio,diagtype,fixedq
#endif

! check if ratio is meaningful

      if (inratio .le. 0.d0 .or. inratio .gt. 1.e10) then
        result=0.d0
        return
      endif

      ndim1t3=3*ndim1
      g=0
      itrana=0
      itranb=0
      itranc=0
      valtest=0.d0
      test=0
      d=0.d0

      read(levu,*) ((ITRANA(LL,KK),LL=1,2),KK=1,ndim2)!150)
      read(levl,*) ((ITRANB(LL,KK),LL=1,2),KK=1,ndim2)!150)


!Transfer atomic data to local variables
      nion = 0
      do i = 1,iion
         if(atomicdata(i)%ion .eq. ion) nion=i
      enddo
      if (nion .eq. 0) then
          print *,"I'm afraid. I'm afraid, Dave."
          nion = 1
      endif
!          print*,nion,atomicdata(nion)%ion,ion
      nlev=atomicdata(nion)%nlevs
      ntemp=atomicdata(nion)%ntemps
          T(1:ntemp)= log10(atomicdata(nion)%Temps(1:ntemp))
      ROOTT(1:ntemp)=atomicdata(nion)%rootT(1:ntemp)
      Label(1:nlev)=atomicdata(nion)%labels(1:nlev)
      QOM(1:ntemp,1:nlev,1:nlev)=atomicdata(nion)%col_str(1:ntemp,1:nlev,1:nlev)
      A(1:nlev,1:nlev)=atomicdata(nion)%A_coeffs(1:nlev,1:nlev)
      E(1:nlev)=atomicdata(nion)%Waveno(1:nlev)
      G(1:nlev)=atomicdata(nion)%G(1:nlev)

      irats=0
      ibig=0

      NLEV1 = NLEV - 1
      ITRANC = 0

      !*****LOOP STARTS HERE*************************
      do LOOP = 1, 9
      if (diagtype .eq. "t" .or. diagtype .eq. "T") then

        if (LOOP .eq. 1) then
                TEMPI=tlower
        else
                TEMPI= valtest(1)
        endif

        INT=4
        TINC=(tupper-tlower)/((INT-1)**(LOOP))

        densi=fixedq
        dinc=0
        ind=1

        allocate(RESULTS(3,INT))

      else

        if (LOOP .eq. 1) then
          densi=0
        else
          densi=valtest(2)
        endif

        IND=4
        DINC=(100000.)/((IND-1)**(LOOP))

        TempI=fixedq
        TINC=0
        INT=1

        allocate(results(3,IND))
      endif

      if (densi .le. 0) densi=1
          if (tempi .lt. tlower) tempi=tlower

                                                               !Start of Te loop
      do JT = 1, INT
        TEMP=TEMPI+(JT-1)*TINC
                                                               !Start of Ne loop

        do JJD = 1, IND
          DENS=DENSI+(JJD-1)*DINC
          if (TEMP.LE.0.D0.OR.DENS.LE.0.D0) then
            write (6,6100)
                print *,"Temp = ", TEMP, ", Dens = ", DENS, ", Ion = ",ion,diagtype
            call exit(200)
          endif
          DLOGD = log10(DENS)
          TLOGT = log10(TEMP)
          TEMP2= sqrt(TEMP)
                                                                 !Form matrices

          X=0.0
          CS=0.0
          QEFF=0.0
          TNIJ=0.0
          Y=0.0

          IOPT=0
          if (NTEMP.EQ.1) then
            write (6,*)
            write (6,*)                                                 &
     &      'Coll. strengths available for 1 Te only - assuming const'
          elseif (NTEMP.EQ.2) then
            write (6,*)
            write (6,*)                                                 &
     &      'Coll. strengths available for 2 Te only - linear interp'
          else
            call SPLMAT(T, NTEMP, IOPT, NDIM1, NDIM1T3, HMH)
            call CFD(TLOGT,T,NTEMP,NDIM1,HMH,D)
          endif
          do I = 2, NLEV
            do J = I, NLEV
                                                           !Negative!
              DELTEK = (E(I-1)-E(J))*1.4388463D0
              EXPE = EXP(DELTEK/TEMP)
              do IT = 1, NTEMP
                if (IRATS.EQ.0.D+00) then
                  QQ(IT) = QOM(IT,I-1,J)
                else
                                                      !Take out the exp. depend.
                  QQ(IT) = QOM(IT,I-1,J) / EXPE
                                                      !before interpolation
                endif
              enddo
              if (NTEMP.EQ.1) then
                 DD = QQ(1)
              elseif (NTEMP.EQ.2) then
                 DD = QQ(1) +                                           &
     &            (QQ(2) - QQ(1))/(T(2) - T(1)) * (TLOGT - T(1))
              else
                call CFY(TLOGT, DD, T, QQ, NTEMP, NDIM1, D)
              endif
              if (IRATS.EQ.0.D+00) then
                CS(I-1,J) = DD
              else
                CS(I-1,J) = DD * EXPE
              endif
              if (IRATS .EQ. 0.D+00) then
                QEFF(I-1,J) = 8.63D-06*CS(I-1,J) * EXPE / (G(I-1)*TEMP2)
                QEFF(J,I-1) = 8.63D-06 * CS(I-1,J) / (G(J)*TEMP2)
              else
                QEFF(I-1,J) = CS(I-1,J) * 10. ** IRATS
                                                                     !Be careful
                QEFF(J,I-1) = G(I-1) * QEFF(I-1,J) / (EXPE * G(J))
                                                                     !G integer!
              endif
            enddo
          enddo
          do I = 2, NLEV
            do J = 1, NLEV
              if (J.NE.I) then
                X(I,J) = X(I,J) + DENS * QEFF(J,I)
                X(I,I) = X(I,I) - DENS * QEFF(I,J)
                if (J.GT.I) then
                  X(I,J) = X(I,J) + A(J,I)
                else
                  X(I,I) = X(I,I) - A(I,J)
                endif
              endif
            enddo
          enddo
          do I = 2, NLEV
            IM1 = I - 1
            VALUE = 0.D0 - X(I,1)
            Y(IM1) = VALUE
            Y2(IM1) = VALUE
            YKEEP(IM1) = VALUE
            do J = 2, NLEV
              JM1 = J - 1
              VALUE = X(I,J)
              X(IM1,JM1) = VALUE
              X2(IM1,JM1) = VALUE
              XKEEP(IM1,JM1) = VALUE
            enddo
          enddo
                                              !Solve matrices for populations
          call LUSLV(X,Y,NLEV1,NDIM2)
          do I = NLEV, 2, -1
            N(I) = Y(I-1)
          enddo
          SUMN = 1.D0
          do I = 2, NLEV
            SUMN = SUMN + N(I)
          enddo
          do I = 2, NLEV
            N(I) = N(I) / SUMN
          enddo
          N(1) = 1.D0 / SUMN
                                                                  !Output data
          TTT=TEMP*1.0D-4
          TTP=TTT**(-0.87D0)
                                       !Eff. recombination coef. of Hb
!          AHB=3.036D-14*TTP
          ahb=10**(gamm4861(temp,dens)+11.38871) ! as used in RL routines.  differs by <<1% from previous simpler approximation unless ne>>1e4
          do I = 1, NLEV1
            IP1 = I + 1
            do J = IP1, NLEV
               if (A(J,I).NE.0.D0) then
                 EJI = E(J) - E(I)
                 WAV = 1.D8 / EJI
                 RLINT = A(J,I) * EJI
                 RLINT = RLINT *N(J)
                 TNIJ(I,J)=RLINT
                 FINT=N(J)*A(J,I)*4861.D0/(DENS*AHB*WAV)
                 FINTIJ(I,J)=FINT
               endif
            enddo
          enddo
                     !Search ITRANA, ITRANB & ITRANC for transitions & sum up
          SUMA=0.D0
          SUMB=0.D0

          IAPR=0
          IBPR=0
          ICPR=0
          do IKT = 1, NDIM2
            IA1=ITRANA(1,IKT)
            IA2=ITRANA(2,IKT)
            if (IA1.NE.0.AND.IA2.NE.0) then
              SUMA=SUMA+TNIJ(IA1,IA2)
              IAPR=IAPR+1
            endif
            IB1=ITRANB(1,IKT)
            IB2=ITRANB(2,IKT)
            if (IB1.NE.0.AND.IB2.NE.0) then
             IBPR=IBPR+1
             SUMB=SUMB+TNIJ(IB1,IB2)
            endif

            IC1=ITRANC(1,IKT)
            IC2=ITRANC(2,IKT)
            if (IC1.NE.0.AND.IC2.NE.0) then
             ICPR=ICPR+1
            endif
          enddo

          if (SUMB.eq.0.d0) then
            FRAT= 0.d0
          else
            FRAT=SUMA/SUMB
          endif

          if (diagtype .eq. "t" .or. diagtype .eq. "T") then
            RESULTS(1, JT) = TEMP
            RESULTS(2, JT) = DENS
            RESULTS(3, JT) = FRAT-inratio
          else
            RESULTS(1, JJD) = TEMP
            RESULTS(2, JJD) = DENS
            RESULTS(3, JJD) = FRAT-inratio
          endif                                   !End of the Ne loop
        enddo

        do IA = 1, IAPR
          I1=ITRANA(1,IA)
          I2=ITRANA(2,IA)
          DEE=E(I2)-E(I1)
          WAVA(IA)=1.D8/DEE
        enddo
        do IB = 1, IBPR
          I1=ITRANB(1,IB)
          I2=ITRANB(2,IB)
          DEE=E(I2)-E(I1)
          WAVB(IB)=1.D8/DEE
        enddo
        do IC = 1, ICPR
          I1=ITRANC(1,IC)
          I2=ITRANC(2,IC)
          DEE=E(I2)-E(I1)
          WAVC(IC)=1.D8/DEE
        enddo

      enddo                                              !End of the Te loop

! here, find the value in RESULTS which is closest to zero
! sort values in results to find two lowest values
!print*,results
!print*," "
      if (diagtype .eq. "D" .or. diagtype .eq. "d") then
        INT = ind
      endif

        ! loop through array and find out where the sign changes.

        do I=2,INT
            test=0
                if (sign(results(3,I),results(3,1)) .ne. results(3,I)) then !when this condition is fulfilled, the values in the array are now a different sign to the first value in the array
                        valtest(:) = (results(:,I-1)) ! return the value before the sign change so that the next loop starts at a sensible value
                        test=1
                        exit
                endif
        enddo

        if(test .eq. 0 .and. loop .lt. 9) then !test fails if no change of sign
                             !this kicks in then, and checks if it should be upper or lower limit
            if(abs(results(3,1)) .lt. abs(results(3,INT))) then
                valtest(:)=results(:,1)
            elseif(abs(results(3,INT)) .lt. abs(results(3,1))) then
                valtest(:)=results(:,INT-1)
            else                   !A simplistic work-around for this problem:
                result = -1d0      !it flags the value as ill-defined so that it
                deallocate(results)!can be dealt with without breaking the code.
                return             !It's set to zero later and excluded from
                                   !the averaging. Turns out long train rides
                                   !are good for sorting such problems.
                !print*,"Valtest failed"
                !print*,ion,levu,levl,loop,inratio,diagtype
                !print*,results
                !call exit(1)
            endif
        elseif(test .eq. 0 .and. loop .eq. 9) then !test fails if no change of sign
                             !this kicks in then, and checks if it should be upper or lower limit
            if(abs(results(3,1)) .lt. abs(results(3,INT))) then
                valtest(:)=results(:,1)
            elseif(abs(results(3,INT)) .lt. abs(results(3,1))) then
                valtest(:)=results(:,INT)
            else
                result = -1d0
                deallocate(results)
                return
                !print*,"Valtest failed"
                !print*,ion,levu,levl,loop,inratio,diagtype
                !print*,results
                !call exit(1)
            endif
        endif

        !LOOP = LOOP + 1
        DEALLOCATE(RESULTS) !thanks Bruce!
      enddo
      !********LOOP WOULD END HERE**********************

      if (diagtype .eq. "D" .or. diagtype .eq. "d") then
        result = valtest(2)
      else
        result = valtest(1)
      endif

      return

 6100 FORMAT (' PROCESSING COMPLETED'/                                  &
     & ' GOODBYE!!'///)


end subroutine get_diagnostic

      subroutine get_abundance(ion,levels,tempi,densi,iobs,abund,ndim2,ndim1,atomicdata,iion)
      use mod_atomicdata
      implicit none

          !INTEGER maxlevs,maxtemps
      integer :: NDIM1, NDIM2, NDIM1T3, MAXND
                                                      !Maximum no of Te & levels
      !PARAMETER (NDIM1=maxtemps, NDIM2=maxlevs)!35, NDIM2=150)
                                             !NDIM1T3 should be at least 3*NDIM1
      !PARAMETER (NDIM1T3 = 3*NDIM1)!105)
                                                   !Maximum no. of Ne increments
      PARAMETER (MAXND=100)
      integer :: G(NDIM2),                                                 &
     &  ITRANA(2,NDIM2),ITRANB(2,NDIM2),ITRANC(2,NDIM2)
      type(atomic_data),dimension(:),intent(in) :: atomicdata
      integer :: nion,iion
      real(kind=dp) :: N(NDIM2)
      real(kind=dp) :: TDRAT(2,MAXND), TNIJ(NDIM2,NDIM2), FINTIJ(NDIM2,NDIM2),    &
     & WAVA(NDIM2), WAVB(NDIM2), WAVC(NDIM2), CS(NDIM2,NDIM2),          &
     & QEFF(NDIM2,NDIM2), QQ(NDIM1),                                    &
     & QOM(NDIM1,NDIM2,NDIM2), A(NDIM2,NDIM2), E(NDIM2), T(NDIM1),      &
     & ROOTT(NDIM1), X(NDIM2,NDIM2), Y(NDIM2),                          &
     & X2(NDIM2,NDIM2), XKEEP(NDIM2,NDIM2), Y2(NDIM2), YKEEP(NDIM2),    &
     & HMH(NDIM1,NDIM1), D(NDIM1)
      character(len=20) :: LABEL(NDIM2)
      character(len=10) :: ION
      integer :: I, I1, I2, J, KK, LL, JT, JJD,                            &
     & NLEV, NTEMP, IBIG, IRATS,                                        &
     & NLEV1, INT, IND, IOPT, IT, IM1, JM1, IP1,                        &
     & IAPR, IBPR, ICPR, IKT, IA, IB, IC, IA1, IA2, IB1, IB2, IC1, IC2
      real(kind=dp) :: TEMPI, TINC, DENSI, DINC, DENS, DLOGD, TEMP, TLOGT,        &
     & TEMP2, DD, DELTEK, EXPE, VALUE, SUMN, TTT, TTP, AHB, EJI, WAV,   &
     & RLINT, FINT, SUMA, SUMB, SUMC, FRAT, DEE

      character(len=20) :: levels
      real(kind=dp) :: iobs, abund

!debugging
#ifdef CO
        print *,"subroutine: get_abundance. ion=",ion
#endif

      ndim1t3=3*ndim1
      g=0
      itrana=0
      itranb=0
      itranc=0
      d=0.d0

      read(levels,*) ((ITRANC(LL,KK),LL=1,2),KK=1,ndim2)!150)

      tinc=0
      dinc=0
      int=1
      ind=1

      nion = 0
      do i = 1,iion
         if(atomicdata(i)%ion .eq. ion(1:10)) nion=i
      enddo
      if (nion .eq. 0) then
          print *, "Dave, my mind is going. I can feel it."
          nion = 1
      endif
!          print*,nion,atomicdata(nion)%ion,ion
      nlev=atomicdata(nion)%nlevs
      ntemp=atomicdata(nion)%ntemps
          T(1:ntemp)=log10(atomicdata(nion)%Temps(1:ntemp))
      ROOTT(1:ntemp)=atomicdata(nion)%rootT(1:ntemp)
      Label(1:nlev)=atomicdata(nion)%labels(1:nlev)
      QOM(1:ntemp,1:nlev,1:nlev)=atomicdata(nion)%col_str(1:ntemp,1:nlev,1:nlev)
      A(1:nlev,1:nlev)=atomicdata(nion)%A_coeffs(1:nlev,1:nlev)
      E(1:nlev)=atomicdata(nion)%Waveno(1:nlev)
      G(1:nlev)=atomicdata(nion)%G(1:nlev)

      irats=0
      ibig=0

      NLEV1 = NLEV - 1
                                                               !Start of Te loop
      do JT = 1, INT
        TEMP=TEMPI+(JT-1)*TINC

                                                               !Start of Ne loop
        do JJD = 1, IND
          DENS=DENSI+(JJD-1)*DINC
          if (TEMP.LE.0.D0.OR.DENS.LE.0.D0) then
            write (6,6100)
            print *,"Temp = ", TEMP, ", Dens = ", DENS, ", Ion = ",ion
            call exit(200)
          endif
          DLOGD = log10(DENS)
          TLOGT = log10(TEMP)
          TEMP2= sqrt(TEMP)
                                                                  !Form matrices
          X = 0.D0
          CS=0.D0
          QEFF=0.D0
          TNIJ = 0.D0
          Y=0.D0

          IOPT=0
          if (NTEMP.EQ.1) then
            write (6,*)
            write (6,*)                                                 &
     &      'Coll. strengths available for 1 Te only - assuming const'
          elseif (NTEMP.EQ.2) then
            write (6,*)
            write (6,*)                                                 &
     &      'Coll. strengths available for 2 Te only - linear interp'
          else
            call SPLMAT(T, NTEMP, IOPT, NDIM1, NDIM1T3, HMH)
            call CFD(TLOGT,T,NTEMP,NDIM1,HMH,D)
          endif
          do I = 2, NLEV
            do J = I, NLEV
                                                           !Negative!
              DELTEK = (E(I-1)-E(J))*1.4388463D0
              EXPE = EXP(DELTEK/TEMP)
              do IT = 1, NTEMP
                if (IRATS.EQ.0.D+00) then
                  QQ(IT) = QOM(IT,I-1,J)
                else
                                                      !Take out the exp. depend.
                  QQ(IT) = QOM(IT,I-1,J) / EXPE
                                                      !before interpolation
                endif
              enddo
              if (NTEMP.EQ.1) then
                 DD = QQ(1)
              elseif (NTEMP.EQ.2) then
                 DD = QQ(1) +                                           &
     &            (QQ(2) - QQ(1))/(T(2) - T(1)) * (TLOGT - T(1))
              else
                call CFY(TLOGT, DD, T, QQ, NTEMP, NDIM1, D)
              endif
              if (IRATS.EQ.0.D+00) then
                CS(I-1,J) = DD
              else
                CS(I-1,J) = DD * EXPE
              endif
              if (IRATS .EQ. 0.D+00) then
                QEFF(I-1,J) = 8.63D-06*CS(I-1,J) * EXPE / (G(I-1)*TEMP2)
                QEFF(J,I-1) = 8.63D-06 * CS(I-1,J) / (G(J)*TEMP2)
              else
                QEFF(I-1,J) = CS(I-1,J) * 10. ** IRATS
                                                                     !Be careful
                QEFF(J,I-1) = G(I-1) * QEFF(I-1,J) / (EXPE * G(J))
                                                                     !G integer!
              endif
            enddo
          enddo
          do I = 2, NLEV
            do J = 1, NLEV
              if (J.NE.I) then
                X(I,J) = X(I,J) + DENS * QEFF(J,I)
                X(I,I) = X(I,I) - DENS * QEFF(I,J)
                if (J.GT.I) then
                  X(I,J) = X(I,J) + A(J,I)
                else
                  X(I,I) = X(I,I) - A(I,J)
                endif
              endif
            enddo
          enddo
          do I = 2, NLEV
            IM1 = I - 1
            VALUE = 0.D0 - X(I,1)
            Y(IM1) = VALUE
            Y2(IM1) = VALUE
            YKEEP(IM1) = VALUE
            do J = 2, NLEV
              JM1 = J - 1
              VALUE = X(I,J)
              X(IM1,JM1) = VALUE
              X2(IM1,JM1) = VALUE
              XKEEP(IM1,JM1) = VALUE
            enddo
          enddo
                                              !Solve matrices for populations
          call LUSLV(X,Y,NLEV1,NDIM2)
          do I = NLEV, 2, -1
            N(I) = Y(I-1)
          enddo
          SUMN = 1.D0
          do I = 2, NLEV
            SUMN = SUMN + N(I)
          enddo
          do I = 2, NLEV
            N(I) = N(I) / SUMN
          enddo
          N(1) = 1.D0 / SUMN
                                                                  !Output data
          TTT=TEMP*1.0D-4
          TTP=TTT**(-0.87D0)
                                       !Eff. recombination coef. of Hb
          AHB=3.036D-14*TTP
          do I = 1, NLEV1
            IP1 = I + 1
            do J = IP1, NLEV
               if (A(J,I).NE.0.D0) then
                 EJI = E(J) - E(I)
                 WAV = 1.D8 / EJI
                 RLINT = A(J,I) * EJI
                 RLINT = RLINT *N(J)
                 TNIJ(I,J)=RLINT
                 FINT=N(J)*A(J,I)*4861.D0/(DENS*AHB*WAV)
                 FINTIJ(I,J)=FINT
               endif
            enddo
          enddo
                     !Search ITRANA, ITRANB & ITRANC for transitions & sum up
          SUMA=0.D0
          SUMB=0.D0
          SUMC=0.D0
          IAPR=0
          IBPR=0
          ICPR=0
          do IKT = 1, NDIM2
            IA1=ITRANA(1,IKT)
            IA2=ITRANA(2,IKT)
            if (IA1.NE.0.AND.IA2.NE.0) then
              SUMA=SUMA+TNIJ(IA1,IA2)
              IAPR=IAPR+1
            endif
            IB1=ITRANB(1,IKT)
            IB2=ITRANB(2,IKT)
            if (IB1.NE.0.AND.IB2.NE.0) then
             IBPR=IBPR+1
             SUMB=SUMB+TNIJ(IB1,IB2)
            endif

            IC1=ITRANC(1,IKT)
            IC2=ITRANC(2,IKT)
            if (IC1.NE.0.AND.IC2.NE.0) then
             ICPR=ICPR+1
             SUMC=SUMC+FINTIJ(IC1,IC2)
            endif
          enddo

          if (sumb .gt. 0) then
            FRAT=SUMA/SUMB
          else
            FRAT=0.d0
          endif

          if (SUMC.ne.0.d0) SUMC = 1./SUMC
          TDRAT(1,JJD)=DENS
          TDRAT(2,JJD)=FRAT
          abund = sumc*iobs/100
                                                       !End of the Ne loop
        enddo
        do IA = 1, IAPR
          I1=ITRANA(1,IA)
          I2=ITRANA(2,IA)
          DEE=E(I2)-E(I1)
          WAVA(IA)=1.D8/DEE
        enddo
        do IB = 1, IBPR
          I1=ITRANB(1,IB)
          I2=ITRANB(2,IB)
          DEE=E(I2)-E(I1)
          WAVB(IB)=1.D8/DEE
        enddo
        do IC = 1, ICPR
          I1=ITRANC(1,IC)
          I2=ITRANC(2,IC)
          DEE=E(I2)-E(I1)
          WAVC(IC)=1.D8/DEE
        enddo
                                                         !End of the Te loop
      enddo

      return

 6100 FORMAT (' PROCESSING COMPLETED'/                                  &
     & ' GOODBYE!!'///)

end subroutine get_abundance

!---- PROC LUSLV
                                                       !Solving linear equations
      subroutine LUSLV(A,B,N,M)
      implicit none

      integer :: M, N
      real(kind=dp) :: A(M,M),B(M)
#ifdef DEBUG
        print *,"subroutine: luslv"
#endif
      call LURED(A,N,M)
      call RESLV(A,B,N,M)
      return

end subroutine luslv
!
!---- PROC LURED
      subroutine LURED(A,N,NR)
      implicit none

      integer :: N, NR, NM1, I, J, K, IP1
      real(kind=dp) :: A(NR,NR), FACT

#ifdef DEBUG
        print *,"subroutine: lured"
#endif

      if (N.EQ.1) return
      NM1=N-1
      do I=1,NM1
        IP1=I+1
        do K=IP1,N
          FACT=A(K,I)/A(I,I)
          do J=IP1,N
            A(K,J)=A(K,J)-A(I,J)*FACT
          enddo
        enddo
      enddo
      return

end subroutine lured
!
!---- PROC RESLV
                                                               !Resolve A with B
      subroutine RESLV(A,B,N,NR)
      implicit none

      integer :: N, NR, NM1, I, J, K, L, IP1
      real(kind=dp) :: A(NR,NR),B(NR)

#ifdef DEBUG
        print *,"subroutine: reslv"
#endif

      if (N.EQ.1) GOTO 1
      NM1=N-1
      do I=1,NM1
        IP1=I+1
        do J=IP1,N
          B(J)=B(J)-B(I)*A(J,I)/A(I,I)
        enddo
      enddo
      B(N)=B(N)/A(N,N)
      do I=1,NM1
        K=N-I
        L=K+1
        do J=L,N
          B(K)=B(K)-B(J)*A(K,J)
        enddo
        B(K)=B(K)/A(K,K)
      enddo
      return
    1 B(N)=B(N)/A(N,N)
      return

end subroutine reslv
!
!---- PROC SPLMAT
      subroutine SPLMAT(XX,NPT,IOPT,NDIM, NDIMT3, HMH)
      implicit none

      integer :: NDIM, NDIMT3, NPT, IOPT, NPM, NELEM
      real(kind=dp) :: XX(NDIM),GH(NDIMT3),Y(NDIM), HMH(NDIM,NDIM)

#ifdef DEBUG
        print *,"subroutine: splmat"
#endif

      NPM=NPT-2
      call GHGEN(GH,XX,NPT,IOPT,NDIM,NDIMT3)
      NELEM=3*NPM-2
      call ELU(GH,NPM,NDIMT3)
      call HGEN(XX,GH,Y,NPT,IOPT,NDIM,NDIMT3,HMH)
      return

end subroutine splmat
!
!---- PROC DERIV
!     Calculate the first derivative of the lagrangian interpolator
!     of a function F, tabulated at the N points XY(I), I=1 to N.
!     The derivative is given as the coefficients of F(I), I=1 to N,
!     in the array D(I), I=1 to N.
      subroutine DERIV(XY,D,X,N,NDIM)
      implicit none

      integer :: N ,NDIM, I, J, K
      real(kind=dp) :: XY(NDIM),D(NDIM), X, P1, P2, S

#ifdef DEBUG
        print *,"subroutine: deriv"
#endif

      do I=1,N
        P1=1.
        S=0.
        do J=1,N
          if (J.NE.I) then
            P1=P1*(XY(I)-XY(J))
            P2=1.
            do K=1,N
              if (K.NE.I.AND.K.NE.J) P2=P2*(X-XY(K))
            enddo
            S=S+P2
          endif
        enddo
        D(I)=S/P1
      enddo
      return

end subroutine deriv
!
!---- PROC HGEN
!     Cubic spline interpolation
!     The equation for the second derivatives at internal points
!     is of the form G*YPP=B, where G has been evaluated and LU
!     decomposed.
!     this routine writes B=HMH*Y and then solves YPP=G**(-1)*HMH*Y,
!     =HMH*Y.
!     Three options are provided for boundary conditions-
!     IOPT = 0  YPP=0 at end points
!     IOPT = 1  YP=0  at end points
!     IOPT = 2  YP at end points from lagarnge interpolant of a set of
!     internal points.
      subroutine HGEN(XX,GH,Y,NPT,IOPT,NDIM,NDIMT3,HMH)
      implicit none

      integer :: NPT, IOPT, NDIM, NDIMT3, NDIM3, NIP, I, J, K, NPM,        &
     & INDX
      real(kind=dp) :: XX(NDIM), GH(NDIMT3), Y(NDIM), HMH(NDIM,NDIM),             &
     & XY(5),D(5),C(2,5), A0, AN1, H1, H2

#ifdef DEBUG
        print *,"subroutine: hgen"
#endif

      if (IOPT.EQ.2) then          !Case of derivative boundary condition, with
        NDIM3=5                    !derivatives from NIP-point Lagrange at
        NIP=3                      !internal points
        do J=1,2
          do I=1,NIP
            K=(NPT-NIP)*(J-1)
            XY(I)=XX(K+I)
          enddo
          K=1+(NPT-1)*(J-1)
          call DERIV(XY,D,XX(K),NIP,NDIM3)
          do I=1,NIP
            C(J,I)=D(I)
          enddo
        enddo
      endif
                                             !Set up matrix equation G*YPP=HMH*Y
      A0=XX(2)-XX(1)
      AN1=XX(NPT)-XX(NPT-1)
      NPM=NPT-2
      hmh=0d0
      do I=1,NPM
        H1=6./(XX(I+1)-XX(I))
        H2=6./(XX(I+2)-XX(I+1))
        HMH(I,I)=H1
        HMH(I,I+1)=-H1-H2
        HMH(I,I+2)=H2
!        do J=1,NPT
!          HMH(I,J)=0.
!          if (J.EQ.I) HMH(I,J)=H1
!          if (J.EQ.I+2) HMH(I,J)=H2
!          if (J.EQ.I+1) HMH(I,J)=-H1-H2
!        enddo
      enddo
                                                 !Correct matrix for case of
      if (IOPT.EQ.1.OR.IOPT.EQ.2) then
                                                 !derivative boundary conditions
        HMH(1,1)=HMH(1,1)+3/A0
        HMH(1,2)=HMH(1,2)-3/A0
        HMH(NPM,NPT-1)=HMH(NPM,NPT-1)-3/AN1
        HMH(NPM,NPT)=HMH(NPM,NPT)+3/AN1
      endif
      if (IOPT.EQ.2) then
        do J=1,NIP
          HMH(1,J)=HMH(1,J)+3*C(1,J)
          K=NPT+J-NIP
          HMH(NPM,K)=HMH(NPM,K)-3*C(2,J)
        enddo
      endif
                                 !Solve matrix equation with results in the form
      do I=1,NPT
                                 !YPP=HMH*Y. matrix g has been LU decomposed
        Y(1)=HMH(1,I)
        INDX=0
        do J=2,NPM
          INDX=INDX+3
          Y(J)=HMH(J,I)-GH(INDX)*Y(J-1)
        enddo
        INDX=INDX+1
        Y(NPM)=Y(NPM)/GH(INDX)
        do J=2,NPM
          K=NPM-J+1
          INDX=INDX-3
          Y(K)=(Y(K)-GH(INDX+1)*Y(K+1))/GH(INDX)
        enddo

        HMH(2:npm+1,I)=Y(1:npm)
                                    !Insert values for second derivative at end
        HMH(1,I)=0.
                                    !points: first and last rows of the matrix
        HMH(NPT,I)=0.
      enddo
                                         !Case of derivative boundary conditions
      if (IOPT.GT.0) then
        do J=1,NPT
          HMH(1,J)=-0.5*HMH(2,J)
          HMH(NPT,J)=-0.5*HMH(NPT-1,J)
        enddo
        HMH(1,1)=HMH(1,1)-3/(A0*A0)
        HMH(1,2)=HMH(1,2)+3/(A0*A0)
        HMH(NPT,NPT-1)=HMH(NPT,NPT-1)+3/(AN1*AN1)
        HMH(NPT,NPT)=HMH(NPT,NPT)-3/(AN1*AN1)
      endif
      if (IOPT.EQ.2) then
        do J=1,NIP
          HMH(1,J)=HMH(1,J)-3*C(1,J)/A0
          K=NPT+J-NIP
          HMH(NPT,K)=HMH(NPT,K)+3*C(2,J)/AN1
        enddo
      endif
      return
end subroutine hgen
!
!---- PROC GHGEN
      subroutine GHGEN(GH,XX,NPT,IOPT,NDIM,NDIMT3)
      implicit none

      integer :: NPT, IOPT, NDIM, NDIMT3, INDX, NPTM, I, J, IP, JP, IK
      real(kind=dp) :: XX(NDIM),GH(NDIMT3)

#ifdef DEBUG
        print *,"subroutine: ghgen"
#endif

      INDX=0
      NPTM=NPT-1
      do I=2,NPTM
        IP=I-1
        do J=1,3
          JP=IP+J-2
          if (JP.GE.1.AND.JP.LE.NPTM-1) then
            INDX=INDX+1
            if (J.EQ.2) then
              GH(INDX)=2*(XX(I+1)-XX(I-1))
            else
              IK=I+(J-1)/2
              GH(INDX)=XX(IK)-XX(IK-1)
            endif
          endif
        enddo
      enddo
      if (IOPT.GE.1) then
        GH(1)=GH(1)-(XX(2)-XX(1))/2.
        GH(INDX)=GH(INDX)-(XX(NPT)-XX(NPT-1))/2.
      endif
      return
end subroutine ghgen
!
!---- PROC ELU
      subroutine ELU(GH,N,NDIM)
      implicit none

      integer :: N, NDIM, INDX, I, J, JP
      real(kind=dp) :: GH(NDIM)

#ifdef DEBUG
        print *,"subroutine: elu"
#endif

      INDX=0
      do I=1,N
        do J=1,3
          JP=I+J-2
          if (JP.GE.1.AND.JP.LE.N) then
            INDX=INDX+1
            if (I.GT.1) then
              if (J.EQ.1) then
                GH(INDX)=GH(INDX)/GH(INDX-2)
              endif
              if (J.EQ.2) then
                GH(INDX)=GH(INDX)-GH(INDX-1)*GH(INDX-2)
              endif
            endif
          endif
        enddo
      enddo
      return
end subroutine elu
!
!---- PROC CFY
      subroutine CFY(X,Y,XX,YY,NPT,NDIM,D)
      implicit none

      integer :: NPT, NDIM, J
      real(kind=dp) :: XX(NDIM),YY(NDIM), D(NDIM), X, Y, TT

#ifdef DEBUG
        print *,"subroutine: cfy"
#endif

      if (X.LT.XX(1)) then
        Y=YY(1)
      endif
      if (X.GT.XX(NPT)) then
        Y=YY(NPT)
      endif
      TT=0.
      do J=1,NPT
        TT=TT+D(J)*YY(J)
      enddo
      Y=TT
      return

end subroutine cfy
!
!---- PROC CFD
      subroutine CFD(X,XX,NPT,NDIM, HMH, D)
      implicit none

      integer :: NPT, NDIM, NPTM, I, J
      real(kind=dp) :: X, XX(NDIM), HMH(NDIM,NDIM), D(NDIM), X1, X2, A1, A2, HI

#ifdef DEBUG
        print *,"subroutine: cfd"
#endif

      if (X.LT.XX(1)) then
        !write(6,400) XX(1)
        return
      endif
      if (X.GT.XX(NPT)) then
        !write(6,401) XX(NPT)
        return
      endif
      NPTM=NPT-1
      do I=1,NPTM
        if (X.LT.XX(I+1)) then
          X1=XX(I+1)-X
          X2=X-XX(I)
          HI=XX(I+1)-XX(I)
          A1=X1*(X1*X1/(6*HI)-HI/6)
          A2=X2*(X2*X2/(6*HI)-HI/6)
          do J=1,NPT
            D(J)=(A1*HMH(I,J)+A2*HMH(I+1,J))
          enddo
          D(I)=D(I)+X1/HI
          D(I+1)=D(I+1)+X2/HI
          return
        endif
      enddo

end subroutine cfd

end module mod_equib