-
Notifications
You must be signed in to change notification settings - Fork 0
/
w3sic5md.ftn
executable file
·2049 lines (2041 loc) · 72.5 KB
/
w3sic5md.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 "w3macros.h"
!/ ------------------------------------------------------------------- /
MODULE W3SIC5MD
!/
!/ +-----------------------------------+
!/ | WAVEWATCH III NOAA/NCEP |
!/ | Q. Liu |
!/ | E. Rogers |
!/ | FORTRAN 90 |
!/ | Last update : 02-Jun-2017 |
!/ +-----------------------------------+
!/
!/ 15-Mar-2016 : Origination. ( version 5.10 )
!/ ( Q. Liu )
!/ 15-Mar-2016 : Started from w3sic1/2/3/4 module ( Q. Liu )
!/
!/ 24-Apr-2017 : Adding more filters ( Q. Liu )
!/
!/ 29-Apr-2017 : Introducing CMPLX_TANH2 ( Q. Liu )
!/
!/ 02-Jun-2017 : Update to version 5.16 ( Q. Liu )
!/
!/ 17-Jun-2017 : Remove some unnecessary lines ( Q. Liu )
!/ (cg_ice, detla function, complx_tanh etc.)
!/
!/ 20-Aug-2018 : Ready to be merged to master (v6.06)( Q. Liu)
!/
!/ 1. Purpose :
! Calculate ice source term S_{ice} according to a viscoelastic sea
! ice model (Mosig et al. 2015)
!
! Reference:
! Mosig, J. E. M., F. Montiel, and V. A. Squire (2015),
! Comparison of viscoelastic-type models for ocean wave attenuation
! in ice-covered seas, J. Geophys. Res. Oceans, 120, 6072–6090,
! doi:10.1002/2015JC010881.
!
! 2. Variables and types :
!
! Name Type Scope Description
! ----------------------------------------------------------------
! KSP Int. Private the kind parameter for single precision
! real variables
! KDP Int. Private Same as KSP but for double precision
! KSPC Int. Private the kind parameter for single precision
! complex variables
! KDPC Int. Private Same as KSPC but for double precision
! ERRTOL Real Private A real parameter used for "==" test
! ----------------------------------------------------------------
!
! 3. Subroutines and functions :
!
! Name Type Scope Description
! ----------------------------------------------------------------
! W3SIC5 Subr. Public Ice source term
! W3IC5WNCG Subr. Public Wavenumber and group velocity of ice-
! coupled waves
! FSDISP Subr. Public Solving the ice-coupled wave dispersion
! BALANCING_MATRIX
! Subr. Private Balancing the matrix before we try to
! find its eigenvalues
! EIG_HQR Subr. Private QR algorithm for real Hessenberg matrix
! (eigenvalues-finding algorithm)
! POLYROOTS Subr. Private Finding roots of a general polynomial
! NR_CORR Func. Private Get the Newton-Raphson correction term
! for iteration
! NR_ROOT Func. Private Newton-Raphson algorithm for solving
! the ice-coupled wave dispersion
! CMPLX_SINH, CMPLX_COSH, CMPLX_TANH2
! Func. Private sinh, cosh, tanh for complex inputs
! INIT_RANDOM_SEED
! Subr. Private Initialize the random seed based on
! the system's time
! ----------------------------------------------------------------
!
! 4. Subroutines and functions used :
! See subroutine documentation
!
! 5. Remarks :
!
! 6. Switches :
! See subroutine documentation
!
! 7. Source code:
!/
!/ ------------------------------------------------------------------- /
IMPLICIT NONE
!/
PUBLIC :: W3SIC5, W3IC5WNCG, FSDISP
PRIVATE :: BALANCING_MATRIX, EIG_HQR, POLYROOTS
PRIVATE :: NR_CORR, NR_ROOT
PRIVATE :: CMPLX_SINH, CMPLX_COSH, CMPLX_TANH2
PRIVATE :: INIT_RANDOM_SEED
!/
PRIVATE :: KSP, KDP, KSPC, KDPC, ERRTOL
!/ ------------------------------------------------------------------- /
!/ Parameter list
! Kind for single- and double-precision real type
INTEGER, PARAMETER :: KSP = KIND(1.0)
INTEGER, PARAMETER :: KDP = KIND(1.0D0)
!
! Kind for single- and double-precision complex type
INTEGER, PARAMETER :: KSPC = KIND((1.0, 1.0))
INTEGER, PARAMETER :: KDPC = KIND((1.0D0, 1.0D0))
REAL, PARAMETER :: ERRTOL = 1.E-12
!/
CONTAINS
!/ ------------------------------------------------------------------- /
!/
SUBROUTINE W3SIC5 (A, DEPTH, CG, WN, IX, IY, S, D)
!/
!/ +-----------------------------------+
!/ | WAVEWATCH III NOAA/NCEP |
!/ | Q. Liu |
!/ | E. Rogers |
!/ | FORTRAN 90 |
!/ | Last update : 25-Apr-2017 |
!/ +-----------------------------------+
!/
!/ 23-Mar-2016 : Origination ( version 5.10 )
!/ ( Q. Liu)
!/ 23-Mar-2016 : Started from w3sic1/2/3/4 subr. ( Q. Liu)
!/ 05-Apr-2016 : Options for Cg_{ice} or Cg ( Q. Liu)
!/ 25-Apr-2017 : Add more filters ( Q. Liu)
!/ 20-Aug-2018 : Ready to be merged to master (v6.06)( Q. Liu)
!/
!/ Copyright 2009 National Weather Service (NWS),
!/ National Oceanic and Atmospheric Administration. All rights
!/ reserved. WAVEWATCH III is a trademark of the NWS.
!/ No unauthorized use without permission.
!/
!/ 1. Purpose :
! Calculate ice source term S_{ice} according to a viscoelastic sea
! ice model (Mosig et al. 2015)
!
! Reference:
! Mosig, J. E. M., F. Montiel, and V. A. Squire (2015),
! Comparison of viscoelastic-type models for ocean wave attenuation
! in ice-covered seas, J. Geophys. Res. Oceans, 120, 6072–6090,
! doi:10.1002/2015JC010881.
!
! 2. Method :
! Regarding i/o (general to all Sice modules): S_{ice} source term
! is calculated using up to 5 parameters read from input files.
! These parameters are allowed to vary in space and time.
! The parameters control the exponential decay rate k_i.
! Since there are 5 parameters, this permits description of
! dependence of k_i on frequency or wavenumber.
!
! Sea ice affects the wavenumber k of wind-generated ocean waves.
! The ice-modified wavenumber can be expressed as a complex number
! k = k_r + i * k_i, with the real part k_r representing impact of
! the sea ice on the physical wavelength and propagation speeds,
! producing something analogous to shoaling and refraction by
! bathymetry, whereas the imaginary part of the complex
! wavenumber, k_i, is an exponential decay coefficient
! k_i(x,y,t,sigma) (depending on location, time and frequency,
! respectively), representing wave attenuation, and can be
! introduced in a wave model such as WW3 as S_ice/E=-2*Cg*k_i,
! where S_ice is one of several dissipation mechanisms, along
! with whitecapping, for example, S_ds=S_wc+S_ice+⋯. The k_r -
! modified by ice would enter the model via the C calculations
! on the left-hand side of the governing equation.The fundamentals
! are straightforward, e.g. Rogers and Holland (2009 and
! subsequent unpublished work) modified a similar model, SWAN
! (Booij et al. 1999) to include the effects of a viscous mud
! layer using the same approach (k = k_r + i*k_i) previously.
!
! General approach is analogous to Rogers and Holland (2009)
! approach for mud.
! See text near their eq. 1 :
! k = k_r + i * k_i
! eta(x,t) = Real( a * exp( i * ( k * x - sigma * t ) ) )
! a = a0 * exp( -k_i * x )
! S / E = -2 * Cg * k_i (see also Komen et al. (1994, pg. 170)
!
! Following W3SBT1 as a guide, equation 1 of W3SBT1 says:
! S = D * E
! However, the code of W3SBT1 has
! S = D * A
! This leads me to believe that the calling routine is
! expecting "S/sigma" not "S"
! Thus we will use D = S/E = -2 * Cg * k_i
!
! The calling routine is expecting "S/sigma" not "S"
! Thus we will use D = S/E = -2 * Cg * k_i
! (see also documentation of W3SIC1)
!
! Notes regarding numerics:
! -------------------------
! Experiments with constant k_i values suggest that results may be
! dependent on resolution if insufficient resolution is used.
! For detailed information, see documentation of W3SIC1.
!
! Note regarding applicability/validity:
! --------------------------------------
! Similar to the Wang and Shen model used in w3sic3md, the FS model
! used here is also a viscoelastic-type model which treats the sea
! ice cover as a continuum and uses two empirical rheological para-
! meters, i.e., the effective shear modulus of ice G and the effec-
! tive viscosity η to characterize sea ices of various type. Please
! see the documentation of w3sic3md for a detailed discussion of
! this kind of model.
!
! 3. Parameters :
!
! Parameter list
! ----------------------------------------------------------------
! A R.A. I Action density spectrum (1-D).
! DEPTH Real I Local water depth.
! CG R.A. I Group velocities.
! WN R.A. I Wavenumbers
! IX,IY I.S. I Grid indices.
! S R.A. O Source term (1-D version).
! D R.A. O Diagonal term of derivative (1-D version).
! ----------------------------------------------------------------
!
! 4. Subroutines used :
!
! Name Type Module Description
! ----------------------------------------------------------------
! STRACE Subr. W3SERVMD Subroutine tracing (!/S switch).
! PRT2DS Subr. W3ARRYMD Print plot output (!/T0 switch).
! OUTMAT Subr. W3ARRYMD Matrix output (!/T1 switch).
! W3IC5WNCG Subr. / Wavenumber and group velocity of ice-
! coupled waves
! ----------------------------------------------------------------
! * / means this module
!
! 5. Called by :
!
! Name Type Module Description
! ----------------------------------------------------------------
! W3SRCE Subr. W3SRCEMD Source term integration.
! W3EXPO Subr. N/A ASCII Point output post-processor.
! W3EXNC Subr. N/A NetCDF Point output post-processor.
! GXEXPO Subr. N/A GrADS point output post-processor.
! ----------------------------------------------------------------
!
! 6. Error messages :
!
! None.
!
! 7. Remarks :
!
! If ice parameter 1 is zero, no calculations are made.
!
! 8. Structure :
!
! See source code.
!
! 9. Switches :
!
! !/S Enable subroutine tracing.
! !/T Enable general test output.
! !/T0 2-D print plot of source term.
! !/T1 Print arrays.
!
! 10. Source code :
!/ ------------------------------------------------------------------- /
!/
!/T USE W3ODATMD, ONLY: NDST
!/S USE W3SERVMD, ONLY: STRACE
!/T0 USE W3ARRYMD, ONLY: PRT2DS
!/T1 USE W3ARRYMD, ONLY: OUTMAT
!/
USE CONSTANTS, ONLY: TPI
USE W3SERVMD, ONLY: EXTCDE
USE W3ODATMD, ONLY: NDSE, IAPROC, NAPROC, NAPERR
USE W3GDATMD, ONLY: NK, NTH, NSPEC, SIG, MAPWN, IC5PARS
USE W3IDATMD, ONLY: INFLAGS2, ICEP1, ICEP2, ICEP3, ICEP4, ICEI
!
IMPLICIT NONE
!
!/
!/ ------------------------------------------------------------------- /
!/ Parameter list
!/
REAL, INTENT(IN) :: CG(NK), WN(NK), A(NSPEC), DEPTH
REAL, INTENT(OUT) :: S(NSPEC), D(NSPEC)
INTEGER, INTENT(IN) :: IX, IY
!/ ------------------------------------------------------------------- /
!/ Local parameters
!/
!/S INTEGER, SAVE :: IENT = 0
!/T0 INTEGER :: ITH
!/T0 REAL :: DOUT(NK,NTH)
!/
REAL :: ICECOEF1, ICECOEF2, ICECOEF3, &
ICECOEF4, ICECONC
REAL, DIMENSION(NK) :: D1D, WN_R, WN_I
REAL :: TWN_R, TWN_I
INTEGER :: IK, IKTH
LOGICAL :: NOICE
!/ ------------------------------------------------------------------- /
!/
!/S CALL STRACE (IENT, 'W3SIC5')
!
! 0. Initializations ------------------------------------------------ /
D = 0.
D1D = 0.
WN_R = 0.
WN_I = 0.
!
ICECOEF1 = 0.
ICECOEF2 = 0.
ICECOEF3 = 0.
ICECOEF4 = 0.
ICECONC = 0.
!
! Set the ice parameters from input
IF (INFLAGS2(-7)) THEN
ICECOEF1 = ICEP1(IX, IY) ! ice thickness h_i
ELSE
IF ( IAPROC .EQ. NAPERR ) &
WRITE (NDSE,1001) 'ICE PARAMETER 1 (HICE)'
CALL EXTCDE(2)
ENDIF
!
IF (INFLAGS2(-6)) THEN
ICECOEF2 = ICEP2(IX, IY) ! effective viscosity of ice η
ELSE
IF ( IAPROC .EQ. NAPERR ) &
WRITE (NDSE,1001) 'ICE PARAMETER 2 (VISC)'
CALL EXTCDE(2)
ENDIF
!
IF (INFLAGS2(-5)) THEN
ICECOEF3 = ICEP3(IX, IY) ! density of ice ρ_i
ELSE
IF ( IAPROC .EQ. NAPERR ) &
WRITE (NDSE,1001) 'ICE PARAMETER 3 (DENS)'
CALL EXTCDE(2)
ENDIF
!
IF (INFLAGS2(-4)) THEN
ICECOEF4 = ICEP4(IX, IY) ! effective shear modulus of ice G
ELSE
IF ( IAPROC .EQ. NAPERR ) &
WRITE (NDSE,1001) 'ICE PARAMETER 4 (ELAS)'
CALL EXTCDE(2)
ENDIF
!
IF (INFLAGS2(4)) ICECONC = ICEI(IX, IY) ! ice concentration
!
! 1. No ice --------------------------------------------------------- /
NOICE = .FALSE.
! Zero ice thickness
! Very small ice thickness may cause problems in POLYROOTS because
! the first coefficient C1 may be very close to zero. So we regard
! cases where hice is less than 0.0001 as no ice.
! IF (ICECOEF1 < ERRTOL) NOICE = .TRUE.
IF (ICECOEF1 < 0.0001) NOICE = .TRUE.
! zero ice concentration
IF (INFLAGS2(4) .AND. ICECONC < ERRTOL) NOICE = .TRUE.
!
! Calculate the decay rate k_i
IF ( NOICE ) THEN
D1D = 0.
!
! 2. Ice ------------------------------------------------------------- /
ELSE
! W3IC5WNCG(WN_R, WN_I, CG, HICE, IVISC, RHOI, ISMODG, HWAT)
CALL W3IC5WNCG(WN_R, WN_I, CG, ICECOEF1, ICECOEF2, &
ICECOEF3, ICECOEF4, DEPTH)
! recall that D=S/E=-2*Cg_{ice}*k_i
! In some cases, the FS model yields very larg Cg_{ice}, which
! subquently may result in numerical failure due to the violation of CFL
! conditions, therefore we still use ice-free group velocity to advect
! wave packets.
!
DO IK = 1, NK
D1D(IK) = -2.0 * CG(IK) * WN_I(IK)
END DO
END IF
!
! 2.1 Fill diagonal matrix
DO IKTH = 1, NSPEC
D(IKTH) = D1D(MAPWN(IKTH))
END DO
S = D * A
!
! ... Test output of arrays
!
!/T0 DO IK=1, NK
!/T0 DO ITH=1, NTH
!/T0 DOUT(IK,ITH) = D(ITH+(IK-1)*NTH)
!/T0 END DO
!/T0 END DO
!
!/T0 CALL PRT2DS (NDST, NK, NK, NTH, DOUT, SIG(1:), ' ', 1., &
!/T0 0.0, 0.001, 'Diag Sice', ' ', 'NONAME')
!
!/T1 CALL OUTMAT (NDST, D, NTH, NTH, NK, 'diag Sice')
!
! Formats
!
1001 FORMAT(/' *** WAVEWATCH III ERROR IN W3SIC5MD : '/ &
' ',A,' IS NOT DEFINED IN ww3_shel.inp.')
!/
!/ End of W3SIC5------------------------------------------------------ /
!/
END SUBROUTINE W3SIC5
!/ ------------------------------------------------------------------- /
!/
SUBROUTINE W3IC5WNCG(WN_R, WN_I, CG, HICE, IVISC, RHOI, ISMODG, &
HWAT)
!/
!/ +-----------------------------------+
!/ | WAVEWATCH III NOAA/NCEP |
!/ | Q. Liu |
!/ | E. Rogers |
!/ | FORTRAN 90 |
!/ | Last update : 25-Apr-2017 |
!/ +-----------------------------------+
!/
!/ 17-Apr-2016 : Origination ( version 5.10)
!/ ( Q. Liu )
!/ 17-Apr-2016 : Start from W3IC3WNCG_CHENG ( Q. Liu )
!/
!/ 1. Purpose:
! Calculation of complex wavenumber arrays for ice-coupled waves.
!
! This also allows us to use Cg_ice in the advection part of the
! radiative transfer energy equation (RTE). --- abandoned in the end
!
! 2. Method:
! Using the Fox-Squire dispersion relations to get (kr, ki) and
! then get cg by cg = dσ / dk (here dk uses kr)
!
! 3. Parameters:
!
! Parameter list:
! ----------------------------------------------------------------
! Name Type Intent Description
! ----------------------------------------------------------------
! WN_R R.A. I/O the real. part of the wave number
! WN_I R.A. I/O the imag. part of the wave number
! CG R.A. I group velocity (m s^{-1})
! HICE Real. I thickness of ice (m)
! IVISC Real. I viscosity parameter of ice (m^2 s^{-1})
! RHOI Real. I the density of ice (kg m^{-3})
! ISMODG Real. I effecitive shear modulus G of ice (Pa)
! HWAT Real. I water depth
! ----------------------------------------------------------------
! * the intent of WN_R/I must be inout
! * CG is unchanged but still kept here because some legacy reasons.
!
! 4. Subroutines used:
!
! Name Type Module Description
! ----------------------------------------------------------------
! FSDISP Subr. / dispersion relations for ice-coupled waves
! CGINICE5 Subr. / group velocity for given (σ, kr) array
! ----------------------------------------------------------------
!
! 5. Called by :
!
! Name Type Module Description
! ----------------------------------------------------------------
! W3SIC5 Subr. Public Ice source term
! W3WAVE Subr. W3WAVEMD WW3 integration
! ----------------------------------------------------------------
!
! 6. Error messages :
!
! 7. Remarks :
!
! 8. Structure :
!
! See source code.
!
! 9. Switches :
!
! !/S Enable subroutine tracing.
!
! 10. Source code :
!
!/ ------------------------------------------------------------------- /
!/S USE W3SERVMD, ONLY: STRACE
USE CONSTANTS, ONLY: TPI
USE W3GDATMD, ONLY: NK, SIG
USE W3ODATMD, ONLY: NDSE, IAPROC, NAPROC, NAPERR
USE W3SERVMD, ONLY: EXTCDE
!/
IMPLICIT NONE
!/
!/ ------------------------------------------------------------------- /
!/ Parameter list
REAL, INTENT(INOUT) :: WN_R(:), WN_I(:)
REAL, INTENT(IN) :: CG(:)
REAL, INTENT(IN) :: HICE, IVISC, RHOI, ISMODG, HWAT
!/
!/ ------------------------------------------------------------------- /
!/ Local parameters
REAL, ALLOCATABLE :: SIGMA(:)
INTEGER :: KL, KU, IK
REAL :: TWN_R, TWN_I
!/
!/S INTEGER, SAVE :: IENT = 0
!/
!/ ------------------------------------------------------------------- /
!/
!/S CALL STRACE (IENT, 'W3IC5WNCG')
!/
! Initialize SIGMA {in w3gdatmd: SIG (0: NK+1)}
IF (ALLOCATED(SIGMA)) DEALLOCATE(SIGMA); ALLOCATE(SIGMA(SIZE(CG)))
SIGMA = 0.
IF (SIZE(WN_R, 1) .EQ. NK) THEN
KL = 1
KU = NK
SIGMA = SIG(1:NK)
ELSE IF (SIZE(WN_R,1) .EQ. NK+2) THEN
KL = 1
KU = NK+2
SIGMA = SIG(0:NK+1)
ELSE
IF ( IAPROC .EQ. NAPERR ) WRITE(NDSE,900) 'W3IC5WNCG'
CALL EXTCDE(3)
END IF
!
! Fox-Squire dispersion
DO IK = KL, KU
! FSDISP(HICE, IVISC, RHOI, ISMODG, HWAT, WT, WNR, WNI)
CALL FSDISP(HICE, IVISC, RHOI, ISMODG, HWAT, TPI/SIGMA(IK), &
TWN_R, TWN_I)
WN_R(IK) = TWN_R
WN_I(IK) = TWN_I
END DO
!
DEALLOCATE(SIGMA)
!
900 FORMAT(/' *** WAVEWATCH III ERROR IN W3SIC5MD : '/ &
' Subr. ', A, ': Cannot determine bounds of&
& wavenumber array.'/)
!/
!/ End of W3IC5WNCG -------------------------------------------------- /
!/
END SUBROUTINE W3IC5WNCG
!/ ------------------------------------------------------------------- /
!/
SUBROUTINE FSDISP(HICE, IVISC, RHOI, ISMODG, HWAT, WT, WNR, WNI)
!/
!/ +-----------------------------------+
!/ | WAVEWATCH III NOAA/NCEP |
!/ | Q. Liu |
!/ | FORTRAN 90 |
!/ | Last update : 25-Apr-2017 |
!/ +-----------------------------------+
!/
!/ 17-Mar-2016 : Origination ( version 5.10)
!/ ( Q. Liu )
!/ 17-Mar-2016 : Start from the Matlab code `FoxSquire.m` (provided
!/ by Prof. Vernon Squire from University of Otago)
!/ ( Q. Liu )
!/ 25-Apr-2017 : Add more filters ( Q. Liu)
!/
! 1. Purpose :
!
! Calculate the complex wavenumber for waves in ice according to
! the viscoelastic model, i.e., FS model in Mosig et al. (2015)
!
! 2. Method :
! Solving the dispersion relations of FS model (Eq. (20) and (24)
! in Mosig et al. (2015))
!
! Reference:
! Mosig, J. E. M., F. Montiel, and V. A. Squire (2015),
! Comparison of viscoelastic-type models for ocean wave attenuation
! in ice-covered seas, J. Geophys. Res. Oceans, 120, 6072–6090,
! doi:10.1002/2015JC010881.
!
! 3. Parameters :
!
! Parameter list
! ----------------------------------------------------------------
! Name Type Intent Description
! ----------------------------------------------------------------
! HICE Real. IN thickness of ice (m)
! IVISC Real. IN viscosity parameter of ice (m^2 s^{-1})
! RHOI Real. IN the density of ice (kg m^{-3})
! ISMODG Real. IN effecitive shear modulus G of ice (Pa)
! HWAT Real. IN water depth
! WT Real. IN wave period (s; 1/freq)
! WNR Real. Out the real. part of the wave number
! WNI Real. Out the imag. part of the wave number
! ----------------------------------------------------------------
!
! 4. Subroutines used :
!
! Name Type Module Description
! ----------------------------------------------------------------
! STRACE Subr. W3SERVMD Subroutine tracing.
! POLYROOTS Subr. / Find the roots of a general polynomial
! NR_ROOT Func. / Newton-Raphson root finding
! ----------------------------------------------------------------
!
! 5. Called by :
!
! Name Type Module Description
! ----------------------------------------------------------------
! W3IC5WNCG Subr. / Wavenumber and group velocity of ice-
! coupled waves
! ----------------------------------------------------------------
!
! 6. Error messages :
!
! See Format 1000, 1001, 1002
!
! 7. Remarks :
!
! 8. Structure :
!
! See source code.
!
! 9. Switches :
!
! !/S Enable subroutine tracing.
!
! 10. Source code :
!
!/ ------------------------------------------------------------------- /
!/S USE W3SERVMD, ONLY: STRACE
!/
USE CONSTANTS, ONLY: GRAV, TPI
USE W3DISPMD, ONLY: WAVNU1
USE W3SERVMD, ONLY: EXTCDE
USE W3ODATMD, ONLY: NDSE, IAPROC, NAPROC, NAPERR
USE W3GDATMD, ONLY: IC5PARS
USE W3GSRUMD, ONLY: W3INAN
!/
IMPLICIT NONE
!/
!/ ------------------------------------------------------------------- /
!/ Parameter list
REAL, INTENT(IN) :: HICE, IVISC, RHOI, ISMODG, HWAT, WT
REAL, INTENT(OUT) :: WNR, WNI
!/
!/ ------------------------------------------------------------------- /
!/ Local parameters
!
REAL :: IC5MINIG, IC5MINWT, IC5MAXKRATIO, &
IC5MAXKI, IC5MINHW
REAL :: TISMODG, TWT, TRATIO, THW
REAL, PARAMETER :: NU = 0.3, RHOW = 1025.
COMPLEX :: GV, C1
REAL :: SIGMA, C2, WNO, CGO, THKH, &
RTRL(5), RTIM(5), RTANG(5)
INTEGER :: IREAL
COMPLEX(KDPC) :: GUESS, CROOT, C1D
REAL(KDP) :: C2D, HWATD
!/
!/S INTEGER, SAVE :: IENT = 0
!/
!/ ------------------------------------------------------------------- /
!/
!/S CALL STRACE (IENT, 'FSDISP')
! Note, same as W3IC3WNCG_xx in w3sic3md :
! HICE → ICE1
! IVISC → ICE2
! RHOI → ICE3
! ISMODG → ICE4
! 0. Initializations ------------------------------------------------ *
! Set limiters
!
! When G = 0, the FS method does not provide a solution. It is not
! unexpected because the FS model is originally devised as a
! thin elastic plate model in which elasticity is necessary.
!
! The FS algorithm may also have issues for very short wave periods,
! shallow waters and low G (e.g., T~3 s, d~10 m, hi~0.5 m, G<10^6 Pa)
!
IC5MINIG = IC5PARS(1) ! Minimum G
IC5MINWT = IC5PARS(2) ! Minimum T
IC5MAXKRATIO = IC5PARS(3) ! Maximum k_{ow}/k_r
IC5MAXKI = IC5PARS(4) ! Maximum k_i
IC5MINHW = IC5PARS(5) ! Minimum d
!
TISMODG = MAX(IC5MINIG, ISMODG)
TWT = MAX(IC5MINWT, WT)
THW = MAX(IC5MINHW, HWAT)
!
! G <= 0. is not allowed
IF (ABS(TISMODG) < ERRTOL) THEN
IF ( IAPROC .EQ. NAPERR ) WRITE(NDSE, 1000) 'FSDISP'
CALL EXTCDE(1)
END IF
!
! σ = 2π / T
SIGMA = TPI / TWT
!
! Complex shear modulus Gv = G - i σ ρ η
GV = CMPLX(TISMODG, -1. * SIGMA * RHOI * IVISC)
!
! -------------------------------------------------------------------- *
! Note that Eq. (24) in Mosig et al. (2015) can be written like below:
! (c1 * k^5 + c2 * k) * tanh(HWAT*k) - 1 = 0
! Most Important part of this module --------------------------------- *
C1 = GV * HICE**3. / (6. * RHOW * SIGMA**2.)
!
! To be divided by (1-NU) or multiplied by (1+NU) ??
! Beam model: then multiplied by (1+ν)
! Plate model: then divided by (1-ν)
! The beam version is more theoretically (J.E.M. Mosig, personal
! communication, 2016), although there is only very marginal difference
! between this two version as (1+NU = 1.3 and 1/(1-NU) ~ 1.4)
C1 = C1 * (1+NU)
! C1 = C1 / (1-NU)
!
! C2
C2 = GRAV / SIGMA**2. - RHOI * HICE / RHOW
!
! Use the dispersion in open water to get an approximation of
! tanh(HWAT * k). We can also roughly use the dispersion in deep
! water case, that is tanh(HWAT*k) ~ 1.
! Wavenumber in the open water
! WAVNU1(SI, H, K, CG)
CALL WAVNU1(SIGMA, THW, WNO, CGO)
THKH = TANH(WNO * THW)
!
! Get the first guess of the complex wavenumber
CALL POLYROOTS(6, &
(/REAL(REAL(C1))*THKH, 0., 0., 0., C2*THKH, -1./),&
RTRL, RTIM)
RTANG = ATAN2(RTIM, RTRL)
!
! There should only be one real root in RTRL + i * RTIM because in
! this case (ivisc=0) the original viscoelastic-type model reduced to
! the thin elastic plate model which has only one real solution.
! Find its index ...
!
IREAL = MINLOC(ABS(RTANG), DIM=1)
IF (RTRL(IREAL) <= 0. .OR. ABS(RTIM(IREAL)) > ERRTOL) THEN
IF ( IAPROC .EQ. NAPERR ) WRITE(NDSE, 1001) 'FSDISP'
CALL EXTCDE(2)
END IF
!
! Get the first guess for iteration
GUESS = RTRL(IREAL) * EXP(CMPLX(0., 1E-6))
!
! Newton-Raphson method
! Turn c1, c2, hwat to be double
C1D = C1; C2D = C2; HWATD = THW
CROOT = NR_ROOT(C1D, C2D, HWATD, GUESS)
WNR = REAL(REAL(CROOT))
WNI = REAL(AIMAG(CROOT))
!
! RATIO Check
! Using the ratio k0 / kr as a basic check for the reliability of
! FSDISP. The FS dispersion relation can give a very different kr from
! k0, especially for small wave periods (k0/kr is as high as 100).
! From my tests, using IC5MAXKRATIO = 1000. can basically detect most
! spurious solutions (although not all of them)
!
! ISNAN Check
! Common ways used are:
! NAN = SQRT(-1.) or
! a /= a then a is NaN or
! ISNAN func (supported by gfortran & ifort)
! --- ISNAN -> W3INAN because ISNAN is not supported by pgi
! For very few cases, we can get nan | negative ki | kr
!
! (N.B.) NaN problem solved by using CMPLX_TANH2
!
TRATIO = WNO / WNR
IF (W3INAN(WNR) .OR. W3INAN(WNI) .OR. WNR <= 0 .OR. WNI <= 0. &
.OR. TRATIO >= IC5MAXKRATIO) THEN
IF ( IAPROC .EQ. NAPERR ) &
WRITE(NDSE, 1002) 'FSDISP', HICE, IVISC, TISMODG, HWAT, TWT, &
WNO, WNR, WNI
CALL EXTCDE(3)
END IF
!
! Filter high ki
WNI = MIN(IC5MAXKI, WNI)
!
! FORMAT
1000 FORMAT(/' *** WAVEWATCH III ERROR IN W3SIC5MD : '/ &
' Subr. ', A, ': Zero shear modulus G is not allowed&
& in the FS viscoelastic model'/)
!
1001 FORMAT(/' *** WAVEWATCH III ERROR IN W3SIC5MD : '/ &
' Subr. ', A, ': get a bad first guess'/)
!
1002 FORMAT(/' *** WAVEWATCH III ERROR IN W3SIC5MD : '/ &
' -----------------------------------------------------'/&
' Subr. ', A,' : get NaN/NeG/Huge kr or ki for' /&
' -----------------------------------------------------'/&
' Ice thickness : ', F9.1, ' m'/ &
' Ice viscosity : ', E9.2, ' m2/s'/ &
' Ice shear modulus : ', E9.2, ' Pa' / &
' Water depth : ', F9.1, ' m'/ &
' Wave period : ', F10.2, ' s'/ &
' Wave number (Ko) : ', F11.3, ' rad/m'/ &
' Wave number (Kr) : ', F11.3, ' rad/m'/ &
' Attenu. Rate (Ki) : ', E9.2, ' /m'/)
!/
!/ End of FSDISP ----------------------------------------------------- /
!/
END SUBROUTINE FSDISP
!/ ------------------------------------------------------------------- /
!/
SUBROUTINE BALANCING_MATRIX(NMAT, MATRIX)
!/
!/ +-----------------------------------+
!/ | WAVEWATCH III NOAA/NCEP |
!/ | Q. Liu |
!/ | FORTRAN 90 |
!/ | Last update : 15-Mar-2016 |
!/ +-----------------------------------+
!/
!/ 15-Mar-2016 : Origination ( version 5.10)
!/ ( Q. Liu )
!/ 15-Mar-2016 : Borrowed from Numerical Recipes in Fortran
!/ ( Q. Liu )
! 1. Purpose :
! Reducing the sensitivity of eigenvalues to rounding errors during
! the execution of some algorithms.
!
! 2. Method :
! The errors in the eigensystem found by a numerical procedure are
! generally proportional to the Euclidean norm of the matrix, that
! is, to the square root of the sum of the squares of the elements
! (sqrt(sum(a_{i, j} ** 2.)). The idea of balancing is to use
! similarity transformations to make corresponding rows and columns
! of the matrix have comparable norms, thus reducing the overall
! norm of the matrix while leaving the eigenvalues unchanged. Note
! that the symmetric matrix is already balanced.
!
! The output is matrix that is balanced in the norm given by
! summing the absolute magnitudes of the matrix elements(
! sum(abs(a_{i, j})) ). This is more efficient than using the
! Euclidean norm, and equally effective: a large reduction in
! one norm implies a large reduction in the other.
!
! For the details of this method, please refer to
! 1) Numerical Recipes in Fortran 77 (Volume 1, 2nd Edition)
! [Chapter 11.5 / subroutine balanc]
! 2) Numerical Recipes in Fortran 90 (Volume 2)
! [Chapter B11 / subroutine balanc]
!
! 3. Parameters :
!
! Parameter list
! ----------------------------------------------------------------
! Name Type Intent Description
! ----------------------------------------------------------------
! NMAT Int. I The size of one dimension of MATRIX
! MATRIX R.A. I/O A matrix with the shape (NMAT, NMAT)
! ----------------------------------------------------------------
!
! 4. Subroutines used :
!
! Name Type Module Description
! ----------------------------------------------------------------
! STRACE Subr. W3SERVMD Subroutine tracing (!/S switch).
! ----------------------------------------------------------------
!
! 5. Called by :
!
! Name Type Module Description
! ----------------------------------------------------------------
! POLYROOTS Subr. / Find the roots of polynomials
! ----------------------------------------------------------------
!
! 6. Error messages :
!
! None.
!
! 7. Remarks:
! Balancing only needs marginal computational efforts but can
! substantially improve the accuracy of the eigenvalues computed
! for a badly balanced matrix. It is therefore recommended that
! you always balance nonsymmetric matrices.
!
! Given a (NMAT, NMAT) MATRIX, this routine replaces it by a
! balanced matrix with identical eigenvalues. A symmetric matrix is
! already balanced and is unaffected by this procedure.
!
! 8. Structure :
!
! See the source code.
!
! 9. Switches :
!
! !/S Enable subroutine tracing.
!
! 10. Source code :
!
!/ ------------------------------------------------------------------- /
!/S USE W3SERVMD, ONLY: STRACE
!
IMPLICIT NONE
!/
!/ ------------------------------------------------------------------- /
!/ Parameter list
INTEGER, INTENT(IN) :: NMAT
REAL, INTENT(INOUT) :: MATRIX(NMAT, NMAT)
!/ ------------------------------------------------------------------- /
!/ Local parameter
!/S INTEGER, SAVE :: IENT = 0
! the parameter radx is the machine's floating-point radix
REAL, PARAMETER :: RADX = RADIX(MATRIX), &
SQRADX = RADX ** 2
INTEGER :: I, LAST
REAL :: C, F, G, R, S
!/ ------------------------------------------------------------------- /
!/S CALL STRACE (IENT, 'BALANCING_MATRIX')
!
DO
LAST = 1
DO I = 1, NMAT
! Calculate row and column norms
C = SUM( ABS(MATRIX(:, I)) ) - MATRIX(I, I)
R = SUM( ABS(MATRIX(I, :)) ) - MATRIX(I, I)
! If both are non-zero
IF (C /= 0.0 .AND. R /= 0.0) THEN
! Find the integer power of the machine radix that comes closest to
! balancing the matrix (get G, F from C, R)
G = R / RADX
F = 1.0
S = C + R
DO
IF (C >= G) EXIT
F = F * RADX
C = C * SQRADX
END DO
!
G = R * RADX
DO
IF (C <= G) EXIT
F = F / RADX
C = C / SQRADX
END DO
!
IF ( (C+R)/F < 0.95*S) THEN
LAST = 0
G = 1.0 / F
! Apply similarity tranformation
MATRIX(I, :) = MATRIX(I, :) * G
MATRIX(:, I) = MATRIX(:, I) * F
END IF
END IF
END DO
IF (LAST /= 0) EXIT
END DO
!/
!/ End of subroutine BALANCING_MATRIX -------------------------------- /
!/
END SUBROUTINE BALANCING_MATRIX
!/ ------------------------------------------------------------------- /
!/
SUBROUTINE EIG_HQR (NMAT, HMAT, EIGR, EIGI)
!/
!/ +-----------------------------------+
!/ | WAVEWATCH III NOAA/NCEP |
!/ | Q. Liu |
!/ | FORTRAN 90 |
!/ | Last update : 17-Mar-2016 |
!/ +-----------------------------------+
!/
!/ 16-Mar-2016 : Origination ( version 5.10)
!/ ( Q. Liu )
!/ 16-Mar-2016 : Borrowed from Numerical Recipes in Fortran
!/ ( Q. Liu )
!/ 17-Mar-2016 : Update the NR code v2.08 to v2.10 ( Q. Liu )
!/
! 1. Purpose :
!
! When we calculate the eigenvalues of a general matrix, we first
! reduce the matrix to a simpler form (e.g., Hessenberg form) and
! then we perform the iterative procedures.
!
! A upper Hessenberg matrix has zeros everywhere below the diagnal
! except for the first subdiagonal row. For example, in the 6x6
! case, the non-zero elements are:
! |x x x x x x|
! |x x x x x x|
! | x x x x x|
! | x x x x|
! | x x x|
! | x x|
!
! This subroutine uses QR algorithm to get the eigenvalues of a
! Hessenberg matrix. So make sure the input array HMAT is a
! Hessenberg-type matrix.
!
! 2. Method :
! QR algorithm for real Hessenberg matrices.
! (I did not understand this algorithm well, so I could not give
! any detailed explanations)
!
! For the details of this HQR method, please refer to
! 1) Numerical Recipes in Fortran 77 (Volume 1, 2nd Edition)
! [Chapter 11.6 / subroutine hqr]
! 2) Numerical Recipes in Fortran 90 (Volume 2)