374   type(verticalgrid_type),   
intent(in)    :: gv
   375   real, 
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
   377   real, 
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
   379   real, 
dimension(SZI_(G),SZJ_(G),SZK_(G)),  &
   381   real, 
dimension(SZI_(G),SZJ_(G),SZK_(G)),  &
   383   real, 
dimension(SZI_(G),SZJ_(G),SZK_(G)),  &
   385   type(thermo_var_ptrs),     
intent(inout) :: tv
   387   type(forcing),             
intent(in)    :: fluxes
   389   type(optics_type),         
pointer       :: optics
   390   type(vertvisc_type),       
intent(inout) :: visc
   393   real,                      
intent(in)    :: dt
   394   type(set_diffusivity_cs),  
pointer       :: cs
   395   real, 
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
   397   real, 
dimension(SZI_(G),SZJ_(G),SZK_(G)+1), &
   398                    optional, 
intent(out)   :: kd_int
   417   real, 
dimension(SZI_(G)) :: &
   420   real, 
dimension(SZI_(G), SZJ_(G)) :: &
   423   type(diffusivity_diags) :: dd 
   425   real, 
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: &
   429   real, 
dimension(SZI_(G),SZK_(G)) :: &
   436   real, 
dimension(SZI_(G),SZK_(G)+1) :: &
   460   type(p3d) :: z_ptrs(6)    
   461   integer   :: kb(szi_(g))  
   463   integer   :: num_z_diags  
   465   logical   :: showcalltree 
   467   integer :: i, j, k, is, ie, js, je, nz
   468   integer :: isd, ied, jsd, jed
   473   is  = g%isc ; ie  = g%iec ; js  = g%jsc ; je  = g%jec ; nz = g%ke
   474   isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
   475   showcalltree = calltree_showquery()
   476   if (showcalltree) 
call calltree_enter(
"set_diffusivity(), MOM_set_diffusivity.F90")
   478   if (.not.
associated(cs)) 
call mom_error(fatal,
"set_diffusivity: "//&
   479          "Module must be initialized before it is used.")
   484   deg_to_rad = atan(1.0)/45.0 
   485   omega2     = cs%Omega*cs%Omega
   486   i_2omega   = 0.5/cs%Omega
   489   use_eos = 
associated(tv%eqn_of_state)
   491   if ((cs%double_diffusion) .and. &
   492       .not.(
associated(visc%Kd_extra_T) .and. 
associated(visc%Kd_extra_S)) ) &
   493     call mom_error(fatal, 
"set_diffusivity: visc%Kd_extra_T and "//&
   494          "visc%Kd_extra_S must be associated when DOUBLE_DIFFUSION is true.")
   498   if ((cs%id_N2 > 0) .or. (cs%id_N2_z > 0)) 
then   499     allocate(dd%N2_3d(isd:ied,jsd:jed,nz+1)) ; dd%N2_3d(:,:,:) = 0.0
   501   if ((cs%id_Kd_itidal > 0) .or. (cs%id_Kd_itidal_z > 0) .or. &
   502       (cs%id_Kd_Itidal_work > 0)) 
then   503     allocate(dd%Kd_itidal(isd:ied,jsd:jed,nz+1)) ; dd%Kd_itidal(:,:,:) = 0.0
   505   if ((cs%id_Kd_lowmode > 0) .or. (cs%id_Kd_lowmode_z > 0) .or. &
   506       (cs%id_Kd_lowmode_work > 0)) 
then   507     allocate(dd%Kd_lowmode(isd:ied,jsd:jed,nz+1)) ; dd%Kd_lowmode(:,:,:) = 0.0
   509   if ( (cs%id_Fl_itidal > 0) ) 
then   510     allocate(dd%Fl_itidal(isd:ied,jsd:jed,nz+1)) ; dd%Fl_itidal(:,:,:) = 0.0
   512   if ( (cs%id_Fl_lowmode > 0) ) 
then   513     allocate(dd%Fl_lowmode(isd:ied,jsd:jed,nz+1)) ; dd%Fl_lowmode(:,:,:) = 0.0
   515   if ( (cs%id_Polzin_decay_scale > 0) ) 
then   516     allocate(dd%Polzin_decay_scale(isd:ied,jsd:jed))
   517     dd%Polzin_decay_scale(:,:) = 0.0
   519   if ( (cs%id_Polzin_decay_scale_scaled > 0) ) 
