/
analytical_solution_viscoelastic_2D_plane_strain_Carcione_correct_with_1_over_L.f90
7824 lines (7658 loc) · 259 KB
/
analytical_solution_viscoelastic_2D_plane_strain_Carcione_correct_with_1_over_L.f90
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
program analytical_solution
!! DK DK we compute the solution for velocity instead of for displacement in this version of the analytical code.
! This program implements the analytical solution for the velocity vector in a 2D plane-strain viscoelastic medium
! with a vertical force source located in (0,0),
! from Appendix B of Carcione et al., Wave propagation simulation in a linear viscoelastic medium, GJI, vol. 95, p. 597-611 (1988)
! (note that that Appendix contains two typos, fixed in this code; I added two comments below to mention them).
! The amplitude of the force is called F and is defined below.
implicit none
!! DK DK May 2018: the missing 1/L factor in older Carcione papers
!! DK DK May 2018: has been added to this code by Quentin Brissaud and by Etienne Bachmann
!! DK DK for the viscoacoustic code in directory EXAMPLES/attenuation/viscoacoustic,
!! DK DK it would be very easy to copy the changes from there to this viscoelastic version;
!! DK DK but then all the values of the tau_epsilon in the code below would need to change.
!! DK DK Dimitri Komatitsch, CNRS Marseille, France, April 2017: added the elastic reference calculation.
! compute the non-viscoacoustic case as a reference if needed, i.e. turn attenuation off
logical, parameter :: TURN_ATTENUATION_OFF = .false.
! to see how small the contribution of the near-field term is,
! here the user can ask not to include it, to then compare with the full result obtained with this flag set to false
logical, parameter :: DO_NOT_COMPUTE_THE_NEAR_FIELD = .false.
integer, parameter :: iratio = 64
integer, parameter :: nfreq = 524288
integer, parameter :: nt = iratio * nfreq
double precision, parameter :: freqmax = 400.d0
!! DK DK to print the velocity if we want to display the curve of how velocity varies with frequency
!! DK DK for instance to compute the unrelaxed velocity in the Zener model
! double precision, parameter :: freqmax = 20000.d0
double precision, parameter :: freqseuil = 0.00005d0
double precision, parameter :: pi = 3.141592653589793d0
! for the solution in time domain
integer it,i
real wsave(4*nt+15)
complex c(nt)
!! DK DK for my slow inverse Discrete Fourier Transform using a double loop
complex :: input(nt), i_imaginary_constant
integer :: j,m
! density of the medium
double precision, parameter :: rho = 2000.d0
! unrelaxed (f = +infinity) values
! these values for the unrelaxed state are computed from the relaxed state values (Vp = 3000, Vs = 2000, rho = 2000)
! given in Carcione et al. 1988 GJI vol 95 p 604 Table 1
double precision, parameter :: Vp = 2000.d0
double precision, parameter :: Vs = Vp / 1.732d0
! unrelaxed (f = +infinity) values, i.e. using the fastest Vp and Vs velocities
double precision, parameter :: M2_unrelaxed = Vs**2 * 2.d0 * rho
double precision, parameter :: M1_unrelaxed = 2.d0 * Vp**2 * rho - M2_unrelaxed
! amplitude of the force source
double precision, parameter :: F = 1.d0
! definition position recepteur Carcione
double precision x1,x2
! Definition source Dimitri
double precision, parameter :: f0 = 35.d0
double precision, parameter :: t0 = 1.2d0 / f0
! Definition source Carcione
! double precision f0,t0,eta,epsil
! parameter(f0 = 50.d0)
! parameter(t0 = 0.075d0)
! parameter(epsil = 1.d0)
! parameter(eta = 0.5d0)
! number of Zener standard linear solids in parallel
integer, parameter :: Lnu = 3
! DK DK I implemented a very simple and slow inverse Discrete Fourier Transform
! DK DK at some point, for verification, using a double loop. I keep it just in case.
! DK DK For large number of points it is extremely slow because of the double loop.
! DK DK Thus there is no reason to turn this flag on.
logical, parameter :: USE_SLOW_FOURIER_TRANSFORM = .false.
!! DK DK March 2018: this missing 1/L factor has been added to this code by Quentin Brissaud
!! DK DK for the viscoacoustic code in directory EXAMPLES/attenuation/viscoacoustic,
!! DK DK it would be very easy to copy the changes from there to this viscoelastic version;
!! DK DK but then all the values of the tau_epsilon below would need to change.
double precision, dimension(Lnu) :: tau_sigma_nu1,tau_sigma_nu2,tau_epsilon_nu1,tau_epsilon_nu2
integer :: ifreq,icalculation
double precision :: deltafreq,freq,omega,omega0,deltat,time,a
double complex :: comparg
! Fourier transform of the Ricker wavelet source
double complex fomega(0:nfreq)
! real and imaginary parts
double precision ra(0:nfreq),rb(0:nfreq)
! spectral amplitude
double precision ampli(0:nfreq)
! analytical solution for the two components
double complex phi1(-nfreq:nfreq)
double complex phi2(-nfreq:nfreq)
! external functions
double complex, external :: u1,u2
! modules elastiques
double complex :: M1C, M2C, E, V1, V2, temp
! ********** end of variable declarations ************
! classical least-squares constants
tau_epsilon_nu1 = (/ 2.408158185753685d-002, 4.699608990861351d-003, 9.567997872435925d-004 /)
tau_sigma_nu1 = (/ 2.256014638636808d-002, 4.508471279712252d-003, 8.937876403768840d-004 /)
tau_epsilon_nu2 = (/ 2.430544480527216d-002, 4.728107829226396d-003, 9.667252695863502d-004 /)
tau_sigma_nu2 = (/ 2.250919779429490d-002, 4.501388007338097d-003, 8.917332095369118d-004 /)
! do the calculation twice, because in the finite-difference code that we want to check using this analytical code
! the Vy component is staggered half a grid cell away from Vx, thus to compute the analytical solution we need
! to slightly change the location at which the calculation is done when computing the second component,
! by half a grid cell
do icalculation = 1,2
! position of the receiver
if (icalculation == 1) then
x1 = +801. - 1.5/2.
x2 = +801. - 1.5/2.
else if (icalculation == 2) then
x1 = +801.
x2 = +801.
else
stop 'wrong value of icalculation'
endif
print *,'Force source located at the origin (0,0)'
print *,'Receiver located in (x,z) = ',x1,x2
if (TURN_ATTENUATION_OFF) then
print *,'BEWARE: computing the elastic reference solution (i.e., without attenuation) instead of the viscoelastic solution'
else
print *,'Computing the viscoelastic solution'
endif
if (DO_NOT_COMPUTE_THE_NEAR_FIELD) then
print *,'BEWARE: computing the far-field solution only, rather than the full Green function'
else
print *,'Computing the full solution, including the near-field term of the Green function'
endif
! step in frequency
deltafreq = freqmax / dble(nfreq)
! define parameters for the Ricker source
omega0 = 2.d0 * pi * f0
a = pi**2 * f0**2
deltat = 1.d0 / (freqmax*dble(iratio))
! define the spectrum of the source
do ifreq=0,nfreq
freq = deltafreq * dble(ifreq)
omega = 2.d0 * pi * freq
! typo in equation (B7) of Carcione et al., Wave propagation simulation in a linear viscoelastic medium,
! Geophysical Journal, vol. 95, p. 597-611 (1988), the exponential should be of -i omega t0,
! fixed here by adding the minus sign
comparg = dcmplx(0.d0,-omega*t0)
! definir le spectre du Ricker de Carcione avec cos()
! equation (B7) of Carcione et al., Wave propagation simulation in a linear viscoelastic medium,
! Geophysical Journal, vol. 95, p. 597-611 (1988)
! fomega(ifreq) = pi * dsqrt(pi/eta) * (1.d0/omega0) * cdexp(comparg) * ( dexp(- (pi*pi/eta) * (epsil/2 - omega/omega0)**2) &
! + dexp(- (pi*pi/eta) * (epsil/2 + omega/omega0)**2) )
! definir le spectre d'un Ricker classique (centre en t0)
fomega(ifreq) = dsqrt(pi) * cdexp(comparg) * omega**2 * dexp(-omega**2/(4.d0*a)) / (2.d0 * dsqrt(a**3))
!! DK DK multiply by i omega in order to get the solution for velocity instead of for displacement
fomega(ifreq) = fomega(ifreq) * dcmplx(0.d0,omega)
ra(ifreq) = dreal(fomega(ifreq))
rb(ifreq) = dimag(fomega(ifreq))
! prendre le module de l'amplitude spectrale
ampli(ifreq) = dsqrt(ra(ifreq)**2 + rb(ifreq)**2)
enddo
! sauvegarde du spectre d'amplitude de la source en Hz au format Gnuplot
open(unit=10,file='spectrum_of_the_source_used.gnu',status='unknown')
do ifreq = 0,nfreq
freq = deltafreq * dble(ifreq)
write(10,*) sngl(freq),sngl(ampli(ifreq))
enddo
close(10)
! ************** calcul solution analytique ****************
! d'apres Carcione GJI vol 95 p 611 (1988)
do ifreq=0,nfreq
freq = deltafreq * dble(ifreq)
omega = 2.d0 * pi * freq
! critere ad-hoc pour eviter singularite en zero
if (freq < freqseuil) omega = 2.d0 * pi * freqseuil
! use standard infinite frequency (unrelaxed) reference,
! in which waves slow down when attenuation is turned on.
temp = dcmplx(0.d0,0.d0)
do i=1,Lnu
temp = temp + dcmplx(1.d0,omega*tau_epsilon_nu1(i)) / dcmplx(1.d0,omega*tau_sigma_nu1(i))
enddo
M1C = (M1_unrelaxed /(sum(tau_epsilon_nu1(:)/tau_sigma_nu1(:)))) * temp
temp = dcmplx(0.d0,0.d0)
do i=1,Lnu
temp = temp + dcmplx(1.d0,omega*tau_epsilon_nu2(i)) / dcmplx(1.d0,omega*tau_sigma_nu2(i))
enddo
M2C = (M2_unrelaxed /(sum(tau_epsilon_nu2(:)/tau_sigma_nu2(:)))) * temp
if (TURN_ATTENUATION_OFF) then
! from Etienne Bachmann, May 2018: pour calculer la solution sans attenuation, il faut donner le Mu_unrelaxed et pas le Mu_relaxed.
! En effet, pour comparer avec SPECFEM, il faut simplement partir de la bonne reference.
! SPECFEM est defini en unrelaxed et les constantes unrelaxed dans Carcione matchent parfaitement les Vp et Vs definis dans SPECFEM.
M1C = M1_unrelaxed
M2C = M2_unrelaxed
endif
E = (M1C + M2C) / 2
V1 = cdsqrt(E / rho) !! DK DK this is Vp
!! DK DK print the velocity if we want to display the curve of how velocity varies with frequency
!! DK DK for instance to compute the unrelaxed velocity in the Zener model
! print *,freq,dsqrt(real(V1)**2 + imag(V1)**2)
V2 = cdsqrt(M2C / (2.d0 * rho)) !! DK DK this is Vs
!! DK DK print the velocity if we want to display the curve of how velocity varies with frequency
!! DK DK for instance to compute the unrelaxed velocity in the Zener model
! print *,freq,dsqrt(real(V2)**2 + imag(V2)**2)
! calcul de la solution analytique en frequence
phi1(ifreq) = u1(omega,V1,V2,x1,x2,rho,F,DO_NOT_COMPUTE_THE_NEAR_FIELD) * fomega(ifreq)
phi2(ifreq) = u2(omega,V1,V2,x1,x2,rho,F,DO_NOT_COMPUTE_THE_NEAR_FIELD) * fomega(ifreq)
enddo
! take the conjugate value for negative frequencies
do ifreq=-nfreq,-1
phi1(ifreq) = dconjg(phi1(-ifreq))
phi2(ifreq) = dconjg(phi2(-ifreq))
enddo
! save the result in the frequency domain
! open(unit=11,file='cmplx_phi',status='unknown')
! do ifreq=-nfreq,nfreq
! freq = deltafreq * dble(ifreq)
! write(11,*) sngl(freq),sngl(dreal(phi1(ifreq))),sngl(dimag(phi1(ifreq))),sngl(dreal(phi2(ifreq))),sngl(dimag(phi2(ifreq)))
! enddo
! close(11)
! ***************************************************************************
! Calculation of the time domain solution (using routine "cfftb" from Netlib)
! ***************************************************************************
! **********
! Compute Vx
! **********
if (icalculation == 1) then
! initialize FFT arrays
call cffti(nt,wsave)
! clear array of Fourier coefficients
do it = 1,nt
c(it) = cmplx(0.,0.)
enddo
! use the Fourier values for Vx
c(1) = cmplx(phi1(0))
do ifreq=1,nfreq-2
c(ifreq+1) = cmplx(phi1(ifreq))
c(nt+1-ifreq) = conjg(cmplx(phi1(ifreq)))
enddo
! perform the inverse FFT for Vx
if (.not. USE_SLOW_FOURIER_TRANSFORM) then
call cfftb(nt,c,wsave)
else
! DK DK I implemented a very simple and slow inverse Discrete Fourier Transform here
! DK DK at some point, for verification, using a double loop. I keep it just in case.
! DK DK For large number of points it is extremely slow because of the double loop.
input(:) = c(:)
! imaginary constant "i"
i_imaginary_constant = (0.,1.)
do it = 1,nt
if (mod(it,1000) == 0) print *,'FFT inverse it = ',it,' out of ',nt
j = it
c(j) = cmplx(0.,0.)
do m = 1,nt
c(j) = c(j) + input(m) * exp(2.d0 * PI * i_imaginary_constant * dble((m-1) * (j-1)) / nt)
enddo
enddo
endif
! in the inverse Discrete Fourier transform one needs to divide by N, the number of samples (number of time steps here)
c(:) = c(:) / nt
! value of a time step
deltat = 1.d0 / (freqmax*dble(iratio))
! to get the amplitude right, we need to divide by the time step
c(:) = c(:) / deltat
! save time result inverse FFT for Vx
if (TURN_ATTENUATION_OFF) then
if (DO_NOT_COMPUTE_THE_NEAR_FIELD) then
open(unit=11,file='Vx_time_analytical_solution_elastic_without_near_field.dat',status='unknown')
else
open(unit=11,file='Vx_time_analytical_solution_elastic.dat',status='unknown')
endif
else
if (DO_NOT_COMPUTE_THE_NEAR_FIELD) then
open(unit=11,file='Vx_time_analytical_solution_viscoelastic_without_near_field.dat',status='unknown')
else
open(unit=11,file='Vx_time_analytical_solution_viscoelastic.dat',status='unknown')
endif
endif
do it=1,nt
! DK DK Dec 2011: subtract t0 to be consistent with the SPECFEM2D code
time = dble(it-1)*deltat - t0
! the seismograms are very long due to the very large number of FFT points used,
! thus keeping the useful part of the signal only (the first six seconds of the seismogram)
if (time >= 0.d0 .and. time <= 6.d0) write(11,*) sngl(time),real(c(it))
enddo
close(11)
endif ! of if (icalculation == 1) then
! **********
! Compute Vz
! **********
if (icalculation == 2) then
! clear array of Fourier coefficients
do it = 1,nt
c(it) = cmplx(0.,0.)
enddo
! use the Fourier values for Vz
c(1) = cmplx(phi2(0))
do ifreq=1,nfreq-2
c(ifreq+1) = cmplx(phi2(ifreq))
c(nt+1-ifreq) = conjg(cmplx(phi2(ifreq)))
enddo
! perform the inverse FFT for Vz
if (.not. USE_SLOW_FOURIER_TRANSFORM) then
call cfftb(nt,c,wsave)
else
! DK DK I implemented a very simple and slow inverse Discrete Fourier Transform here
! DK DK at some point, for verification, using a double loop. I keep it just in case.
! DK DK For large number of points it is extremely slow because of the double loop.
input(:) = c(:)
! imaginary constant "i"
i_imaginary_constant = (0.,1.)
do it = 1,nt
if (mod(it,1000) == 0) print *,'FFT inverse it = ',it,' out of ',nt
j = it
c(j) = cmplx(0.,0.)
do m = 1,nt
c(j) = c(j) + input(m) * exp(2.d0 * PI * i_imaginary_constant * dble((m-1) * (j-1)) / nt)
enddo
enddo
endif
! in the inverse Discrete Fourier transform one needs to divide by N, the number of samples (number of time steps here)
c(:) = c(:) / nt
! value of a time step
deltat = 1.d0 / (freqmax*dble(iratio))
! to get the amplitude right, we need to divide by the time step
c(:) = c(:) / deltat
! save time result inverse FFT for Vz
if (TURN_ATTENUATION_OFF) then
if (DO_NOT_COMPUTE_THE_NEAR_FIELD) then
open(unit=11,file='Vz_time_analytical_solution_elastic_without_near_field.dat',status='unknown')
else
open(unit=11,file='Vz_time_analytical_solution_elastic.dat',status='unknown')
endif
else
if (DO_NOT_COMPUTE_THE_NEAR_FIELD) then
open(unit=11,file='Vz_time_analytical_solution_viscoelastic_without_near_field.dat',status='unknown')
else
open(unit=11,file='Vz_time_analytical_solution_viscoelastic.dat',status='unknown')
endif
endif
do it=1,nt
! DK DK Dec 2011: subtract t0 to be consistent with the SPECFEM2D code
time = dble(it-1)*deltat - t0
! the seismograms are very long due to the very large number of FFT points used,
! thus keeping the useful part of the signal only (the first six seconds of the seismogram)
if (time >= 0.d0 .and. time <= 6.d0) write(11,*) sngl(time),real(c(it))
enddo
close(11)
endif ! of if (icalculation == 2) then
enddo ! of do icalculation = 1,2
end
! -----------
double complex function u1(omega,v1,v2,x1,x2,rho,F,DO_NOT_COMPUTE_THE_NEAR_FIELD)
implicit none
double precision omega
double complex v1,v2
logical :: DO_NOT_COMPUTE_THE_NEAR_FIELD
double complex G1,G2
external G1,G2
double precision, parameter :: pi = 3.141592653589793d0
! amplitude of the force
double precision F
double precision x1,x2,r,rho
! source-receiver distance
r = dsqrt(x1**2 + x2**2)
u1 = F * x1 * x2 * (G1(r,omega,v1,v2,DO_NOT_COMPUTE_THE_NEAR_FIELD) + &
G2(r,omega,v1,v2,DO_NOT_COMPUTE_THE_NEAR_FIELD)) / (2.d0 * pi * rho * r**2)
end
! -----------
double complex function u2(omega,v1,v2,x1,x2,rho,F,DO_NOT_COMPUTE_THE_NEAR_FIELD)
implicit none
double precision omega
double complex v1,v2
logical :: DO_NOT_COMPUTE_THE_NEAR_FIELD
double complex G1,G2
external G1,G2
double precision, parameter :: pi = 3.141592653589793d0
! amplitude of the force
double precision F
double precision x1,x2,r,rho
! source-receiver distance
r = dsqrt(x1**2 + x2**2)
u2 = F * (x2*x2*G1(r,omega,v1,v2,DO_NOT_COMPUTE_THE_NEAR_FIELD) - &
x1*x1*G2(r,omega,v1,v2,DO_NOT_COMPUTE_THE_NEAR_FIELD)) / (2.d0 * pi * rho * r**2)
end
! -----------
double complex function G1(r,omega,v1,v2,DO_NOT_COMPUTE_THE_NEAR_FIELD)
implicit none
double precision r,omega
double complex v1,v2
logical :: DO_NOT_COMPUTE_THE_NEAR_FIELD
double complex hankel0,hankel1
external hankel0,hankel1
double precision, parameter :: pi = 3.141592653589793d0
! typo in equations (B4a) and (B4b) of Carcione et al., Wave propagation simulation in a linear viscoelastic medium,
! Geophysical Journal, vol. 95, p. 597-611 (1988), fixed here: omega/(r*v) -> omega*r/v
if (DO_NOT_COMPUTE_THE_NEAR_FIELD) then
G1 = (hankel0(omega*r/v1)/(v1**2)) * dcmplx(0.d0,-pi/2.d0)
else
G1 = (hankel0(omega*r/v1)/(v1**2) + hankel1(omega*r/v2)/(omega*r*v2) - hankel1(omega*r/v1)/(omega*r*v1)) * dcmplx(0.d0,-pi/2.d0)
endif
end
! -----------
double complex function G2(r,omega,v1,v2,DO_NOT_COMPUTE_THE_NEAR_FIELD)
implicit none
double precision r,omega
double complex v1,v2
logical :: DO_NOT_COMPUTE_THE_NEAR_FIELD
double complex hankel0,hankel1
external hankel0,hankel1
double precision, parameter :: pi = 3.141592653589793d0
! typo in equations (B4a) and (B4b) of Carcione et al., Wave propagation simulation in a linear viscoelastic medium,
! Geophysical Journal, vol. 95, p. 597-611 (1988), fixed here: omega/(r*v) -> omega*r/v
if (DO_NOT_COMPUTE_THE_NEAR_FIELD) then
G2 = (hankel0(omega*r/v2)/(v2**2)) * dcmplx(0.d0,+pi/2.d0)
else
G2 = (hankel0(omega*r/v2)/(v2**2) - hankel1(omega*r/v2)/(omega*r*v2) + hankel1(omega*r/v1)/(omega*r*v1)) * dcmplx(0.d0,+pi/2.d0)
endif
end
! -----------
double complex function hankel0(z)
implicit none
double complex z
! on utilise la routine NAG appelee S17DLE (simple precision)
integer ifail,nz
complex result
ifail = -1
call S17DLE(2,0.0,cmplx(z),1,'U',result,nz,ifail)
if (ifail /= 0) stop 'S17DLE failed in hankel0'
if (nz > 0) print *,nz,' termes mis a zero par underflow'
hankel0 = dcmplx(result)
end
! -----------
double complex function hankel1(z)
implicit none
double complex z
! on utilise la routine NAG appelee S17DLE (simple precision)
integer ifail,nz
complex result
ifail = -1
call S17DLE(2,1.0,cmplx(z),1,'U',result,nz,ifail)
if (ifail /= 0) stop 'S17DLE failed in hankel1'
if (nz > 0) print *,nz,' termes mis a zero par underflow'
hankel1 = dcmplx(result)
end
! ***************** routine de FFT pour signal en temps ****************
! FFT routine taken from Netlib
subroutine CFFTB (N,C,WSAVE)
DIMENSION C(1) ,WSAVE(1)
if (N == 1) return
IW1 = N+N+1
IW2 = IW1+N+N
CALL CFFTB1 (N,C,WSAVE,WSAVE(IW1),WSAVE(IW2))
return
END
subroutine CFFTB1 (N,C,CH,WA,IFAC)
DIMENSION CH(1) ,C(1) ,WA(1) ,IFAC(1)
NF = IFAC(2)
NA = 0
L1 = 1
IW = 1
DO 116 K1=1,NF
IP = IFAC(K1+2)
L2 = IP*L1
IDO = N/L2
IDOT = IDO+IDO
IDL1 = IDOT*L1
if (IP /= 4) goto 103
IX2 = IW+IDOT
IX3 = IX2+IDOT
if (NA /= 0) goto 101
CALL PASSB4 (IDOT,L1,C,CH,WA(IW),WA(IX2),WA(IX3))
goto 102
101 CALL PASSB4 (IDOT,L1,CH,C,WA(IW),WA(IX2),WA(IX3))
102 NA = 1-NA
goto 115
103 if (IP /= 2) goto 106
if (NA /= 0) goto 104
CALL PASSB2 (IDOT,L1,C,CH,WA(IW))
goto 105
104 CALL PASSB2 (IDOT,L1,CH,C,WA(IW))
105 NA = 1-NA
goto 115
106 if (IP /= 3) goto 109
IX2 = IW+IDOT
if (NA /= 0) goto 107
CALL PASSB3 (IDOT,L1,C,CH,WA(IW),WA(IX2))
goto 108
107 CALL PASSB3 (IDOT,L1,CH,C,WA(IW),WA(IX2))
108 NA = 1-NA
goto 115
109 if (IP /= 5) goto 112
IX2 = IW+IDOT
IX3 = IX2+IDOT
IX4 = IX3+IDOT
if (NA /= 0) goto 110
CALL PASSB5 (IDOT,L1,C,CH,WA(IW),WA(IX2),WA(IX3),WA(IX4))
goto 111
110 CALL PASSB5 (IDOT,L1,CH,C,WA(IW),WA(IX2),WA(IX3),WA(IX4))
111 NA = 1-NA
goto 115
112 if (NA /= 0) goto 113
CALL PASSB (NAC,IDOT,IP,L1,IDL1,C,C,C,CH,CH,WA(IW))
goto 114
113 CALL PASSB (NAC,IDOT,IP,L1,IDL1,CH,CH,CH,C,C,WA(IW))
114 if (NAC /= 0) NA = 1-NA
115 L1 = L2
IW = IW+(IP-1)*IDOT
116 continue
if (NA == 0) return
N2 = N+N
DO 117 I=1,N2
C(I) = CH(I)
117 continue
return
END
subroutine PASSB (NAC,IDO,IP,L1,IDL1,CC,C1,C2,CH,CH2,WA)
DIMENSION CH(IDO,L1,IP) ,CC(IDO,IP,L1), &
C1(IDO,L1,IP) ,WA(1) ,C2(IDL1,IP), &
CH2(IDL1,IP)
IDOT = IDO/2
NT = IP*IDL1
IPP2 = IP+2
IPPH = (IP+1)/2
IDP = IP*IDO
!
if (IDO < L1) goto 106
DO 103 J=2,IPPH
JC = IPP2-J
DO 102 K=1,L1
DO 101 I=1,IDO
CH(I,K,J) = CC(I,J,K)+CC(I,JC,K)
CH(I,K,JC) = CC(I,J,K)-CC(I,JC,K)
101 continue
102 continue
103 continue
DO 105 K=1,L1
DO 104 I=1,IDO
CH(I,K,1) = CC(I,1,K)
104 continue
105 continue
goto 112
106 DO 109 J=2,IPPH
JC = IPP2-J
DO 108 I=1,IDO
DO 107 K=1,L1
CH(I,K,J) = CC(I,J,K)+CC(I,JC,K)
CH(I,K,JC) = CC(I,J,K)-CC(I,JC,K)
107 continue
108 continue
109 continue
DO 111 I=1,IDO
DO 110 K=1,L1
CH(I,K,1) = CC(I,1,K)
110 continue
111 continue
112 IDL = 2-IDO
INC = 0
DO 116 L=2,IPPH
LC = IPP2-L
IDL = IDL+IDO
DO 113 IK=1,IDL1
C2(IK,L) = CH2(IK,1)+WA(IDL-1)*CH2(IK,2)
C2(IK,LC) = WA(IDL)*CH2(IK,IP)
113 continue
IDLJ = IDL
INC = INC+IDO
DO 115 J=3,IPPH
JC = IPP2-J
IDLJ = IDLJ+INC
if (IDLJ > IDP) IDLJ = IDLJ-IDP
WAR = WA(IDLJ-1)
WAI = WA(IDLJ)
DO 114 IK=1,IDL1
C2(IK,L) = C2(IK,L)+WAR*CH2(IK,J)
C2(IK,LC) = C2(IK,LC)+WAI*CH2(IK,JC)
114 continue
115 continue
116 continue
DO 118 J=2,IPPH
DO 117 IK=1,IDL1
CH2(IK,1) = CH2(IK,1)+CH2(IK,J)
117 continue
118 continue
DO 120 J=2,IPPH
JC = IPP2-J
DO 119 IK=2,IDL1,2
CH2(IK-1,J) = C2(IK-1,J)-C2(IK,JC)
CH2(IK-1,JC) = C2(IK-1,J)+C2(IK,JC)
CH2(IK,J) = C2(IK,J)+C2(IK-1,JC)
CH2(IK,JC) = C2(IK,J)-C2(IK-1,JC)
119 continue
120 continue
NAC = 1
if (IDO == 2) return
NAC = 0
DO 121 IK=1,IDL1
C2(IK,1) = CH2(IK,1)
121 continue
DO 123 J=2,IP
DO 122 K=1,L1
C1(1,K,J) = CH(1,K,J)
C1(2,K,J) = CH(2,K,J)
122 continue
123 continue
if (IDOT > L1) goto 127
IDIJ = 0
DO 126 J=2,IP
IDIJ = IDIJ+2
DO 125 I=4,IDO,2
IDIJ = IDIJ+2
DO 124 K=1,L1
C1(I-1,K,J) = WA(IDIJ-1)*CH(I-1,K,J)-WA(IDIJ)*CH(I,K,J)
C1(I,K,J) = WA(IDIJ-1)*CH(I,K,J)+WA(IDIJ)*CH(I-1,K,J)
124 continue
125 continue
126 continue
return
127 IDJ = 2-IDO
DO 130 J=2,IP
IDJ = IDJ+IDO
DO 129 K=1,L1
IDIJ = IDJ
DO 128 I=4,IDO,2
IDIJ = IDIJ+2
C1(I-1,K,J) = WA(IDIJ-1)*CH(I-1,K,J)-WA(IDIJ)*CH(I,K,J)
C1(I,K,J) = WA(IDIJ-1)*CH(I,K,J)+WA(IDIJ)*CH(I-1,K,J)
128 continue
129 continue
130 continue
return
END
subroutine PASSB2 (IDO,L1,CC,CH,WA1)
DIMENSION CC(IDO,2,L1) ,CH(IDO,L1,2), &
WA1(1)
if (IDO > 2) goto 102
DO 101 K=1,L1
CH(1,K,1) = CC(1,1,K)+CC(1,2,K)
CH(1,K,2) = CC(1,1,K)-CC(1,2,K)
CH(2,K,1) = CC(2,1,K)+CC(2,2,K)
CH(2,K,2) = CC(2,1,K)-CC(2,2,K)
101 continue
return
102 DO 104 K=1,L1
DO 103 I=2,IDO,2
CH(I-1,K,1) = CC(I-1,1,K)+CC(I-1,2,K)
TR2 = CC(I-1,1,K)-CC(I-1,2,K)
CH(I,K,1) = CC(I,1,K)+CC(I,2,K)
TI2 = CC(I,1,K)-CC(I,2,K)
CH(I,K,2) = WA1(I-1)*TI2+WA1(I)*TR2
CH(I-1,K,2) = WA1(I-1)*TR2-WA1(I)*TI2
103 continue
104 continue
return
END
subroutine PASSB3 (IDO,L1,CC,CH,WA1,WA2)
DIMENSION CC(IDO,3,L1) ,CH(IDO,L1,3), &
WA1(1) ,WA2(1)
DATA TAUR,TAUI /-.5,.866025403784439/
if (IDO /= 2) goto 102
DO 101 K=1,L1
TR2 = CC(1,2,K)+CC(1,3,K)
CR2 = CC(1,1,K)+TAUR*TR2
CH(1,K,1) = CC(1,1,K)+TR2
TI2 = CC(2,2,K)+CC(2,3,K)
CI2 = CC(2,1,K)+TAUR*TI2
CH(2,K,1) = CC(2,1,K)+TI2
CR3 = TAUI*(CC(1,2,K)-CC(1,3,K))
CI3 = TAUI*(CC(2,2,K)-CC(2,3,K))
CH(1,K,2) = CR2-CI3
CH(1,K,3) = CR2+CI3
CH(2,K,2) = CI2+CR3
CH(2,K,3) = CI2-CR3
101 continue
return
102 DO 104 K=1,L1
DO 103 I=2,IDO,2
TR2 = CC(I-1,2,K)+CC(I-1,3,K)
CR2 = CC(I-1,1,K)+TAUR*TR2
CH(I-1,K,1) = CC(I-1,1,K)+TR2
TI2 = CC(I,2,K)+CC(I,3,K)
CI2 = CC(I,1,K)+TAUR*TI2
CH(I,K,1) = CC(I,1,K)+TI2
CR3 = TAUI*(CC(I-1,2,K)-CC(I-1,3,K))
CI3 = TAUI*(CC(I,2,K)-CC(I,3,K))
DR2 = CR2-CI3
DR3 = CR2+CI3
DI2 = CI2+CR3
DI3 = CI2-CR3
CH(I,K,2) = WA1(I-1)*DI2+WA1(I)*DR2
CH(I-1,K,2) = WA1(I-1)*DR2-WA1(I)*DI2
CH(I,K,3) = WA2(I-1)*DI3+WA2(I)*DR3
CH(I-1,K,3) = WA2(I-1)*DR3-WA2(I)*DI3
103 continue
104 continue
return
END
subroutine PASSB4 (IDO,L1,CC,CH,WA1,WA2,WA3)
DIMENSION CC(IDO,4,L1) ,CH(IDO,L1,4), &
WA1(1) ,WA2(1) ,WA3(1)
if (IDO /= 2) goto 102
DO 101 K=1,L1
TI1 = CC(2,1,K)-CC(2,3,K)
TI2 = CC(2,1,K)+CC(2,3,K)
TR4 = CC(2,4,K)-CC(2,2,K)
TI3 = CC(2,2,K)+CC(2,4,K)
TR1 = CC(1,1,K)-CC(1,3,K)
TR2 = CC(1,1,K)+CC(1,3,K)
TI4 = CC(1,2,K)-CC(1,4,K)
TR3 = CC(1,2,K)+CC(1,4,K)
CH(1,K,1) = TR2+TR3
CH(1,K,3) = TR2-TR3
CH(2,K,1) = TI2+TI3
CH(2,K,3) = TI2-TI3
CH(1,K,2) = TR1+TR4
CH(1,K,4) = TR1-TR4
CH(2,K,2) = TI1+TI4
CH(2,K,4) = TI1-TI4
101 continue
return
102 DO 104 K=1,L1
DO 103 I=2,IDO,2
TI1 = CC(I,1,K)-CC(I,3,K)
TI2 = CC(I,1,K)+CC(I,3,K)
TI3 = CC(I,2,K)+CC(I,4,K)
TR4 = CC(I,4,K)-CC(I,2,K)
TR1 = CC(I-1,1,K)-CC(I-1,3,K)
TR2 = CC(I-1,1,K)+CC(I-1,3,K)
TI4 = CC(I-1,2,K)-CC(I-1,4,K)
TR3 = CC(I-1,2,K)+CC(I-1,4,K)
CH(I-1,K,1) = TR2+TR3
CR3 = TR2-TR3
CH(I,K,1) = TI2+TI3
CI3 = TI2-TI3
CR2 = TR1+TR4
CR4 = TR1-TR4
CI2 = TI1+TI4
CI4 = TI1-TI4
CH(I-1,K,2) = WA1(I-1)*CR2-WA1(I)*CI2
CH(I,K,2) = WA1(I-1)*CI2+WA1(I)*CR2
CH(I-1,K,3) = WA2(I-1)*CR3-WA2(I)*CI3
CH(I,K,3) = WA2(I-1)*CI3+WA2(I)*CR3
CH(I-1,K,4) = WA3(I-1)*CR4-WA3(I)*CI4
CH(I,K,4) = WA3(I-1)*CI4+WA3(I)*CR4
103 continue
104 continue
return
END
subroutine PASSB5 (IDO,L1,CC,CH,WA1,WA2,WA3,WA4)
DIMENSION CC(IDO,5,L1) ,CH(IDO,L1,5), &
WA1(1) ,WA2(1) ,WA3(1) ,WA4(1)
DATA TR11,TI11,TR12,TI12 /.309016994374947,.951056516295154, &
-.809016994374947,.587785252292473/
if (IDO /= 2) goto 102
DO 101 K=1,L1
TI5 = CC(2,2,K)-CC(2,5,K)
TI2 = CC(2,2,K)+CC(2,5,K)
TI4 = CC(2,3,K)-CC(2,4,K)
TI3 = CC(2,3,K)+CC(2,4,K)
TR5 = CC(1,2,K)-CC(1,5,K)
TR2 = CC(1,2,K)+CC(1,5,K)
TR4 = CC(1,3,K)-CC(1,4,K)
TR3 = CC(1,3,K)+CC(1,4,K)
CH(1,K,1) = CC(1,1,K)+TR2+TR3
CH(2,K,1) = CC(2,1,K)+TI2+TI3
CR2 = CC(1,1,K)+TR11*TR2+TR12*TR3
CI2 = CC(2,1,K)+TR11*TI2+TR12*TI3
CR3 = CC(1,1,K)+TR12*TR2+TR11*TR3
CI3 = CC(2,1,K)+TR12*TI2+TR11*TI3
CR5 = TI11*TR5+TI12*TR4
CI5 = TI11*TI5+TI12*TI4
CR4 = TI12*TR5-TI11*TR4
CI4 = TI12*TI5-TI11*TI4
CH(1,K,2) = CR2-CI5
CH(1,K,5) = CR2+CI5
CH(2,K,2) = CI2+CR5
CH(2,K,3) = CI3+CR4
CH(1,K,3) = CR3-CI4
CH(1,K,4) = CR3+CI4
CH(2,K,4) = CI3-CR4
CH(2,K,5) = CI2-CR5
101 continue
return
102 DO 104 K=1,L1
DO 103 I=2,IDO,2
TI5 = CC(I,2,K)-CC(I,5,K)
TI2 = CC(I,2,K)+CC(I,5,K)
TI4 = CC(I,3,K)-CC(I,4,K)
TI3 = CC(I,3,K)+CC(I,4,K)
TR5 = CC(I-1,2,K)-CC(I-1,5,K)
TR2 = CC(I-1,2,K)+CC(I-1,5,K)
TR4 = CC(I-1,3,K)-CC(I-1,4,K)
TR3 = CC(I-1,3,K)+CC(I-1,4,K)
CH(I-1,K,1) = CC(I-1,1,K)+TR2+TR3
CH(I,K,1) = CC(I,1,K)+TI2+TI3
CR2 = CC(I-1,1,K)+TR11*TR2+TR12*TR3
CI2 = CC(I,1,K)+TR11*TI2+TR12*TI3
CR3 = CC(I-1,1,K)+TR12*TR2+TR11*TR3
CI3 = CC(I,1,K)+TR12*TI2+TR11*TI3
CR5 = TI11*TR5+TI12*TR4
CI5 = TI11*TI5+TI12*TI4
CR4 = TI12*TR5-TI11*TR4
CI4 = TI12*TI5-TI11*TI4
DR3 = CR3-CI4
DR4 = CR3+CI4
DI3 = CI3+CR4
DI4 = CI3-CR4
DR5 = CR2+CI5
DR2 = CR2-CI5
DI5 = CI2-CR5
DI2 = CI2+CR5
CH(I-1,K,2) = WA1(I-1)*DR2-WA1(I)*DI2
CH(I,K,2) = WA1(I-1)*DI2+WA1(I)*DR2
CH(I-1,K,3) = WA2(I-1)*DR3-WA2(I)*DI3
CH(I,K,3) = WA2(I-1)*DI3+WA2(I)*DR3
CH(I-1,K,4) = WA3(I-1)*DR4-WA3(I)*DI4
CH(I,K,4) = WA3(I-1)*DI4+WA3(I)*DR4
CH(I-1,K,5) = WA4(I-1)*DR5-WA4(I)*DI5
CH(I,K,5) = WA4(I-1)*DI5+WA4(I)*DR5
103 continue
104 continue
return
END
subroutine CFFTI (N,WSAVE)
DIMENSION WSAVE(1)
if (N == 1) return
IW1 = N+N+1
IW2 = IW1+N+N
CALL CFFTI1 (N,WSAVE(IW1),WSAVE(IW2))
return
END
subroutine CFFTI1 (N,WA,IFAC)
DIMENSION WA(1) ,IFAC(1) ,NTRYH(4)
DATA NTRYH(1),NTRYH(2),NTRYH(3),NTRYH(4)/3,4,2,5/
NL = N
NF = 0
J = 0
101 J = J+1
if (J-4) 102,102,103
102 NTRY = NTRYH(J)
goto 104
103 NTRY = NTRY+2
104 NQ = NL/NTRY
NR = NL-NTRY*NQ
if (NR) 101,105,101
105 NF = NF+1
IFAC(NF+2) = NTRY
NL = NQ
if (NTRY /= 2) goto 107
if (NF == 1) goto 107
DO 106 I=2,NF
IB = NF-I+2
IFAC(IB+2) = IFAC(IB+1)
106 continue
IFAC(3) = 2
107 if (NL /= 1) goto 104
IFAC(1) = N
IFAC(2) = NF
TPI = 6.28318530717959