forked from mom-ocean/MOM6
-
Notifications
You must be signed in to change notification settings - Fork 60
/
MOM_MEKE.F90
1703 lines (1579 loc) · 86.3 KB
/
MOM_MEKE.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
!> Implements the Mesoscale Eddy Kinetic Energy framework
!! with topographic beta effect included in computing beta in Rhines scale
module MOM_MEKE
! This file is part of MOM6. See LICENSE.md for the license.
use MOM_debugging, only : hchksum, uvchksum
use MOM_cpu_clock, only : cpu_clock_id, cpu_clock_begin, cpu_clock_end, CLOCK_ROUTINE
use MOM_diag_mediator, only : post_data, register_diag_field, safe_alloc_ptr
use MOM_diag_mediator, only : diag_ctrl, time_type
use MOM_domains, only : create_group_pass, do_group_pass, group_pass_type
use MOM_error_handler, only : MOM_error, FATAL, WARNING, NOTE, MOM_mesg
use MOM_file_parser, only : read_param, get_param, log_version, param_file_type
use MOM_grid, only : ocean_grid_type
use MOM_hor_index, only : hor_index_type
use MOM_io, only : vardesc, var_desc
use MOM_restart, only : MOM_restart_CS, register_restart_field, query_initialized
use MOM_unit_scaling, only : unit_scale_type
use MOM_variables, only : vertvisc_type
use MOM_verticalGrid, only : verticalGrid_type
use MOM_MEKE_types, only : MEKE_type
implicit none ; private
#include <MOM_memory.h>
public step_forward_MEKE, MEKE_init, MEKE_alloc_register_restart, MEKE_end
!> Control structure that contains MEKE parameters and diagnostics handles
type, public :: MEKE_CS ; private
! Parameters
real, dimension(:,:), pointer :: equilibrium_value => NULL() !< The equilbrium value
!! of MEKE to be calculated at each time step [L2 T-2 ~> m2 s-2]
real :: MEKE_FrCoeff !< Efficiency of conversion of ME into MEKE [nondim]
real :: MEKE_GMcoeff !< Efficiency of conversion of PE into MEKE [nondim]
real :: MEKE_GMECoeff !< Efficiency of conversion of MEKE into ME by GME [nondim]
real :: MEKE_damping !< Local depth-independent MEKE dissipation rate [T-1 ~> s-1].
real :: MEKE_Cd_scale !< The ratio of the bottom eddy velocity to the column mean
!! eddy velocity, i.e. sqrt(2*MEKE). This should be less than 1
!! to account for the surface intensification of MEKE.
real :: MEKE_Cb !< Coefficient in the \f$\gamma_{bot}\f$ expression [nondim]
real :: MEKE_min_gamma!< Minimum value of gamma_b^2 allowed [nondim]
real :: MEKE_Ct !< Coefficient in the \f$\gamma_{bt}\f$ expression [nondim]
logical :: visc_drag !< If true use the vertvisc_type to calculate bottom drag.
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 thickness diffusion.
logical :: MEKE_equilibrium_alt !< If true, use an alternative calculation for the
!! equilibrium value of MEKE.
logical :: MEKE_equilibrium_restoring !< If true, restore MEKE back to its equilibrium value,
!! which is calculated at each time step.
logical :: GM_src_alt !< If true, use the GM energy conversion form S^2*N^2*kappa rather
!! than the streamfunction for the MEKE GM source term.
logical :: Rd_as_max_scale !< If true the length scale can not exceed the
!! first baroclinic deformation radius.
logical :: use_old_lscale !< Use the old formula for mixing length scale.
logical :: use_min_lscale !< Use simple minimum for mixing length scale.
real :: cdrag !< The bottom drag coefficient for MEKE [nondim].
real :: MEKE_BGsrc !< Background energy source for MEKE [L2 T-3 ~> W kg-1] (= m2 s-3).
real :: MEKE_dtScale !< Scale factor to accelerate time-stepping [nondim]
real :: MEKE_KhCoeff !< Scaling factor to convert MEKE into Kh [nondim]
real :: MEKE_Uscale !< MEKE velocity scale for bottom drag [L T-1 ~> m s-1]
real :: MEKE_KH !< Background lateral diffusion of MEKE [L2 T-1 ~> m2 s-1]
real :: MEKE_K4 !< Background bi-harmonic diffusivity (of MEKE) [L4 T-1 ~> m4 s-1]
real :: KhMEKE_Fac !< A factor relating MEKE%Kh to the diffusivity used for
!! MEKE itself [nondim].
real :: viscosity_coeff_Ku !< The scaling coefficient in the expression for
!! viscosity used to parameterize lateral harmonic momentum mixing
!! by unresolved eddies represented by MEKE.
real :: viscosity_coeff_Au !< The scaling coefficient in the expression for
!! viscosity used to parameterize lateral biharmonic momentum mixing
!! by unresolved eddies represented by MEKE.
real :: Lfixed !< Fixed mixing length scale [L ~> m].
real :: aDeform !< Weighting towards deformation scale of mixing length [nondim]
real :: aRhines !< Weighting towards Rhines scale of mixing length [nondim]
real :: aFrict !< Weighting towards frictional arrest scale of mixing length [nondim]
real :: aEady !< Weighting towards Eady scale of mixing length [nondim]
real :: aGrid !< Weighting towards grid scale of mixing length [nondim]
real :: MEKE_advection_factor !< A scaling in front of the advection of MEKE [nondim]
real :: MEKE_topographic_beta !< Weight for how much topographic beta is considered
!! when computing beta in Rhines scale [nondim]
real :: MEKE_restoring_rate !< Inverse of the timescale used to nudge MEKE toward its equilibrium value [s-1].
logical :: fixed_total_depth !< If true, use the nominal bathymetric depth as the estimate of
!! the time-varying ocean depth. Otherwise base the depth on the total
!! ocean mass per unit area.
logical :: kh_flux_enabled !< If true, lateral diffusive MEKE flux is enabled.
logical :: initialize !< If True, invokes a steady state solver to calculate MEKE.
logical :: debug !< If true, write out checksums of data for debugging
type(diag_ctrl), pointer :: diag => NULL() !< A type that regulates diagnostics output
!>@{ Diagnostic handles
integer :: id_MEKE = -1, id_Ue = -1, id_Kh = -1, id_src = -1
integer :: id_Ub = -1, id_Ut = -1
integer :: id_GM_src = -1, id_mom_src = -1, id_GME_snk = -1, id_decay = -1
integer :: id_KhMEKE_u = -1, id_KhMEKE_v = -1, id_Ku = -1, id_Au = -1
integer :: id_Le = -1, id_gamma_b = -1, id_gamma_t = -1
integer :: id_Lrhines = -1, id_Leady = -1
integer :: id_MEKE_equilibrium = -1
!>@}
! Infrastructure
integer :: id_clock_pass !< Clock for group pass calls
type(group_pass_type) :: pass_MEKE !< Group halo pass handle for MEKE%MEKE and maybe MEKE%Kh_diff
type(group_pass_type) :: pass_Kh !< Group halo pass handle for MEKE%Kh, MEKE%Ku, and/or MEKE%Au
end type MEKE_CS
contains
!> Integrates forward-in-time the MEKE eddy energy equation.
!! See \ref section_MEKE_equations.
subroutine step_forward_MEKE(MEKE, h, SN_u, SN_v, visc, dt, G, GV, US, CS, hu, hv)
type(MEKE_type), pointer :: MEKE !< MEKE data.
type(ocean_grid_type), intent(inout) :: G !< Ocean grid.
type(verticalGrid_type), intent(in) :: GV !< Ocean 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(SZIB_(G),SZJ_(G)), intent(in) :: SN_u !< Eady growth rate at u-points [T-1 ~> s-1].
real, dimension(SZI_(G),SZJB_(G)), intent(in) :: SN_v !< Eady growth rate at v-points [T-1 ~> s-1].
type(vertvisc_type), intent(in) :: visc !< The vertical viscosity type.
real, intent(in) :: dt !< Model(baroclinic) time-step [T ~> s].
type(MEKE_CS), pointer :: CS !< MEKE control structure.
real, dimension(SZIB_(G),SZJ_(G),SZK_(GV)), intent(in) :: hu !< Accumlated zonal mass flux [H L2 ~> m3 or kg]
real, dimension(SZI_(G),SZJB_(G),SZK_(GV)), intent(in) :: hv !< Accumlated meridional mass flux [H L2 ~> m3 or kg]
! Local variables
real, dimension(SZI_(G),SZJ_(G)) :: &
mass, & ! The total mass of the water column [R Z ~> kg m-2].
I_mass, & ! The inverse of mass [R-1 Z-1 ~> m2 kg-1].
depth_tot, & ! The depth of the water column [Z ~> m].
src, & ! The sum of all MEKE sources [L2 T-3 ~> W kg-1] (= m2 s-3).
MEKE_decay, & ! A diagnostic of the MEKE decay timescale [T-1 ~> s-1].
drag_rate_visc, & ! Near-bottom velocity contribution to bottom dratg [L T-1 ~> m s-1]
drag_rate, & ! The MEKE spindown timescale due to bottom drag [T-1 ~> s-1].
drag_rate_J15, & ! The MEKE spindown timescale due to bottom drag with the Jansen 2015 scheme.
! Unfortunately, as written the units seem inconsistent. [T-1 ~> s-1].
del2MEKE, & ! Laplacian of MEKE, used for bi-harmonic diffusion [T-2 ~> s-2].
del4MEKE, & ! Time-integrated MEKE tendency arising from the biharmonic of MEKE [L2 T-2 ~> m2 s-2].
LmixScale, & ! Eddy mixing length [L ~> m].
barotrFac2, & ! Ratio of EKE_barotropic / EKE [nondim]
bottomFac2, & ! Ratio of EKE_bottom / EKE [nondim]
tmp ! Temporary variable for diagnostic computation
real, dimension(SZIB_(G),SZJ_(G)) :: &
MEKE_uflux, & ! The zonal advective and diffusive flux of MEKE with units of [R Z L4 T-3 ~> kg m-2 s-3].
! In one place, MEKE_uflux is used as temporary work space with units of [L2 T-2 ~> m2 s-2].
Kh_u, & ! The zonal diffusivity that is actually used [L2 T-1 ~> m2 s-1].
baroHu, & ! Depth integrated accumulated zonal mass flux [R Z L2 ~> kg].
drag_vel_u ! A (vertical) viscosity associated with bottom drag at
! u-points [Z T-1 ~> m s-1].
real, dimension(SZI_(G),SZJB_(G)) :: &
MEKE_vflux, & ! The meridional advective and diffusive flux of MEKE with units of [R Z L4 T-3 ~> kg m-2 s-3].
! In one place, MEKE_vflux is used as temporary work space with units of [L2 T-2 ~> m2 s-2].
Kh_v, & ! The meridional diffusivity that is actually used [L2 T-1 ~> m2 s-1].
baroHv, & ! Depth integrated accumulated meridional mass flux [R Z L2 ~> kg].
drag_vel_v ! A (vertical) viscosity associated with bottom drag at
! v-points [Z T-1 ~> m s-1].
real :: Kh_here ! The local horizontal viscosity [L2 T-1 ~> m2 s-1]
real :: Inv_Kh_max ! The inverse of the local horizontal viscosity [T L-2 ~> s m-2]
real :: K4_here ! The local horizontal biharmonic viscosity [L4 T-1 ~> m4 s-1]
real :: Inv_K4_max ! The inverse of the local horizontal biharmonic viscosity [T L-4 ~> s m-4]
real :: cdrag2
real :: advFac ! The product of the advection scaling factor and 1/dt [T-1 ~> s-1]
real :: mass_neglect ! A negligible mass [R Z ~> kg m-2].
real :: ldamping ! The MEKE damping rate [T-1 ~> s-1].
real :: Rho0 ! A density used to convert mass to distance [R ~> kg m-3]
real :: I_Rho0 ! The inverse of the density used to convert mass to distance [R-1 ~> m3 kg-1]
real :: sdt ! dt to use locally [T ~> s] (could be scaled to accelerate)
real :: sdt_damp ! dt for damping [T ~> s] (sdt could be split).
logical :: use_drag_rate ! Flag to indicate drag_rate is finite
integer :: i, j, k, is, ie, js, je, Isq, Ieq, Jsq, Jeq, nz
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
if (.not.associated(CS)) call MOM_error(FATAL, &
"MOM_MEKE: Module must be initialized before it is used.")
if (.not.associated(MEKE)) call MOM_error(FATAL, &
"MOM_MEKE: MEKE must be initialized before it is used.")
if ((CS%MEKE_damping > 0.0) .or. (CS%MEKE_Cd_scale > 0.0) .or. (CS%MEKE_Cb>0.) &
.or. CS%visc_drag) then
use_drag_rate = .true.
else
use_drag_rate = .false.
endif
! Only integrate the MEKE equations if MEKE is required.
if (.not.associated(MEKE%MEKE)) then
! call MOM_error(FATAL, "MOM_MEKE: MEKE%MEKE is not associated!")
return
endif
if (CS%debug) then
if (associated(MEKE%mom_src)) &
call hchksum(MEKE%mom_src, 'MEKE mom_src', G%HI, scale=US%RZ3_T3_to_W_m2*US%L_to_Z**2)
if (associated(MEKE%GME_snk)) &
call hchksum(MEKE%GME_snk, 'MEKE GME_snk', G%HI, scale=US%RZ3_T3_to_W_m2*US%L_to_Z**2)
if (associated(MEKE%GM_src)) &
call hchksum(MEKE%GM_src, 'MEKE GM_src', G%HI, scale=US%RZ3_T3_to_W_m2*US%L_to_Z**2)
if (associated(MEKE%MEKE)) call hchksum(MEKE%MEKE, 'MEKE MEKE', G%HI, scale=US%L_T_to_m_s**2)
call uvchksum("MEKE SN_[uv]", SN_u, SN_v, G%HI, scale=US%s_to_T, &
scalar_pair=.true.)
call uvchksum("MEKE h[uv]", hu, hv, G%HI, haloshift=1, &
scale=GV%H_to_m*(US%L_to_m**2))
endif
sdt = dt*CS%MEKE_dtScale ! Scaled dt to use for time-stepping
Rho0 = GV%Rho0
I_Rho0 = 1.0 / GV%Rho0
mass_neglect = GV%H_to_RZ * GV%H_subroundoff
cdrag2 = CS%cdrag**2
! With a depth-dependent (and possibly strong) damping, it seems
! advisable to use Strang splitting between the damping and diffusion.
sdt_damp = sdt ; if (CS%MEKE_KH >= 0.0 .or. CS%MEKE_K4 >= 0.) sdt_damp = 0.5*sdt
! Calculate depth integrated mass exchange if doing advection [R Z L2 ~> kg]
if (CS%MEKE_advection_factor>0.) then
do j=js,je ; do I=is-1,ie
baroHu(I,j) = 0.
enddo ; enddo
do k=1,nz
do j=js,je ; do I=is-1,ie
baroHu(I,j) = hu(I,j,k) * GV%H_to_RZ
enddo ; enddo
enddo
do J=js-1,je ; do i=is,ie
baroHv(i,J) = 0.
enddo ; enddo
do k=1,nz
do J=js-1,je ; do i=is,ie
baroHv(i,J) = hv(i,J,k) * GV%H_to_RZ
enddo ; enddo
enddo
endif
if (CS%MEKE_Cd_scale == 0.0 .and. .not. CS%visc_drag) then
!$OMP parallel do default(shared) private(ldamping)
do j=js,je ; do i=is,ie
drag_rate(i,j) = 0. ; drag_rate_J15(i,j) = 0.
enddo ; enddo
endif
! Calculate drag_rate_visc(i,j) which accounts for the model bottom mean flow
if (CS%visc_drag) then
!$OMP parallel do default(shared)
do j=js,je ; do I=is-1,ie
drag_vel_u(I,j) = 0.0
if ((G%mask2dCu(I,j) > 0.0) .and. (visc%bbl_thick_u(I,j) > 0.0)) &
drag_vel_u(I,j) = visc%Kv_bbl_u(I,j) / visc%bbl_thick_u(I,j)
enddo ; enddo
!$OMP parallel do default(shared)
do J=js-1,je ; do i=is,ie
drag_vel_v(i,J) = 0.0
if ((G%mask2dCv(i,J) > 0.0) .and. (visc%bbl_thick_v(i,J) > 0.0)) &
drag_vel_v(i,J) = visc%Kv_bbl_v(i,J) / visc%bbl_thick_v(i,J)
enddo ; enddo
!$OMP parallel do default(shared)
do j=js,je ; do i=is,ie
drag_rate_visc(i,j) = (0.25*G%IareaT(i,j) * US%Z_to_L * &
((G%areaCu(I-1,j)*drag_vel_u(I-1,j) + &
G%areaCu(I,j)*drag_vel_u(I,j)) + &
(G%areaCv(i,J-1)*drag_vel_v(i,J-1) + &
G%areaCv(i,J)*drag_vel_v(i,J)) ) )
enddo ; enddo
else
!$OMP parallel do default(shared)
do j=js,je ; do i=is,ie
drag_rate_visc(i,j) = 0.
enddo ; enddo
endif
!$OMP parallel do default(shared)
do j=js-1,je+1
do i=is-1,ie+1 ; mass(i,j) = 0.0 ; enddo
do k=1,nz ; do i=is-1,ie+1
mass(i,j) = mass(i,j) + G%mask2dT(i,j) * (GV%H_to_RZ * h(i,j,k)) ! [R Z ~> kg m-2]
enddo ; enddo
do i=is-1,ie+1
I_mass(i,j) = 0.0
if (mass(i,j) > 0.0) I_mass(i,j) = 1.0 / mass(i,j) ! [R-1 Z-1 ~> m2 kg-1]
enddo
enddo
if (CS%fixed_total_depth) then
!$OMP parallel do default(shared)
do j=js-1,je+1 ; do i=is-1,ie+1
depth_tot(i,j) = G%bathyT(i,j)
enddo ; enddo
else
!$OMP parallel do default(shared)
do j=js-1,je+1 ; do i=is-1,ie+1
depth_tot(i,j) = mass(i,j) * I_Rho0
enddo ; enddo
endif
if (CS%initialize) then
call MEKE_equilibrium(CS, MEKE, G, GV, US, SN_u, SN_v, drag_rate_visc, I_mass, depth_tot)
CS%initialize = .false.
endif
! Calculates bottomFac2, barotrFac2 and LmixScale
call MEKE_lengthScales(CS, MEKE, G, GV, US, SN_u, SN_v, MEKE%MEKE, depth_tot, bottomFac2, barotrFac2, LmixScale)
if (CS%debug) then
if (CS%visc_drag) &
call uvchksum("MEKE drag_vel_[uv]", drag_vel_u, drag_vel_v, G%HI, &
scale=US%Z_to_m*US%s_to_T, scalar_pair=.true.)
call hchksum(mass, 'MEKE mass',G%HI,haloshift=1, scale=US%RZ_to_kg_m2)
call hchksum(drag_rate_visc, 'MEKE drag_rate_visc', G%HI, scale=US%L_T_to_m_s)
call hchksum(bottomFac2, 'MEKE bottomFac2', G%HI)
call hchksum(barotrFac2, 'MEKE barotrFac2', G%HI)
call hchksum(LmixScale, 'MEKE LmixScale', G%HI,scale=US%L_to_m)
endif
! Aggregate sources of MEKE (background, frictional and GM)
!$OMP parallel do default(shared)
do j=js,je ; do i=is,ie
src(i,j) = CS%MEKE_BGsrc
enddo ; enddo
if (associated(MEKE%mom_src)) then
!$OMP parallel do default(shared)
do j=js,je ; do i=is,ie
src(i,j) = src(i,j) - CS%MEKE_FrCoeff*I_mass(i,j)*MEKE%mom_src(i,j)
enddo ; enddo
endif
if (associated(MEKE%GME_snk)) then
!$OMP parallel do default(shared)
do j=js,je ; do i=is,ie
src(i,j) = src(i,j) - CS%MEKE_GMECoeff*I_mass(i,j)*MEKE%GME_snk(i,j)
enddo ; enddo
endif
if (associated(MEKE%GM_src)) then
if (CS%GM_src_alt) then
!$OMP parallel do default(shared)
do j=js,je ; do i=is,ie
src(i,j) = src(i,j) - CS%MEKE_GMcoeff*MEKE%GM_src(i,j) / &
(GV%Rho0 * MAX(1.0*US%m_to_Z, depth_tot(i,j)))
enddo ; enddo
else
!$OMP parallel do default(shared)
do j=js,je ; do i=is,ie
src(i,j) = src(i,j) - CS%MEKE_GMcoeff*I_mass(i,j)*MEKE%GM_src(i,j)
enddo ; enddo
endif
endif
if (CS%MEKE_equilibrium_restoring) then
call MEKE_equilibrium_restoring(CS, G, US, SN_u, SN_v, depth_tot)
do j=js,je ; do i=is,ie
src(i,j) = src(i,j) - CS%MEKE_restoring_rate*(MEKE%MEKE(i,j) - CS%equilibrium_value(i,j))
enddo ; enddo
endif
! Increase EKE by a full time-steps worth of source
!$OMP parallel do default(shared)
do j=js,je ; do i=is,ie
MEKE%MEKE(i,j) = (MEKE%MEKE(i,j) + sdt*src(i,j))*G%mask2dT(i,j)
enddo ; enddo
if (use_drag_rate) then
! Calculate a viscous drag rate (includes BBL contributions from mean flow and eddies)
!$OMP parallel do default(shared)
do j=js,je ; do i=is,ie
drag_rate(i,j) = (US%L_to_Z*Rho0 * I_mass(i,j)) * sqrt( drag_rate_visc(i,j)**2 + &
cdrag2 * ( max(0.0, 2.0*bottomFac2(i,j)*MEKE%MEKE(i,j)) + CS%MEKE_Uscale**2 ) )
enddo ; enddo
endif
! First stage of Strang splitting
!$OMP parallel do default(shared)
do j=js,je ; do i=is,ie
ldamping = CS%MEKE_damping + drag_rate(i,j) * bottomFac2(i,j)
if (MEKE%MEKE(i,j) < 0.) ldamping = 0.
! notice that the above line ensures a damping only if MEKE is positive,
! while leaving MEKE unchanged if it is negative
MEKE%MEKE(i,j) = MEKE%MEKE(i,j) / (1.0 + sdt_damp*ldamping)
MEKE_decay(i,j) = ldamping*G%mask2dT(i,j)
enddo ; enddo
if (CS%kh_flux_enabled .or. CS%MEKE_K4 >= 0.0) then
! Update MEKE in the halos for lateral or bi-harmonic diffusion
call cpu_clock_begin(CS%id_clock_pass)
call do_group_pass(CS%pass_MEKE, G%Domain)
call cpu_clock_end(CS%id_clock_pass)
endif
if (CS%MEKE_K4 >= 0.0) then
! Calculate Laplacian of MEKE using MEKE_uflux and MEKE_vflux as temporary work space.
!$OMP parallel do default(shared)
do j=js-1,je+1 ; do I=is-2,ie+1
! MEKE_uflux is used here as workspace with units of [L2 T-2 ~> m2 s-2].
MEKE_uflux(I,j) = ((G%dy_Cu(I,j)*G%IdxCu(I,j)) * G%mask2dCu(I,j)) * &
(MEKE%MEKE(i+1,j) - MEKE%MEKE(i,j))
! This would have units of [R Z L2 T-2 ~> kg s-2]
! MEKE_uflux(I,j) = ((G%dy_Cu(I,j)*G%IdxCu(I,j)) * &
! ((2.0*mass(i,j)*mass(i+1,j)) / ((mass(i,j)+mass(i+1,j)) + mass_neglect)) ) * &
! (MEKE%MEKE(i+1,j) - MEKE%MEKE(i,j))
enddo ; enddo
!$OMP parallel do default(shared)
do J=js-2,je+1 ; do i=is-1,ie+1
! MEKE_vflux is used here as workspace with units of [L2 T-2 ~> m2 s-2].
MEKE_vflux(i,J) = ((G%dx_Cv(i,J)*G%IdyCv(i,J)) * G%mask2dCv(i,J)) * &
(MEKE%MEKE(i,j+1) - MEKE%MEKE(i,j))
! This would have units of [R Z L2 T-2 ~> kg s-2]
! MEKE_vflux(i,J) = ((G%dx_Cv(i,J)*G%IdyCv(i,J)) * &
! ((2.0*mass(i,j)*mass(i,j+1)) / ((mass(i,j)+mass(i,j+1)) + mass_neglect)) ) * &
! (MEKE%MEKE(i,j+1) - MEKE%MEKE(i,j))
enddo ; enddo
!$OMP parallel do default(shared)
do j=js-1,je+1 ; do i=is-1,ie+1 ! del2MEKE has units [T-2 ~> s-2].
del2MEKE(i,j) = G%IareaT(i,j) * &
((MEKE_uflux(I,j) - MEKE_uflux(I-1,j)) + (MEKE_vflux(i,J) - MEKE_vflux(i,J-1)))
enddo ; enddo
! Bi-harmonic diffusion of MEKE
!$OMP parallel do default(shared) private(K4_here,Inv_K4_max)
do j=js,je ; do I=is-1,ie
K4_here = CS%MEKE_K4 ! [L4 T-1 ~> m4 s-1]
! Limit Kh to avoid CFL violations.
Inv_K4_max = 64.0 * sdt * ((G%dy_Cu(I,j)*G%IdxCu(I,j)) * &
max(G%IareaT(i,j), G%IareaT(i+1,j)))**2
if (K4_here*Inv_K4_max > 0.3) K4_here = 0.3 / Inv_K4_max
! Here the units of MEKE_uflux are [R Z L4 T-3 ~> kg m2 s-3].
MEKE_uflux(I,j) = ((K4_here * (G%dy_Cu(I,j)*G%IdxCu(I,j))) * &
((2.0*mass(i,j)*mass(i+1,j)) / ((mass(i,j)+mass(i+1,j)) + mass_neglect)) ) * &
(del2MEKE(i+1,j) - del2MEKE(i,j))
enddo ; enddo
!$OMP parallel do default(shared) private(K4_here,Inv_K4_max)
do J=js-1,je ; do i=is,ie
K4_here = CS%MEKE_K4 ! [L4 T-1 ~> m4 s-1]
Inv_K4_max = 64.0 * sdt * ((G%dx_Cv(i,J)*G%IdyCv(i,J)) * max(G%IareaT(i,j), G%IareaT(i,j+1)))**2
if (K4_here*Inv_K4_max > 0.3) K4_here = 0.3 / Inv_K4_max
! Here the units of MEKE_vflux are [R Z L4 T-3 ~> kg m2 s-3].
MEKE_vflux(i,J) = ((K4_here * (G%dx_Cv(i,J)*G%IdyCv(i,J))) * &
((2.0*mass(i,j)*mass(i,j+1)) / ((mass(i,j)+mass(i,j+1)) + mass_neglect)) ) * &
(del2MEKE(i,j+1) - del2MEKE(i,j))
enddo ; enddo
! Store change in MEKE arising from the bi-harmonic in del4MEKE [L2 T-2 ~> m2 s-2].
!$OMP parallel do default(shared)
do j=js,je ; do i=is,ie
del4MEKE(i,j) = (sdt*(G%IareaT(i,j)*I_mass(i,j))) * &
((MEKE_uflux(I-1,j) - MEKE_uflux(I,j)) + &
(MEKE_vflux(i,J-1) - MEKE_vflux(i,J)))
enddo ; enddo
endif !
if (CS%kh_flux_enabled) then
! Lateral diffusion of MEKE
Kh_here = max(0., CS%MEKE_Kh)
!$OMP parallel do default(shared) firstprivate(Kh_here) private(Inv_Kh_max)
do j=js,je ; do I=is-1,ie
! Limit Kh to avoid CFL violations.
if (associated(MEKE%Kh)) &
Kh_here = max(0., CS%MEKE_Kh) + &
CS%KhMEKE_Fac*0.5*(MEKE%Kh(i,j)+MEKE%Kh(i+1,j))
if (associated(MEKE%Kh_diff)) &
Kh_here = max(0.,CS%MEKE_Kh) + &
CS%KhMEKE_Fac*0.5*(MEKE%Kh_diff(i,j)+MEKE%Kh_diff(i+1,j))
Inv_Kh_max = 2.0*sdt * ((G%dy_Cu(I,j)*G%IdxCu(I,j)) * &
max(G%IareaT(i,j),G%IareaT(i+1,j)))
if (Kh_here*Inv_Kh_max > 0.25) Kh_here = 0.25 / Inv_Kh_max
Kh_u(I,j) = Kh_here
! Here the units of MEKE_uflux and MEKE_vflux are [R Z L4 T-3 ~> kg m2 s-3].
MEKE_uflux(I,j) = ((Kh_here * (G%dy_Cu(I,j)*G%IdxCu(I,j))) * &
((2.0*mass(i,j)*mass(i+1,j)) / ((mass(i,j)+mass(i+1,j)) + mass_neglect)) ) * &
(MEKE%MEKE(i,j) - MEKE%MEKE(i+1,j))
enddo ; enddo
!$OMP parallel do default(shared) firstprivate(Kh_here) private(Inv_Kh_max)
do J=js-1,je ; do i=is,ie
if (associated(MEKE%Kh)) &
Kh_here = max(0.,CS%MEKE_Kh) + CS%KhMEKE_Fac * 0.5*(MEKE%Kh(i,j)+MEKE%Kh(i,j+1))
if (associated(MEKE%Kh_diff)) &
Kh_here = max(0.,CS%MEKE_Kh) + CS%KhMEKE_Fac * 0.5*(MEKE%Kh_diff(i,j)+MEKE%Kh_diff(i,j+1))
Inv_Kh_max = 2.0*sdt * ((G%dx_Cv(i,J)*G%IdyCv(i,J)) * max(G%IareaT(i,j),G%IareaT(i,j+1)))
if (Kh_here*Inv_Kh_max > 0.25) Kh_here = 0.25 / Inv_Kh_max
Kh_v(i,J) = Kh_here
! Here the units of MEKE_uflux and MEKE_vflux are [R Z L4 T-3 ~> kg m2 s-3].
MEKE_vflux(i,J) = ((Kh_here * (G%dx_Cv(i,J)*G%IdyCv(i,J))) * &
((2.0*mass(i,j)*mass(i,j+1)) / ((mass(i,j)+mass(i,j+1)) + mass_neglect)) ) * &
(MEKE%MEKE(i,j) - MEKE%MEKE(i,j+1))
enddo ; enddo
if (CS%MEKE_advection_factor>0.) then
advFac = CS%MEKE_advection_factor / sdt ! [T-1 ~> s-1]
!$OMP parallel do default(shared)
do j=js,je ; do I=is-1,ie
! Here the units of the quantities added to MEKE_uflux are [R Z L4 T-3 ~> kg m2 s-3].
if (baroHu(I,j)>0.) then
MEKE_uflux(I,j) = MEKE_uflux(I,j) + baroHu(I,j)*MEKE%MEKE(i,j)*advFac
elseif (baroHu(I,j)<0.) then
MEKE_uflux(I,j) = MEKE_uflux(I,j) + baroHu(I,j)*MEKE%MEKE(i+1,j)*advFac
endif
enddo ; enddo
!$OMP parallel do default(shared)
do J=js-1,je ; do i=is,ie
! Here the units of the quantities added to MEKE_vflux are [R Z L4 T-3 ~> kg m2 s-3].
if (baroHv(i,J)>0.) then
MEKE_vflux(i,J) = MEKE_vflux(i,J) + baroHv(i,J)*MEKE%MEKE(i,j)*advFac
elseif (baroHv(i,J)<0.) then
MEKE_vflux(i,J) = MEKE_vflux(i,J) + baroHv(i,J)*MEKE%MEKE(i,j+1)*advFac
endif
enddo ; enddo
endif
!$OMP parallel do default(shared)
do j=js,je ; do i=is,ie
MEKE%MEKE(i,j) = MEKE%MEKE(i,j) + (sdt*(G%IareaT(i,j)*I_mass(i,j))) * &
((MEKE_uflux(I-1,j) - MEKE_uflux(I,j)) + &
(MEKE_vflux(i,J-1) - MEKE_vflux(i,J)))
enddo ; enddo
endif ! MEKE_KH>0
! Add on bi-harmonic tendency
if (CS%MEKE_K4 >= 0.0) then
!$OMP parallel do default(shared)
do j=js,je ; do i=is,ie
MEKE%MEKE(i,j) = MEKE%MEKE(i,j) + del4MEKE(i,j)
enddo ; enddo
endif
! Second stage of Strang splitting
if (CS%MEKE_KH >= 0.0 .or. CS%MEKE_K4 >= 0.0) then
if (sdt>sdt_damp) then
! Recalculate the drag rate, since MEKE has changed.
if (use_drag_rate) then
!$OMP parallel do default(shared)
do j=js,je ; do i=is,ie
drag_rate(i,j) = (US%L_to_Z*Rho0 * I_mass(i,j)) * sqrt( drag_rate_visc(i,j)**2 + &
cdrag2 * ( max(0.0, 2.0*bottomFac2(i,j)*MEKE%MEKE(i,j)) + CS%MEKE_Uscale**2 ) )
enddo ; enddo
!$OMP parallel do default(shared)
do j=js,je ; do i=is,ie
ldamping = CS%MEKE_damping + drag_rate(i,j) * bottomFac2(i,j)
if (MEKE%MEKE(i,j) < 0.) ldamping = 0.
! notice that the above line ensures a damping only if MEKE is positive,
! while leaving MEKE unchanged if it is negative
MEKE%MEKE(i,j) = MEKE%MEKE(i,j) / (1.0 + sdt_damp*ldamping)
MEKE_decay(i,j) = ldamping*G%mask2dT(i,j)
enddo ; enddo
endif
endif
endif ! MEKE_KH>=0
! do j=js,je ; do i=is,ie
! MEKE%MEKE(i,j) = MAX(MEKE%MEKE(i,j),0.0)
! enddo ; enddo
call cpu_clock_begin(CS%id_clock_pass)
call do_group_pass(CS%pass_MEKE, G%Domain)
call cpu_clock_end(CS%id_clock_pass)
! Calculate diffusivity for main model to use
if (CS%MEKE_KhCoeff>0.) then
if (.not.CS%MEKE_GEOMETRIC) then
if (CS%use_old_lscale) then
if (CS%Rd_as_max_scale) then
!$OMP parallel do default(shared)
do j=js,je ; do i=is,ie
MEKE%Kh(i,j) = (CS%MEKE_KhCoeff * &
sqrt(2.*max(0.,barotrFac2(i,j)*MEKE%MEKE(i,j))*G%areaT(i,j)) ) * &
min(MEKE%Rd_dx_h(i,j), 1.0)
enddo ; enddo
else
!$OMP parallel do default(shared)
do j=js,je ; do i=is,ie
MEKE%Kh(i,j) = CS%MEKE_KhCoeff * &
sqrt(2.*max(0., barotrFac2(i,j)*MEKE%MEKE(i,j))*G%areaT(i,j))
enddo ; enddo
endif
else
!$OMP parallel do default(shared)
do j=js,je ; do i=is,ie
MEKE%Kh(i,j) = CS%MEKE_KhCoeff * &
sqrt(2.*max(0., barotrFac2(i,j)*MEKE%MEKE(i,j))) * LmixScale(i,j)
enddo ; enddo
endif
endif
endif
! Calculate viscosity for the main model to use
if (CS%viscosity_coeff_Ku /=0.) then
do j=js,je ; do i=is,ie
MEKE%Ku(i,j) = CS%viscosity_coeff_Ku * sqrt(2.*max(0.,MEKE%MEKE(i,j))) * LmixScale(i,j)
enddo ; enddo
endif
if (CS%viscosity_coeff_Au /=0.) then
do j=js,je ; do i=is,ie
MEKE%Au(i,j) = CS%viscosity_coeff_Au * sqrt(2.*max(0.,MEKE%MEKE(i,j))) * LmixScale(i,j)**3
enddo ; enddo
endif
if (associated(MEKE%Kh) .or. associated(MEKE%Ku) .or. associated(MEKE%Au)) then
call cpu_clock_begin(CS%id_clock_pass)
call do_group_pass(CS%pass_Kh, G%Domain)
call cpu_clock_end(CS%id_clock_pass)
endif
! Offer fields for averaging.
if (any([CS%id_Ue, CS%id_Ub, CS%id_Ut] > 0)) &
tmp(:,:) = 0.
if (CS%id_MEKE>0) call post_data(CS%id_MEKE, MEKE%MEKE, CS%diag)
if (CS%id_Ue>0) then
do j=js,je ; do i=is,ie
tmp(i,j) = sqrt(max(0., 2. * MEKE%MEKE(i,j)))
enddo ; enddo
call post_data(CS%id_Ue, tmp, CS%diag)
endif
if (CS%id_Ub>0) then
do j=js,je ; do i=is,ie
tmp(i,j) = sqrt(max(0., 2. * MEKE%MEKE(i,j) * bottomFac2(i,j)))
enddo ; enddo
call post_data(CS%id_Ub, tmp, CS%diag)
endif
if (CS%id_Ut>0) then
do j=js,je ; do i=is,ie
tmp(i,j) = sqrt(max(0., 2. * MEKE%MEKE(i,j) * barotrFac2(i,j)))
enddo ; enddo
call post_data(CS%id_Ut, tmp, CS%diag)
endif
if (CS%id_Kh>0) call post_data(CS%id_Kh, MEKE%Kh, CS%diag)
if (CS%id_Ku>0) call post_data(CS%id_Ku, MEKE%Ku, CS%diag)
if (CS%id_Au>0) call post_data(CS%id_Au, MEKE%Au, CS%diag)
if (CS%id_KhMEKE_u>0) call post_data(CS%id_KhMEKE_u, Kh_u, CS%diag)
if (CS%id_KhMEKE_v>0) call post_data(CS%id_KhMEKE_v, Kh_v, CS%diag)
if (CS%id_src>0) call post_data(CS%id_src, src, CS%diag)
if (CS%id_decay>0) call post_data(CS%id_decay, MEKE_decay, CS%diag)
if (CS%id_GM_src>0) call post_data(CS%id_GM_src, MEKE%GM_src, CS%diag)
if (CS%id_mom_src>0) call post_data(CS%id_mom_src, MEKE%mom_src, CS%diag)
if (CS%id_GME_snk>0) call post_data(CS%id_GME_snk, MEKE%GME_snk, CS%diag)
if (CS%id_Le>0) call post_data(CS%id_Le, LmixScale, CS%diag)
if (CS%id_gamma_b>0) then
do j=js,je ; do i=is,ie
bottomFac2(i,j) = sqrt(bottomFac2(i,j))
enddo ; enddo
call post_data(CS%id_gamma_b, bottomFac2, CS%diag)
endif
if (CS%id_gamma_t>0) then
do j=js,je ; do i=is,ie
barotrFac2(i,j) = sqrt(barotrFac2(i,j))
enddo ; enddo
call post_data(CS%id_gamma_t, barotrFac2, CS%diag)
endif
end subroutine step_forward_MEKE
!> Calculates the equilibrium solution where the source depends only on MEKE diffusivity
!! and there is no lateral diffusion of MEKE.
!! Results is in MEKE%MEKE.
subroutine MEKE_equilibrium(CS, MEKE, G, GV, US, SN_u, SN_v, drag_rate_visc, I_mass, depth_tot)
type(ocean_grid_type), intent(inout) :: G !< Ocean grid.
type(verticalGrid_type), intent(in) :: GV !< Ocean vertical grid structure.
type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
type(MEKE_CS), pointer :: CS !< MEKE control structure.
type(MEKE_type), pointer :: MEKE !< A structure with MEKE data.
real, dimension(SZIB_(G),SZJ_(G)), intent(in) :: SN_u !< Eady growth rate at u-points [T-1 ~> s-1].
real, dimension(SZI_(G),SZJB_(G)), intent(in) :: SN_v !< Eady growth rate at v-points [T-1 ~> s-1].
real, dimension(SZI_(G),SZJ_(G)), intent(in) :: drag_rate_visc !< Mean flow velocity contribution
!! to the MEKE drag rate [L T-1 ~> m s-1]
real, dimension(SZI_(G),SZJ_(G)), intent(in) :: I_mass !< Inverse of column mass [R-1 Z-1 ~> m2 kg-1].
real, dimension(SZI_(G),SZJ_(G)), intent(in) :: depth_tot !< The depth of the water column [Z ~> m].
! Local variables
real :: beta ! Combined topograpic and planetary vorticity gradient [T-1 L-1 ~> s-1 m-1]
real :: SN ! The local Eady growth rate [T-1 ~> s-1]
real :: bottomFac2, barotrFac2 ! Vertical structure factors [nondim]
real :: LmixScale, LRhines, LEady ! Various mixing length scales [L ~> m]
real :: I_H, KhCoeff
real :: Kh ! A lateral diffusivity [L2 T-1 ~> m2 s-1]
real :: Ubg2 ! Background (tidal?) velocity squared [L2 T-2 ~> m2 s-2]
real :: cd2
real :: drag_rate ! The MEKE spindown timescale due to bottom drag [T-1 ~> s-1].
real :: src ! The sum of MEKE sources [L2 T-3 ~> W kg-1]
real :: ldamping ! The MEKE damping rate [T-1 ~> s-1].
real :: EKE, EKEmin, EKEmax, EKEerr ! [L2 T-2 ~> m2 s-2]
real :: resid, ResMin, ResMax ! Residuals [L2 T-3 ~> W kg-1]
real :: FatH ! Coriolis parameter at h points; to compute topographic beta [T-1 ~> s-1]
real :: beta_topo_x, beta_topo_y ! Topographic PV gradients in x and y [T-1 L-1 ~> s-1 m-1]
real :: dZ_neglect ! A negligible change in height [Z ~> m]
integer :: i, j, is, ie, js, je, n1, n2
real :: tolerance ! Width of EKE bracket [L2 T-2 ~> m2 s-2].
logical :: useSecant, debugIteration
real :: Lgrid, Ldeform, Lfrict
is = G%isc ; ie = G%iec ; js = G%jsc ; je = G%jec
debugIteration = .false.
KhCoeff = CS%MEKE_KhCoeff
Ubg2 = CS%MEKE_Uscale**2
cd2 = CS%cdrag**2
tolerance = 1.0e-12*US%m_s_to_L_T**2
dZ_neglect = GV%H_to_Z*GV%H_subroundoff
!$OMP do
do j=js,je ; do i=is,ie
! SN = 0.25*max( (SN_u(I,j) + SN_u(I-1,j)) + (SN_v(i,J) + SN_v(i,J-1)), 0.)
! This avoids extremes values in equilibrium solution due to bad values in SN_u, SN_v
SN = min(SN_u(I,j), SN_u(I-1,j), SN_v(i,J), SN_v(i,J-1))
if (CS%MEKE_equilibrium_alt) then
MEKE%MEKE(i,j) = (CS%MEKE_GEOMETRIC_alpha * SN * US%Z_to_m*depth_tot(i,j))**2 / cd2
else
FatH = 0.25*((G%CoriolisBu(I,J) + G%CoriolisBu(I-1,J-1)) + &
(G%CoriolisBu(I-1,J) + G%CoriolisBu(I,J-1))) ! Coriolis parameter at h points
! Since zero-bathymetry cells are masked, this avoids calculations on land
if (CS%MEKE_topographic_beta == 0. .or. (depth_tot(i,j) == 0.0)) then
beta_topo_x = 0. ; beta_topo_y = 0.
else
!### Consider different combinations of these estimates of topographic beta.
beta_topo_x = -CS%MEKE_topographic_beta * FatH * 0.5 * ( &
(depth_tot(i+1,j)-depth_tot(i,j)) * G%IdxCu(I,j) &
/ max(depth_tot(i+1,j), depth_tot(i,j), dZ_neglect) &
+ (depth_tot(i,j)-depth_tot(i-1,j)) * G%IdxCu(I-1,j) &
/ max(depth_tot(i,j), depth_tot(i-1,j), dZ_neglect) )
beta_topo_y = -CS%MEKE_topographic_beta * FatH * 0.5 * ( &
(depth_tot(i,j+1)-depth_tot(i,j)) * G%IdyCv(i,J) &
/ max(depth_tot(i,j+1), depth_tot(i,j), dZ_neglect) + &
(depth_tot(i,j)-depth_tot(i,j-1)) * G%IdyCv(i,J-1) &
/ max(depth_tot(i,j), depth_tot(i,j-1), dZ_neglect) )
endif
beta = sqrt((G%dF_dx(i,j) + beta_topo_x)**2 + &
(G%dF_dy(i,j) + beta_topo_y)**2 )
I_H = US%L_to_Z*GV%Rho0 * I_mass(i,j)
if (KhCoeff*SN*I_H>0.) then
! Solve resid(E) = 0, where resid = Kh(E) * (SN)^2 - damp_rate(E) E
EKEmin = 0. ! Use the trivial root as the left bracket
ResMin = 0. ! Need to detect direction of left residual
EKEmax = 0.01*US%m_s_to_L_T**2 ! First guess at right bracket
useSecant = .false. ! Start using a bisection method
! First find right bracket for which resid<0
resid = 1.0*US%m_to_L**2*US%T_to_s**3 ; n1 = 0
do while (resid>0.)
n1 = n1 + 1
EKE = EKEmax
call MEKE_lengthScales_0d(CS, US, G%areaT(i,j), beta, depth_tot(i,j), &
MEKE%Rd_dx_h(i,j), SN, EKE, &
bottomFac2, barotrFac2, LmixScale, LRhines, LEady)
! TODO: Should include resolution function in Kh
Kh = (KhCoeff * sqrt(2.*barotrFac2*EKE) * LmixScale)
src = Kh * (SN * SN)
drag_rate = I_H * sqrt(drag_rate_visc(i,j)**2 + cd2 * ( 2.0*bottomFac2*EKE + Ubg2 ) )
ldamping = CS%MEKE_damping + drag_rate * bottomFac2
resid = src - ldamping * EKE
! if (debugIteration) then
! write(0,*) n1, 'EKE=',EKE,'resid=',resid
! write(0,*) 'EKEmin=',EKEmin,'ResMin=',ResMin
! write(0,*) 'src=',src,'ldamping=',ldamping
! write(0,*) 'gamma-b=',bottomFac2,'gamma-t=',barotrFac2
! write(0,*) 'drag_visc=',drag_rate_visc(i,j),'Ubg2=',Ubg2
! endif
if (resid>0.) then ! EKE is to the left of the root
EKEmin = EKE ! so we move the left bracket here
EKEmax = 10. * EKE ! and guess again for the right bracket
if (resid<ResMin) useSecant = .true.
ResMin = resid
if (EKEmax > 2.e17*US%m_s_to_L_T**2) then
if (debugIteration) stop 'Something has gone very wrong'
debugIteration = .true.
resid = 1. ; n1 = 0
EKEmin = 0. ; ResMin = 0.
EKEmax = 0.01*US%m_s_to_L_T**2
useSecant = .false.
endif
endif
enddo ! while(resid>0.) searching for right bracket
ResMax = resid
! Bisect the bracket
n2 = 0 ; EKEerr = EKEmax - EKEmin
do while (EKEerr > tolerance)
n2 = n2 + 1
if (useSecant) then
EKE = EKEmin + (EKEmax - EKEmin) * (ResMin / (ResMin - ResMax))
else
EKE = 0.5 * (EKEmin + EKEmax)
endif
EKEerr = min( EKE-EKEmin, EKEmax-EKE )
! TODO: Should include resolution function in Kh
Kh = (KhCoeff * sqrt(2.*barotrFac2*EKE) * LmixScale)
src = Kh * (SN * SN)
drag_rate = I_H * sqrt( drag_rate_visc(i,j)**2 + cd2 * ( 2.0*bottomFac2*EKE + Ubg2 ) )
ldamping = CS%MEKE_damping + drag_rate * bottomFac2
resid = src - ldamping * EKE
if (useSecant .and. resid>ResMin) useSecant = .false.
if (resid>0.) then ! EKE is to the left of the root
EKEmin = EKE ! so we move the left bracket here
if (resid<ResMin) useSecant = .true.
ResMin = resid ! Save this for the secant method
elseif (resid<0.) then ! EKE is to the right of the root
EKEmax = EKE ! so we move the right bracket here
ResMax = resid ! Save this for the secant method
else
exit ! resid=0 => EKE is exactly at the root
endif
if (n2>200) stop 'Failing to converge?'
enddo ! while(EKEmax-EKEmin>tolerance)
else
EKE = 0.
endif
MEKE%MEKE(i,j) = EKE
endif
enddo ; enddo
end subroutine MEKE_equilibrium
!< This subroutine calculates a new equilibrium value for MEKE at each time step. This is not copied into
!! MEKE%MEKE; rather, it is used as a restoring term to nudge MEKE%MEKE back to an equilibrium value
subroutine MEKE_equilibrium_restoring(CS, G, US, SN_u, SN_v, depth_tot)
type(ocean_grid_type), intent(inout) :: G !< Ocean grid.
type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type.
type(MEKE_CS), pointer :: CS !< MEKE control structure.
real, dimension(SZIB_(G),SZJ_(G)), intent(in) :: SN_u !< Eady growth rate at u-points [T-1 ~> s-1].
real, dimension(SZI_(G),SZJB_(G)), intent(in) :: SN_v !< Eady growth rate at v-points [T-1 ~> s-1].
real, dimension(SZI_(G),SZJ_(G)), intent(in) :: depth_tot !< The depth of the water column [Z ~> m].
! Local variables
real :: SN ! The local Eady growth rate [T-1 ~> s-1]
integer :: i, j, is, ie, js, je ! local indices
real :: cd2 ! bottom drag
is = G%isc ; ie = G%iec ; js = G%jsc ; je = G%jec
cd2 = CS%cdrag**2
if (.not. associated(CS%equilibrium_value)) allocate(CS%equilibrium_value(SZI_(G),SZJ_(G)))
CS%equilibrium_value(:,:) = 0.0
!$OMP do
do j=js,je ; do i=is,ie
! SN = 0.25*max( (SN_u(I,j) + SN_u(I-1,j)) + (SN_v(i,J) + SN_v(i,J-1)), 0.)
! This avoids extremes values in equilibrium solution due to bad values in SN_u, SN_v
SN = min(SN_u(I,j), SN_u(I-1,j), SN_v(i,J), SN_v(i,J-1))
CS%equilibrium_value(i,j) = (CS%MEKE_GEOMETRIC_alpha * SN * US%Z_to_L*depth_tot(i,j))**2 / cd2
enddo ; enddo
if (CS%id_MEKE_equilibrium>0) call post_data(CS%id_MEKE_equilibrium, CS%equilibrium_value, CS%diag)
end subroutine MEKE_equilibrium_restoring
!> Calculates the eddy mixing length scale and \f$\gamma_b\f$ and \f$\gamma_t\f$
!! functions that are ratios of either bottom or barotropic eddy energy to the
!! column eddy energy, respectively. See \ref section_MEKE_equations.
subroutine MEKE_lengthScales(CS, MEKE, G, GV, US, SN_u, SN_v, EKE, depth_tot, &
bottomFac2, barotrFac2, LmixScale)
type(MEKE_CS), pointer :: CS !< MEKE control structure.
type(MEKE_type), pointer :: MEKE !< MEKE data.
type(ocean_grid_type), intent(inout) :: G !< Ocean grid.
type(verticalGrid_type), intent(in) :: GV !< Ocean vertical grid structure.
type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
real, dimension(SZIB_(G),SZJ_(G)), intent(in) :: SN_u !< Eady growth rate at u-points [T-1 ~> s-1].
real, dimension(SZI_(G),SZJB_(G)), intent(in) :: SN_v !< Eady growth rate at v-points [T-1 ~> s-1].
real, dimension(SZI_(G),SZJ_(G)), intent(in) :: EKE !< Eddy kinetic energy [L2 T-2 ~> m2 s-2].
real, dimension(SZI_(G),SZJ_(G)), intent(in) :: depth_tot !< The depth of the water column [Z ~> m].
real, dimension(SZI_(G),SZJ_(G)), intent(out) :: bottomFac2 !< gamma_b^2 [nondim]
real, dimension(SZI_(G),SZJ_(G)), intent(out) :: barotrFac2 !< gamma_t^2 [nondim]
real, dimension(SZI_(G),SZJ_(G)), intent(out) :: LmixScale !< Eddy mixing length [L ~> m].
! Local variables
real, dimension(SZI_(G),SZJ_(G)) :: LRhines, LEady ! Possible mixing length scales [L ~> m]
real :: beta ! Combined topograpic and planetary vorticity gradient [T-1 L-1 ~> s-1 m-1]
real :: SN ! The local Eady growth rate [T-1 ~> s-1]
real :: FatH ! Coriolis parameter at h points [T-1 ~> s-1]
real :: beta_topo_x, beta_topo_y ! Topographic PV gradients in x and y [T-1 L-1 ~> s-1 m-1]
real :: dZ_neglect ! A negligible change in height [Z ~> m]
integer :: i, j, is, ie, js, je
is = G%isc ; ie = G%iec ; js = G%jsc ; je = G%jec
dZ_neglect = GV%H_to_Z*GV%H_subroundoff
!$OMP do
do j=js,je ; do i=is,ie
if (.not.CS%use_old_lscale) then
if (CS%aEady > 0.) then
SN = 0.25 * ( (SN_u(I,j) + SN_u(I-1,j)) + (SN_v(i,J) + SN_v(i,J-1)) )
else
SN = 0.
endif
FatH = 0.25* ( ( G%CoriolisBu(I,J) + G%CoriolisBu(I-1,J-1) ) + &
( G%CoriolisBu(I-1,J) + G%CoriolisBu(I,J-1) ) ) ! Coriolis parameter at h points
! If depth_tot is zero, then a division by zero FPE will be raised. In this
! case, we apply Adcroft's rule of reciprocals and set the term to zero.
! Since zero-bathymetry cells are masked, this should not affect values.
if (CS%MEKE_topographic_beta == 0. .or. (depth_tot(i,j) == 0.0)) then
beta_topo_x = 0. ; beta_topo_y = 0.
else
!### Consider different combinations of these estimates of topographic beta.
beta_topo_x = -CS%MEKE_topographic_beta * FatH * 0.5 * ( &
(depth_tot(i+1,j)-depth_tot(i,j)) * G%IdxCu(I,j) &
/ max(depth_tot(i+1,j), depth_tot(i,j), dZ_neglect) &
+ (depth_tot(i,j)-depth_tot(i-1,j)) * G%IdxCu(I-1,j) &
/ max(depth_tot(i,j), depth_tot(i-1,j), dZ_neglect) )
beta_topo_y = -CS%MEKE_topographic_beta * FatH * 0.5 * ( &
(depth_tot(i,j+1)-depth_tot(i,j)) * G%IdyCv(i,J) &
/ max(depth_tot(i,j+1), depth_tot(i,j), dZ_neglect) + &
(depth_tot(i,j)-depth_tot(i,j-1)) * G%IdyCv(i,J-1) &
/ max(depth_tot(i,j), depth_tot(i,j-1), dZ_neglect) )
endif
beta = sqrt((G%dF_dx(i,j) + beta_topo_x)**2 + &
(G%dF_dy(i,j) + beta_topo_y)**2 )
else
beta = 0.
endif
! Returns bottomFac2, barotrFac2 and LmixScale
call MEKE_lengthScales_0d(CS, US, G%areaT(i,j), beta, depth_tot(i,j), &
MEKE%Rd_dx_h(i,j), SN, MEKE%MEKE(i,j), &
bottomFac2(i,j), barotrFac2(i,j), LmixScale(i,j), &
LRhines(i,j), LEady(i,j))
enddo ; enddo
if (CS%id_Lrhines>0) call post_data(CS%id_LRhines, LRhines, CS%diag)
if (CS%id_Leady>0) call post_data(CS%id_LEady, LEady, CS%diag)
end subroutine MEKE_lengthScales
!> Calculates the eddy mixing length scale and \f$\gamma_b\f$ and \f$\gamma_t\f$
!! functions that are ratios of either bottom or barotropic eddy energy to the
!! column eddy energy, respectively. See \ref section_MEKE_equations.
subroutine MEKE_lengthScales_0d(CS, US, area, beta, depth, Rd_dx, SN, EKE, & ! Z_to_L, &
bottomFac2, barotrFac2, LmixScale, Lrhines, Leady)
type(MEKE_CS), pointer :: CS !< MEKE control structure.
type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
real, intent(in) :: area !< Grid cell area [L2 ~> m2]
real, intent(in) :: beta !< Planetary beta = \f$ \nabla f\f$ [T-1 L-1 ~> s-1 m-1]
real, intent(in) :: depth !< Ocean depth [Z ~> m]
real, intent(in) :: Rd_dx !< Resolution Ld/dx [nondim].
real, intent(in) :: SN !< Eady growth rate [T-1 ~> s-1].
real, intent(in) :: EKE !< Eddy kinetic energy [L2 T-2 ~> m2 s-2].
! real, intent(in) :: Z_to_L !< A conversion factor from depth units (Z) to
! !! the units for lateral distances (L).
real, intent(out) :: bottomFac2 !< gamma_b^2
real, intent(out) :: barotrFac2 !< gamma_t^2
real, intent(out) :: LmixScale !< Eddy mixing length [L ~> m].
real, intent(out) :: Lrhines !< Rhines length scale [L ~> m].
real, intent(out) :: Leady !< Eady length scale [L ~> m].
! Local variables
real :: Lgrid, Ldeform, Lfrict ! Length scales [L ~> m]
real :: Ue ! An eddy velocity [L T-1 ~> m s-1]
! Length scale for MEKE derived diffusivity
Lgrid = sqrt(area) ! Grid scale
Ldeform = Lgrid * Rd_dx ! Deformation scale
Lfrict = (US%Z_to_L * depth) / CS%cdrag ! Frictional arrest scale
! gamma_b^2 is the ratio of bottom eddy energy to mean column eddy energy
! used in calculating bottom drag
bottomFac2 = CS%MEKE_CD_SCALE**2
if (Lfrict*CS%MEKE_Cb>0.) bottomFac2 = bottomFac2 + 1./( 1. + CS%MEKE_Cb*(Ldeform/Lfrict) )**0.8
bottomFac2 = max(bottomFac2, CS%MEKE_min_gamma)
! gamma_t^2 is the ratio of barotropic eddy energy to mean column eddy energy
! used in the velocity scale for diffusivity
barotrFac2 = 1.
if (Lfrict*CS%MEKE_Ct>0.) barotrFac2 = 1. / ( 1. + CS%MEKE_Ct*(Ldeform/Lfrict) )**0.25
barotrFac2 = max(barotrFac2, CS%MEKE_min_gamma)
if (CS%use_old_lscale) then
if (CS%Rd_as_max_scale) then
LmixScale = min(Ldeform, Lgrid) ! The smaller of Ld or dx
else
LmixScale = Lgrid
endif
else
Ue = sqrt( 2.0 * max( 0., barotrFac2*EKE ) ) ! Barotropic eddy flow scale
Lrhines = sqrt( Ue / max( beta, 1.e-30*US%T_to_s*US%L_to_m ) ) ! Rhines scale
if (CS%aEady > 0.) then
Leady = Ue / max( SN, 1.e-15*US%T_to_s ) ! Bound Eady time-scale < 1e15 seconds
else
Leady = 0.
endif
if (CS%use_min_lscale) then
LmixScale = 1.e7
if (CS%aDeform*Ldeform > 0.) LmixScale = min(LmixScale,CS%aDeform*Ldeform)
if (CS%aFrict *Lfrict > 0.) LmixScale = min(LmixScale,CS%aFrict *Lfrict)
if (CS%aRhines*Lrhines > 0.) LmixScale = min(LmixScale,CS%aRhines*Lrhines)
if (CS%aEady *Leady > 0.) LmixScale = min(LmixScale,CS%aEady *Leady)
if (CS%aGrid *Lgrid > 0.) LmixScale = min(LmixScale,CS%aGrid *Lgrid)
if (CS%Lfixed > 0.) LmixScale = min(LmixScale,CS%Lfixed)
else
LmixScale = 0.
if (CS%aDeform*Ldeform > 0.) LmixScale = LmixScale + 1./(CS%aDeform*Ldeform)
if (CS%aFrict *Lfrict > 0.) LmixScale = LmixScale + 1./(CS%aFrict *Lfrict)
if (CS%aRhines*Lrhines > 0.) LmixScale = LmixScale + 1./(CS%aRhines*Lrhines)
if (CS%aEady *Leady > 0.) LmixScale = LmixScale + 1./(CS%aEady *Leady)
if (CS%aGrid *Lgrid > 0.) LmixScale = LmixScale + 1./(CS%aGrid *Lgrid)
if (CS%Lfixed > 0.) LmixScale = LmixScale + 1./CS%Lfixed
if (LmixScale > 0.) LmixScale = 1. / LmixScale
endif