then   520     allocate(dd%Polzin_decay_scale_scaled(isd:ied,jsd:jed))
   521     dd%Polzin_decay_scale_scaled(:,:) = 0.0
   523   if ( (cs%id_N2_bot > 0) ) 
then   524     allocate(dd%N2_bot(isd:ied,jsd:jed)) ; dd%N2_bot(:,:) = 0.0
   526   if ( (cs%id_N2_meanz > 0) ) 
then   527     allocate(dd%N2_meanz(isd:ied,jsd:jed)) ; dd%N2_meanz(:,:) = 0.0
   529   if ((cs%id_Kd_Niku > 0) .or. (cs%id_Kd_Niku_z > 0) .or. &
   530       (cs%id_Kd_Niku_work > 0)) 
then   531     allocate(dd%Kd_Niku(isd:ied,jsd:jed,nz+1)) ; dd%Kd_Niku(:,:,:) = 0.0
   533   if ((cs%id_Kd_user > 0) .or. (cs%id_Kd_user_z > 0)) 
then   534     allocate(dd%Kd_user(isd:ied,jsd:jed,nz+1)) ; dd%Kd_user(:,:,:) = 0.0
   536   if (cs%id_Kd_work > 0) 
then   537     allocate(dd%Kd_work(isd:ied,jsd:jed,nz)) ; dd%Kd_work(:,:,:) = 0.0
   539   if (cs%id_Kd_Niku_work > 0) 
then   540     allocate(dd%Kd_Niku_work(isd:ied,jsd:jed,nz)) ; dd%Kd_Niku_work(:,:,:) = 0.0
   542   if (cs%id_Kd_Itidal_work > 0) 
then   543     allocate(dd%Kd_Itidal_work(isd:ied,jsd:jed,nz))
   544     dd%Kd_Itidal_work(:,:,:) = 0.0
   546   if (cs%id_Kd_Lowmode_Work > 0) 
then   547     allocate(dd%Kd_Lowmode_Work(isd:ied,jsd:jed,nz))
   548     dd%Kd_Lowmode_Work(:,:,:) = 0.0
   550   if (cs%id_TKE_itidal > 0) 
then   551     allocate(dd%TKE_Itidal_used(isd:ied,jsd:jed)) ; dd%TKE_Itidal_used(:,:) = 0.
   553   if (cs%id_maxTKE > 0) 
then   554     allocate(dd%maxTKE(isd:ied,jsd:jed,nz)) ; dd%maxTKE(:,:,:) = 0.0
   556   if (cs%id_TKE_to_Kd > 0) 
then   557     allocate(dd%TKE_to_Kd(isd:ied,jsd:jed,nz)) ; dd%TKE_to_Kd(:,:,:) = 0.0
   559   if ((cs%id_KT_extra > 0) .or. (cs%id_KT_extra_z > 0)) 
then   560     allocate(dd%KT_extra(isd:ied,jsd:jed,nz+1)) ; dd%KT_extra(:,:,:) = 0.0
   562   if ((cs%id_KS_extra > 0) .or. (cs%id_KS_extra_z > 0)) 
then   563     allocate(dd%KS_extra(isd:ied,jsd:jed,nz+1)) ; dd%KS_extra(:,:,:) = 0.0
   565   if ((cs%id_Kd_BBL > 0) .or. (cs%id_Kd_BBL_z > 0)) 
then   566     allocate(dd%Kd_BBL(isd:ied,jsd:jed,nz+1)) ; dd%Kd_BBL(:,:,:) = 0.0
   572       call hchksum(tv%T, 
"before vert_fill_TS tv%T",g%HI)
   573       call hchksum(tv%S, 
"before vert_fill_TS tv%S",g%HI)
   574       call hchksum(h, 
"before vert_fill_TS h",g%HI, scale=gv%H_to_m)
   576     call vert_fill_ts(h, tv%T, tv%S, kappa_fill, dt_fill, t_f, s_f, g, gv)
   578       call hchksum(tv%T, 
"after vert_fill_TS tv%T",g%HI)
   579       call hchksum(tv%S, 
"after vert_fill_TS tv%S",g%HI)
   580       call hchksum(h, 
"after vert_fill_TS h",g%HI, scale=gv%H_to_m)
   584   if (cs%useKappaShear) 
