-
Notifications
You must be signed in to change notification settings - Fork 1
/
libtracks_tools.ftn
3963 lines (3440 loc) · 121 KB
/
libtracks_tools.ftn
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
#include <tracks_cte.h>
c***********************************************************************
c---- Sous-routine GLOBAL_OR_LAM ---
c***********************************************************************
subroutine global_or_lam(global, lon, ni, nj)
implicit none
integer ni, nj
real lon(ni,nj)
integer global
real dx, next_lon
integer j
j = nj/2
dx = lon(2,j) - lon(1,j)
next_lon = lon(ni,j) + dx
if ( next_lon .ge. 360. ) next_lon = next_lon - 360.
if ( ( next_lon .ge. lon(1,j)-dx/4 .and.
$ next_lon .le. lon(1,j)+dx/4 ) .or.
$ ( next_lon .ge. lon(2,j)-dx/4 .and.
$ next_lon .le. lon(2,j)+dx/4 )) then
if ( lon(1,j) .eq. lon(ni,j) ) then
global = 2
print*
print *,' *** Grille globale ou i = ni (repetition) '
else
global = 1
print*
print *,' *** Grille globale ou i /= ni '
end if
else
global = 0
print*
print *,' *** Grille a aire limitee '
end if
end subroutine global_or_lam
c***********************************************************************
c---- Sous-routine GRAD_VORT ---
c***********************************************************************
#ifndef NETCDF
subroutine grad_vort(VORT, h, delta_x, delta_y, slat, gdid,
$ ni, nj)
#else
subroutine grad_vort(VORT, h, delta_x, delta_y, slat, slon,
$ ni, nj)
#endif
implicit none
real vort(ni,nj), h(ni,nj), slat(ni,nj)
real ug(ni,nj), vg(ni,nj), wind_spd(ni,nj)
real work1(ni,nj), work2(ni,nj), delta_x(ni,nj), delta_y(ni,nj)
real fc(ni,nj)
real slatt
real ome
integer ni, nj
#ifndef NETCDF
integer gdid
#else
real slon(ni,nj)
#endif
integer i, j
ome=K_OME
! H en entree doit etre en DAM !
DO J=1,NJ
DO I=1,NI
SLATT=SLAT(I,J)
IF(ABS(SLAT(I,J)).LT.20.) THEN
IF(SLAT(I,J).LE.0.)SLATT=-20.
IF(SLAT(I,J).GT.0.)SLATT=20.
ENDIF
FC(I,J)=2.*OME*SIN(SLATT*K_DEGTOR) ! coriolis parameter
END DO
END DO
call derhn_gr(VG,H,'x',10.0*K_GRAVITY,delta_x,delta_y,ni,nj)
call derhn_gr(UG,H,'y',10.0*K_GRAVITY,delta_x,delta_y,ni,nj)
#ifndef NETCDF
CALL GRADIENT_WIND_FACTOR(H,WORK1,delta_x,delta_y,gdid,ni,nj)
#else
CALL GRADIENT_WIND_FACTOR(H,WORK1,delta_x,delta_y,slat,slon,ni,nj)
#endif
VG(:,:)= VG(:,:)*WORK1(:,:)/FC(:,:)
UG(:,:)=-UG(:,:)*WORK1(:,:)/FC(:,:)
WIND_SPD(:,:)=SQRT(UG(:,:)**2+VG(:,:)**2)
call derhn_gr(WORK1,UG,'y',1.0,delta_x,delta_y,ni,nj)
call derhn_gr(WORK2,VG,'x',1.0,delta_x,delta_y,ni,nj)
VORT(:,:)=WORK2(:,:)-WORK1(:,:)
end subroutine grad_vort
c***********************************************************************
c---- Sous-routine VORT ---
c***********************************************************************
subroutine vorticity(VORT, u, v, delta_x, delta_y, ni, nj)
implicit none
real vort(ni,nj)
real u(ni,nj), v(ni,nj)
real wk1(ni,nj), wk2(ni,nj), delta_x(ni,nj), delta_y(ni,nj)
integer ni, nj
call derhn_gr(wk1, u, 'y', 1.0, delta_x, delta_y, ni, nj)
call derhn_gr(wk2, v, 'x', 1.0, delta_x, delta_y, ni, nj)
call gdadgd(vort, wk2, wk1, 1.0, -1.0, ni, nj, 0)
end subroutine vorticity
c***********************************************************************
c---- Sous-routine THERM_WIND ---
c***********************************************************************
subroutine therm_wind(UT, VT, DZ, delta_x, delta_y, lat, ni, nj)
implicit none
real ut(ni,nj), vt(ni,nj)
real dz(ni,nj), lat(ni,nj)
real wk1(ni,nj), wk2(ni,nj), delta_x(ni,nj), delta_y(ni,nj)
real fc, ome, latt
integer ni, nj, i, j
! DZ en entree doit etre en DAM !
ome=K_OME
call derhn_gr(wk1, dz, 'y', 10.0*K_GRAVITY, delta_x, delta_y,
$ ni, nj)
call derhn_gr(wk2, dz, 'x', 10.0*K_GRAVITY, delta_x, delta_y,
$ ni, nj)
DO J=1,NJ
DO I=1,NI
LATT=LAT(I,J)
IF(ABS(LAT(I,J)).LT.20.) THEN
IF(LAT(I,J).LE.0.)LATT=-20.
IF(LAT(I,J).GT.0.)LATT=20.
ENDIF
FC=2.*OME*SIN(LATT*K_DEGTOR) ! coriolis parameter
ut(i,j) = -1.0 * wk1(i,j) / fc
vt(i,j) = wk2(i,j) / fc
END DO
END DO
end subroutine therm_wind
c***********************************************************************
c---- Sous-routine CLOSEST_CENTRE ---
c***********************************************************************
subroutine closest_centre(DIST_MIN, F_LAT, F_LON, c_lat, c_lon,
$ x_lat, x_lon, nx, nc)
implicit none
integer nc, nx
real dist_min, c_lat, c_lon
real x_lat(nc), x_lon(nc)
real daz, dist, f_lat, f_lon
integer x
dist_min = 1.E30 ! metres
f_lat = 0.0
f_lon = 0.0
do x = 1, nx
call gcpazd(x_lat(x), x_lon(x), c_lat, c_lon, daz, dist)
if (dist .lt. dist_min) then
dist_min = dist
f_lat = x_lat(x)
f_lon = x_lon(x)
end if
end do
end subroutine closest_centre
c***********************************************************************
c---- Sous-routine SURROUND_MINMAX ---
c***********************************************************************
subroutine surround_minmax(MIN, MAX, c_i, c_j, c_lat, c_lon,
$ field, radius, delta_x, delta_y, lat, lon, ni, nj)
implicit none
integer ni, nj, c_i, c_j, c_lat, c_lon
real field(ni,nj), delta_x(ni,nj), delta_y(ni,nj)
real lat(ni,nj), lon(ni,nj)
real min, max, radius, daz, dist
integer i, j, point_num_x, point_num_y
max = -9.E19
min = +9.E19
point_num_x = nint( radius * 1.E3 / delta_x(c_i,c_j) ) + 1
point_num_y = nint( radius * 1.E3 / delta_y(c_i,c_j) ) + 1
! print*, c_i, c_j, radius, delta_x(c_i,c_j), delta_y(c_i,c_j),
! $ point_num_x, point_num_y
do i = c_i - point_num_x, c_i + point_num_x
do j = c_j - point_num_y, c_j + point_num_y
if (i .lt. 1 ) cycle
if (i .gt. ni) cycle
if (j .lt. 1 ) cycle
if (j .gt. nj) cycle
call GCPAZD(float(c_lat), float(c_lon), lat(i,j), lon(i,j),
$ daz, DIST)
if ( dist .le. radius * 1.E3 ) then
if ( field(i,j) .gt. max ) then
max = field(i,j)
end if
if ( field(i,j) .lt. min ) then
min = field(i,j)
end if
end if
end do
end do
end subroutine surround_minmax
c***********************************************************************
c---- Sous-routine SURROUND_MEAN ---
c***********************************************************************
subroutine surround_mean(MEAN, c_i, c_j, c_lat, c_lon, field,
$ radius, delta_x, delta_y, lat, lon, ni, nj)
implicit none
integer ni, nj, c_i, c_j
real c_lat, c_lon
real field(ni,nj), delta_x(ni,nj), delta_y(ni,nj)
real lat(ni,nj), lon(ni,nj)
real mean, radius
real dist, daz, count_inv
integer i, j, point_num_x, point_num_y
integer count
mean = 0.
point_num_x = nint( radius * 1.E3 / delta_x(c_i,c_j) ) + 1
point_num_y = nint( radius * 1.E3 / delta_y(c_i,c_j) ) + 1
count = 0
do i = c_i - point_num_x, c_i + point_num_x
do j = c_j - point_num_y, c_j + point_num_y
if (i .lt. 1 ) cycle
if (i .gt. ni) cycle
if (j .lt. 1 ) cycle
if (j .gt. nj) cycle
call GCPAZD(c_lat, c_lon, lat(i,j), lon(i,j), daz, DIST)
if ( dist .le. radius * 1.E3 ) then
count = count + 1
count_inv = 1.0 / real(count)
mean = ( real(count-1) * mean + field(i,j) ) * count_inv
end if
end do
end do
end subroutine surround_mean
c***********************************************************************
c---- Sous-routine SURROUND_WC_MEAN ---
c***********************************************************************
subroutine surround_wc_mean(WARM, COLD, c_i, c_j, c_lat, c_lon,
$ field, radius, tuvdir, delta_x, delta_y, lat, lon, ni, nj)
implicit none
integer ni, nj, c_i, c_j
real c_lat, c_lon
real field(ni,nj), delta_x(ni,nj), delta_y(ni,nj)
real lat(ni,nj), lon(ni,nj)
real warm, cold, radius, tuvdir
real dist, daz, icount_w, icount_c, dummy, dir, ddir
integer i, j, point_num_x, point_num_y
integer count_w, count_c
warm = 0.
cold = 0.
point_num_x = nint( radius * 1.E3 / delta_x(c_i,c_j) ) + 1
point_num_y = nint( radius * 1.E3 / delta_y(c_i,c_j) ) + 1
count_w = 0
count_c = 0
do i = c_i - point_num_x, c_i + point_num_x
do j = c_j - point_num_y, c_j + point_num_y
if (i .lt. 1 ) cycle
if (i .gt. ni) cycle
if (j .lt. 1 ) cycle
if (j .gt. nj) cycle
! Distance entre le centre du cyclone et le point de grille
call GCPAZD(c_lat, c_lon, lat(i,j), lon(i,j), daz, DIST)
if ( dist .le. radius * 1.E3 ) then
! Angle entre le point (i,j) et le point le plus pres du
! centre du cyclone (c_i,c_j)
call UVTOSD(float(i-c_i), float(j-c_j), dummy, DIR)
ddir = dir - tuvdir
if (ddir .lt. 0.0) ddir = ddir + 360.0
if (ddir .gt. 360.0) ddir = ddir - 360.0
if ( (ddir .ge. 0.0) .and.
$ (ddir .lt. 180.0) ) then ! Cote droit
count_w = count_w + 1
icount_w = 1.0 / real(count_w)
warm = ( real(count_w-1) * warm + field(i,j) ) *
$ icount_w
else
count_c = count_c + 1
icount_c = 1.0 / real(count_c)
cold = ( real(count_c-1) * cold + field(i,j) ) *
$ icount_c
end if
end if
end do
end do
end subroutine surround_wc_mean
#ifndef NETCDF
c***********************************************************************
c---- Sous-routine ECRIT_FST ---
c***********************************************************************
subroutine ecrit_fst(iun,buffer,nomvar,etiket,ip1,ip2,npas,
$ dateo,ni,nj,ig1,ig2,ig3,ig4,ip3,grtyp,typvar,deet)
implicit none
character*2 nomvar
character*1 typvar, grtyp
character*12 etiket
integer dateo, deet, npas, ni, nj, nk, npak, datyp
integer ip1, ip2, ip3
integer ig1, ig2, ig3, ig4, iun
integer fstecr, ier
real buffer(ni,nj)
real wke (ni,nj)
npak = -16
nk = 1
ier = fstecr(buffer,wke,npak,iun,dateo,deet,npas,ni,nj,nk,
$ ip1,ip2,ip3,typvar,nomvar,etiket,grtyp,ig1,ig2,ig3,
$ ig4,1,.false.)
return
end
#else
c***********************************************************************
c---- Sous-routine ECRIT_NC ---
c***********************************************************************
subroutine ecrit_nc(iun,buffer,nomvar,t,nt,londimid,latdimid,
$ timedimid,ni,nj)
implicit none
include 'netcdf.inc'
integer ni, nj, nt
character*(*) nomvar
integer londimid, latdimid, timedimid
integer iun, t
real buffer(ni,nj)
integer ier
integer rhdims(3), start(3), count(3)
integer varid, ndims
character*256 strerr
c Define variable in header if needed
ier = nf_inq_varid(iun, nomvar, varid)
if (ier .ne. nf_noerr) then
c Variable does not already exists, must create it
ier = nf_redef(iun)
if (ier .ne. nf_noerr) then
write(strerr,*) 'ecrit_nc: Cannot go into NetCDF redefine mode for output file ID:',iun
call handle_err(ier,strerr)
endif
rhdims(1) = londimid
rhdims(2) = latdimid
if (nt .gt. 1) then
rhdims(3) = timedimid
ndims = 3
else
ndims = 2
endif
ier = nf_def_var(iun, nomvar, nf_float, ndims, rhdims, varid)
if (ier .ne. nf_noerr) then
write(strerr,*) 'ecrit_nc: Cannot define NetCDF variable:',nomvar,' number of dimensions:',ndims
call handle_err(ier,strerr)
endif
ier = nf_enddef(iun)
if (ier .ne. nf_noerr) then
write(strerr,*) 'ecrit_nc: Cannot close NetCDF redefine mode for output file ID:',iun
call handle_err(ier,strerr)
endif
endif
c Write variable onto disk
start(1) = 1
start(2) = 1
count(1) = ni
count(2) = nj
if (nt .gt. 1) then
start(3) = t
count(3) = 1
endif
ier = nf_put_vara_real(iun, varid, start, count, buffer)
if (ier .ne. nf_noerr) then
write(strerr,*) 'ecrit_nc: Cannot write NetCDF variable ',nomvar,' in output file ID:',iun
call handle_err(ier,strerr)
endif
return
end
#endif
c***********************************************************************
c---- Sous-routine DIFFDIR ---
c***********************************************************************
subroutine diffdir(ddir,dist1,dist2,clat,clon,nn,it,nt,nc,tlat,
$ tlon,tstep)
implicit none
c ----------------------------------------------------------
c Ecrite par JF Caron - Septembre 2010 -
c ----------------------------------------------------------
integer nt,nc,nn,it,tstep
real clat(nt,nc),clon(nt,nc),dist
real tlat,tlon,dlat,dlon,dir1,dir2,dist1,dist2,ddir
logical pass
c Direction entre le centre a l'etude et la derniere position
c de la trajectoire
call distance(dist1,dir1,dlon,dlat,tlat,tlon,clat(nn,it),
$ clon(nn,it),1.0/(tstep*1000.0))
c Derniere direction de la trajecoire
call distance(dist2,dir2,dlon,dlat,clat(nn,it),clon(nn,it),
$ clat(nn-1,it),clon(nn-1,it),1.0/(tstep*1000.0))
c Changement de direction implique par le centre a l'etude
ddir=abs(dir1-dir2)
if(ddir.gt.180.0)then
ddir=360.-ddir
endif
return
end
c***********************************************************************
c---- Sous-routine DISTANCE ---
c***********************************************************************
subroutine distance(DT, DIR, DX, DY, lat2, lon2, lat1, lon1,
$ fact)
implicit none
c ----------------------------------------------------------
c Ecrite par JF Caron - Septembre 2010 -
c ----------------------------------------------------------
c Calcul utilisant la Formule d'Haversine
c Reference: R.W. Sinnott,'Virtues of Haversine',Sky and Telescope,
c vol.68, no.2, 1984, p.159)
c Les distance sont en m
real lat1, lon1, lat2, lon2, dx, dy ,dt, dir, fact, dt2
real dlat, dlon, a, c, R
c Rayon moyen de la Terre
R = K_EARTHR
*
dlon = lon2 - lon1
dlat = 0.0
a = (sin(dlat/2.0*K_DEGTOR))**2 +
$ cos(lat1*K_DEGTOR)*cos(lat2*K_DEGTOR)*
$ (sin(dlon/2.0*K_DEGTOR))**2
c = 2.0 * atan2(sqrt(a),sqrt(1-a))
dx= R * c * fact
! On suppose que les longitudes sont ordonnes de 0 a 359 !!!
if( (dlon.lt.0.0 .and. dlon.gt.-180.0) .or.
$ dlon.gt.180.0 ) then
dx=-1.0*dx
endif
*
dlon = 0.0
dlat = lat2 - lat1
a = (sin(dlat/2.0*K_DEGTOR))**2 +
$ cos(lat1*K_DEGTOR)*cos(lat2*K_DEGTOR)*
$ (sin(dlon/2.0*K_DEGTOR))**2
c = 2.0 * atan2(sqrt(a),sqrt(1-a))
dy= R * c * fact
if(dlat.lt.0.0)then
dy=-1.0*dy
endif
*
dlon = lon2 - lon1
dlat = lat2 - lat1
a = (sin(dlat/2.0*K_DEGTOR))**2 +
$ cos(lat1*K_DEGTOR)*cos(lat2*K_DEGTOR)*
$ (sin(dlon/2.0*K_DEGTOR))**2
c = 2.0 * atan2(sqrt(a),sqrt(1-a))
dt= R * c * fact
*
dt2=sqrt(dx**2+dy**2)
if(dt2.ne.0.0)then
dir=acos(dy/dt2)/K_DEGTOR
if(dx.lt.0.0)then
dir=360.-dir
endif
else
dir=0.0
endif
*
return
end
c***********************************************************************
c---- Sous-routine GCPAZD ---
c***********************************************************************
SUBROUTINE GCPAZD(DLAT1,DLON1,DLAT2,DLON2,DAZ12,DIST12)
implicit none
! DESCRIPTION
!
C TO COMPUTE THE GREAT-CIRCLE PATH AZIMUTH (AZ12) AND
C DISTANCE (DIST12) BETWEEN THE POINTS (DLAT1,DLON1)
C AND (DLAT2,DLON2)
C
C ALL ANGLES AND LAT/LON ARE IN DEGREES IF THEY START WITH D,
C IN RADIANS IF THEY START WITH R.
c Arguments
REAL DLAT1 ! Input latitude of first point
REAL DLON1 ! Input longitude of first point
REAL DLAT2 ! Input latitude of second point
REAL DLON2 ! Input longitude of second point
REAL DAZ12 ! Output azimuth from 1 to 2 (degrees true)
REAL DIST12 ! Output distance from 1 to 2 (m)
c Local Variables
real EARTHR, DEGTOR, RLAT1, RLAT2, RLON1, RLON2, COS12
real SRLAT1, SRLAT2, CRLAT1, CLON12, COSBET, BETA, RAZ12
! ----------------------------------------------------------
! AUTHORS: Neil Gordon
! CREATION DATE: 10-Jun-1988 (operational version)
! ----------------------------------------------------------
EARTHR=K_EARTHR
DEGTOR=K_DEGTOR
RLAT1=DLAT1*DEGTOR
RLAT2=DLAT2*DEGTOR
RLON1=DLON1*DEGTOR
RLON2=DLON2*DEGTOR
C
IF (DLAT1.EQ.DLAT2 .AND. DLON1.EQ.DLON2) then
DAZ12 = 0.0
DIST12 = 0.0
else
C
SRLAT1=SIN(RLAT1)
SRLAT2=SIN(RLAT2)
CRLAT1=COS(RLAT1)
C
CLON12=COS(RLON1-RLON2)
C
C CALCULATE ANGLE SUBTENDED BETWEEN POINTS
C
COSBET=SRLAT1*SRLAT2+CRLAT1*COS(RLAT2)*CLON12
IF(COSBET.LT.-1.) COSBET=-1. !IN CASE OF ROUNDOFF
IF(COSBET.GT.+1.) COSBET=+1. !DITTO
BETA=ABS(ACOS(COSBET))
DIST12=EARTHR*BETA
C
if ( DIST12 == 0.0 ) then
! stop pour eviter une division par zero
DAZ12 = 0.0
return
end if
COS12=(SRLAT2-COSBET*SRLAT1)/(SIN(BETA)*COS(RLAT1))
IF(COS12.LT.-1.) COS12=-1.
IF(COS12.GT.+1.) COS12=+1. !IN CASE OF ROUNDOFF
RAZ12=ACOS(COS12)
C
DAZ12=RAZ12/DEGTOR
C
C TEST FOR AZIMUTH BEING 180-360, RATHER THAN 0-180
C
IF(SIN(RLON2-RLON1).LT.0.) DAZ12=360.-DAZ12
endif
return
END
c***********************************************************************
c---- Sous-routine GCPMOV ---
c***********************************************************************
SUBROUTINE GCPMOV(DLAT1,DLON1,DLAT2,DLON2,DAZ12,DIST12)
implicit none
! DESCRIPTION
C
C TO RETURN THE POINT DLAT2,DLON2 WHICH IS AT THE GIVEN
C BEARING AND DISTANCE FROM DLAT1,DLON1
C
C ALL ANGLES AND LAT/LON ARE IN DEGREES IF THEY START WITH D,
C IN RADIANS IF THEY START WITH R.
C
C ALL ANGLES AND LAT/LON ARE IN DEGREES IF THEY START WITH D,
C IN RADIANS IF THEY START WITH R.
! ARGUMENTS
REAL*4 DLAT1 ! Input latitude of first point
REAL*4 DLON1 ! Input longitude of first point
REAL*4 DLAT2 ! Output latitude of second point
REAL*4 DLON2 ! Output longitude of second point
REAL*4 DAZ12 ! Input azimuth from 1 to 2 (degrees true)
REAL*4 DIST12 ! Input distance from 1 to 2 (m)
c Local Variables
real EARTHR, DEGTOR, RLAT1, RLON1, R12, CR12, SR12
real SRLAT1, CRLAT1, RAZ12, CRAZ12, DLON12
! ----------------------------------------------------------
! AUTHORS: Neil Gordon
! CREATION DATE: 10-Jun-1988 (operational version)
! ----------------------------------------------------------
EARTHR=K_EARTHR
DEGTOR=K_DEGTOR
RLAT1=DLAT1*DEGTOR
SRLAT1=SIN(RLAT1)
CRLAT1=COS(RLAT1)
RLON1=DLON1*DEGTOR
R12=DIST12/EARTHR
CR12=COS(R12)
SR12=SIN(R12)
RAZ12=DAZ12*DEGTOR
CRAZ12=COS(RAZ12)
DLAT2=ASIN(SRLAT1*CR12+CRLAT1*SR12*CRAZ12)/DEGTOR
DLON12=ATAN2(SR12*SIN(RAZ12),CR12*CRLAT1-SR12*SRLAT1*CRAZ12)/
$ DEGTOR
DLON2=MOD(DLON1+DLON12+720.0,360.)
IF (DLON2.GT.180.0) DLON2=DLON2-360.
RETURN
END
c***********************************************************************
c---- Sous-routine GET_TIMES ---
c***********************************************************************
#ifndef NETCDF
subroutine get_times(year, month, day, hour, nstamp,
$ ntimes, nincr, hi, hf, iun, nomvar, ip1, ntw)
implicit none
integer ntw, iun, naj, hi, hf
integer year(ntw), month(ntw), day(ntw), hour(ntw)
character*2 nomvar
integer nstamp(ntw), ntimes, nincr
integer ier
integer dateo,deet,npas,ni,nj,nbits,datyp,nk,datev
integer ip1,ip2,ip3,swa,lng,dltf,ubc,ip1d
integer ig1,ig2,ig3,ig4,extra1,extra2,extra3
integer liste(ntw),infon,vip2(ntw)
character*1 typvar, grtyp
character*12 etiket
character*2 nomvar_f
integer,allocatable,dimension (:) :: ibwork
integer fstprm, fstinl, newdate
integer i, k, cmcstamp, curdate, curhour, nnincr, t
ier = fstinl(iun, ni, nj, nk, -1, ' ', ip1, -1, -1, ' ', nomvar,
$ liste, infon, ntw)
k = 1
t = 1
do i=1,infon
ier = fstprm(liste(i),dateo,deet,npas,ni,nj,nk,nbits,
$ datyp,ip1d,ip2,ip3,typvar,nomvar_f,etiket,grtyp,ig1,ig2,
$ ig3,ig4,swa,lng,dltf,ubc,extra1,extra2,extra3)
call incdat(cmcstamp, dateo, ip2)
call ins_unique_vect(nstamp, cmcstamp, k, ntw)
call ins_unique_vect(vip2, ip2, t, ntw)
enddo
nk = k - 1
if (nk.gt.0) then
allocate(ibwork(nk))
call c_sort(nstamp, ibwork, nk)
call c_sort(vip2, ibwork, nk)
deallocate (ibwork)
hi = vip2(1)
hf = vip2(nk)
else
print*
print*,'Nombre de temps egal a zero.'
print*,'Impossible de continuer.'
print*
stop
endif
call difdat(nstamp(1), nstamp(2), nincr)
nincr = abs(nincr)
ntimes = 0
do k=1,nk
ntimes = ntimes + 1
ier = newdate(nstamp(k), curdate, curhour, -3)
year(ntimes) = curdate / 10000
month(ntimes) = (curdate - (year(ntimes)*10000)) / 100
day(ntimes) = curdate - (year(ntimes)*10000) -
$ (month(ntimes)*100)
hour(ntimes) = curhour / 1000000
if (k.gt.1) then
call difdat(nstamp(k-1), nstamp(k), nnincr)
nnincr = abs(nnincr)
if (nnincr.ne.nincr) then
print*
print*,'Time increment is not constant in the datafile.'
print*,'It must be for make_tracks to work.'
print*,'Aborting...'
stop
endif
endif
enddo
IF(NINCR.LE.0.OR.NINCR.GT.24) then
print*
print*,'GET_TIMES: Invalid time increment'
stop
else
print*
print*,'Les donnees sont a toutes les :',NINCR,'heures'
print*
end if
IF(NTIMES.EQ.0.OR.NTIMES.GT.NTW) then
print*
print*,'GET_TIMES: NTIMES greater then NTW !!! Stop'
stop
end if
return
end
#else
subroutine get_times(yyyy, mm, dd, hh, nstamp,
$ ntimes, nincr, iun, nomvar, ntw)
USE datetime_module
implicit none
include 'netcdf.inc'
TYPE(datetime) :: dat_a, dat_b
TYPE(timedelta) :: dif_dat
real*8 prev_dif_dat_sec
integer ntw, iun, naj
integer yyyy(ntw), mm(ntw), dd(ntw), hh(ntw), mins(ntw)
real*8 secs(ntw)
character*(*) nomvar
integer ntimes, nincr
integer ier
integer liste(ntw),infon,vip2(ntw)
real*8 nstamp(ntw)
integer real_len
integer nomvarid, ndims, dimids(NF_MAX_VAR_DIMS)
character*256 tunits
character*256 strerr
integer i, k, cmcstamp, curdate, curhour, nnincr, t
integer,external :: get_calendar
c Get ID of nomvar time
ier = nf_inq_varid(iun, trim(nomvar), nomvarid)
if (ier .ne. nf_noerr) then
write(strerr,*) 'get_times: Cannot get NetCDF variable ID: ',
$ trim(nomvar),' NetCDF variable from input file ID:',iun
call handle_err(ier,strerr)
endif
c Get Dimension ID of time
ier = nf_inq_vardimid(iun, nomvarid, dimids)
if (ier .ne. nf_noerr) then
write(strerr,*) 'get_times: Cannot get NetCDF dimensions ID of variable: ',
$ trim(nomvar),' NetCDF variable from input file ID:',iun
call handle_err(ier,strerr)
endif
c Get Time Units
ier = nf_get_att_text(iun, nomvarid, "units", tunits)
if (ier .ne. nf_noerr) then
write(strerr,*) 'get_times: Cannot get NetCDF attribute units of variable: ',
$ trim(nomvar),' NetCDF variable from input file ID:',iun
call handle_err(ier,strerr)
endif
c Get number of time steps
ier = nf_inq_dimlen(iun, dimids(1), ntimes)
if (ier .ne. nf_noerr) then
write(strerr,*) 'get_times: Cannot dimension length of variable: ',trim(nomvar)
call handle_err(ier,strerr)
endif
c Read times
ier = nf_get_vara_double(iun, nomvarid, 1, ntimes, nstamp)
if (ier .ne. nf_noerr) then
write(strerr,*) 'get_times: Cannot read ',trim(nomvar),' NetCDF variable from input file ID:',iun
call handle_err(ier,strerr)
endif
print*, 'get_times: READ ',nomvar,' successful'
tunits = tunits(1:real_len(tunits))//char(0)
ier = get_calendar(yyyy, mm, dd, hh, mins, secs, tunits, nstamp, ntimes)
do i=1,ntimes
if (i.ge.2) then
dat_b = datetime(yyyy(i), mm(i), dd(i), hh(i), mins(i))
dif_dat = dat_b - dat_a
if (i.ge.3) then
if (dif_dat%total_seconds() .ne. prev_dif_dat_sec) then
print*
print*,'Time increment is not constant in the datafile: ',
$ dif_dat%total_seconds(),' vs ',prev_dif_dat_sec
print*,'It must be for make_tracks to work.'
print*,'Aborting...'
stop
endif
endif
prev_dif_dat_sec = dif_dat%total_seconds()
endif
dat_a = datetime(yyyy(i), mm(i), dd(i), hh(i), mins(i))
enddo
nincr = (dif_dat%total_seconds())/3600
if (ntimes.le.0) then
print*
print*,'Nombre de temps egal a zero.'
print*,'Impossible de continuer.'
print*
stop
endif
IF(NINCR.LE.0.OR.NINCR.GT.24) then
print*
print*,'GET_TIMES: Invalid time increment'
stop
else
print*
print*,'Les donnees sont a toutes les :',NINCR,'heures'
print*
end if
IF(NTIMES.EQ.0.OR.NTIMES.GT.NTW) then
print*
print*,'GET_TIMES: NTIMES greater then NTW !!! Stop'
stop
end if
return
end
#endif
c***********************************************************************
c---- Sous-routine GRADIENT_WIND_FACTOR ---
c***********************************************************************
#ifndef NETCDF
SUBROUTINE GRADIENT_WIND_FACTOR(H,F_GRAD,delta_x,delta_y,gdid,
$ ni,nj)
#else
SUBROUTINE GRADIENT_WIND_FACTOR(H,F_GRAD,delta_x,delta_y,lat,lon,
$ ni,nj)
#endif
implicit none
!
! DESCRIPTION
!
! Computes F_GRAD, the factor that geostrophic winds are multiplied by
! to get gradient wind
!
! Valid both hemispheres
!
! ARGUMENTS
INTEGER NI,NJ !defined below in COMMON /GRID/
REAL H(NI*NJ) !.rf.ra input geopotential field
REAL F_GRAD(NI*NJ) !.wf.ra ouput gradient wind factor
REAL delta_x(NI*NJ), delta_y(NI*NJ)
#ifndef NETCDF
integer gdid
#else
real lat(ni,nj), lon(ni,nj)
#endif
! Local variables
integer ier, nij, i, j, k, ihem, ij
#ifndef NETCDF
integer gdllfxy
#endif
real pi, ome, rd, slat, slon, fij, term1, term2, term3, denom
real vg, fk, v, ro
real dHdx(NI*NJ),dHdy(NI*NJ),d2Hdxdy(NI*NJ),d2Hdx2(NI*NJ)
real d2Hdy2(NI*NJ)
! ----------------------------------------------------------
! AUTHORS: M.R.SINCLAIR, based on Trenberth's T.I.C. 158 pp. 10-12
! CREATION DATE: 09-MAR-1989
! ----------------------------------------------------------
RD=K_DEGTOR !PI/180.
OME=K_OME
NIJ=NI*NJ
IF(NIJ.EQ.0) then
STOP 'GRADIENT_WIND_FACTOR: CHECK NI,NJ'
endif
! H en entree doit etre en DAM !
call derhn_gr(dHdx,H,'x',10.0*K_GRAVITY,delta_x,delta_y,ni,nj)
call derhn_gr(dHdy,H,'y',10.0*K_GRAVITY,delta_x,delta_y,ni,nj)
call derhn_gr(d2Hdx2,dHdx,'x',1.0,delta_x,delta_y,ni,nj)
call derhn_gr(d2Hdxdy,dHdx,'y',1.0,delta_x,delta_y,ni,nj)
call derhn_gr(d2Hdxdy,dHdy,'x',1.0,delta_x,delta_y,ni,nj)
call derhn_gr(d2Hdy2,dHdy,'y',1.0,delta_x,delta_y,ni,nj)
DO J=1,NJ
DO I=1,NI
IJ=I+(J-1)*NI
#ifndef NETCDF
ier = gdllfxy(gdid, SLAT, SLON, I*1., J*1., 1)
#else
slat = lat(i,j)
slon = lon(i,j)
#endif
IF(SLAT.LE.0.0)THEN
FIJ=2.*OME*SIN(MIN(SLAT,-20.)*K_DEGTOR)
IHEM=-1
ELSE
FIJ=2.*OME*SIN(MAX(SLAT,20.)*K_DEGTOR)
IHEM=1
ENDIF
C
C USE FORMULAE FROM TRENBERTH'S TECH NOTE TIC 158 P. 11
C