-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathwwm_main.F90
914 lines (828 loc) · 32.8 KB
/
wwm_main.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
#include "wwm_functions.h"
!
! __ __ __ __ _____ ____ ____
!/ \ / \/ \ / \/ \ \ \ / /
!\ \/\/ /\ \/\/ / \ / \ ______ \ Y /
! \ / \ / Y \ /_____/ \ /
! \__/\ / \__/\ /\____|__ / \___/
! \/ \/ \/
!
! WWM-V (Wind Wave Model)
!
! The 1st version of the WWM-I code was written by Jian-Ming Liau in his thesis supervised by Tai-Wen Hsu (Liau et al. 2002).
! The source code served as the basis for my thesis that was as well supervised by Tai-Wen Hsu and Ulrich Zanke. In my thesis work
! new numerics and source terms have beend developed (Roland, 2008) and resulted in the WWM-II version of the code. Following this
! the code has served from than as a basis for a 10 year development. In this time the source code was significantly rewritten and
! enhanced with various capabilities. The numerics have been completely revised (Roland, 2008) and parallelized.
! The source term package of Ardhuin et al. 2009, 2010 and from ECMWF (courtesy Jean-Bidlot) was implemented in the WWM-III.
! The code has served as a basis for a 10 year development. In this time the source code was significantly rewritten and
! enhanced with various capabilities. The numerics have been completely revised (Roland, 2008), which lead to the version of WWM-II.
! In WWM-III the model was fully parallelized using Domain Decomposition and coupled to SCHISM. Moreover,
! the source term package of Ardhuin et al. 2009, 2010 and from ECMWF (courtesy Jean-Bidlot) was implemented in the WWM-III.
! The I/O was completely rewritten in NETCDF and various common wind fields can be read such as CFRS, ECMWF, NCEP or others.
! Parallelization is done using the PDLIB decomposition library developed by BGS IT&E GmbH.
!
! Since WWM-II was never released and WWM-III is only given within the SCHISM repository, we introduce finally, the stand alone github
! repository of WWM-V, which is now the official version.
!
!
! Developers
!
! Lead:
!
! Aron Roland (Roland & Partner, Darmstadt)
! Mathieu Dutour Sikiric (IRB, Zagreb)
! Lorenzo Mentaschi (XXXX, Bologna, Italy)
! Yinglong Joseph Zhang (VIMS, USA)
! Tai-Wen Hsu (NTOU, Taiwan)
! Jian-Ming Liau (XXXX, Taiwan)
!
! Contributors:
!
! Christian Ferrarin (ISMAR-CNR, Venice, Italy),
! Fabrice Ardhuin (IFREMER, Brest, France),
! Thomas Huxhorn (BGS IT&E, Darmstadt, Germany),
! Andrea Fortunato (LNEC, Lissabon, Portugal),
! Guillaume Dodet (IFREMER, Brest, France),
! Jean Bidlot (ECMWF, Reading, U.K.)
! Xavier Bertin (La Rochelle, France)
! Alexander Babanin (XXXX, Australia)
! Xavier Bertin (UNR, La Rochelle, France),
! Jean Bidlot (ECMWF, Reading, U.K.)
! Guillaume Dodet (IFREMER, Brest, France),
! Christian Ferrarin (ISMAR-CNR, Venice, Italy),
! Andrea Fortunato (LNEC, Lissabon, Portugal),
! Thomas Huxhorn (BGS IT&E, Darmstadt, Germany),
! Ivica Janekovic, (UWA, Perth, Australia)
! Kai Li (XXXX, PR China),
! Kevin Maarten (UNR, La Rochele, France),
! Peter Janssen (ECWMF, Reading, U.K.)
! Stefan Zieger (XXXX, Australia)
!
! Copyright: 2008 - 2021 (Aron Roland, Jian-Ming Liau, Tai-Wen Hsu)
! All Rights Reserved
!
! http://www.apache.org/licenses/LICENSE-2.0
!
! Unless required by applicable law or agreed to in writing, software
! distributed under the License is distributed on an "AS IS" BASIS,
! WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
! See the License for the specific language governing permissions and
! limitations under the License.
!
!**********************************************************************
!* *
!**********************************************************************
#ifdef SCHISM
!!! SUBROUTINE WWM_II(IT_SCHISM,icou_elfe_wwm,DT_SELFE0,NSTEP_WWM0)
SUBROUTINE WWM_II(IT_SCHISM,icou_elfe_wwm,DT_SCHISM0,NSTEP_WWM0,RADFLAG2)
USE DATAPOOL
use schism_msgp !, only : myrank,parallel_abort,itype,comm,ierr
use schism_glbl, only : iplg,ielg,wwave_force
IMPLICIT NONE
INTEGER, INTENT(IN) :: NSTEP_WWM0, icou_elfe_wwm
REAL(rkind), INTENT(IN) :: DT_SCHISM0
CHARACTER(LEN=3), INTENT(OUT) :: RADFLAG2
! REAL(rkind), INTENT(OUT) :: STOKES_X,STOKES_Y,JPRESS,SBR,SBF
REAL(rkind), SAVE :: SIMUTIME
REAL(rkind) :: T1, T2
REAL(rkind) :: TIME1, TIME2, TIME3, TIME4, TIME5, TIME6, TIME7
INTEGER :: I, IP, IT_SCHISM, K
REAL(rkind) :: DT_PROVIDED
REAL(rkind) :: OUTPAR(OUTVARS), OUTWINDPAR(WINDVARS), WALOC(NUMSIG,NUMDIR)
character(LEN=15) :: CALLFROM
IF (PrintLOG) THEN
WRITE(STAT%FHNDL,'("+TRACE...",A)') 'ENTERING WWM_II'
FLUSH(STAT%FHNDL)
END IF
#ifdef TIMINGS
TIME1 = mpi_wtime()
#endif
IF (LNANINFCHK) THEN
WRITE(DBG%FHNDL,*) ' STARTING WWM FROM SCHISM ', SUM(AC2)
IF (SUM(AC2) .NE. SUM(AC2)) call wwm_abort('NAN IN MAIN 1')
ENDIF
if(WINDVARS/=size(WIND_INTPAR,2)) call wwm_abort('Dimension mismatch: OUTWINDPAR and out_wwm_windpar')
if(OUTVARS/=size(OUTT_INTPAR,2)) call wwm_abort('Dimension mismatch: OUTPAR and out_wwm')
NSTEPWWM = NSTEP_WWM0
DT_SCHISM = DT_SCHISM0
#ifdef TIMINGS
T1 = MyREAL(IT_SCHISM-NSTEPWWM)*DT_SELFE0 ! Beginn time step ...STOKES_VEL
T2 = MyREAL(IT_SCHISM)*DT_SELFE0 ! End of time time step ...
#endif
DT_PROVIDED=NSTEPWWM*DT_SCHISM
IF (abs(MAIN%DELT - DT_PROVIDED).gt.THR) THEN
WRITE(DBG%FHNDL,*) 'MAIN%DELT=', MAIN%DELT, ' in wwminput.nml'
WRITE(DBG%FHNDL,*) 'But nstep_wwm*dt=', DT_PROVIDED
WRITE(DBG%FHNDL,*) 'nstep_wwm=', NSTEPWWM
WRITE(DBG%FHNDL,*) ' dt=', DT_SCHISM
CALL WWM_ABORT('Correct coupled model time-steppings')
ENDIF
IF (LNANINFCHK) THEN
WRITE(DBG%FHNDL,*) ' FIRST SUM IN MAIN ', SUM(AC2)
IF (SUM(AC2) .NE. SUM(AC2)) CALL WWM_ABORT('NAN IN MAIN 2')
ENDIF
IF (PrintLOG) THEN
WRITE(STAT%FHNDL,'("+TRACE...",A)') ' ---- ALL CHECKS DONE'
FLUSH(STAT%FHNDL)
END IF
SIMUTIME = SIMUTIME + MAIN%DELT
IF (.NOT. LWINDFROMWWM) THEN
WINDXY(:,1) = WINDX0
WINDXY(:,2) = WINDY0
END IF
IF (icou_elfe_wwm == 1) THEN ! Full coupling
WLDEP = DEP8
WATLEV = ETA2
WATLEVOLD = ETA1
DEP = MAX(ZERO,WLDEP + WATLEV)
CURTXY(:,1) = UU2(NVRT,:)
CURTXY(:,2) = VV2(NVRT,:)
LSECU = .TRUE.
LSEWL = .TRUE.
ELSE IF (icou_elfe_wwm == 0) THEN ! No interaction at all
WLDEP = DEP8
WATLEV = ZERO
WATLEVOLD = ZERO
DEP = WLDEP
CURTXY(:,1) = ZERO !REAL(rkind)(UU2(NVRT,:))
CURTXY(:,2) = ZERO !REAL(rkind)(VV2(NVRT,:))
LSECU = .FALSE.
LSEWL = .FALSE.
ELSE IF (icou_elfe_wwm == 2) THEN ! Currents and water levels in wwm but no radiation stress in SCHISM
WLDEP = DEP8
WATLEV = ETA2
WATLEVOLD = ETA1
DEP = MAX(ZERO, WLDEP + WATLEV)
CURTXY(:,1) = UU2(NVRT,:)
CURTXY(:,2) = VV2(NVRT,:)
LSECU = .TRUE.
LSEWL = .TRUE.
ELSE IF (icou_elfe_wwm == 3) THEN ! No current and no water levels in wwm but radiation stress in SCHISM
WLDEP = DEP8
WATLEV = ZERO
WATLEVOLD = ZERO
DEP = WLDEP
CURTXY(:,1) = ZERO !REAL(rkind)(UU2(NVRT,:))
CURTXY(:,2) = ZERO !REAL(rkind)(VV2(NVRT,:))
LSECU = .FALSE.
LSEWL = .FALSE.
ELSE IF (icou_elfe_wwm == 4) THEN ! No current but water levels in wwm and radiation stresss in SCHISM
WLDEP = DEP8
WATLEV = ETA2
WATLEVOLD = ETA1
DEP = WLDEP
CURTXY(:,1) = 0.!UU2(NVRT,:)
CURTXY(:,2) = 0.!UU2(NVRT,:)
LSECU = .FALSE.
LSEWL = .TRUE.
ELSE IF (icou_elfe_wwm == 5) THEN ! No current but water levels in wwm and no radiation stress in SCHISM
WLDEP = DEP
WATLEV = ETA2
WATLEVOLD = ETA1
DEP = WLDEP
CURTXY(:,1) = 0.!UU2(NVRT,:)
CURTXY(:,2) = 0.!UU2(NVRT,:)
LSECU = .FALSE.
LSEWL = .TRUE.
ELSE IF (icou_elfe_wwm == 6) THEN ! Currents but no water levels in wwm and radiation stress in SCHISM
WLDEP = DEP
WATLEV = ZERO
WATLEVOLD = ZERO
DEP = WLDEP
CURTXY(:,1) = UU2(NVRT,:)
CURTXY(:,2) = VV2(NVRT,:)
LSECU = .TRUE.
LSEWL = .FALSE.
ELSE IF (icou_elfe_wwm == 7) THEN ! Currents but no water levels in wwm and no radiation stress in SCHISM
WLDEP = DEP
WATLEV = ZERO
WATLEVOLD = ZERO
DEP = WLDEP
CURTXY(:,1) = UU2(NVRT,:)
CURTXY(:,2) = UU2(NVRT,:)
LSECU = .TRUE.
LSEWL = .FALSE.
END IF
LCALC = .TRUE.
DEPDT = (WATLEV - WATLEVOLD) / DT_SCHISM0
IF (LNANINFCHK) THEN
CALL SCHISM_NANCHECK_INPUT_A
END IF
IF (LBCSE) THEN
CALL SET_WAVE_BOUNDARY_CONDITION
END IF
IF (LNANINFCHK) THEN
WRITE(DBG%FHNDL,*) ' AFTER SETTING BOUNDARY CONDITION IN MAIN ', SUM(AC2)
IF (SUM(AC2) .NE. SUM(AC2)) CALL WWM_ABORT('NAN IN MAIN 3')
ENDIF
IF (LFIRSTSTEP) THEN
IF (INITSTYLE == 1) CALL INITIAL_CONDITION!We need to call for the case of wind dependent intiial guess this call since before we have no wind from SCHISM
LFIRSTSTEP = .FALSE.
LCALC = .TRUE.
END IF
#ifdef TIMINGS
TIME2 = mpi_wtime()
#endif
IF (LNANINFCHK) THEN
WRITE(DBG%FHNDL,*) ' BEFORE COMPUTE ', SUM(AC2)
IF (SUM(AC2) .NE. SUM(AC2)) CALL WWM_ABORT('NAN IN MAIN 4')
ENDIF
IF (PrintLOG) THEN
WRITE(STAT%FHNDL,'("+TRACE...",A)') 'ENTERING COMPUTE'
FLUSH(STAT%FHNDL)
END IF
CALLFROM='SCHISM'
IF (LQSTEA) THEN
CALL QUASI_STEADY(KKK)
ELSE
CALL UN_STEADY(KKK)
END IF
IF (LNANINFCHK) THEN
WRITE(DBG%FHNDL,*) ' AFTER COMPUTE ', SUM(AC2)
IF (SUM(AC2) .NE. SUM(AC2)) CALL WWM_ABORT('NAN IN MAIN 5')
ENDIF
#ifdef TIMINGS
TIME3 = mpi_wtime()
#endif
IF (PrintLOG) THEN
WRITE(STAT%FHNDL,'("+TRACE...",A,F15.4)') 'FINISHED COMPUTE nth call to WWM', SIMUTIME
FLUSH(STAT%FHNDL)
END IF
DO IP = 1, MNP
WALOC = AC2(:,:,IP)
IF (DEP(IP) .GT. DMIN) THEN
CALL INTPAR(IP, NUMSIG, WALOC, OUTPAR)
OUTT_INTPAR(IP,:) = OUTPAR
CALL WINDPAR(IP,OUTWINDPAR)
WIND_INTPAR(IP,:) = OUTWINDPAR
IF (LMONO_OUT) THEN
OUTT_INTPAR(IP,1) = OUTT_INTPAR(IP,1) / SQRT(TWO)
END IF
ELSE
OUTT_INTPAR(IP,:) = ZERO
WIND_INTPAR(IP,:) = ZERO
END IF
END DO
! menta
!removing the line below, otherwise with icou_elfe_wwm the the output will be always 0
!IF (icou_elfe_wwm .eq. 0) OUTT_INTPAR = ZERO
IF (PrintLOG) THEN
WRITE(STAT%FHNDL,'("+TRACE...",A,F15.4)') 'FINISHED FILLING RESULTS', SIMUTIME
FLUSH(STAT%FHNDL)
END IF
#ifdef TIMINGS
TIME4 = mpi_wtime()
#endif
!
! Compute radiation stress ...
!
! RADFLAG=VOR , then coupling with SCHISM will gives stokes_velocity (Eq. 17 from Bennis 2011), Wave-induced pressure (Eq. 20) and source momentums (Eq.21)
RADFLAG2=RADFLAG !for output into SCHISM
IF (icou_elfe_wwm == 0 .OR. icou_elfe_wwm == 2 .OR. icou_elfe_wwm == 5 .OR. icou_elfe_wwm == 7) THEN
WWAVE_FORCE = ZERO
!STOKES_X=ZERO
!STOKES_Y=ZERO
!STOKES_VEL=0
!JPRESS=ZERO
!SBR=ZERO
!YJZ: I changed the following to SBF
!SBF=ZERO
ELSE
!MENTA: updated to the calls implemented btw schism and WWMIII
IF (RADFLAG == 'VOR') THEN ! Vortex force formalism as described in Bennis (2011)
CALL STOKES_STRESS_INTEGRAL_SCHISM ! Compute Stokes drift velocities and pressure terms
CALL COMPUTE_CONSERVATIVE_VF_TERMS_SCHISM ! Conservative terms (relative to Stokes drift advection, Coriolis and pressure head: Eq. 17, 19 and 20 from Bennis 2011)
IF (fwvor_breaking == 1) THEN ! BM
CALL COMPUTE_BREAKING_VF_TERMS_SCHISM ! Sink of momentum due to wave breaking and update wwave_force
END IF
ELSE ! Radiation stress formalism (Longuet-Higgins and Stewart, 1962 and 1964) as described in Battjes (1974)
CALL RADIATION_STRESS_SCHISM
ENDIF
IF(SHOREWAFO == 1) CALL SHORELINE_WAVE_FORCES
! Apply ramp on wave forces if wafo_obcramp == 1 (in
! param.nml)
IF (wafo_obcramp == 1) CALL APPLY_WAFO_OPBND_RAMP
IF (.false.) THEN
WRITE(STAT%FHNDL,'("+TRACE...",A,F15.4)') 'FINISHED FILLING VORTEX', SIMUTIME
CALL FLUSH(STAT%FHNDL)
END IF
END IF
! end modif AD
IF (PrintLOG) THEN
WRITE(STAT%FHNDL,'("+TRACE...",A,F15.4)') 'FINISHED FILLING VORTEX', SIMUTIME
FLUSH(STAT%FHNDL)
END IF
#ifdef TIMINGS
TIME5 = mpi_wtime()
#endif
IF (LNANINFCHK) THEN
CALL SCHISM_NANCHECK_INPUT_B
END IF
KKK = KKK + 1
#ifdef TIMINGS
TIME6 = mpi_wtime()
#endif
IF (LNANINFCHK) THEN
WRITE(DBG%FHNDL,*) ' END OF MAIN ', SUM(AC2)
IF (SUM(AC2) .NE. SUM(AC2)) CALL WWM_ABORT ('NAN IN MAIN 5')
ENDIF
IF (PrintLOG) THEN
WRITE(STAT%FHNDL,'("+TRACE...",A,F15.4)') 'END OF COMPUTATIONS NOW RETURN TO SCHISM', SIMUTIME
FLUSH(STAT%FHNDL)
END IF
#ifdef TIMINGS
IF (PrintLOG) THEN
WRITE(STAT%FHNDL,'("+TRACE...",A,F15.6)') '-----TOTAL TIMINGS-----'
WRITE(STAT%FHNDL,'("+TRACE...",A,F15.6)') 'PREPARATION ', TIME2-TIME1
WRITE(STAT%FHNDL,'("+TRACE...",A,F15.6)') 'INTEGRATION ', TIME3-TIME2
WRITE(STAT%FHNDL,'("+TRACE...",A,F15.6)') 'OUTPUT TO SCHISM ', TIME4-TIME3
WRITE(STAT%FHNDL,'("+TRACE...",A,F15.6)') 'RADIATION STRESSES ', TIME5-TIME4
WRITE(STAT%FHNDL,'("+TRACE...",A,F15.6)') 'NAN CHECK ', TIME6-TIME5
WRITE(STAT%FHNDL,'("+TRACE...",A,F15.6)') 'TOTAL TIME ', TIME6-TIME1
WRITE(STAT%FHNDL,'("+TRACE...",A,F15.6)') '------END-TIMINGS- ---'
WRITE(STAT%FHNDL,'("+TRACE...",A,F15.4)') 'FINISHED WITH WWM', SIMUTIME
FLUSH(STAT%FHNDL)
END IF
#endif
! menta: to avoid 'nan from wwm' errors
WHERE(WWAVE_FORCE /= WWAVE_FORCE)
WWAVE_FORCE = 0
END WHERE
END SUBROUTINE WWM_II
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE SCHISM_NANCHECK_INPUT_A
USE DATAPOOL
implicit none
integer IP
DO IP = 1, MNP
IF (WINDXY(IP,1) .NE. WINDXY(IP,1)) THEN
WRITE(DBG%FHNDL,*) 'NaN in WINDX', IP, WINDXY(IP,1)
FLUSH(DBG%FHNDL)
END IF
IF (WINDXY(IP,2) .NE. WINDXY(IP,2)) THEN
WRITE(DBG%FHNDL,*) 'NaN in WINDY', IP, WINDXY(IP,2)
FLUSH(DBG%FHNDL)
END IF
IF (WATLEV(IP) .NE. WATLEV(IP)) THEN
WRITE(DBG%FHNDL,*) 'NaN in WATLEV', IP, WATLEV(IP)
FLUSH(DBG%FHNDL)
END IF
IF (WATLEVOLD(IP) .NE. WATLEVOLD(IP)) THEN
WRITE(DBG%FHNDL,*) 'NaN in WATLEV', IP, WATLEV(IP)
FLUSH(DBG%FHNDL)
END IF
IF (CURTXY(IP,1) .NE. CURTXY(IP,1)) THEN
WRITE(DBG%FHNDL,*) 'NaN in CURTX', IP, CURTXY(IP,1)
FLUSH(DBG%FHNDL)
END IF
IF (CURTXY(IP,2) .NE. CURTXY(IP,2)) THEN
WRITE(DBG%FHNDL,*) 'NaN in CURTY', IP, CURTXY(IP,2)
FLUSH(DBG%FHNDL)
END IF
END DO
END SUBROUTINE SCHISM_NANCHECK_INPUT_A
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE SCHISM_NANCHECK_INPUT_B
USE DATAPOOL
implicit none
integer IP, I
DO IP = 1, MNP
IF (SUM(OUTT_INTPAR(IP,:)) .NE. SUM(OUTT_INTPAR(IP,:))) THEN
DO I = 1, SIZE(OUTT_INTPAR(IP,:))
WRITE(DBG%FHNDL,*) 'NaN in OUTT_INTPAR', IP, I, OUTT_INTPAR(IP,I)
FLUSH(DBG%FHNDL)
END DO
END IF
IF (SUM(WIND_INTPAR(IP,:)) .NE. SUM(WIND_INTPAR(IP,:))) THEN
DO I = 1, SIZE(WIND_INTPAR(IP,:))
WRITE(DBG%FHNDL,*) 'NaN in WIND_INTPAR', IP, I, WIND_INTPAR(IP,I)
FLUSH(DBG%FHNDL)
END DO
END IF
!Error: force defined at side centers
IF (SUM(WWAVE_FORCE(:,IP,:)) .NE. SUM(WWAVE_FORCE(:,IP,:))) THEN
DO I = 1, SIZE(WWAVE_FORCE(1,:,IP)) ! loop over layers ...
WRITE(DBG%FHNDL,*) 'NaN in WWAVE_FORCE', IP, I, WWAVE_FORCE(1,I,IP), WWAVE_FORCE(2,I,IP)
FLUSH(DBG%FHNDL)
END DO
END IF
END DO
END SUBROUTINE SCHISM_NANCHECK_INPUT_B
#endif
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE UN_STEADY(K)
USE DATAPOOL
IMPLICIT NONE
INTEGER, INTENT(IN) :: K
REAL(rkind) :: CONV1, CONV2, CONV3, CONV4, CONV5
REAL(rkind) :: TIME1, TIME2, TIME3, TIME4, TIME5, TIME6
CHARACTER(LEN=15) :: CTIME
#ifdef TIMINGS
CALL WAV_MY_WTIME(TIME1)
#endif
CALL IO_1(K)
#ifdef TIMINGS
CALL WAV_MY_WTIME(TIME2)
#endif
IF (LCFL_CASD) THEN
CALL CFLSPEC()
ENDIF
! CALL Print_SumAC2("Before the advection")
IF (ICOMP .EQ. 0) THEN
CALL COMPUTE_EXPLICIT
ELSE IF (ICOMP .EQ. 1) THEN
CALL COMPUTE_SEMI_IMPLICIT
ELSE IF (ICOMP .EQ. 2) THEN
CALL COMPUTE_SEMI_IMPLICIT
ELSE IF (ICOMP .EQ. 3) THEN
CALL COMPUTE_IMPLICIT
END IF
! CALL Print_SumAC2("After the advection")
#ifdef TIMINGS
CALL WAV_MY_WTIME(TIME3)
#endif
#ifdef WWM_SETUP
IF (LZETA_SETUP) THEN
CALL WAVE_SETUP_COMPUTATION
END IF
#endif
#ifdef TIMINGS
CALL WAV_MY_WTIME(TIME4)
#endif
MAIN%TMJD = MAIN%BMJD + MyREAL(K)*MAIN%DELT*SEC2DAY
RTIME = MAIN%TMJD - MAIN%BMJD
#ifndef SCHISM
# if defined WWM_MPI
IF (myrank.eq.0) WRITE(*,101) K, MAIN%ISTP, RTIME
# else
WRITE(*,101) K, MAIN%ISTP, RTIME
# endif
#endif
CALL IO_2(K)
! CALL Print_SumAC2("After IO_2")
#ifdef TIMINGS
CALL WAV_MY_WTIME(TIME5)
#endif
IF (LCONV .AND. LQSTEA) THEN
CALL CHECK_STEADY(RTIME,CONV1,CONV2,CONV3,CONV4,CONV5)
END IF
IF (ABORT_BLOWUP) THEN
CALL CHECK_FOR_BLOW
END IF
#ifdef TIMINGS
CALL WAV_MY_WTIME(TIME6)
#endif
CALL MJD2CT(MAIN%TMJD, CTIME)
#ifdef TIMINGS
IF (PrintLOG) THEN
WRITE(STAT%FHNDL,'("+TRACE...",A,F15.6,A20)') '-----SIMULATION TIME----- ', MAIN%TMJD, CTIME
WRITE(STAT%FHNDL,'("+TRACE...",A,F15.6)') '-----TOTAL RUN TIMES----- '
WRITE(STAT%FHNDL,'("+TRACE...",A,F15.6)') 'PREPROCESSING ', TIME2-TIME1
WRITE(STAT%FHNDL,'("+TRACE...",A,F15.6)') 'INTEGRATION ', TIME3-TIME2
WRITE(STAT%FHNDL,'("+TRACE...",A,F15.6)') 'WAVE SETUP ', TIME4-TIME3
WRITE(STAT%FHNDL,'("+TRACE...",A,F15.6)') 'POSTPROCESSING ', TIME5-TIME4
WRITE(STAT%FHNDL,'("+TRACE...",A,F15.6)') 'CHECK STEADY ', TIME6-TIME5
WRITE(STAT%FHNDL,'("+TRACE...",A,F15.6)') 'TOTAL TIME ', TIME6-TIME1
WRITE(STAT%FHNDL,'("+TRACE...",A,F15.6)') '-------------TIMINGS-------------'
FLUSH(STAT%FHNDL)
ENDIF
#endif
IF (PrintLOG) THEN
WRITE(STAT%FHNDL,'("+TRACE......",A)') 'LEAVING UN_STEADY'
END IF
101 FORMAT ('+STEP = ',I10,'/',I10,' ( TIME = ',F15.4,' DAYS)')
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE QUASI_STEADY(K)
USE DATAPOOL
IMPLICIT NONE
INTEGER, INTENT(IN) :: K
INTEGER :: IT
REAL(rkind) :: ITERTIME
REAL(rkind) :: CONV1, CONV2, CONV3, CONV4, CONV5
CALL IO_1(K)
IF (LCFL_CASD) CALL CFLSPEC()
IF (LCHKCONV) IP_IS_STEADY = 0 ! Reset local convergence indicators ...
IF (LCHKCONV) IE_IS_STEADY = 0
#ifdef MPI_PARALL_GRID
! NQSITER = NSTEPWWM ! this is not very flexible!
#endif
DO IT = 1, NQSITER
DT_ITER = MAIN%DELT/MyREAL(NQSITER)
IF (ICOMP .EQ. 0) THEN
CALL COMPUTE_EXPLICIT
ELSE IF (ICOMP .EQ. 1) THEN
CALL COMPUTE_SEMI_IMPLICIT
ELSE IF (ICOMP .EQ. 2) THEN
CALL COMPUTE_SEMI_IMPLICIT
ELSE IF (ICOMP .EQ. 3) THEN
CALL COMPUTE_IMPLICIT
END IF
ITERTIME = RTIME*DAY2SEC + IT*DT_ITER
IF (LCHKCONV) THEN
CALL CHECK_STEADY(ITERTIME,CONV1,CONV2,CONV3,CONV4,CONV5)
! DO IP = 1, MNP
! IF (IP_IS_STEADY(IP) .GE. 1) AC2(:,:,IP) = AC1(IP,:,:)
! ENDDO
IF ( (CONV1 .GT. 100._rkind*QSCONV1 .AND. &
& CONV2 .GT. 100._rkind*QSCONV2 .AND. &
& CONV3 .GT. 100._rkind*QSCONV3 .AND. &
& CONV4 .GT. 100._rkind*QSCONV4 .AND. &
& CONV5 .GT. 100._rkind*QSCONV5 .AND. &
& K .NE. 1) .OR. &
& IT .EQ. NQSITER ) THEN
#ifndef SCHISM
WRITE(QSTEA%FHNDL,'(3I10,5F15.8)') K, IT, NQSITER, CONV1, CONV2, CONV3, CONV4, CONV5
#else
if (myrank == 0) WRITE(QSTEA%FHNDL,'(3I10,5F15.8)') K, IT, NQSITER, CONV1, CONV2, CONV3, CONV4, CONV5
#endif
FLUSH(QSTEA%FHNDL)
EXIT
END IF
END IF
IF (LOUTITER) CALL GENERAL_OUTPUT(ITERTIME)
END DO
#ifdef MPI_PARALL_GRID
MAIN%TMJD = MAIN%BMJD + MyREAL(K)*MAIN%DELT*SEC2DAY
RTIME = MAIN%TMJD - MAIN%BMJD
IF (myrank == 0) WRITE(STAT%FHNDL,101) K, MAIN%ISTP, RTIME*DAY2SEC
#else
MAIN%TMJD = MAIN%BMJD + MyREAL(K)*MAIN%DELT*SEC2DAY
RTIME = MAIN%TMJD - MAIN%BMJD
WRITE(STAT%FHNDL,101) K, MAIN%ISTP, RTIME*DAY2SEC
#endif
FLUSH(STAT%FHNDL)
CALL IO_2(K)
101 FORMAT ('+STEP = ',I5,'/',I5,' ( TIME = ',F15.4,'HR )')
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE IO_1(K)
#if defined ROMS_WWM_PGMCL_COUPLING || defined MODEL_COUPLING_ATM_WAV || defined MODEL_COUPLING_OCN_WAV
USE WWMaOCN_PGMCL
#endif
#if defined MODEL_COUPLING_ATM_WAV || defined MODEL_COUPLING_OCN_WAV
USE pgmcl_lib_WWM, only : WAV_all_import_export
#endif
USE DATAPOOL
IMPLICIT NONE
INTEGER, INTENT(IN) :: K
IF (LWINDFROMWWM) THEN
CALL UPDATE_WIND(K)
END IF
#ifndef SCHISM
IF (.NOT. LCPL) THEN
IF (LSECU) THEN
CALL UPDATE_CURRENT(K)
END IF
IF (LSEWL) THEN
CALL UPDATE_WATLEV(K)
END IF
END IF
#endif
IF (MESIC.EQ.1 .AND. LICEFROMWWM) THEN
CALL UPDATE_ICE(K)
END IF
#ifndef SCHISM
IF (LBCSE) THEN
CALL SET_WAVE_BOUNDARY_CONDITION
END IF
#endif
!
! *** coupling via pipe *** read pipe
!
#if !defined SCHISM && !defined ROMS_WWM_PGMCL_COUPLING && !defined MODEL_COUPLING_ATM_WAV && !defined MODEL_COUPLING_OCN_WAV
IF (LCPL .AND. LTIMOR) THEN
CALL PIPE_TIMOR_IN(K)
# ifdef SHYFEM_COUPLING
ELSE IF (LCPL .AND. LSHYFEM) THEN
CALL PIPE_SHYFEM_IN(K)
# endif
ELSE IF (LCPL .AND. LROMS) THEN
CALL PIPE_ROMS_IN(K)
END IF
#endif
#ifdef ROMS_WWM_PGMCL_COUPLING
IF ( K-INT(K/MAIN%ICPLT)*MAIN%ICPLT .EQ. 0 ) THEN
CALL WAV_ocnAwav_import(K)
END IF
IF (K == 1) CALL INITIAL_CONDITION
#endif
#if defined MODEL_COUPLING_ATM_WAV || defined MODEL_COUPLING_OCN_WAV
CALL WAV_all_import_export(K, IFILE, IT)
#endif
!
! *** recalculate water level and current related values
!
IF (LSEWL .OR. LSECU .OR. LCPL) THEN ! LCPL makes sure that when the model is coupled it gets into this part for 100%
DEP = MAX(ZERO,WLDEP + WATLEV) ! d = a + h if -h .gt. a set d to zero
CALL SETSHALLOW
CALL GRADDEP
IF (MESTR == 6) CALL GRAD_CG_K
CALL WAVE_K_C_CG
CALL GRADCURT
CALL SET_IOBPD
CALL SET_IOBPD_BY_DEP
IF (LCFL_CASD) THEN
CALL CFLSPEC
ENDIF
IF (LMAXETOT .AND. MESBR == 0) THEN ! Needed for certain tests without sources terms, where only the depth limiter is acting
CALL SET_HMAX
ENDIF
END IF
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE IO_2(K)
#if defined ROMS_WWM_PGMCL_COUPLING || defined MODEL_COUPLING_ATM_WAV || defined MODEL_COUPLING_OCN_WAV
USE WWMaOCN_PGMCL
#endif
USE DATAPOOL
IMPLICIT NONE
INTEGER, INTENT(IN) :: K
CALL GENERAL_OUTPUT(RTIME*DAY2SEC)
#ifndef SCHISM
# if !defined ROMS_WWM_PGMCL_COUPLING && !defined MODEL_COUPLING_ATM_WAV && !defined MODEL_COUPLING_OCN_WAV
IF (LCPL .AND. LTIMOR) THEN
CALL PIPE_TIMOR_OUT(K)
# ifdef SHYFEM_COUPLING
ELSE IF (LCPL .AND. LSHYFEM) THEN
CALL PIPE_SHYFEM_OUT(K)
# endif
ELSE IF (LCPL .AND. LROMS) THEN
CALL PIPE_ROMS_OUT(K)
END IF
# endif
# ifdef ROMS_WWM_PGMCL_COUPLING
IF ( K-INT(K/MAIN%ICPLT)*MAIN%ICPLT .EQ. 0 ) THEN
CALL WAV_ocnAwav_export(K)
END IF
# endif
#endif
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
#if !defined SCHISM && !defined PDLIB && defined MPI_PARALL_GRID
SUBROUTINE SIMPLE_PRE_READ
USE DATAPOOL
!MENTA
!USE schism_glbl, only : NUMSIG2, NUMDIR2, ics
USE schism_glbl, only : MSC2, MDC2, ics
IMPLICIT NONE
CHARACTER(LEN=20) :: BEGTC, UNITC, ENDTC
REAL(rkind) DELTC
NAMELIST /PROC/ PROCNAME, DIMMODE, LSTEA, LQSTEA, LSPHE, &
& LNAUTIN, LNAUTOUT, LMONO_OUT, LMONO_IN, &
& BEGTC, DELTC, UNITC, ENDTC, DMIN
NAMELIST /GRID/ LCIRD, LSTAG, MINDIR, MAXDIR, NUMDIR, FRLOW, &
& FRHIGH, NUMSIG, FILEGRID, IGRIDTYPE, LSLOP, SLMAX, LVAR1D, &
& LOPTSIG, CART2LATLON, LATLON2CART
INTEGER FHNDL
!
FHNDL=12
CALL TEST_FILE_EXIST_DIE("Missing input file : ", TRIM(INP%FNAME))
OPEN(FHNDL, FILE = TRIM(INP%FNAME))
READ(FHNDL, NML = PROC)
READ(FHNDL, NML = GRID)
CLOSE(FHNDL)
IF (LSPHE) THEN
ics=2
ELSE
ics=1
ENDIF
IF (CART2LATLON) THEN
ics=1
END IF
IF (LATLON2CART) THEN
ics=2
END IF
IF (CART2LATLON .and. LATLON2CART) THEN
CALL WWM_ABORT('You cannot have both CART2LATLON and CART2LONLAT')
END IF
!MENTA
!NUMSIG2=NUMSIG
!NUMDIR2=NUMDIR
MSC2=NUMSIG
MDC2=NUMDIR
END SUBROUTINE
#endif
!**********************************************************************
!* *
!**********************************************************************
#if !defined SCHISM
# if defined ROMS_WWM_PGMCL_COUPLING || defined MODEL_COUPLING_ATM_WAV || defined MODEL_COUPLING_OCN_WAV
SUBROUTINE WWMIII_MPI(MyCOMM)
# else
PROGRAM WWMIII_MPI
# endif
# ifdef ROMS_WWM_PGMCL_COUPLING
USE mod_coupler, only : WAV_COMM_WORLD, MyRankGlobal
# endif
# if defined MODEL_COUPLING_ATM_WAV || defined MODEL_COUPLING_OCN_WAV
USE coupling_var, only : WAV_COMM_WORLD, MyRankGlobal
# endif
USE DATAPOOL, only: MAIN, STAT, LQSTEA, RKIND, LCALC, PrintLOG
# ifdef MPI_PARALL_GRID
USE datapool, only: comm, myrank, ierr, nproc, &
& parallel_finalize
# endif
implicit none
# if defined MPI_PARALL_GRID || defined ROMS_WWM_PGMCL_COUPLING || defined MODEL_COUPLING_ATM_WAV || defined MODEL_COUPLING_OCN_WAV
include 'mpif.h'
# endif
# if defined ROMS_WWM_PGMCL_COUPLING || defined MODEL_COUPLING_ATM_WAV || defined MODEL_COUPLING_OCN_WAV
integer, intent(in) :: MyCOMM
# endif
# ifdef TIMINGS
REAL(rkind) :: TIME1, TIME2
# endif
integer :: k
# if defined DEBUG && (defined MODEL_COUPLING_ATM_WAV || defined MODEL_COUPLING_OCN_WAV)
write(740+MyRankGlobal,*) 'WWMIII_MPI, before mpi_init'
FLUSH(740+MyRankGlobal)
# endif
# if defined WWM_MPI && !defined ROMS_WWM_PGMCL_COUPLING && !defined MODEL_COUPLING_ATM_WAV && !defined MODEL_COUPLING_OCN_WAV
call mpi_init(ierr)
if(ierr/=MPI_SUCCESS) call wwm_abort('Error at mpi_init')
# endif
# if defined DEBUG && (defined MODEL_COUPLING_ATM_WAV || defined MODEL_COUPLING_OCN_WAV)
write(740+MyRankGlobal,*) 'WWMIII_MPI, after mpi_init'
FLUSH(740+MyRankGlobal)
# endif
# ifdef TIMINGS
CALL WAV_MY_WTIME(TIME1)
# endif
# if defined DEBUG && (defined MODEL_COUPLING_ATM_WAV || defined MODEL_COUPLING_OCN_WAV)
write(740+MyRankGlobal,*) 'WWMIII_MPI, after WAV_MY_WTIME'
FLUSH(740+MyRankGlobal)
# endif
# if defined ROMS_WWM_PGMCL_COUPLING || defined MODEL_COUPLING_ATM_WAV || defined MODEL_COUPLING_OCN_WAV
comm=MyCOMM
WAV_COMM_WORLD=MyCOMM
# else
# ifdef WWM_MPI
call mpi_comm_dup(MPI_COMM_WORLD,comm,ierr)
if(ierr/=MPI_SUCCESS) call wwm_abort('Error at mpi_comm_dup')
# endif
# endif
# if defined DEBUG && (defined MODEL_COUPLING_ATM_WAV || defined MODEL_COUPLING_OCN_WAV)
write(740+MyRankGlobal,*) 'WWMIII_MPI, after mpi_comm_dup and WAV_COMM_WORLD'
FLUSH(740+MyRankGlobal)
# endif
# ifdef MPI_PARALL_GRID
call mpi_comm_size(comm,nproc,ierr)
if(ierr/=MPI_SUCCESS) call wwm_abort('Error at mpi_comm_size')
call mpi_comm_rank(comm,myrank,ierr)
if(ierr/=MPI_SUCCESS) call wwm_abort('Error at mpi_comm_rank')
# ifndef PDLIB
CALL SIMPLE_PRE_READ
# endif
# endif
# if defined DEBUG && (defined MODEL_COUPLING_ATM_WAV || defined MODEL_COUPLING_OCN_WAV)
write(740+MyRankGlobal,*) 'WWMIII_MPI, after mpi_comm_size/rank'
FLUSH(740+MyRankGlobal)
# endif
CALL INITIALIZE_WWM
# if defined DEBUG && (defined MODEL_COUPLING_ATM_WAV || defined MODEL_COUPLING_OCN_WAV)
write(740+MyRankGlobal,*) 'WWMIII_MPI, after INITIALIZE_WWM'
FLUSH(740+MyRankGlobal)
# endif
DO K = 1, MAIN%ISTP
CALL Print_SumAC2("In the time loop")
IF (LQSTEA) THEN
CALL QUASI_STEADY(K)
ELSE
CALL UN_STEADY(K)
END IF
LCALC=.FALSE.
END DO
# ifdef TIMINGS
CALL WAV_MY_WTIME(TIME2)
IF (PrintLOG) THEN
WRITE(STAT%FHNDL,'("+TRACE...",A,F15.6)') '-----TOTAL TIME IN PROG-----', TIME2-TIME1
END IF
# endif
# if defined MPI_PARALL_GRID && !defined ROMS_WWM_PGMCL_COUPLING && !defined MODEL_COUPLING_ATM_WAV && !defined MODEL_COUPLING_OCN_WAV
call parallel_finalize
# endif
# if defined ROMS_WWM_PGMCL_COUPLING || defined MODEL_COUPLING_ATM_WAV || defined MODEL_COUPLING_OCN_WAV
END SUBROUTINE
# else
END PROGRAM
# endif
#endif
!**********************************************************************
!* *
!**********************************************************************