then   586       call hchksum(u_h, 
"before calc_KS u_h",g%HI)
   587       call hchksum(v_h, 
"before calc_KS v_h",g%HI)
   589     call cpu_clock_begin(id_clock_kappashear)
   592     call calculate_kappa_shear(u_h, v_h, h, tv, fluxes%p_surf, visc%Kd_turb, visc%TKE_turb, &
   593                                visc%Kv_turb, dt, g, gv, cs%kappaShear_CSp)
   594     call cpu_clock_end(id_clock_kappashear)
   596       call hchksum(visc%Kd_turb, 
"after calc_KS visc%Kd_turb",g%HI)
   597       call hchksum(visc%Kv_turb, 
"after calc_KS visc%Kv_turb",g%HI)
   598       call hchksum(visc%TKE_turb, 
"after calc_KS visc%TKE_turb",g%HI)
   600     if (showcalltree) 
call calltree_waypoint(
"done with calculate_kappa_shear (set_diffusivity)")
   601  elseif (cs%useCVMix) 
then   603     call calculate_cvmix_shear(u_h, v_h, h, tv, visc%Kd_turb, visc%Kv_turb,g,gv,cs%CVMix_shear_CSp)
   604   elseif (
associated(visc%Kv_turb)) 
then   605     visc%Kv_turb(:,:,:) = 0. 
   612   if (cs%Bryan_Lewis_diffusivity) 
then   614     do j=js,je ; 
do i=is,ie
   615       kd_sfc(i,j) = cs%Kd_Bryan_Lewis_surface
   619     do j=js,je ; 
do i=is,ie
   623   if (cs%Henyey_IGW_background) 
then   624     i_x30 = 2.0 / invcosh(cs%N0_2Omega*2.0) 
   627     do j=js,je ; 
do i=is,ie
   628       abs_sin = abs(sin(g%geoLatT(i,j)*deg_to_rad))
   629       kd_sfc(i,j) = max(cs%Kd_min, kd_sfc(i,j) * &
   630            ((abs_sin * invcosh(cs%N0_2Omega/max(epsilon,abs_sin))) * i_x30) )
   632   elseif (cs%Kd_tanh_lat_fn) 
then   634     do j=js,je ; 
do i=is,ie
   638       kd_sfc(i,j) = max(cs%Kd_min, kd_sfc(i,j) * (1.0 + &
   639           cs%Kd_tanh_lat_scale * 0.5*tanh((abs(g%geoLatT(i,j)) - 35.0)/5.0) ))
   643   if (cs%debug) 
call hchksum(kd_sfc,
"Kd_sfc",g%HI,haloshift=0)
   653     call find_n2(h, tv, t_f, s_f, fluxes, j, g, gv, cs, drho_int, n2_lay, n2_int, n2_bot)
   654     if (
associated(dd%N2_3d)) 
then   655       do k=1,nz+1 ; 
do i=is,ie ; dd%N2_3d(i,j,k) = n2_int(i,k) ;
 enddo ;
 enddo   659     if (cs%Bryan_Lewis_diffusivity) 
then   660       i_trans = 1.0 / cs%Bryan_Lewis_width_trans
   661       atan_fn_sfc = atan(cs%Bryan_Lewis_depth_cent*i_trans)
   662       i_atan_fn = 1.0 / (2.0*atan(1.0) + atan_fn_sfc)
   663       do i=is,ie ; depth(i) = 0.0 ;
 enddo   664       do k=1,nz ; 
do i=is,ie
   665         atan_fn_lay = atan((cs%Bryan_Lewis_depth_cent - &
   666                             (depth(i)+0.5*gv%H_to_m*h(i,j,k)))*i_trans)
   667         kd(i,j,k) = kd_sfc(i,j) + (cs%Kd_Bryan_Lewis_deep - kd_sfc(i,j)) * &
   668                                   (atan_fn_sfc - atan_fn_lay) * i_atan_fn
   669         depth(i) = depth(i) + gv%H_to_m*h(i,j,k)
   671     elseif ((.not.cs%bulkmixedlayer) .and. (cs%Kd /= cs%Kdml)) 
then   672       i_hmix = 1.0 / cs%Hmix
   673       do i=is,ie ; depth(i) = 0.0 ;
 enddo   674       do k=1,nz ; 
do i=is,ie
   675         depth_c = depth(i) + 0.5*gv%H_to_m*h(i,j,k)
   677         if (depth_c <= cs%Hmix) 
