forked from mom-ocean/MOM6
-
Notifications
You must be signed in to change notification settings - Fork 60
/
MOM_tidal_mixing.F90
1747 lines (1537 loc) · 95.9 KB
/
MOM_tidal_mixing.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
!> Interface to vertical tidal mixing schemes including CVMix tidal mixing.
module MOM_tidal_mixing
! This file is part of MOM6. See LICENSE.md for the license.
use MOM_diag_mediator, only : diag_ctrl, time_type, register_diag_field
use MOM_diag_mediator, only : safe_alloc_ptr, post_data
use MOM_debugging, only : hchksum
use MOM_error_handler, only : MOM_error, is_root_pe, FATAL, WARNING, NOTE
use MOM_file_parser, only : openParameterBlock, closeParameterBlock
use MOM_file_parser, only : get_param, log_param, log_version, param_file_type
use MOM_grid, only : ocean_grid_type
use MOM_io, only : slasher, MOM_read_data, field_size
use MOM_io, only : read_netCDF_data
use MOM_internal_tides, only : int_tide_CS, get_lowmode_loss
use MOM_remapping, only : remapping_CS, initialize_remapping, remapping_core_h
use MOM_string_functions, only : uppercase, lowercase
use MOM_unit_scaling, only : unit_scale_type
use MOM_variables, only : thermo_var_ptrs, p3d
use MOM_verticalGrid, only : verticalGrid_type
use CVMix_tidal, only : CVMix_init_tidal, CVMix_compute_Simmons_invariant
use CVMix_tidal, only : CVMix_coeffs_tidal, CVMix_tidal_params_type
use CVMix_tidal, only : CVMix_compute_Schmittner_invariant, CVMix_compute_SchmittnerCoeff
use CVMix_tidal, only : CVMix_coeffs_tidal_schmittner
use CVMix_kinds_and_types, only : CVMix_global_params_type
use CVMix_put_get, only : CVMix_put
implicit none ; private
#include <MOM_memory.h>
public tidal_mixing_init
public setup_tidal_diagnostics
public calculate_tidal_mixing
public post_tidal_diagnostics
public tidal_mixing_h_amp
public tidal_mixing_end
! 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.
!> Containers for tidal mixing diagnostics
type, public :: tidal_mixing_diags ; private
real, allocatable :: Kd_itidal(:,:,:) !< internal tide diffusivity at interfaces [Z2 T-1 ~> m2 s-1].
real, allocatable :: Fl_itidal(:,:,:) !< vertical flux of tidal turbulent dissipation [Z3 T-3 ~> m3 s-3]
real, allocatable :: Kd_Niku(:,:,:) !< lee-wave diffusivity at interfaces [Z2 T-1 ~> m2 s-1].
real, allocatable :: Kd_Niku_work(:,:,:) !< layer integrated work by lee-wave driven mixing [R Z3 T-3 ~> W m-2]
real, allocatable :: Kd_Itidal_Work(:,:,:) !< layer integrated work by int tide driven mixing [R Z3 T-3 ~> W m-2]
real, allocatable :: Kd_Lowmode_Work(:,:,:) !< layer integrated work by low mode driven mixing [R Z3 T-3 ~> W m-2]
real, allocatable :: N2_int(:,:,:) !< Buoyancy frequency squared at interfaces [T-2 ~> s-2]
real, allocatable :: vert_dep_3d(:,:,:) !< The 3-d mixing energy deposition vertical fraction [nondim]?
real, allocatable :: Schmittner_coeff_3d(:,:,:) !< The coefficient in the Schmittner et al mixing scheme [nondim]
real, allocatable :: tidal_qe_md(:,:,:) !< Input tidal energy dissipated locally,
!! interpolated to model vertical coordinate [R Z3 T-3 ~> W m-2]
real, allocatable :: Kd_lowmode(:,:,:) !< internal tide diffusivity at interfaces
!! due to propagating low modes [Z2 T-1 ~> m2 s-1].
real, allocatable :: Fl_lowmode(:,:,:) !< vertical flux of tidal turbulent
!! dissipation due to propagating low modes [Z3 T-3 ~> m3 s-3]
real, allocatable :: TKE_itidal_used(:,:) !< internal tide TKE input at ocean bottom [R Z3 T-3 ~> W m-2]
real, allocatable :: N2_bot(:,:) !< bottom squared buoyancy frequency [T-2 ~> s-2]
real, allocatable :: N2_meanz(:,:) !< vertically averaged buoyancy frequency [T-2 ~> s-2]
real, allocatable :: Polzin_decay_scale_scaled(:,:) !< vertical scale of decay for tidal dissipation [Z ~> m]
real, allocatable :: Polzin_decay_scale(:,:) !< vertical decay scale for tidal dissipation with Polzin [Z ~> m]
real, allocatable :: Simmons_coeff_2d(:,:) !< The Simmons et al mixing coefficient [nondim]
end type
!> Control structure with parameters for the tidal mixing module.
type, public :: tidal_mixing_cs ; private
logical :: debug = .true. !< If true, do more extensive debugging checks. This is hard-coded.
! Parameters
logical :: int_tide_dissipation = .false. !< Internal tide conversion (from barotropic)
!! with the schemes of St Laurent et al (2002) & Simmons et al (2004)
integer :: Int_tide_profile !< A coded integer indicating the vertical profile
!! for dissipation of the internal waves. Schemes that are
!! currently encoded are St Laurent et al (2002) and Polzin (2009).
logical :: Lee_wave_dissipation = .false. !< Enable lee-wave driven mixing, following
!! Nikurashin (2010), with a vertical energy
!! deposition profile specified by Lee_wave_profile to be
!! St Laurent et al (2002) or Simmons et al (2004) scheme
integer :: Lee_wave_profile !< A coded integer indicating the vertical profile
!! for dissipation of the lee waves. Schemes that are
!! currently encoded are St Laurent et al (2002) and
!! Polzin (2009).
real :: Int_tide_decay_scale !< decay scale for internal wave TKE [Z ~> m].
real :: Mu_itides !< efficiency for conversion of dissipation
!! to potential energy [nondim]
real :: Gamma_itides !< fraction of local dissipation [nondim]
real :: Gamma_lee !< fraction of local dissipation for lee waves
!! (Nikurashin's energy input) [nondim]
real :: Decay_scale_factor_lee !< Scaling factor for the decay scale of lee
!! wave energy dissipation [nondim]
real :: min_zbot_itides !< minimum depth for internal tide conversion [Z ~> m].
logical :: Lowmode_itidal_dissipation = .false. !< If true, consider mixing due to breaking low
!! modes that have been remotely generated using an internal tidal
!! dissipation scheme to specify the vertical profile of the energy
!! input to drive diapycnal mixing, along the lines of St. Laurent
!! et al. (2002) and Simmons et al. (2004).
real :: Nu_Polzin !< The non-dimensional constant used in Polzin form of
!! the vertical scale of decay of tidal dissipation [nondim]
real :: Nbotref_Polzin !< Reference value for the buoyancy frequency at the
!! ocean bottom used in Polzin formulation of the
!! vertical scale of decay of tidal dissipation [T-1 ~> s-1]
real :: Polzin_decay_scale_factor !< Scaling factor for the decay length scale
!! of the tidal dissipation profile in Polzin [nondim]
real :: Polzin_decay_scale_max_factor !< The decay length scale of tidal dissipation
!! profile in Polzin formulation should not exceed
!! Polzin_decay_scale_max_factor * depth of the ocean [nondim].
real :: Polzin_min_decay_scale !< minimum decay scale of the tidal dissipation
!! profile in Polzin formulation [Z ~> m].
real :: TKE_itide_max !< maximum internal tide conversion [R Z3 T-3 ~> W m-2]
!! available to mix above the BBL
real :: utide !< constant tidal amplitude [Z T-1 ~> m s-1] if READ_TIDEAMP is false.
real :: kappa_itides !< topographic wavenumber and non-dimensional scaling [Z-1 ~> m-1].
real :: kappa_h2_factor !< factor for the product of wavenumber * rms sgs height [nondim]
character(len=200) :: inputdir !< The directory in which to find input files
logical :: use_CVMix_tidal = .false. !< true if CVMix is to be used for determining
!! diffusivity due to tidal mixing
real :: min_thickness !< Minimum thickness allowed [Z ~> m]
! CVMix-specific parameters
integer :: CVMix_tidal_scheme = -1 !< 1 for Simmons, 2 for Schmittner
type(CVMix_tidal_params_type) :: CVMix_tidal_params !< A CVMix-specific type with parameters for tidal mixing
type(CVMix_global_params_type) :: CVMix_glb_params !< CVMix-specific for Prandtl number only
real :: tidal_max_coef !< CVMix-specific maximum allowable tidal
!! diffusivity. [Z2 T-1 ~> m2 s-1]
real :: tidal_diss_lim_tc !< CVMix-specific dissipation limit depth for
!! tidal-energy-constituent data [Z ~> m].
type(remapping_CS) :: remap_CS !< The control structure for remapping
integer :: remap_answer_date !< The vintage of the order of arithmetic and expressions to use
!! for remapping. Values below 20190101 recover the remapping
!! answers from 2018, while higher values use more robust
!! forms of the same remapping expressions.
integer :: tidal_answer_date !< The vintage of the order of arithmetic and expressions in the tidal
!! mixing calculations. Values below 20190101 recover the answers
!! from the end of 2018, while higher values use updated and more robust
!! forms of the same expressions.
type(int_tide_CS), pointer :: int_tide_CSp=> NULL() !< Control structure for a child module
! Data containers
real, allocatable :: TKE_Niku(:,:) !< Lee wave driven Turbulent Kinetic Energy input
!! [R Z3 T-3 ~> W m-2]
real, allocatable :: TKE_itidal(:,:) !< The internal Turbulent Kinetic Energy input divided
!! by the bottom stratification [R Z3 T-2 ~> J m-2].
real, allocatable :: Nb(:,:) !< The near bottom buoyancy frequency [T-1 ~> s-1].
real, allocatable :: mask_itidal(:,:) !< A mask of where internal tide energy is input [nondim]
real, allocatable :: h2(:,:) !< Squared bottom depth variance [Z2 ~> m2].
real, allocatable :: tideamp(:,:) !< RMS tidal amplitude [Z T-1 ~> m s-1]
real, allocatable :: h_src(:) !< tidal constituent input layer thickness [m]
real, allocatable :: tidal_qe_2d(:,:) !< Tidal energy input times the local dissipation
!! fraction, q*E(x,y), with the CVMix implementation
!! of Jayne et al tidal mixing [R Z3 T-3 ~> W m-2].
!! TODO: make this E(x,y) only
real, allocatable :: tidal_qe_3d_in(:,:,:) !< q*E(x,y,z) with the Schmittner parameterization [R Z3 T-3 ~> W m-2]
! Diagnostics
type(diag_ctrl), pointer :: diag => NULL() !< structure to regulate diagnostic output timing
type(tidal_mixing_diags) :: dd !< Tidal mixing diagnostic arrays
!>@{ Diagnostic identifiers
integer :: id_TKE_itidal = -1
integer :: id_TKE_leewave = -1
integer :: id_Kd_itidal = -1
integer :: id_Kd_Niku = -1
integer :: id_Kd_lowmode = -1
integer :: id_Kd_Itidal_Work = -1
integer :: id_Kd_Niku_Work = -1
integer :: id_Kd_Lowmode_Work = -1
integer :: id_Nb = -1
integer :: id_N2_bot = -1
integer :: id_N2_meanz = -1
integer :: id_Fl_itidal = -1
integer :: id_Fl_lowmode = -1
integer :: id_Polzin_decay_scale = -1
integer :: id_Polzin_decay_scale_scaled = -1
integer :: id_N2_int = -1
integer :: id_Simmons_coeff = -1
integer :: id_Schmittner_coeff = -1
integer :: id_tidal_qe_md = -1
integer :: id_vert_dep = -1
!>@}
end type tidal_mixing_cs
!>@{ Coded parmameters for specifying mixing schemes
character*(20), parameter :: STLAURENT_PROFILE_STRING = "STLAURENT_02"
character*(20), parameter :: POLZIN_PROFILE_STRING = "POLZIN_09"
integer, parameter :: STLAURENT_02 = 1
integer, parameter :: POLZIN_09 = 2
character*(20), parameter :: SIMMONS_SCHEME_STRING = "SIMMONS"
character*(20), parameter :: SCHMITTNER_SCHEME_STRING = "SCHMITTNER"
integer, parameter :: SIMMONS = 1
integer, parameter :: SCHMITTNER = 2
!>@}
contains
!> Initializes internal tidal dissipation scheme for diapycnal mixing
logical function tidal_mixing_init(Time, G, GV, US, param_file, int_tide_CSp, diag, CS)
type(time_type), intent(in) :: Time !< The current time.
type(ocean_grid_type), intent(in) :: G !< Grid structure.
type(verticalGrid_type), intent(in) :: GV !< Vertical grid structure.
type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
type(param_file_type), intent(in) :: param_file !< Run-time parameter file handle
type(int_tide_CS),target, intent(in) :: int_tide_CSp !< A pointer to the internal tides control structure
type(diag_ctrl), target, intent(inout) :: diag !< Diagnostics control structure.
type(tidal_mixing_cs), intent(inout) :: CS !< This module's control structure.
! Local variables
logical :: use_CVMix_tidal
logical :: int_tide_dissipation
logical :: read_tideamp
integer :: default_answer_date ! The default setting for the various ANSWER_DATE flags.
logical :: default_2018_answers ! The default setting for the various 2018_ANSWERS flags.
logical :: remap_answers_2018 ! If true, use the order of arithmetic and expressions that
! recover the remapping answers from 2018. If false, use more
! robust forms of the same remapping expressions.
integer :: default_remap_ans_date ! The default setting for remap_answer_date
integer :: default_tide_ans_date ! The default setting for tides_answer_date
logical :: tide_answers_2018 ! If true, use the order of arithmetic and expressions that recover the
! answers from the end of 2018. Otherwise, use updated and more robust
! forms of the same expressions.
character(len=20) :: tmpstr, int_tide_profile_str
character(len=20) :: CVMix_tidal_scheme_str, tidal_energy_type
character(len=200) :: filename, h2_file, Niku_TKE_input_file ! Input file names
character(len=200) :: tideamp_file ! Input file names or paths
character(len=80) :: tideamp_var, rough_var, TKE_input_var ! Input file variable names
real :: hamp ! The magnitude of the sub-gridscale bottom depth variance [Z ~> m]
real :: utide ! The RMS tidal amplitude [Z T-1 ~> m s-1]
real :: max_frac_rough ! A limit on the depth variance as a fraction of the total depth [nondim]
real :: prandtl_tidal ! Prandtl number used by CVMix tidal mixing schemes to convert vertical
! diffusivities into viscosities [nondim]
real :: Niku_scale ! local variable for scaling the Nikurashin TKE flux data [nondim]
integer :: i, j, is, ie, js, je
integer :: isd, ied, jsd, jed
! This include declares and sets the variable "version".
# include "version_variable.h"
character(len=40) :: mdl = "MOM_tidal_mixing" !< This module's name.
is = G%isc ; ie = G%iec ; js = G%jsc ; je = G%jec
isd = G%isd ; ied = G%ied ; jsd = G%jsd ; jed = G%jed
! Read parameters
! NOTE: These are read twice because logfile output is streamed and we want
! to preserve the ordering of module header before parameters.
call get_param(param_file, mdl, "USE_CVMix_TIDAL", use_CVMix_tidal, &
default=.false., do_not_log=.true.)
call get_param(param_file, mdl, "INT_TIDE_DISSIPATION", int_tide_dissipation, &
default=use_CVMix_tidal, do_not_log=.true.)
call log_version(param_file, mdl, version, &
"Vertical Tidal Mixing Parameterization", &
all_default=.not.(use_CVMix_tidal .or. int_tide_dissipation))
call get_param(param_file, mdl, "USE_CVMix_TIDAL", use_CVMix_tidal, &
"If true, turns on tidal mixing via CVMix", &
default=.false.)
call get_param(param_file, mdl, "INT_TIDE_DISSIPATION", int_tide_dissipation, &
"If true, use an internal tidal dissipation scheme to "//&
"drive diapycnal mixing, along the lines of St. Laurent "//&
"et al. (2002) and Simmons et al. (2004).", default=use_CVMix_tidal)
! return if tidal mixing is inactive
tidal_mixing_init = int_tide_dissipation
if (.not. tidal_mixing_init) return
CS%debug = CS%debug.and.is_root_pe()
CS%diag => diag
CS%int_tide_CSp => int_tide_CSp
CS%use_CVmix_tidal = use_CVmix_tidal
CS%int_tide_dissipation = int_tide_dissipation
call get_param(param_file, mdl, "INPUTDIR", CS%inputdir, default=".",do_not_log=.true.)
CS%inputdir = slasher(CS%inputdir)
call get_param(param_file, mdl, "DEFAULT_ANSWER_DATE", default_answer_date, &
"This sets the default value for the various _ANSWER_DATE parameters.", &
default=99991231)
call get_param(param_file, mdl, "DEFAULT_2018_ANSWERS", default_2018_answers, &
"This sets the default value for the various _2018_ANSWERS parameters.", &
default=(default_answer_date<20190101))
call get_param(param_file, mdl, "TIDAL_MIXING_2018_ANSWERS", tide_answers_2018, &
"If true, use the order of arithmetic and expressions that recover the "//&
"answers from the end of 2018. Otherwise, use updated and more robust "//&
"forms of the same expressions.", default=default_2018_answers)
! Revise inconsistent default answer dates for the tidal mixing.
default_tide_ans_date = default_answer_date
if (tide_answers_2018 .and. (default_tide_ans_date >= 20190101)) default_tide_ans_date = 20181231
if (.not.tide_answers_2018 .and. (default_tide_ans_date < 20190101)) default_tide_ans_date = 20190101
call get_param(param_file, mdl, "TIDAL_MIXING_ANSWER_DATE", CS%tidal_answer_date, &
"The vintage of the order of arithmetic and expressions in the tidal mixing "//&
"calculations. Values below 20190101 recover the answers from the end of 2018, "//&
"while higher values use updated and more robust forms of the same expressions. "//&
"If both TIDAL_MIXING_2018_ANSWERS and TIDAL_MIXING_ANSWER_DATE are specified, "//&
"the latter takes precedence.", default=default_tide_ans_date)
call get_param(param_file, mdl, "REMAPPING_2018_ANSWERS", remap_answers_2018, &
"If true, use the order of arithmetic and expressions that recover the "//&
"answers from the end of 2018. Otherwise, use updated and more robust "//&
"forms of the same expressions.", default=default_2018_answers)
! Revise inconsistent default answer dates for remapping.
default_remap_ans_date = default_answer_date
if (remap_answers_2018 .and. (default_remap_ans_date >= 20190101)) default_remap_ans_date = 20181231
if (.not.remap_answers_2018 .and. (default_remap_ans_date < 20190101)) default_remap_ans_date = 20190101
call get_param(param_file, mdl, "REMAPPING_ANSWER_DATE", CS%remap_answer_date, &
"The vintage of the expressions and order of arithmetic to use for remapping. "//&
"Values below 20190101 result in the use of older, less accurate expressions "//&
"that were in use at the end of 2018. Higher values result in the use of more "//&
"robust and accurate forms of mathematically equivalent expressions. "//&
"If both REMAPPING_2018_ANSWERS and REMAPPING_ANSWER_DATE are specified, the "//&
"latter takes precedence.", default=default_remap_ans_date)
if (CS%int_tide_dissipation) then
! Read in CVMix tidal scheme if CVMix tidal mixing is on
if (CS%use_CVMix_tidal) then
call get_param(param_file, mdl, "CVMIX_TIDAL_SCHEME", CVMix_tidal_scheme_str, &
"CVMIX_TIDAL_SCHEME selects the CVMix tidal mixing "//&
"scheme with INT_TIDE_DISSIPATION. Valid values are:\n"//&
"\t SIMMONS - Use the Simmons et al (2004) tidal \n"//&
"\t mixing scheme.\n"//&
"\t SCHMITTNER - Use the Schmittner et al (2014) tidal \n"//&
"\t mixing scheme.", &
default=SIMMONS_SCHEME_STRING)
CVMix_tidal_scheme_str = uppercase(CVMix_tidal_scheme_str)
select case (CVMix_tidal_scheme_str)
case (SIMMONS_SCHEME_STRING) ; CS%CVMix_tidal_scheme = SIMMONS
case (SCHMITTNER_SCHEME_STRING) ; CS%CVMix_tidal_scheme = SCHMITTNER
case default
call MOM_error(FATAL, "tidal_mixing_init: Unrecognized setting "// &
"#define CVMIX_TIDAL_SCHEME "//trim(CVMix_tidal_scheme_str)//" found in input file.")
end select
endif ! CS%use_CVMix_tidal
! Read in vertical profile of tidal energy dissipation
if ( CS%CVMix_tidal_scheme == SCHMITTNER .or. .not. CS%use_CVMix_tidal) then
call get_param(param_file, mdl, "INT_TIDE_PROFILE", int_tide_profile_str, &
"INT_TIDE_PROFILE selects the vertical profile of energy "//&
"dissipation with INT_TIDE_DISSIPATION. Valid values are:\n"//&
"\t STLAURENT_02 - Use the St. Laurent et al exponential \n"//&
"\t decay profile.\n"//&
"\t POLZIN_09 - Use the Polzin WKB-stretched algebraic \n"//&
"\t decay profile.", &
default=STLAURENT_PROFILE_STRING)
int_tide_profile_str = uppercase(int_tide_profile_str)
select case (int_tide_profile_str)
case (STLAURENT_PROFILE_STRING) ; CS%int_tide_profile = STLAURENT_02
case (POLZIN_PROFILE_STRING) ; CS%int_tide_profile = POLZIN_09
case default
call MOM_error(FATAL, "tidal_mixing_init: Unrecognized setting "// &
"#define INT_TIDE_PROFILE "//trim(int_tide_profile_str)//" found in input file.")
end select
endif
elseif (CS%use_CVMix_tidal) then
call MOM_error(FATAL, "tidal_mixing_init: Cannot set INT_TIDE_DISSIPATION to False "// &
"when USE_CVMix_TIDAL is set to True.")
endif
call get_param(param_file, mdl, "LEE_WAVE_DISSIPATION", CS%Lee_wave_dissipation, &
"If true, use an lee wave driven dissipation scheme to "//&
"drive diapycnal mixing, along the lines of Nikurashin "//&
"(2010) and using the St. Laurent et al. (2002) "//&
"and Simmons et al. (2004) vertical profile", default=.false.)
if (CS%lee_wave_dissipation) then
if (CS%use_CVMix_tidal) then
call MOM_error(FATAL, "tidal_mixing_init: Lee wave driven dissipation scheme cannot "// &
"be used when CVMix tidal mixing scheme is active.")
endif
call get_param(param_file, mdl, "LEE_WAVE_PROFILE", tmpstr, &
"LEE_WAVE_PROFILE selects the vertical profile of energy "//&
"dissipation with LEE_WAVE_DISSIPATION. Valid values are:\n"//&
"\t STLAURENT_02 - Use the St. Laurent et al exponential \n"//&
"\t decay profile.\n"//&
"\t POLZIN_09 - Use the Polzin WKB-stretched algebraic \n"//&
"\t decay profile.", &
default=STLAURENT_PROFILE_STRING)
tmpstr = uppercase(tmpstr)
select case (tmpstr)
case (STLAURENT_PROFILE_STRING) ; CS%lee_wave_profile = STLAURENT_02
case (POLZIN_PROFILE_STRING) ; CS%lee_wave_profile = POLZIN_09
case default
call MOM_error(FATAL, "tidal_mixing_init: Unrecognized setting "// &
"#define LEE_WAVE_PROFILE "//trim(tmpstr)//" found in input file.")
end select
endif
call get_param(param_file, mdl, "INT_TIDE_LOWMODE_DISSIPATION", CS%Lowmode_itidal_dissipation, &
"If true, consider mixing due to breaking low modes that "//&
"have been remotely generated; as with itidal drag on the "//&
"barotropic tide, use an internal tidal dissipation scheme to "//&
"drive diapycnal mixing, along the lines of St. Laurent "//&
"et al. (2002) and Simmons et al. (2004).", default=.false.)
if ((CS%Int_tide_dissipation .and. (CS%int_tide_profile == POLZIN_09)) .or. &
(CS%lee_wave_dissipation .and. (CS%lee_wave_profile == POLZIN_09))) then
if (CS%use_CVMix_tidal) then
call MOM_error(FATAL, "tidal_mixing_init: Polzin scheme cannot "// &
"be used when CVMix tidal mixing scheme is active.")
endif
call get_param(param_file, mdl, "NU_POLZIN", CS%Nu_Polzin, &
"When the Polzin decay profile is used, this is a "//&
"non-dimensional constant in the expression for the "//&
"vertical scale of decay for the tidal energy dissipation.", &
units="nondim", default=0.0697)
call get_param(param_file, mdl, "NBOTREF_POLZIN", CS%Nbotref_Polzin, &
"When the Polzin decay profile is used, this is the "//&
"reference value of the buoyancy frequency at the ocean "//&
"bottom in the Polzin formulation for the vertical "//&
"scale of decay for the tidal energy dissipation.", &
units="s-1", default=9.61e-4, scale=US%T_to_s)
call get_param(param_file, mdl, "POLZIN_DECAY_SCALE_FACTOR", &
CS%Polzin_decay_scale_factor, &
"When the Polzin decay profile is used, this is a "//&
"scale factor for the vertical scale of decay of the tidal "//&
"energy dissipation.", default=1.0, units="nondim")
call get_param(param_file, mdl, "POLZIN_SCALE_MAX_FACTOR", &
CS%Polzin_decay_scale_max_factor, &
"When the Polzin decay profile is used, this is a factor "//&
"to limit the vertical scale of decay of the tidal "//&
"energy dissipation to POLZIN_DECAY_SCALE_MAX_FACTOR "//&
"times the depth of the ocean.", units="nondim", default=1.0)
call get_param(param_file, mdl, "POLZIN_MIN_DECAY_SCALE", CS%Polzin_min_decay_scale, &
"When the Polzin decay profile is used, this is the "//&
"minimum vertical decay scale for the vertical profile\n"//&
"of internal tide dissipation with the Polzin (2009) formulation", &
units="m", default=0.0, scale=US%m_to_Z)
endif
if (CS%Int_tide_dissipation .or. CS%Lee_wave_dissipation) then
call get_param(param_file, mdl, "INT_TIDE_DECAY_SCALE", CS%Int_tide_decay_scale, &
"The decay scale away from the bottom for tidal TKE with "//&
"the new coding when INT_TIDE_DISSIPATION is used.", &
units="m", default=500.0, scale=US%m_to_Z)
call get_param(param_file, mdl, "MU_ITIDES", CS%Mu_itides, &
"A dimensionless turbulent mixing efficiency used with "//&
"INT_TIDE_DISSIPATION, often 0.2.", units="nondim", default=0.2)
call get_param(param_file, mdl, "GAMMA_ITIDES", CS%Gamma_itides, &
"The fraction of the internal tidal energy that is "//&
"dissipated locally with INT_TIDE_DISSIPATION. "//&
"THIS NAME COULD BE BETTER.", &
units="nondim", default=0.3333)
call get_param(param_file, mdl, "MIN_ZBOT_ITIDES", CS%min_zbot_itides, &
"Turn off internal tidal dissipation when the total "//&
"ocean depth is less than this value.", units="m", default=0.0, scale=US%m_to_Z)
endif
if ( (CS%Int_tide_dissipation .or. CS%Lee_wave_dissipation) .and. &
.not. CS%use_CVMix_tidal) then
allocate(CS%Nb(isd:ied,jsd:jed), source=0.)
allocate(CS%h2(isd:ied,jsd:jed), source=0.)
allocate(CS%TKE_itidal(isd:ied,jsd:jed), source=0.)
allocate(CS%mask_itidal(isd:ied,jsd:jed), source=1.)
call get_param(param_file, mdl, "KAPPA_ITIDES", CS%kappa_itides, &
"A topographic wavenumber used with INT_TIDE_DISSIPATION. "//&
"The default is 2pi/10 km, as in St.Laurent et al. 2002.", &
units="m-1", default=8.e-4*atan(1.0), scale=US%Z_to_m)
call get_param(param_file, mdl, "UTIDE", CS%utide, &
"The constant tidal amplitude used with INT_TIDE_DISSIPATION.", &
units="m s-1", default=0.0, scale=US%m_to_Z*US%T_to_s)
allocate(CS%tideamp(is:ie,js:je), source=CS%utide)
call get_param(param_file, mdl, "KAPPA_H2_FACTOR", CS%kappa_h2_factor, &
"A scaling factor for the roughness amplitude with "//&
"INT_TIDE_DISSIPATION.", units="nondim", default=1.0)
call get_param(param_file, mdl, "TKE_ITIDE_MAX", CS%TKE_itide_max, &
"The maximum internal tide energy source available to mix "//&
"above the bottom boundary layer with INT_TIDE_DISSIPATION.", &
units="W m-2", default=1.0e3, scale=US%W_m2_to_RZ3_T3)
call get_param(param_file, mdl, "READ_TIDEAMP", read_tideamp, &
"If true, read a file (given by TIDEAMP_FILE) containing "//&
"the tidal amplitude with INT_TIDE_DISSIPATION.", default=.false.)
if (read_tideamp) then
if (CS%use_CVMix_tidal) then
call MOM_error(FATAL, "tidal_mixing_init: Tidal amplitude files are "// &
"not compatible with CVMix tidal mixing. ")
endif
call get_param(param_file, mdl, "TIDEAMP_FILE", tideamp_file, &
"The path to the file containing the spatially varying "//&
"tidal amplitudes with INT_TIDE_DISSIPATION.", default="tideamp.nc")
filename = trim(CS%inputdir) // trim(tideamp_file)
call log_param(param_file, mdl, "INPUTDIR/TIDEAMP_FILE", filename)
call get_param(param_file, mdl, "TIDEAMP_VARNAME", tideamp_var, &
"The name of the tidal amplitude variable in the input file.", &
default="tideamp")
! NOTE: There are certain cases where FMS is unable to read this file, so
! we use read_netCDF_data in place of MOM_read_data.
call read_netCDF_data(filename, tideamp_var, CS%tideamp, G%domain, &
rescale=US%m_to_Z*US%T_to_s)
endif
call get_param(param_file, mdl, "H2_FILE", h2_file, &
"The path to the file containing the sub-grid-scale "//&
"topographic roughness amplitude with INT_TIDE_DISSIPATION.", &
fail_if_missing=(.not.CS%use_CVMix_tidal))
filename = trim(CS%inputdir) // trim(h2_file)
call log_param(param_file, mdl, "INPUTDIR/H2_FILE", filename)
call get_param(param_file, mdl, "ROUGHNESS_VARNAME", rough_var, &
"The name in the input file of the squared sub-grid-scale "//&
"topographic roughness amplitude variable.", default="h2")
! NOTE: There are certain cases where FMS is unable to read this file, so
! we use read_netCDF_data in place of MOM_read_data.
call read_netCDF_data(filename, rough_var, CS%h2, G%domain, &
rescale=US%m_to_Z**2)
call get_param(param_file, mdl, "FRACTIONAL_ROUGHNESS_MAX", max_frac_rough, &
"The maximum topographic roughness amplitude as a fraction of the mean depth, "//&
"or a negative value for no limitations on roughness.", &
units="nondim", default=0.1)
do j=js,je ; do i=is,ie
if (G%bathyT(i,j)+G%Z_ref < CS%min_zbot_itides) CS%mask_itidal(i,j) = 0.0
CS%tideamp(i,j) = CS%tideamp(i,j) * CS%mask_itidal(i,j) * G%mask2dT(i,j)
! Restrict rms topo to a fraction (often 10 percent) of the column depth.
if ((CS%tidal_answer_date < 20190101) .and. (max_frac_rough >= 0.0)) then
hamp = min(max_frac_rough*(G%bathyT(i,j)+G%Z_ref), sqrt(CS%h2(i,j)))
CS%h2(i,j) = hamp*hamp
else
if (max_frac_rough >= 0.0) &
CS%h2(i,j) = min((max_frac_rough*(G%bathyT(i,j)+G%Z_ref))**2, CS%h2(i,j))
endif
utide = CS%tideamp(i,j)
! Compute the fixed part of internal tidal forcing.
! The units here are [R Z3 T-2 ~> J m-2 = kg s-2] here.
CS%TKE_itidal(i,j) = 0.5 * CS%kappa_h2_factor * GV%Rho0 * &
CS%kappa_itides * CS%h2(i,j) * utide*utide
enddo ; enddo
endif
if (CS%Lee_wave_dissipation) then
call get_param(param_file, mdl, "NIKURASHIN_TKE_INPUT_FILE", Niku_TKE_input_file, &
"The path to the file containing the TKE input from lee "//&
"wave driven mixing. Used with LEE_WAVE_DISSIPATION.", &
fail_if_missing=.true.)
call get_param(param_file, mdl, "NIKURASHIN_SCALE", Niku_scale, &
"A non-dimensional factor by which to scale the lee-wave "//&
"driven TKE input. Used with LEE_WAVE_DISSIPATION.", &
units="nondim", default=1.0)
filename = trim(CS%inputdir) // trim(Niku_TKE_input_file)
call log_param(param_file, mdl, "INPUTDIR/NIKURASHIN_TKE_INPUT_FILE", filename)
call get_param(param_file, mdl, "TKE_INPUT_VAR", TKE_input_var, &
"The name in the input file of the turbulent kinetic energy input variable.", &
default="TKE_input")
allocate(CS%TKE_Niku(is:ie,js:je), source=0.)
call MOM_read_data(filename, TKE_input_var, CS%TKE_Niku, G%domain, timelevel=1, & ! ??? timelevel -aja
scale=Niku_scale*US%W_m2_to_RZ3_T3)
call get_param(param_file, mdl, "GAMMA_NIKURASHIN",CS%Gamma_lee, &
"The fraction of the lee wave energy that is dissipated "//&
"locally with LEE_WAVE_DISSIPATION.", units="nondim", default=0.3333)
call get_param(param_file, mdl, "DECAY_SCALE_FACTOR_LEE",CS%Decay_scale_factor_lee, &
"Scaling for the vertical decay scale of the local "//&
"dissipation of lee wave dissipation.", units="nondim", default=1.0)
else
CS%Decay_scale_factor_lee = -9.e99 ! This should never be used if CS%Lee_wave_dissipation = False
endif
! Configure CVMix
if (CS%use_CVMix_tidal) then
! Read in CVMix params
!call openParameterBlock(param_file,'CVMix_TIDAL')
call get_param(param_file, mdl, "TIDAL_MAX_COEF", CS%tidal_max_coef, &
"largest acceptable value for tidal diffusivity", &
units="m^2/s", default=50e-4, scale=US%m2_s_to_Z2_T) ! the default is 50e-4 in CVMix, 100e-4 in POP.
call get_param(param_file, mdl, "TIDAL_DISS_LIM_TC", CS%tidal_diss_lim_tc, &
"Min allowable depth for dissipation for tidal-energy-constituent data. "//&
"No dissipation contribution is applied above TIDAL_DISS_LIM_TC.", &
units="m", default=0.0, scale=US%m_to_Z)
call get_param(param_file, mdl, 'MIN_THICKNESS', CS%min_thickness, &
units="m", default=0.001, scale=US%m_to_Z, do_not_log=.True.)
call get_param(param_file, mdl, "PRANDTL_TIDAL", prandtl_tidal, &
"Prandtl number used by CVMix tidal mixing schemes "//&
"to convert vertical diffusivities into viscosities.", &
units="nondim", default=1.0, do_not_log=.true.)
call CVMix_put(CS%CVMix_glb_params, 'Prandtl', prandtl_tidal)
call get_param(param_file, mdl, "TIDAL_ENERGY_TYPE",tidal_energy_type, &
"The type of input tidal energy flux dataset. Valid values are"//&
"\t Jayne\n"//&
"\t ER03 \n",&
fail_if_missing=.true.)
! Check whether tidal energy input format and CVMix tidal mixing scheme are consistent
if ( .not. ( &
(uppercase(tidal_energy_type(1:4)) == 'JAYN' .and. CS%CVMix_tidal_scheme == SIMMONS).or. &
(uppercase(tidal_energy_type(1:4)) == 'ER03' .and. CS%CVMix_tidal_scheme == SCHMITTNER) ) )then
call MOM_error(FATAL, "tidal_mixing_init: Tidal energy file type ("//&
trim(tidal_energy_type)//") is incompatible with CVMix tidal "//&
" mixing scheme: "//trim(CVMix_tidal_scheme_str) )
endif
CVMix_tidal_scheme_str = lowercase(CVMix_tidal_scheme_str)
! Set up CVMix
call CVMix_init_tidal(CVmix_tidal_params_user = CS%CVMix_tidal_params, &
mix_scheme = CVMix_tidal_scheme_str, &
efficiency = CS%Mu_itides, &
vertical_decay_scale = CS%int_tide_decay_scale*US%Z_to_m, &
max_coefficient = CS%tidal_max_coef*US%Z2_T_to_m2_s, &
local_mixing_frac = CS%Gamma_itides, &
depth_cutoff = CS%min_zbot_itides*US%Z_to_m)
call read_tidal_energy(G, US, tidal_energy_type, param_file, CS)
!call closeParameterBlock(param_file)
endif ! CVMix on
! Register Diagnostics fields
if (CS%Int_tide_dissipation .or. CS%Lee_wave_dissipation .or. &
CS%Lowmode_itidal_dissipation) then
CS%id_Kd_itidal = register_diag_field('ocean_model','Kd_itides',diag%axesTi,Time, &
'Internal Tide Driven Diffusivity', 'm2 s-1', conversion=US%Z2_T_to_m2_s)
if (CS%use_CVMix_tidal) then
CS%id_N2_int = register_diag_field('ocean_model','N2_int',diag%axesTi,Time, &
'Bouyancy frequency squared, at interfaces', 's-2', conversion=US%s_to_T**2)
!> TODO: add units
if (CS%CVMix_tidal_scheme .eq. SIMMONS) then
CS%id_Simmons_coeff = register_diag_field('ocean_model','Simmons_coeff',diag%axesT1,Time, &
'time-invariant portion of the tidal mixing coefficient using the Simmons', '')
else if (CS%CVMix_tidal_scheme .eq. SCHMITTNER) then
CS%id_Schmittner_coeff = register_diag_field('ocean_model','Schmittner_coeff',diag%axesTL,Time, &
'time-invariant portion of the tidal mixing coefficient using the Schmittner', '')
CS%id_tidal_qe_md = register_diag_field('ocean_model','tidal_qe_md',diag%axesTL,Time, &
'input tidal energy dissipated locally interpolated to model vertical coordinates', &
'W m-2', conversion=US%RZ3_T3_to_W_m2)
endif
CS%id_vert_dep = register_diag_field('ocean_model','vert_dep',diag%axesTi,Time, &
'vertical deposition function needed for Simmons et al tidal mixing', '')
else
CS%id_TKE_itidal = register_diag_field('ocean_model','TKE_itidal',diag%axesT1,Time, &
'Internal Tide Driven Turbulent Kinetic Energy', &
'W m-2', conversion=US%RZ3_T3_to_W_m2)
CS%id_Nb = register_diag_field('ocean_model','Nb',diag%axesT1,Time, &
'Bottom Buoyancy Frequency', 's-1', conversion=US%s_to_T)
CS%id_Kd_lowmode = register_diag_field('ocean_model','Kd_lowmode',diag%axesTi,Time, &
'Internal Tide Driven Diffusivity (from propagating low modes)', &
'm2 s-1', conversion=US%Z2_T_to_m2_s)
CS%id_Fl_itidal = register_diag_field('ocean_model','Fl_itides',diag%axesTi,Time, &
'Vertical flux of tidal turbulent dissipation', &
'm3 s-3', conversion=(US%Z_to_m**3*US%s_to_T**3))
CS%id_Fl_lowmode = register_diag_field('ocean_model','Fl_lowmode',diag%axesTi,Time, &
'Vertical flux of tidal turbulent dissipation (from propagating low modes)', &
'm3 s-3', conversion=(US%Z_to_m**3*US%s_to_T**3))
CS%id_Polzin_decay_scale = register_diag_field('ocean_model','Polzin_decay_scale',diag%axesT1,Time, &
'Vertical decay scale for the tidal turbulent dissipation with Polzin scheme', &
'm', conversion=US%Z_to_m)
CS%id_Polzin_decay_scale_scaled = register_diag_field('ocean_model', &
'Polzin_decay_scale_scaled', diag%axesT1, Time, &
'Vertical decay scale for the tidal turbulent dissipation with Polzin scheme, '// &
'scaled by N2_bot/N2_meanz', 'm', conversion=US%Z_to_m)
CS%id_N2_bot = register_diag_field('ocean_model','N2_b',diag%axesT1,Time, &
'Bottom Buoyancy frequency squared', 's-2', conversion=US%s_to_T**2)
CS%id_N2_meanz = register_diag_field('ocean_model','N2_meanz', diag%axesT1, Time, &
'Buoyancy frequency squared averaged over the water column', 's-2', conversion=US%s_to_T**2)
CS%id_Kd_Itidal_Work = register_diag_field('ocean_model','Kd_Itidal_Work',diag%axesTL,Time, &
'Work done by Internal Tide Diapycnal Mixing', &
'W m-2', conversion=US%RZ3_T3_to_W_m2)
CS%id_Kd_Niku_Work = register_diag_field('ocean_model','Kd_Nikurashin_Work',diag%axesTL,Time, &
'Work done by Nikurashin Lee Wave Drag Scheme', &
'W m-2', conversion=US%RZ3_T3_to_W_m2)
CS%id_Kd_Lowmode_Work = register_diag_field('ocean_model','Kd_Lowmode_Work',diag%axesTL,Time, &
'Work done by Internal Tide Diapycnal Mixing (low modes)', &
'W m-2', conversion=US%RZ3_T3_to_W_m2)
if (CS%Lee_wave_dissipation) then
CS%id_TKE_leewave = register_diag_field('ocean_model','TKE_leewave',diag%axesT1,Time, &
'Lee wave Driven Turbulent Kinetic Energy', &
'W m-2', conversion=US%RZ3_T3_to_W_m2)
CS%id_Kd_Niku = register_diag_field('ocean_model','Kd_Nikurashin',diag%axesTi,Time, &
'Lee Wave Driven Diffusivity', 'm2 s-1', conversion=US%Z2_T_to_m2_s)
endif
endif ! S%use_CVMix_tidal
endif
end function tidal_mixing_init
!> Depending on whether or not CVMix is active, calls the associated subroutine to compute internal
!! tidal dissipation and to add the effect of internal-tide-driven mixing to the layer or interface
!! diffusivities.
subroutine calculate_tidal_mixing(h, j, N2_bot, N2_lay, N2_int, TKE_to_Kd, max_TKE, &
G, GV, US, CS, Kd_max, Kv, Kd_lay, Kd_int)
type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure
type(verticalGrid_type), intent(in) :: GV !< The ocean's 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 thicknesses [H ~> m or kg m-2]
integer, intent(in) :: j !< The j-index to work on
real, dimension(SZI_(G)), intent(in) :: N2_bot !< The near-bottom squared buoyancy
!! frequency [T-2 ~> s-2].
real, dimension(SZI_(G),SZK_(GV)), intent(in) :: N2_lay !< The squared buoyancy frequency of the
!! layers [T-2 ~> s-2].
real, dimension(SZI_(G),SZK_(GV)+1), intent(in) :: N2_int !< The squared buoyancy frequency at the
!! interfaces [T-2 ~> s-2].
real, dimension(SZI_(G),SZK_(GV)), intent(in) :: TKE_to_Kd !< The conversion rate between the TKE
!! dissipated within a layer and the
!! diapycnal diffusivity within that layer,
!! usually (~Rho_0 / (G_Earth * dRho_lay))
!! [Z2 T-1 / Z3 T-3 = T2 Z-1 ~> s2 m-1]
real, dimension(SZI_(G),SZK_(GV)), intent(in) :: max_TKE !< The energy required to for a layer to entrain
!! to its maximum realizable thickness [Z3 T-3 ~> m3 s-3]
type(tidal_mixing_cs), intent(inout) :: CS !< The control structure for this module
real, intent(in) :: Kd_max !< The maximum increment for diapycnal
!! diffusivity due to TKE-based processes,
!! [Z2 T-1 ~> m2 s-1].
!! Set this to a negative value to have no limit.
real, dimension(:,:,:), pointer :: Kv !< The "slow" vertical viscosity at each interface
!! (not layer!) [Z2 T-1 ~> m2 s-1].
real, dimension(SZI_(G),SZK_(GV)), &
optional, intent(inout) :: Kd_lay !< The diapycnal diffusivity in layers [Z2 T-1 ~> m2 s-1].
real, dimension(SZI_(G),SZK_(GV)+1), &
optional, intent(inout) :: Kd_int !< The diapycnal diffusivity at interfaces,
!! [Z2 T-1 ~> m2 s-1].
if (CS%Int_tide_dissipation .or. CS%Lee_wave_dissipation .or. CS%Lowmode_itidal_dissipation) then
if (CS%use_CVMix_tidal) then
call calculate_CVMix_tidal(h, j, N2_int, G, GV, US, CS, Kv, Kd_lay, Kd_int)
else
call add_int_tide_diffusivity(h, j, N2_bot, N2_lay, TKE_to_Kd, max_TKE, &
G, GV, US, CS, Kd_max, Kd_lay, Kd_int)
endif
endif
end subroutine calculate_tidal_mixing
!> Calls the CVMix routines to compute tidal dissipation and to add the effect of internal-tide-driven
!! mixing to the interface diffusivities.
subroutine calculate_CVMix_tidal(h, j, N2_int, G, GV, US, CS, Kv, Kd_lay, Kd_int)
type(ocean_grid_type), intent(in) :: G !< Grid structure.
type(verticalGrid_type), intent(in) :: GV !< ocean vertical grid structure
type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
type(tidal_mixing_cs), intent(inout) :: CS !< This module's control structure.
real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2].
integer, intent(in) :: j !< The j-index to work on
real, dimension(SZI_(G),SZK_(GV)+1), intent(in) :: N2_int !< The squared buoyancy
!! frequency at the interfaces [T-2 ~> s-2].
real, dimension(:,:,:), pointer :: Kv !< The "slow" vertical viscosity at each interface
!! (not layer!) [Z2 T-1 ~> m2 s-1].
real, dimension(SZI_(G),SZK_(GV)), &
optional, intent(inout) :: Kd_lay!< The diapycnal diffusivity in the layers [Z2 T-1 ~> m2 s-1].
real, dimension(SZI_(G),SZK_(GV)+1), &
optional, intent(inout) :: Kd_int!< The diapycnal diffusivity at interfaces [Z2 T-1 ~> m2 s-1].
! Local variables
real, dimension(SZK_(GV)+1) :: Kd_tidal ! tidal diffusivity [m2 s-1]
real, dimension(SZK_(GV)+1) :: Kv_tidal ! tidal viscosity [m2 s-1]
real, dimension(SZK_(GV)+1) :: vert_dep ! vertical deposition [nondim]
real, dimension(SZK_(GV)+1) :: iFaceHeight ! Height of interfaces [m]
real, dimension(SZK_(GV)+1) :: SchmittnerSocn ! A larger value of the Schmittner coefficint to
! use in the Southern Ocean [nondim]. If this is smaller
! than Schmittner_coeff, that standard value is used.
real, dimension(SZK_(GV)) :: cellHeight ! Height of cell centers [m]
real, dimension(SZK_(GV)) :: tidal_qe_md ! Tidal dissipation energy interpolated from 3d input
! to model coordinates [R Z3 T-3 ~> W m-2]
real, dimension(SZK_(GV)+1) :: N2_int_i ! De-scaled interface buoyancy frequency [s-2]
real, dimension(SZK_(GV)) :: Schmittner_coeff ! A coefficient in the Schmittner et al (2014) mixing
! parameterization [nondim]
real, dimension(SZK_(GV)) :: h_m ! Cell thickness [m]
real, allocatable, dimension(:,:) :: exp_hab_zetar ! A badly documented array that appears to be
! related to the distribution of tidal mixing energy, with unusual array
! extents that are not explained, that is set and used by the CVMix
! tidal mixing schemes, perhaps in [m3 kg-1]?
real :: dh, hcorr ! Limited thicknesses and a cumulative correction [Z ~> m]
real :: Simmons_coeff ! A coefficient in the Simmons et al (2004) mixing parameterization [nondim]
integer :: i, k, is, ie
real, parameter :: rho_fw = 1000.0 ! fresh water density [kg m-3]
! TODO: when coupled, get this from CESM (SHR_CONST_RHOFW)
is = G%isc ; ie = G%iec
select case (CS%CVMix_tidal_scheme)
case (SIMMONS)
do i=is,ie
if (G%mask2dT(i,j)<1) cycle
iFaceHeight = 0.0 ! BBL is all relative to the surface
hcorr = 0.0
! Compute cell center depth and cell bottom in meters (negative values in the ocean)
do k=1,GV%ke
dh = h(i,j,k) * GV%H_to_Z ! Nominal thickness to use for increment, in the units of heights
dh = dh + hcorr ! Take away the accumulated error (could temporarily make dh<0)
hcorr = min( dh - CS%min_thickness, 0. ) ! If inflating then hcorr<0
dh = max(dh, CS%min_thickness) ! Limited increment dh>=min_thickness
cellHeight(k) = iFaceHeight(k) - 0.5 * US%Z_to_m*dh
iFaceHeight(k+1) = iFaceHeight(k) - US%Z_to_m*dh
enddo
call CVMix_compute_Simmons_invariant( nlev = GV%ke, &
energy_flux = US%RZ3_T3_to_W_m2*CS%tidal_qe_2d(i,j), &
rho = rho_fw, &
SimmonsCoeff = Simmons_coeff, &
VertDep = vert_dep, &
zw = iFaceHeight, &
zt = cellHeight, &
CVMix_tidal_params_user = CS%CVMix_tidal_params)
! Since we pass tidal_qe_2d=(CS%Gamma_itides)*tidal_energy_flux_2d, and not tidal_energy_flux_2d in
! above subroutine call, we divide Simmons_coeff by CS%Gamma_itides as a corrective step:
! TODO: (CS%Gamma_itides)*tidal_energy_flux_2d is unnecessary, directly use tidal_energy_flux_2d
Simmons_coeff = Simmons_coeff / CS%Gamma_itides
! XXX: Temporary de-scaling of N2_int(i,:) into a temporary variable
do K=1,GV%ke+1
N2_int_i(K) = US%s_to_T**2 * N2_int(i,K)
enddo
call CVMix_coeffs_tidal( Mdiff_out = Kv_tidal, &
Tdiff_out = Kd_tidal, &
Nsqr = N2_int_i, &
OceanDepth = -iFaceHeight(GV%ke+1),&
SimmonsCoeff = Simmons_coeff, &
vert_dep = vert_dep, &
nlev = GV%ke, &
max_nlev = GV%ke, &
CVMix_params = CS%CVMix_glb_params, &
CVMix_tidal_params_user = CS%CVMix_tidal_params)
! Update diffusivity
if (present(Kd_lay)) then
do k=1,GV%ke
Kd_lay(i,k) = Kd_lay(i,k) + 0.5 * US%m2_s_to_Z2_T * (Kd_tidal(k) + Kd_tidal(k+1))
enddo
endif
if (present(Kd_int)) then
do K=1,GV%ke+1
Kd_int(i,K) = Kd_int(i,K) + (US%m2_s_to_Z2_T * Kd_tidal(K))
enddo
endif
! Update viscosity with the proper unit conversion.
if (associated(Kv)) then
do K=1,GV%ke+1
Kv(i,j,K) = Kv(i,j,K) + US%m2_s_to_Z2_T * Kv_tidal(K) ! Rescale from m2 s-1 to Z2 T-1.
enddo
endif
! diagnostics
if (allocated(CS%dd%Kd_itidal)) then
CS%dd%Kd_itidal(i,j,:) = US%m2_s_to_Z2_T * Kd_tidal(:)
endif
if (allocated(CS%dd%N2_int)) then
CS%dd%N2_int(i,j,:) = N2_int(i,:)
endif
if (allocated(CS%dd%Simmons_coeff_2d)) then
CS%dd%Simmons_coeff_2d(i,j) = Simmons_coeff
endif
if (allocated(CS%dd%vert_dep_3d)) then
CS%dd%vert_dep_3d(i,j,:) = vert_dep(:)
endif
enddo ! i=is,ie
case (SCHMITTNER)
! TODO: correct exp_hab_zetar shapes in CVMix_compute_Schmittner_invariant
! and CVMix_compute_SchmittnerCoeff low subroutines
allocate(exp_hab_zetar(GV%ke+1,GV%ke+1))
do i=is,ie
if (G%mask2dT(i,j)<1) cycle
iFaceHeight(:) = 0.0 ! BBL is all relative to the surface
hcorr = 0.0
! Compute heights at cell center and interfaces, and rescale layer thicknesses
do k=1,GV%ke
h_m(k) = h(i,j,k)*GV%H_to_m ! Rescale thicknesses to m for use by CVmix.
dh = h(i,j,k) * GV%H_to_Z ! Nominal thickness to use for increment, in the units of heights
dh = dh + hcorr ! Take away the accumulated error (could temporarily make dh<0)
hcorr = min( dh - CS%min_thickness, 0. ) ! If inflating then hcorr<0
dh = max(dh, CS%min_thickness) ! Limited increment dh>=min_thickness
cellHeight(k) = iFaceHeight(k) - 0.5 * US%Z_to_m*dh
iFaceHeight(k+1) = iFaceHeight(k) - US%Z_to_m*dh
enddo
SchmittnerSocn = 0.0 ! TODO: compute this
! form the time-invariant part of Schmittner coefficient term
call CVMix_compute_Schmittner_invariant(nlev = GV%ke, &
VertDep = vert_dep, &
efficiency = CS%Mu_itides, &
rho = rho_fw, &
exp_hab_zetar = exp_hab_zetar, &
zw = iFaceHeight, &
CVmix_tidal_params_user = CS%CVMix_tidal_params)
!TODO: in above call, there is no need to pass efficiency, since it gets
! passed via CVMix_init_tidal and stored in CVMix_tidal_params. Change
! CVMix API to prevent this redundancy.
! remap from input z coordinate to model coordinate:
tidal_qe_md(:) = 0.0
call remapping_core_h(CS%remap_cs, size(CS%h_src), CS%h_src, CS%tidal_qe_3d_in(i,j,:), &
GV%ke, h_m, tidal_qe_md, GV%H_subroundoff, GV%H_subroundoff)
! form the Schmittner coefficient that is based on 3D q*E, which is formed from
! summing q_i*TidalConstituent_i over the number of constituents.
call CVMix_compute_SchmittnerCoeff( nlev = GV%ke, &
energy_flux = US%RZ3_T3_to_W_m2*tidal_qe_md(:), &
SchmittnerCoeff = Schmittner_coeff, &
exp_hab_zetar = exp_hab_zetar, &
CVmix_tidal_params_user = CS%CVMix_tidal_params)
! XXX: Temporary de-scaling of N2_int(i,:) into a temporary variable
do k=1,GV%ke+1
N2_int_i(k) = US%s_to_T**2 * N2_int(i,k)
enddo
call CVMix_coeffs_tidal_schmittner( Mdiff_out = Kv_tidal, &
Tdiff_out = Kd_tidal, &
Nsqr = N2_int_i, &
OceanDepth = -iFaceHeight(GV%ke+1), &
nlev = GV%ke, &
max_nlev = GV%ke, &
SchmittnerCoeff = Schmittner_coeff, &
SchmittnerSouthernOcean = SchmittnerSocn, &
CVmix_params = CS%CVMix_glb_params, &
CVmix_tidal_params_user = CS%CVMix_tidal_params)
! Update diffusivity
if (present(Kd_lay)) then
do k=1,GV%ke
Kd_lay(i,k) = Kd_lay(i,k) + 0.5 * US%m2_s_to_Z2_T * (Kd_tidal(k) + Kd_tidal(k+1))
enddo
endif
if (present(Kd_int)) then
do K=1,GV%ke+1
Kd_int(i,K) = Kd_int(i,K) + (US%m2_s_to_Z2_T * Kd_tidal(K))
enddo
endif
! Update viscosity
if (associated(Kv)) then
do K=1,GV%ke+1
Kv(i,j,K) = Kv(i,j,K) + US%m2_s_to_Z2_T * Kv_tidal(K) ! Rescale from m2 s-1 to Z2 T-1.
enddo
endif
! diagnostics
if (allocated(CS%dd%Kd_itidal)) then
CS%dd%Kd_itidal(i,j,:) = US%m2_s_to_Z2_T*Kd_tidal(:)
endif
if (allocated(CS%dd%N2_int)) then
CS%dd%N2_int(i,j,:) = N2_int(i,:)
endif
if (allocated(CS%dd%Schmittner_coeff_3d)) then
CS%dd%Schmittner_coeff_3d(i,j,:) = Schmittner_coeff(:)
endif
if (allocated(CS%dd%tidal_qe_md)) then
CS%dd%tidal_qe_md(i,j,:) = tidal_qe_md(:)
endif
if (allocated(CS%dd%vert_dep_3d)) then
CS%dd%vert_dep_3d(i,j,:) = vert_dep(:)
endif
enddo ! i=is,ie
deallocate(exp_hab_zetar)