-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathampld.lp.f
2074 lines (1904 loc) · 67.3 KB
/
ampld.lp.f
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
C This is NOT the original code by the authors listed below!
C Modified 2009-2013 by Jussi Leinonen ([email protected])
C to be compatible with the Python extension interface.
C The requirement of non-profit use only has been dropped from
C this code with the permission of M. I. Mishchenko; thus the
C license has been made open source compatible.
C New release including the LAPACK matrix inversion procedure.
C We thank Cory Davis (University of Edinburgh) for pointing
C out the possibility of replacing the proprietary NAG matrix
C inversion routine by the public-domain LAPACK equivalent.
C CALCULATION OF THE AMPLITUDE AND PHASE MATRICES FOR
C A PARTICLE WITH AN AXIALLY SYMMETRIC SHAPE
C This version of the code uses DOUBLE PRECISION variables,
C is applicable to spheroids, finite circular cylinders,
C Chebyshev particles, and generalized Chebyshev particles
C (distorted water drops), and must be used along with the
C accompanying files lpd.f and ampld.par.f.
C Last update 08/06/2005
C The code has been developed by Michael Mishchenko at the NASA
C Goddard Institute for Space Studies, New York. The development
C of the code was supported by the NASA Radiation Sciences Program.
C We only request that in any publication using the code the source of
C the code be acknowledged and relevant references (see below) be
C made.
C The computational method is based on the T-matrix approach
C [P. C. Waterman, Phys. Rev. D 3, 825 (1971)], also known as
C the extended boundary condition method. The method was
C improved in the following papers:
C 1. M. I. Mishchenko and L. D. Travis, T-matrix computations
C of light scattering by large spheroidal particles,
C Opt. Commun., vol. 109, 16-21 (1994).
C
C 2. M. I. Mishchenko, L. D. Travis, and A. Macke, Scattering
C of light by polydisperse, randomly oriented, finite
C circular cylinders, Appl. Opt., vol. 35, 4927-4940 (1996).
C
C 3. D. J. Wielaard, M. I. Mishchenko, A. Macke, and B. E. Carlson,
C Improved T-matrix computations for large, nonabsorbing and
C weakly absorbing nonspherical particles and comparison
C with geometrical optics approximation, Appl. Opt., vol. 36,
C 4305-4313 (1997).
C
C A general review of the T-matrix approach can be found in
C
C 4. M. I. Mishchenko, L. D. Travis, and D. W. Mackowski,
C T-matrix computations of light scattering by nonspherical
C particles: a review, J. Quant. Spectrosc. Radiat.
C Transfer, vol. 55, 535-575 (1996).
C
C Additional useful information is contained in the paper
C
C 5. M. I. Mishchenko and L. D. Travis, Capabilities and
C limitations of a current FORTRAN implementation of the
C T-matrix method for randomly oriented, rotationally
C symmetric scatterers, J. Quant. Spectrosc. Radiat. Transfer,
C vol. 60, 309-324 (1998).
C
C The definitions and notation used can also be found in
C
C 6. M. I. Mishchenko, Calculation of the amplitude matrix
C for a nonspherical particle in a fixed orientation,
C Appl. Opt. vol. 39, 1026-1031 (2000).
C These papers are available in the .pdf format at the web site
C
C http://www.giss.nasa.gov/~crmim/publications/
C
C or in hardcopy upon request from Michael Mishchenko
C Please e-mail your request to [email protected].
C
C A comprehensive book "Scattering, Absorption, and Emission of
C Light by Small Particles" (Cambridge University Press, Cambridge,
C 2002) is also available in the .pdf format at the web site
C
C http://www.giss.nasa.gov/~crmim/books.html
C One of the main advantages of the T-matrix method is that the
C T-matrix for a given nonspherical particle needs to be computed
C only once and then can be used any number of times for computing
C the amplitude and phase matrices for any directions of the incident
C and scattered beams. This makes the T-matrix method
C extremely convenient and efficient in averaging over particle
C orientations and/or directions of incidence (or scattering).
C The use of extended precision variables (Ref. 1) can significantly
C increase the maximum convergent equivalent-sphere size parameter
C and make it well larger than 100 (depending on refractive index
C and aspect ratio). The extended-precision code is also available.
C However, the use of extended precision varibales results in a
C greater consumption of CPU time.
C On IBM RISC workstations, that code is approximately
C five times slower than this double-precision code. The
C CPU time difference between the double-precision and extended-
C precision codes can be larger on supercomputers.
C This is the first part of the full T-matrix code. The second part,
C lpq.f, is completely independent of the first part. It contains no
C T-matrix-specific subroutines and can be compiled separately.
C The second part of the code replaces the previously implemented
C standard matrix inversion scheme based on Gaussian elimination
C by a scheme based on the LU factorization technique.
C As described in Ref. 3 above, the use of the LU factorization is
C especially beneficial for nonabsorbing or weakly absorbing particles.
C In this code we use the LAPACK implementation of the LU factorization
C scheme. LAPACK stands for Linear Algebra PACKage. The latter is
C publicly available at the following internet site:
C
C http://www.netlib.org/lapack/
C INPUT PARAMETERS:
C
C AXI - equivalent-sphere radius
C RAT = 1 - particle size is specified in terms of the
C equal-volume-sphere radius
C RAT.ne.1 - particle size is specified in terms of the
C equal-surface-area-sphere radius
C LAM - WAVELENGTH OF INCIDENT LIGHT
C MRR and MRI - real and imaginary parts of the refractive
C index
C EPS and NP - specify the shape of the particles.
C For spheroids NP=-1 and EPS is the ratio of the
C horizontal to rotational axes. EPS is larger than
C 1 for oblate spheroids and smaller than 1 for
C prolate spheroids.
C For cylinders NP=-2 and EPS is the ratio of the
C diameter to the length.
C For Chebyshev particles NP must be positive and
C is the degree of the Chebyshev polynomial, while
C EPS is the deformation parameter (Ref. 5).
C For generalized Chebyshev particles (describing the shape
C of distorted water drops) NP=-3. The coefficients
C of the Chebyshev polynomial expansion of the particle
C shape (Ref. 7) are specified in subroutine DROP.
C DDELT - accuracy of the computations
C NDGS - parameter controlling the number of division points
C in computing integrals over the particle surface (Ref. 5).
C For compact particles, the recommended value is 2.
C For highly aspherical particles larger values (3, 4,...)
C may be necessary to obtain convergence.
C The code does not check convergence over this parameter.
C Therefore, control comparisons of results obtained with
C different NDGS-values are recommended.
C ALPHA and BETA - Euler angles (in degrees) specifying the orientation
C of the scattering particle relative to the laboratory reference
C frame (Refs. 6 and 7).
C THET0 - zenith angle of the incident beam in degrees
C THET - zenith angle of the scattered beam in degrees
C PHI0 - azimuth angle of the incident beam in degrees
C PHI - azimuth angle of the scattered beam in degrees
C (Refs. 6 and 7)
C ALPHA, BETA, THET0, THET, PHI0, and PHI are specified at the end of
C the main program before the line
C
C "CALL AMPL (NMAX,...)"
C
C The part of the main program following the line
C
C "COMPUTATION OF THE AMPLITUDE AND PHASE MATRICES"
C
C can be repeated any number of times. At this point the T-matrix
C for the given scattering particle has already
C been fully computed and can be repeatedly used in computations
C for any directions of illumination and scattering and any particle
C orientations.
C OUTPUT PARAMETERS:
C
C Elements of the 2x2 amplitude matrix
C Elements of the 4x4 phase matrix
C Note that LAM and AXI must be given in the same units of length
C (e.g., microns).
C The convergence of the T-matrix method for particles with
C different sizes, refractive indices, and aspect ratios can be
C dramatically different. Usually, large sizes and large aspect
C ratios cause problems. The user of this code
C should "play" a little with different input parameters in
C order to get an idea of the range of applicability of this
C technique. Sometimes decreasing the aspect ratio
C from 3 to 2 can increase the maximum convergent equivalent-
C sphere size parameter by a factor of several.
C The CPU time required rapidly increases with increasing ratio
C radius/wavelength and/or with increasing particle asphericity.
C This should be taken into account in planning massive computations.
C Execution can be automatically terminated if dimensions of certain
C arrays are not big enough or if the convergence procedure decides
C that the accuracy of double-precision variables is insufficient
C to obtain a converged T-matrix solution for given particles.
C In all cases, a message appears explaining
C the cause of termination.
C The message
C "WARNING: W IS GREATER THAN 1"
C means that the single-scattering albedo exceeds the maximum
C possible value 1. If W is greater than 1 by more than
C DDELT, this message can be an indication of numerical
C instability caused by extreme values of particle parameters.
C The message "WARNING: NGAUSS=NPNG1" means that convergence over
C the parameter NG (see Ref. 2) cannot be obtained for the NPNG1
C value specified in the PARAMETER statement in the file ampld.par.f.
C Often this is not a serious problem, especially for compact
C particles.
C Larger and/or more aspherical particles may require larger
C values of the parameters NPN1, NPN4, and NPNG1 in the file
C ampld.par.f. It is recommended to keep NPN1=NPN4+25 and
C NPNG1=3*NPN1. Note that the memory requirement increases
C as the third power of NPN4. If the memory of
C a computer is too small to accomodate the code in its current
C setting, the parameters NPN1, NPN4, and NPNG1 should be
C decreased. However, this will decrease the maximum size parameter
C that can be handled by the code.
C In some cases any increases of NPN1 will not make the T-matrix
C computations convergent. This means that the particle is just
C too "bad" (extreme size parameter and/or extreme aspect ratio
C and/or extreme refractive index).
C The main program contains several PRINT statements which are
C currently commentd out. If uncommented, these statements will
C produce numbers which show the convergence rate and can be
C used to determine whether T-matrix computations for given particle
C parameters will converge at all, whatever the parameter NPN1 is.
C Some of the common blocks are used to save memory rather than
C to transfer data. Therefore, if a compiler produces a warning
C message that the lengths of a common block are different in
C different subroutines, this is not a real problem.
C The recommended value of DDELT is 0.001. For bigger values,
C false convergence can be obtained.
C In computations for spheres, use EPS=1.000001 instead of EPS=1.
C EPS=1 can cause overflows in some rare cases.
C For some compilers, DACOS must be raplaced by DARCOS and DASIN
C by DARSIN.
C I would highly appreciate informing me of any problems encountered
C with this code. Please send your message to the following
C e-mail address: [email protected].
C WHILE THE COMPUTER PROGRAM HAS BEEN TESTED FOR A VARIETY OF CASES,
C IT IS NOT INCONCEIVABLE THAT IT CONTAINS UNDETECTED ERRORS. ALSO,
C INPUT PARAMETERS CAN BE USED WHICH ARE OUTSIDE THE ENVELOPE OF
C VALUES FOR WHICH RESULTS ARE COMPUTED ACCURATELY. FOR THIS REASON,
C THE AUTHORS AND THEIR ORGANIZATION DISCLAIM ALL LIABILITY FOR
C ANY DAMAGES THAT MAY RESULT FROM THE USE OF THE PROGRAM.
SUBROUTINE CALCTMAT(AXI,RAT,LAM,MRR,MRI,EPS,NP,DDELT,NDGS,NMAX)
IMPLICIT REAL*8 (A-H,O-Z)
Cf2py intent(in) axi
Cf2py intent(in) rat
Cf2py intent(in) lam
Cf2py intent(in) mrr
Cf2py intent(in) mri
Cf2py intent(in) eps
Cf2py intent(in) np
Cf2py intent(in) ddelt
Cf2py intent(in) ndgs
Cf2py intent(out) nmax
INCLUDE 'ampld.par.f'
REAL*8 LAM,MRR,MRI,X(NPNG2),W(NPNG2),S(NPNG2),SS(NPNG2),
* AN(NPN1),R(NPNG2),DR(NPNG2),
* DDR(NPNG2),DRR(NPNG2),DRI(NPNG2),ANN(NPN1,NPN1)
REAL*8 TR1(NPN2,NPN2),TI1(NPN2,NPN2)
REAL*8 XALPHA(300),XBETA(300),WALPHA(300),WBETA(300)
REAL*4
& RT11(NPN6,NPN4,NPN4),RT12(NPN6,NPN4,NPN4),
& RT21(NPN6,NPN4,NPN4),RT22(NPN6,NPN4,NPN4),
& IT11(NPN6,NPN4,NPN4),IT12(NPN6,NPN4,NPN4),
& IT21(NPN6,NPN4,NPN4),IT22(NPN6,NPN4,NPN4)
COMMON /CT/ TR1,TI1
COMMON /TMAT/ RT11,RT12,RT21,RT22,IT11,IT12,IT21,IT22
COMMON /CHOICE/ ICHOICE
C OPEN FILES *******************************************************
C OPEN (6,FILE='test')
C INPUT DATA ********************************************************
C AXI=10D0
C RAT=0.1 D0
C LAM=DACOS(-1D0)*2D0
C MRR=1.5 D0
C MRI=0.02 D0
C EPS=0.5 D0
C NP=-1
C DDELT=0.001D0
C NDGS=2
P=DACOS(-1D0)
NCHECK=0
IF (NP.EQ.-1.OR.NP.EQ.-2) NCHECK=1
IF (NP.GT.0.AND.(-1)**NP.EQ.1) NCHECK=1
C WRITE (6,5454) NCHECK
C 5454 FORMAT ('NCHECK=',I1)
IF (DABS(RAT-1D0).GT.1D-8.AND.NP.EQ.-1) CALL SAREA (EPS,RAT)
IF (DABS(RAT-1D0).GT.1D-8.AND.NP.GE.0) CALL SURFCH(NP,EPS,RAT)
IF (DABS(RAT-1D0).GT.1D-8.AND.NP.EQ.-2) CALL SAREAC (EPS,RAT)
IF (NP.EQ.-3) CALL DROP (RAT)
C PRINT 8000, RAT
C 8000 FORMAT ('RAT=',F8.6)
C IF(NP.EQ.-1.AND.EPS.GE.1D0) PRINT 7000,EPS
C IF(NP.EQ.-1.AND.EPS.LT.1D0) PRINT 7001,EPS
C IF(NP.GE.0) PRINT 7100,NP,EPS
C IF(NP.EQ.-2.AND.EPS.GE.1D0) PRINT 7150,EPS
C IF(NP.EQ.-2.AND.EPS.LT.1D0) PRINT 7151,EPS
C IF(NP.EQ.-3) PRINT 7160
C PRINT 7400, LAM,MRR,MRI
C PRINT 7200,DDELT
C 7000 FORMAT('OBLATE SPHEROIDS, A/B=',F11.7)
C 7001 FORMAT('PROLATE SPHEROIDS, A/B=',F11.7)
C 7100 FORMAT('CHEBYSHEV PARTICLES, T',
C & I1,'(',F5.2,')')
C 7150 FORMAT('OBLATE CYLINDERS, D/L=',F11.7)
C 7151 FORMAT('PROLATE CYLINDERS, D/L=',F11.7)
C 7160 FORMAT('GENERALIZED CHEBYSHEV PARTICLES')
C 7200 FORMAT ('ACCURACY OF COMPUTATIONS DDELT = ',D8.2)
C 7400 FORMAT('LAM=',F10.6,3X,'MRR=',D10.4,3X,'MRI=',D10.4)
C DDELT=0.1D0*DDELT
C IF (DABS(RAT-1D0).LE.1D-6) PRINT 8003, AXI
C IF (DABS(RAT-1D0).GT.1D-6) PRINT 8004, AXI
C 8003 FORMAT('EQUAL-VOLUME-SPHERE RADIUS=',F8.4)
C 8004 FORMAT('EQUAL-SURFACE-AREA-SPHERE RADIUS=',F8.4)
A=RAT*AXI
XEV=2D0*P*A/LAM
IXXX=XEV+4.05D0*XEV**0.333333D0
INM1=MAX0(4,IXXX)
C IF (INM1.GE.NPN1) PRINT 7333, NPN1
IF (INM1.GE.NPN1) STOP
C 7333 FORMAT('CONVERGENCE IS NOT OBTAINED FOR NPN1=',I3,
C & '. EXECUTION TERMINATED')
QEXT1=0D0
QSCA1=0D0
DO 50 NMA=INM1,NPN1
NMAX=NMA
MMAX=1
NGAUSS=NMAX*NDGS
C IF (NGAUSS.GT.NPNG1) PRINT 7340, NGAUSS
IF (NGAUSS.GT.NPNG1) STOP
C 7340 FORMAT('NGAUSS =',I3,' I.E. IS GREATER THAN NPNG1.',
C & ' EXECUTION TERMINATED')
C 7334 FORMAT(' NMAX =', I3,' DC2=',D8.2,' DC1=',D8.2)
CALL CONST(NGAUSS,NMAX,MMAX,P,X,W,AN,ANN,S,SS,NP,EPS)
CALL VARY(LAM,MRR,MRI,A,EPS,NP,NGAUSS,X,P,PPI,PIR,PII,R,
& DR,DDR,DRR,DRI,NMAX)
CALL TMATR0 (NGAUSS,X,W,AN,ANN,S,SS,PPI,PIR,PII,R,DR,
& DDR,DRR,DRI,NMAX,NCHECK)
QEXT=0D0
QSCA=0D0
DO 4 N=1,NMAX
N1=N+NMAX
TR1NN=TR1(N,N)
TI1NN=TI1(N,N)
TR1NN1=TR1(N1,N1)
TI1NN1=TI1(N1,N1)
DN1=DFLOAT(2*N+1)
QSCA=QSCA+DN1*(TR1NN*TR1NN+TI1NN*TI1NN
& +TR1NN1*TR1NN1+TI1NN1*TI1NN1)
QEXT=QEXT+(TR1NN+TR1NN1)*DN1
4 CONTINUE
DSCA=DABS((QSCA1-QSCA)/QSCA)
DEXT=DABS((QEXT1-QEXT)/QEXT)
QEXT1=QEXT
QSCA1=QSCA
C PRINT 7334, NMAX,DSCA,DEXT
IF(DSCA.LE.DDELT.AND.DEXT.LE.DDELT) GO TO 55
C IF (NMA.EQ.NPN1) PRINT 7333, NPN1
IF (NMA.EQ.NPN1) STOP
50 CONTINUE
55 NNNGGG=NGAUSS+1
MMAX=NMAX
C IF (NGAUSS.EQ.NPNG1) PRINT 7336
IF (NGAUSS.EQ.NPNG1) GO TO 155
DO 150 NGAUS=NNNGGG,NPNG1
NGAUSS=NGAUS
NGGG=2*NGAUSS
C 7336 FORMAT('WARNING: NGAUSS=NPNG1')
C 7337 FORMAT(' NG=',I3,' DC2=',D8.2,' DC1=',D8.2)
CALL CONST(NGAUSS,NMAX,MMAX,P,X,W,AN,ANN,S,SS,NP,EPS)
CALL VARY(LAM,MRR,MRI,A,EPS,NP,NGAUSS,X,P,PPI,PIR,PII,R,
& DR,DDR,DRR,DRI,NMAX)
CALL TMATR0 (NGAUSS,X,W,AN,ANN,S,SS,PPI,PIR,PII,R,DR,
& DDR,DRR,DRI,NMAX,NCHECK)
QEXT=0D0
QSCA=0D0
DO 104 N=1,NMAX
N1=N+NMAX
TR1NN=TR1(N,N)
TI1NN=TI1(N,N)
TR1NN1=TR1(N1,N1)
TI1NN1=TI1(N1,N1)
DN1=DFLOAT(2*N+1)
QSCA=QSCA+DN1*(TR1NN*TR1NN+TI1NN*TI1NN
& +TR1NN1*TR1NN1+TI1NN1*TI1NN1)
QEXT=QEXT+(TR1NN+TR1NN1)*DN1
104 CONTINUE
DSCA=DABS((QSCA1-QSCA)/QSCA)
DEXT=DABS((QEXT1-QEXT)/QEXT)
C PRINT 7337, NGGG,DSCA,DEXT
QEXT1=QEXT
QSCA1=QSCA
IF(DSCA.LE.DDELT.AND.DEXT.LE.DDELT) GO TO 155
C IF (NGAUS.EQ.NPNG1) PRINT 7336
150 CONTINUE
155 CONTINUE
QSCA=0D0
QEXT=0D0
NNM=NMAX*2
DO 204 N=1,NNM
QEXT=QEXT+TR1(N,N)
204 CONTINUE
DO 213 N2=1,NMAX
NN2=N2+NMAX
DO 213 N1=1,NMAX
NN1=N1+NMAX
ZZ1=TR1(N1,N2)
RT11(1,N1,N2)=ZZ1
ZZ2=TI1(N1,N2)
IT11(1,N1,N2)=ZZ2
ZZ3=TR1(N1,NN2)
RT12(1,N1,N2)=ZZ3
ZZ4=TI1(N1,NN2)
IT12(1,N1,N2)=ZZ4
ZZ5=TR1(NN1,N2)
RT21(1,N1,N2)=ZZ5
ZZ6=TI1(NN1,N2)
IT21(1,N1,N2)=ZZ6
ZZ7=TR1(NN1,NN2)
RT22(1,N1,N2)=ZZ7
ZZ8=TI1(NN1,NN2)
IT22(1,N1,N2)=ZZ8
QSCA=QSCA+ZZ1*ZZ1+ZZ2*ZZ2+ZZ3*ZZ3+ZZ4*ZZ4
& +ZZ5*ZZ5+ZZ6*ZZ6+ZZ7*ZZ7+ZZ8*ZZ8
213 CONTINUE
C PRINT 7800,0,DABS(QEXT),QSCA,NMAX
DO 220 M=1,NMAX
CALL TMATR(M,NGAUSS,X,W,AN,ANN,S,SS,PPI,PIR,PII,R,DR,
& DDR,DRR,DRI,NMAX,NCHECK)
NM=NMAX-M+1
M1=M+1
QSC=0D0
DO 214 N2=1,NM
NN2=N2+M-1
N22=N2+NM
DO 214 N1=1,NM
NN1=N1+M-1
N11=N1+NM
ZZ1=TR1(N1,N2)
RT11(M1,NN1,NN2)=ZZ1
ZZ2=TI1(N1,N2)
IT11(M1,NN1,NN2)=ZZ2
ZZ3=TR1(N1,N22)
RT12(M1,NN1,NN2)=ZZ3
ZZ4=TI1(N1,N22)
IT12(M1,NN1,NN2)=ZZ4
ZZ5=TR1(N11,N2)
RT21(M1,NN1,NN2)=ZZ5
ZZ6=TI1(N11,N2)
IT21(M1,NN1,NN2)=ZZ6
ZZ7=TR1(N11,N22)
RT22(M1,NN1,NN2)=ZZ7
ZZ8=TI1(N11,N22)
IT22(M1,NN1,NN2)=ZZ8
QSC=QSC+(ZZ1*ZZ1+ZZ2*ZZ2+ZZ3*ZZ3+ZZ4*ZZ4
& +ZZ5*ZZ5+ZZ6*ZZ6+ZZ7*ZZ7+ZZ8*ZZ8)*2D0
214 CONTINUE
NNM=2*NM
QXT=0D0
DO 215 N=1,NNM
QXT=QXT+TR1(N,N)*2D0
215 CONTINUE
QSCA=QSCA+QSC
QEXT=QEXT+QXT
C PRINT 7800,M,DABS(QXT),QSC,NMAX
C 7800 FORMAT(' m=',I3,' qxt=',D12.6,' qsc=',D12.6,
C & ' nmax=',I3)
220 CONTINUE
WALB=-QSCA/QEXT
C IF (WALB.GT.1D0+DDELT) PRINT 9111
C 9111 FORMAT ('WARNING: W IS GREATER THAN 1')
RETURN
END
C COMPUTATION OF THE AMPLITUDE AND PHASE MATRICES
C ALPHA=145D0
C BETA=52D0
C THET0=56D0
C THET=65D0
C PHI0=114D0
C PHI=128D0
SUBROUTINE CALCAMPL(NMAX,LAM,THET0,THET,PHI0,PHI,ALPHA,BETA,S,Z)
IMPLICIT REAL*8 (A-H,O-Z)
REAL*8 LAM
COMPLEX*16 S(2,2)
REAL*8 Z(4,4)
Cf2py intent(in) nmax
Cf2py intent(in) lam
Cf2py intent(in) thet0
Cf2py intent(in) thet
Cf2py intent(in) phi0
Cf2py intent(in) phi
Cf2py intent(in) alpha
Cf2py intent(in) beta
Cf2py intent(out) s
Cf2py intent(out) z
C AMPLITUDE MATRIX [Eqs. (2)-(4) of Ref. 6]
CALL AMPL (NMAX,LAM,THET0,THET,PHI0,PHI,ALPHA,BETA,
& S(1,1),S(1,2),S(2,1),S(2,2))
C PHASE MATRIX [Eqs. (13)-(29) of Ref. 6]
Z(1,1)=0.5D0*(S(1,1)*DCONJG(S(1,1))+S(1,2)*DCONJG(S(1,2))
& +S(2,1)*DCONJG(S(2,1))+S(2,2)*DCONJG(S(2,2)))
Z(1,2)=0.5D0*(S(1,1)*DCONJG(S(1,1))-S(1,2)*DCONJG(S(1,2))
& +S(2,1)*DCONJG(S(2,1))-S(2,2)*DCONJG(S(2,2)))
Z(1,3)=-S(1,1)*DCONJG(S(1,2))-S(2,2)*DCONJG(S(2,1))
Z(1,4)=(0D0,1D0)*(S(1,1)*DCONJG(S(1,2))-S(2,2)*DCONJG(S(2,1)))
Z(2,1)=0.5D0*(S(1,1)*DCONJG(S(1,1))+S(1,2)*DCONJG(S(1,2))
& -S(2,1)*DCONJG(S(2,1))-S(2,2)*DCONJG(S(2,2)))
Z(2,2)=0.5D0*(S(1,1)*DCONJG(S(1,1))-S(1,2)*DCONJG(S(1,2))
& -S(2,1)*DCONJG(S(2,1))+S(2,2)*DCONJG(S(2,2)))
Z(2,3)=-S(1,1)*DCONJG(S(1,2))+S(2,2)*DCONJG(S(2,1))
Z(2,4)=(0D0,1D0)*(S(1,1)*DCONJG(S(1,2))+S(2,2)*DCONJG(S(2,1)))
Z(3,1)=-S(1,1)*DCONJG(S(2,1))-S(2,2)*DCONJG(S(1,2))
Z(3,2)=-S(1,1)*DCONJG(S(2,1))+S(2,2)*DCONJG(S(1,2))
Z(3,3)=S(1,1)*DCONJG(S(2,2))+S(1,2)*DCONJG(S(2,1))
Z(3,4)=(0D0,-1D0)*(S(1,1)*DCONJG(S(2,2))+S(2,1)*DCONJG(S(1,2)))
Z(4,1)=(0D0,1D0)*(S(2,1)*DCONJG(S(1,1))+S(2,2)*DCONJG(S(1,2)))
Z(4,2)=(0D0,1D0)*(S(2,1)*DCONJG(S(1,1))-S(2,2)*DCONJG(S(1,2)))
Z(4,3)=(0D0,-1D0)*(S(2,2)*DCONJG(S(1,1))-S(1,2)*DCONJG(S(2,1)))
Z(4,4)=S(2,2)*DCONJG(S(1,1))-S(1,2)*DCONJG(S(2,1))
C WRITE (6,5000)
C WRITE (6,5001) Z11,Z12,Z13,Z14
C WRITE (6,5001) Z21,Z22,Z23,Z24
C WRITE (6,5001) Z31,Z32,Z33,Z34
C WRITE (6,5001) Z41,Z42,Z43,Z44
C 5000 FORMAT ('PHASE MATRIX')
C 5001 FORMAT (4F10.4)
C ITIME=MCLOCK()
C TIME=DFLOAT(ITIME)/6000D0
C PRINT 1001,TIME
C 1001 FORMAT (' time =',F8.2,' min')
C STOP
RETURN
END
C********************************************************************
C CALCULATION OF THE AMPLITUDE MATRIX
SUBROUTINE AMPL (NMAX,DLAM,TL,TL1,PL,PL1,ALPHA,BETA,
& VV,VH,HV,HH)
INCLUDE 'ampld.par.f'
IMPLICIT REAL*8 (A-B,D-H,O-Z), COMPLEX*16 (C)
REAL*8 AL(3,2),AL1(3,2),AP(2,3),AP1(2,3),B(3,3),
* R(2,2),R1(2,2),C(3,2),CA,CB,CT,CP,CTP,CPP,CT1,CP1,
* CTP1,CPP1
REAL*8 DV1(NPN6),DV2(NPN6),DV01(NPN6),DV02(NPN6)
REAL*4
& TR11(NPN6,NPN4,NPN4),TR12(NPN6,NPN4,NPN4),
& TR21(NPN6,NPN4,NPN4),TR22(NPN6,NPN4,NPN4),
& TI11(NPN6,NPN4,NPN4),TI12(NPN6,NPN4,NPN4),
& TI21(NPN6,NPN4,NPN4),TI22(NPN6,NPN4,NPN4)
COMPLEX*16 CAL(NPN4,NPN4),VV,VH,HV,HH
COMMON /TMAT/ TR11,TR12,TR21,TR22,TI11,TI12,TI21,TI22
IF (ALPHA.LT.0D0.OR.ALPHA.GT.360D0.OR.
& BETA.LT.0D0.OR.BETA.GT.180D0.OR.
& TL.LT.0D0.OR.TL.GT.180D0.OR.
& TL1.LT.0D0.OR.TL1.GT.180D0.OR.
& PL.LT.0D0.OR.PL.GT.360D0.OR.
& PL1.LT.0D0.OR.PL1.GT.360D0) THEN
C WRITE (6,2000)
STOP
ELSE
CONTINUE
ENDIF
C 2000 FORMAT ('AN ANGULAR PARAMETER IS OUTSIDE ITS',
C & ' ALLOWABLE RANGE')
PIN=DACOS(-1D0)
PIN2=PIN*0.5D0
PI=PIN/180D0
ALPH=ALPHA*PI
BET=BETA*PI
THETL=TL*PI
PHIL=PL*PI
THETL1=TL1*PI
PHIL1=PL1*PI
EPS=1D-7
IF (THETL.LT.PIN2) THETL=THETL+EPS
IF (THETL.GT.PIN2) THETL=THETL-EPS
IF (THETL1.LT.PIN2) THETL1=THETL1+EPS
IF (THETL1.GT.PIN2) THETL1=THETL1-EPS
IF (PHIL.LT.PIN) PHIL=PHIL+EPS
IF (PHIL.GT.PIN) PHIL=PHIL-EPS
IF (PHIL1.LT.PIN) PHIL1=PHIL1+EPS
IF (PHIL1.GT.PIN) PHIL1=PHIL1-EPS
IF (BET.LE.PIN2.AND.PIN2-BET.LE.EPS) BET=BET-EPS
IF (BET.GT.PIN2.AND.BET-PIN2.LE.EPS) BET=BET+EPS
C_____________COMPUTE THETP, PHIP, THETP1, AND PHIP1, EQS. (8), (19), AND (20)
CB=DCOS(BET)
SB=DSIN(BET)
CT=DCOS(THETL)
ST=DSIN(THETL)
CP=DCOS(PHIL-ALPH)
SP=DSIN(PHIL-ALPH)
CTP=CT*CB+ST*SB*CP
THETP=DACOS(CTP)
CPP=CB*ST*CP-SB*CT
SPP=ST*SP
PHIP=DATAN(SPP/CPP)
IF (PHIP.GT.0D0.AND.SP.LT.0D0) PHIP=PHIP+PIN
IF (PHIP.LT.0D0.AND.SP.GT.0D0) PHIP=PHIP+PIN
IF (PHIP.LT.0D0) PHIP=PHIP+2D0*PIN
CT1=DCOS(THETL1)
ST1=DSIN(THETL1)
CP1=DCOS(PHIL1-ALPH)
SP1=DSIN(PHIL1-ALPH)
CTP1=CT1*CB+ST1*SB*CP1
THETP1=DACOS(CTP1)
CPP1=CB*ST1*CP1-SB*CT1
SPP1=ST1*SP1
PHIP1=DATAN(SPP1/CPP1)
IF (PHIP1.GT.0D0.AND.SP1.LT.0D0) PHIP1=PHIP1+PIN
IF (PHIP1.LT.0D0.AND.SP1.GT.0D0) PHIP1=PHIP1+PIN
IF (PHIP1.LT.0D0) PHIP1=PHIP1+2D0*PIN
C____________COMPUTE MATRIX BETA, EQ. (21)
CA=DCOS(ALPH)
SA=DSIN(ALPH)
B(1,1)=CA*CB
B(1,2)=SA*CB
B(1,3)=-SB
B(2,1)=-SA
B(2,2)=CA
B(2,3)=0D0
B(3,1)=CA*SB
B(3,2)=SA*SB
B(3,3)=CB
C____________COMPUTE MATRICES AL AND AL1, EQ. (14)
CP=DCOS(PHIL)
SP=DSIN(PHIL)
CP1=DCOS(PHIL1)
SP1=DSIN(PHIL1)
AL(1,1)=CT*CP
AL(1,2)=-SP
AL(2,1)=CT*SP
AL(2,2)=CP
AL(3,1)=-ST
AL(3,2)=0D0
AL1(1,1)=CT1*CP1
AL1(1,2)=-SP1
AL1(2,1)=CT1*SP1
AL1(2,2)=CP1
AL1(3,1)=-ST1
AL1(3,2)=0D0
C____________COMPUTE MATRICES AP^(-1) AND AP1^(-1), EQ. (15)
CT=CTP
ST=DSIN(THETP)
CP=DCOS(PHIP)
SP=DSIN(PHIP)
CT1=CTP1
ST1=DSIN(THETP1)
CP1=DCOS(PHIP1)
SP1=DSIN(PHIP1)
AP(1,1)=CT*CP
AP(1,2)=CT*SP
AP(1,3)=-ST
AP(2,1)=-SP
AP(2,2)=CP
AP(2,3)=0D0
AP1(1,1)=CT1*CP1
AP1(1,2)=CT1*SP1
AP1(1,3)=-ST1
AP1(2,1)=-SP1
AP1(2,2)=CP1
AP1(2,3)=0D0
C____________COMPUTE MATRICES R AND R^(-1), EQ. (13)
DO I=1,3
DO J=1,2
X=0D0
DO K=1,3
X=X+B(I,K)*AL(K,J)
ENDDO
C(I,J)=X
ENDDO
ENDDO
DO I=1,2
DO J=1,2
X=0D0
DO K=1,3
X=X+AP(I,K)*C(K,J)
ENDDO
R(I,J)=X
ENDDO
ENDDO
DO I=1,3
DO J=1,2
X=0D0
DO K=1,3
X=X+B(I,K)*AL1(K,J)
ENDDO
C(I,J)=X
ENDDO
ENDDO
DO I=1,2
DO J=1,2
X=0D0
DO K=1,3
X=X+AP1(I,K)*C(K,J)
ENDDO
R1(I,J)=X
ENDDO
ENDDO
D=1D0/(R1(1,1)*R1(2,2)-R1(1,2)*R1(2,1))
X=R1(1,1)
R1(1,1)=R1(2,2)*D
R1(1,2)=-R1(1,2)*D
R1(2,1)=-R1(2,1)*D
R1(2,2)=X*D
CI=(0D0,1D0)
DO 5 NN=1,NMAX
DO 5 N=1,NMAX
CN=CI**(NN-N-1)
DNN=DFLOAT((2*N+1)*(2*NN+1))
DNN=DNN/DFLOAT( N*NN*(N+1)*(NN+1) )
RN=DSQRT(DNN)
CAL(N,NN)=CN*RN
5 CONTINUE
DCTH0=CTP
DCTH=CTP1
PH=PHIP1-PHIP
VV=(0D0,0D0)
VH=(0D0,0D0)
HV=(0D0,0D0)
HH=(0D0,0D0)
DO 500 M=0,NMAX
M1=M+1
NMIN=MAX(M,1)
CALL VIGAMPL (DCTH, NMAX, M, DV1, DV2)
CALL VIGAMPL (DCTH0, NMAX, M, DV01, DV02)
FC=2D0*DCOS(M*PH)
FS=2D0*DSIN(M*PH)
DO 400 NN=NMIN,NMAX
DV1NN=M*DV01(NN)
DV2NN=DV02(NN)
DO 400 N=NMIN,NMAX
DV1N=M*DV1(N)
DV2N=DV2(N)
CT11=DCMPLX(TR11(M1,N,NN),TI11(M1,N,NN))
CT22=DCMPLX(TR22(M1,N,NN),TI22(M1,N,NN))
IF (M.EQ.0) THEN
CN=CAL(N,NN)*DV2N*DV2NN
VV=VV+CN*CT22
HH=HH+CN*CT11
ELSE
CT12=DCMPLX(TR12(M1,N,NN),TI12(M1,N,NN))
CT21=DCMPLX(TR21(M1,N,NN),TI21(M1,N,NN))
CN1=CAL(N,NN)*FC
CN2=CAL(N,NN)*FS
D11=DV1N*DV1NN
D12=DV1N*DV2NN
D21=DV2N*DV1NN
D22=DV2N*DV2NN
VV=VV+(CT11*D11+CT21*D21
& +CT12*D12+CT22*D22)*CN1
VH=VH+(CT11*D12+CT21*D22
& +CT12*D11+CT22*D21)*CN2
HV=HV-(CT11*D21+CT21*D11
& +CT12*D22+CT22*D12)*CN2
HH=HH+(CT11*D22+CT21*D12
& +CT12*D21+CT22*D11)*CN1
ENDIF
400 CONTINUE
500 CONTINUE
DK=2D0*PIN/DLAM
VV=VV/DK
VH=VH/DK
HV=HV/DK
HH=HH/DK
CVV=VV*R(1,1)+VH*R(2,1)
CVH=VV*R(1,2)+VH*R(2,2)
CHV=HV*R(1,1)+HH*R(2,1)
CHH=HV*R(1,2)+HH*R(2,2)
VV=R1(1,1)*CVV+R1(1,2)*CHV
VH=R1(1,1)*CVH+R1(1,2)*CHH
HV=R1(2,1)*CVV+R1(2,2)*CHV
HH=R1(2,1)*CVH+R1(2,2)*CHH
C WRITE (6,1005) TL,TL1,PL,PL1,ALPHA,BETA
C WRITE (6,1006)
C PRINT 1101, VV
C PRINT 1102, VH
C PRINT 1103, HV
C PRINT 1104, HH
C 1101 FORMAT ('S11=',D11.5,' + i*',D11.5)
C 1102 FORMAT ('S12=',D11.5,' + i*',D11.5)
C 1103 FORMAT ('S21=',D11.5,' + i*',D11.5)
C 1104 FORMAT ('S22=',D11.5,' + i*',D11.5)
C 1005 FORMAT ('thet0=',F6.2,' thet=',F6.2,' phi0=',F6.2,
C & ' phi=',F6.2,' alpha=',F6.2,' beta=',F6.2)
C 1006 FORMAT ('AMPLITUDE MATRIX')
RETURN
END
C*****************************************************************
C
C Calculation of the functions
C DV1(N)=dvig(0,m,n,arccos x)/sin(arccos x)
C and
C DV2(N)=[d/d(arccos x)] dvig(0,m,n,arccos x)
C 1.LE.N.LE.NMAX
C 0.LE.X.LE.1
SUBROUTINE VIGAMPL (X, NMAX, M, DV1, DV2)
INCLUDE 'ampld.par.f'
IMPLICIT REAL*8 (A-H,O-Z)
REAL*8 DV1(NPN6), DV2(NPN6)
DO 1 N=1,NMAX
DV1(N)=0D0
DV2(N)=0D0
1 CONTINUE
DX=DABS(X)
IF (DABS(1D0-DX).LE.1D-10) GO TO 100
A=1D0
QS=DSQRT(1D0-X*X)
QS1=1D0/QS
DSI=QS1
IF (M.NE.0) GO TO 20
D1=1D0
D2=X
DO 5 N=1,NMAX
QN=DFLOAT(N)
QN1=DFLOAT(N+1)
QN2=DFLOAT(2*N+1)
D3=(QN2*X*D2-QN*D1)/QN1
DER=QS1*(QN1*QN/QN2)*(-D1+D3)
DV1(N)=D2*DSI
DV2(N)=DER
D1=D2
D2=D3
5 CONTINUE
RETURN
20 QMM=DFLOAT(M*M)
DO 25 I=1,M
I2=I*2
A=A*DSQRT(DFLOAT(I2-1)/DFLOAT(I2))*QS
25 CONTINUE
D1=0D0
D2=A
DO 30 N=M,NMAX
QN=DFLOAT(N)
QN2=DFLOAT(2*N+1)
QN1=DFLOAT(N+1)
QNM=DSQRT(QN*QN-QMM)
QNM1=DSQRT(QN1*QN1-QMM)
D3=(QN2*X*D2-QNM*D1)/QNM1
DER=QS1*(-QN1*QNM*D1+QN*QNM1*D3)/QN2
DV1(N)=D2*DSI
DV2(N)=DER
D1=D2
D2=D3
30 CONTINUE
RETURN
100 IF (M.NE.1) RETURN
DO 110 N=1,NMAX
DN=DFLOAT(N*(N+1))
DN=0.5D0*DSQRT(DN)
IF (X.LT.0D0) DN=DN*(-1)**(N+1)
DV1(N)=DN
IF (X.LT.0D0) DN=-DN
DV2(N)=DN
110 CONTINUE
RETURN
END
C**********************************************************************
SUBROUTINE CONST (NGAUSS,NMAX,MMAX,P,X,W,AN,ANN,S,SS,NP,EPS)
IMPLICIT REAL*8 (A-H,O-Z)
INCLUDE 'ampld.par.f'
REAL*8 X(NPNG2),W(NPNG2),X1(NPNG1),W1(NPNG1),
* X2(NPNG1),W2(NPNG1),
* S(NPNG2),SS(NPNG2),
* AN(NPN1),ANN(NPN1,NPN1),DD(NPN1)
DO 10 N=1,NMAX
NN=N*(N+1)
AN(N)=DFLOAT(NN)
D=DSQRT(DFLOAT(2*N+1)/DFLOAT(NN))
DD(N)=D
DO 10 N1=1,N
DDD=D*DD(N1)*0.5D0
ANN(N,N1)=DDD
ANN(N1,N)=DDD
10 CONTINUE
NG=2*NGAUSS
IF (NP.EQ.-2) GO TO 11
CALL GAUSS(NG,0,0,X,W)
GO TO 19
11 NG1=DFLOAT(NGAUSS)/2D0
NG2=NGAUSS-NG1
XX=-DCOS(DATAN(EPS))
CALL GAUSS(NG1,0,0,X1,W1)
CALL GAUSS(NG2,0,0,X2,W2)
DO 12 I=1,NG1
W(I)=0.5D0*(XX+1D0)*W1(I)
X(I)=0.5D0*(XX+1D0)*X1(I)+0.5D0*(XX-1D0)
12 CONTINUE
DO 14 I=1,NG2
W(I+NG1)=-0.5D0*XX*W2(I)
X(I+NG1)=-0.5D0*XX*X2(I)+0.5D0*XX
14 CONTINUE
DO 16 I=1,NGAUSS
W(NG-I+1)=W(I)
X(NG-I+1)=-X(I)
16 CONTINUE
19 DO 20 I=1,NGAUSS
Y=X(I)
Y=1D0/(1D0-Y*Y)
SS(I)=Y
SS(NG-I+1)=Y
Y=DSQRT(Y)
S(I)=Y
S(NG-I+1)=Y
20 CONTINUE
RETURN
END
C**********************************************************************
SUBROUTINE VARY (LAM,MRR,MRI,A,EPS,NP,NGAUSS,X,P,PPI,PIR,PII,
* R,DR,DDR,DRR,DRI,NMAX)
INCLUDE 'ampld.par.f'
IMPLICIT REAL*8 (A-H,O-Z)
REAL*8 X(NPNG2),R(NPNG2),DR(NPNG2),MRR,MRI,LAM,
* Z(NPNG2),ZR(NPNG2),ZI(NPNG2),
* J(NPNG2,NPN1),Y(NPNG2,NPN1),JR(NPNG2,NPN1),
* JI(NPNG2,NPN1),DJ(NPNG2,NPN1),
* DJR(NPNG2,NPN1),DJI(NPNG2,NPN1),DDR(NPNG2),
* DRR(NPNG2),DRI(NPNG2),
* DY(NPNG2,NPN1)
COMMON /CBESS/ J,Y,JR,JI,DJ,DY,DJR,DJI
NG=NGAUSS*2
IF (NP.GT.0) CALL RSP2(X,NG,A,EPS,NP,R,DR)
IF (NP.EQ.-1) CALL RSP1(X,NG,NGAUSS,A,EPS,NP,R,DR)
IF (NP.EQ.-2) CALL RSP3(X,NG,NGAUSS,A,EPS,R,DR)
IF (NP.EQ.-3) CALL RSP4(X,NG,A,R,DR)
PI=P*2D0/LAM
PPI=PI*PI
PIR=PPI*MRR
PII=PPI*MRI
V=1D0/(MRR*MRR+MRI*MRI)
PRR=MRR*V
PRI=-MRI*V
TA=0D0