then ; kd(i,j,k) = cs%Kdml
   678         elseif (depth_c >= 2.0*cs%Hmix) 
then ; kd(i,j,k) = kd_sfc(i,j)
   680           kd(i,j,k) = ((kd_sfc(i,j) - cs%Kdml) * i_hmix) * depth_c + &
   681                       (2.0*cs%Kdml - kd_sfc(i,j))
   684         depth(i) = depth(i) + gv%H_to_m*h(i,j,k)
   686     elseif (cs%Henyey_IGW_background_new) 
then   687       i_x30 = 2.0 / invcosh(cs%N0_2Omega*2.0) 
   688       do k=1,nz ; 
do i=is,ie
   689         abs_sin = max(epsilon,abs(sin(g%geoLatT(i,j)*deg_to_rad)))
   690         n_2omega = max(abs_sin,sqrt(n2_lay(i,k))*i_2omega)
   691         n02_n2 = (cs%N0_2Omega/n_2omega)**2
   692         kd(i,j,k) = max(cs%Kd_min, kd_sfc(i,j) * &
   693              ((abs_sin * invcosh(n_2omega/abs_sin)) * i_x30)*n02_n2)
   696       do k=1,nz ; 
do i=is,ie
   697         kd(i,j,k) = kd_sfc(i,j)
   701     if (cs%double_diffusion) 
then   702       call double_diffusion(tv, h, t_f, s_f, j, g, gv, cs, kt_extra, ks_extra)
   703       do k=2,nz ; 
do i=is,ie
   704         if (ks_extra(i,k) > kt_extra(i,k)) 
then    705           kd(i,j,k-1) = kd(i,j,k-1) + 0.5*kt_extra(i,k)
   706           kd(i,j,k)   = kd(i,j,k)   + 0.5*kt_extra(i,k)
   707           visc%Kd_extra_S(i,j,k) = ks_extra(i,k)-kt_extra(i,k)
   708           visc%Kd_extra_T(i,j,k) = 0.0
   709         elseif (kt_extra(i,k) > 0.0) 
then    710           kd(i,j,k-1) = kd(i,j,k-1) + 0.5*ks_extra(i,k)
   711           kd(i,j,k)   = kd(i,j,k)   + 0.5*ks_extra(i,k)
   712           visc%Kd_extra_T(i,j,k) = kt_extra(i,k)-ks_extra(i,k)
   713           visc%Kd_extra_S(i,j,k) = 0.0
   715           visc%Kd_extra_T(i,j,k) = 0.0
   716           visc%Kd_extra_S(i,j,k) = 0.0
   719       if (
associated(dd%KT_extra)) 
then ; 
do k=1,nz+1 ; 
do i=is,ie
   720         dd%KT_extra(i,j,k) = kt_extra(i,k)
   721       enddo ;
 enddo ;
 endif   723       if (
associated(dd%KS_extra)) 
then ; 
do k=1,nz+1 ; 
do i=is,ie
   724         dd%KS_extra(i,j,k) = ks_extra(i,k)
   725       enddo ;
 enddo ;
 endif   729     if (cs%useKappaShear .or. cs%useCVMix) 
then   730       if (
present(kd_int)) 
then   731         do k=2,nz ; 
do i=is,ie
   732           kd_int(i,j,k) = visc%Kd_turb(i,j,k) + 0.5*(kd(i,j,k-1) + kd(i,j,k))
   735           kd_int(i,j,1) = visc%Kd_turb(i,j,1) 
   736           kd_int(i,j,nz+1) = 0.0
   739       do k=1,nz ; 
do i=is,ie
   740         kd(i,j,k) = kd(i,j,k) + 0.5*(visc%Kd_turb(i,j,k) + visc%Kd_turb(i,j,k+1))
   743       if (
present(kd_int)) 
then   745           kd_int(i,j,1) = kd(i,j,1) ; kd_int(i,j,nz+1) = 0.0
   747         do k=2,nz ; 
