forked from mom-ocean/MOM6
-
Notifications
You must be signed in to change notification settings - Fork 62
/
Copy pathMOM_hor_visc.F90
2981 lines (2756 loc) · 147 KB
/
MOM_hor_visc.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
!> Calculates horizontal viscosity and viscous stresses
module MOM_hor_visc
! This file is part of MOM6. See LICENSE.md for the license.
use MOM_checksums, only : hchksum, Bchksum, uvchksum
use MOM_coms, only : min_across_PEs
use MOM_diag_mediator, only : post_data, register_diag_field, safe_alloc_ptr
use MOM_diag_mediator, only : post_product_u, post_product_sum_u
use MOM_diag_mediator, only : post_product_v, post_product_sum_v
use MOM_diag_mediator, only : diag_ctrl, time_type
use MOM_domains, only : pass_var, CORNER, pass_vector, AGRID, BGRID_NE
use MOM_domains, only : To_All, Scalar_Pair
use MOM_error_handler, only : MOM_error, FATAL, WARNING, is_root_pe
use MOM_file_parser, only : get_param, log_version, param_file_type
use MOM_grid, only : ocean_grid_type
use MOM_lateral_mixing_coeffs, only : VarMix_CS, calc_QG_Leith_viscosity
use MOM_barotropic, only : barotropic_CS, barotropic_get_tav
use MOM_thickness_diffuse, only : thickness_diffuse_CS, thickness_diffuse_get_KH
use MOM_io, only : MOM_read_data, slasher
use MOM_MEKE_types, only : MEKE_type
use MOM_open_boundary, only : ocean_OBC_type, OBC_DIRECTION_E, OBC_DIRECTION_W
use MOM_open_boundary, only : OBC_DIRECTION_N, OBC_DIRECTION_S, OBC_NONE
use MOM_unit_scaling, only : unit_scale_type
use MOM_verticalGrid, only : verticalGrid_type
use MOM_variables, only : accel_diag_ptrs
implicit none ; private
#include <MOM_memory.h>
public horizontal_viscosity, hor_visc_init, hor_visc_end
!> Control structure for horizontal viscosity
type, public :: hor_visc_CS ; private
logical :: initialized = .false. !< True if this control structure has been initialized.
logical :: Laplacian !< Use a Laplacian horizontal viscosity if true.
logical :: biharmonic !< Use a biharmonic horizontal viscosity if true.
logical :: debug !< If true, write verbose checksums for debugging purposes.
logical :: no_slip !< If true, no slip boundary conditions are used.
!! Otherwise free slip boundary conditions are assumed.
!! The implementation of the free slip boundary
!! conditions on a C-grid is much cleaner than the
!! no slip boundary conditions. The use of free slip
!! b.c.s is strongly encouraged. The no slip b.c.s
!! are not implemented with the biharmonic viscosity.
logical :: bound_Kh !< If true, the Laplacian coefficient is locally
!! limited to guarantee stability.
logical :: better_bound_Kh !< If true, use a more careful bounding of the
!! Laplacian viscosity to guarantee stability.
logical :: bound_Ah !< If true, the biharmonic coefficient is locally
!! limited to guarantee stability.
logical :: better_bound_Ah !< If true, use a more careful bounding of the
!! biharmonic viscosity to guarantee stability.
real :: Re_Ah !! If nonzero, the biharmonic coefficient is scaled
!< so that the biharmonic Reynolds number is equal to this.
real :: bound_coef !< The nondimensional coefficient of the ratio of
!! the viscosity bounds to the theoretical maximum
!! for stability without considering other terms [nondim].
!! The default is 0.8.
logical :: Smagorinsky_Kh !< If true, use Smagorinsky nonlinear eddy
!! viscosity. KH is the background value.
logical :: Smagorinsky_Ah !< If true, use a biharmonic form of Smagorinsky
!! nonlinear eddy viscosity. AH is the background.
logical :: Leith_Kh !< If true, use 2D Leith nonlinear eddy
!! viscosity. KH is the background value.
logical :: Modified_Leith !< If true, use extra component of Leith viscosity
!! to damp divergent flow. To use, still set Leith_Kh=.TRUE.
logical :: use_beta_in_Leith !< If true, includes the beta term in the Leith viscosity
logical :: Leith_Ah !< If true, use a biharmonic form of 2D Leith
!! nonlinear eddy viscosity. AH is the background.
logical :: use_QG_Leith_visc !< If true, use QG Leith nonlinear eddy viscosity.
!! KH is the background value.
logical :: bound_Coriolis !< If true & SMAGORINSKY_AH is used, the biharmonic
!! viscosity is modified to include a term that
!! scales quadratically with the velocity shears.
logical :: use_Kh_bg_2d !< Read 2d background viscosity from a file.
logical :: Kh_bg_2d_bug !< If true, retain an answer-changing horizontal indexing bug
!! in setting the corner-point viscosities when USE_KH_BG_2D=True.
real :: Kh_bg_min !< The minimum value allowed for Laplacian horizontal
!! viscosity [L2 T-1 ~> m2 s-1]. The default is 0.0.
logical :: use_land_mask !< Use the land mask for the computation of thicknesses
!! at velocity locations. This eliminates the dependence on
!! arbitrary values over land or outside of the domain.
!! Default is False to maintain answers with legacy experiments
!! but should be changed to True for new experiments.
logical :: anisotropic !< If true, allow anisotropic component to the viscosity.
logical :: add_LES_viscosity!< If true, adds the viscosity from Smagorinsky and Leith to
!! the background viscosity instead of taking the maximum.
real :: Kh_aniso !< The anisotropic viscosity [L2 T-1 ~> m2 s-1].
logical :: dynamic_aniso !< If true, the anisotropic viscosity is recomputed as a function
!! of state. This is set depending on ANISOTROPIC_MODE.
logical :: res_scale_MEKE !< If true, the viscosity contribution from MEKE is scaled by
!! the resolution function.
logical :: use_GME !< If true, use GME backscatter scheme.
logical :: answers_2018 !< If true, use the order of arithmetic and expressions that recover the
!! answers from the end of 2018. Otherwise, use updated and more robust
!! forms of the same expressions.
real :: GME_h0 !< The strength of GME tapers quadratically to zero when the bathymetric
!! depth is shallower than GME_H0 [Z ~> m]
real :: GME_efficiency !< The nondimensional prefactor multiplying the GME coefficient [nondim]
real :: GME_limiter !< The absolute maximum value the GME coefficient is allowed to take [L2 T-1 ~> m2 s-1].
real :: min_grid_Kh !< Minimum horizontal Laplacian viscosity used to
!! limit the grid Reynolds number [L2 T-1 ~> m2 s-1]
real :: min_grid_Ah !< Minimun horizontal biharmonic viscosity used to
!! limit grid Reynolds number [L4 T-1 ~> m4 s-1]
real ALLOCABLE_, dimension(NIMEM_,NJMEM_) :: Kh_bg_xx
!< The background Laplacian viscosity at h points [L2 T-1 ~> m2 s-1].
!! The actual viscosity may be the larger of this
!! viscosity and the Smagorinsky and Leith viscosities.
real ALLOCABLE_, dimension(NIMEM_,NJMEM_) :: Kh_bg_2d
!< The background Laplacian viscosity at h points [L2 T-1 ~> m2 s-1].
!! The actual viscosity may be the larger of this
!! viscosity and the Smagorinsky and Leith viscosities.
real ALLOCABLE_, dimension(NIMEM_,NJMEM_) :: Ah_bg_xx
!< The background biharmonic viscosity at h points [L4 T-1 ~> m4 s-1].
!! The actual viscosity may be the larger of this
!! viscosity and the Smagorinsky and Leith viscosities.
real ALLOCABLE_, dimension(NIMEM_,NJMEM_) :: reduction_xx
!< The amount by which stresses through h points are reduced
!! due to partial barriers [nondim].
real ALLOCABLE_, dimension(NIMEM_,NJMEM_) :: &
Kh_Max_xx, & !< The maximum permitted Laplacian viscosity [L2 T-1 ~> m2 s-1].
Ah_Max_xx, & !< The maximum permitted biharmonic viscosity [L4 T-1 ~> m4 s-1].
n1n2_h, & !< Factor n1*n2 in the anisotropic direction tensor at h-points
n1n1_m_n2n2_h, & !< Factor n1**2-n2**2 in the anisotropic direction tensor at h-points
grid_sp_h2, & !< Harmonic mean of the squares of the grid [L2 ~> m2]
grid_sp_h3 !< Harmonic mean of the squares of the grid^(3/2) [L3 ~> m3]
real ALLOCABLE_, dimension(NIMEMB_PTR_,NJMEMB_PTR_) :: Kh_bg_xy
!< The background Laplacian viscosity at q points [L2 T-1 ~> m2 s-1].
!! The actual viscosity may be the larger of this
!! viscosity and the Smagorinsky and Leith viscosities.
real ALLOCABLE_, dimension(NIMEMB_PTR_,NJMEMB_PTR_) :: Ah_bg_xy
!< The background biharmonic viscosity at q points [L4 T-1 ~> m4 s-1].
!! The actual viscosity may be the larger of this
!! viscosity and the Smagorinsky and Leith viscosities.
real ALLOCABLE_, dimension(NIMEMB_PTR_,NJMEMB_PTR_) :: reduction_xy
!< The amount by which stresses through q points are reduced
!! due to partial barriers [nondim].
real ALLOCABLE_, dimension(NIMEMB_PTR_,NJMEMB_PTR_) :: &
Kh_Max_xy, & !< The maximum permitted Laplacian viscosity [L2 T-1 ~> m2 s-1].
Ah_Max_xy, & !< The maximum permitted biharmonic viscosity [L4 T-1 ~> m4 s-1].
n1n2_q, & !< Factor n1*n2 in the anisotropic direction tensor at q-points
n1n1_m_n2n2_q !< Factor n1**2-n2**2 in the anisotropic direction tensor at q-points
real ALLOCABLE_, dimension(NIMEM_,NJMEM_) :: &
dx2h, & !< Pre-calculated dx^2 at h points [L2 ~> m2]
dy2h, & !< Pre-calculated dy^2 at h points [L2 ~> m2]
dx_dyT, & !< Pre-calculated dx/dy at h points [nondim]
dy_dxT !< Pre-calculated dy/dx at h points [nondim]
real ALLOCABLE_, dimension(NIMEMB_PTR_,NJMEMB_PTR_) :: &
dx2q, & !< Pre-calculated dx^2 at q points [L2 ~> m2]
dy2q, & !< Pre-calculated dy^2 at q points [L2 ~> m2]
dx_dyBu, & !< Pre-calculated dx/dy at q points [nondim]
dy_dxBu !< Pre-calculated dy/dx at q points [nondim]
real ALLOCABLE_, dimension(NIMEMB_PTR_,NJMEM_) :: &
Idx2dyCu, & !< 1/(dx^2 dy) at u points [L-3 ~> m-3]
Idxdy2u !< 1/(dx dy^2) at u points [L-3 ~> m-3]
real ALLOCABLE_, dimension(NIMEM_,NJMEMB_PTR_) :: &
Idx2dyCv, & !< 1/(dx^2 dy) at v points [L-3 ~> m-3]
Idxdy2v !< 1/(dx dy^2) at v points [L-3 ~> m-3]
! The following variables are precalculated time-invariant combinations of
! parameters and metric terms.
real ALLOCABLE_, dimension(NIMEM_,NJMEM_) :: &
Laplac2_const_xx, & !< Laplacian metric-dependent constants [L2 ~> m2]
Biharm6_const_xx, & !< Biharmonic metric-dependent constants [L6 ~> m6]
Laplac3_const_xx, & !< Laplacian metric-dependent constants [L3 ~> m3]
Biharm_const_xx, & !< Biharmonic metric-dependent constants [L4 ~> m4]
Biharm_const2_xx, & !< Biharmonic metric-dependent constants [T L4 ~> s m4]
Re_Ah_const_xx !< Biharmonic metric-dependent constants [L3 ~> m3]
real ALLOCABLE_, dimension(NIMEMB_PTR_,NJMEMB_PTR_) :: &
Laplac2_const_xy, & !< Laplacian metric-dependent constants [L2 ~> m2]
Biharm6_const_xy, & !< Biharmonic metric-dependent constants [L6 ~> m6]
Laplac3_const_xy, & !< Laplacian metric-dependent constants [L3 ~> m3]
Biharm_const_xy, & !< Biharmonic metric-dependent constants [L4 ~> m4]
Biharm_const2_xy, & !< Biharmonic metric-dependent constants [T L4 ~> s m4]
Re_Ah_const_xy !< Biharmonic metric-dependent constants [L3 ~> m3]
type(diag_ctrl), pointer :: diag => NULL() !< structure to regulate diagnostics
! real, allocatable :: hf_diffu(:,:,:) ! Zonal hor. visc. accel. x fract. thickness [L T-2 ~> m s-2].
! real, allocatable :: hf_diffv(:,:,:) ! Merdional hor. visc. accel. x fract. thickness [L T-2 ~> m s-2].
! 3D diagnostics hf_diffu(diffv) are commented because there is no clarity on proper remapping grid option.
! The code is retained for degugging purposes in the future.
integer :: num_smooth_gme !< number of smoothing passes for the GME fluxes.
!>@{
!! Diagnostic id
integer :: id_grid_Re_Ah = -1, id_grid_Re_Kh = -1
integer :: id_diffu = -1, id_diffv = -1
! integer :: id_hf_diffu = -1, id_hf_diffv = -1
integer :: id_h_diffu = -1, id_h_diffv = -1
integer :: id_hf_diffu_2d = -1, id_hf_diffv_2d = -1
integer :: id_intz_diffu_2d = -1, id_intz_diffv_2d = -1
integer :: id_diffu_visc_rem = -1, id_diffv_visc_rem = -1
integer :: id_Ah_h = -1, id_Ah_q = -1
integer :: id_Kh_h = -1, id_Kh_q = -1
integer :: id_GME_coeff_h = -1, id_GME_coeff_q = -1
integer :: id_dudx_bt = -1, id_dvdy_bt = -1
integer :: id_dudy_bt = -1, id_dvdx_bt = -1
integer :: id_vort_xy_q = -1, id_div_xx_h = -1
integer :: id_sh_xy_q = -1, id_sh_xx_h = -1
integer :: id_FrictWork = -1, id_FrictWorkIntz = -1
integer :: id_FrictWork_GME = -1
integer :: id_normstress = -1, id_shearstress = -1
!>@}
end type hor_visc_CS
contains
!> Calculates the acceleration due to the horizontal viscosity.
!!
!! A combination of biharmonic and Laplacian forms can be used. The coefficient
!! may either be a constant or a shear-dependent form. The biharmonic is
!! determined by twice taking the divergence of an appropriately defined stress
!! tensor. The Laplacian is determined by doing so once.
!!
!! To work, the following fields must be set outside of the usual
!! is:ie range before this subroutine is called:
!! u[is-2:ie+2,js-2:je+2]
!! v[is-2:ie+2,js-2:je+2]
!! h[is-1:ie+1,js-1:je+1]
subroutine horizontal_viscosity(u, v, h, diffu, diffv, MEKE, VarMix, G, GV, US, &
CS, OBC, BT, TD, ADp)
type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure.
type(verticalGrid_type), intent(in) :: GV !< The ocean's vertical grid structure.
real, dimension(SZIB_(G),SZJ_(G),SZK_(GV)), &
intent(in) :: u !< The zonal velocity [L T-1 ~> m s-1].
real, dimension(SZI_(G),SZJB_(G),SZK_(GV)), &
intent(in) :: v !< The meridional velocity [L T-1 ~> m s-1].
real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
intent(inout) :: h !< Layer thicknesses [H ~> m or kg m-2].
real, dimension(SZIB_(G),SZJ_(G),SZK_(GV)), &
intent(out) :: diffu !< Zonal acceleration due to convergence of
!! along-coordinate stress tensor [L T-2 ~> m s-2]
real, dimension(SZI_(G),SZJB_(G),SZK_(GV)), &
intent(out) :: diffv !< Meridional acceleration due to convergence
!! of along-coordinate stress tensor [L T-2 ~> m s-2].
type(MEKE_type), intent(inout) :: MEKE !< MEKE fields
!! related to Mesoscale Eddy Kinetic Energy.
type(VarMix_CS), intent(inout) :: VarMix !< Variable mixing control struct
type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
type(hor_visc_CS), intent(in) :: CS !< Horizontal viscosity control struct
type(ocean_OBC_type), optional, pointer :: OBC !< Pointer to an open boundary condition type
type(barotropic_CS), intent(in), optional :: BT !< Barotropic control struct
type(thickness_diffuse_CS), intent(in), optional :: TD !< Thickness diffusion control struct
type(accel_diag_ptrs), intent(in), optional :: ADp !< Acceleration diagnostics
! Local variables
real, dimension(SZIB_(G),SZJ_(G)) :: &
Del2u, & ! The u-compontent of the Laplacian of velocity [L-1 T-1 ~> m-1 s-1]
h_u, & ! Thickness interpolated to u points [H ~> m or kg m-2].
vort_xy_dy, & ! y-derivative of vertical vorticity (d/dy(dv/dx - du/dy)) [L-1 T-1 ~> m-1 s-1]
div_xx_dx, & ! x-derivative of horizontal divergence (d/dx(du/dx + dv/dy)) [L-1 T-1 ~> m-1 s-1]
ubtav ! zonal barotropic vel. ave. over baroclinic time-step [L T-1 ~> m s-1]
real, dimension(SZI_(G),SZJB_(G)) :: &
Del2v, & ! The v-compontent of the Laplacian of velocity [L-1 T-1 ~> m-1 s-1]
h_v, & ! Thickness interpolated to v points [H ~> m or kg m-2].
vort_xy_dx, & ! x-derivative of vertical vorticity (d/dx(dv/dx - du/dy)) [L-1 T-1 ~> m-1 s-1]
div_xx_dy, & ! y-derivative of horizontal divergence (d/dy(du/dx + dv/dy)) [L-1 T-1 ~> m-1 s-1]
vbtav ! meridional barotropic vel. ave. over baroclinic time-step [L T-1 ~> m s-1]
real, dimension(SZI_(G),SZJ_(G)) :: &
dudx_bt, dvdy_bt, & ! components in the barotropic horizontal tension [T-1 ~> s-1]
div_xx, & ! Estimate of horizontal divergence at h-points [T-1 ~> s-1]
sh_xx, & ! horizontal tension (du/dx - dv/dy) including metric terms [T-1 ~> s-1]
sh_xx_bt, & ! barotropic horizontal tension (du/dx - dv/dy) including metric terms [T-1 ~> s-1]
str_xx,& ! str_xx is the diagonal term in the stress tensor [H L2 T-2 ~> m3 s-2 or kg s-2]
str_xx_GME,& ! smoothed diagonal term in the stress tensor from GME [H L2 T-2 ~> m3 s-2 or kg s-2]
bhstr_xx, & ! A copy of str_xx that only contains the biharmonic contribution [H L2 T-2 ~> m3 s-2 or kg s-2]
FrictWorkIntz, & ! depth integrated energy dissipated by lateral friction [R L2 T-3 ~> W m-2]
grad_vort_mag_h, & ! Magnitude of vorticity gradient at h-points [L-1 T-1 ~> m-1 s-1]
grad_vort_mag_h_2d, & ! Magnitude of 2d vorticity gradient at h-points [L-1 T-1 ~> m-1 s-1]
Del2vort_h, & ! Laplacian of vorticity at h-points [L-2 T-1 ~> m-2 s-1]
grad_div_mag_h, & ! Magnitude of divergence gradient at h-points [L-1 T-1 ~> m-1 s-1]
dudx, dvdy, & ! components in the horizontal tension [T-1 ~> s-1]
grad_vel_mag_h, & ! Magnitude of the velocity gradient tensor squared at h-points [T-2 ~> s-2]
grad_vel_mag_bt_h, & ! Magnitude of the barotropic velocity gradient tensor squared at h-points [T-2 ~> s-2]
grad_d2vel_mag_h, & ! Magnitude of the Laplacian of the velocity vector, squared [L-2 T-2 ~> m-2 s-2]
GME_effic_h, & ! The filtered efficiency of the GME terms at h points [nondim]
htot, & ! The total thickness of all layers [Z ~> m]
boundary_mask_h ! A mask that zeroes out cells with at least one land edge [nondim]
real, allocatable, dimension(:,:,:) :: h_diffu ! h x diffu [H L T-2 ~> m2 s-2]
real, allocatable, dimension(:,:,:) :: h_diffv ! h x diffv [H L T-2 ~> m2 s-2]
real, allocatable, dimension(:,:,:) :: diffu_visc_rem ! diffu x visc_rem_u [L T-2 ~> m s-2]
real, allocatable, dimension(:,:,:) :: diffv_visc_rem ! diffv x visc_rem_v [L T-2 ~> m s-2]
real, dimension(SZIB_(G),SZJB_(G)) :: &
dvdx, dudy, & ! components in the shearing strain [T-1 ~> s-1]
dDel2vdx, dDel2udy, & ! Components in the biharmonic equivalent of the shearing strain [L-2 T-1 ~> m-2 s-1]
dvdx_bt, dudy_bt, & ! components in the barotropic shearing strain [T-1 ~> s-1]
sh_xy, & ! horizontal shearing strain (du/dy + dv/dx) including metric terms [T-1 ~> s-1]
sh_xy_bt, & ! barotropic horizontal shearing strain (du/dy + dv/dx) inc. metric terms [T-1 ~> s-1]
str_xy, & ! str_xy is the cross term in the stress tensor [H L2 T-2 ~> m3 s-2 or kg s-2]
str_xy_GME, & ! smoothed cross term in the stress tensor from GME [H L2 T-2 ~> m3 s-2 or kg s-2]
bhstr_xy, & ! A copy of str_xy that only contains the biharmonic contribution [H L2 T-2 ~> m3 s-2 or kg s-2]
vort_xy, & ! Vertical vorticity (dv/dx - du/dy) including metric terms [T-1 ~> s-1]
Leith_Kh_q, & ! Leith Laplacian viscosity at q-points [L2 T-1 ~> m2 s-1]
Leith_Ah_q, & ! Leith bi-harmonic viscosity at q-points [L4 T-1 ~> m4 s-1]
grad_vort_mag_q, & ! Magnitude of vorticity gradient at q-points [L-1 T-1 ~> m-1 s-1]
grad_vort_mag_q_2d, & ! Magnitude of 2d vorticity gradient at q-points [L-1 T-1 ~> m-1 s-1]
Del2vort_q, & ! Laplacian of vorticity at q-points [L-2 T-1 ~> m-2 s-1]
grad_div_mag_q, & ! Magnitude of divergence gradient at q-points [L-1 T-1 ~> m-1 s-1]
grad_vel_mag_q, & ! Magnitude of the velocity gradient tensor squared at q-points [T-2 ~> s-2]
hq, & ! harmonic mean of the harmonic means of the u- & v point thicknesses [H ~> m or kg m-2]
! This form guarantees that hq/hu < 4.
grad_vel_mag_bt_q, & ! Magnitude of the barotropic velocity gradient tensor squared at q-points [T-2 ~> s-2]
GME_effic_q, & ! The filtered efficiency of the GME terms at q points [nondim]
boundary_mask_q ! A mask that zeroes out cells with at least one land edge [nondim]
real, dimension(SZIB_(G),SZJB_(G),SZK_(GV)) :: &
Ah_q, & ! biharmonic viscosity at corner points [L4 T-1 ~> m4 s-1]
Kh_q, & ! Laplacian viscosity at corner points [L2 T-1 ~> m2 s-1]
vort_xy_q, & ! vertical vorticity at corner points [T-1 ~> s-1]
sh_xy_q, & ! horizontal shearing strain at corner points [T-1 ~> s-1]
GME_coeff_q, & !< GME coeff. at q-points [L2 T-1 ~> m2 s-1]
ShSt ! A diagnostic array of shear stress [T-1 ~> s-1].
real, dimension(SZIB_(G),SZJ_(G),SZK_(GV)+1) :: &
KH_u_GME !< interface height diffusivities in u-columns [L2 T-1 ~> m2 s-1]
real, dimension(SZI_(G),SZJB_(G),SZK_(GV)+1) :: &
KH_v_GME !< interface height diffusivities in v-columns [L2 T-1 ~> m2 s-1]
real, dimension(SZI_(G),SZJ_(G),SZK_(GV)) :: &
Ah_h, & ! biharmonic viscosity at thickness points [L4 T-1 ~> m4 s-1]
Kh_h, & ! Laplacian viscosity at thickness points [L2 T-1 ~> m2 s-1]
FrictWork, & ! work done by MKE dissipation mechanisms [R L2 T-3 ~> W m-2]
FrictWork_GME, & ! work done by GME [R L2 T-3 ~> W m-2]
div_xx_h, & ! horizontal divergence [T-1 ~> s-1]
sh_xx_h, & ! horizontal tension (du/dx - dv/dy) including metric terms [T-1 ~> s-1]
NoSt ! A diagnostic array of normal stress [T-1 ~> s-1].
real, dimension(SZI_(G),SZJ_(G),SZK_(G)) :: &
grid_Re_Kh, & ! Grid Reynolds number for Laplacian horizontal viscosity at h points [nondim]
grid_Re_Ah, & ! Grid Reynolds number for Biharmonic horizontal viscosity at h points [nondim]
GME_coeff_h ! GME coeff. at h-points [L2 T-1 ~> m2 s-1]
real :: AhSm ! Smagorinsky biharmonic viscosity [L4 T-1 ~> m4 s-1]
real :: AhLth ! 2D Leith biharmonic viscosity [L4 T-1 ~> m4 s-1]
real :: mod_Leith ! nondimensional coefficient for divergence part of modified Leith
! viscosity. Here set equal to nondimensional Laplacian Leith constant.
! This is set equal to zero if modified Leith is not used.
real :: Shear_mag_bc ! Shear_mag value in backscatter [T-1 ~> s-1]
real :: sh_xx_sq ! Square of tension (sh_xx) [T-2 ~> s-2]
real :: sh_xy_sq ! Square of shearing strain (sh_xy) [T-2 ~> s-2]
real :: h2uq, h2vq ! temporary variables [H2 ~> m2 or kg2 m-4].
real :: hu, hv ! Thicknesses interpolated by arithmetic means to corner
! points; these are first interpolated to u or v velocity
! points where masks are applied [H ~> m or kg m-2].
real :: h_arith_q ! The arithmetic mean total thickness at q points [Z ~> m]
real :: h_harm_q ! The harmonic mean total thickness at q points [Z ~> m]
real :: I_hq ! The inverse of the arithmetic mean total thickness at q points [Z-1 ~> m-1]
real :: I_GME_h0 ! The inverse of GME tapering scale [Z-1 ~> m-1]
real :: h_neglect ! thickness so small it can be lost in roundoff and so neglected [H ~> m or kg m-2]
real :: h_neglect3 ! h_neglect^3 [H3 ~> m3 or kg3 m-6]
real :: h_min ! Minimum h at the 4 neighboring velocity points [H ~> m]
real :: Kh_scale ! A factor between 0 and 1 by which the horizontal
! Laplacian viscosity is rescaled [nondim]
real :: RoScl ! The scaling function for MEKE source term [nondim]
real :: FatH ! abs(f) at h-point for MEKE source term [T-1 ~> s-1]
real :: local_strain ! Local variable for interpolating computed strain rates [T-1 ~> s-1].
real :: meke_res_fn ! A copy of the resolution scaling factor if being applied to MEKE. Otherwise =1.
real :: GME_coeff ! The GME (negative) viscosity coefficient [L2 T-1 ~> m2 s-1]
real :: GME_coeff_limiter ! Maximum permitted value of the GME coefficient [L2 T-1 ~> m2 s-1]
real :: FWfrac ! Fraction of maximum theoretical energy transfer to use when scaling GME coefficient [nondim]
real :: DY_dxBu ! Ratio of meridional over zonal grid spacing at vertices [nondim]
real :: DX_dyBu ! Ratio of zonal over meridiononal grid spacing at vertices [nondim]
real :: DY_dxCv ! Ratio of meridional over zonal grid spacing at faces [nondim]
real :: DX_dyCu ! Ratio of zonal over meridional grid spacing at faces [nondim]
real :: Sh_F_pow ! The ratio of shear over the absolute value of f raised to some power and rescaled [nondim]
real :: backscat_subround ! The ratio of f over Shear_mag that is so small that the backscatter
! calculation gives the same value as if f were 0 [nondim].
real :: H0_GME ! Depth used to scale down GME coefficient in shallow areas [Z ~> m]
real :: KE ! Local kinetic energy [L2 T-2 ~> m2 s-2]
real :: d_del2u ! dy-weighted Laplacian(u) diff in x [L-2 T-1 ~> m-2 s-1]
real :: d_del2v ! dx-weighted Laplacian(v) diff in y [L-2 T-1 ~> m-2 s-1]
real :: d_str ! Stress tensor update [H L2 T-2 ~> m3 s-2 or kg s-2]
real :: grad_vort ! Vorticity gradient magnitude [L-1 T-1 ~> m-1 s-1]
real :: grad_vort_qg ! QG-based vorticity gradient magnitude [L-1 T-1 ~> m-1 s-1]
real :: grid_Kh ! Laplacian viscosity bound by grid [L2 T-1 ~> m2 s-1]
real :: grid_Ah ! Biharmonic viscosity bound by grid [L4 T-1 ~> m4 s-1]
logical :: rescale_Kh, legacy_bound
logical :: find_FrictWork
logical :: apply_OBC = .false.
logical :: use_MEKE_Ku
logical :: use_MEKE_Au
integer :: is, ie, js, je, Isq, Ieq, Jsq, Jeq, nz
integer :: i, j, k, n
real :: inv_PI3, inv_PI2, inv_PI6
! Fields evaluated on active layers, used for constructing 3D stress fields
! NOTE: The position of these declarations can impact performance, due to the
! very large number of stack arrays in this function. Move with caution!
real, dimension(SZIB_(G),SZJB_(G)) :: &
Ah, & ! biharmonic viscosity (h or q) [L4 T-1 ~> m4 s-1]
Kh, & ! Laplacian viscosity [L2 T-1 ~> m2 s-1]
Shear_mag, & ! magnitude of the shear [T-1 ~> s-1]
vert_vort_mag, & ! magnitude of the vertical vorticity gradient [L-1 T-1 ~> m-1 s-1]
hrat_min, & ! h_min divided by the thickness at the stress point (h or q) [nondim]
visc_bound_rem ! fraction of overall viscous bounds that remain to be applied [nondim]
real, dimension(SZIB_(G),SZJ_(G)) :: &
hf_diffu_2d, & ! Depth sum of hf_diffu, hf_diffv [L T-2 ~> m s-2]
intz_diffu_2d ! Depth-integral of diffu [H L T-2 ~> m2 s-2]
real, dimension(SZI_(G),SZJB_(G)) :: &
hf_diffv_2d, & ! Depth sum of hf_diffu, hf_diffv [L T-2 ~> m s-2]
intz_diffv_2d ! Depth-integral of diffv [H L T-2 ~> m2 s-2]
is = G%isc ; ie = G%iec ; js = G%jsc ; je = G%jec ; nz = GV%ke
Isq = G%IscB ; Ieq = G%IecB ; Jsq = G%JscB ; Jeq = G%JecB
h_neglect = GV%H_subroundoff
h_neglect3 = h_neglect**3
inv_PI3 = 1.0/((4.0*atan(1.0))**3)
inv_PI2 = 1.0/((4.0*atan(1.0))**2)
inv_PI6 = inv_PI3 * inv_PI3
if (present(OBC)) then ; if (associated(OBC)) then ; if (OBC%OBC_pe) then
apply_OBC = OBC%Flather_u_BCs_exist_globally .or. OBC%Flather_v_BCs_exist_globally
apply_OBC = .true.
endif ; endif ; endif
if (.not.CS%initialized) call MOM_error(FATAL, &
"MOM_hor_visc: Module must be initialized before it is used.")
if (.not.(CS%Laplacian .or. CS%biharmonic)) return
find_FrictWork = (CS%id_FrictWork > 0)
if (CS%id_FrictWorkIntz > 0) find_FrictWork = .true.
if (allocated(MEKE%mom_src)) find_FrictWork = .true.
backscat_subround = 0.0
if (find_FrictWork .and. allocated(MEKE%mom_src) .and. (MEKE%backscatter_Ro_c > 0.0) .and. &
(MEKE%backscatter_Ro_Pow /= 0.0)) &
backscat_subround = (1.0e-16/MEKE%backscatter_Ro_c)**(1.0/MEKE%backscatter_Ro_Pow)
! Toggle whether to use a Laplacian viscosity derived from MEKE
use_MEKE_Ku = allocated(MEKE%Ku)
use_MEKE_Au = allocated(MEKE%Au)
rescale_Kh = .false.
if (VarMix%use_variable_mixing) then
rescale_Kh = VarMix%Resoln_scaled_Kh
if ((rescale_Kh .or. CS%res_scale_MEKE) &
.and. (.not. allocated(VarMix%Res_fn_h) .or. .not. allocated(VarMix%Res_fn_q))) &
call MOM_error(FATAL, "MOM_hor_visc: VarMix%Res_fn_h and VarMix%Res_fn_q "//&
"both need to be associated with Resoln_scaled_Kh or RES_SCALE_MEKE_VISC.")
elseif (CS%res_scale_MEKE) then
call MOM_error(FATAL, "MOM_hor_visc: VarMix needs to be associated if "//&
"RES_SCALE_MEKE_VISC is True.")
endif
legacy_bound = (CS%Smagorinsky_Kh .or. CS%Leith_Kh) .and. &
(CS%bound_Kh .and. .not.CS%better_bound_Kh)
if (CS%use_GME) then
do j=Jsq-1,Jeq+2 ; do i=Isq-1,Ieq+2
boundary_mask_h(i,j) = (G%mask2dCu(I,j) * G%mask2dCv(i,J) * G%mask2dCu(I-1,j) * G%mask2dCv(i,J-1))
enddo ; enddo
do J=Jsq-2,Jeq+1 ; do I=Isq-2,Ieq+1
boundary_mask_q(I,J) = (G%mask2dCv(i,J) * G%mask2dCv(i+1,J) * G%mask2dCu(I,j) * G%mask2dCu(I,j+1))
enddo ; enddo
! initialize diag. array with zeros
GME_coeff_h(:,:,:) = 0.0
GME_coeff_q(:,:,:) = 0.0
str_xx_GME(:,:) = 0.0
str_xy_GME(:,:) = 0.0
! Get barotropic velocities and their gradients
call barotropic_get_tav(BT, ubtav, vbtav, G, US)
call pass_vector(ubtav, vbtav, G%Domain)
! Calculate the barotropic horizontal tension
do j=Jsq-2,Jeq+2 ; do i=Isq-2,Ieq+2
dudx_bt(i,j) = CS%DY_dxT(i,j)*(G%IdyCu(I,j) * ubtav(I,j) - &
G%IdyCu(I-1,j) * ubtav(I-1,j))
dvdy_bt(i,j) = CS%DX_dyT(i,j)*(G%IdxCv(i,J) * vbtav(i,J) - &
G%IdxCv(i,J-1) * vbtav(i,J-1))
sh_xx_bt(i,j) = dudx_bt(i,j) - dvdy_bt(i,j)
enddo ; enddo
! Components for the barotropic shearing strain
do J=Jsq-2,Jeq+2 ; do I=Isq-2,Ieq+2
dvdx_bt(I,J) = CS%DY_dxBu(I,J)*(vbtav(i+1,J)*G%IdyCv(i+1,J) &
- vbtav(i,J)*G%IdyCv(i,J))
dudy_bt(I,J) = CS%DX_dyBu(I,J)*(ubtav(I,j+1)*G%IdxCu(I,j+1) &
- ubtav(I,j)*G%IdxCu(I,j))
enddo ; enddo
! post barotropic tension and strain
if (CS%id_dudx_bt > 0) call post_data(CS%id_dudx_bt, dudx_bt, CS%diag)
if (CS%id_dvdy_bt > 0) call post_data(CS%id_dvdy_bt, dvdy_bt, CS%diag)
if (CS%id_dudy_bt > 0) call post_data(CS%id_dudy_bt, dudy_bt, CS%diag)
if (CS%id_dvdx_bt > 0) call post_data(CS%id_dvdx_bt, dvdx_bt, CS%diag)
if (CS%no_slip) then
do J=Jsq-2,Jeq+2 ; do I=Isq-2,Ieq+2
sh_xy_bt(I,J) = (2.0-G%mask2dBu(I,J)) * ( dvdx_bt(I,J) + dudy_bt(I,J) )
enddo ; enddo
else
do J=Jsq-2,Jeq+2 ; do I=Isq-2,Ieq+2
sh_xy_bt(I,J) = G%mask2dBu(I,J) * ( dvdx_bt(I,J) + dudy_bt(I,J) )
enddo ; enddo
endif
do j=Jsq-1,Jeq+2 ; do i=Isq-1,Ieq+2
grad_vel_mag_bt_h(i,j) = G%mask2dT(I,J) * boundary_mask_h(i,j) * (dudx_bt(i,j)**2 + dvdy_bt(i,j)**2 + &
(0.25*((dvdx_bt(I,J)+dvdx_bt(I-1,J-1))+(dvdx_bt(I,J-1)+dvdx_bt(I-1,J))))**2 + &
(0.25*((dudy_bt(I,J)+dudy_bt(I-1,J-1))+(dudy_bt(I,J-1)+dudy_bt(I-1,J))))**2)
enddo ; enddo
do J=Jsq-2,Jeq+1 ; do I=Isq-2,Ieq+1
grad_vel_mag_bt_q(I,J) = G%mask2dBu(I,J) * boundary_mask_q(I,J) * (dvdx_bt(I,J)**2 + dudy_bt(I,J)**2 + &
(0.25*((dudx_bt(i,j)+dudx_bt(i+1,j+1))+(dudx_bt(i,j+1)+dudx_bt(i+1,j))))**2 + &
(0.25*((dvdy_bt(i,j)+dvdy_bt(i+1,j+1))+(dvdy_bt(i,j+1)+dvdy_bt(i+1,j))))**2)
enddo ; enddo
call pass_var(h, G%domain, halo=2)
do j=js-2,je+2 ; do i=is-2,ie+2
htot(i,j) = 0.0
enddo ; enddo
do k=1,nz ; do j=js-2,je+2 ; do i=is-2,ie+2
htot(i,j) = htot(i,j) + GV%H_to_Z*h(i,j,k)
enddo ; enddo ; enddo
I_GME_h0 = 1.0 / CS%GME_h0
do j=Jsq-1,Jeq+2 ; do i=Isq-1,Ieq+2
if (grad_vel_mag_bt_h(i,j)>0) then
GME_effic_h(i,j) = CS%GME_efficiency * G%mask2dT(I,J) * &
(MIN(htot(i,j) * I_GME_h0, 1.0)**2)
else
GME_effic_h(i,j) = 0.0
endif
enddo ; enddo
do J=Jsq-2,Jeq+1 ; do I=Isq-2,Ieq+1
if (grad_vel_mag_bt_q(I,J)>0) then
h_arith_q = 0.25 * ((htot(i,j) + htot(i+1,j+1)) + (htot(i+1,j) + htot(i,j+1)))
I_hq = 1.0 / h_arith_q
h_harm_q = 0.25 * h_arith_q * ((htot(i,j)*I_hq + htot(i+1,j+1)*I_hq) + &
(htot(i+1,j)*I_hq + htot(i,j+1)*I_hq))
GME_effic_q(I,J) = CS%GME_efficiency * G%mask2dBu(I,J) * (MIN(h_harm_q * I_GME_h0, 1.0)**2)
else
GME_effic_q(I,J) = 0.0
endif
enddo ; enddo
call thickness_diffuse_get_KH(TD, KH_u_GME, KH_v_GME, G, GV)
call pass_vector(KH_u_GME, KH_v_GME, G%domain, To_All+Scalar_Pair)
if (CS%debug) &
call uvchksum("GME KH[u,v]_GME", KH_u_GME, KH_v_GME, G%HI, haloshift=2, scale=US%L_to_m**2*US%s_to_T)
endif ! use_GME
!$OMP parallel do default(none) &
!$OMP shared( &
!$OMP CS, G, GV, US, OBC, VarMix, MEKE, u, v, h, &
!$OMP is, ie, js, je, Isq, Ieq, Jsq, Jeq, nz, &
!$OMP apply_OBC, rescale_Kh, legacy_bound, find_FrictWork, &
!$OMP use_MEKE_Ku, use_MEKE_Au, boundary_mask_h, boundary_mask_q, &
!$OMP backscat_subround, GME_coeff_limiter, GME_effic_h, GME_effic_q, &
!$OMP h_neglect, h_neglect3, FWfrac, inv_PI3, inv_PI6, H0_GME, &
!$OMP diffu, diffv, Kh_h, Kh_q, Ah_h, Ah_q, FrictWork, FrictWork_GME, &
!$OMP div_xx_h, sh_xx_h, vort_xy_q, sh_xy_q, GME_coeff_h, GME_coeff_q, &
!$OMP KH_u_GME, KH_v_GME, grid_Re_Kh, grid_Re_Ah, NoSt, ShSt &
!$OMP ) &
!$OMP private( &
!$OMP i, j, k, n, &
!$OMP dudx, dudy, dvdx, dvdy, sh_xx, sh_xy, h_u, h_v, &
!$OMP Del2u, Del2v, DY_dxBu, DX_dyBu, sh_xx_bt, sh_xy_bt, &
!$OMP str_xx, str_xy, bhstr_xx, bhstr_xy, str_xx_GME, str_xy_GME, &
!$OMP vort_xy, vort_xy_dx, vort_xy_dy, div_xx, div_xx_dx, div_xx_dy, &
!$OMP grad_div_mag_h, grad_div_mag_q, grad_vort_mag_h, grad_vort_mag_q, &
!$OMP grad_vort, grad_vort_qg, grad_vort_mag_h_2d, grad_vort_mag_q_2d, &
!$OMP grad_vel_mag_h, grad_vel_mag_q, sh_xx_sq, sh_xy_sq, &
!$OMP grad_vel_mag_bt_h, grad_vel_mag_bt_q, grad_d2vel_mag_h, &
!$OMP meke_res_fn, Shear_mag, Shear_mag_bc, vert_vort_mag, h_min, hrat_min, visc_bound_rem, &
!$OMP grid_Ah, grid_Kh, d_Del2u, d_Del2v, d_str, &
!$OMP Kh, Ah, AhSm, AhLth, local_strain, Sh_F_pow, &
!$OMP dDel2vdx, dDel2udy, DY_dxCv, DX_dyCu, Del2vort_q, Del2vort_h, KE, &
!$OMP h2uq, h2vq, hu, hv, hq, FatH, RoScl, GME_coeff &
!$OMP )
do k=1,nz
! The following are the forms of the horizontal tension and horizontal
! shearing strain advocated by Smagorinsky (1993) and discussed in
! Griffies and Hallberg (2000).
! NOTE: There is a ~1% speedup when the tension and shearing loops below
! are fused (presumably due to shared access of Id[xy]C[uv]). However,
! this breaks the center/vertex index case convention, and also evaluates
! the dudx and dvdy terms beyond their valid bounds.
! TODO: Explore methods for retaining both the syntax and speedup.
! Calculate horizontal tension
do j=Jsq-1,Jeq+2 ; do i=Isq-1,Ieq+2
dudx(i,j) = CS%DY_dxT(i,j)*(G%IdyCu(I,j) * u(I,j,k) - &
G%IdyCu(I-1,j) * u(I-1,j,k))
dvdy(i,j) = CS%DX_dyT(i,j)*(G%IdxCv(i,J) * v(i,J,k) - &
G%IdxCv(i,J-1) * v(i,J-1,k))
sh_xx(i,j) = dudx(i,j) - dvdy(i,j)
enddo ; enddo
! Components for the shearing strain
do J=Jsq-2,Jeq+2 ; do I=Isq-2,Ieq+2
dvdx(I,J) = CS%DY_dxBu(I,J)*(v(i+1,J,k)*G%IdyCv(i+1,J) - v(i,J,k)*G%IdyCv(i,J))
dudy(I,J) = CS%DX_dyBu(I,J)*(u(I,j+1,k)*G%IdxCu(I,j+1) - u(I,j,k)*G%IdxCu(I,j))
enddo ; enddo
if (CS%id_normstress > 0) then
do j=Jsq-1,Jeq+2 ; do i=Isq-1,Ieq+2
NoSt(i,j,k) = sh_xx(i,j)
enddo ; enddo
endif
! Interpolate the thicknesses to velocity points.
! The extra wide halos are to accommodate the cross-corner-point projections
! in OBCs, which are not ordinarily be necessary, and might not be necessary
! even with OBCs if the accelerations are zeroed at OBC points, in which
! case the j-loop for h_u could collapse to j=js=1,je+1. -RWH
if (CS%use_land_mask) then
do j=js-2,je+2 ; do I=Isq-1,Ieq+1
h_u(I,j) = 0.5 * (G%mask2dT(i,j)*h(i,j,k) + G%mask2dT(i+1,j)*h(i+1,j,k))
enddo ; enddo
do J=Jsq-1,Jeq+1 ; do i=is-2,ie+2
h_v(i,J) = 0.5 * (G%mask2dT(i,j)*h(i,j,k) + G%mask2dT(i,j+1)*h(i,j+1,k))
enddo ; enddo
else
do j=js-2,je+2 ; do I=Isq-1,Ieq+1
h_u(I,j) = 0.5 * (h(i,j,k) + h(i+1,j,k))
enddo ; enddo
do J=Jsq-1,Jeq+1 ; do i=is-2,ie+2
h_v(i,J) = 0.5 * (h(i,j,k) + h(i,j+1,k))
enddo ; enddo
endif
! Adjust contributions to shearing strain and interpolated values of
! thicknesses on open boundaries.
if (apply_OBC) then ; do n=1,OBC%number_of_segments
J = OBC%segment(n)%HI%JsdB ; I = OBC%segment(n)%HI%IsdB
if (OBC%zero_strain .or. OBC%freeslip_strain .or. OBC%computed_strain) then
if (OBC%segment(n)%is_N_or_S .and. (J >= js-2) .and. (J <= Jeq+1)) then
do I=OBC%segment(n)%HI%IsdB,OBC%segment(n)%HI%IedB
if (OBC%zero_strain) then
dvdx(I,J) = 0. ; dudy(I,J) = 0.
elseif (OBC%freeslip_strain) then
dudy(I,J) = 0.
elseif (OBC%computed_strain) then
if (OBC%segment(n)%direction == OBC_DIRECTION_N) then
dudy(I,J) = 2.0*CS%DX_dyBu(I,J)* &
(OBC%segment(n)%tangential_vel(I,J,k) - u(I,j,k))*G%IdxCu(I,j)
else
dudy(I,J) = 2.0*CS%DX_dyBu(I,J)* &
(u(I,j+1,k) - OBC%segment(n)%tangential_vel(I,J,k))*G%IdxCu(I,j+1)
endif
elseif (OBC%specified_strain) then
if (OBC%segment(n)%direction == OBC_DIRECTION_N) then
dudy(I,J) = CS%DX_dyBu(I,J)*OBC%segment(n)%tangential_grad(I,J,k)*G%IdxCu(I,j)*G%dxBu(I,J)
else
dudy(I,J) = CS%DX_dyBu(I,J)*OBC%segment(n)%tangential_grad(I,J,k)*G%IdxCu(I,j+1)*G%dxBu(I,J)
endif
endif
enddo
elseif (OBC%segment(n)%is_E_or_W .and. (I >= is-2) .and. (I <= Ieq+1)) then
do J=OBC%segment(n)%HI%JsdB,OBC%segment(n)%HI%JedB
if (OBC%zero_strain) then
dvdx(I,J) = 0. ; dudy(I,J) = 0.
elseif (OBC%freeslip_strain) then
dvdx(I,J) = 0.
elseif (OBC%computed_strain) then
if (OBC%segment(n)%direction == OBC_DIRECTION_E) then
dvdx(I,J) = 2.0*CS%DY_dxBu(I,J)* &
(OBC%segment(n)%tangential_vel(I,J,k) - v(i,J,k))*G%IdyCv(i,J)
else
dvdx(I,J) = 2.0*CS%DY_dxBu(I,J)* &
(v(i+1,J,k) - OBC%segment(n)%tangential_vel(I,J,k))*G%IdyCv(i+1,J)
endif
elseif (OBC%specified_strain) then
if (OBC%segment(n)%direction == OBC_DIRECTION_E) then
dvdx(I,J) = CS%DY_dxBu(I,J)*OBC%segment(n)%tangential_grad(I,J,k)*G%IdyCv(i,J)*G%dxBu(I,J)
else
dvdx(I,J) = CS%DY_dxBu(I,J)*OBC%segment(n)%tangential_grad(I,J,k)*G%IdyCv(i+1,J)*G%dxBu(I,J)
endif
endif
enddo
endif
endif
if (OBC%segment(n)%direction == OBC_DIRECTION_N) then
! There are extra wide halos here to accommodate the cross-corner-point
! OBC projections, but they might not be necessary if the accelerations
! are always zeroed out at OBC points, in which case the i-loop below
! becomes do i=is-1,ie+1. -RWH
if ((J >= Jsq-1) .and. (J <= Jeq+1)) then
do i = max(is-2,OBC%segment(n)%HI%isd), min(ie+2,OBC%segment(n)%HI%ied)
h_v(i,J) = h(i,j,k)
enddo
endif
elseif (OBC%segment(n)%direction == OBC_DIRECTION_S) then
if ((J >= Jsq-1) .and. (J <= Jeq+1)) then
do i = max(is-2,OBC%segment(n)%HI%isd), min(ie+2,OBC%segment(n)%HI%ied)
h_v(i,J) = h(i,j+1,k)
enddo
endif
elseif (OBC%segment(n)%direction == OBC_DIRECTION_E) then
if ((I >= Isq-1) .and. (I <= Ieq+1)) then
do j = max(js-2,OBC%segment(n)%HI%jsd), min(je+2,OBC%segment(n)%HI%jed)
h_u(I,j) = h(i,j,k)
enddo
endif
elseif (OBC%segment(n)%direction == OBC_DIRECTION_W) then
if ((I >= Isq-1) .and. (I <= Ieq+1)) then
do j = max(js-2,OBC%segment(n)%HI%jsd), min(je+2,OBC%segment(n)%HI%jed)
h_u(I,j) = h(i+1,j,k)
enddo
endif
endif
enddo ; endif
! Now project thicknesses across corner points on OBCs.
if (apply_OBC) then ; do n=1,OBC%number_of_segments
J = OBC%segment(n)%HI%JsdB ; I = OBC%segment(n)%HI%IsdB
if (OBC%segment(n)%direction == OBC_DIRECTION_N) then
if ((J >= js-2) .and. (J <= je)) then
do I = max(Isq-1,OBC%segment(n)%HI%IsdB), min(Ieq+1,OBC%segment(n)%HI%IedB)
h_u(I,j+1) = h_u(I,j)
enddo
endif
elseif (OBC%segment(n)%direction == OBC_DIRECTION_S) then
if ((J >= js-1) .and. (J <= je+1)) then
do I = max(Isq-1,OBC%segment(n)%HI%isd), min(Ieq+1,OBC%segment(n)%HI%ied)
h_u(I,j) = h_u(I,j+1)
enddo
endif
elseif (OBC%segment(n)%direction == OBC_DIRECTION_E) then
if ((I >= is-2) .and. (I <= ie)) then
do J = max(Jsq-1,OBC%segment(n)%HI%jsd), min(Jeq+1,OBC%segment(n)%HI%jed)
h_v(i+1,J) = h_v(i,J)
enddo
endif
elseif (OBC%segment(n)%direction == OBC_DIRECTION_W) then
if ((I >= is-1) .and. (I <= ie+1)) then
do J = max(Jsq-1,OBC%segment(n)%HI%jsd), min(Jeq+1,OBC%segment(n)%HI%jed)
h_v(i,J) = h_v(i+1,J)
enddo
endif
endif
enddo ; endif
! Shearing strain (including no-slip boundary conditions at the 2-D land-sea mask).
! dudy and dvdx include modifications at OBCs from above.
if (CS%no_slip) then
do J=js-2,Jeq+1 ; do I=is-2,Ieq+1
sh_xy(I,J) = (2.0-G%mask2dBu(I,J)) * ( dvdx(I,J) + dudy(I,J) )
if (CS%id_shearstress > 0) ShSt(I,J,k) = sh_xy(I,J)
enddo ; enddo
else
do J=js-2,Jeq+1 ; do I=is-2,Ieq+1
sh_xy(I,J) = G%mask2dBu(I,J) * ( dvdx(I,J) + dudy(I,J) )
if (CS%id_shearstress > 0) ShSt(I,J,k) = sh_xy(I,J)
enddo ; enddo
endif
! Evaluate Del2u = x.Div(Grad u) and Del2v = y.Div( Grad u)
if (CS%biharmonic) then
do j=js-1,Jeq+1 ; do I=Isq-1,Ieq+1
Del2u(I,j) = CS%Idxdy2u(I,j)*(CS%dy2h(i+1,j)*sh_xx(i+1,j) - CS%dy2h(i,j)*sh_xx(i,j)) + &
CS%Idx2dyCu(I,j)*(CS%dx2q(I,J)*sh_xy(I,J) - CS%dx2q(I,J-1)*sh_xy(I,J-1))
enddo ; enddo
do J=Jsq-1,Jeq+1 ; do i=is-1,Ieq+1
Del2v(i,J) = CS%Idxdy2v(i,J)*(CS%dy2q(I,J)*sh_xy(I,J) - CS%dy2q(I-1,J)*sh_xy(I-1,J)) - &
CS%Idx2dyCv(i,J)*(CS%dx2h(i,j+1)*sh_xx(i,j+1) - CS%dx2h(i,j)*sh_xx(i,j))
enddo ; enddo
if (apply_OBC) then ; if (OBC%zero_biharmonic) then
do n=1,OBC%number_of_segments
I = OBC%segment(n)%HI%IsdB ; J = OBC%segment(n)%HI%JsdB
if (OBC%segment(n)%is_N_or_S .and. (J >= Jsq-1) .and. (J <= Jeq+1)) then
do I=OBC%segment(n)%HI%isd,OBC%segment(n)%HI%ied
Del2v(i,J) = 0.
enddo
elseif (OBC%segment(n)%is_E_or_W .and. (I >= Isq-1) .and. (I <= Ieq+1)) then
do j=OBC%segment(n)%HI%jsd,OBC%segment(n)%HI%jed
Del2u(I,j) = 0.
enddo
endif
enddo
endif ; endif
endif
! Vorticity
if (CS%no_slip) then
do J=Jsq-2,Jeq+2 ; do I=Isq-2,Ieq+2
vort_xy(I,J) = (2.0-G%mask2dBu(I,J)) * ( dvdx(I,J) - dudy(I,J) )
enddo ; enddo
else
do J=Jsq-2,Jeq+2 ; do I=Isq-2,Ieq+2
vort_xy(I,J) = G%mask2dBu(I,J) * ( dvdx(I,J) - dudy(I,J) )
enddo ; enddo
endif
! Divergence
do j=Jsq-1,Jeq+2 ; do i=Isq-1,Ieq+2
div_xx(i,j) = dudx(i,j) + dvdy(i,j)
enddo ; enddo
if ((CS%Leith_Kh) .or. (CS%Leith_Ah)) then
! Vorticity gradient
do J=Jsq-1,Jeq+1 ; do i=Isq-1,Ieq+2
DY_dxBu = G%dyBu(I,J) * G%IdxBu(I,J)
vort_xy_dx(i,J) = DY_dxBu * (vort_xy(I,J) * G%IdyCu(I,j) - vort_xy(I-1,J) * G%IdyCu(I-1,j))
enddo ; enddo
do j=Jsq-1,Jeq+2 ; do I=Isq-1,Ieq+1
DX_dyBu = G%dxBu(I,J) * G%IdyBu(I,J)
vort_xy_dy(I,j) = DX_dyBu * (vort_xy(I,J) * G%IdxCv(i,J) - vort_xy(I,J-1) * G%IdxCv(i,J-1))
enddo ; enddo
! Laplacian of vorticity
do J=Jsq-1,Jeq+1 ; do I=Isq-1,Ieq+1
DY_dxBu = G%dyBu(I,J) * G%IdxBu(I,J)
DX_dyBu = G%dxBu(I,J) * G%IdyBu(I,J)
Del2vort_q(I,J) = DY_dxBu * (vort_xy_dx(i+1,J) * G%IdyCv(i+1,J) - vort_xy_dx(i,J) * G%IdyCv(i,J)) + &
DX_dyBu * (vort_xy_dy(I,j+1) * G%IdyCu(I,j+1) - vort_xy_dy(I,j) * G%IdyCu(I,j))
enddo ; enddo
do J=Jsq,Jeq+1 ; do I=Isq,Ieq+1
Del2vort_h(i,j) = 0.25*(Del2vort_q(I,J) + Del2vort_q(I-1,J) + Del2vort_q(I,J-1) + Del2vort_q(I-1,J-1))
enddo ; enddo
if (CS%modified_Leith) then
! Divergence gradient
do j=Jsq-1,Jeq+2 ; do I=Isq-1,Ieq+1
div_xx_dx(I,j) = G%IdxCu(I,j)*(div_xx(i+1,j) - div_xx(i,j))
enddo ; enddo
do J=Jsq-1,Jeq+1 ; do i=Isq-1,Ieq+2
div_xx_dy(i,J) = G%IdyCv(i,J)*(div_xx(i,j+1) - div_xx(i,j))
enddo ; enddo
! Magnitude of divergence gradient
do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
grad_div_mag_h(i,j) = sqrt((0.5*(div_xx_dx(I,j) + div_xx_dx(I-1,j)))**2 + &
(0.5*(div_xx_dy(i,J) + div_xx_dy(i,J-1)))**2)
enddo ; enddo
do j=Jsq-1,Jeq+1 ; do i=Isq-1,Ieq+1
grad_div_mag_q(I,J) = sqrt((0.5*(div_xx_dx(I,j) + div_xx_dx(I,j+1)))**2 + &
(0.5*(div_xx_dy(i,J) + div_xx_dy(i+1,J)))**2)
enddo ; enddo
else
do j=Jsq-1,Jeq+2 ; do I=is-2,Ieq+1
div_xx_dx(I,j) = 0.0
enddo ; enddo
do J=Jsq-1,Jeq+1 ; do i=Isq-1,Ieq+2
div_xx_dy(i,J) = 0.0
enddo ; enddo
do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
grad_div_mag_h(i,j) = 0.0
enddo ; enddo
do J=Jsq-1,Jeq+1 ; do I=Isq-1,Ieq+1
grad_div_mag_q(I,J) = 0.0
enddo ; enddo
endif ! CS%modified_Leith
! Add in beta for the Leith viscosity
if (CS%use_beta_in_Leith) then
do J=js-2,Jeq+1 ; do i=is-1,Ieq+1
vort_xy_dx(i,J) = vort_xy_dx(i,J) + 0.5 * ( G%dF_dx(i,j) + G%dF_dx(i,j+1))
enddo ; enddo
do j=js-1,Jeq+1 ; do I=is-2,Ieq+1
vort_xy_dy(I,j) = vort_xy_dy(I,j) + 0.5 * ( G%dF_dy(i,j) + G%dF_dy(i+1,j))
enddo ; enddo
endif ! CS%use_beta_in_Leith
if (CS%use_QG_Leith_visc) then
do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
grad_vort_mag_h_2d(i,j) = SQRT((0.5*(vort_xy_dx(i,J) + vort_xy_dx(i,J-1)))**2 + &
(0.5*(vort_xy_dy(I,j) + vort_xy_dy(I-1,j)))**2 )
enddo ; enddo
do J=js-1,Jeq ; do I=is-1,Ieq
grad_vort_mag_q_2d(I,J) = SQRT((0.5*(vort_xy_dx(i,J) + vort_xy_dx(i+1,J)))**2 + &
(0.5*(vort_xy_dy(I,j) + vort_xy_dy(I,j+1)))**2 )
enddo ; enddo
! This accumulates terms, some of which are in VarMix, so rescaling can not be done here.
call calc_QG_Leith_viscosity(VarMix, G, GV, US, h, k, div_xx_dx, div_xx_dy, &
vort_xy_dx, vort_xy_dy)
endif
do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
grad_vort_mag_h(i,j) = SQRT((0.5*(vort_xy_dx(i,J) + vort_xy_dx(i,J-1)))**2 + &
(0.5*(vort_xy_dy(I,j) + vort_xy_dy(I-1,j)))**2 )
enddo ; enddo
do J=js-1,Jeq ; do I=is-1,Ieq
grad_vort_mag_q(I,J) = SQRT((0.5*(vort_xy_dx(i,J) + vort_xy_dx(i+1,J)))**2 + &
(0.5*(vort_xy_dy(I,j) + vort_xy_dy(I,j+1)))**2 )
enddo ; enddo
endif ! CS%Leith_Kh
if ((CS%Smagorinsky_Kh) .or. (CS%Smagorinsky_Ah)) then
do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
sh_xx_sq = sh_xx(i,j)**2
sh_xy_sq = 0.25 * ( (sh_xy(I-1,J-1)**2 + sh_xy(I,J)**2) &
+ (sh_xy(I-1,J)**2 + sh_xy(I,J-1)**2) )
Shear_mag(i,j) = sqrt(sh_xx_sq + sh_xy_sq)
enddo ; enddo
endif
if (CS%better_bound_Ah .or. CS%better_bound_Kh) then
do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
h_min = min(h_u(I,j), h_u(I-1,j), h_v(i,J), h_v(i,J-1))
hrat_min(i,j) = min(1.0, h_min / (h(i,j,k) + h_neglect))
enddo ; enddo
if (CS%better_bound_Kh) then
do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
visc_bound_rem(i,j) = 1.0
enddo ; enddo
endif
endif
if (CS%Laplacian) then
if ((CS%Leith_Kh) .or. (CS%Leith_Ah)) then
if (CS%use_QG_Leith_visc) then
do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
grad_vort = grad_vort_mag_h(i,j) + grad_div_mag_h(i,j)
grad_vort_qg = 3. * grad_vort_mag_h_2d(i,j)
vert_vort_mag(i,j) = min(grad_vort, grad_vort_qg)
enddo ; enddo
else
do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
vert_vort_mag(i,j) = grad_vort_mag_h(i,j) + grad_div_mag_h(i,j)
enddo ; enddo
endif
endif
! Determine the Laplacian viscosity at h points, using the
! largest value from several parameterizations.
! Static (pre-computed) background viscosity
do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
Kh(i,j) = CS%Kh_bg_xx(i,j)
enddo ; enddo
! NOTE: The following do-block can be decomposed and vectorized after the
! stack size has been reduced.
do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
if (CS%add_LES_viscosity) then
if (CS%Smagorinsky_Kh) &
Kh(i,j) = Kh(i,j) + CS%Laplac2_const_xx(i,j) * Shear_mag(i,j)
if (CS%Leith_Kh) &
Kh(i,j) = Kh(i,j) + CS%Laplac3_const_xx(i,j) * vert_vort_mag(i,j) * inv_PI3
else
if (CS%Smagorinsky_Kh) &
Kh(i,j) = max(Kh(i,j), CS%Laplac2_const_xx(i,j) * Shear_mag(i,j))
if (CS%Leith_Kh) &
Kh(i,j) = max(Kh(i,j), CS%Laplac3_const_xx(i,j) * vert_vort_mag(i,j) * inv_PI3)
endif
enddo ; enddo
! All viscosity contributions above are subject to resolution scaling
if (rescale_Kh) then
do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
Kh(i,j) = VarMix%Res_fn_h(i,j) * Kh(i,j)
enddo ; enddo
endif
if (legacy_bound) then
! Older method of bounding for stability
do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
Kh(i,j) = min(Kh(i,j), CS%Kh_Max_xx(i,j))
enddo ; enddo
endif
! Place a floor on the viscosity, if desired.
do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
Kh(i,j) = max(Kh(i,j), CS%Kh_bg_min)
enddo ; enddo
if (use_MEKE_Ku) then
! *Add* the MEKE contribution (which might be negative)
if (CS%res_scale_MEKE) then
do j=Jsq,Jeq+1 ; do i=Isq,Ieq+1
Kh(i,j) = Kh(i,j) + MEKE%Ku(i,j) * VarMix%Res_fn_h(i,j)
enddo ; enddo
else