source: trunk/LMDZ.TITAN/libf/phytitan/dmiess.F @ 306

Last change on this file since 306 was 3, checked in by slebonnois, 14 years ago

Creation de repertoires:

  • chantiers : pour communiquer sur nos projets de modifs
  • documentation : pour stocker les docs

Ajout de:

  • libf/phytitan : physique de Titan
  • libf/chimtitan: chimie de Titan
  • libf/phyvenus : physique de Venus
File size: 16.1 KB
Line 
1      SUBROUTINE DMIESS(  RO,      RFR,     RFI,     THETD,     JX,
2     2                    QEXT,    QSCAT,   CTBRQS,  ELTRMX,    PI,
3     3                    TAU,     CSTHT,   SI2THT,  ACAP,      IT,
4     4                    LL,      R,       RE2,     TMAG2,     WVNO  )
5C                                                         
6C THIS VERSION OF THE FAMOUS DMIESS CODE IS ADAPTED FOR THE VAX
7C ITS ACCURACY HAS BEEN CHECKED BY COMPARISON TO CRAY RUNS.
8C              -C.P. MCKAY
9C
10      IMPLICIT REAL  (A-H,O-Z)
11C
12      COMPLEX    FNAP,      FNBP,      ACAP(LL),     W,
13     2              FNA,       FNB,       RF,           RRF,
14     3              RRFX,      WM1,       FN1,          FN2,
15     4              TC1,       TC2,       WFN(2),       Z(4),
16     5              K1,        K2,        K3,           SIN,
17     6              COS,       RC,        U(8),         DH1,
18     7              DH2,       DH4,       P24H24,       P24H21,
19     8              PSTORE,    HSTORE,    DUMMY,        DUMSQ
20C
21      DIMENSION         W(3,9000)
22      DIMENSION     T(5),      TA(4),     TB(2),        TC(2),
23     2              TD(2),     TE(2),     PI( 3,IT ),   TAU( 3,IT ),
24     3              CSTHT(IT), THETD(IT), SI2THT(IT),   ELTRMX( 4,IT,2 )
25C
26      EQUIVALENCE   ( WFN(1),TA(1) ),     (  FNA,TB(1) ),
27     2              (    FNB,TC(1) ),     ( FNAP,TD(1) ),
28     3              (   FNBP,TE(1) )
29C
30C
31C    THIS SUBROUTINE COMPUTES MIE SCATTERING BY A STRATIFIED SPHERE
32C    I.E. A PARTICLE CONSISTING OF A SPHERICAL CORE SURROUNDED BY A
33C    SPHERICAL SHELL.  THE BASIC CODE USED WAS THAT DESCRIBED IN THE
34C    REPORT: " SUBROUTINES FOR COMPUTING THE PARAMETERS OF THE
35C    ELECTROMAGNETIC RADIATION SCATTERED BY A SPHERE " J.V. DAVE,
36C    I B M SCIENTIFIC CENTER, PALO ALTO , CALIFORNIA.
37C    REPORT NO. 320 - 3236 .. MAY 1968 .
38C
39C    THE MODIFICATIONS FOR STRATIFIED SPHERES ARE DESCRIBED IN
40C        TOON AND ACKERMAN, APPL. OPTICS, IN PRESS, 1981
41C
42C
43C    THE PARAMETERS IN THE CALLING STATEMENT ARE DEFINED AS FOLLOWS :
44C      RO IS THE OUTER (SHELL) RADIUS;
45C      R  IS THE CORE RADIUS;
46C      RFR,   RFI  ARE THE REAL AND IMAGINARY PARTS OF THE SHELL INDEX
47C                  OF REFRACTION IN THE FORM (RFR - I* RFI);
48C      RE2, TMAG2  ARE THE INDEX PARTS FOR THE CORE;
49C          ( WE ASSUME SPACE HAS UNIT INDEX. )
50C      THETD(J): ANGLE IN DEGREES BETWEEN THE DIRECTIONS OF THE INCIDENT
51C          AND THE SCATTERED RADIATION.  THETD(J) IS< OR= 90.0
52C          IF THETD(J) SHOULD HAPPEN TO BE GREATER THAN 90.0, ENTER WITH
53C          SUPPLEMENTARY VALUE, SEE COMMENTS BELOW ON ELTRMX;
54C      JX: TOTAL NUMBER OF THETD FOR WHICH THE COMPUTATIONS ARE
55C          REQUIRED.  JX SHOULD NOT EXCEED IT UNLESS THE DIMENSIONS
56C          STATEMENTS ARE APPROPRIATEDLY MODIFIED;
57C
58C      THE DEFINITIONS FOR THE FOLLOWING SYMBOLS CAN BE FOUND IN"LIGHT
59C          SCATTERING BY SMALL PARTICLES,H.C.VAN DE HULST, JOHN WILEY '
60C          SONS, INC., NEW YORK, 1957" .
61C      QEXT: EFFIECIENCY FACTOR FOR EXTINCTION,VAN DE HULST,P.14 ' 127.
62C      QSCAT: EFFIECINCY FACTOR FOR SCATTERING,V.D. HULST,P.14 ' 127.
63C      CTBRQS: AVERAGE(COSINE THETA) * QSCAT,VAN DE HULST,P.128
64C      ELTRMX(I,J,K): ELEMENTS OF THE TRANSFORMATION MATRIX F,V.D.HULST
65C          ,P.34,45 ' 125. I=1: ELEMENT M SUB 2..I=2: ELEMENT M SUB 1..
66C          I = 3: ELEMENT S SUB 21.. I = 4: ELEMENT D SUB 21..
67C      ELTRMX(I,J,1) REPRESENTS THE ITH ELEMENT OF THE MATRIX FOR
68C          THE ANGLE THETD(J).. ELTRMX(I,J,2) REPRESENTS THE ITH ELEMENT
69C          OF THE MATRIX FOR THE ANGLE 180.0 - THETD(J) ..
70C
71C      IT: IS THE DIMENSION OF THETD, ELTRMX, CSTHT, PI, TAU, SI2THT,
72C          IT MUST CORRESPOND EXACTLY TO THE SECOND DIMENSION OF ELTRMX.
73C      LL: IS THE DIMENSION OF ACAP
74C          IN THE ORIGINAL PROGRAM THE DIMENSION OF ACAP WAS 7000.
75C          FOR CONSERVING SPACE THIS SHOULD BE NOT MUCH HIGHER THAN
76C          THE VALUE, N=1.1*(NREAL**2 + NIMAG**2)**.5 * X + 1
77C      WVNO: 2*PI / WAVELENGTH
78C
79C      THIS SUBROUTINE COMPUTES THE CAPITAL A FUNCTION BY MAKING USE OF
80C      DOWNWARD RECURRENCE RELATIONSHIP.
81C
82C      TA(1): REAL PART OF WFN(1).  TA(2): IMAGINARY PART OF WFN(1).
83C      TA(3): REAL PART OF WFN(2).  TA(4): IMAGINARY PART OF WFN(2).
84C      TB(1): REAL PART OF FNA.     TB(2): IMAGINARY PART OF FNA.
85C      TC(1): REAL PART OF FNB.     TC(2): IMAGINARY PART OF FNB.
86C      TD(1): REAL PART OF FNAP.    TD(2): IMAGINARY PART OF FNAP.
87C      TE(1): REAL PART OF FNBP.    TE(2): IMAGINARY PART OF FNBP.
88C      FNAP, FNBP  ARE THE PRECEDING VALUES OF FNA, FNB RESPECTIVELY.
89C
90C   IF THE CORE IS SMALL SCATTERING IS COMPUTED FOR THE SHELL ONLY
91
92c     print*,'debut dmiess ',second(0.)
93
94      IFLAG = 1
95      IF ( R/RO .LT. 1.0E-06 )   IFLAG = 2
96C
97      IF ( JX .LE. IT )   GO TO 20
98         WRITE( 6,7 )
99         WRITE( 6,6 )
100         STOP 30
101   20 RF =  CMPLX( RFR,  -RFI )
102      RC =  CMPLX( RE2,-TMAG2 )
103      X  =  RO * WVNO
104      K1 =  RC * WVNO
105      K2 =  RF * WVNO
106      ZET=0.0
107      K3 =  CMPLX( WVNO, ZET )
108      Z(1) =  K2 * RO
109      Z(2) =  K3 * RO
110      Z(3) =  K1 * R
111      Z(4) =  K2 * R
112      X1   = REAL( Z(1) )
113      Y1   = AIMAG( Z(1) )
114      X4   = REAL( Z(4) )
115      Y4   = AIMAG( Z(4) )
116c      X1   = REAL( Z(1) )
117c      Y1   = AIMAG( Z(1) )
118c      X4   = REAL( Z(4) )
119c      Y4   = AIMAG( Z(4) )
120      RRF  =  1.0 / RF
121      RX   =  1.0 / X
122      RRFX =  RRF * RX
123      T(1) =  ( X**2 ) * ( RFR**2 + RFI**2 )
124      T(1) =  SQRT( T(1) )
125      NMX1 =  1.10 * T(1)
126
127C
128      IF ( NMX1 .LE. LL-1 )   GO TO 21
129         WRITE(6,8)
130         PRINT*,'LL  =',LL
131         STOP 32
132   21 NMX2 = T(1)
133      IF ( NMX1 .GT.  150 )   GO TO 22
134         NMX1 = 150
135         NMX2 = 135
136C
137   22 ACAP( NMX1+1 )  =  ( 0.0,0.0 )
138      IF ( IFLAG .EQ. 2 )   GO TO 26
139         DO 29   N = 1,3
140   29    W( N,NMX1+1 )  =  ( 0.0,0.0 )
141   26 CONTINUE
142      DO 23   N = 1,NMX1
143         NN = NMX1 - N + 1
144         ACAP(NN) = (NN+1) * RRFX - 1.0 / ( (NN+1) * RRFX + ACAP(NN+1) )
145         IF ( IFLAG .EQ. 2 )   GO TO 23
146            DO 31   M = 1,3
147   31       W( M,NN ) = (NN+1) / Z(M+1)  -
148     1                       1.0 / (  (NN+1) / Z(M+1)  +  W( M,NN+1 )  )
149   23 CONTINUE
150      DO 30   J = 1,JX
151      IF ( THETD(J) .LT. 0.0 )  THETD(J) =  ABS( THETD(J) )
152      IF ( THETD(J) .GT. 0.0 )  GO TO 24
153      CSTHT(J)  = 1.0
154      SI2THT(J) = 0.0
155      GO TO 30
156   24 IF ( THETD(J) .GE. 90.0 )  GO TO 25
157      T(1)      =  ( 3.14159265359 * THETD(J) ) / 180.0
158      CSTHT(J)  =  COS( T(1) )
159      SI2THT(J) =  1.0 - CSTHT(J)**2
160      GO TO 30
161   25 IF ( THETD(J) .GT. 90.0 )  GO TO 28
162      CSTHT(J)  =  0.0
163      SI2THT(J) =  1.0
164      GO TO 30
165   28 WRITE( 6,5 )  THETD(J)
166      WRITE( 6,6 )
167      STOP 34
168   30 CONTINUE
169C
170      DO 35  J = 1,JX
171      PI(1,J)  =  0.0
172      PI(2,J)  =  1.0
173      TAU(1,J) =  0.0
174      TAU(2,J) =  CSTHT(J)
175   35 CONTINUE
176C
177C   INITIALIZATION OF HOMOGENEOUS SPHERE
178      T(1)   =  COS(X)
179      T(2)   =  SIN(X)
180      WM1    =  CMPLX( T(1),-T(2) )
181      WFN(1) =  CMPLX( T(2), T(1) )
182      WFN(2) =  RX * WFN(1) - WM1
183      IF ( IFLAG .EQ. 2 )   GO TO 560
184      N = 1
185C
186C INITIALIZATION PROCEDURE FOR STRATIFIED SPHERE BEGINS HERE
187C
188      SINX1   =  SIN( X1 )
189      SINX4   =  SIN( X4 )
190      COSX1   =  COS( X1 )
191      COSX4   =  COS( X4 )
192      EY1     =  EXP( Y1 )
193      E2Y1    =  EY1 * EY1
194      EY4     =  EXP( Y4 )
195      EY1MY4  =  EXP( Y1 - Y4 )
196      EY1PY4  =  EY1 * EY4
197      EY1MY4  =  EXP( Y1 - Y4 )
198      AA  =  SINX4 * ( EY1PY4 + EY1MY4 )
199      BB  =  COSX4 * ( EY1PY4 - EY1MY4 )
200      CC  =  SINX1 * ( E2Y1 + 1.0 )
201      DD  =  COSX1 * ( E2Y1 - 1.0 )
202      DENOM   =  1.0  +  E2Y1 * ( 4.0 * SINX1 * SINX1 - 2.0 + E2Y1 )
203      REALP   =  ( AA * CC  +  BB * DD ) / DENOM
204      AMAGP   =  ( BB * CC  -  AA * DD ) / DENOM
205      DUMMY   =  CMPLX( REALP, AMAGP )
206      AA  =  SINX4 * SINX4 - 0.5
207      BB  =  COSX4 * SINX4
208      P24H24  =  0.5 + CMPLX( AA,BB ) * EY4 * EY4
209      AA  =  SINX1 * SINX4  -  COSX1 * COSX4
210      BB  =  SINX1 * COSX4  +  COSX1 * SINX4
211      CC  =  SINX1 * SINX4  +  COSX1 * COSX4
212      DD  = -SINX1 * COSX4  +  COSX1 * SINX4
213      P24H21  =  0.5 * CMPLX( AA,BB ) * EY1 * EY4  +
214     2           0.5 * CMPLX( CC,DD ) * EY1MY4
215      DH4  =  Z(4) / ( 1.0 + ( 0.0,1.0 ) * Z(4) )  -  1.0 / Z(4)
216      DH1  =  Z(1) / ( 1.0 + ( 0.0,1.0 ) * Z(1) )  -  1.0 / Z(1)
217      DH2  =  Z(2) / ( 1.0 + ( 0.0,1.0 ) * Z(2) )  -  1.0 / Z(2)
218      PSTORE  =  ( DH4 + N / Z(4) )  *  ( W(3,N) + N / Z(4) )
219      P24H24  =  P24H24 / PSTORE
220      HSTORE  =  ( DH1 + N / Z(1) )  *  ( W(3,N) + N / Z(4) )
221      P24H21  =  P24H21 / HSTORE
222      PSTORE  =  ( ACAP(N) + N / Z(1) )  /  ( W(3,N) + N / Z(4) )
223      DUMMY   =  DUMMY * PSTORE
224      DUMSQ   =  DUMMY * DUMMY
225C
226C NOTE:  THE DEFINITIONS OF U(I) IN THIS PROGRAM ARE NOT THE SAME AS
227C
228C          USUB1 = U(1)                       USUB2 = U(5)
229C          USUB3 = U(7)                       USUB4 = DUMSQ
230C          USUB5 = U(2)                       USUB6 = U(3)
231C          USUB7 = U(6)                       USUB8 = U(4)
232C          RATIO OF SPHERICAL BESSEL FTN TO SPHERICAL HENKAL FTN = U(8)
233C
234      U(1) =  K3 * ACAP(N)  -  K2 * W(1,N)
235      U(2) =  K3 * ACAP(N)  -  K2 * DH2
236      U(3) =  K2 * ACAP(N)  -  K3 * W(1,N)
237      U(4) =  K2 * ACAP(N)  -  K3 * DH2
238      U(5) =  K1 *  W(3,N)  -  K2 * W(2,N)
239      U(6) =  K2 *  W(3,N)  -  K1 * W(2,N)
240      U(7) =  ( 0.0,-1.0 )  *  ( DUMMY * P24H21 - P24H24 )
241      U(8) =  TA(3) / WFN(2)
242C
243      FNA  =  U(8) * ( U(1)*U(5)*U(7)  +  K1*U(1)  -  DUMSQ*K3*U(5) ) /
244     2               ( U(2)*U(5)*U(7)  +  K1*U(2)  -  DUMSQ*K3*U(5) )
245      FNB  =  U(8) * ( U(3)*U(6)*U(7)  +  K2*U(3)  -  DUMSQ*K2*U(6) ) /
246     2               ( U(4)*U(6)*U(7)  +  K2*U(4)  -  DUMSQ*K2*U(6) )
247      GO TO 561
248  560 TC1  =  ACAP(1) * RRF  +  RX
249      TC2  =  ACAP(1) * RF   +  RX
250      FNA  =  ( TC1 * TA(3)  -  TA(1) ) / ( TC1 * WFN(2)  -  WFN(1) )
251      FNB  =  ( TC2 * TA(3)  -  TA(1) ) / ( TC2 * WFN(2)  -  WFN(1) )
252  561 CONTINUE
253      FNAP = FNA
254      FNBP = FNB
255      T(1) = 1.50
256C
257C    FROM HERE TO THE STATMENT NUMBER 90, ELTRMX(I,J,K) HAS
258C    FOLLOWING MEANING:
259C    ELTRMX(1,J,K): REAL PART OF THE FIRST COMPLEX AMPLITUDE.
260C    ELTRMX(2,J,K): IMAGINARY PART OF THE FIRST COMPLEX AMPLITUDE.
261C    ELTRMX(3,J,K): REAL PART OF THE SECOND COMPLEX AMPLITUDE.
262C    ELTRMX(4,J,K): IMAGINARY PART OF THE SECOND COMPLEX AMPLITUDE.
263C    K = 1 : FOR THETD(J) AND K = 2 : FOR 180.0 - THETD(J)
264C    DEFINITION OF THE COMPLEX AMPLITUDE: VAN DE HULST,P.125.
265      TB(1) = T(1) * TB(1)
266      TB(2) = T(1) * TB(2)
267      TC(1) = T(1) * TC(1)
268      TC(2) = T(1) * TC(2)
269      DO 60 J = 1,JX
270          ELTRMX(1,J,1) = TB(1) * PI(2,J) + TC(1) * TAU(2,J)
271          ELTRMX(2,J,1) = TB(2) * PI(2,J) + TC(2) * TAU(2,J)
272          ELTRMX(3,J,1) = TC(1) * PI(2,J) + TB(1) * TAU(2,J)
273          ELTRMX(4,J,1) = TC(2) * PI(2,J) + TB(2) * TAU(2,J)
274          ELTRMX(1,J,2) = TB(1) * PI(2,J) - TC(1) * TAU(2,J)
275          ELTRMX(2,J,2) = TB(2) * PI(2,J) - TC(2) * TAU(2,J)
276          ELTRMX(3,J,2) = TC(1) * PI(2,J) - TB(1) * TAU(2,J)
277          ELTRMX(4,J,2) = TC(2) * PI(2,J) - TB(2) * TAU(2,J)
278   60 CONTINUE
279C
280      QEXT   = 2.0 * ( TB(1) + TC(1))
281      QSCAT  = ( TB(1)**2 + TB(2)**2 + TC(1)**2 + TC(2)**2 ) / 0.75
282      CTBRQS = 0.0
283      N = 2
284   65 T(1) = 2*N - 1
285      T(2) =   N - 1
286      T(3) = 2*N + 1
287      DO 70  J = 1,JX
288          PI(3,J)  = ( T(1) * PI(2,J) * CSTHT(J) - N * PI(1,J) ) / T(2)
289          TAU(3,J) = CSTHT(J) * ( PI(3,J) - PI(1,J) )  -
290     1                          T(1) * SI2THT(J) * PI(2,J)  +  TAU(1,J)
291   70 CONTINUE
292C
293C   HERE SET UP HOMOGENEOUS SPHERE
294      WM1    =  WFN(1)
295      WFN(1) =  WFN(2)
296      WFN(2) =  T(1) * RX * WFN(1)  -  WM1
297      IF ( IFLAG .EQ. 2 )   GO TO 1000
298C
299C   HERE SET UP STRATIFIED SPHERE
300C
301      DH2  =  - N / Z(2)  +  1.0 / ( N / Z(2) - DH2 )
302      DH4  =  - N / Z(4)  +  1.0 / ( N / Z(4) - DH4 )
303      DH1  =  - N / Z(1)  +  1.0 / ( N / Z(1) - DH1 )
304      PSTORE  =  ( DH4 + N / Z(4) )  *  ( W(3,N) + N / Z(4) )
305      P24H24  =  P24H24 / PSTORE
306      HSTORE  =  ( DH1 + N / Z(1) )  *  ( W(3,N) + N / Z(4) )
307      P24H21  =  P24H21 / HSTORE
308      PSTORE  =  ( ACAP(N) + N / Z(1) )  /  ( W(3,N) + N / Z(4) )
309      DUMMY   =  DUMMY * PSTORE
310      DUMSQ   =  DUMMY * DUMMY
311C
312      U(1) =  K3 * ACAP(N)  -  K2 * W(1,N)
313      U(2) =  K3 * ACAP(N)  -  K2 * DH2
314      U(3) =  K2 * ACAP(N)  -  K3 * W(1,N)
315      U(4) =  K2 * ACAP(N)  -  K3 * DH2
316      U(5) =  K1 *  W(3,N)  -  K2 * W(2,N)
317      U(6) =  K2 *  W(3,N)  -  K1 * W(2,N)
318      U(7) =  ( 0.0,-1.0 )  *  ( DUMMY * P24H21 - P24H24 )
319      U(8) =  TA(3) / WFN(2)
320C
321      FNA  =  U(8) * ( U(1)*U(5)*U(7)  +  K1*U(1)  -  DUMSQ*K3*U(5) ) /
322     2               ( U(2)*U(5)*U(7)  +  K1*U(2)  -  DUMSQ*K3*U(5) )
323      FNB  =  U(8) * ( U(3)*U(6)*U(7)  +  K2*U(3)  -  DUMSQ*K2*U(6) ) /
324     2               ( U(4)*U(6)*U(7)  +  K2*U(4)  -  DUMSQ*K2*U(6) )
325C
326 1000 CONTINUE
327      TC1  =  ACAP(N) * RRF  +  N * RX
328      TC2  =  ACAP(N) * RF   +  N * RX
329      FN1  =  ( TC1 * TA(3)  -  TA(1) ) / ( TC1 * WFN(2) - WFN(1) )
330      FN2  =  ( TC2 * TA(3)  -  TA(1) ) / ( TC2 * WFN(2) - WFN(1) )
331      M    =  WVNO * R
332      IF ( N .LT. M )   GO TO 1002
333      IF ( IFLAG .EQ. 2 )   GO TO 1001
334C!!!!!!!!!!!! WARNING MODIF PERSO
335C     IF (  ABS(  ( FN1-FNA ) / FN1  )  .LT.  1.0E-09   .AND.
336C    1      ABS(  ( FN2-FNB ) / FN2  )  .LT . 1.0E-09  )       IFLAG = 2
337
338      IF (  ABS(  ( FN1-FNA ) / FN1  )  .LT.  1.0E-4   .AND.
339     1      ABS(  ( FN2-FNB ) / FN2  )  .LT . 1.0E-4  )       IFLAG = 2
340      IF ( IFLAG .EQ. 1 )   GO TO 1002
341 1001 FNA  =  FN1
342      FNB  =  FN2
343 1002 CONTINUE
344      T(5)  =  N
345      T(4)  =  T(1) / ( T(5) * T(2) )
346      T(2)  =  (  T(2) * ( T(5) + 1.0 )  ) / T(5)
347C
348      CTBRQS  =  CTBRQS  +  T(2) * ( TD(1) * TB(1)  +  TD(2) * TB(2)  +
349     1                               TE(1) * TC(1)  +  TE(2) * TC(2) )
350     2                   +  T(4) * ( TD(1) * TE(1)  +  TD(2) * TE(2) )
351      QEXT    =   QEXT  +  T(3) * ( TB(1) + TC(1) )
352      T(4)    =  TB(1)**2 + TB(2)**2 + TC(1)**2 + TC(2)**2
353      QSCAT   =  QSCAT  +  T(3) * T(4)
354      T(2)    =  N * (N+1)
355      T(1)    =  T(3) / T(2)
356      K = (N/2)*2
357      DO 80 J = 1,JX
358       ELTRMX(1,J,1) = ELTRMX(1,J,1)+T(1)*(TB(1)*PI(3,J)+TC(1)*TAU(3,J))
359       ELTRMX(2,J,1) = ELTRMX(2,J,1)+T(1)*(TB(2)*PI(3,J)+TC(2)*TAU(3,J))
360       ELTRMX(3,J,1) = ELTRMX(3,J,1)+T(1)*(TC(1)*PI(3,J)+TB(1)*TAU(3,J))
361       ELTRMX(4,J,1) = ELTRMX(4,J,1)+T(1)*(TC(2)*PI(3,J)+TB(2)*TAU(3,J))
362      IF ( K .EQ. N )   GO TO 75
363       ELTRMX(1,J,2) = ELTRMX(1,J,2)+T(1)*(TB(1)*PI(3,J)-TC(1)*TAU(3,J))
364       ELTRMX(2,J,2) = ELTRMX(2,J,2)+T(1)*(TB(2)*PI(3,J)-TC(2)*TAU(3,J))
365       ELTRMX(3,J,2) = ELTRMX(3,J,2)+T(1)*(TC(1)*PI(3,J)-TB(1)*TAU(3,J))
366       ELTRMX(4,J,2) = ELTRMX(4,J,2)+T(1)*(TC(2)*PI(3,J)-TB(2)*TAU(3,J))
367      GO TO 80
368   75  ELTRMX(1,J,2) =ELTRMX(1,J,2)+T(1)*(-TB(1)*PI(3,J)+TC(1)*TAU(3,J))
369       ELTRMX(2,J,2) =ELTRMX(2,J,2)+T(1)*(-TB(2)*PI(3,J)+TC(2)*TAU(3,J))
370       ELTRMX(3,J,2) =ELTRMX(3,J,2)+T(1)*(-TC(1)*PI(3,J)+TB(1)*TAU(3,J))
371       ELTRMX(4,J,2) =ELTRMX(4,J,2)+T(1)*(-TC(2)*PI(3,J)+TB(2)*TAU(3,J))
372   80 CONTINUE
373C
374C!!!!!!!!!!!! WARNING MODIF PERSO
375C      IF ( T(4) .LT. 1.0E-14 )   GO TO 100
376      IF ( T(4) .LT. 1.0E-4 )   GO TO 100
377      N = N + 1
378      DO 90 J = 1,JX
379         PI(1,J)   =   PI(2,J)
380         PI(2,J)   =   PI(3,J)
381         TAU(1,J)  =  TAU(2,J)
382         TAU(2,J)  =  TAU(3,J)
383   90 CONTINUE
384      FNAP  =  FNA
385      FNBP  =  FNB
386c     print*,'NMX2 =',nmx2
387      IF ( N .LE. NMX2 )   GO TO 65
388         WRITE( 6,8 )
389         STOP 36
390  100 DO 120 J = 1,JX
391      DO 120 K = 1,2
392         DO  115  I= 1,4
393         T(I)  =  ELTRMX(I,J,K)
394  115    CONTINUE
395         ELTRMX(2,J,K)  =      T(1)**2  +  T(2)**2
396         ELTRMX(1,J,K)  =      T(3)**2  +  T(4)**2
397         ELTRMX(3,J,K)  =  T(1) * T(3)  +  T(2) * T(4)
398         ELTRMX(4,J,K)  =  T(2) * T(3)  -  T(4) * T(1)
399  120 CONTINUE
400      T(1)    =    2.0 * RX**2
401      QEXT    =   QEXT * T(1)
402      QSCAT   =  QSCAT * T(1)
403      CTBRQS  =  2.0 * CTBRQS * T(1)
404C
405c     print*,'NMX1= ',nmx1,'  LL=',ll
406c     print*,'fin dmiess ',second(0.)
407      RETURN
408C
409    5 FORMAT( 10X,' THE VALUE OF THE SCATTERING ANGLE IS GREATER THAN
410     1 90.0 DEGREES. IT IS ', E15.4 )
411    6 FORMAT( // 10X, 'PLEASE READ COMMENTS.' // )
412    7 FORMAT( // 10X, 'THE VALUE OF THE ARGUMENT JX IS GREATER THAN IT')
413    8 FORMAT( // 10X, 'THE UPPER LIMIT FOR ACAP IS NOT ENOUGH. SUGGEST
414     1 GET DETAILED OUTPUT AND MODIFY SUBROUTINE' // )
415C
416      END
Note: See TracBrowser for help on using the repository browser.