do i=is,ie
   748           kd_int(i,j,k) = 0.5*(kd(i,j,k-1) + kd(i,j,k))
   753     call find_tke_to_kd(h, tv, drho_int, n2_lay, j, dt, g, gv, cs, tke_to_kd, &
   755     if (
associated(dd%maxTKE)) 
then ; 
do k=1,nz ; 
do i=is,ie
   756       dd%maxTKE(i,j,k) = maxtke(i,k)
   757     enddo ;
 enddo ;
 endif   758     if (
associated(dd%TKE_to_Kd)) 
then ; 
do k=1,nz ; 
do i=is,ie
   759       dd%TKE_to_Kd(i,j,k) = tke_to_kd(i,k)
   760     enddo ;
 enddo ;
 endif   763     if (cs%ML_radiation) &
   764       call add_mlrad_diffusivity(h, fluxes, j, g, gv, cs, kd, tke_to_kd, kd_int)
   767     if (cs%Int_tide_dissipation .or. cs%Lee_wave_dissipation .or. cs%Lowmode_itidal_dissipation) &
   768       call add_int_tide_diffusivity(h, n2_bot, j, tke_to_kd, maxtke, g, gv, cs, &
   769                                     dd, n2_lay, kd, kd_int)
   773     if (cs%bottomdraglaw .and. (cs%BBL_effic>0.0)) 
then   774       if (cs%use_LOTW_BBL_diffusivity) 
then   775         call add_lotw_bbl_diffusivity(h, u, v, tv, fluxes, visc, j, n2_int, g, gv, cs,  &
   776                                       kd, kd_int, dd%Kd_BBL)
   778         call add_drag_diffusivity(h, u, v, tv, fluxes, visc, j, tke_to_kd, &
   779                                   maxtke, kb, g, gv, cs, kd, kd_int, dd%Kd_BBL)
   783     if (cs%limit_dissipation) 
then   784       do k=2,nz-1 ; 
do i=is,ie
   790         dissip = max( cs%dissip_min, &   
   791                       cs%dissip_N0 + cs%dissip_N1 * sqrt(n2_lay(i,k)), & 
   792                       cs%dissip_N2 * n2_lay(i,k) ) 
   793         kd(i,j,k) = max( kd(i,j,k) , &  
   794            dissip * (cs%FluxRi_max / (gv%Rho0 * (n2_lay(i,k) + omega2))) )
   797       if (
present(kd_int)) 
then ; 
do k=2,nz ; 
do i=is,ie
   803         dissip = max( cs%dissip_min, &   
   804                       cs%dissip_N0 + cs%dissip_N1 * sqrt(n2_int(i,k)), & 
   805                       cs%dissip_N2 * n2_int(i,k) ) 
   806         kd_int(i,j,k) = max( kd_int(i,j,k) , &  
   807            dissip * (cs%FluxRi_max / (gv%Rho0 * (n2_int(i,k) + omega2))) )
   808       enddo ;
 enddo ;
 endif   811     if (
associated(dd%Kd_work)) 
then   812       do k=1,nz ; 
do i=is,ie
   813         dd%Kd_Work(i,j,k)  = gv%Rho0 * kd(i,j,k) * n2_lay(i,k) * &
   820     call hchksum(kd,
"BBL Kd",g%HI,haloshift=0)
   821     if (cs%useKappaShear) 
call hchksum(visc%Kd_turb,
"Turbulent Kd",g%HI,haloshift=0)
   822     if (
associated(visc%kv_bbl_u) .and. 
associated(visc%kv_bbl_v)) 
then   823       call uvchksum(
"BBL Kv_bbl_[uv]", visc%kv_bbl_u, visc%kv_bbl_v, &
   824                     g%HI, 0, symmetric=.true.)
   826     if (
associated(visc%bbl_thick_u) .and. 
associated(visc%bbl_thick_v)) 
then   827       call uvchksum(
"BBL bbl_thick_[uv]", visc%bbl_thick_u, &
   828                     visc%bbl_thick_v, g%HI, 0, symmetric=.true.)
   830     if (
associated(visc%Ray_u) .and. 
associated(visc%Ray_v)) 
then   831       call uvchksum(
"Ray_[uv]", visc%Ray_u, visc%Ray_v, g%HI, 0, symmetric=.true.)
   835   if (cs%Kd_add > 0.0) 
then   836     if (
present(kd_int)) 
then   838       do k=1,nz ; 
do j=js,je ; 
do i=is,ie
   839         kd_int(i,j,k) = kd_int(i,j,k) + cs%Kd_add
   840         kd(i,j,k) = kd(i,j,k) + cs%Kd_add
   841       enddo ;
 enddo ;
 enddo   844       do k=1,nz ; 
do j=js,je ; 
do i=is,ie
   845         kd(i,j,k) = kd(i,j,k) + cs%Kd_add
   846       enddo ;
 enddo ;
 enddo   850   if (cs%user_change_diff) 
then   851     call user_change_diff(h, tv, g, cs%user_change_diff_CSp, kd, kd_int, &
   852                           t_f, s_f, dd%Kd_user)
   855   if (cs%id_Kd_layer > 0) 
call post_data(cs%id_Kd_layer, kd, cs%diag)
   858   if (cs%Int_tide_dissipation .or. cs%Lee_wave_dissipation .or. cs%Lowmode_itidal_dissipation) 
then   859     if (cs%id_TKE_itidal  > 0) 
call post_data(cs%id_TKE_itidal,  dd%TKE_itidal_used, cs%diag)
   860     if (cs%id_TKE_leewave > 0) 
call post_data(cs%id_TKE_leewave, cs%TKE_Niku,        cs%diag)
   861     if (cs%id_Nb          > 0) 
call post_data(cs%id_Nb,      cs%Nb,      cs%diag)
   862     if (cs%id_N2          > 0) 
call post_data(cs%id_N2,      dd%N2_3d,   cs%diag)
   863     if (cs%id_N2_bot      > 0) 
call post_data(cs%id_N2_bot,  dd%N2_bot,  cs%diag)
   864     if (cs%id_N2_meanz    > 0) 
call post_data(cs%id_N2_meanz,dd%N2_meanz,cs%diag)
   866     if (cs%id_Fl_itidal > 0) 
call post_data(cs%id_Fl_itidal, dd%Fl_itidal, cs%diag)
   867     if (cs%id_Kd_itidal > 0) 
call post_data(cs%id_Kd_itidal, dd%Kd_itidal, cs%diag)
   868     if (cs%id_Kd_Niku   > 0) 
call post_data(cs%id_Kd_Niku,   dd%Kd_Niku,   cs%diag)
   869     if (cs%id_Kd_lowmode> 0) 
call post_data(cs%id_Kd_lowmode, dd%Kd_lowmode, cs%diag)
   870     if (cs%id_Fl_lowmode> 0) 
call post_data(cs%id_Fl_lowmode, dd%Fl_lowmode, cs%diag)
   871     if (cs%id_Kd_user   > 0) 
call post_data(cs%id_Kd_user,   dd%Kd_user,   cs%diag)
   872     if (cs%id_Kd_Work   > 0) 
call post_data(cs%id_Kd_Work,   dd%Kd_Work,   cs%diag)
   873     if (cs%id_Kd_Itidal_Work > 0) &
   874       call post_data(cs%id_Kd_Itidal_Work, dd%Kd_Itidal_Work, cs%diag)
   875     if (cs%id_Kd_Niku_Work > 0) 
call post_data(cs%id_Kd_Niku_Work, dd%Kd_Niku_Work, cs%diag)
   876     if (cs%id_Kd_Lowmode_Work > 0) &
   877       call post_data(cs%id_Kd_Lowmode_Work, dd%Kd_Lowmode_Work, cs%diag)
   878     if (cs%id_maxTKE > 0) 
call post_data(cs%id_maxTKE, dd%maxTKE, cs%diag)
   879     if (cs%id_TKE_to_Kd > 0) 
call post_data(cs%id_TKE_to_Kd, dd%TKE_to_Kd, cs%diag)
   881     if (cs%id_Polzin_decay_scale > 0 ) &
   882       call post_data(cs%id_Polzin_decay_scale, dd%Polzin_decay_scale, cs%diag)
   883     if (cs%id_Polzin_decay_scale_scaled > 0 ) &
   884       call post_data(cs%id_Polzin_decay_scale_scaled, dd%Polzin_decay_scale_scaled, cs%diag)
   886     if (cs%id_Kd_itidal_z > 0) 
then   887       num_z_diags        = num_z_diags + 1
   888       z_ids(num_z_diags) = cs%id_Kd_itidal_z
   889       z_ptrs(num_z_diags)%p => dd%Kd_itidal
   892     if (cs%id_Kd_Niku_z > 0) 
then   893       num_z_diags        = num_z_diags + 1
   894       z_ids(num_z_diags) = cs%id_Kd_Niku_z
   895       z_ptrs(num_z_diags)%p => dd%Kd_Niku
   898     if (cs%id_Kd_lowmode_z > 0) 
then   899       num_z_diags        = num_z_diags + 1
   900       z_ids(num_z_diags) = cs%id_Kd_lowmode_z
   901       z_ptrs(num_z_diags)%p => dd%Kd_lowmode
   904     if (cs%id_N2_z > 0) 
then   905       num_z_diags = num_z_diags + 1
   906       z_ids(num_z_diags) = cs%id_N2_z
   907       z_ptrs(num_z_diags)%p => dd%N2_3d
   910     if (cs%id_Kd_user_z > 0) 
then   911       num_z_diags        = num_z_diags + 1
   912       z_ids(num_z_diags) = cs%id_Kd_user_z
   913       z_ptrs(num_z_diags)%p => dd%Kd_user
   918   if (cs%id_KT_extra > 0) 
call post_data(cs%id_KT_extra, dd%KT_extra, cs%diag)
   919   if (cs%id_KS_extra > 0) 
call post_data(cs%id_KS_extra, dd%KS_extra, cs%diag)
   920   if (cs%id_Kd_BBL > 0)   
call post_data(cs%id_Kd_BBL, dd%Kd_BBL, cs%diag)
   922   if (cs%id_KT_extra_z > 0) 
then   923       num_z_diags = num_z_diags + 1
   924       z_ids(num_z_diags) = cs%id_KT_extra_z
   925       z_ptrs(num_z_diags)%p => dd%KT_extra
   928   if (cs%id_KS_extra_z > 0) 
then   929       num_z_diags = num_z_diags + 1
   930       z_ids(num_z_diags) = cs%id_KS_extra_z
   931       z_ptrs(num_z_diags)%p => dd%KS_extra
   934   if (cs%id_Kd_BBL_z > 0) 
then   935       num_z_diags = num_z_diags + 1
   936       z_ids(num_z_diags) = cs%id_Kd_BBL_z
   937       z_ptrs(num_z_diags)%p => dd%KS_extra
   940   if (num_z_diags > 0) &
   941     call calc_zint_diags(h, z_ptrs, z_ids, num_z_diags, g, gv, cs%diag_to_Z_CSp)
   943   if (
associated(dd%N2_3d)) 
deallocate(dd%N2_3d)
   944   if (
associated(dd%Kd_itidal)) 
deallocate(dd%Kd_itidal)
   945   if (
associated(dd%Kd_lowmode)) 
deallocate(dd%Kd_lowmode)
   946   if (
associated(dd%Fl_itidal)) 
deallocate(dd%Fl_itidal)
   947   if (
associated(dd%Fl_lowmode)) 
deallocate(dd%Fl_lowmode)
   948   if (
associated(dd%Polzin_decay_scale)) 
deallocate(dd%Polzin_decay_scale)
   949   if (
associated(dd%Polzin_decay_scale_scaled)) 
deallocate(dd%Polzin_decay_scale_scaled)
   950   if (
associated(dd%N2_bot)) 
deallocate(dd%N2_bot)
   951   if (
associated(dd%N2_meanz)) 
deallocate(dd%N2_meanz)
   952   if (
associated(dd%Kd_work)) 
deallocate(dd%Kd_work)
   953   if (
associated(dd%Kd_user)) 
deallocate(dd%Kd_user)
   954   if (
associated(dd%Kd_Niku)) 
deallocate(dd%Kd_Niku)
   955   if (
associated(dd%Kd_Niku_work)) 
deallocate(dd%Kd_Niku_work)
   956   if (
associated(dd%Kd_Itidal_Work))  
deallocate(dd%Kd_Itidal_Work)
   957   if (
associated(dd%Kd_Lowmode_Work)) 
deallocate(dd%Kd_Lowmode_Work)
   958   if (
associated(dd%TKE_itidal_used)) 
deallocate(dd%TKE_itidal_used)
   959   if (
associated(dd%maxTKE)) 
deallocate(dd%maxTKE)
   960   if (
associated(dd%TKE_to_Kd)) 
deallocate(dd%TKE_to_Kd)
   961   if (
associated(dd%KT_extra)) 
deallocate(dd%KT_extra)
   962   if (
associated(dd%KS_extra)) 
deallocate(dd%KS_extra)
   963   if (
associated(dd%Kd_BBL)) 
deallocate(dd%Kd_BBL)
   965   if (showcalltree) 
call calltree_leave(
"set_diffusivity()")
 Ocean grid type. See mom_grid for details.