forked from mom-ocean/MOM6
-
Notifications
You must be signed in to change notification settings - Fork 63
/
Copy pathMOM_thickness_diffuse.F90
2529 lines (2297 loc) · 131 KB
/
MOM_thickness_diffuse.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
!> Isopycnal height diffusion (or Gent McWilliams diffusion)
module MOM_thickness_diffuse
! This file is part of MOM6. See LICENSE.md for the license.
use MOM_debugging, only : hchksum, uvchksum
use MOM_diag_mediator, only : post_data, query_averaging_enabled, diag_ctrl
use MOM_diag_mediator, only : register_diag_field, safe_alloc_ptr, time_type
use MOM_diag_mediator, only : diag_update_remap_grids
use MOM_domains, only : pass_var, CORNER, pass_vector
use MOM_error_handler, only : MOM_error, FATAL, WARNING, is_root_pe
use MOM_EOS, only : calculate_density, calculate_density_derivs, EOS_domain
use MOM_EOS, only : calculate_density_second_derivs
use MOM_file_parser, only : get_param, log_version, param_file_type
use MOM_grid, only : ocean_grid_type
use MOM_io, only : MOM_read_data, slasher
use MOM_interface_heights, only : find_eta, thickness_to_dz
use MOM_isopycnal_slopes, only : vert_fill_TS
use MOM_lateral_mixing_coeffs, only : VarMix_CS
use MOM_MEKE_types, only : MEKE_type
use MOM_unit_scaling, only : unit_scale_type
use MOM_variables, only : thermo_var_ptrs, cont_diag_ptrs
use MOM_verticalGrid, only : verticalGrid_type
implicit none ; private
#include <MOM_memory.h>
public thickness_diffuse, thickness_diffuse_init, thickness_diffuse_end
public thickness_diffuse_get_KH
! A note on unit descriptions in comments: MOM6 uses units that can be rescaled for dimensional
! consistency testing. These are noted in comments with units like Z, H, L, and T, along with
! their mks counterparts with notation like "a velocity [Z T-1 ~> m s-1]". If the units
! vary with the Boussinesq approximation, the Boussinesq variant is given first.
!> Control structure for thickness_diffuse
type, public :: thickness_diffuse_CS ; private
logical :: initialized = .false. !< True if this control structure has been initialized.
real :: Khth !< Background isopycnal depth diffusivity [L2 T-1 ~> m2 s-1]
real :: Khth_Slope_Cff !< Slope dependence coefficient of Khth [nondim]
real :: max_Khth_CFL !< Maximum value of the diffusive CFL for isopycnal height diffusion [nondim]
real :: Khth_Min !< Minimum value of Khth [L2 T-1 ~> m2 s-1]
real :: Khth_Max !< Maximum value of Khth [L2 T-1 ~> m2 s-1], or 0 for no max
real :: Kh_eta_bg !< Background isopycnal height diffusivity [L2 T-1 ~> m2 s-1]
real :: Kh_eta_vel !< Velocity scale that is multiplied by the grid spacing to give
!! the isopycnal height diffusivity [L T-1 ~> m s-1]
real :: slope_max !< Slopes steeper than slope_max are limited in some way [Z L-1 ~> nondim]
real :: kappa_smooth !< Vertical diffusivity used to interpolate more sensible values
!! of T & S into thin layers [H Z T-1 ~> m2 s-1 or kg m-1 s-1]
logical :: thickness_diffuse !< If true, interfaces heights are diffused.
logical :: use_FGNV_streamfn !< If true, use the streamfunction formulation of
!! Ferrari et al., 2010, which effectively emphasizes
!! graver vertical modes by smoothing in the vertical.
real :: FGNV_scale !< A coefficient scaling the vertical smoothing term in the
!! Ferrari et al., 2010, streamfunction formulation [nondim].
real :: FGNV_c_min !< A minimum wave speed used in the Ferrari et al., 2010,
!! streamfunction formulation [L T-1 ~> m s-1].
real :: N2_floor !< A floor for squared buoyancy frequency in the Ferrari et al., 2010,
!! streamfunction formulation [T-2 ~> s-2].
logical :: detangle_interfaces !< If true, add 3-d structured interface height
!! diffusivities to horizontally smooth jagged layers.
real :: detangle_time !< If detangle_interfaces is true, this is the
!! timescale over which maximally jagged grid-scale
!! thickness variations are suppressed [T ~> s]. This must be
!! longer than DT, or 0 (the default) to use DT.
integer :: nkml !< number of layers within mixed layer
logical :: debug !< write verbose checksums for debugging purposes
logical :: use_GME_thickness_diffuse !< If true, passes GM coefficients to MOM_hor_visc for use
!! with GME closure.
logical :: MEKE_GEOMETRIC !< If true, uses the GM coefficient formulation from the GEOMETRIC
!! framework (Marshall et al., 2012)
real :: MEKE_GEOMETRIC_alpha!< The nondimensional coefficient governing the efficiency of
!! the GEOMETRIC isopycnal height diffusion [nondim]
real :: MEKE_GEOMETRIC_epsilon !< Minimum Eady growth rate for the GEOMETRIC thickness
!! diffusivity [T-1 ~> s-1].
integer :: MEKE_GEOM_answer_date !< The vintage of the expressions in the MEKE_GEOMETRIC
!! calculation. Values below 20190101 recover the answers from the
!! original implementation, while higher values use expressions that
!! satisfy rotational symmetry.
logical :: Use_KH_in_MEKE !< If true, uses the isopycnal height diffusivity calculated here to diffuse MEKE.
real :: MEKE_min_depth_diff !< The minimum total depth over which to average the diffusivity
!! used for MEKE [H ~> m or kg m-2]. When the total depth is less
!! than this, the diffusivity is scaled away.
logical :: GM_src_alt !< If true, use the GM energy conversion form S^2*N^2*kappa rather
!! than the streamfunction for the GM source term for MEKE.
integer :: MEKE_src_answer_date !< The vintage of the expressions in the GM energy conversion
!! calculation when MEKE_GM_SRC_ALT is true. Values below 20240601
!! recover the answers from the original implementation, while higher
!! values use expressions that satisfy rotational symmetry.
logical :: MEKE_src_slope_bug !< If true, use a bug that limits the positive values, but not the
!! negative values, of the slopes used when MEKE_GM_SRC_ALT is true.
!! When this is true, it breaks rotational symmetry.
logical :: use_GM_work_bug !< If true, use the incorrect sign for the
!! top-level work tendency on the top layer.
real :: Stanley_det_coeff !< The coefficient correlating SGS temperature variance with the mean
!! temperature gradient in the deterministic part of the Stanley parameterization.
!! Negative values disable the scheme. [nondim]
logical :: read_khth !< If true, read a file containing the spatially varying horizontal
!! isopycnal height diffusivity
logical :: use_stanley_gm !< If true, also use the Stanley parameterization in MOM_thickness_diffuse
type(diag_ctrl), pointer :: diag => NULL() !< structure used to regulate timing of diagnostics
real, allocatable :: GMwork(:,:) !< Work by isopycnal height diffusion [R Z L2 T-3 ~> W m-2]
real, allocatable :: diagSlopeX(:,:,:) !< Diagnostic: zonal neutral slope [Z L-1 ~> nondim]
real, allocatable :: diagSlopeY(:,:,:) !< Diagnostic: zonal neutral slope [Z L-1 ~> nondim]
real, allocatable :: Kh_eta_u(:,:) !< Isopycnal height diffusivities at u points [L2 T-1 ~> m2 s-1]
real, allocatable :: Kh_eta_v(:,:) !< Isopycnal height diffusivities in v points [L2 T-1 ~> m2 s-1]
real, allocatable :: KH_u_GME(:,:,:) !< Isopycnal height diffusivities in u-columns [L2 T-1 ~> m2 s-1]
real, allocatable :: KH_v_GME(:,:,:) !< Isopycnal height diffusivities in v-columns [L2 T-1 ~> m2 s-1]
real, allocatable :: khth2d(:,:) !< 2D isopycnal height diffusivity at h-points [L2 T-1 ~> m2 s-1]
!>@{
!! Diagnostic identifier
integer :: id_uhGM = -1, id_vhGM = -1, id_GMwork = -1
integer :: id_KH_u = -1, id_KH_v = -1, id_KH_t = -1
integer :: id_KH_u1 = -1, id_KH_v1 = -1, id_KH_t1 = -1
integer :: id_slope_x = -1, id_slope_y = -1
integer :: id_sfn_unlim_x = -1, id_sfn_unlim_y = -1, id_sfn_x = -1, id_sfn_y = -1
!>@}
end type thickness_diffuse_CS
contains
!> Calculates isopycnal height diffusion coefficients and applies isopycnal height diffusion
!! by modifying to the layer thicknesses, h. Diffusivities are limited to ensure stability.
!! Also returns along-layer mass fluxes used in the continuity equation.
subroutine thickness_diffuse(h, uhtr, vhtr, tv, dt, G, GV, US, MEKE, VarMix, CDp, CS)
type(ocean_grid_type), intent(in) :: G !< Ocean grid structure
type(verticalGrid_type), intent(in) :: GV !< Vertical grid structure
type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(inout) :: h !< Layer thickness [H ~> m or kg m-2]
real, dimension(SZIB_(G),SZJ_(G),SZK_(GV)), intent(inout) :: uhtr !< Accumulated zonal mass flux
!! [L2 H ~> m3 or kg]
real, dimension(SZI_(G),SZJB_(G),SZK_(GV)), intent(inout) :: vhtr !< Accumulated meridional mass flux
!! [L2 H ~> m3 or kg]
type(thermo_var_ptrs), intent(in) :: tv !< Thermodynamics structure
real, intent(in) :: dt !< Time increment [T ~> s]
type(MEKE_type), intent(inout) :: MEKE !< MEKE fields
type(VarMix_CS), target, intent(in) :: VarMix !< Variable mixing coefficients
type(cont_diag_ptrs), intent(inout) :: CDp !< Diagnostics for the continuity equation
type(thickness_diffuse_CS), intent(inout) :: CS !< Control structure for thickness_diffuse
! Local variables
real :: e(SZI_(G),SZJ_(G),SZK_(GV)+1) ! heights of interfaces, relative to mean
! sea level [Z ~> m], positive up.
real :: uhD(SZIB_(G),SZJ_(G),SZK_(GV)) ! Diffusive u*h fluxes [L2 H T-1 ~> m3 s-1 or kg s-1]
real :: vhD(SZI_(G),SZJB_(G),SZK_(GV)) ! Diffusive v*h fluxes [L2 H T-1 ~> m3 s-1 or kg s-1]
real, dimension(SZIB_(G),SZJ_(G),SZK_(GV)+1) :: &
KH_u, & ! Isopycnal height diffusivities in u-columns [L2 T-1 ~> m2 s-1]
int_slope_u ! A nondimensional ratio from 0 to 1 that gives the relative
! weighting of the interface slopes to that calculated also
! using density gradients at u points. The physically correct
! slopes occur at 0, while 1 is used for numerical closures [nondim].
real, dimension(SZI_(G),SZJB_(G),SZK_(GV)+1) :: &
KH_v, & ! Isopycnal height diffusivities in v-columns [L2 T-1 ~> m2 s-1]
int_slope_v ! A nondimensional ratio from 0 to 1 that gives the relative
! weighting of the interface slopes to that calculated also
! using density gradients at v points. The physically correct
! slopes occur at 0, while 1 is used for numerical closures [nondim].
real, dimension(SZI_(G),SZJ_(G),SZK_(GV)) :: &
KH_t ! diagnosed diffusivity at tracer points [L2 T-1 ~> m2 s-1]
real, dimension(SZIB_(G),SZJ_(G)) :: &
KH_u_CFL ! The maximum stable isopycnal height diffusivity at u grid points [L2 T-1 ~> m2 s-1]
real, dimension(SZI_(G),SZJB_(G)) :: &
KH_v_CFL ! The maximum stable isopycnal height diffusivity at v grid points [L2 T-1 ~> m2 s-1]
real, dimension(SZI_(G),SZJ_(G)) :: &
htot ! The sum of the total layer thicknesses [H ~> m or kg m-2]
real :: Khth_Loc_u(SZIB_(G),SZJ_(G)) ! The isopycnal height diffusivity at u points [L2 T-1 ~> m2 s-1]
real :: Khth_Loc_v(SZI_(G),SZJB_(G)) ! The isopycnal height diffusivity at v points [L2 T-1 ~> m2 s-1]
real :: h_neglect ! A thickness that is so small it is usually lost
! in roundoff and can be neglected [H ~> m or kg m-2].
real, dimension(:,:), pointer :: cg1 => null() !< Wave speed [L T-1 ~> m s-1]
real :: hu(SZI_(G),SZJ_(G)) ! A thickness-based mask at u points, used for diagnostics [nondim]
real :: hv(SZI_(G),SZJ_(G)) ! A thickness-based mask at v points, used for diagnostics [nondim]
real :: KH_u_lay(SZI_(G),SZJ_(G)) ! Diagnostic of isopycnal height diffusivities at u-points averaged
! to layer centers [L2 T-1 ~> m2 s-1]
real :: KH_v_lay(SZI_(G),SZJ_(G)) ! Diagnostic of isopycnal height diffusivities at v-points averaged
! to layer centers [L2 T-1 ~> m2 s-1]
logical :: use_VarMix, Resoln_scaled, Depth_scaled, use_stored_slopes, khth_use_ebt_struct, use_Visbeck
logical :: use_QG_Leith
integer :: i, j, k, is, ie, js, je, nz
if (.not. CS%initialized) call MOM_error(FATAL, "MOM_thickness_diffuse: "//&
"Module must be initialized before it is used.")
if ((.not.CS%thickness_diffuse) &
.or. .not. (CS%Khth > 0.0 .or. CS%read_khth &
.or. VarMix%use_variable_mixing)) return
is = G%isc ; ie = G%iec ; js = G%jsc ; je = G%jec ; nz = GV%ke
h_neglect = GV%H_subroundoff
if (allocated(MEKE%GM_src)) then
do j=js,je ; do i=is,ie ; MEKE%GM_src(i,j) = 0. ; enddo ; enddo
endif
use_VarMix = .false. ; Resoln_scaled = .false. ; use_stored_slopes = .false.
khth_use_ebt_struct = .false. ; use_Visbeck = .false. ; use_QG_Leith = .false.
Depth_scaled = .false.
if (VarMix%use_variable_mixing) then
use_VarMix = VarMix%use_variable_mixing .and. (CS%KHTH_Slope_Cff > 0.)
Resoln_scaled = VarMix%Resoln_scaled_KhTh
Depth_scaled = VarMix%Depth_scaled_KhTh
use_stored_slopes = VarMix%use_stored_slopes
khth_use_ebt_struct = VarMix%khth_use_ebt_struct
use_Visbeck = VarMix%use_Visbeck
use_QG_Leith = VarMix%use_QG_Leith_GM
if (allocated(VarMix%cg1)) cg1 => VarMix%cg1
else
cg1 => null()
endif
!$OMP parallel do default(shared)
do j=js,je ; do I=is-1,ie
KH_u_CFL(I,j) = (0.25*CS%max_Khth_CFL) / &
(dt * (G%IdxCu(I,j)*G%IdxCu(I,j) + G%IdyCu(I,j)*G%IdyCu(I,j)))
enddo ; enddo
!$OMP parallel do default(shared)
do j=js-1,je ; do I=is,ie
KH_v_CFL(i,J) = (0.25*CS%max_Khth_CFL) / &
(dt * (G%IdxCv(i,J)*G%IdxCv(i,J) + G%IdyCv(i,J)*G%IdyCv(i,J)))
enddo ; enddo
! Calculates interface heights, e, in [Z ~> m].
call find_eta(h, tv, G, GV, US, e, halo_size=1)
! Set the diffusivities.
!$OMP parallel default(shared)
if (.not. CS%read_khth) then
!$OMP do
do j=js,je ; do I=is-1,ie
Khth_loc_u(I,j) = CS%Khth
enddo ; enddo
else ! use 2d KHTH that was read in from file
!$OMP do
do j=js,je ; do I=is-1,ie
Khth_loc_u(I,j) = 0.5 * (CS%khth2d(i,j) + CS%khth2d(i+1,j))
enddo ; enddo
endif
if (use_VarMix) then
if (use_Visbeck) then
!$OMP do
do j=js,je ; do I=is-1,ie
Khth_loc_u(I,j) = Khth_loc_u(I,j) + &
CS%KHTH_Slope_Cff*VarMix%L2u(I,j) * VarMix%SN_u(I,j)
enddo ; enddo
endif
endif
if (allocated(MEKE%Kh)) then
if (CS%MEKE_GEOMETRIC) then
!$OMP do
do j=js,je ; do I=is-1,ie
Khth_loc_u(I,j) = Khth_loc_u(I,j) + G%OBCmaskCu(I,j) * CS%MEKE_GEOMETRIC_alpha * &
0.5*(MEKE%MEKE(i,j)+MEKE%MEKE(i+1,j)) / &
(VarMix%SN_u(I,j) + CS%MEKE_GEOMETRIC_epsilon)
enddo ; enddo
else
do j=js,je ; do I=is-1,ie
Khth_loc_u(I,j) = Khth_loc_u(I,j) + MEKE%KhTh_fac*sqrt(MEKE%Kh(i,j)*MEKE%Kh(i+1,j))
enddo ; enddo
endif
endif
if (Resoln_scaled) then
!$OMP do
do j=js,je ; do I=is-1,ie
Khth_loc_u(I,j) = Khth_loc_u(I,j) * VarMix%Res_fn_u(I,j)
enddo ; enddo
endif
if (Depth_scaled) then
!$OMP do
do j=js,je ; do I=is-1,ie
Khth_loc_u(I,j) = Khth_loc_u(I,j) * VarMix%Depth_fn_u(I,j)
enddo ; enddo
endif
if (CS%Khth_Max > 0) then
!$OMP do
do j=js,je ; do I=is-1,ie
Khth_loc_u(I,j) = max(CS%Khth_Min, min(Khth_loc_u(I,j), CS%Khth_Max))
enddo ; enddo
else
!$OMP do
do j=js,je ; do I=is-1,ie
Khth_loc_u(I,j) = max(CS%Khth_Min, Khth_loc_u(I,j))
enddo ; enddo
endif
!$OMP do
do j=js,je ; do I=is-1,ie
KH_u(I,j,1) = min(KH_u_CFL(I,j), Khth_loc_u(I,j))
enddo ; enddo
if (khth_use_ebt_struct) then
!$OMP do
do K=2,nz+1 ; do j=js,je ; do I=is-1,ie
KH_u(I,j,K) = KH_u(I,j,1) * 0.5 * ( VarMix%ebt_struct(i,j,k-1) + VarMix%ebt_struct(i+1,j,k-1) )
enddo ; enddo ; enddo
else
!$OMP do
do K=2,nz+1 ; do j=js,je ; do I=is-1,ie
KH_u(I,j,K) = KH_u(I,j,1)
enddo ; enddo ; enddo
endif
if (use_VarMix) then
if (use_QG_Leith) then
!$OMP do
do k=1,nz ; do j=js,je ; do I=is-1,ie
KH_u(I,j,k) = VarMix%KH_u_QG(I,j,k)
enddo ; enddo ; enddo
endif
endif
if (CS%use_GME_thickness_diffuse) then
!$OMP do
do k=1,nz+1 ; do j=js,je ; do I=is-1,ie
CS%KH_u_GME(I,j,k) = KH_u(I,j,k)
enddo ; enddo ; enddo
endif
if (.not. CS%read_khth) then
!$OMP do
do J=js-1,je ; do i=is,ie
Khth_loc_v(i,J) = CS%Khth
enddo ; enddo
else ! read KHTH from file
!$OMP do
do J=js-1,je ; do i=is,ie
Khth_loc_v(i,J) = 0.5 * (CS%khth2d(i,j) + CS%khth2d(i,j+1))
enddo ; enddo
endif
if (use_VarMix) then
if (use_Visbeck) then
!$OMP do
do J=js-1,je ; do i=is,ie
Khth_loc_v(i,J) = Khth_loc_v(i,J) + CS%KHTH_Slope_Cff*VarMix%L2v(i,J)*VarMix%SN_v(i,J)
enddo ; enddo
endif
endif
if (allocated(MEKE%Kh)) then
if (CS%MEKE_GEOMETRIC) then
!$OMP do
do J=js-1,je ; do i=is,ie
Khth_loc_v(i,J) = Khth_loc_v(i,J) + G%OBCmaskCv(i,J) * CS%MEKE_GEOMETRIC_alpha * &
0.5*(MEKE%MEKE(i,j)+MEKE%MEKE(i,j+1)) / &
(VarMix%SN_v(i,J) + CS%MEKE_GEOMETRIC_epsilon)
enddo ; enddo
else
do J=js-1,je ; do i=is,ie
Khth_loc_v(i,J) = Khth_loc_v(i,J) + MEKE%KhTh_fac*sqrt(MEKE%Kh(i,j)*MEKE%Kh(i,j+1))
enddo ; enddo
endif
endif
if (Resoln_scaled) then
!$OMP do
do J=js-1,je ; do i=is,ie
Khth_loc_v(i,J) = Khth_loc_v(i,J) * VarMix%Res_fn_v(i,J)
enddo ; enddo
endif
if (Depth_scaled) then
!$OMP do
do J=js-1,je ; do i=is,ie
Khth_loc_v(i,J) = Khth_loc_v(i,J) * VarMix%Depth_fn_v(i,J)
enddo ; enddo
endif
if (CS%Khth_Max > 0) then
!$OMP do
do J=js-1,je ; do i=is,ie
Khth_loc_v(i,J) = max(CS%Khth_Min, min(Khth_loc_v(i,J), CS%Khth_Max))
enddo ; enddo
else
!$OMP do
do J=js-1,je ; do i=is,ie
Khth_loc_v(i,J) = max(CS%Khth_Min, Khth_loc_v(i,J))
enddo ; enddo
endif
if (CS%max_Khth_CFL > 0.0) then
!$OMP do
do J=js-1,je ; do i=is,ie
KH_v(i,J,1) = min(KH_v_CFL(i,J), Khth_loc_v(i,J))
enddo ; enddo
endif
if (khth_use_ebt_struct) then
!$OMP do
do K=2,nz+1 ; do J=js-1,je ; do i=is,ie
KH_v(i,J,K) = KH_v(i,J,1) * 0.5 * ( VarMix%ebt_struct(i,j,k-1) + VarMix%ebt_struct(i,j+1,k-1) )
enddo ; enddo ; enddo
else
!$OMP do
do K=2,nz+1 ; do J=js-1,je ; do i=is,ie
KH_v(i,J,K) = KH_v(i,J,1)
enddo ; enddo ; enddo
endif
if (use_VarMix) then
if (use_QG_Leith) then
!$OMP do
do k=1,nz ; do J=js-1,je ; do i=is,ie
KH_v(i,J,k) = VarMix%KH_v_QG(i,J,k)
enddo ; enddo ; enddo
endif
endif
if (CS%use_GME_thickness_diffuse) then
!$OMP do
do k=1,nz+1 ; do J=js-1,je ; do i=is,ie
CS%KH_v_GME(i,J,k) = KH_v(i,J,k)
enddo ; enddo ; enddo
endif
if (allocated(MEKE%Kh)) then
if (CS%MEKE_GEOMETRIC) then
if (CS%MEKE_GEOM_answer_date < 20190101) then
!$OMP do
do j=js,je ; do I=is,ie
! This does not give bitwise rotational symmetry.
MEKE%Kh(i,j) = CS%MEKE_GEOMETRIC_alpha * MEKE%MEKE(i,j) / &
(0.25*(VarMix%SN_u(I,j)+VarMix%SN_u(I-1,j) + &
VarMix%SN_v(i,J)+VarMix%SN_v(i,J-1)) + &
CS%MEKE_GEOMETRIC_epsilon)
enddo ; enddo
else
!$OMP do
do j=js,je ; do I=is,ie
! With the additional parentheses this gives bitwise rotational symmetry.
MEKE%Kh(i,j) = CS%MEKE_GEOMETRIC_alpha * MEKE%MEKE(i,j) / &
(0.25*((VarMix%SN_u(I,j)+VarMix%SN_u(I-1,j)) + &
(VarMix%SN_v(i,J)+VarMix%SN_v(i,J-1))) + &
CS%MEKE_GEOMETRIC_epsilon)
enddo ; enddo
endif
endif
endif
!$OMP do
do K=1,nz+1 ; do j=js,je ; do I=is-1,ie ; int_slope_u(I,j,K) = 0.0 ; enddo ; enddo ; enddo
!$OMP do
do K=1,nz+1 ; do J=js-1,je ; do i=is,ie ; int_slope_v(i,J,K) = 0.0 ; enddo ; enddo ; enddo
!$OMP end parallel
if (CS%detangle_interfaces) then
call add_detangling_Kh(h, e, Kh_u, Kh_v, KH_u_CFL, KH_v_CFL, tv, dt, G, GV, US, &
CS, int_slope_u, int_slope_v)
endif
if ((CS%Kh_eta_bg > 0.0) .or. (CS%Kh_eta_vel > 0.0)) then
call add_interface_Kh(G, GV, US, CS, Kh_u, Kh_v, KH_u_CFL, KH_v_CFL, int_slope_u, int_slope_v)
endif
if (CS%debug) then
call uvchksum("Kh_[uv]", Kh_u, Kh_v, G%HI, haloshift=0, &
scale=(US%L_to_m**2)*US%s_to_T, scalar_pair=.true.)
call uvchksum("Kh_[uv]_CFL", Kh_u_CFL, Kh_v_CFL, G%HI, haloshift=0, &
scale=(US%L_to_m**2)*US%s_to_T, scalar_pair=.true.)
if (Resoln_scaled) then
call uvchksum("Res_fn_[uv]", VarMix%Res_fn_u, VarMix%Res_fn_v, G%HI, haloshift=0, &
scale=1.0, scalar_pair=.true.)
endif
call uvchksum("int_slope_[uv]", int_slope_u, int_slope_v, G%HI, haloshift=0)
call hchksum(h, "thickness_diffuse_1 h", G%HI, haloshift=1, scale=GV%H_to_m)
call hchksum(e, "thickness_diffuse_1 e", G%HI, haloshift=1, scale=US%Z_to_m)
if (use_stored_slopes) then
call uvchksum("VarMix%slope_[xy]", VarMix%slope_x, VarMix%slope_y, &
G%HI, haloshift=0, scale=US%Z_to_L)
endif
if (associated(tv%eqn_of_state)) then
call hchksum(tv%T, "thickness_diffuse T", G%HI, haloshift=1, scale=US%C_to_degC)
call hchksum(tv%S, "thickness_diffuse S", G%HI, haloshift=1, scale=US%S_to_ppt)
endif
endif
! Calculate uhD, vhD from h, e, KH_u, KH_v, tv%T/S
if (use_stored_slopes) then
call thickness_diffuse_full(h, e, Kh_u, Kh_v, tv, uhD, vhD, cg1, dt, G, GV, US, MEKE, CS, &
int_slope_u, int_slope_v, VarMix%slope_x, VarMix%slope_y)
else
call thickness_diffuse_full(h, e, Kh_u, Kh_v, tv, uhD, vhD, cg1, dt, G, GV, US, MEKE, CS, &
int_slope_u, int_slope_v)
endif
if (VarMix%use_variable_mixing) then
if (allocated(MEKE%Rd_dx_h) .and. allocated(VarMix%Rd_dx_h)) then
!$OMP parallel do default(shared)
do j=js,je ; do i=is,ie
MEKE%Rd_dx_h(i,j) = VarMix%Rd_dx_h(i,j)
enddo ; enddo
endif
endif
! offer diagnostic fields for averaging
if (query_averaging_enabled(CS%diag)) then
if (CS%id_uhGM > 0) call post_data(CS%id_uhGM, uhD, CS%diag)
if (CS%id_vhGM > 0) call post_data(CS%id_vhGM, vhD, CS%diag)
if (CS%id_GMwork > 0) call post_data(CS%id_GMwork, CS%GMwork, CS%diag)
if (CS%id_KH_u > 0) call post_data(CS%id_KH_u, KH_u, CS%diag)
if (CS%id_KH_v > 0) call post_data(CS%id_KH_v, KH_v, CS%diag)
if (CS%id_KH_u1 > 0) call post_data(CS%id_KH_u1, KH_u(:,:,1), CS%diag)
if (CS%id_KH_v1 > 0) call post_data(CS%id_KH_v1, KH_v(:,:,1), CS%diag)
! Diagnose diffusivity at T-cell point. Do a simple average, rather than a
! thickness-weighted average, so that KH_t is depth-independent when KH_u and KH_v
! are depth independent. If a thickness-weighted average were used, the variations
! of thickness could give a spurious depth dependence to the diagnosed KH_t.
if (CS%id_KH_t > 0 .or. CS%id_KH_t1 > 0 .or. CS%Use_KH_in_MEKE) then
do k=1,nz
! thicknesses across u and v faces, converted to 0/1 mask
! layer average of the interface diffusivities KH_u and KH_v
do j=js,je ; do I=is-1,ie
! This expression uses harmonic mean thicknesses:
! hu(I,j) = 2.0*h(i,j,k)*h(i+1,j,k) / (h(i,j,k)+h(i+1,j,k)+h_neglect)
! This expression is a 0/1 mask based on depths where there are thick layers:
hu(I,j) = 0.0 ; if (h(i,j,k)*h(i+1,j,k) /= 0.0) hu(I,j) = 1.0
KH_u_lay(I,j) = 0.5*(KH_u(I,j,k)+KH_u(I,j,k+1))
enddo ; enddo
do J=js-1,je ; do i=is,ie
! This expression uses harmonic mean thicknesses:
! hv(i,J) = 2.0*h(i,j,k)*h(i,j+1,k)/(h(i,j,k)+h(i,j+1,k)+h_neglect)
! This expression is a 0/1 mask based on depths where there are thick layers:
hv(i,J) = 0.0 ; if (h(i,j,k)*h(i,j+1,k) /= 0.0) hv(i,J) = 1.0
KH_v_lay(i,J) = 0.5*(KH_v(i,J,k)+KH_v(i,J,k+1))
enddo ; enddo
! diagnose diffusivity at T-points
do j=js,je ; do i=is,ie
Kh_t(i,j,k) = ((hu(I-1,j)*KH_u_lay(i-1,j) + hu(I,j)*KH_u_lay(I,j)) + &
(hv(i,J-1)*KH_v_lay(i,J-1) + hv(i,J)*KH_v_lay(i,J))) / &
((hu(I-1,j)+hu(I,j)) + (hv(i,J-1)+hv(i,J)) + 1.0e-20)
! Use this denominator instead if hu and hv are actual thicknesses rather than a 0/1 mask:
! ((hu(I-1,j)+hu(I,j)) + (hv(i,J-1)+hv(i,J)) + h_neglect)
enddo ; enddo
enddo
if (CS%Use_KH_in_MEKE) then
MEKE%Kh_diff(:,:) = 0.0
htot(:,:) = 0.0
do k=1,nz
do j=js,je ; do i=is,ie
MEKE%Kh_diff(i,j) = MEKE%Kh_diff(i,j) + Kh_t(i,j,k) * h(i,j,k)
htot(i,j) = htot(i,j) + h(i,j,k)
enddo ; enddo
enddo
do j=js,je ; do i=is,ie
MEKE%Kh_diff(i,j) = MEKE%Kh_diff(i,j) / MAX(CS%MEKE_min_depth_diff, htot(i,j))
enddo ; enddo
endif
if (CS%id_KH_t > 0) call post_data(CS%id_KH_t, KH_t, CS%diag)
if (CS%id_KH_t1 > 0) call post_data(CS%id_KH_t1, KH_t(:,:,1), CS%diag)
endif
endif
!$OMP parallel do default(shared)
do k=1,nz
do j=js,je ; do I=is-1,ie
uhtr(I,j,k) = uhtr(I,j,k) + uhD(I,j,k) * dt
if (associated(CDp%uhGM)) CDp%uhGM(I,j,k) = uhD(I,j,k)
enddo ; enddo
do J=js-1,je ; do i=is,ie
vhtr(i,J,k) = vhtr(i,J,k) + vhD(i,J,k) * dt
if (associated(CDp%vhGM)) CDp%vhGM(i,J,k) = vhD(i,J,k)
enddo ; enddo
do j=js,je ; do i=is,ie
h(i,j,k) = h(i,j,k) - dt * G%IareaT(i,j) * &
((uhD(I,j,k) - uhD(I-1,j,k)) + (vhD(i,J,k) - vhD(i,J-1,k)))
if (h(i,j,k) < GV%Angstrom_H) h(i,j,k) = GV%Angstrom_H
enddo ; enddo
enddo
! Whenever thickness changes let the diag manager know, target grids
! for vertical remapping may need to be regenerated.
! This needs to happen after the H update and before the next post_data.
call diag_update_remap_grids(CS%diag)
if (CS%debug) then
call uvchksum("thickness_diffuse [uv]hD", uhD, vhD, &
G%HI, haloshift=0, scale=GV%H_to_m*US%L_to_m**2*US%s_to_T)
call uvchksum("thickness_diffuse [uv]htr", uhtr, vhtr, &
G%HI, haloshift=0, scale=US%L_to_m**2*GV%H_to_m)
call hchksum(h, "thickness_diffuse h", G%HI, haloshift=0, scale=GV%H_to_m)
endif
end subroutine thickness_diffuse
!> Calculates parameterized layer transports for use in the continuity equation.
!! Fluxes are limited to give positive definite thicknesses.
!! Called by thickness_diffuse().
subroutine thickness_diffuse_full(h, e, Kh_u, Kh_v, tv, uhD, vhD, cg1, dt, G, GV, US, MEKE, &
CS, int_slope_u, int_slope_v, slope_x, slope_y)
type(ocean_grid_type), intent(in) :: G !< Ocean grid structure
type(verticalGrid_type), intent(in) :: GV !< Vertical grid structure
type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), intent(in) :: h !< Layer thickness [H ~> m or kg m-2]
real, dimension(SZI_(G),SZJ_(G),SZK_(GV)+1), intent(in) :: e !< Interface positions [Z ~> m]
real, dimension(SZIB_(G),SZJ_(G),SZK_(GV)+1), intent(in) :: Kh_u !< Isopycnal height diffusivity
!! at u points [L2 T-1 ~> m2 s-1]
real, dimension(SZI_(G),SZJB_(G),SZK_(GV)+1), intent(in) :: Kh_v !< Isopycnal height diffusivity
!! at v points [L2 T-1 ~> m2 s-1]
type(thermo_var_ptrs), intent(in) :: tv !< Thermodynamics structure
real, dimension(SZIB_(G),SZJ_(G),SZK_(GV)), intent(out) :: uhD !< Zonal mass fluxes
!! [H L2 T-1 ~> m3 s-1 or kg s-1]
real, dimension(SZI_(G),SZJB_(G),SZK_(GV)), intent(out) :: vhD !< Meridional mass fluxes
!! [H L2 T-1 ~> m3 s-1 or kg s-1]
real, dimension(:,:), pointer :: cg1 !< Wave speed [L T-1 ~> m s-1]
real, intent(in) :: dt !< Time increment [T ~> s]
type(MEKE_type), intent(inout) :: MEKE !< MEKE fields
type(thickness_diffuse_CS), intent(inout) :: CS !< Control structure for thickness_diffuse
real, dimension(SZIB_(G),SZJ_(G),SZK_(GV)+1), intent(in) :: int_slope_u !< Ratio that determine how much of
!! the isopycnal slopes are taken directly from
!! the interface slopes without consideration of
!! density gradients [nondim].
real, dimension(SZI_(G),SZJB_(G),SZK_(GV)+1), intent(in) :: int_slope_v !< Ratio that determine how much of
!! the isopycnal slopes are taken directly from
!! the interface slopes without consideration of
!! density gradients [nondim].
real, dimension(SZIB_(G),SZJ_(G),SZK_(GV)+1), optional, intent(in) :: slope_x !< Isopyc. slope at u [Z L-1 ~> nondim]
real, dimension(SZI_(G),SZJB_(G),SZK_(GV)+1), optional, intent(in) :: slope_y !< Isopyc. slope at v [Z L-1 ~> nondim]
! Local variables
real, dimension(SZI_(G),SZJ_(G),SZK_(GV)) :: &
T, & ! The temperature [C ~> degC], with the values in
! in massless layers filled vertically by diffusion.
S, & ! The filled salinity [S ~> ppt], with the values in
! in massless layers filled vertically by diffusion.
h_avail, & ! The mass available for diffusion out of each face, divided
! by dt [H L2 T-1 ~> m3 s-1 or kg s-1].
h_frac ! The fraction of the mass in the column above the bottom
! interface of a layer that is within a layer [nondim]. 0<h_frac<=1
real :: dz(SZI_(G),SZJ_(G),SZK_(GV)) ! Height change across layers [Z ~> m]
real, dimension(SZI_(G),SZJB_(G),SZK_(GV)+1) :: &
Slope_y_PE, & ! 3D array of neutral slopes at v-points, set equal to Slope (below) [Z L-1 ~> nondim]
hN2_y_PE ! Harmonic mean of thicknesses around the interfaces times the buoyancy frequency
! at v-points with unit conversion factors [H L2 Z-2 T-2 ~> m s-2 or kg m-2 s-2],
! used for calculating the potential energy release
real, dimension(SZIB_(G),SZJ_(G),SZK_(GV)+1) :: &
Slope_x_PE, & ! 3D array of neutral slopes at u-points, set equal to Slope (below) [Z L-1 ~> nondim]
hN2_x_PE ! Harmonic mean of thicknesses around the interfaces times the buoyancy frequency
! at u-points with unit conversion factors [H L2 Z-2 T-2 ~> m s-2 or kg m-2 s-2],
! used for calculating the potential energy release
real, dimension(SZI_(G),SZJ_(G),SZK_(GV)+1) :: &
pres, & ! The pressure at an interface [R L2 T-2 ~> Pa].
h_avail_rsum ! The running sum of h_avail above an interface [H L2 T-1 ~> m3 s-1 or kg s-1].
real, dimension(SZIB_(G)) :: &
drho_dT_u, & ! The derivative of density with temperature at u points [R C-1 ~> kg m-3 degC-1]
drho_dS_u ! The derivative of density with salinity at u points [R S-1 ~> kg m-3 ppt-1].
real, dimension(SZIB_(G)) :: scrap ! An array to pass to calculate_density_second_derivs()
! with various units that will be ignored [various]
real, dimension(SZI_(G)) :: &
drho_dT_v, & ! The derivative of density with temperature at v points [R C-1 ~> kg m-3 degC-1]
drho_dS_v, & ! The derivative of density with salinity at v points [R S-1 ~> kg m-3 ppt-1].
drho_dT_dT_h, & ! The second derivative of density with temperature at h points [R C-2 ~> kg m-3 degC-2]
drho_dT_dT_hr ! The second derivative of density with temperature at h (+1) points [R C-2 ~> kg m-3 degC-2]
real :: uhtot(SZIB_(G),SZJ_(G)) ! The vertical sum of uhD [H L2 T-1 ~> m3 s-1 or kg s-1].
real :: vhtot(SZI_(G),SZJB_(G)) ! The vertical sum of vhD [H L2 T-1 ~> m3 s-1 or kg s-1].
real, dimension(SZIB_(G)) :: &
T_u, & ! Temperature on the interface at the u-point [C ~> degC].
S_u, & ! Salinity on the interface at the u-point [S ~> ppt].
pres_u ! Pressure on the interface at the u-point [R L2 T-2 ~> Pa].
real, dimension(SZI_(G)) :: &
T_v, & ! Temperature on the interface at the v-point [C ~> degC].
S_v, & ! Salinity on the interface at the v-point [S ~> ppt].
pres_v, & ! Pressure on the interface at the v-point [R L2 T-2 ~> Pa].
T_h, & ! Temperature on the interface at the h-point [C ~> degC].
S_h, & ! Salinity on the interface at the h-point [S ~> ppt].
pres_h, & ! Pressure on the interface at the h-point [R L2 T-2 ~> Pa].
T_hr, & ! Temperature on the interface at the h (+1) point [C ~> degC].
S_hr, & ! Salinity on the interface at the h (+1) point [S ~> ppt].
pres_hr ! Pressure on the interface at the h (+1) point [R L2 T-2 ~> Pa].
real :: Work_u(SZIB_(G),SZJ_(G)) ! The work done by the isopycnal height diffusion
! integrated over u-point water columns [R Z L4 T-3 ~> W]
real :: Work_v(SZI_(G),SZJB_(G)) ! The work done by the isopycnal height diffusion
! integrated over v-point water columns [R Z L4 T-3 ~> W]
real :: Work_h ! The work averaged over an h-cell [R Z L2 T-3 ~> W m-2].
real :: PE_release_h ! The amount of potential energy released by GM averaged over an h-cell
! [R Z L2 T-3 ~> W m-2]. The calculation equals rho0 * h * S^2 * N^2 * kappa_GM.
real :: I4dt ! 1 / 4 dt [T-1 ~> s-1].
real :: drdiA, drdiB ! Along layer zonal potential density gradients in the layers above (A)
! and below (B) the interface times the grid spacing [R ~> kg m-3].
real :: drdjA, drdjB ! Along layer meridional potential density gradients in the layers above (A)
! and below (B) the interface times the grid spacing [R ~> kg m-3].
real :: drdkL, drdkR ! Vertical density differences across an interface [R ~> kg m-3].
real :: drdi_u(SZIB_(G),SZK_(GV)) ! Copy of drdi at u-points [R ~> kg m-3].
real :: drdj_v(SZI_(G),SZK_(GV)) ! Copy of drdj at v-points [R ~> kg m-3].
real :: drdkDe_u(SZIB_(G),SZK_(GV)+1) ! Lateral difference of product of drdk and e at u-points
! [Z R ~> kg m-2].
real :: drdkDe_v(SZI_(G),SZK_(GV)+1) ! Lateral difference of product of drdk and e at v-points
! [Z R ~> kg m-2].
real :: hg2A, hg2B, hg2L, hg2R ! Squares of geometric mean thicknesses [H2 ~> m2 or kg2 m-4].
real :: haA, haB, haL, haR ! Arithmetic mean thicknesses [H ~> m or kg m-2].
real :: dzg2A, dzg2B ! Squares of geometric mean vertical layer extents [Z2 ~> m2].
real :: dzaA, dzaB ! Arithmetic mean vertical layer extents [Z ~> m].
real :: dzaL, dzaR ! Temporary vertical layer extents [Z ~> m]
real :: wtA, wtB ! Unnormalized weights of the slopes above and below [H3 ~> m3 or kg3 m-6]
real :: wtL, wtR ! Unnormalized weights of the slopes to the left and right [H3 Z ~> m4 or kg3 m-5]
real :: drdx, drdy ! Zonal and meridional density gradients [R L-1 ~> kg m-4].
real :: drdz ! Vertical density gradient [R Z-1 ~> kg m-4].
real :: dz_harm ! Harmonic mean layer vertical extent [Z ~> m].
real :: c2_dz_u(SZIB_(G),SZK_(GV)+1) ! Wave speed squared divided by dz at u-points times rescaling
! factors from depths to thicknesses [H2 L2 Z-3 T-2 ~> m s-2 or kg m-2 s-2]
real :: c2_dz_v(SZI_(G),SZK_(GV)+1) ! Wave speed squared divided by dz at v-points times rescaling
! factors from depths to thicknesses [H L2 Z-2 T-2 ~> m s-2 or kg m-2 s-2]
real :: dzN2_u(SZIB_(G),SZK_(GV)+1) ! Vertical extent times N2 at interfaces above u-points times
! rescaling factors from vertical to horizontal distances [L2 Z-1 T-2 ~> m s-2]
real :: dzN2_v(SZI_(G),SZK_(GV)+1) ! Vertical extent times N2 at interfaces above v-points times
! rescaling factors from vertical to horizontal distances [L2 Z-1 T-2 ~> m s-2]
real :: Sfn_est ! A preliminary estimate (before limiting) of the overturning
! streamfunction [H L2 T-1 ~> m3 s-1 or kg s-1].
real :: Sfn_unlim_u(SZIB_(G),SZK_(GV)+1) ! Volume streamfunction for u-points [Z L2 T-1 ~> m3 s-1]
real :: Sfn_unlim_v(SZI_(G),SZK_(GV)+1) ! Volume streamfunction for v-points [Z L2 T-1 ~> m3 s-1]
real :: slope2_Ratio_u(SZIB_(G),SZK_(GV)+1) ! The ratio of the slope squared to slope_max squared [nondim]
real :: slope2_Ratio_v(SZI_(G),SZK_(GV)+1) ! The ratio of the slope squared to slope_max squared [nondim]
real :: Sfn_in_h ! The overturning streamfunction [H L2 T-1 ~> m3 s-1 or kg s-1] (note that
! the units are different from other Sfn vars).
real :: Sfn_safe ! The streamfunction that goes linearly back to 0 at the surface
! [H L2 T-1 ~> m3 s-1 or kg s-1]. This is a good value to use when the
! slope is so large as to be meaningless, usually due to weak stratification.
real :: Slope ! The slope of density surfaces, calculated in a way that is always
! between -1 and 1 after undoing dimensional scaling, [Z L-1 ~> nondim]
real :: mag_grad2 ! The squared magnitude of the 3-d density gradient [R2 L-2 ~> kg2 m-8].
real :: I_slope_max2 ! The inverse of slope_max squared [L2 Z-2 ~> nondim].
real :: h_neglect ! A thickness that is so small it is usually lost
! in roundoff and can be neglected [H ~> m or kg m-2].
real :: hn_2 ! Half of h_neglect [H ~> m or kg m-2].
real :: h_neglect2 ! h_neglect^2 [H2 ~> m2 or kg2 m-4].
real :: dz_neglect ! A thickness [Z ~> m], that is so small it is usually lost
! in roundoff and can be neglected [Z ~> m].
real :: dz_neglect2 ! dz_neglect^2 [Z2 ~> m2]
real :: G_scale ! The gravitational acceleration times a unit conversion
! factor [L2 H-1 T-2 ~> m s-2 or m4 kg-1 s-2].
logical :: use_EOS ! If true, density is calculated from T & S using an equation of state.
logical :: find_work ! If true, find the change in energy due to the fluxes.
integer :: nk_linear ! The number of layers over which the streamfunction goes to 0.
real :: G_rho0 ! g/Rho0 [L2 R-1 Z-1 T-2 ~> m4 kg-1 s-2].
real :: Rho_avg ! The in situ density averaged to an interface [R ~> kg m-3]
real :: N2_floor ! A floor for N2 to avoid degeneracy in the elliptic solver
! times unit conversion factors [L2 Z-2 T-2 ~> s-2]
real :: N2_unlim ! An unlimited estimate of the buoyancy frequency
! times unit conversion factors [L2 Z-2 T-2 ~> s-2]
real :: Tl(5) ! copy of T in local stencil [C ~> degC]
real :: mn_T ! mean of T in local stencil [C ~> degC]
real :: mn_T2 ! mean of T**2 in local stencil [C2 ~> degC2]
real :: hl(5) ! Copy of local stencil of H [H ~> m]
real :: r_sm_H ! Reciprocal of sum of H in local stencil [H-1 ~> m-1]
real :: Z_to_H ! A conversion factor from heights to thicknesses, perhaps based on
! a spatially variable local density [H Z-1 ~> nondim or kg m-3]
real :: Tsgs2(SZI_(G),SZJ_(G),SZK_(GV)) ! Sub-grid temperature variance [C2 ~> degC2]
real :: diag_sfn_x(SZIB_(G),SZJ_(G),SZK_(GV)+1) ! Diagnostic of the x-face streamfunction
! [H L2 T-1 ~> m3 s-1 or kg s-1]
real :: diag_sfn_unlim_x(SZIB_(G),SZJ_(G),SZK_(GV)+1) ! Diagnostic of the x-face streamfunction before
! applying limiters [Z L2 T-1 ~> m3 s-1]
real :: diag_sfn_y(SZI_(G),SZJB_(G),SZK_(GV)+1) ! Diagnostic of the y-face streamfunction
! [H L2 T-1 ~> m3 s-1 or kg s-1]
real :: diag_sfn_unlim_y(SZI_(G),SZJB_(G),SZK_(GV)+1) ! Diagnostic of the y-face streamfunction before
! applying limiters [Z L2 T-1 ~> m3 s-1]
logical :: present_slope_x, present_slope_y, calc_derivatives
integer, dimension(2) :: EOSdom_u ! The shifted I-computational domain to use for equation of
! state calculations at u-points.
integer, dimension(2) :: EOSdom_v ! The shifted i-computational domain to use for equation of
! state calculations at v-points.
integer, dimension(2) :: EOSdom_h1 ! The shifted i-computational domain to use for equation of
! state calculations at h points with 1 extra halo point
logical :: use_stanley
integer :: is, ie, js, je, nz, IsdB, halo
integer :: i, j, k
is = G%isc ; ie = G%iec ; js = G%jsc ; je = G%jec ; nz = GV%ke ; IsdB = G%IsdB
I4dt = 0.25 / dt
I_slope_max2 = 1.0 / (CS%slope_max**2)
h_neglect = GV%H_subroundoff ; h_neglect2 = h_neglect**2 ; hn_2 = 0.5*h_neglect
dz_neglect = GV%dZ_subroundoff ; dz_neglect2 = dz_neglect**2
if (GV%Boussinesq) G_rho0 = GV%g_Earth / GV%Rho0
N2_floor = CS%N2_floor * US%Z_to_L**2
use_EOS = associated(tv%eqn_of_state)
present_slope_x = PRESENT(slope_x)
present_slope_y = PRESENT(slope_y)
use_stanley = CS%use_stanley_gm
nk_linear = max(GV%nkml, 1)
Slope_x_PE(:,:,:) = 0.0
Slope_y_PE(:,:,:) = 0.0
hN2_x_PE(:,:,:) = 0.0
hN2_y_PE(:,:,:) = 0.0
find_work = allocated(MEKE%GM_src)
find_work = (allocated(CS%GMwork) .or. find_work)
if (use_EOS) then
halo = 1 ! Default halo to fill is 1
call vert_fill_TS(h, tv%T, tv%S, CS%kappa_smooth*dt, T, S, G, GV, US, halo, larger_h_denom=.true.)
endif
! Rescale the thicknesses, perhaps using the specific volume.
call thickness_to_dz(h, tv, dz, G, GV, US, halo_size=1)
if (CS%use_FGNV_streamfn .and. .not. associated(cg1)) call MOM_error(FATAL, &
"cg1 must be associated when using FGNV streamfunction.")
!$OMP parallel default(shared) private(hl,r_sm_H,Tl,mn_T,mn_T2)
! Find the maximum and minimum permitted streamfunction.
!$OMP do
do j=js-1,je+1 ; do i=is-1,ie+1
h_avail_rsum(i,j,1) = 0.0
pres(i,j,1) = 0.0
if (associated(tv%p_surf)) then ; pres(i,j,1) = tv%p_surf(i,j) ; endif
h_avail(i,j,1) = max(I4dt*G%areaT(i,j)*(h(i,j,1)-GV%Angstrom_H),0.0)
h_avail_rsum(i,j,2) = h_avail(i,j,1)
h_frac(i,j,1) = 1.0
pres(i,j,2) = pres(i,j,1) + (GV%g_Earth*GV%H_to_RZ) * h(i,j,1)
enddo ; enddo
do j=js-1,je+1
do k=2,nz ; do i=is-1,ie+1
h_avail(i,j,k) = max(I4dt*G%areaT(i,j)*(h(i,j,k)-GV%Angstrom_H),0.0)
h_avail_rsum(i,j,k+1) = h_avail_rsum(i,j,k) + h_avail(i,j,k)
h_frac(i,j,k) = 0.0 ; if (h_avail(i,j,k) > 0.0) &
h_frac(i,j,k) = h_avail(i,j,k) / h_avail_rsum(i,j,k+1)
pres(i,j,K+1) = pres(i,j,K) + (GV%g_Earth*GV%H_to_RZ) * h(i,j,k)
enddo ; enddo
enddo
!$OMP do
do j=js,je ; do I=is-1,ie
uhtot(I,j) = 0.0 ; Work_u(I,j) = 0.0
enddo ; enddo
!$OMP do
do J=js-1,je ; do i=is,ie
vhtot(i,J) = 0.0 ; Work_v(i,J) = 0.0
enddo ; enddo
!$OMP end parallel
if (CS%id_sfn_x > 0) then ; diag_sfn_x(:,:,1) = 0.0 ; diag_sfn_x(:,:,nz+1) = 0.0 ; endif
if (CS%id_sfn_y > 0) then ; diag_sfn_y(:,:,1) = 0.0 ; diag_sfn_y(:,:,nz+1) = 0.0 ; endif
if (CS%id_sfn_unlim_x > 0) then ; diag_sfn_unlim_x(:,:,1) = 0.0 ; diag_sfn_unlim_x(:,:,nz+1) = 0.0 ; endif
if (CS%id_sfn_unlim_y > 0) then ; diag_sfn_unlim_y(:,:,1) = 0.0 ; diag_sfn_unlim_y(:,:,nz+1) = 0.0 ; endif
EOSdom_u(1) = (is-1) - (G%IsdB-1) ; EOSdom_u(2) = ie - (G%IsdB-1)
EOSdom_v(:) = EOS_domain(G%HI)
EOSdom_h1(:) = EOS_domain(G%HI, halo=1)
!$OMP parallel do default(none) shared(nz,is,ie,js,je,find_work,use_EOS,G,GV,US,pres,T,S, &
!$OMP nk_linear,IsdB,tv,h,h_neglect,e,dz,dz_neglect,dz_neglect2, &
!$OMP h_neglect2,hn_2,I_slope_max2,int_slope_u,KH_u,uhtot, &
!$OMP h_frac,h_avail_rsum,uhD,h_avail,Work_u,CS,slope_x,cg1, &
!$OMP diag_sfn_x,diag_sfn_unlim_x,N2_floor,EOSdom_u,EOSdom_h1, &
!$OMP use_stanley,Tsgs2,present_slope_x,G_rho0,Slope_x_PE,hN2_x_PE) &
!$OMP private(drdiA,drdiB,drdkL,drdkR,pres_u,T_u,S_u,G_scale, &
!$OMP drho_dT_u,drho_dS_u,hg2A,hg2B,hg2L,hg2R,haA, &
!$OMP drho_dT_dT_h,scrap,pres_h,T_h,S_h,N2_unlim, &
!$OMP haB,haL,haR,dzaL,dzaR,wtA,wtB,wtL,wtR,drdz, &
!$OMP dzg2A,dzg2B,dzaA,dzaB,dz_harm,Z_to_H, &
!$OMP drdx,mag_grad2,Slope,slope2_Ratio_u,dzN2_u, &
!$OMP Sfn_unlim_u,Rho_avg,drdi_u,drdkDe_u,c2_dz_u, &
!$OMP Sfn_safe,Sfn_est,Sfn_in_h,calc_derivatives)
do j=js,je
do I=is-1,ie ; dzN2_u(I,1) = 0. ; dzN2_u(I,nz+1) = 0. ; enddo
do K=nz,2,-1
if (find_work .and. .not.(use_EOS)) then
drdiA = 0.0 ; drdiB = 0.0
drdkL = GV%Rlay(k) - GV%Rlay(k-1) ; drdkR = drdkL
endif
calc_derivatives = use_EOS .and. (k >= nk_linear) .and. &
(find_work .or. .not. present_slope_x .or. CS%use_FGNV_streamfn .or. use_stanley)
! Calculate the zonal fluxes and gradients.
if (calc_derivatives) then
do I=is-1,ie
pres_u(I) = 0.5*(pres(i,j,K) + pres(i+1,j,K))
T_u(I) = 0.25*((T(i,j,k) + T(i+1,j,k)) + (T(i,j,k-1) + T(i+1,j,k-1)))
S_u(I) = 0.25*((S(i,j,k) + S(i+1,j,k)) + (S(i,j,k-1) + S(i+1,j,k-1)))
enddo
call calculate_density_derivs(T_u, S_u, pres_u, drho_dT_u, drho_dS_u, &
tv%eqn_of_state, EOSdom_u)
endif
if (use_stanley) then
do i=is-1,ie+1
pres_h(i) = pres(i,j,K)
T_h(i) = 0.5*(T(i,j,k) + T(i,j,k-1))
S_h(i) = 0.5*(S(i,j,k) + S(i,j,k-1))
enddo
! The second line below would correspond to arguments
! drho_dS_dS, drho_dS_dT, drho_dT_dT, drho_dS_dP, drho_dT_dP, &
call calculate_density_second_derivs(T_h, S_h, pres_h, &
scrap, scrap, drho_dT_dT_h, scrap, scrap, &
tv%eqn_of_state, EOSdom_h1)
endif
do I=is-1,ie
if (calc_derivatives) then
! Estimate the horizontal density gradients along layers.
drdiA = drho_dT_u(I) * (T(i+1,j,k-1)-T(i,j,k-1)) + &
drho_dS_u(I) * (S(i+1,j,k-1)-S(i,j,k-1))
drdiB = drho_dT_u(I) * (T(i+1,j,k)-T(i,j,k)) + &
drho_dS_u(I) * (S(i+1,j,k)-S(i,j,k))
! Estimate the vertical density gradients times the grid spacing.
drdkL = (drho_dT_u(I) * (T(i,j,k)-T(i,j,k-1)) + &
drho_dS_u(I) * (S(i,j,k)-S(i,j,k-1)))
drdkR = (drho_dT_u(I) * (T(i+1,j,k)-T(i+1,j,k-1)) + &
drho_dS_u(I) * (S(i+1,j,k)-S(i+1,j,k-1)))
drdkDe_u(I,K) = drdkR * e(i+1,j,K) - drdkL * e(i,j,K)
elseif (find_work) then ! This is used in pure stacked SW mode
drdkDe_u(I,K) = drdkR * e(i+1,j,K) - drdkL * e(i,j,K)
endif
if (use_stanley) then
! Correction to the horizontal density gradient due to nonlinearity in
! the EOS rectifying SGS temperature anomalies
drdiA = drdiA + 0.5 * ((drho_dT_dT_h(i+1) * tv%varT(i+1,j,k-1)) - &
(drho_dT_dT_h(i) * tv%varT(i,j,k-1)) )
drdiB = drdiB + 0.5 * ((drho_dT_dT_h(i+1) * tv%varT(i+1,j,k)) - &
(drho_dT_dT_h(i) * tv%varT(i,j,k)) )
endif
if (find_work) drdi_u(I,k) = drdiB
if (k > nk_linear) then
if (use_EOS) then
if (CS%use_FGNV_streamfn .or. find_work .or. .not.present_slope_x) then
hg2L = h(i,j,k-1)*h(i,j,k) + h_neglect2
hg2R = h(i+1,j,k-1)*h(i+1,j,k) + h_neglect2
haL = 0.5*(h(i,j,k-1) + h(i,j,k)) + h_neglect
haR = 0.5*(h(i+1,j,k-1) + h(i+1,j,k)) + h_neglect
if (GV%Boussinesq) then
dzaL = haL * GV%H_to_Z ; dzaR = haR * GV%H_to_Z
elseif (GV%semi_Boussinesq) then
dzaL = 0.5*(e(i,j,K-1) - e(i,j,K+1)) + dz_neglect
dzaR = 0.5*(e(i+1,j,K-1) - e(i+1,j,K+1)) + dz_neglect
else
dzaL = 0.5*(dz(i,j,k-1) + dz(i,j,k)) + dz_neglect
dzaR = 0.5*(dz(i+1,j,k-1) + dz(i+1,j,k)) + dz_neglect
endif
! Use the harmonic mean thicknesses to weight the horizontal gradients.
! These unnormalized weights have been rearranged to minimize divisions.
wtL = hg2L*(haR*dzaR) ; wtR = hg2R*(haL*dzaL)
drdz = (wtL * drdkL + wtR * drdkR) / (dzaL*wtL + dzaR*wtR)
! The expression for drdz above is mathematically equivalent to:
! drdz = ((hg2L/haL) * drdkL/dzaL + (hg2R/haR) * drdkR/dzaR) / &
! ((hg2L/haL) + (hg2R/haR))
hg2A = h(i,j,k-1)*h(i+1,j,k-1) + h_neglect2
hg2B = h(i,j,k)*h(i+1,j,k) + h_neglect2
haA = 0.5*(h(i,j,k-1) + h(i+1,j,k-1)) + h_neglect
haB = 0.5*(h(i,j,k) + h(i+1,j,k)) + h_neglect
if (GV%Boussinesq) then
N2_unlim = drdz*G_rho0
else
N2_unlim = (GV%g_Earth*GV%RZ_to_H) * &
((wtL * drdkL + wtR * drdkR) / (haL*wtL + haR*wtR))
endif
dzg2A = dz(i,j,k-1)*dz(i+1,j,k-1) + dz_neglect2
dzg2B = dz(i,j,k)*dz(i+1,j,k) + dz_neglect2
dzaA = 0.5*(dz(i,j,k-1) + dz(i+1,j,k-1)) + dz_neglect
dzaB = 0.5*(dz(i,j,k) + dz(i+1,j,k)) + dz_neglect
! dzN2_u is used with the FGNV streamfunction formulation
dzN2_u(I,K) = (0.5 * ( dzg2A / dzaA + dzg2B / dzaB )) * max(N2_unlim, N2_floor)
if (find_work .and. CS%GM_src_alt) &
hN2_x_PE(I,j,k) = (0.5 * ( hg2A / haA + hg2B / haB )) * max(N2_unlim, N2_floor)
endif
if (present_slope_x) then
Slope = slope_x(I,j,k)
slope2_Ratio_u(I,K) = Slope**2 * I_slope_max2
else
! Use the harmonic mean thicknesses to weight the horizontal gradients.
! These unnormalized weights have been rearranged to minimize divisions.
wtA = hg2A*haB ; wtB = hg2B*haA
! This is the gradient of density along geopotentials.
drdx = ((wtA * drdiA + wtB * drdiB) / (wtA + wtB) - &
drdz * (e(i,j,K)-e(i+1,j,K))) * G%IdxCu(I,j)
! This estimate of slope is accurate for small slopes, but bounded
! to be between -1 and 1.
mag_grad2 = (US%Z_to_L*drdx)**2 + drdz**2
if (mag_grad2 > 0.0) then
Slope = drdx / sqrt(mag_grad2)
slope2_Ratio_u(I,K) = Slope**2 * I_slope_max2
else ! Just in case mag_grad2 = 0 ever.
Slope = 0.0
slope2_Ratio_u(I,K) = 1.0e20 ! Force the use of the safe streamfunction.
endif
endif