forked from mom-ocean/MOM6
-
Notifications
You must be signed in to change notification settings - Fork 62
/
Copy pathMOM_set_viscosity.F90
2568 lines (2344 loc) · 131 KB
/
MOM_set_viscosity.F90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
!> Calculates various values related to the bottom boundary layer, such as the viscosity and
!! thickness of the BBL (set_viscous_BBL).
module MOM_set_visc
! This file is part of MOM6. See LICENSE.md for the license.
use MOM_ALE, only : ALE_CS, ALE_remap_velocities, ALE_remap_interface_vals, ALE_remap_vertex_vals
use MOM_cpu_clock, only : cpu_clock_id, cpu_clock_begin, cpu_clock_end, CLOCK_ROUTINE
use MOM_cvmix_conv, only : cvmix_conv_is_used
use MOM_CVMix_ddiff, only : CVMix_ddiff_is_used
use MOM_cvmix_shear, only : cvmix_shear_is_used
use MOM_debugging, only : uvchksum, hchksum
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 : pass_var, CORNER
use MOM_EOS, only : calculate_density, calculate_density_derivs, calculate_specific_vol_derivs
use MOM_error_handler, only : MOM_error, FATAL, WARNING
use MOM_file_parser, only : get_param, log_param, log_version, param_file_type
use MOM_forcing_type, only : forcing, mech_forcing, find_ustar
use MOM_grid, only : ocean_grid_type
use MOM_hor_index, only : hor_index_type
use MOM_interface_heights, only : thickness_to_dz
use MOM_io, only : slasher, MOM_read_data
use MOM_kappa_shear, only : kappa_shear_is_used, kappa_shear_at_vertex
use MOM_open_boundary, only : ocean_OBC_type, OBC_segment_type, OBC_NONE, OBC_DIRECTION_E
use MOM_open_boundary, only : OBC_DIRECTION_W, OBC_DIRECTION_N, OBC_DIRECTION_S
use MOM_restart, only : register_restart_field, query_initialized, MOM_restart_CS
use MOM_restart, only : register_restart_field_as_obsolete
use MOM_safe_alloc, only : safe_alloc_ptr, safe_alloc_alloc
use MOM_unit_scaling, only : unit_scale_type
use MOM_variables, only : thermo_var_ptrs, vertvisc_type, porous_barrier_type
use MOM_verticalGrid, only : verticalGrid_type
implicit none ; private
#include <MOM_memory.h>
public set_viscous_BBL, set_viscous_ML, set_visc_init, set_visc_end
public set_visc_register_restarts, remap_vertvisc_aux_vars
! A note on unit descriptions in comments: MOM6 uses units that can be rescaled for dimensional
! consistency testing. These are noted in comments with units like Z, H, L, and T, along with
! their mks counterparts with notation like "a velocity [Z T-1 ~> m s-1]". If the units
! vary with the Boussinesq approximation, the Boussinesq variant is given first.
!> Control structure for MOM_set_visc
type, public :: set_visc_CS ; private
logical :: initialized = .false. !< True if this control structure has been initialized.
real :: Hbbl !< The static bottom boundary layer thickness [H ~> m or kg m-2].
!! Runtime parameter `HBBL`.
real :: dz_bbl !< The static bottom boundary layer thickness in height units [Z ~> m].
!! Runtime parameter `HBBL`.
real :: cdrag !< The quadratic drag coefficient [nondim].
!! Runtime parameter `CDRAG`.
real :: c_Smag !< The Laplacian Smagorinsky coefficient for
!! calculating the drag in channels [nondim].
real :: drag_bg_vel !< An assumed unresolved background velocity for
!! calculating the bottom drag [L T-1 ~> m s-1].
!! Runtime parameter `DRAG_BG_VEL`.
real :: BBL_thick_min !< The minimum bottom boundary layer thickness [Z ~> m].
!! This might be Kv / (cdrag * drag_bg_vel) to give
!! Kv as the minimum near-bottom viscosity.
real :: Htbl_shelf !< A nominal thickness of the surface boundary layer for use
!! in calculating the near-surface velocity [H ~> m or kg m-2].
real :: Htbl_shelf_min !< The minimum surface boundary layer thickness [Z ~> m].
real :: KV_BBL_min !< The minimum viscosity in the bottom boundary layer [H Z T-1 ~> m2 s-1 or Pa s]
real :: KV_TBL_min !< The minimum viscosity in the top boundary layer [H Z T-1 ~> m2 s-1 or Pa s]
logical :: bottomdraglaw !< If true, the bottom stress is calculated with a
!! drag law c_drag*|u|*u. The velocity magnitude
!! may be an assumed value or it may be based on the
!! actual velocity in the bottommost `HBBL`, depending
!! on whether linear_drag is true.
!! Runtime parameter `BOTTOMDRAGLAW`.
logical :: body_force_drag !< If true, the bottom stress is imposed as an explicit body force
!! applied over a fixed distance from the bottom, rather than as an
!! implicit calculation based on an enhanced near-bottom viscosity.
logical :: BBL_use_EOS !< If true, use the equation of state in determining
!! the properties of the bottom boundary layer.
logical :: linear_drag !< If true, the drag law is cdrag*`DRAG_BG_VEL`*u.
!! Runtime parameter `LINEAR_DRAG`.
logical :: Channel_drag !< If true, the drag is exerted directly on each layer
!! according to what fraction of the bottom they overlie.
real :: Chan_drag_max_vol !< The maximum bottom boundary layer volume within which the
!! channel drag is applied, normalized by the full cell area,
!! or a negative value to apply no maximum [Z ~> m].
logical :: correct_BBL_bounds !< If true, uses the correct bounds on the BBL thickness and
!! viscosity so that the bottom layer feels the intended drag.
logical :: RiNo_mix !< If true, use Richardson number dependent mixing.
logical :: dynamic_viscous_ML !< If true, use a bulk Richardson number criterion to
!! determine the mixed layer thickness for viscosity.
real :: bulk_Ri_ML !< The bulk mixed layer used to determine the
!! thickness of the viscous mixed layer [nondim]
real :: omega !< The Earth's rotation rate [T-1 ~> s-1].
real :: ustar_min !< A minimum value of ustar to avoid numerical
!! problems [H T-1 ~> m s-1 or kg m-2 s-1]. If the value is
!! small enough, this should not affect the solution.
real :: TKE_decay !< The ratio of the natural Ekman depth to the TKE
!! decay scale [nondim]
real :: omega_frac !< When setting the decay scale for turbulence, use this
!! fraction of the absolute rotation rate blended with the local
!! value of f, as sqrt((1-of)*f^2 + of*4*omega^2) [nondim]
integer :: answer_date !< The vintage of the order of arithmetic and expressions in the set
!! viscosity 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.
logical :: debug !< If true, write verbose checksums for debugging purposes.
logical :: BBL_use_tidal_bg !< If true, use a tidal background amplitude for the bottom velocity
!! when computing the bottom stress.
character(len=200) :: inputdir !< The directory for input files.
type(ocean_OBC_type), pointer :: OBC => NULL() !< Open boundaries control structure
type(diag_ctrl), pointer :: diag => NULL() !< A structure that is used to
!! regulate the timing of diagnostic output.
! Allocatable data arrays
real, allocatable, dimension(:,:) :: tideamp !< RMS tidal amplitude at h points [Z T-1 ~> m s-1]
! Diagnostic arrays
real, allocatable, dimension(:,:) :: bbl_u !< BBL mean U current [L T-1 ~> m s-1]
real, allocatable, dimension(:,:) :: bbl_v !< BBL mean V current [L T-1 ~> m s-1]
!>@{ Diagnostics handles
integer :: id_bbl_thick_u = -1, id_kv_bbl_u = -1, id_bbl_u = -1
integer :: id_bbl_thick_v = -1, id_kv_bbl_v = -1, id_bbl_v = -1
integer :: id_Ray_u = -1, id_Ray_v = -1
integer :: id_nkml_visc_u = -1, id_nkml_visc_v = -1
!>@}
end type set_visc_CS
contains
!> Calculates the thickness of the bottom boundary layer and the viscosity within that layer.
subroutine set_viscous_BBL(u, v, h, tv, visc, G, GV, US, CS, pbv)
type(ocean_grid_type), intent(inout) :: 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(SZIB_(G),SZJ_(G),SZK_(GV)), &
intent(in) :: u !< The zonal velocity [L T-1 ~> m s-1].
real, dimension(SZI_(G),SZJB_(G),SZK_(GV)), &
intent(in) :: v !< The meridional velocity [L T-1 ~> m s-1].
real, dimension(SZI_(G),SZJ_(G),SZK_(GV)), &
intent(in) :: h !< Layer thicknesses [H ~> m or kg m-2].
type(thermo_var_ptrs), intent(in) :: tv !< A structure containing pointers to any
!! available thermodynamic fields. Absent fields
!! have NULL pointers.
type(vertvisc_type), intent(inout) :: visc !< A structure containing vertical viscosities and
!! related fields.
type(set_visc_CS), intent(inout) :: CS !< The control structure returned by a previous
!! call to set_visc_init.
type(porous_barrier_type),intent(in) :: pbv !< porous barrier fractional cell metrics
! Local variables
real, dimension(SZIB_(G)) :: &
ustar, & ! The bottom friction velocity [H T-1 ~> m s-1 or kg m-2 s-1].
T_EOS, & ! The temperature used to calculate the partial derivatives
! of density with T and S [C ~> degC].
S_EOS, & ! The salinity used to calculate the partial derivatives
! of density with T and S [S ~> ppt].
dR_dT, & ! Partial derivative of the density in the bottom boundary
! layer with temperature [R C-1 ~> kg m-3 degC-1].
dR_dS, & ! Partial derivative of the density in the bottom boundary
! layer with salinity [R S-1 ~> kg m-3 ppt-1].
press, & ! The pressure at which dR_dT and dR_dS are evaluated [R L2 T-2 ~> Pa].
umag_avg, & ! The average magnitude of velocities in the bottom boundary layer [L T-1 ~> m s-1].
h_bbl_drag, & ! The thickness over which to apply drag as a body force [H ~> m or kg m-2].
dz_bbl_drag ! The vertical height over which to apply drag as a body force [Z ~> m].
real :: htot ! Sum of the layer thicknesses up to some point [H ~> m or kg m-2].
real :: dztot ! Distance from the bottom up to some point [Z ~> m].
real :: htot_vel ! Sum of the layer thicknesses up to some point [H ~> m or kg m-2].
real :: dztot_vel ! Distance from the bottom up to some point [Z ~> m].
real :: Rhtot ! Running sum of thicknesses times the layer potential
! densities [H R ~> kg m-2 or kg2 m-5].
real, dimension(SZIB_(G),SZJ_(G)) :: &
D_u, & ! Bottom depth linearly interpolated to u points [Z ~> m].
mask_u ! A mask that disables any contributions from u points that
! are land or past open boundary conditions [nondim], 0 or 1.
real, dimension(SZI_(G),SZJB_(G)) :: &
D_v, & ! Bottom depth linearly interpolated to v points [Z ~> m].
mask_v ! A mask that disables any contributions from v points that
! are land or past open boundary conditions [nondim], 0 or 1.
real, dimension(SZIB_(G),SZK_(GV)) :: &
h_at_vel, & ! Layer thickness at a velocity point, using an upwind-biased
! second order accurate estimate based on the previous velocity
! direction [H ~> m or kg m-2].
h_vel, & ! Arithmetic mean of the layer thicknesses adjacent to a
! velocity point [H ~> m or kg m-2].
dz_at_vel, & ! Vertical extent of a layer, using an upwind-biased
! second order accurate estimate based on the previous velocity
! direction [Z ~> m].
dz_vel, & ! Arithmetic mean of the difference in across the layers adjacent
! to a velocity point [Z ~> m].
T_vel, & ! Arithmetic mean of the layer temperatures adjacent to a
! velocity point [C ~> degC].
S_vel, & ! Arithmetic mean of the layer salinities adjacent to a
! velocity point [S ~> ppt].
SpV_vel, & ! Arithmetic mean of the layer averaged specific volumes adjacent to a
! velocity point [R-1 ~> kg m-3].
Rml_vel ! Arithmetic mean of the layer coordinate densities adjacent
! to a velocity point [R ~> kg m-3].
real :: dz(SZI_(G),SZJ_(G),SZK_(GV)) ! Height change across layers [Z ~> m]
real :: h_vel_pos ! The arithmetic mean thickness at a velocity point
! plus H_neglect to avoid 0 values [H ~> m or kg m-2].
real :: ustarsq ! 400 times the square of ustar, times
! Rho0 divided by G_Earth and the conversion
! from m to thickness units [H R ~> kg m-2 or kg2 m-5].
real :: cdrag_sqrt ! Square root of the drag coefficient [nondim].
real :: cdrag_sqrt_H ! Square root of the drag coefficient, times a unit conversion factor
! from lateral lengths to layer thicknesses [H L-1 ~> nondim or kg m-3].
real :: cdrag_sqrt_H_RL ! Square root of the drag coefficient, times a unit conversion factor from
! density times lateral lengths to layer thicknesses [H L-1 R-1 ~> m3 kg-1 or nondim]
real :: cdrag_L_to_H ! The drag coeffient times conversion factors from lateral
! distance to thickness units [H L-1 ~> nondim or kg m-3]
real :: cdrag_RL_to_H ! The drag coeffient times conversion factors from density times lateral
! distance to thickness units [H L-1 R-1 ~> m3 kg-1 or nondim]
real :: cdrag_conv ! The drag coeffient times a combination of static conversion factors and in
! situ density or Boussinesq reference density [H L-1 ~> nondim or kg m-3]
real :: oldfn ! The integrated energy required to
! entrain up to the bottom of the layer,
! divided by G_Earth [H R ~> kg m-2 or kg2 m-5].
real :: Dfn ! The increment in oldfn for entraining
! the layer [H R ~> kg m-2 or kg2 m-5].
real :: frac_used ! The fraction of the present layer that contributes to Dh and Ddz [nondim]
real :: Dh ! The increment in layer thickness from
! the present layer [H ~> m or kg m-2].
real :: Ddz ! The increment in height change from the present layer [Z ~> m].
real :: bbl_thick ! The thickness of the bottom boundary layer [Z ~> m].
real :: BBL_thick_max ! A huge upper bound on the boundary layer thickness [Z ~> m].
real :: kv_bbl ! The bottom boundary layer viscosity [H Z T-1 ~> m2 s-1 or Pa s]
real :: C2f ! C2f = 2*f at velocity points [T-1 ~> s-1].
real :: U_bg_sq ! The square of an assumed background
! velocity, for calculating the mean
! magnitude near the bottom for use in the
! quadratic bottom drag [L2 T-2 ~> m2 s-2].
real :: hwtot ! Sum of the thicknesses used to calculate
! the near-bottom velocity magnitude [H ~> m or kg m-2].
real :: I_hwtot ! The Adcroft reciprocal of hwtot [H-1 ~> m-1 or m2 kg-1].
real :: dzwtot ! The vertical extent of the region used to calculate
! the near-bottom velocity magnitude [Z ~> m].
real :: hutot ! Running sum of thicknesses times the velocity
! magnitudes [H L T-1 ~> m2 s-1 or kg m-1 s-1].
real :: Thtot ! Running sum of thickness times temperature [C H ~> degC m or degC kg m-2].
real :: Shtot ! Running sum of thickness times salinity [S H ~> ppt m or ppt kg m-2].
real :: SpV_htot ! Running sum of thickness times specific volume [R-1 H ~> m4 kg-1 or m]
real :: hweight ! The thickness of a layer that is within Hbbl
! of the bottom [H ~> m or kg m-2].
real :: dzweight ! The counterpart of hweight in height units [Z ~> m].
real :: v_at_u, u_at_v ! v at a u point or vice versa [L T-1 ~> m s-1].
real :: Rho0x400_G ! 400*Rho0/G_Earth, times unit conversion factors
! [R T2 H-1 ~> kg s2 m-4 or s2 m-1].
! The 400 is a constant proposed by Killworth and Edwards, 1999.
real, dimension(SZI_(G),SZJ_(G),max(GV%nk_rho_varies,1)) :: &
Rml ! The mixed layer coordinate density [R ~> kg m-3].
real :: p_ref(SZI_(G)) ! The pressure used to calculate the coordinate
! density [R L2 T-2 ~> Pa] (usually set to 2e7 Pa = 2000 dbar).
real :: D_vel ! The bottom depth at a velocity point [Z ~> m].
real :: Dp, Dm ! The depths at the edges of a velocity cell [Z ~> m].
real :: crv ! crv is the curvature of the bottom depth across a
! cell, times the cell width squared [Z ~> m].
real :: crv_3 ! crv/3 [Z ~> m].
real :: C24_crv ! 24/crv [Z-1 ~> m-1].
real :: slope ! The absolute value of the bottom depth slope across
! a cell times the cell width [Z ~> m].
real :: apb_4a, ax2_3apb ! Various nondimensional ratios of crv and slope [nondim].
real :: a2x48_apb3, Iapb, Ibma_2 ! Combinations of crv (a) and slope (b) [Z-1 ~> m-1]
! All of the following "volumes" have units of vertical heights because they are normalized
! by the full horizontal area of a velocity cell.
real :: Vol_bbl_chan ! The volume of the bottom boundary layer as used in the channel
! drag parameterization, normalized by the full horizontal area
! of the velocity cell [Z ~> m].
real :: Vol_open ! The cell volume above which it is open [Z ~> m].
real :: Vol_direct ! With less than Vol_direct [Z ~> m], there is a direct
! solution of a cubic equation for L.
real :: Vol_2_reg ! The cell volume above which there are two separate
! open areas that must be integrated [Z ~> m].
real :: vol ! The volume below the interface whose normalized
! width is being sought [Z ~> m].
real :: vol_below ! The volume below the interface below the one that
! is currently under consideration [Z ~> m].
real :: Vol_err ! The error in the volume with the latest estimate of
! L, or the error for the interface below [Z ~> m].
real :: Vol_quit ! The volume error below which to quit iterating [Z ~> m].
real :: Vol_tol ! A volume error tolerance [Z ~> m].
real :: L(SZK_(GV)+1) ! The fraction of the full cell width that is open at
! the depth of each interface [nondim].
real :: L_direct ! The value of L above volume Vol_direct [nondim].
real :: L_max, L_min ! Upper and lower bounds on the correct value for L [nondim].
real :: Vol_err_max ! The volume error for the upper bound on the correct value for L [Z ~> m]
real :: Vol_err_min ! The volume error for the lower bound on the correct value for L [Z ~> m]
real :: Vol_0 ! A deeper volume with known width L0 [Z ~> m].
real :: L0 ! The value of L above volume Vol_0 [nondim].
real :: dVol ! vol - Vol_0 [Z ~> m].
real :: dV_dL2 ! The partial derivative of volume with L squared
! evaluated at L=L0 [Z ~> m].
real :: h_neglect ! A thickness that is so small it is usually lost
! in roundoff and can be neglected [H ~> m or kg m-2].
real :: dz_neglect ! A vertical distance that is so small it is usually lost
! in roundoff and can be neglected [Z ~> m].
real :: ustH ! ustar converted to units of H T-1 [H T-1 ~> m s-1 or kg m-2 s-1].
real :: root ! A temporary variable [H T-1 ~> m s-1 or kg m-2 s-1].
real :: Cell_width ! The transverse width of the velocity cell [L ~> m].
real :: Rayleigh ! A factor that is multiplied by the layer's velocity magnitude
! to give the Rayleigh drag velocity, times a lateral distance to
! thickness conversion factor [H L-1 ~> nondim or kg m-3].
real :: gam ! The ratio of the change in the open interface width
! to the open interface width atop a cell [nondim].
real :: BBL_frac ! The fraction of a layer's drag that goes into the
! viscous bottom boundary layer [nondim].
real :: BBL_visc_frac ! The fraction of all the drag that is expressed as
! a viscous bottom boundary layer [nondim].
real :: h_bbl_fr ! The fraction of the bottom boundary layer in a layer [nondim].
real :: h_sum ! The sum of the thicknesses of the layers below the one being
! worked on [H ~> m or kg m-2].
real, parameter :: C1_3 = 1.0/3.0, C1_6 = 1.0/6.0, C1_12 = 1.0/12.0 ! Rational constants [nondim]
real :: C2pi_3 ! An irrational constant, 2/3 pi. [nondim]
real :: tmp ! A temporary variable, sometimes in [Z ~> m]
real :: tmp_val_m1_to_p1 ! A temporary variable [nondim]
real :: curv_tol ! Numerator of curvature cubed, used to estimate
! accuracy of a single L(:) Newton iteration [Z5 ~> m5]
logical :: use_L0, do_one_L_iter ! Control flags for L(:) Newton iteration
logical :: use_BBL_EOS, do_i(SZIB_(G))
integer, dimension(2) :: EOSdom ! The computational domain for the equation of state
integer :: i, j, k, is, ie, js, je, Isq, Ieq, Jsq, Jeq, nz, m, n, K2, nkmb, nkml
integer :: itt, maxitt=20
type(ocean_OBC_type), pointer :: OBC => NULL()
is = G%isc ; ie = G%iec ; js = G%jsc ; je = G%jec ; nz = GV%ke
Isq = G%isc-1 ; Ieq = G%IecB ; Jsq = G%jsc-1 ; Jeq = G%JecB
nkmb = GV%nk_rho_varies ; nkml = GV%nkml
h_neglect = GV%H_subroundoff
dz_neglect = GV%dZ_subroundoff
Rho0x400_G = 400.0*(GV%H_to_RZ / (US%L_to_Z**2 * GV%g_Earth))
Vol_quit = (0.9*GV%Angstrom_Z + dz_neglect)
C2pi_3 = 8.0*atan(1.0)/3.0
if (.not.CS%initialized) call MOM_error(FATAL,"MOM_set_viscosity(BBL): "//&
"Module must be initialized before it is used.")
if (.not.CS%bottomdraglaw) return
if (CS%debug) then
call uvchksum("Start set_viscous_BBL [uv]", u, v, G%HI, haloshift=1, scale=US%L_T_to_m_s)
call hchksum(h,"Start set_viscous_BBL h", G%HI, haloshift=1, scale=GV%H_to_m)
if (associated(tv%T)) call hchksum(tv%T, "Start set_viscous_BBL T", G%HI, haloshift=1, scale=US%C_to_degC)
if (associated(tv%S)) call hchksum(tv%S, "Start set_viscous_BBL S", G%HI, haloshift=1, scale=US%S_to_ppt)
if (allocated(tv%SpV_avg)) &
call hchksum(tv%SpV_avg, "Start set_viscous_BBL SpV_avg", G%HI, haloshift=1, scale=US%kg_m3_to_R)
if (allocated(tv%SpV_avg)) call hchksum(tv%SpV_avg, "Cornerless SpV_avg", G%HI, &
haloshift=1, omit_corners=.true., scale=US%kg_m3_to_R)
if (associated(tv%T)) call hchksum(tv%T, "Cornerless T", G%HI, haloshift=1, omit_corners=.true., scale=US%C_to_degC)
if (associated(tv%S)) call hchksum(tv%S, "Cornerless S", G%HI, haloshift=1, omit_corners=.true., scale=US%S_to_ppt)
endif
use_BBL_EOS = associated(tv%eqn_of_state) .and. CS%BBL_use_EOS
OBC => CS%OBC
U_bg_sq = CS%drag_bg_vel * CS%drag_bg_vel
cdrag_sqrt = sqrt(CS%cdrag)
cdrag_sqrt_H = cdrag_sqrt * US%L_to_m * GV%m_to_H
cdrag_sqrt_H_RL = cdrag_sqrt * US%L_to_Z * GV%RZ_to_H
cdrag_L_to_H = CS%cdrag * US%L_to_m * GV%m_to_H
cdrag_RL_to_H = CS%cdrag * US%L_to_Z * GV%RZ_to_H
BBL_thick_max = G%Rad_Earth_L * US%L_to_Z
K2 = max(nkmb+1, 2)
! Find the vertical distances across layers.
call thickness_to_dz(h, tv, dz, G, GV, US, halo_size=1)
! With a linear drag law, the friction velocity is already known.
! if (CS%linear_drag) ustar(:) = cdrag_sqrt_H*CS%drag_bg_vel
if ((nkml>0) .and. .not.use_BBL_EOS) then
EOSdom(1) = Isq - (G%isd-1) ; EOSdom(2) = G%iec+1 - (G%isd-1)
do i=Isq,Ieq+1 ; p_ref(i) = tv%P_Ref ; enddo
!$OMP parallel do default(shared)
do k=1,nkmb ; do j=Jsq,Jeq+1
call calculate_density(tv%T(:,j,k), tv%S(:,j,k), p_ref, Rml(:,j,k), tv%eqn_of_state, &
EOSdom)
enddo ; enddo
endif
!$OMP parallel do default(shared)
do J=js-1,je ; do i=is-1,ie+1
D_v(i,J) = 0.5*(G%bathyT(i,j) + G%bathyT(i,j+1)) + G%Z_ref
mask_v(i,J) = G%mask2dCv(i,J)
enddo ; enddo
!$OMP parallel do default(shared)
do j=js-1,je+1 ; do I=is-1,ie
D_u(I,j) = 0.5*(G%bathyT(i,j) + G%bathyT(i+1,j)) + G%Z_ref
mask_u(I,j) = G%mask2dCu(I,j)
enddo ; enddo
if (associated(OBC)) then ; do n=1,OBC%number_of_segments
if (.not. OBC%segment(n)%on_pe) cycle
! Use a one-sided projection of bottom depths at OBC points.
I = OBC%segment(n)%HI%IsdB ; J = OBC%segment(n)%HI%JsdB
if (OBC%segment(n)%is_N_or_S .and. (J >= js-1) .and. (J <= je)) then
do i = max(is-1,OBC%segment(n)%HI%isd), min(ie+1,OBC%segment(n)%HI%ied)
if (OBC%segment(n)%direction == OBC_DIRECTION_N) D_v(i,J) = G%bathyT(i,j) + G%Z_ref
if (OBC%segment(n)%direction == OBC_DIRECTION_S) D_v(i,J) = G%bathyT(i,j+1) + G%Z_ref
enddo
elseif (OBC%segment(n)%is_E_or_W .and. (I >= is-1) .and. (I <= ie)) then
do j = max(js-1,OBC%segment(n)%HI%jsd), min(je+1,OBC%segment(n)%HI%jed)
if (OBC%segment(n)%direction == OBC_DIRECTION_E) D_u(I,j) = G%bathyT(i,j) + G%Z_ref
if (OBC%segment(n)%direction == OBC_DIRECTION_W) D_u(I,j) = G%bathyT(i+1,j) + G%Z_ref
enddo
endif
enddo ; endif
if (associated(OBC)) then ; do n=1,OBC%number_of_segments
! Now project bottom depths across cell-corner points in the OBCs. The two
! projections have to occur in sequence and can not be combined easily.
if (.not. OBC%segment(n)%on_pe) cycle
! Use a one-sided projection of bottom depths at OBC points.
I = OBC%segment(n)%HI%IsdB ; J = OBC%segment(n)%HI%JsdB
if (OBC%segment(n)%is_N_or_S .and. (J >= js-1) .and. (J <= je)) then
do I = max(is-1,OBC%segment(n)%HI%IsdB), min(ie,OBC%segment(n)%HI%IedB)
if (OBC%segment(n)%direction == OBC_DIRECTION_N) then
D_u(I,j+1) = D_u(I,j) ; mask_u(I,j+1) = 0.0
elseif (OBC%segment(n)%direction == OBC_DIRECTION_S) then
D_u(I,j) = D_u(I,j+1) ; mask_u(I,j) = 0.0
endif
enddo
elseif (OBC%segment(n)%is_E_or_W .and. (I >= is-1) .and. (I <= ie)) then
do J = max(js-1,OBC%segment(n)%HI%JsdB), min(je,OBC%segment(n)%HI%JedB)
if (OBC%segment(n)%direction == OBC_DIRECTION_E) then
D_v(i+1,J) = D_v(i,J) ; mask_v(i+1,J) = 0.0
elseif (OBC%segment(n)%direction == OBC_DIRECTION_W) then
D_v(i,J) = D_v(i+1,J) ; mask_v(i,J) = 0.0
endif
enddo
endif
enddo ; endif
if (.not.use_BBL_EOS) Rml_vel(:,:) = 0.0
if (allocated(visc%Ray_u)) visc%Ray_u(:,:,:) = 0.0
if (allocated(visc%Ray_v)) visc%Ray_v(:,:,:) = 0.0
!$OMP parallel do default(private) shared(u,v,h,dz,tv,visc,G,GV,US,CS,Rml,nz,nkmb,nkml,K2, &
!$OMP Isq,Ieq,Jsq,Jeq,h_neglect,dz_neglect,Rho0x400_G,C2pi_3, &
!$OMP U_bg_sq,cdrag_sqrt,cdrag_sqrt_H,cdrag_sqrt_H_RL, &
!$OMP cdrag_L_to_H,cdrag_RL_to_H,use_BBL_EOS,BBL_thick_max, &
!$OMP OBC,maxitt,D_u,D_v,mask_u,mask_v,pbv) &
!$OMP firstprivate(Vol_quit)
do j=Jsq,Jeq ; do m=1,2
if (m==1) then
! m=1 refers to u-points
if (j<G%Jsc) cycle
is = Isq ; ie = Ieq
do i=is,ie
do_i(i) = (G%mask2dCu(I,j) > 0.0)
enddo
else
! m=2 refers to v-points
is = G%isc ; ie = G%iec
do i=is,ie
do_i(i) = (G%mask2dCv(i,J) > 0.0)
enddo
endif
! Calculate thickness at velocity points (u or v depending on value of m).
! Also interpolate the ML density or T/S properties.
if (m==1) then ! u-points
do k=1,nz ; do I=is,ie
if (do_i(I)) then
if (u(I,j,k) * (h(i+1,j,k) - h(i,j,k)) >= 0) then
! If the flow is from thin to thick then bias towards the thinner thickness
h_at_vel(I,k) = 2.0*h(i,j,k)*h(i+1,j,k) / &
(h(i,j,k) + h(i+1,j,k) + h_neglect)
dz_at_vel(I,k) = 2.0*dz(i,j,k)*dz(i+1,j,k) / &
(dz(i,j,k) + dz(i+1,j,k) + dz_neglect)
else
! If the flow is from thick to thin then use the simple average thickness
h_at_vel(I,k) = 0.5 * (h(i,j,k) + h(i+1,j,k))
dz_at_vel(I,k) = 0.5 * (dz(i,j,k) + dz(i+1,j,k))
endif
endif
h_vel(I,k) = 0.5 * (h(i,j,k) + h(i+1,j,k))
dz_vel(I,k) = 0.5 * (dz(i,j,k) + dz(i+1,j,k))
enddo ; enddo
if (use_BBL_EOS) then ; do k=1,nz ; do I=is,ie
! Perhaps these should be thickness weighted.
T_vel(I,k) = 0.5 * (tv%T(i,j,k) + tv%T(i+1,j,k))
S_vel(I,k) = 0.5 * (tv%S(i,j,k) + tv%S(i+1,j,k))
enddo ; enddo ; else ; do k=1,nkmb ; do I=is,ie
Rml_vel(I,k) = 0.5 * (Rml(i,j,k) + Rml(i+1,j,k))
enddo ; enddo ; endif
if (allocated(tv%SpV_avg)) then ; do k=1,nz ; do I=is,ie
SpV_vel(I,k) = 0.5 * (tv%SpV_avg(i,j,k) + tv%SpV_avg(i+1,j,k))
enddo ; enddo ; endif
else ! v-points
do k=1,nz ; do i=is,ie
if (do_i(i)) then
if (v(i,J,k) * (h(i,j+1,k) - h(i,j,k)) >= 0) then
! If the flow is from thin to thick then bias towards the thinner thickness
h_at_vel(i,k) = 2.0*h(i,j,k)*h(i,j+1,k) / &
(h(i,j,k) + h(i,j+1,k) + h_neglect)
dz_at_vel(i,k) = 2.0*dz(i,j,k)*dz(i,j+1,k) / &
(dz(i,j,k) + dz(i,j+1,k) + dz_neglect)
else
! If the flow is from thick to thin then use the simple average thickness
h_at_vel(i,k) = 0.5 * (h(i,j,k) + h(i,j+1,k))
dz_at_vel(i,k) = 0.5 * (dz(i,j,k) + dz(i,j+1,k))
endif
endif
h_vel(i,k) = 0.5 * (h(i,j,k) + h(i,j+1,k))
dz_vel(i,k) = 0.5 * (dz(i,j,k) + dz(i,j+1,k))
enddo ; enddo
if (use_BBL_EOS) then ; do k=1,nz ; do i=is,ie
! Perhaps these should be thickness weighted.
T_vel(i,k) = 0.5 * (tv%T(i,j,k) + tv%T(i,j+1,k))
S_vel(i,k) = 0.5 * (tv%S(i,j,k) + tv%S(i,j+1,k))
enddo ; enddo ; else ; do k=1,nkmb ; do i=is,ie
Rml_vel(i,k) = 0.5 * (Rml(i,j,k) + Rml(i,j+1,k))
enddo ; enddo ; endif
if (allocated(tv%SpV_avg)) then ; do k=1,nz ; do i=is,ie
SpV_vel(i,k) = 0.5 * (tv%SpV_avg(i,j,k) + tv%SpV_avg(i,j+1,k))
enddo ; enddo ; endif
endif
if (associated(OBC)) then ; if (OBC%number_of_segments > 0) then
! Apply a zero gradient projection of thickness across OBC points.
if (m==1) then
do I=is,ie ; if (do_i(I) .and. (OBC%segnum_u(I,j) /= OBC_NONE)) then
if (OBC%segment(OBC%segnum_u(I,j))%direction == OBC_DIRECTION_E) then
do k=1,nz
h_at_vel(I,k) = h(i,j,k) ; h_vel(I,k) = h(i,j,k)
dz_at_vel(I,k) = dz(i,j,k) ; dz_vel(I,k) = dz(i,j,k)
enddo
if (use_BBL_EOS) then
do k=1,nz
T_vel(I,k) = tv%T(i,j,k) ; S_vel(I,k) = tv%S(i,j,k)
enddo
else
do k=1,nkmb
Rml_vel(I,k) = Rml(i,j,k)
enddo
endif
if (allocated(tv%SpV_avg)) then ; do k=1,nz
SpV_vel(I,k) = tv%SpV_avg(i,j,k)
enddo ; endif
elseif (OBC%segment(OBC%segnum_u(I,j))%direction == OBC_DIRECTION_W) then
do k=1,nz
h_at_vel(I,k) = h(i+1,j,k) ; h_vel(I,k) = h(i+1,j,k)
dz_at_vel(I,k) = dz(i+1,j,k) ; dz_vel(I,k) = dz(i+1,j,k)
enddo
if (use_BBL_EOS) then
do k=1,nz
T_vel(I,k) = tv%T(i+1,j,k) ; S_vel(I,k) = tv%S(i+1,j,k)
enddo
else
do k=1,nkmb
Rml_vel(I,k) = Rml(i+1,j,k)
enddo
endif
if (allocated(tv%SpV_avg)) then ; do k=1,nz
SpV_vel(I,k) = tv%SpV_avg(i+1,j,k)
enddo ; endif
endif
endif ; enddo
else
do i=is,ie ; if (do_i(i) .and. (OBC%segnum_v(i,J) /= OBC_NONE)) then
if (OBC%segment(OBC%segnum_v(i,J))%direction == OBC_DIRECTION_N) then
do k=1,nz
h_at_vel(i,k) = h(i,j,k) ; h_vel(i,k) = h(i,j,k)
dz_at_vel(i,k) = dz(i,j,k) ; dz_vel(i,k) = dz(i,j,k)
enddo
if (use_BBL_EOS) then
do k=1,nz
T_vel(i,k) = tv%T(i,j,k) ; S_vel(i,k) = tv%S(i,j,k)
enddo
else
do k=1,nkmb
Rml_vel(i,k) = Rml(i,j,k)
enddo
endif
if (allocated(tv%SpV_avg)) then ; do k=1,nz
SpV_vel(i,k) = tv%SpV_avg(i,j,k)
enddo ; endif
elseif (OBC%segment(OBC%segnum_v(i,J))%direction == OBC_DIRECTION_S) then
do k=1,nz
h_at_vel(i,k) = h(i,j+1,k) ; h_vel(i,k) = h(i,j+1,k)
dz_at_vel(i,k) = dz(i,j+1,k) ; dz_vel(i,k) = dz(i,j+1,k)
enddo
if (use_BBL_EOS) then
do k=1,nz
T_vel(i,k) = tv%T(i,j+1,k) ; S_vel(i,k) = tv%S(i,j+1,k)
enddo
else
do k=1,nkmb
Rml_vel(i,k) = Rml(i,j+1,k)
enddo
endif
if (allocated(tv%SpV_avg)) then ; do k=1,nz
SpV_vel(i,k) = tv%SpV_avg(i,j+1,k)
enddo ; endif
endif
endif ; enddo
endif
endif ; endif
if (use_BBL_EOS .or. CS%body_force_drag .or. .not.CS%linear_drag) then
! Calculate the mean velocity magnitude over the bottommost CS%Hbbl of
! the water column for determining the quadratic bottom drag.
! Used in ustar(i)
do i=is,ie ; if (do_i(i)) then
htot_vel = 0.0 ; hwtot = 0.0 ; hutot = 0.0
dztot_vel = 0.0 ; dzwtot = 0.0
Thtot = 0.0 ; Shtot = 0.0 ; SpV_htot = 0.0
do k=nz,1,-1
if (htot_vel>=CS%Hbbl) exit ! terminate the k loop
hweight = MIN(CS%Hbbl - htot_vel, h_at_vel(i,k))
if (hweight < 1.5*GV%Angstrom_H + h_neglect) cycle
dzweight = MIN(CS%dz_bbl - dztot_vel, dz_at_vel(i,k))
htot_vel = htot_vel + h_at_vel(i,k)
hwtot = hwtot + hweight
dztot_vel = dztot_vel + dz_at_vel(i,k)
dzwtot = dzwtot + dzweight
if ((.not.CS%linear_drag) .and. (hweight >= 0.0)) then ; if (m==1) then
v_at_u = set_v_at_u(v, h, G, GV, i, j, k, mask_v, OBC)
if (CS%BBL_use_tidal_bg) then
U_bg_sq = 0.5*( G%mask2dT(i,j)*(CS%tideamp(i,j)*CS%tideamp(i,j))+ &
G%mask2dT(i+1,j)*(CS%tideamp(i+1,j)*CS%tideamp(i+1,j)) )
endif
hutot = hutot + hweight * sqrt(u(I,j,k)*u(I,j,k) + v_at_u*v_at_u + U_bg_sq)
else
u_at_v = set_u_at_v(u, h, G, GV, i, j, k, mask_u, OBC)
if (CS%BBL_use_tidal_bg) then
U_bg_sq = 0.5*( G%mask2dT(i,j)*(CS%tideamp(i,j)*CS%tideamp(i,j))+ &
G%mask2dT(i,j+1)*(CS%tideamp(i,j+1)*CS%tideamp(i,j+1)) )
endif
hutot = hutot + hweight * sqrt(v(i,J,k)*v(i,J,k) + u_at_v*u_at_v + U_bg_sq)
endif ; endif
if (use_BBL_EOS .and. (hweight >= 0.0)) then
Thtot = Thtot + hweight * T_vel(i,k)
Shtot = Shtot + hweight * S_vel(i,k)
endif
if (allocated(tv%SpV_avg) .and. (hweight >= 0.0)) then
SpV_htot = SpV_htot + hweight * SpV_vel(i,k)
endif
enddo ! end of k loop
! Find the Adcroft reciprocal of the total thickness weights
I_hwtot = 0.0 ; if (hwtot > 0.0) I_hwtot = 1.0 / hwtot
! Set u* based on u*^2 = Cdrag u_bbl^2
if ((hwtot <= 0.0) .or. (CS%linear_drag .and. .not.allocated(tv%SpV_avg))) then
ustar(i) = cdrag_sqrt_H * CS%drag_bg_vel
elseif (CS%linear_drag .and. allocated(tv%SpV_avg)) then
ustar(i) = cdrag_sqrt_H_RL * CS%drag_bg_vel * (hwtot / SpV_htot)
elseif (allocated(tv%SpV_avg)) then ! (.not.CS%linear_drag)
ustar(i) = cdrag_sqrt_H_RL * hutot / SpV_htot
else ! (.not.CS%linear_drag .and. .not.allocated(tv%SpV_avg))
ustar(i) = cdrag_sqrt_H * hutot / hwtot
endif
umag_avg(i) = hutot * I_hwtot
h_bbl_drag(i) = hwtot
dz_bbl_drag(i) = dzwtot
if (use_BBL_EOS) then ; if (hwtot > 0.0) then
T_EOS(i) = Thtot/hwtot ; S_EOS(i) = Shtot/hwtot
else
T_EOS(i) = 0.0 ; S_EOS(i) = 0.0
endif ; endif
! Diagnostic BBL flow speed at u- and v-points.
if (CS%id_bbl_u>0 .and. m==1) then
if (hwtot > 0.0) CS%bbl_u(I,j) = hutot/hwtot
elseif (CS%id_bbl_v>0 .and. m==2) then
if (hwtot > 0.0) CS%bbl_v(i,J) = hutot/hwtot
endif
endif ; enddo
else
do i=is,ie ; ustar(i) = cdrag_sqrt_H*CS%drag_bg_vel ; enddo
endif ! Not linear_drag
if (use_BBL_EOS) then
if (associated(tv%p_surf)) then
if (m==1) then ; do i=is,ie ; press(I) = 0.5*(tv%p_surf(i,j) + tv%p_surf(i+1,j)) ; enddo
else ; do i=is,ie ; press(i) = 0.5*(tv%p_surf(i,j) + tv%p_surf(i,j+1)) ; enddo ; endif
else
do i=is,ie ; press(i) = 0.0 ; enddo
endif
do i=is,ie ; if (.not.do_i(i)) then ; T_EOS(i) = 0.0 ; S_EOS(i) = 0.0 ; endif ; enddo
do k=1,nz ; do i=is,ie
press(i) = press(i) + (GV%H_to_RZ*GV%g_Earth) * h_vel(i,k)
enddo ; enddo
call calculate_density_derivs(T_EOS, S_EOS, press, dR_dT, dR_dS, tv%eqn_of_state, &
(/is-G%IsdB+1,ie-G%IsdB+1/) )
endif
! Find a BBL thickness given by equation 2.20 of Killworth and Edwards, 1999:
! ( f h / Cn u* )^2 + ( N h / Ci u* ) = 1
! where Cn=0.5 and Ci=20 (constants suggested by Zilitinkevich and Mironov, 1996).
! Eq. 2.20 can be expressed in terms of boundary layer thicknesses limited by
! rotation (h_f) and stratification (h_N):
! ( h / h_f )^2 + ( h / h_N ) = 1
! When stratification dominates h_N<<h_f, and vice versa.
do i=is,ie ; if (do_i(i)) then
! The 400.0 in this expression is the square of a Ci introduced in KW99, eq. 2.22.
ustarsq = Rho0x400_G * ustar(i)**2 ! Note not in units of u*^2 but [H R ~> kg m-2 or kg2 m-5]
htot = 0.0
dztot = 0.0
! Calculate the thickness of a stratification limited BBL ignoring rotation:
! h_N = Ci u* / N (limit of KW99 eq. 2.20 for |f|->0)
! For layer mode, N^2 = g'/h. Since (Ci u*)^2 = (h_N N)^2 = h_N g' then
! h_N = (Ci u*)^2 / g' (KW99, eq, 2.22)
! Starting from the bottom, integrate the stratification upward until h_N N balances Ci u*
! or in layer mode
! h_N Delta rho ~ (Ci u*)^2 rho0 / g
! where the rhs is stored in variable ustarsq.
! The method was described in Stephens and Hallberg 2000 (unpublished and lost manuscript).
if (use_BBL_EOS) then
Thtot = 0.0 ; Shtot = 0.0 ; oldfn = 0.0
do k=nz,2,-1
if (h_at_vel(i,k) <= 0.0) cycle
! Delta rho * h_bbl assuming everything below is homogenized
oldfn = dR_dT(i)*(Thtot - T_vel(i,k)*htot) + &
dR_dS(i)*(Shtot - S_vel(i,k)*htot)
if (oldfn >= ustarsq) exit
! Local Delta rho * h_bbl at interface
Dfn = (dR_dT(i)*(T_vel(i,k) - T_vel(i,k-1)) + &
dR_dS(i)*(S_vel(i,k) - S_vel(i,k-1))) * &
(h_at_vel(i,k) + htot)
if ((oldfn + Dfn) <= ustarsq) then
! Use whole layer
Dh = h_at_vel(i,k)
Ddz = dz_at_vel(i,k)
else
! Use only part of the layer
frac_used = sqrt((ustarsq-oldfn) / (Dfn))
Dh = h_at_vel(i,k) * frac_used
Ddz = dz_at_vel(i,k) * frac_used
endif
! Increment total BBL thickness and cumulative T and S
htot = htot + Dh
dztot = dztot + Ddz
Thtot = Thtot + T_vel(i,k)*Dh ; Shtot = Shtot + S_vel(i,k)*Dh
enddo
if ((oldfn < ustarsq) .and. h_at_vel(i,1) > 0.0) then
! Layer 1 might be part of the BBL.
if (dR_dT(i) * (Thtot - T_vel(i,1)*htot) + &
dR_dS(i) * (Shtot - S_vel(i,1)*htot) < ustarsq) then
htot = htot + h_at_vel(i,1)
dztot = dztot + dz_at_vel(i,1)
endif
endif ! Examination of layer 1.
else ! Use Rlay and/or the coordinate density as density variables.
Rhtot = 0.0
do k=nz,K2,-1
oldfn = Rhtot - GV%Rlay(k)*htot
Dfn = (GV%Rlay(k) - GV%Rlay(k-1))*(h_at_vel(i,k)+htot)
if (oldfn >= ustarsq) then
cycle
elseif ((oldfn + Dfn) <= ustarsq) then
Dh = h_at_vel(i,k)
Ddz = dz_at_vel(i,k)
else
frac_used = sqrt((ustarsq-oldfn) / (Dfn))
Dh = h_at_vel(i,k) * frac_used
Ddz = dz_at_vel(i,k) * frac_used
endif
htot = htot + Dh
dztot = dztot + Ddz
Rhtot = Rhtot + GV%Rlay(k)*Dh
enddo
if (nkml>0) then
do k=nkmb,2,-1
oldfn = Rhtot - Rml_vel(i,k)*htot
Dfn = (Rml_vel(i,k) - Rml_vel(i,k-1)) * (h_at_vel(i,k)+htot)
if (oldfn >= ustarsq) then
cycle
elseif ((oldfn + Dfn) <= ustarsq) then
Dh = h_at_vel(i,k)
Ddz = dz_at_vel(i,k)
else
frac_used = sqrt((ustarsq-oldfn) / (Dfn))
Dh = h_at_vel(i,k) * frac_used
Ddz = dz_at_vel(i,k) * frac_used
endif
htot = htot + Dh
dztot = dztot + Ddz
Rhtot = Rhtot + Rml_vel(i,k)*Dh
enddo
if (Rhtot - Rml_vel(i,1)*htot < ustarsq) then
htot = htot + h_at_vel(i,1)
dztot = dztot + dz_at_vel(i,1)
endif
else
if (Rhtot - GV%Rlay(1)*htot < ustarsq) then
htot = htot + h_at_vel(i,1)
dztot = dztot + dz_at_vel(i,1)
endif
endif
endif ! use_BBL_EOS
! Value of 2*f at u- or v-points.
if (m==1) then ; C2f = G%CoriolisBu(I,J-1) + G%CoriolisBu(I,J)
else ; C2f = G%CoriolisBu(I-1,J) + G%CoriolisBu(I,J) ; endif
! The thickness of a rotation limited BBL ignoring stratification is
! h_f ~ Cn u* / f (limit of KW99 eq. 2.20 for N->0).
! The buoyancy limit of BBL thickness (h_N) is already in the variable htot from above.
! Substituting x = h_N/h into KW99 eq. 2.20 yields the quadratic
! x^2 - x = (h_N / h_f)^2
! for which the positive root is
! xp = 1/2 + sqrt( 1/4 + (h_N/h_f)^2 )
! and thus h_bbl = h_N / xp . Since h_f = Cn u*/f and Cn=0.5
! xp = 1/2 + sqrt( 1/4 + (2 f h_N/u*)^2 )
! To avoid dividing by zero if u*=0 then
! xp u* = 1/2 u* + sqrt( 1/4 u*^2 + (2 f h_N)^2 )
if (CS%cdrag * U_bg_sq <= 0.0) then
! This avoids NaNs and overflows, and could be used in all cases,
! but is not bitwise identical to the current code.
ustH = ustar(i) ; root = sqrt(0.25*ustH**2 + (htot*C2f)**2)
if (dztot*ustH <= (CS%BBL_thick_min+dz_neglect) * (0.5*ustH + root)) then
bbl_thick = CS%BBL_thick_min
else
! The following expression reads
! h_bbl = h_N u* / ( 1/2 u* + sqrt( 1/4 u*^2 + ( 2 f h_N )^2 ) )
! which is h_bbl = h_N u*/(xp u*) as described above.
bbl_thick = (dztot * ustH) / (0.5*ustH + root)
endif
else
! The following expression reads
! h_bbl = h_N / ( 1/2 + sqrt( 1/4 + ( 2 f h_N / u* )^2 ) )
! which is h_bbl = h_N/xp as described above.
bbl_thick = dztot / (0.5 + sqrt(0.25 + htot*htot*C2f*C2f / (ustar(i)*ustar(i)) ) )
if (bbl_thick < CS%BBL_thick_min) bbl_thick = CS%BBL_thick_min
endif
! Store the normalized bottom boundary layer volume.
if (CS%Channel_drag) Vol_bbl_chan = bbl_thick
! If there is Richardson number dependent mixing, that determines
! the vertical extent of the bottom boundary layer, and there is no
! need to set that scale here. In fact, viscously reducing the
! shears over an excessively large region reduces the efficacy of
! the Richardson number dependent mixing.
! In other words, if using RiNo_mix then CS%dz_bbl acts as an upper bound on
! bbl_thick.
if ((bbl_thick > 0.5*CS%dz_bbl) .and. (CS%RiNo_mix)) bbl_thick = 0.5*CS%dz_bbl
! If drag is a body force, bbl_thick is HBBL
if (CS%body_force_drag) bbl_thick = dz_bbl_drag(i)
if (CS%Channel_drag) then
! The drag within the bottommost Vol_bbl_chan is applied as a part of
! an enhanced bottom viscosity, while above this the drag is applied
! directly to the layers in question as a Rayleigh drag term.
! Restrict the volume over which the channel drag is applied.
if (CS%Chan_drag_max_vol >= 0.0) Vol_bbl_chan = min(Vol_bbl_chan, CS%Chan_drag_max_vol)
!### The harmonic mean edge depths here are not invariant to offsets!
if (m==1) then
D_vel = D_u(I,j)
tmp = G%mask2dCu(I,j+1) * D_u(I,j+1)
Dp = 2.0 * D_vel * tmp / (D_vel + tmp)
tmp = G%mask2dCu(I,j-1) * D_u(I,j-1)
Dm = 2.0 * D_vel * tmp / (D_vel + tmp)
else
D_vel = D_v(i,J)
tmp = G%mask2dCv(i+1,J) * D_v(i+1,J)
Dp = 2.0 * D_vel * tmp / (D_vel + tmp)
tmp = G%mask2dCv(i-1,J) * D_v(i-1,J)
Dm = 2.0 * D_vel * tmp / (D_vel + tmp)
endif
if (Dm > Dp) then ; tmp = Dp ; Dp = Dm ; Dm = tmp ; endif
crv_3 = (Dp + Dm - 2.0*D_vel) ; crv = 3.0*crv_3
slope = Dp - Dm
! If the curvature is small enough, there is no reason not to assume
! a uniformly sloping or flat bottom.
if (abs(crv) < 1e-2*(slope + CS%BBL_thick_min)) crv = 0.0
! Each cell extends from x=-1/2 to 1/2, and has a topography
! given by D(x) = crv*x^2 + slope*x + D - crv/12.
! Calculate the volume above which the entire cell is open and the
! other volumes at which the equation that is solved for L changes.
if (crv > 0.0) then
if (slope >= crv) then
Vol_open = D_vel - Dm ; Vol_2_reg = Vol_open
else
tmp = slope/crv
Vol_open = 0.25*slope*tmp + C1_12*crv
Vol_2_reg = 0.5*tmp**2 * (crv - C1_3*slope)
endif
! Define some combinations of crv & slope for later use.
C24_crv = 24.0/crv ; Iapb = 1.0/(crv+slope)
apb_4a = (slope+crv)/(4.0*crv) ; a2x48_apb3 = (48.0*(crv*crv))*(Iapb**3)
ax2_3apb = 2.0*C1_3*crv*Iapb
elseif (crv == 0.0) then
Vol_open = 0.5*slope
if (slope > 0) Iapb = 1.0/slope
else ! crv < 0.0
Vol_open = D_vel - Dm
if (slope >= -crv) then
Iapb = 1.0e30*US%Z_to_m ; if (slope+crv /= 0.0) Iapb = 1.0/(crv+slope)
Vol_direct = 0.0 ; L_direct = 0.0 ; C24_crv = 0.0
else
C24_crv = 24.0/crv ; Iapb = 1.0/(crv+slope)
L_direct = 1.0 + slope/crv ! L_direct < 1 because crv < 0
Vol_direct = -C1_6*crv*L_direct**3
endif
Ibma_2 = 2.0 / (slope - crv)
endif
L(nz+1) = 0.0 ; vol = 0.0 ; Vol_err = 0.0 ; BBL_visc_frac = 0.0
! Determine the normalized open length at each interface.
do K=nz,1,-1
vol_below = vol
vol = vol + dz_vel(i,k)
h_vel_pos = h_vel(i,k) + h_neglect
if (vol >= Vol_open) then ; L(K) = 1.0
elseif (crv == 0) then ! The bottom has no curvature.
L(K) = sqrt(2.0*vol*Iapb)
elseif (crv > 0) then
! There may be a minimum depth, and there are
! analytic expressions for L for all cases.
if (vol < Vol_2_reg) then
! In this case, there is a contiguous open region and
! vol = 0.5*L^2*(slope + crv/3*(3-4L)).
if (a2x48_apb3*vol < 1e-8) then ! Could be 1e-7?
! There is a very good approximation here for massless layers.
L0 = sqrt(2.0*vol*Iapb) ; L(K) = L0*(1.0 + ax2_3apb*L0)
else
L(K) = apb_4a * (1.0 - &
2.0 * cos(C1_3*acos(a2x48_apb3*vol - 1.0) - C2pi_3))
endif
! To check the answers.
! Vol_err = 0.5*(L(K)*L(K))*(slope + crv_3*(3.0-4.0*L(K))) - vol
else ! There are two separate open regions.
! vol = slope^2/4crv + crv/12 - (crv/12)*(1-L)^2*(1+2L)
! At the deepest volume, L = slope/crv, at the top L = 1.
!L(K) = 0.5 - cos(C1_3*acos(1.0 - C24_crv*(Vol_open - vol)) - C2pi_3)
tmp_val_m1_to_p1 = 1.0 - C24_crv*(Vol_open - vol)
tmp_val_m1_to_p1 = max(-1., min(1., tmp_val_m1_to_p1))
L(K) = 0.5 - cos(C1_3*acos(tmp_val_m1_to_p1) - C2pi_3)
! To check the answers.
! Vol_err = Vol_open - 0.25*crv_3*(1.0+2.0*L(K)) * (1.0-L(K))**2 - vol
endif
else ! a < 0.
if (vol <= Vol_direct) then
! Both edges of the cell are bounded by walls.
L(K) = (-0.25*C24_crv*vol)**C1_3
else
! x_R is at 1/2 but x_L is in the interior & L is found by solving
! vol = 0.5*L^2*(slope + crv/3*(3-4L))
! Vol_err = 0.5*(L(K+1)*L(K+1))*(slope + crv_3*(3.0-4.0*L(K+1))) - vol_below
! Change to ...
! if (min(vol_below + Vol_err, vol) <= Vol_direct) then ?
if (vol_below + Vol_err <= Vol_direct) then
L0 = L_direct ; Vol_0 = Vol_direct
else
L0 = L(K+1) ; Vol_0 = vol_below + Vol_err
! Change to Vol_0 = min(vol_below + Vol_err, vol) ?
endif
! Try a relatively simple solution that usually works well
! for massless layers.
dV_dL2 = 0.5*(slope+crv) - crv*L0 ; dVol = (vol-Vol_0)
! dV_dL2 = 0.5*(slope+crv) - crv*L0 ; dVol = max(vol-Vol_0, 0.0)
use_L0 = .false.
do_one_L_iter = .false.
if (CS%answer_date < 20190101) then
curv_tol = GV%Angstrom_Z*dV_dL2**2 &
* (0.25 * dV_dL2 * GV%Angstrom_Z - crv * L0 * dVol)
do_one_L_iter = (crv * crv * dVol**3) < curv_tol
else
! The following code is more robust when GV%Angstrom_H=0, but
! it changes answers.
use_L0 = (dVol <= 0.)
Vol_tol = max(0.5 * GV%Angstrom_Z + dz_neglect, 1e-14 * vol)
Vol_quit = max(0.9 * GV%Angstrom_Z + dz_neglect, 1e-14 * vol)