-
Notifications
You must be signed in to change notification settings - Fork 24
/
Copy pathvmix_kpp.F90
4206 lines (3426 loc) · 143 KB
/
vmix_kpp.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
!|||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||
module vmix_kpp
!BOP
! !MODULE: vmix_kpp
!
! !DESCRIPTION:
! This module contains routines for initializing and computing
! vertical mixing coefficients for the KPP parameterization
! (see Large, McWilliams and Doney, Reviews of Geophysics, 32, 363
! November 1994).
!
! !REVISION HISTORY:
! SVN:$Id$
! !USES:
use kinds_mod
use blocks
use distribution
use domain_size
use domain
use constants
use grid
use broadcast
use io
use state_mod
use exit_mod
use sw_absorption
use tavg, only: define_tavg_field, accumulate_tavg_field, accumulate_tavg_now
use io_types, only: stdout
use communicate, only: my_task, master_task
use niw_mixing
use tidal_mixing
use registry
use prognostic
use time_management
use cvmix_kinds_and_types
use cvmix_put_get
use cvmix_convection
use cvmix_ddiff
use cvmix_tidal
use cvmix_shear
use cvmix_kpp
use cvmix_math
use shr_sys_mod
use forcing_fields, only: LAMULT, LASL, IFRAC
implicit none
private
save
! !PUBLIC MEMBER FUNCTIONS:
public :: init_vmix_kpp, &
vmix_coeffs_kpp, &
add_kpp_sources, &
smooth_hblt, &
linertial
! !PUBLIC DATA MEMBERS:
real (r8), dimension(:,:,:), allocatable, public :: &
HMXL, &! mixed layer depth
HMXL_DR, &! mixed layer depth with density criterion, QL, 150526
KPP_HBLT, &! boundary layer depth
BOLUS_SP ! scaled eddy-induced (bolus) speed used in inertial
! mixing parameterization
real (r8), public :: &
bckgrnd_vdc2 ! variation in diffusivity
real (r8), public :: cvmix2popL2, pop2cvmixL2
real (r8), public :: cvmix2popL , pop2cvmixL
!EOP
!BOC
!-----------------------------------------------------------------------
!
! mixing constants
!
!-----------------------------------------------------------------------
real (r8) :: &
rich_mix ! coefficient for rich number term
real (r8), dimension(:,:,:,:), allocatable :: &
bckgrnd_vvc, &! background value for viscosity
bckgrnd_vdc ! background value for diffusivity
logical (log_kind) :: &
lrich, &! flag for computing Ri-dependent mixing
ldbl_diff, &! flag for computing double-diffusive mixing
lshort_wave, &! flag for computing short-wave forcing
lcheckekmo, &! check Ekman, Monin-Obhukov depth limit
linertial, &! flag for using inertial mixing parameterization
lcvmix, &! flag for using CVMix KPP routine instead of POP
lccsm_control_compatible !flag for backwards compatibility with ccsm4 control
character (char_len) :: &
langmuir_opt ! lanmguir turbulence parameterization option
! valid: null, vr12-ma, lf17
integer (int_kind) :: &
num_v_smooth_Ri ! num of times to vertically smooth Ri
real (r8), parameter :: &
epssfc = 0.1_r8 ! non-dimensional extent of sfc layer
real (r8) :: &
Prandtl ! Prandtl number
!-----------------------------------------------------------------------
!
! non-local mixing (counter-gradient mixing), treated as source term
!
!-----------------------------------------------------------------------
real (r8), dimension(:,:,:,:,:), allocatable :: &
KPP_SRC ! non-local mixing (treated as source term)
!-----------------------------------------------------------------------
!
! parameters for subroutine bldepth: computes bndy layer depth
!
! concv = ratio of interior buoyancy frequency to
! buoyancy frequency at entrainment depth
! parameter statement sets the minimum value.
!
!-----------------------------------------------------------------------
real (r8), dimension(km) :: &
Ricr ! crit Rich no. for bndy layer depth
! as a function of vertical resolution
real (r8), parameter :: &
cekman = 0.7_r8, &! coefficient for Ekman depth
cmonob = 1.0_r8, &! coefficient for Monin-Obukhov depth
concv = 1.7_r8, &! minimum allowed value
hbf = 1.0_r8 ! frac of layer affected by sol forcing
!-----------------------------------------------------------------------
!
! parameters for subroutine ri_iwmix which computes
! vertical mixing coefficients below boundary layer due to shear
! instabilities, internal waves and convection
!
!-----------------------------------------------------------------------
real (r8), parameter :: &
Riinfty = 0.8_r8, &! Rich. no. limit for shear instab.
BVSQcon = c0 ! Brunt-Vaisala square cutoff(s**-2)
!-----------------------------------------------------------------------
!
! parameters for subroutine ddmix (double-diffusive mixing)
!
!-----------------------------------------------------------------------
real (r8), parameter :: &
Rrho0 = 2.55_r8, &! limit for double-diff density ratio
dsfmax = 1.0_r8 ! max diffusivity for salt fingering
!-----------------------------------------------------------------------
!
! parameters for subroutine blmix: mixing within boundary layer
!
! cstar = proportionality coefficient for nonlocal transport
! cg = non-dimensional coefficient for counter-gradient term
!
!-----------------------------------------------------------------------
real (r8), parameter :: &
cstar = 10.0_r8 ! coeff for nonlocal transport
real (r8) :: &
cg, &! coefficient for counter-gradient term
Vtc ! resolution and buoyancy independent part of the
! turbulent velocity shear coefficient (for bulk Ri no)
!-----------------------------------------------------------------------
!
! parameters for velocity scale function (from Large et al.)
!
!-----------------------------------------------------------------------
real (r8), parameter :: &
zeta_m = -0.2_r8, &
zeta_s = -1.0_r8, &
c_m = 8.38_r8, &
c_s = 98.96_r8, &
a_m = 1.26_r8, &
a_s = -28.86_r8
!-----------------------------------------------------------------------
!
! common vertical grid arrays used by KPP mixing
!
!-----------------------------------------------------------------------
real (r8), dimension(:), allocatable :: &
zgrid, &! depth at cell interfaces
hwide ! layer thickness at interfaces
!-----------------------------------------------------------------------
!
! module variables to indicate which tidal-mixing method is in use
! defined here to improve efficiency in subroutine ri_iwmix
!
!-----------------------------------------------------------------------
logical (log_kind) :: luse_simmons, luse_schmittner, luse_polzin
!-----------------------------------------------------------------------
!
! ids for tavg diagnostics computed in this module
!
!-----------------------------------------------------------------------
integer (int_kind) :: &
tavg_QSW_HBL, &! tavg id for solar short-wave heat flux in bndry layer
tavg_VDC_BCK, &! tavg id for bckgrnd vertical tracer diffusivity
tavg_VVC_BCK, &! tavg id for bckgrnd vertical momentum viscosity
tavg_KVMIX, &! tavg id for tidal+bckgrnd vertical tracer diffusivity
tavg_KVMIX_M, &! tavg id for tidal+bckgrnd vertical momentum viscosity
tavg_TPOWER
integer (int_kind), dimension(nt) :: &
tavg_KPP_SRC ! tavg id for KPP_SRC for each tracer
! Parameters for CVMix
type(cvmix_global_params_type) :: CVmix_params
!EOC
!***********************************************************************
contains
!***********************************************************************
!BOP
! !IROUTINE: init_vmix_kpp
! !INTERFACE:
subroutine init_vmix_kpp(CVmix_vars, VDC, VVC, convect_diff, convect_visc)
! !DESCRIPTION:
! Initializes constants and reads options for the KPP parameterization.
!
! !REVISION HISTORY:
! same as module
! !INPUT PARAMETERS:
real (r8), intent(in) :: &
convect_diff, &! diffusivity to mimic convection
convect_visc ! viscosity to mimic convection
! !INPUT/OUTPUT PARAMETERS:
type(cvmix_data_type), dimension(:,:), intent(inout) :: CVmix_vars
real (r8), dimension(:,:,:,:), intent(inout) :: &
VVC ! viscosity for momentum diffusion
real (r8), dimension(:,:,0:,:,:),intent(inout) :: &
VDC ! diffusivity for tracer diffusion
!EOP
!BOC
!-----------------------------------------------------------------------
!
! local variables and namelist inputs
!
!-----------------------------------------------------------------------
integer (int_kind) :: &
n, &! local dummy index for tracer
k, &! local dummy index for vertical lvl
i, j, iblock, &! local dummy indexes
nml_error ! namelist i/o error flag
integer (int_kind) :: bid, jny, inx, ic, nlev
real (r8), dimension(:), allocatable :: tmp_array
real (r8), dimension(:,:), allocatable :: tmp_array2
real (r8) :: &
bckgrnd_vdc1, &! background diffusivity (Ledwell)
bckgrnd_vdc_eq, &! equatorial diffusivity (Gregg)
bckgrnd_vdc_psim, &! Max. PSI induced diffusivity (MacKinnon)
bckgrnd_vdc_psin, &! PSI diffusivity in northern hemisphere
bckgrnd_vdc_psis, &! PSI diffusivity in southern hemisphere
bckgrnd_vdc_ban, &! Banda Sea diffusivity (Gordon)
bckgrnd_vdc_dpth, &! depth at which diff equals vdc1
bckgrnd_vdc_linv ! inverse length for transition region
logical (log_kind) :: &
lhoriz_varying_bckgrnd, &
larctic_bckgrnd_vdc ! historically only used as a suboption of lniw_mixing
namelist /vmix_kpp_nml/bckgrnd_vdc1, bckgrnd_vdc2, &
bckgrnd_vdc_eq, bckgrnd_vdc_psim, &
bckgrnd_vdc_ban, &
bckgrnd_vdc_dpth, bckgrnd_vdc_linv, &
Prandtl, rich_mix, &
num_v_smooth_Ri, lrich, ldbl_diff, &
lshort_wave, lcheckekmo, &
larctic_bckgrnd_vdc, &
lhoriz_varying_bckgrnd, &
linertial, lcvmix, langmuir_opt
character (16), parameter :: &
fmt_real = '(a30,2x,1pe12.5)'
character (11), parameter :: &
fmt_log = '(a30,2x,l7)'
character (11), parameter :: &
fmt_int = '(a30,2x,i5)'
character (char_len) :: &
string, string2
character (char_len) :: &
langmuir_mixing_opt, langmuir_entrainment_opt
!-----------------------------------------------------------------------
!
! set defaults for mixing coefficients, then read them from namelist
!
!-----------------------------------------------------------------------
bckgrnd_vdc1 = 0.1_r8
bckgrnd_vdc2 = c0
bckgrnd_vdc_eq = 0.01_r8
bckgrnd_vdc_psim = 0.13_r8
bckgrnd_vdc_ban = c1
bckgrnd_vdc_dpth = 2500.0e02_r8
bckgrnd_vdc_linv = 4.5e-05_r8
Prandtl = 10.0_r8
rich_mix = 50.0_r8
lrich = .true.
ldbl_diff = .false.
lshort_wave = .false.
lcheckekmo = .false.
larctic_bckgrnd_vdc = .false.
lhoriz_varying_bckgrnd = .false.
linertial = .false.
lcvmix = .false.
langmuir_opt = 'null'
num_v_smooth_Ri = 1
if (my_task == master_task) then
open (nml_in, file=nml_filename, status='old',iostat=nml_error)
if (nml_error /= 0) then
nml_error = -1
else
nml_error = 1
endif
do while (nml_error > 0)
read(nml_in, nml=vmix_kpp_nml, iostat=nml_error)
end do
if (nml_error == 0) close(nml_in)
endif
call broadcast_scalar(nml_error, master_task)
if (nml_error /= 0) then
call exit_POP(sigAbort,'ERROR (init_vmix_kpp) reading vmix_kpp_nml')
endif
if (my_task == master_task) then
write(stdout,delim_fmt)
write(stdout,blank_fmt)
write(stdout,*) ' vmix_kpp_nml namelist settings:'
write(stdout,blank_fmt)
write(stdout,vmix_kpp_nml)
write(stdout,blank_fmt)
write(stdout,fmt_real) ' bckgrnd_vdc1 =', bckgrnd_vdc1
write(stdout,fmt_real) ' bckgrnd_vdc2 =', bckgrnd_vdc2
write(stdout,fmt_real) ' bckgrnd_vdc_dpth =', bckgrnd_vdc_dpth
write(stdout,fmt_real) ' bckgrnd_vdc_linv =', bckgrnd_vdc_linv
write(stdout,fmt_real) ' bckgrnd_vdc_eq =', bckgrnd_vdc_eq
write(stdout,fmt_real) ' bckgrnd_vdc_psim =', bckgrnd_vdc_psim
write(stdout,fmt_real) ' bckgrnd_vdc_ban =', bckgrnd_vdc_ban
write(stdout,fmt_real) ' Prandtl =', Prandtl
write(stdout,fmt_real) ' rich_mix =', rich_mix
write(stdout,fmt_log ) ' Ri mixing =', lrich
write(stdout,fmt_log ) ' double-diff =', ldbl_diff
write(stdout,fmt_log ) ' short_wave =', lshort_wave
write(stdout,fmt_log ) ' lcheckekmo =', lcheckekmo
write(stdout,fmt_int ) ' num_smooth_Ri =', num_v_smooth_Ri
write(stdout,fmt_log ) ' larctic_bckgrnd_vdc =', larctic_bckgrnd_vdc
write(stdout,fmt_log ) ' lhoriz_varying_bckgrnd =', lhoriz_varying_bckgrnd
write(stdout,fmt_log ) ' inertial mixing param. =', linertial
write(stdout,fmt_log ) ' use CVMix KPP routines =', lcvmix
select case (langmuir_opt)
case ('null')
write(stdout,'(a38)') ' no Lanmguir turbulence parameterization'
langmuir_mixing_opt = 'NONE'
langmuir_entrainment_opt = 'NONE'
case ('vr12-ma')
write(stdout,'(a38)') ' Langmuir param. option: vr12-ma '
langmuir_mixing_opt = 'LWF16'
langmuir_entrainment_opt = 'LWF16'
case ('lf17')
write(stdout,'(a38)') ' Langmuir param. option: lf17 '
langmuir_mixing_opt = 'LWF16'
langmuir_entrainment_opt = 'LF17'
case default
call exit_POP(sigAbort, &
'ERROR: Unknown option for Langmuir turbulence parameterization')
end select
endif
call broadcast_scalar(bckgrnd_vdc1, master_task)
call broadcast_scalar(bckgrnd_vdc2, master_task)
call broadcast_scalar(bckgrnd_vdc_dpth, master_task)
call broadcast_scalar(bckgrnd_vdc_linv, master_task)
call broadcast_scalar(bckgrnd_vdc_eq , master_task)
call broadcast_scalar(bckgrnd_vdc_psim, master_task)
call broadcast_scalar(bckgrnd_vdc_ban , master_task)
call broadcast_scalar(Prandtl, master_task)
call broadcast_scalar(rich_mix, master_task)
call broadcast_scalar(num_v_smooth_Ri, master_task)
call broadcast_scalar(lrich, master_task)
call broadcast_scalar(ldbl_diff, master_task)
call broadcast_scalar(lshort_wave, master_task)
call broadcast_scalar(lcheckekmo, master_task)
call broadcast_scalar(larctic_bckgrnd_vdc, master_task)
call broadcast_scalar(lhoriz_varying_bckgrnd,master_task)
call broadcast_scalar(linertial, master_task)
call broadcast_scalar(lcvmix, master_task)
call broadcast_scalar(langmuir_opt, master_task)
call broadcast_scalar(langmuir_mixing_opt, master_task)
call broadcast_scalar(langmuir_entrainment_opt, master_task)
if (lcvmix) call register_string('lcvmix')
!-----------------------------------------------------------------------
!
! determine if this case must be backwards compatible with ccsm4 control
!
!-----------------------------------------------------------------------
lccsm_control_compatible = registry_match('lccsm_control_compatible')
!-----------------------------------------------------------------------
!
! define some non-dimensional constants
!
!-----------------------------------------------------------------------
Vtc = sqrt(0.2_r8/c_s/epssfc)/vonkar**2
cg = cstar*vonkar*(c_s*vonkar*epssfc)**p33
pop2cvmixL = 1.0e-2_r8
if (pop2cvmixL /= c0) then
cvmix2popL = c1/pop2cvmixL
else
call exit_POP (sigAbort, 'ERROR (init_vmix_kpp): pop2cvmixL = 0')
endif
pop2cvmixL2 = 1.0e-4_r8
if (pop2cvmixL2 /= c0) then
cvmix2popL2 = c1/pop2cvmixL2
else
call exit_POP (sigAbort, 'ERROR (init_vmix_kpp): pop2cvmixL2 = 0')
endif
!-----------------------------------------------------------------------
!
! define vertical grid coordinates and cell widths
! compute horizontally or vertically varying background
! (internal wave) diffusivity and viscosity
!
! the vertical profile has the functional form
!
! BCKGRND_VDC(z) = vdc1 + vdc2*atan((|z| - dpth)/L) or
! = vdc1 + vdc2*atan((|z| - dpth)*linv)
!
! where
!
! vdc1 = vertical diffusivity at |z|=D (cm^2/s)
! vdc2 = amplitude of variation (cm^2/s)
! linv = inverse length scale for transition region (1/cm)
! dpth = depth where diffusivity is vdc1 (cm)
!
! the viscosity has the same form but multiplied by a constant
! Prandtl number
!
!-----------------------------------------------------------------------
!-----------------------------------------------------------------------
!
! if tidal mixing is on, use vertically constant (internal wave)
! background diffusivity and viscosity (see consistency testing in POP_check)
!
!-----------------------------------------------------------------------
if (bckgrnd_vdc2 /= c0 .and. my_task == master_task) then
write (stdout,blank_fmt)
write (stdout,'(a43)') &
'Vertical Profile for background diffusivity'
endif
!-----------------------------------------------------------------------
!
! error checking
!
!-----------------------------------------------------------------------
if (lhoriz_varying_bckgrnd .and. bckgrnd_vdc2 /= c0) then
call exit_POP (sigAbort, &
'ERROR (init_vmix_kpp): lhoriz_varying_bckgrnd .and. bckgrnd_vdc2 /= c0')
endif
!-----------------------------------------------------------------------
!
! initialize grid info (only need one block since the vertical grid is
! the same across blocks)
!
!-----------------------------------------------------------------------
allocate (zgrid(0:km+1), hwide(0:km+1))
zgrid(0) = eps
hwide(0) = eps
do k=1,km
zgrid(k) = -zt(k)
hwide(k) = dz(k)
enddo
zgrid(km+1) = -zw(km)
hwide(km+1) = eps
allocate (bckgrnd_vvc(nx_block,ny_block,km,nblocks_clinic))
allocate (bckgrnd_vdc(nx_block,ny_block,km,nblocks_clinic))
if (lhoriz_varying_bckgrnd) then
k = 1
do iblock=1,nblocks_clinic
do j=1,ny_block
do i=1,nx_block
bckgrnd_vdc_psis= bckgrnd_vdc_psim*exp(-(0.4_r8*(TLATD(i,j,iblock)+28.9_r8))**c2)
bckgrnd_vdc_psin= bckgrnd_vdc_psim*exp(-(0.4_r8*(TLATD(i,j,iblock)-28.9_r8))**c2)
bckgrnd_vdc(i,j,k,iblock)=bckgrnd_vdc_eq+bckgrnd_vdc_psin+bckgrnd_vdc_psis
if ( TLATD(i,j,iblock) .lt. -10.0_r8 ) then
bckgrnd_vdc(i,j,k,iblock) = bckgrnd_vdc(i,j,k,iblock) + bckgrnd_vdc1
elseif ( TLATD(i,j,iblock) .le. 10.0_r8 ) then
bckgrnd_vdc(i,j,k,iblock) = bckgrnd_vdc(i,j,k,iblock) + &
bckgrnd_vdc1 * (TLATD(i,j,iblock)/10.0_r8)**c2
else
bckgrnd_vdc(i,j,k,iblock)=bckgrnd_vdc(i,j,k,iblock) + bckgrnd_vdc1
endif
!----------------
! North Banda Sea
!----------------
if ( (TLATD(i,j,iblock) .lt. -1.0_r8) .and. (TLATD(i,j,iblock) .gt. -4.0_r8) .and. &
(TLOND(i,j,iblock) .gt. 103.0_r8) .and. (TLOND(i,j,iblock) .lt. 134.0_r8)) then
bckgrnd_vdc(i,j,k,iblock) = bckgrnd_vdc_ban
endif
!-----------------
! Middle Banda Sea
!-----------------
if ( (TLATD(i,j,iblock) .le. -4.0_r8) .and. (TLATD(i,j,iblock) .gt. -7.0_r8) .and. &
(TLOND(i,j,iblock) .gt. 106.0_r8) .and. (TLOND(i,j,iblock) .lt. 140.0_r8)) then
bckgrnd_vdc(i,j,k,iblock) = bckgrnd_vdc_ban
endif
!----------------
! South Banda Sea
!----------------
if ( (TLATD(i,j,iblock) .le. -7.0_r8) .and. (TLATD(i,j,iblock) .gt. -8.3_r8) .and. &
(TLOND(i,j,iblock) .gt. 111.0_r8) .and. (TLOND(i,j,iblock) .lt. 142.0_r8)) then
bckgrnd_vdc(i,j,k,iblock) = bckgrnd_vdc_ban
endif
!----------------
! Arctic
!----------------
if (larctic_bckgrnd_vdc) then ! historically only used as a suboption of lniw_mixing
if (TLATD(i,j,iblock) .ge. 70.0_r8) then
bckgrnd_vdc(i,j,k,iblock) = bckgrnd_vdc_eq
endif
endif
bckgrnd_vvc(i,j,k,iblock) = Prandtl*bckgrnd_vdc(i,j,k,iblock)
end do ! i
end do ! j
end do ! iblock
do k=2,km
bckgrnd_vdc(:,:,k,:) = bckgrnd_vdc(:,:,1,:)
bckgrnd_vvc(:,:,k,:) = bckgrnd_vvc(:,:,1,:)
enddo
else
!-----------------------------------------------------------------------
!
! only need one block since the vertical grid is the same across blocks
!
!-----------------------------------------------------------------------
do k=1,km
bckgrnd_vdc(:,:,k,:) = bckgrnd_vdc1 + bckgrnd_vdc2* &
atan(bckgrnd_vdc_linv* &
(zw(k)-bckgrnd_vdc_dpth))
bckgrnd_vvc(:,:,k,:) = Prandtl*bckgrnd_vdc(:,:,k,:)
if (bckgrnd_vdc2 /= c0 .and. my_task == master_task) then
write (stdout,'(2x,e13.6)') bckgrnd_vdc(1,1,k,1)
endif
end do
endif ! lhoriz_varying_bckgrnd
!-----------------------------------------------------------------------
!
! compute crit Rich number for bndy layer depth as a function
! of model vertical resolution. note that hwide is in cm.
!
!-----------------------------------------------------------------------
Ricr(1:km) = 0.3_r8
!-----------------------------------------------------------------------
!
! allocate and initialize KPP-specific arrays
!
!-----------------------------------------------------------------------
allocate (HMXL (nx_block,ny_block,nblocks_clinic), &
HMXL_DR (nx_block,ny_block,nblocks_clinic), & ! QL, 150526
KPP_HBLT (nx_block,ny_block,nblocks_clinic), &
KPP_SRC (nx_block,ny_block,km,nt,nblocks_clinic))
HMXL = c0
HMXL_DR = c0 ! QL, 150526
KPP_HBLT = c0
KPP_SRC = c0
VDC = c0
VVC = c0
!-----------------------------------------------------------------------
!
! Set up CVMix variables
!
!-----------------------------------------------------------------------
if ( lcvmix ) then
call cvmix_init_ddiff(strat_param_max=2.55_r8, &
kappa_ddiff_s=1e-4_r8, &
ddiff_exp1=1.0_r8, &
ddiff_exp2=3.0_r8, &
mol_diff=1.5e-6_r8, &
kappa_ddiff_param1=0.909_r8, &
kappa_ddiff_param2=4.6_r8, &
kappa_ddiff_param3=-0.54_r8, &
diff_conv_type="MC76")
call cvmix_init_conv(convect_diff=convect_diff*1e-4_r8, &
convect_visc=convect_visc*1e-4_r8, &
lBruntVaisala=.true., &
BVsqr_convect=BVSQcon)
call cvmix_init_shear(mix_scheme="KPP", &
KPP_nu_zero=rich_mix*1e-4_r8, &
KPP_Ri_zero=Riinfty, &
KPP_exp=real(3,r8))
call cvmix_init_kpp(lEkman=lcheckekmo, &
lMonOb=lcheckekmo, &
surf_layer_ext = epssfc, &
minVtsqr=c0, &
lnoDGat1=.false., &
langmuir_mixing_str=langmuir_mixing_opt, &
langmuir_entrainment_str=langmuir_entrainment_opt, &
MatchTechnique="MatchBoth")
call cvmix_put_kpp("a_m", a_m)
call cvmix_put_kpp("a_s", a_s)
call cvmix_put_kpp("c_m", c_m)
call cvmix_put_kpp("c_s", c_s)
call cvmix_put(CVmix_params, 'prandtl', Prandtl)
! convert rho from g/cm^3 -> kg / m^3
call cvmix_put(CVmix_params, 'FreshWaterDensity', rho_fw*c1000)
do bid=1,nblocks_clinic
do jny=1,ny_block
do inx=1,nx_block
ic = (jny-1)*nx_block + inx
nlev = KMT(inx, jny, bid)
call cvmix_put(CVmix_vars(ic,bid), 'nlev', nlev)
call cvmix_put(CVmix_vars(ic,bid), 'max_nlev', nlev)
if (nlev.gt.0) then
call cvmix_put(CVmix_vars(ic,bid), 'Mdiff', 0._r8)
call cvmix_put(CVmix_vars(ic,bid), 'Tdiff', 0._r8)
call cvmix_put(CVmix_vars(ic,bid), 'Sdiff', 0._r8)
call cvmix_put(CVmix_vars(ic,bid), 'strat_param_num', 0._r8)
call cvmix_put(CVmix_vars(ic,bid), 'strat_param_denom', 0._r8)
call cvmix_put(CVmix_vars(ic,bid), 'Nsqr', 0._r8)
call cvmix_put(CVmix_vars(ic,bid), 'rho', 0._r8)
call cvmix_put(CVmix_vars(ic,bid), 'rho_lwr', 0._r8)
call cvmix_put(CVmix_vars(ic,bid), 'lat', TLATD(inx,jny,bid))
call cvmix_put(CVmix_vars(ic,bid), 'lon', TLOND(inx,jny,bid))
allocate(tmp_array(nlev+1))
tmp_array = 0._r8
allocate(tmp_array2(nlev+1,nlev+1))
tmp_array2 = 0._r8
! Shear Richardson
call cvmix_put(CVmix_vars(ic,bid), 'Ri_iface', tmp_array)
! Bulk Richardson
call cvmix_put(CVmix_vars(ic,bid), 'Ri_bulk', tmp_array(1:nlev))
! vertical levels (cell center and interfaces)
call cvmix_put(CVmix_vars(ic,bid), 'zt', -zt(1:nlev)*1e-2_r8)
call cvmix_put(CVmix_vars(ic,bid), 'zw', tmp_array)
CVmix_vars(ic,bid)%zw_iface(2:nlev+1) = -zw(1:nlev)*1e-2_r8
! Depth at bottom of ocean
call cvmix_put(CVmix_vars(ic,bid), 'ocn_depth', zw(nlev)*1e-2_r8)
if (partial_bottom_cells) then
CVMix_vars(ic,bid)%zt_cntr(nlev) = (-zw(nlev-1)-p5*DZT(inx, jny, nlev, bid))*1e-2_r8
CVMix_vars(ic,bid)%zw_iface(nlev+1) = (-zw(nlev-1)-DZT(inx, jny, nlev, bid))*1e-2_r8
CVMix_vars(ic,bid)%OceanDepth = CVMix_vars(ic,bid)%zw_iface(nlev+1)
endif
! Initialize SchmittnerCoeff
call cvmix_put(CVmix_vars(ic,bid), 'SchmittnerCoeff', tmp_array)
! Initialize SchmittnerSouthernOcean
call cvmix_put(CVmix_vars(ic,bid), 'SchmittnerSouthernOcean', tmp_array)
! Initialize Schmittner time-invariant term
call cvmix_put(CVmix_vars(ic,bid), 'exp_hab_zetar', tmp_array2)
deallocate(tmp_array,tmp_array2)
! NOTE: cvmix_init_tidal is called from tidal_mixing.F90 but we
! still need to compute time-invariant tidal quantites
! (requires nlev, rho, zt, and zw)
if (ltidal_mixing) then
select case (tidal_mixing_method_itype)
case (tidal_mixing_method_jayne)
!-------------------------------------------------------------
! cvmix_compute_Simmons forms SimmonsCoefficient from
! (input)q*(input)E. However, when using ER03 and
! GN13 options with a 2D method such as Simmons, only
!(qE) is available, not E decoupled from q. So in POP,
! we will always use as input to cvmix_compute_Simmons_invariant
! the quantity TIDAL_QE_2D = (qE), regardless of whether the energy
! source is 2D or 3D, and then always divide SimmonsCoeff by q
! afterwards. Note that the operation
! (f(q*E)*q)/q is round-off level different from E*q
!--------------------------------------------------------------
call cvmix_compute_Simmons_invariant(CVmix_vars(ic,bid), &
CVmix_params, &
TIDAL_QE_2D(inx,jny,bid)/c1000)
call cvmix_put(CVmix_vars(ic,bid), 'SimmonsCoeff', &
CVmix_vars(ic,bid)%SimmonsCoeff/tidal_local_mixing_fraction)
!--------------------------------------------------------------
! form the time-invariant part of Schmittner coefficient term
!--------------------------------------------------------------
call cvmix_compute_socn_tidal_invariant(CVmix_vars(ic,bid))
case (tidal_mixing_method_schmittner)
!--------------------------------------------------------------
! form the time-invariant part of Schmittner coefficient term
!--------------------------------------------------------------
call cvmix_compute_Schmittner_invariant(CVmix_vars(ic,bid), &
CVmix_params)
!--------------------------------------------------------------
! 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(CVmix_vars(ic,bid), &
nlev, &
TIDAL_QE_3D(inx,jny,:,bid)/c1000)
!-----------------------------------------------------------------------
! form the time-invariant part of the Schmittner coefficient term
!-----------------------------------------------------------------------
call cvmix_compute_socn_tidal_invariant(CVmix_vars(ic,bid))
end select
end if
end if
end do
end do
end do
end if
!-----------------------------------------------------------------------
!
! define module-wide logical values to indicate which tidal-mixing
! method is active. This should improve efficiency in
! subroutine ri_iwmix.
!
!-----------------------------------------------------------------------
luse_polzin = .false.
luse_simmons = .false.
luse_schmittner = .false.
select case (tidal_mixing_method_itype)
case (tidal_mixing_method_jayne)
luse_simmons = .true.
case (tidal_mixing_method_schmittner)
luse_schmittner = .true.
case (tidal_mixing_method_polzin)
luse_polzin = .true.
end select
!-----------------------------------------------------------------------
!
! allocate and initialize the eddy-induced speed array
!
!-----------------------------------------------------------------------
if ( linertial ) then
allocate ( BOLUS_SP(nx_block,ny_block,nblocks_clinic) )
BOLUS_SP = c0
endif
!-----------------------------------------------------------------------
!
! define tavg fields computed from vmix_kpp module routines
!
!-----------------------------------------------------------------------
string = 'Solar Short-Wave Heat Flux in bndry layer'
call define_tavg_field(tavg_QSW_HBL,'QSW_HBL',2, &
long_name=trim(string), &
units='watt/m^2', grid_loc='2110', &
coordinates='TLONG TLAT time')
if (ltidal_mixing) then
string = 'Vertical diabatic diffusivity due to Tidal Mixing + background'
else
string = 'Vertical diabatic diffusivity due to background'
endif
call define_tavg_field(tavg_KVMIX,'KVMIX',3, &
long_name=trim(string), &
units='centimeter^2/s', &
grid_loc='3113', &
coordinates ='TLONG TLAT z_w_bot time' )
string = 'Vertical diabatic diffusivity due to background'
call define_tavg_field(tavg_VDC_BCK,'VDC_BCK',3, &
long_name=trim(string), &
units='centimeter^2/s', &
grid_loc='3113', &
coordinates ='TLONG TLAT z_w_bot time' )
if (ltidal_mixing) then
string = 'Vertical viscosity due to Tidal Mixing + background'
else
string = 'Vertical viscosity due to background'
endif
call define_tavg_field(tavg_KVMIX_M,'KVMIX_M',3, &
long_name=trim(string), &
units='centimeter^2/s', &
grid_loc='3113', &
coordinates ='TLONG TLAT z_w_bot time' )
string = 'Vertical viscosity due to background'
call define_tavg_field(tavg_VVC_BCK,'VVC_BCK',3, &
long_name=trim(string), &
units='centimeter^2/s', &
grid_loc='3113', &
coordinates ='TLONG TLAT z_w_bot time' )
string = 'Energy Used by Vertical Mixing'
call define_tavg_field(tavg_TPOWER,'TPOWER',3, &
long_name=trim(string), &
units='erg/s/cm^3', &
grid_loc='3113', &
coordinates ='TLONG TLAT z_w_bot time' )
do n = 1,nt
string = 'KPP_SRC_' /&
&/ trim(tracer_d(n)%short_name)
string2 = trim(tracer_d(n)%short_name) /&
&/ ' tendency from KPP non local mixing term'
call define_tavg_field(tavg_KPP_SRC(n),trim(string),3, &
long_name=trim(string2), &
units=trim(tracer_d(n)%tend_units), &
scale_factor=tracer_d(n)%scale_factor,&
grid_loc='3111', &
coordinates ='TLONG TLAT z_t time' )
enddo
!-----------------------------------------------------------------------
!EOC
call flushm (stdout)
end subroutine init_vmix_kpp
!***********************************************************************
!BOP
! !IROUTINE: vmix_coeffs_kpp
! !INTERFACE:
subroutine vmix_coeffs_kpp(CVmix_vars, VDC, VVC, TRCR, UUU, VVV, UCUR, VCUR, &
RHOMIX, STF, SHF_QSW, this_block, convect_diff, &
convect_visc, SMF, SMFT)
! !DESCRIPTION:
! This is the main driver routine which calculates the vertical
! mixing coefficients for the KPP mixing scheme as outlined in
! Large, McWilliams and Doney, Reviews of Geophysics, 32, 363
! (November 1994). The non-local mixing is also computed here, but
! is treated as a source term in baroclinic.
!
!----------------------------------------------------------------------------------
! Updated late 2010/early 2011 to include near inertial wave parameterization.
! Final code description is as follows. We use the diffusivity k to describe
! the order of computation.
!
! k is diffusivity, k_w = background, k_n = near inertial wave, k_t = tidal,
! k_s = shear (Richardson), k_d = double diffusion, and k_c = convection.
! HBLT is boundary layer depth. The diffusivities in the description below
! show total combination and limits after each routine is called. k_n_max is
! the maximum near inertial wave diffusivity allowed, and k_t_max is the
! corresponding maximum for the tidal diffusivity.
!
! routine description
!
! buoydiff (computes buoyancy difference with surface)
! bldepth (computes boundary layer depth HBLT)
! iw_reset initialize k: k = k_w
! niw_mix (compute niwm diffusivity k_n below HBLT,
! and extend upward the top value into the BL for the
! matching slope condition on k required by blmix)
! k_w' = min(max(k_w,k_n),k_n_max)
! ri_iwmix k = min((k_w' + k_t),k_t_max) + k_s
! if( dbl_diff) ddmix k = min((k_w' + k_t),k_t_max) + k_s + k_d
! blmix (computes k in boundary layer using present interior k)
! ..... (computes interior convective mixing)
! ..... k = min((min(max(k_w,k_n),k_n_max) + k_t),k_t_max) + k_s + k_d + k_c
!
! Thus, in words, the background diffusivity k_w is initialized first. Then the
! near inertial wave diffusivity is set, which is then combined with the background
! such that when k_n is greater than the background, it is chosen; otherwise, the
! background is unchanged. Then, the tidal is evaluated and combined with the
! modified background (including near inertial wave) in a similar manner. Finally,
! the shear, double diffusion (if present) and the convection diffusivities are
! added to the modified background/near inertial wave/tidal diffusivity.
!
! Bruce P. Briegleb April 2011
!----------------------------------------------------------------------------------
!
! !REVISION HISTORY:
! same as module
! !INPUT PARAMETERS:
real (r8), dimension(nx_block,ny_block,km,nt), intent(in) :: &
TRCR ! tracers at current time
real (r8), dimension(nx_block,ny_block,km), intent(in) :: &
UUU,VVV, &! velocities at mix time
UCUR,VCUR ! velocities at current time
real (r8), dimension(nx_block,ny_block,km), intent(in) :: &
RHOMIX ! density at mix time
real (r8), dimension(nx_block,ny_block,nt), intent(in) :: &
STF ! surface forcing for all tracers
real (r8), dimension(nx_block,ny_block), intent(in) :: &
SHF_QSW ! short-wave forcing
real (r8), dimension(nx_block,ny_block,2), intent(in), optional :: &
SMF, &! surface momentum forcing at U points
SMFT ! surface momentum forcing at T points
! *** either one or the other (not
! *** both) should be passed
real (r8), intent(in) :: &
convect_diff, &! diffusivity to mimic convection
convect_visc ! viscosity to mimic convection
type (block), intent(in) :: &
this_block ! block information for current block