-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathwwm_wind.F90
3137 lines (2952 loc) · 125 KB
/
wwm_wind.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
#include "wwm_functions.h"
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE INIT_WIND_INPUT
#ifdef NCDF
USE NETCDF
#endif
USE DATAPOOL
IMPLICIT NONE
INTEGER :: IT, IFILE
REAL(rkind) :: WDIRT
REAL(rkind) :: cf_w1, cf_w2
#ifdef MPI_PARALL_GRID
IF (MULTIPLE_IN_WIND) THEN
MNP_WIND=MNP
allocate(XP_WIND(MNP_WIND), YP_WIND(MNP_WIND), stat=istat)
IF (istat/=0) CALL WWM_ABORT('wwm_wind, allocate error 1')
XP_WIND=XP
YP_WIND=YP
ELSE
MNP_WIND=np_total
IF (myrank .eq. 0) THEN
allocate(XP_WIND(MNP_WIND), YP_WIND(MNP_WIND), stat=istat)
IF (istat/=0) CALL WWM_ABORT('wwm_wind, allocate error 1')
XP_WIND=XPtotal
YP_WIND=YPtotal
END IF
END IF
#else
MNP_WIND=MNP
allocate(XP_WIND(MNP_WIND), YP_WIND(MNP_WIND), stat=istat)
IF (istat/=0) CALL WWM_ABORT('wwm_wind, allocate error 1')
XP_WIND=XP
YP_WIND=YP
#endif
IF (LSTWD) THEN
IF (LCWIN) THEN
IF (PrintLOG) THEN
WRITE(WINDBG%FHNDL,'("+TRACE...",A)') 'HOMOGENOUS STEADY WIND FIELD IS USED'
WRITE(WINDBG%FHNDL,'("+TRACE...",A,I10)') 'WIND IS COMING FROM WWM - WINDFORMAT', IWINDFORMAT, LWDIR
FLUSH(WINDBG%FHNDL)
END IF
IF (LWDIR) THEN
CALL DEG2NAUT(WDIR, WDIRT, LNAUTIN)
WINDXY(:,1) = WVEL * COS(WDIRT * DEGRAD)
WINDXY(:,2) = WVEL * SIN(WDIRT * DEGRAD)
ELSE
WINDXY(:,1) = CWINDX
WINDXY(:,2) = CWINDY
END IF
ELSE ! LCWIN
IF (PrintLOG) THEN
WRITE(WINDBG%FHNDL,'("+TRACE...",A,I10)') 'WIND IS COMING FROM WWM - WINDFORMAT', IWINDFORMAT
WRITE(WINDBG%FHNDL,'("+TRACE...",A)') 'SPATIAL VARIABLE WIND FIELD IS USED'
FLUSH(WINDBG%FHNDL)
END IF
IF (IWINDFORMAT == 1) THEN
CALL CSEVAL( WIN%FHNDL, TRIM(WIN%FNAME), .FALSE., 2, WINDXY, MULTIPLE_IN_WIND)
#ifdef NCDF
ELSE IF (IWINDFORMAT == 2) THEN ! NETCDF created using ncl_convert2nc using DWD grib
CALL INIT_NETCDF_DWD
CALL FIND_WIND_NEAREST_LOWER_IDX(SEWI%BMJD, idxWind)
IFILE=WIND_TIME_IFILE(idxWind)
IT=WIND_TIME_IT(idxWind)
CALL READ_NETCDF_DWD(IFILE, IT, WINDXY)
ELSE IF (IWINDFORMAT == 3) THEN ! NETCDF created using cdo -f nc copy file.grb file.nc this is CFRS
CALL INIT_NETCDF_CRFS
CALL FIND_WIND_NEAREST_LOWER_IDX(SEWI%BMJD, idxWind)
IFILE=WIND_TIME_IFILE(idxWind)
IT=WIND_TIME_IT(idxWind)
CALL READ_NETCDF_CRFS(IFILE, IT, WINDXY)
ELSE IF (IWINDFORMAT == 4) THEN ! NETCDF NARR downloaded from NOMAD
CALL INIT_NETCDF_NARR
CALL FIND_WIND_NEAREST_LOWER_IDX(SEWI%BMJD, idxWind)
IFILE=WIND_TIME_IFILE(idxWind)
IT=WIND_TIME_IT(idxWind)
CALL READ_NETCDF_NARR(IFILE, IT, WINDXY)
ELSE IF (IWINDFORMAT == 5) THEN ! NETCDF CF_COMPLIANT STATIONARY FIELD
IF (PrintLOG) THEN
WRITE(WINDBG%FHNDL,'("+TRACE...",A)') 'COMPUTING CF INTERPOLATION COEFS AND LOADING WIND_TIME_MJD'
FLUSH(WINDBG%FHNDL)
END IF
CALL INIT_NETCDF_CF_WWM_WIND(eVAR_WIND)
ALLOCATE(tmp_wind1(MNP,2),tmp_wind2(MNP,2), stat=istat)
IF (istat/=0) CALL WWM_ABORT('wwm_wind, allocate error 1')
CALL GET_CF_TIME_INDEX(eVAR_WIND, REC1_wind_new,REC2_wind_new,cf_w1,cf_w2)
CALL READ_INTERP_NETCDF_CF_WWM_WIND(REC1_wind_new,tmp_wind1)
IF (cf_w1.NE.1) THEN
CALL READ_INTERP_NETCDF_CF_WWM_WIND(REC2_wind_new,tmp_wind2)
WINDXY(:,:) = cf_w1*tmp_wind1(:,:)+cf_w2*tmp_wind2(:,:)
ELSE
WINDXY(:,:) = cf_w1*tmp_wind1(:,:)
END IF
ELSE IF (IWINDFORMAT == 6) THEN ! DIRECT WWM forcing (no interp)
CALL INIT_DIRECT_NETCDF_CF(eVAR_WIND, MULTIPLE_IN_WIND, WIN%FNAME, "Uwind")
ALLOCATE(tmp_wind1(MNP,2),tmp_wind2(MNP,2), stat=istat)
IF (istat/=0) CALL WWM_ABORT('wwm_wind, allocate error 1')
CALL GET_CF_TIME_INDEX(eVAR_WIND, REC1_wind_new,REC2_wind_new,cf_w1,cf_w2)
CALL READ_DIRECT_NETCDF_CF(eVAR_wind, REC1_wind_new,tmp_wind1)
IF (cf_w1.NE.1) THEN
CALL READ_DIRECT_NETCDF_CF(eVAR_wind, REC2_wind_new,tmp_wind2)
WINDXY(:,:) = cf_w1*tmp_wind1(:,:)+cf_w2*tmp_wind2(:,:)
ELSE
WINDXY(:,:) = cf_w1*tmp_wind1(:,:)
END IF
#endif
#ifdef GRIB_API_ECMWF
ELSE IF (IWINDFORMAT == 7) THEN ! GRIB forcing from ecmwf
ALLOCATE(tmp_wind1(MNP,2),tmp_wind2(MNP,2), stat=istat)
IF (istat/=0) CALL WWM_ABORT('wwm_wind, allocate error 1')
CALL GET_CF_TIME_INDEX(eVAR_WIND, REC1_wind_new,REC2_wind_new,cf_w1,cf_w2)
CALL READ_GRIB_WIND(REC1_wind_new,tmp_wind1)
IF (cf_w1.NE.1) THEN
CALL READ_GRIB_WIND(REC2_wind_new,tmp_wind2)
WINDXY(:,:) = cf_w1*tmp_wind1(:,:)+cf_w2*tmp_wind2(:,:)
ELSE
WINDXY(:,:) = cf_w1*tmp_wind1(:,:)
END IF
#endif
ELSE
CALL wwm_abort('Wrong choice of IWINDFORMAT (maybe need NETCDF or GRIB)')
END IF
ENDIF
ELSE IF (LSEWD) THEN
IF (LCWIN) THEN
CALL wwm_abort('LSEWD + LCWIN NOT READY')
! CALL READ_WIND_TIME_SERIES(IT) ! set time according to wwminput.nml and get initial time step
! CALL SET_INITIAL_WIND(IT) !
ELSE
IF (PrintLOG) THEN
WRITE(WINDBG%FHNDL,'("+TRACE...",A,I10)') 'WIND IS COMING FROM WWM - WINDFORMAT', IWINDFORMAT
WRITE(WINDBG%FHNDL,'("+TRACE...",A)') 'NONSTATIONARY WIND FIELD IS USED '
FLUSH(WINDBG%FHNDL)
END IF
SEWI%TOTL = (SEWI%EMJD - SEWI%BMJD) * DAY2SEC
SEWI%ISTP = NINT( SEWI%TOTL / SEWI%DELT ) + 1
SEWI%TMJD = SEWI%BMJD
IF (PrintLOG) THEN
WRITE(WINDBG%FHNDL,*) SEWI%BEGT, SEWI%ENDT, SEWI%ISTP, SEWI%TOTL/3600.0, SEWI%DELT
WRITE(WINDBG%FHNDL,'("+TRACE...",A)') 'SPATIAL VARIABLE WIND FIELD IS USED'
WRITE(WINDBG%FHNDL,*) 'IWINDFORMAT=', IWINDFORMAT
FLUSH(WINDBG%FHNDL)
END IF
IF (IWINDFORMAT == 1) THEN
OPEN(WIN%FHNDL, FILE = TRIM(WIN%FNAME), STATUS = 'OLD')
CALL CSEVAL( WIN%FHNDL, TRIM(WIN%FNAME), .TRUE., 2, WINDXY, MULTIPLE_IN_WIND)
#ifdef NCDF
ELSE IF (IWINDFORMAT == 2) THEN ! NETCDF created using ncl_convert2nc using DWD grib
CALL INIT_NETCDF_DWD
CALL FIND_WIND_NEAREST_LOWER_IDX(MAIN%TMJD, idxWind)
IFILE=WIND_TIME_IFILE(idxWind)
IT=WIND_TIME_IT(idxWind)
CALL READ_NETCDF_DWD(IFILE, IT, WINDXY)
ELSE IF (IWINDFORMAT == 3) THEN ! NETCDF created using cdo -f nc copy file.grb file.nc
CALL INIT_NETCDF_CRFS
CALL FIND_WIND_NEAREST_LOWER_IDX(MAIN%TMJD, idxWind)
IFILE=WIND_TIME_IFILE(idxWind)
IT=WIND_TIME_IT(idxWind)
CALL READ_NETCDF_CRFS(IFILE, IT, WINDXY)
ELSE IF (IWINDFORMAT == 4) THEN ! NETCDF created using cdo -f nc copy file.grb file.nc
CALL INIT_NETCDF_NARR
CALL FIND_WIND_NEAREST_LOWER_IDX(MAIN%TMJD, idxWind)
IFILE=WIND_TIME_IFILE(idxWind)
IT=WIND_TIME_IT(idxWind)
CALL READ_NETCDF_NARR(IFILE, IT, WINDXY)
ELSE IF (IWINDFORMAT == 5) THEN
IF (PrintLOG) THEN
WRITE(WINDBG%FHNDL,'("+TRACE...",A)') 'SPATIAL/TEMPORAL VARIABLE WIND FIELD IS USED CF NETCDF'
WRITE(WINDBG%FHNDL,'("+TRACE...",A)') 'COMPUTING CF INTERPOLATION COEFS AND LOADING WIND_TIME_MJD'
FLUSH(WINDBG%FHNDL)
END IF
CALL INIT_NETCDF_CF_WWM_WIND(eVAR_WIND)
ALLOCATE(tmp_wind1(MNP,2), tmp_wind2(MNP,2), stat=istat)
IF (istat/=0) CALL WWM_ABORT('wwm_wind, allocate error 2')
CALL GET_CF_TIME_INDEX(eVAR_WIND, REC1_wind_new,REC2_wind_new,cf_w1,cf_w2)
CALL READ_INTERP_NETCDF_CF_WWM_WIND(REC1_wind_new,tmp_wind1)
IF (cf_w1.NE.1) THEN
CALL READ_INTERP_NETCDF_CF_WWM_WIND(REC2_wind_new,tmp_wind2)
WINDXY(:,:) = cf_w1*tmp_wind1(:,:)+cf_w2*tmp_wind2(:,:)
ELSE
WINDXY(:,:) = cf_w1*tmp_wind1(:,:)
END IF
ELSE IF (IWINDFORMAT == 6) THEN
IF (PrintLOG) THEN
WRITE(WINDBG%FHNDL,'("+TRACE...",A)') 'SPATIAL/TEMPORAL VARIABLE WIND FIELD IS USED CF NETCDF'
WRITE(WINDBG%FHNDL,'("+TRACE...",A)') 'COMPUTING CF INTERPOLATION COEFS AND LOADING WIND_TIME_MJD'
FLUSH(WINDBG%FHNDL)
END IF
CALL INIT_DIRECT_NETCDF_CF(eVAR_WIND, MULTIPLE_IN_WIND, WIN%FNAME, "Uwind")
ALLOCATE(tmp_wind1(MNP,2), tmp_wind2(MNP,2), stat=istat)
IF (istat/=0) CALL WWM_ABORT('wwm_wind, allocate error 2')
CALL GET_CF_TIME_INDEX(eVAR_WIND, REC1_wind_new,REC2_wind_new,cf_w1,cf_w2)
CALL READ_DIRECT_NETCDF_CF(eVAR_wind, REC1_wind_new,tmp_wind1)
IF (cf_w1.NE.1) THEN
CALL READ_DIRECT_NETCDF_CF(eVAR_wind, REC2_wind_new,tmp_wind2)
WINDXY(:,:) = cf_w1*tmp_wind1(:,:)+cf_w2*tmp_wind2(:,:)
ELSE
WINDXY(:,:) = cf_w1*tmp_wind1(:,:)
END IF
#endif
#ifdef GRIB_API_ECMWF
ELSE IF (IWINDFORMAT == 7) THEN
CALL INIT_GRIB_WIND
ALLOCATE(tmp_wind1(MNP,2), tmp_wind2(MNP,2), stat=istat)
IF (istat/=0) CALL WWM_ABORT('wwm_wind, allocate error 2')
CALL GET_CF_TIME_INDEX(eVAR_WIND, REC1_wind_new,REC2_wind_new,cf_w1,cf_w2)
CALL READ_GRIB_WIND(REC1_wind_new,tmp_wind1)
IF (cf_w1.NE.1) THEN
CALL READ_GRIB_WIND(REC2_wind_new,tmp_wind2)
WINDXY(:,:) = cf_w1*tmp_wind1(:,:)+cf_w2*tmp_wind2(:,:)
ELSE
WINDXY(:,:) = cf_w1*tmp_wind1(:,:)
END IF
#endif
ENDIF
ENDIF
ENDIF
IF (PrintLOG) THEN
write(WINDBG%FHNDL,'("+TRACE... Done with CF init, Uwind ",F7.2,2x,F7.2)')minval(WINDXY(:,1)),maxval(WINDXY(:,1))
write(WINDBG%FHNDL,'("+TRACE... Done with CF init, Vwind ",F7.2,2x,F7.2)')minval(WINDXY(:,2)),maxval(WINDXY(:,2))
FLUSH(WINDBG%FHNDL)
END IF
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE UPDATE_WIND(K)
#ifdef NCDF
USE NETCDF
#endif
USE DATAPOOL
IMPLICIT NONE
REAL(rkind) :: TMP(MNP,2)
#if defined NCDF || defined GRIB_API_ECMWF
REAL(rkind) :: cf_w1, cf_w2
INTEGER :: IT, IFILE
#endif
INTEGER, intent(in) :: K
!AR: All crap ... defining K without using means that nobody has ever checked the results or anything else, so why coding at all?
!AR: Mathieu can you please fix this !!!
IF (PrintLOG) THEN
WRITE(WINDBG%FHNDL,*) 'MAIN%TMJD=', MAIN%TMJD
WRITE(WINDBG%FHNDL,*) 'SEWI(TMJD,EMJD)=', SEWI%TMJD, SEWI%EMJD
END IF
IF ( LSEWD .AND. (MAIN%TMJD .ge. SEWI%TMJD-1.E-8) .AND. (MAIN%TMJD .le. SEWI%EMJD+1.e-8) ) THEN
IF (IWINDFORMAT == 1) THEN
!NDM: Need to add the facility for LINTERWD
CALL CSEVAL( WIN%FHNDL, WIN%FNAME, .TRUE., 2, TMP, MULTIPLE_IN_WIND)
DVWIND = (TMP-WINDXY)/SEWI%DELT*MAIN%DELT
#ifdef NCDF
ELSE IF (IWINDFORMAT == 2) THEN ! DWD_NETCDF
CALL MOVE_BY_ONE_INDEX(IFILE, IT)
CALL READ_NETCDF_DWD(IFILE, IT, TMP)
DVWIND = (TMP-WINDXY)/SEWI%DELT*MAIN%DELT
ELSE IF (IWINDFORMAT == 3) THEN ! NOAA CFRS ... the 1st step is analysis and then we have 5 + 1 forecasts, which give one the option to use either only the 6 forecast's after the analysis or use the analysis with 5 forecast's
CALL MOVE_BY_ONE_INDEX(IFILE, IT)
CALL READ_NETCDF_CRFS(IFILE, IT, TMP)
DVWIND = (TMP-WINDXY)/SEWI%DELT*MAIN%DELT
ELSE IF (IWINDFORMAT == 4) THEN ! NOAA NARR ...
CALL MOVE_BY_ONE_INDEX(IFILE, IT)
CALL READ_NETCDF_NARR(IFILE, IT, TMP)
DVWIND = (TMP-WINDXY)/SEWI%DELT*MAIN%DELT
ELSE IF (IWINDFORMAT == 5) THEN
IF (K.EQ.1) THEN
REC1_wind_old = 0
REC2_wind_old = 0
END IF
CALL GET_CF_TIME_INDEX(eVAR_WIND, REC1_wind_new,REC2_wind_new,cf_w1,cf_w2)
IF (REC1_wind_new.NE.REC1_wind_old) THEN
CALL READ_INTERP_NETCDF_CF_WWM_WIND(REC1_wind_new,tmp_wind1)
END IF
IF (REC2_wind_new.NE.REC2_wind_old) THEN
CALL READ_INTERP_NETCDF_CF_WWM_WIND(REC2_wind_new,tmp_wind2)
END IF
IF (cf_w1.NE.1) THEN
WINDXY(:,:) = cf_w1*tmp_wind1(:,:)+cf_w2*tmp_wind2(:,:)
ELSE
WINDXY(:,:) = cf_w1*tmp_wind1(:,:)
END IF
REC1_wind_old = REC1_wind_new
REC2_wind_old = REC2_wind_new
ELSE IF (IWINDFORMAT == 6) THEN
IF (K.EQ.1) THEN
REC1_wind_old = 0
REC2_wind_old = 0
END IF
CALL GET_CF_TIME_INDEX(eVAR_WIND, REC1_wind_new,REC2_wind_new,cf_w1,cf_w2)
IF (REC1_wind_new.NE.REC1_wind_old) THEN
CALL READ_DIRECT_NETCDF_CF(eVAR_wind, REC1_wind_new,tmp_wind1)
END IF
IF (REC2_wind_new.NE.REC2_wind_old) THEN
CALL READ_DIRECT_NETCDF_CF(eVAR_wind, REC2_wind_new,tmp_wind2)
END IF
IF (cf_w1.NE.1) THEN
WINDXY(:,:) = cf_w1*tmp_wind1(:,:)+cf_w2*tmp_wind2(:,:)
ELSE
WINDXY(:,:) = cf_w1*tmp_wind1(:,:)
END IF
REC1_wind_old = REC1_wind_new
REC2_wind_old = REC2_wind_new
#endif
#ifdef GRIB_API_ECMWF
ELSE IF (IWINDFORMAT == 7) THEN
IF (K.EQ.1) THEN
REC1_wind_old = 0
REC2_wind_old = 0
END IF
CALL GET_CF_TIME_INDEX(eVAR_WIND, REC1_wind_new,REC2_wind_new,cf_w1,cf_w2)
IF (REC1_wind_new.NE.REC1_wind_old) THEN
CALL READ_GRIB_WIND(REC1_wind_new,tmp_wind1)
END IF
IF (REC2_wind_new.NE.REC2_wind_old) THEN
CALL READ_GRIB_WIND(REC2_wind_new,tmp_wind2)
END IF
IF (cf_w1.NE.1) THEN
WINDXY(:,:) = cf_w1*tmp_wind1(:,:)+cf_w2*tmp_wind2(:,:)
ELSE
WINDXY(:,:) = cf_w1*tmp_wind1(:,:)
END IF
REC1_wind_old = REC1_wind_new
REC2_wind_old = REC2_wind_new
#endif
END IF
SEWI%TMJD = SEWI%TMJD + SEWI%DELT*SEC2DAY
END IF
IF (PrintLOG) THEN
write(WINDBG%FHNDL,'("max WINDXY:",2F7.2)')maxval(WINDXY(:,1)),maxval(WINDXY(:,2))
write(WINDBG%FHNDL,'("min WINDXY:",2F7.2)')minval(WINDXY(:,1)),minval(WINDXY(:,2))
END IF
IF (LWINDSWAN) THEN
WRITE(3333,*) SEWI%TMJD
WRITE(3333,*) WINDXY(:,1)
WRITE(3333,*) WINDXY(:,2)
END IF
IF (LSEWD.AND.(IWINDFORMAT.NE.5).AND.(IWINDFORMAT.NE.6) ) THEN
WINDXY = WINDXY + DVWIND
END IF
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE MOVE_BY_ONE_INDEX(IFILE, IT)
USE DATAPOOL
implicit none
integer, intent(out) :: IFILE, IT
idxWind =idxWind+1
IF (idxWind .gt. NDT_WIND_ALL_FILES) THEN
CALL WWM_ABORT('Need wind after the time')
END IF
IFILE=WIND_TIME_IFILE(idxWind)
IT=WIND_TIME_IT(idxWind)
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE FIND_WIND_NEAREST_LOWER_IDX(eTime, idx)
USE DATAPOOL
implicit none
real(rkind), intent(in) :: eTime
integer, intent(out) :: idx
CHARACTER(LEN=15) :: eTimeStr
integer eIdxF, eIdx
eIdxF=-1
DO eIdx=1,NDT_WIND_ALL_FILES
IF (WIND_TIME_ALL_FILES(eIdx) .le. eTime + THR8) THEN
eIdxF=eIdx
ENDIF
END DO
IF (eIdxF .eq. -1) THEN
WRITE(WINDBG%FHNDL,*) 'NDT_WIND_ALL_FILES=', NDT_WIND_ALL_FILES
DO eIdx=1,NDT_WIND_ALL_FILES
CALL MJD2CT(WIND_TIME_ALL_FILES(eIdx),eTimeStr)
WRITE(WINDBG%FHNDL,*) ' eIdx=', eIdx
WRITE(WINDBG%FHNDL,*) ' eTime=', WIND_TIME_ALL_FILES(eIdx)
WRITE(WINDBG%FHNDL,*) ' eTimeStr=', eTimeStr
END DO
FLUSH(WINDBG%FHNDL)
CALL WWM_ABORT('We failed to find the wind index')
END IF
idx=eIdxF
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE KERNEL_INTERP_UV_WINDFD(outwind)
USE DATAPOOL
IMPLICIT NONE
INTEGER I, J
REAL(rkind), INTENT(out) :: outwind(MNP_WIND,2)
REAL(rkind) :: Uw, Vw
INTEGER IX, IY
REAL(rkind) :: cf_scale_factor, cf_add_offset
INTEGER SHIFTXY(4,2)
SHIFTXY(1,1)=0
SHIFTXY(1,2)=0
SHIFTXY(2,1)=1
SHIFTXY(2,2)=0
SHIFTXY(3,1)=0
SHIFTXY(3,2)=1
SHIFTXY(4,1)=1
SHIFTXY(4,2)=1
cf_scale_factor = eVAR_WIND % cf_scale_factor
cf_add_offset = eVAR_WIND % cf_add_offset
DO I = 1, MNP_WIND
Uw=ZERO
Vw=ZERO
IX=CF_IX(I)
IY=CF_IY(I)
DO J=1,4
Uw=Uw + CF_COEFF(J,I)*UWIND_FD(IX+SHIFTXY(J,1),IY+SHIFTXY(J,2))
Vw=Vw + CF_COEFF(J,I)*VWIND_FD(IX+SHIFTXY(J,1),IY+SHIFTXY(J,2))
END DO
outwind(I,1)=Uw*cf_scale_factor + cf_add_offset
outwind(I,2)=Vw*cf_scale_factor + cf_add_offset
END DO
IF (PrintLOG) THEN
WRITE(WINDBG%FHNDL,*) 'KERNEL_INTERP_UV_WINDFD'
WRITE(WINDBG%FHNDL,*) 'UWIND_FD, min/max=', minval(UWIND_FD), maxval(UWIND_FD)
WRITE(WINDBG%FHNDL,*) 'VWIND_FD, min/max=', minval(VWIND_FD), maxval(VWIND_FD)
WRITE(WINDBG%FHNDL,*) 'UWIND_FE, min/max=', minval(outwind(:,1)), maxval(outwind(:,1))
WRITE(WINDBG%FHNDL,*) 'VWIND_FE, min/max=', minval(outwind(:,2)), maxval(outwind(:,2))
! WRITE(WINDBG%FHNDL,*) 'max(CF_COEFF)=', maxval(abs(CF_COEFF))
! WRITE(WINDBG%FHNDL,*) 'cf_scale_factor=', cf_scale_factor
! WRITE(WINDBG%FHNDL,*) 'cf_add_offset=', cf_add_offset
FLUSH(WINDBG%FHNDL)
END IF
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
#ifdef NCDF
SUBROUTINE LOAD_INTERP_ARRAY(FileSave, success)
USE DATAPOOL
USE NETCDF
IMPLICIT NONE
logical, intent(out) :: success
character(len=256), intent(in) :: FileSave
character (len = *), parameter :: CallFct = "LOAD_INTERP_ARRAY"
integer, allocatable :: CF_IX_GLOBAL(:), CF_IY_GLOBAL(:)
real(rkind), allocatable :: CF_COEFF_GLOBAL(:,:)
integer, allocatable :: ListFirstMNP(:)
integer, allocatable :: CF_IX_loc(:), CF_IY_loc(:)
real(rkind), allocatable :: CF_COEFF_loc(:,:)
integer iret, ncid, varid
integer IP, IPglob, iPROC, NPloc, IPloc
INQUIRE(FILE=TRIM(FileSave), EXIST=LPRECOMP_EXIST)
IF (LPRECOMP_EXIST .eqv. .FALSE.) THEN
success=.FALSE.
RETURN
END IF
success=.TRUE.
# ifdef MPI_PARALL_GRID
IF (myrank .eq. 0) THEN
# endif
allocate(CF_IX_GLOBAL(np_total), CF_IY_GLOBAL(np_total), CF_COEFF_GLOBAL(4,np_total), stat=istat)
IF (istat/=0) CALL WWM_ABORT('wwm_wind, allocate error 52')
!
iret=nf90_open(TRIM(FileSave), NF90_NOWRITE, ncid)
CALL GENERIC_NETCDF_ERROR_WWM(CallFct, 1, iret)
!
iret=nf90_inq_varid(ncid, "CF_IX", varid)
CALL GENERIC_NETCDF_ERROR_WWM(CallFct, 5, ISTAT)
iret=NF90_GET_VAR(ncid, varid, CF_IX_GLOBAL)
CALL GENERIC_NETCDF_ERROR_WWM(CallFct, 11, ISTAT)
!
iret=nf90_inq_varid(ncid, "CF_IY", varid)
CALL GENERIC_NETCDF_ERROR_WWM(CallFct, 5, ISTAT)
iret=NF90_GET_VAR(ncid, varid, CF_IY_GLOBAL)
CALL GENERIC_NETCDF_ERROR_WWM(CallFct, 11, ISTAT)
!
iret=nf90_inq_varid(ncid, "CF_COEFF", varid)
CALL GENERIC_NETCDF_ERROR_WWM(CallFct, 5, ISTAT)
iret=NF90_GET_VAR(ncid, varid, CF_COEFF_GLOBAL)
CALL GENERIC_NETCDF_ERROR_WWM(CallFct, 11, ISTAT)
!
iret=nf90_close(ncid)
CALL GENERIC_NETCDF_ERROR_WWM(CallFct, 27, iret)
!
# ifdef MPI_PARALL_GRID
END IF
# endif
!
# ifdef MPI_PARALL_GRID
IF (MULTIPLE_IN_WIND .eqv. .FALSE.) THEN
CF_IX=CF_IX_GLOBAL
CF_IY=CF_IY_GLOBAL
CF_COEFF=CF_COEFF_GLOBAL
deallocate(CF_IX_GLOBAL, CF_IY_GLOBAL, CF_COEFF_GLOBAL)
ELSE
IF (myrank .eq. 0) THEN
allocate(ListFirstMNP(nproc), stat=istat)
IF (istat/=0) CALL WWM_ABORT('wwm_wind, allocate error 52')
ListFirstMNP=0
DO iProc=2,nproc
ListFirstMNP(iProc)=ListFirstMNP(iProc-1) + ListMNP(iProc-1)
END DO
DO IP=1,MNP
IPglob=iplg(IP)
CF_IX(IP)=CF_IX_GLOBAL(IPglob)
CF_IY(IP)=CF_IY_GLOBAL(IPglob)
CF_COEFF(:,IP)=CF_COEFF_GLOBAL(:,IPglob)
END DO
!
DO iPROC=2,nproc
NPloc=ListMNP(iPROC)
allocate(CF_IX_loc(NPloc), CF_IY_loc(NPloc), CF_COEFF_loc(4,NPloc), stat=istat)
IF (istat/=0) CALL WWM_ABORT('wwm_wind, allocate error 52')
DO IPloc=1,NPloc
IPglob=ListIPLG(IPloc + ListFirstMNP(iPROC))
CF_IX_loc(IPloc)=CF_IX_GLOBAL(IPglob)
CF_IY_loc(IPloc)=CF_IY_GLOBAL(IPglob)
CF_COEFF_loc(:, IPloc)=CF_COEFF_GLOBAL(:, IPglob)
END DO
CALL MPI_SEND(CF_IX_loc, NPloc, itype, iPROC-1, 711, comm, ierr)
CALL MPI_SEND(CF_IY_loc, NPloc, itype, iPROC-1, 712, comm, ierr)
CALL MPI_SEND(CF_COEFF_loc, 4*NPloc, rtype, iPROC-1, 713, comm, ierr)
deallocate(CF_IX_loc, CF_IY_loc, CF_COEFF_loc)
END DO
deallocate(ListFirstMNP)
deallocate(CF_IX_GLOBAL, CF_IY_GLOBAL, CF_COEFF_GLOBAL)
ELSE
CALL MPI_RECV(CF_IX, MNP, itype, 0, 711, comm, istatus, ierr)
CALL MPI_RECV(CF_IY, MNP, itype, 0, 712, comm, istatus, ierr)
CALL MPI_RECV(CF_COEFF, 4*MNP, rtype, 0, 713, comm, istatus, ierr)
END IF
END IF
# else
CF_IX=CF_IX_GLOBAL
CF_IY=CF_IY_GLOBAL
CF_COEFF=CF_COEFF_GLOBAL
deallocate(CF_IX_GLOBAL, CF_IY_GLOBAL, CF_COEFF_GLOBAL)
# endif
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE SAVE_INTERP_ARRAY(FileSave)
USE DATAPOOL
USE NETCDF
IMPLICIT NONE
character(len=256), intent(in) :: FileSave
character (len = *), parameter :: CallFct = "SAVE_INTERP_ARRAY"
integer, allocatable :: CF_IX_GLOBAL(:), CF_IY_GLOBAL(:)
real(rkind), allocatable :: CF_COEFF_GLOBAL(:,:)
integer, allocatable :: CF_IX_loc(:), CF_IY_loc(:)
real(rkind), allocatable :: CF_COEFF_loc(:,:)
integer, allocatable :: ListFirstMNP(:)
integer ncid, iret, var_id
integer mnp_dims, four_dims
integer IP, IPglob, iPROC, NP_RESloc, IPloc
IF (PrintLOG) THEN
WRITE(STAT%FHNDL,*) 'minval(CF_IX)=', minval(CF_IX)
WRITE(STAT%FHNDL,*) 'minval(CF_IY)=', minval(CF_IY)
WRITE(STAT%FHNDL,*) 'minval(CF_COEFF)=', minval(CF_COEFF)
END IF
# ifdef MPI_PARALL_GRID
IF (MULTIPLE_IN_WIND .eqv. .TRUE.) THEN
IF (myrank .eq. 0) THEN
allocate(ListFirstMNP(nproc), stat=istat)
IF (istat/=0) CALL WWM_ABORT('wwm_wind, allocate error 52')
ListFirstMNP=0
DO iProc=2,nproc
ListFirstMNP(iProc)=ListFirstMNP(iProc-1) + ListMNP(iProc-1)
END DO
!
allocate(CF_IX_GLOBAL(np_total), CF_IY_GLOBAL(np_total), CF_COEFF_GLOBAL(4,np_total), stat=istat)
IF (istat/=0) CALL WWM_ABORT('wwm_wind, allocate error 52')
!
DO IP=1,NP_RES
IPglob=iplg(IP)
CF_IX_GLOBAL(IPglob)=CF_IX(IP)
CF_IY_GLOBAL(IPglob)=CF_IY(IP)
CF_COEFF_GLOBAL(:, IPglob)=CF_COEFF(:,IP)
END DO
DO iPROC=2,nproc
NP_RESloc=ListNP_RES(iPROC)
IF (PrintLOG) THEN
WRITE(STAT%FHNDL,*) ' iPROC=', iPROC, ' NP_RES_loc=', NP_RESloc
END IF
allocate(CF_IX_loc(NP_RESloc), CF_IY_loc(NP_RESloc), CF_COEFF_loc(4,NP_RESloc), stat=istat)
IF (istat/=0) CALL WWM_ABORT('wwm_wind, allocate error 52')
CALL MPI_RECV(CF_IX_loc, NP_RESloc, itype, iProc-1, 611, comm, istatus, ierr)
CALL MPI_RECV(CF_IY_loc, NP_RESloc, itype, iProc-1, 612, comm, istatus, ierr)
CALL MPI_RECV(CF_COEFF_loc, 4*NP_RESloc, rtype, iProc-1, 613, comm, istatus, ierr)
DO IPloc=1,NP_RESloc
IPglob=ListIPLG(IPloc + ListFirstMNP(iProc))
CF_IX_GLOBAL(IPglob)=CF_IX_loc(IPloc)
CF_IY_GLOBAL(IPglob)=CF_IY_loc(IPloc)
CF_COEFF_GLOBAL(:, IPglob)=CF_COEFF_loc(:, IPloc)
END DO
deallocate(CF_IX_loc, CF_IY_loc, CF_COEFF_loc)
END DO
deallocate(ListFirstMNP)
ELSE
CALL MPI_SEND(CF_IX, NP_RES, itype, 0, 611, comm, ierr)
CALL MPI_SEND(CF_IY, NP_RES, itype, 0, 612, comm, ierr)
CALL MPI_SEND(CF_COEFF, 4*NP_RES, rtype, 0, 613, comm, ierr)
END IF
ELSE
allocate(CF_IX_GLOBAL(np_total), CF_IY_GLOBAL(np_total), CF_COEFF_GLOBAL(4,np_total), stat=istat)
IF (istat/=0) CALL WWM_ABORT('wwm_wind, allocate error 52')
CF_IX_GLOBAL=CF_IX
CF_IY_GLOBAL=CF_IY
CF_COEFF_GLOBAL=CF_COEFF
END IF
# else
allocate(CF_IX_GLOBAL(np_total), CF_IY_GLOBAL(np_total), CF_COEFF_GLOBAL(4,np_total), stat=istat)
IF (istat/=0) CALL WWM_ABORT('wwm_wind, allocate error 52')
CF_IX_GLOBAL=CF_IX
CF_IY_GLOBAL=CF_IY
CF_COEFF_GLOBAL=CF_COEFF
# endif
!
! Now writing up
!
# ifdef MPI_PARALL_GRID
IF (myrank .eq. 0) THEN
# endif
iret=nf90_create(TRIM(FileSave), nf90_CLOBBER, ncid)
CALL GENERIC_NETCDF_ERROR_WWM(CallFct, 1, iret)
!
iret = nf90_def_dim(ncid, 'mnp', np_total, mnp_dims)
CALL GENERIC_NETCDF_ERROR_WWM(CallFct, 2, iret)
!
iret = nf90_def_dim(ncid, 'four', 4, four_dims)
CALL GENERIC_NETCDF_ERROR_WWM(CallFct, 3, iret)
!
iret=nf90_def_var(ncid,'CF_IX',NF90_INT,(/ mnp_dims/),var_id)
CALL GENERIC_NETCDF_ERROR_WWM(CallFct, 4, iret)
!
iret=nf90_def_var(ncid,'CF_IY',NF90_INT,(/ mnp_dims/),var_id)
CALL GENERIC_NETCDF_ERROR_WWM(CallFct, 5, iret)
!
iret=nf90_def_var(ncid,'CF_COEFF',NF90_RUNTYPE,(/ four_dims, mnp_dims/),var_id)
CALL GENERIC_NETCDF_ERROR_WWM(CallFct, 6, iret)
!
iret=nf90_close(ncid)
CALL GENERIC_NETCDF_ERROR_WWM(CallFct, 7, iret)
!
! Now writing the data
!
IF (PrintLOG) THEN
WRITE(STAT%FHNDL,*) 'minval(CF_IX_GLOBAL)=', minval(CF_IX_GLOBAL)
WRITE(STAT%FHNDL,*) 'minval(CF_IY_GLOBAL)=', minval(CF_IY_GLOBAL)
WRITE(STAT%FHNDL,*) 'minval(CF_COEFF_GLOBAL)=', minval(CF_COEFF_GLOBAL)
END IF
!
iret=nf90_open(TRIM(FileSave), NF90_WRITE, ncid)
CALL GENERIC_NETCDF_ERROR_WWM(CallFct, 8, iret)
!
iret=nf90_inq_varid(ncid,'CF_IX',var_id)
CALL GENERIC_NETCDF_ERROR_WWM(CallFct, 9, iret)
!
iret=nf90_put_var(ncid,var_id,CF_IX_GLOBAL,start = (/ 1 /), count=(/np_total/))
CALL GENERIC_NETCDF_ERROR_WWM(CallFct, 10, iret)
!
iret=nf90_inq_varid(ncid,'CF_IY',var_id)
CALL GENERIC_NETCDF_ERROR_WWM(CallFct, 11, iret)
!
iret=nf90_put_var(ncid,var_id,CF_IY_GLOBAL,start = (/ 1 /), count=(/np_total/))
CALL GENERIC_NETCDF_ERROR_WWM(CallFct, 12, iret)
!
iret=nf90_inq_varid(ncid,'CF_COEFF',var_id)
CALL GENERIC_NETCDF_ERROR_WWM(CallFct, 13, iret)
!
iret=nf90_put_var(ncid,var_id,CF_COEFF_GLOBAL,start = (/ 1, 1 /), count=(/4, np_total/))
CALL GENERIC_NETCDF_ERROR_WWM(CallFct, 14, iret)
!
iret=nf90_close(ncid)
CALL GENERIC_NETCDF_ERROR_WWM(CallFct, 15, iret)
!
deallocate(CF_IX_GLOBAL, CF_IY_GLOBAL, CF_COEFF_GLOBAL)
# ifdef MPI_PARALL_GRID
END IF
# endif
END SUBROUTINE SAVE_INTERP_ARRAY
#endif
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE COMPUTE_SINGLE_INTERPOLATION_INFO(TheInfo, EXTRAPO_IN, eX, eY, eCF_IX, eCF_IY, eCF_COEFF, EXTRAPO_OUT)
USE DATAPOOL
IMPLICIT NONE
type(FD_FORCING_GRID), intent(in) :: TheInfo
logical, intent(in) :: EXTRAPO_IN
real(rkind), intent(in) :: eX, eY
integer, intent(out) :: eCF_IX, eCF_IY
real(rkind), intent(out) :: eCF_COEFF(4)
logical, intent(out) :: EXTRAPO_OUT
!
integer IX, IY
integer IXs, IYs
integer IXmin, IYmin, IXmax, IYmax
integer nx, ny
integer aShift
REAL(rkind) :: WI(3), X(3), Y(3), a, b
REAL(rkind) :: MinDist, eDist
nx = TheInfo % nx_dim
ny = TheInfo % ny_dim
MinDist=LARGE
EXTRAPO_OUT=.FALSE.
IXs=-1
IYs=-1
DO IX=1,nx-1
DO IY=1,ny-1
eDist=(eX-TheInfo % LON(IX,IY))**2 + (eY-TheInfo % LAT(IX,IY))**2
IF (eDist .lt. MinDist) THEN
MinDist=eDist
IXs=IX
IYs=IY
END IF
END DO
END DO
aShift=1
DO
IXmin=max(1, IXs - aShift)
IYmin=max(1, IYs - aShift)
IXmax=min(nx-1, IXs+aShift)
IYmax=min(ny-1, IYs+aShift)
DO IX=IXmin,IXmax
DO IY=IYmin,IYmax
!
! First triangle
!
X(1)=TheInfo % LON(IX, IY)
X(2)=TheInfo % LON(IX+1, IY)
X(3)=TheInfo % LON(IX, IY+1)
Y(1)=TheInfo % LAT(IX, IY)
Y(2)=TheInfo % LAT(IX+1, IY)
Y(3)=TheInfo % LAT(IX, IY+1)
CALL INTELEMENT_COEF(X,Y,eX,eY,WI)
IF (minval(WI) .ge. -THR) THEN
eCF_IX=IX
eCF_IY=IY
a=WI(2)
b=WI(3)
eCF_COEFF(1)=(1-a)*(1-b)
eCF_COEFF(2)=a*(1-b)
eCF_COEFF(3)=(1-a)*b
eCF_COEFF(4)=a*b
RETURN
END IF
!
! Second triangle
!
X(1)=TheInfo % LON(IX+1, IY+1)
X(2)=TheInfo % LON(IX+1, IY)
X(3)=TheInfo % LON(IX, IY+1)
Y(1)=TheInfo % LAT(IX+1, IY+1)
Y(2)=TheInfo % LAT(IX+1, IY)
Y(3)=TheInfo % LAT(IX, IY+1)
CALL INTELEMENT_COEF(X,Y,eX,eY,WI)
IF (minval(WI) .ge. -THR) THEN
eCF_IX=IX
eCF_IY=IY
a=1 - WI(3)
b=1 - WI(2)
eCF_COEFF(1)=(1-a)*(1-b)
eCF_COEFF(2)=a*(1-b)
eCF_COEFF(3)=(1-a)*b
eCF_COEFF(4)=a*b
RETURN
END IF
END DO
END DO
IF ((IXmin .eq. 1).and.(IYmin .eq. 1).and.(IXmax .eq. nx-1).and.(IYmax .eq. ny-1)) THEN
EXIT
END IF
aShift=aShift + 1
END DO
IF (EXTRAPO_IN) THEN
EXTRAPO_OUT=.TRUE.
eCF_IX = IXs
eCF_IY = IYs
eCF_COEFF(1)=1
eCF_COEFF(2)=0
eCF_COEFF(3)=0
eCF_COEFF(4)=0
WRITE(STAT % FHNDL,*) 'Point ', eX, '/', eY, ' outside grid'
WRITE(STAT % FHNDL,*) 'MinDist=', MinDist
ELSE
WRITE(STAT % FHNDL,*) 'aShift=', aShift
WRITE(STAT % FHNDL,*) 'eX=', eX, 'eY=', eY
FLUSH(STAT % FHNDL)
CALL WWM_ABORT('We find a model point outside of the available forcing grid')
ENDIF
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE COMPUTE_CF_COEFFICIENTS(TheInfo)
USE DATAPOOL
IMPLICIT NONE
type(FD_FORCING_GRID), intent(in) :: TheInfo
integer I
integer eCF_IX, eCF_IY
real(rkind) eCF_COEFF(4)
integer :: nbExtrapolation = 0
character(len=256) :: FileSave = "wwm_filesave_interp_array.nc"
logical success
logical EXTRAPO_OUT
real(rkind) eX, eY
IF (PrintLOG) THEN
WRITE(WINDBG%FHNDL,*) 'Starting node loop for calcs of coefs'
END IF
allocate(CF_IX(MNP_WIND), CF_IY(MNP_WIND), CF_COEFF(4,MNP_WIND), stat=istat)
IF (istat/=0) CALL WWM_ABORT('wwm_wind, allocate error 52')
#ifdef NCDF
IF (PrintLOG) THEN
WRITE(WINDBG%FHNDL,*) 'LSAVE_INTERP_ARRAY=', LSAVE_INTERP_ARRAY
END IF
IF (LSAVE_INTERP_ARRAY) THEN
CALL LOAD_INTERP_ARRAY(FileSave, success)
IF (PrintLOG) THEN
WRITE(WINDBG%FHNDL,*) 'success=', success
END IF
IF (success .eqv. .TRUE.) RETURN
END IF
#endif
CF_IX=0
CF_IY=0
CF_COEFF=0
IF (PrintLOG) THEN
WRITE(WINDBG%FHNDL,*) 'min(lon)=', minval(TheInfo % LON)
WRITE(WINDBG%FHNDL,*) 'max(lon)=', maxval(TheInfo % LON)
WRITE(WINDBG%FHNDL,*) 'min(lat)=', minval(TheInfo % LAT)
WRITE(WINDBG%FHNDL,*) 'max(lat)=', maxval(TheInfo % LAT)
END IF
DO I = 1, MNP_WIND
eX=XP_WIND(I)
eY=YP_WIND(I)
CALL COMPUTE_SINGLE_INTERPOLATION_INFO(TheInfo, EXTRAPOLATION_ALLOWED_WIND, eX, eY, eCF_IX, eCF_IY, eCF_COEFF, EXTRAPO_OUT)
IF (PrintLOG) THEN
WRITE(WINDBG%FHNDL,'(4I10,10F20.10)') I, MNP_WIND, eCF_IX, eCF_IY, eCF_COEFF
END IF
CF_IX(I) = eCF_IX
CF_IY(I) = eCF_IY
CF_COEFF(:,I) = eCF_COEFF
IF (EXTRAPO_OUT) nbExtrapolation = nbExtrapolation + 1
END DO
#ifdef NCDF
IF (LSAVE_INTERP_ARRAY) THEN
CALL SAVE_INTERP_ARRAY(FileSave)
END IF
#endif
IF (PrintLOG) THEN
IF (EXTRAPOLATION_ALLOWED_WIND .eqv. .TRUE.) THEN
WRITE(WINDBG%FHNDL,*) ' nbExtrapolation=', nbExtrapolation
END IF
WRITE(WINDBG%FHNDL,*) ' done interp calcs'
END IF
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE GET_CF_TIME_INDEX(eVAR, REC1, REC2, w1, w2)
! For given wwm_time and wind_time return records to get and weights for time
! interpolation F(wwm_time)=F(rec1)*w1 + F(rec2)*w2
!
USE DATAPOOL
IMPLICIT NONE
TYPE(VAR_NETCDF_CF), intent(in) :: eVAR
REAL(rkind), INTENT(OUT) :: w1, w2
INTEGER, INTENT(OUT) :: REC1, REC2
REAL(rkind) :: eTime1, eTime2
INTEGER :: iTime
DO iTime=2,eVAR % nbTime
eTime1=eVAR % ListTime(iTime-1)
eTime2=eVAR % ListTime(iTime)
IF ((eTime1 .le. MAIN%TMJD).and.(MAIN%TMJD .le. eTime2)) THEN
REC2=iTime
REC1=iTime-1
w2=(MAIN % TMJD - eTime1)/(eTime2-eTime1)
w1=(eTime2 - MAIN % TMJD)/(eTime2-eTime1)
RETURN
END IF
END DO
IF (PrintLOG) THEN
WRITE(WINDBG%FHNDL,*) 'Time error in wind for CF'
WRITE(WINDBG%FHNDL,*) 'MAIN % TMJD=', MAIN%TMJD
WRITE(WINDBG%FHNDL,*) 'min(wind_time_mjd)=', minval(eVAR % ListTime)
WRITE(WINDBG%FHNDL,*) 'max(wind_time_mjd)=', maxval(eVAR % ListTime)
FLUSH(WINDBG%FHNDL)
END IF
CALL WWM_ABORT('Error in CF wind forcing time setup')
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE SYNCHRONIZE_WIND_TIME_IFILE_IT
USE DATAPOOL
IMPLICIT NONE
integer IPROC, eInt(1)
#ifdef MPI_PARALL_GRID
IF (.NOT. MULTIPLE_IN_WIND) THEN
IF (myrank .eq. 0) THEN
eInt(1)=NDT_WIND_ALL_FILES
DO IPROC=2,nproc
CALL MPI_SEND(eInt,1,itype, iProc-1, 811, comm, ierr)
END DO
DO IPROC=2,nproc
CALL MPI_SEND(WIND_TIME_ALL_FILES,NDT_WIND_ALL_FILES,rtype, iProc-1, 812, comm, ierr)
CALL MPI_SEND(WIND_TIME_IFILE,NDT_WIND_ALL_FILES,itype, iProc-1, 813, comm, ierr)
CALL MPI_SEND(WIND_TIME_IT,NDT_WIND_ALL_FILES,itype, iProc-1, 814, comm, ierr)
END DO
ELSE
CALL MPI_RECV(eInt,1,itype, 0, 811, comm, istatus, ierr)
NDT_WIND_ALL_FILES=eInt(1)
ALLOCATE(WIND_TIME_ALL_FILES(NDT_WIND_ALL_FILES), WIND_TIME_IFILE(NDT_WIND_ALL_FILES), WIND_TIME_IT(NDT_WIND_ALL_FILES), stat=istat)
IF (istat/=0) CALL WWM_ABORT('wwm_wind, allocate error 3')
CALL MPI_RECV(WIND_TIME_ALL_FILES,NDT_WIND_ALL_FILES,rtype, 0, 812, comm, istatus, ierr)
CALL MPI_RECV(WIND_TIME_IFILE,NDT_WIND_ALL_FILES,rtype, 0, 813, comm, istatus, ierr)
CALL MPI_RECV(WIND_TIME_IT,NDT_WIND_ALL_FILES,rtype, 0, 814, comm, istatus, ierr)
END IF
END IF
#endif
SEWI%DELT = ( WIND_TIME_ALL_FILES(2) - WIND_TIME_ALL_FILES(1) ) * DAY2SEC
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
#ifdef NCDF
SUBROUTINE INIT_NETCDF_DWD
USE DATAPOOL
USE NETCDF
IMPLICIT NONE
INTEGER :: IT, IFILE
INTEGER :: ILON_ID, ILAT_ID, ITIME_ID, I, J, COUNTER
REAL(rkind) :: DTMP
REAL(rkind), ALLOCATABLE :: WIND_TIME(:)
character ( len = 20 ) chrtmp
character ( len = 15 ) chrdate
integer, dimension(nf90_max_var_dims) :: dimIDs
character (len = *), parameter :: CallFct="INIT_NETCDF_DWD"
# ifdef MPI_PARALL_GRID
IF (MULTIPLE_IN_WIND .or. (myrank .eq. 0)) THEN
# endif
OPEN(WIN%FHNDL,FILE=WIN%FNAME,STATUS='OLD')
!
! count number of netcdf files in list ...
!
NUM_NETCDF_FILES = 0
DO
READ( WIN%FHNDL, *, IOSTAT = ISTAT )
IF ( ISTAT /= 0 ) EXIT
NUM_NETCDF_FILES = NUM_NETCDF_FILES + 1
END DO
REWIND (WIN%FHNDL)
IF (NUM_NETCDF_FILES .eq. 0) THEN
WRITE(WINDBG%FHNDL,*) 'We have NUM_NETCDF_FILES = 0'
WRITE(WINDBG%FHNDL,*) 'In routine INIT_NETCDF_DWD'
WRITE(WINDBG%FHNDL,*) 'Wrong file is file=', TRIM(WIN%FNAME)
FLUSH(WINDBG%FHNDL)
CALL WWM_ABORT('Please correct your setup')
END IF
ALLOCATE(NETCDF_FILE_NAMES(NUM_NETCDF_FILES), stat=istat)
IF (istat/=0) CALL WWM_ABORT('wwm_wind, allocate error 3')
DO IT = 1, NUM_NETCDF_FILES
READ( WIN%FHNDL, *) NETCDF_FILE_NAMES(IT)
! WRITE(WINDBG%FHNDL,*) IT, NETCDF_FILE_NAMES(IT)
END DO
CLOSE (WIN%FHNDL)
!
! check number of time steps in netcdf file ... it is assumed that all files have the same ammount of time steps ...
!
CALL TEST_FILE_EXIST_DIE("Missing wind file : ", TRIM(NETCDF_FILE_NAMES(1)))
ISTAT = NF90_OPEN(NETCDF_FILE_NAMES(1), NF90_NOWRITE, WIND_NCID)
CALL GENERIC_NETCDF_ERROR_WWM(CallFct, 1, ISTAT)
ISTAT = nf90_inq_varid(WIND_NCID, 'initial_time0_encoded', ITIME_ID)
CALL GENERIC_NETCDF_ERROR_WWM(CallFct, 2, ISTAT)
ISTAT = NF90_INQUIRE_VARIABLE(WIND_NCID, ITIME_ID, dimids = dimids)
CALL GENERIC_NETCDF_ERROR_WWM(CallFct, 3, ISTAT)
ISTAT = nf90_inquire_dimension(WIND_NCID, dimIDs(1), len = NDT_WIND_FILE)
CALL GENERIC_NETCDF_ERROR_WWM(CallFct, 4, ISTAT)
!
! check dimensions in the netcdf ... again it is assumed that this is not changing for all files ...
!
ISTAT = nf90_inq_varid(WIND_NCID, 'g0_lon_2', ILON_ID)
CALL GENERIC_NETCDF_ERROR_WWM(CallFct, 5, ISTAT)
ISTAT = NF90_INQUIRE_VARIABLE(WIND_NCID, ILON_ID, dimids = dimIDs)
CALL GENERIC_NETCDF_ERROR_WWM(CallFct, 6, ISTAT)
ISTAT = nf90_inquire_dimension(WIND_NCID, dimIDs(1), len = NDX_WIND)
CALL GENERIC_NETCDF_ERROR_WWM(CallFct, 7, ISTAT)
ISTAT = nf90_inq_varid(WIND_NCID, 'g0_lat_1', ILAT_ID)