<p><b>mpetersen@lanl.gov</b> 2012-10-03 09:56:44 -0600 (Wed, 03 Oct 2012)</p><p>BRANCH COMMIT<br>
<br>
Add Leith closure to acc branch. This work was done by Todd for the big runs for the first mpas-o paper.<br>
<br>
>From Todd:<br>
<br>
Here are the changes:<br>
M mpas_ocn_vel_hmix_del2.F<br>
This is the implementation of the Leith closure<br>
<br>
M Registry<br>
This adds the Leith config parameters, writes viscosity to output (and turns off tracer1)<br>
<br>
M mpas_ocn_vel_hmix.F<br>
This pushed viscosity into the subroutine arguments<br>
<br>
M mpas_ocn_tendency.F<br>
This pushed viscosity into the subroutine arguments<br>
<br>
At present, the del2 coefficient has to be non-zero so that the Leith code is called.<br>
</p><hr noshade><pre><font color="gray">Modified: branches/ocean_projects/acc/src/core_ocean/Registry
===================================================================
--- branches/ocean_projects/acc/src/core_ocean/Registry        2012-10-03 15:31:41 UTC (rev 2172)
+++ branches/ocean_projects/acc/src/core_ocean/Registry        2012-10-03 15:56:44 UTC (rev 2173)
@@ -58,6 +58,10 @@
namelist logical hmix config_rayleigh_friction false
namelist real hmix config_rayleigh_damping_coeff 0.0
namelist real hmix config_apvm_scale_factor 0.0
+namelist logical hmix config_use_leith_del2 false
+namelist real hmix config_leith_parameter 0.0
+namelist real hmix config_leith_dx 0.0
+namelist real hmix config_leith_visc2_max 1000000.0
namelist character vmix config_vert_visc_type const
namelist character vmix config_vert_diff_type const
namelist logical vmix config_implicit_vertical_mix .true.
@@ -219,7 +223,6 @@
var persistent real rho ( nVertLevels nCells Time ) 2 iro rho state - -
var persistent real temperature ( nVertLevels nCells Time ) 2 iro temperature state tracers dynamics
var persistent real salinity ( nVertLevels nCells Time ) 2 iro salinity state tracers dynamics
-var persistent real tracer1 ( nVertLevels nCells Time ) 2 iro tracer1 state tracers testing
% Tendency variables: neither read nor written to any files
var persistent real tend_u ( nVertLevels nEdges Time ) 1 - u tend - -
@@ -227,7 +230,6 @@
var persistent real tend_ssh ( nCells Time ) 1 - ssh tend - -
var persistent real tend_temperature ( nVertLevels nCells Time ) 1 - temperature tend tracers dynamics
var persistent real tend_salinity ( nVertLevels nCells Time ) 1 - salinity tend tracers dynamics
-var persistent real tend_tracer1 ( nVertLevels nCells Time ) 1 - tracer1 tend tracers testing
% state variables for Split Explicit timesplitting
var persistent real uBtr ( nEdges Time ) 2 - uBtr state - -
@@ -276,6 +278,7 @@
var persistent real pressure ( nVertLevels nCells Time ) 2 - pressure state - -
var persistent real wTop ( nVertLevelsP1 nCells Time ) 2 - wTop state - -
var persistent real rhoDisplaced ( nVertLevels nCells Time ) 2 - rhoDisplaced state - -
+var persistent real viscosity ( nVertLevels nCells Time ) 2 o viscosity state - -
% Other diagnostic variables: neither read nor written to any files
var persistent real vh ( nVertLevels nEdges Time ) 2 - vh state - -
Modified: branches/ocean_projects/acc/src/core_ocean/mpas_ocn_tendency.F
===================================================================
--- branches/ocean_projects/acc/src/core_ocean/mpas_ocn_tendency.F        2012-10-03 15:31:41 UTC (rev 2172)
+++ branches/ocean_projects/acc/src/core_ocean/mpas_ocn_tendency.F        2012-10-03 15:56:44 UTC (rev 2173)
@@ -172,7 +172,7 @@
real (kind=RKIND), dimension(:,:), pointer :: &
h_edge, h, u, rho, zMid, pressure, &
- tend_u, circulation, vorticity, ke, ke_edge, Vor_edge, &
+ tend_u, circulation, vorticity, viscosity, ke, ke_edge, Vor_edge, &
MontPot, wTop, divergence, vertViscTopOfEdge
real (kind=RKIND), dimension(:,:), pointer :: u_src
@@ -186,6 +186,7 @@
wTop => s % wTop % array
zMid => s % zMid % array
h_edge => s % h_edge % array
+ viscosity => s % viscosity % array
vorticity => s % vorticity % array
divergence => s % divergence % array
ke => s % ke % array
@@ -237,7 +238,7 @@
! strictly only valid for config_h_mom_eddy_visc2 == constant
!
call mpas_timer_start("hmix", .false., velHmixTimer)
- call ocn_vel_hmix_tend(grid, divergence, vorticity, tend_u, err)
+ call ocn_vel_hmix_tend(grid, divergence, vorticity, viscosity, tend_u, err)
call mpas_timer_stop("hmix", velHmixTimer)
!
Modified: branches/ocean_projects/acc/src/core_ocean/mpas_ocn_vel_hmix.F
===================================================================
--- branches/ocean_projects/acc/src/core_ocean/mpas_ocn_vel_hmix.F        2012-10-03 15:31:41 UTC (rev 2172)
+++ branches/ocean_projects/acc/src/core_ocean/mpas_ocn_vel_hmix.F        2012-10-03 15:56:44 UTC (rev 2173)
@@ -72,7 +72,7 @@
!
!-----------------------------------------------------------------------
- subroutine ocn_vel_hmix_tend(grid, divergence, vorticity, tend, err)!{{{
+ subroutine ocn_vel_hmix_tend(grid, divergence, vorticity, viscosity, tend, err)!{{{
!-----------------------------------------------------------------
!
@@ -98,6 +98,10 @@
real (kind=RKIND), dimension(:,:), intent(inout) :: &
tend !< Input/Output: velocity tendency
+ real (kind=RKIND), dimension(:,:), intent(inout) :: &
+ viscosity !< Input: viscosity
+
+
!-----------------------------------------------------------------
!
! output variables
@@ -123,7 +127,7 @@
!-----------------------------------------------------------------
call mpas_timer_start("del2", .false., del2Timer)
- call ocn_vel_hmix_del2_tend(grid, divergence, vorticity, tend, err1)
+ call ocn_vel_hmix_del2_tend(grid, divergence, vorticity, viscosity, tend, err1)
call mpas_timer_stop("del2", del2Timer)
call mpas_timer_start("del4", .false., del4Timer)
call ocn_vel_hmix_del4_tend(grid, divergence, vorticity, tend, err2)
Modified: branches/ocean_projects/acc/src/core_ocean/mpas_ocn_vel_hmix_del2.F
===================================================================
--- branches/ocean_projects/acc/src/core_ocean/mpas_ocn_vel_hmix_del2.F        2012-10-03 15:31:41 UTC (rev 2172)
+++ branches/ocean_projects/acc/src/core_ocean/mpas_ocn_vel_hmix_del2.F        2012-10-03 15:56:44 UTC (rev 2173)
@@ -70,7 +70,7 @@
!
!-----------------------------------------------------------------------
- subroutine ocn_vel_hmix_del2_tend(grid, divergence, vorticity, tend, err)!{{{
+ subroutine ocn_vel_hmix_del2_tend(grid, divergence, vorticity, viscosity, tend, err)!{{{
!-----------------------------------------------------------------
!
@@ -96,6 +96,8 @@
real (kind=RKIND), dimension(:,:), intent(inout) :: &
tend !< Input/Output: velocity tendency
+ real (kind=RKIND), dimension(:,:), intent(inout) :: &
+ viscosity !< Input: viscosity
!-----------------------------------------------------------------
!
@@ -116,8 +118,8 @@
integer, dimension(:), pointer :: maxLevelEdgeTop
integer, dimension(:,:), pointer :: cellsOnEdge, verticesOnEdge, edgeMask
- real (kind=RKIND) :: u_diffusion, invLength1, invLength2
- real (kind=RKIND), dimension(:), pointer :: meshScalingDel2, &
+ real (kind=RKIND) :: u_diffusion, invLength1, invLength2, visc2
+ real (kind=RKIND), dimension(:), pointer :: meshScalingDel2, meshScaling, &
dcEdge, dvEdge
!-----------------------------------------------------------------
@@ -135,10 +137,13 @@
cellsOnEdge => grid % cellsOnEdge % array
verticesOnEdge => grid % verticesOnEdge % array
meshScalingDel2 => grid % meshScalingDel2 % array
+ meshScaling => grid % meshScaling % array
edgeMask => grid % edgeMask % array
dcEdge => grid % dcEdge % array
dvEdge => grid % dvEdge % array
+ viscosity = 1.0
+
do iEdge=1,nEdgesSolve
cell1 = cellsOnEdge(1,iEdge)
cell2 = cellsOnEdge(2,iEdge)
@@ -158,8 +163,20 @@
-viscVortCoef &
*( vorticity(k,vertex2) - vorticity(k,vertex1) ) * invLength2
- u_diffusion = meshScalingDel2(iEdge) * eddyVisc2 * u_diffusion
+ if(config_use_leith_del2) then
+ visc2 = ( config_leith_parameter * config_leith_dx * meshScaling(iEdge) / 3.14)**3 &
+ * abs( vorticity(k,vertex2) - vorticity(k,vertex1) ) * invLength1 * sqrt(3.0)
+ visc2 = min(visc2, config_leith_visc2_max)
+ u_diffusion = visc2 * u_diffusion
+ else
+ visc2 = eddyVisc2
+ u_diffusion = meshScalingDel2(iEdge) * visc2 * u_diffusion
+ endif
+ viscosity(k,cell1) = visc2
+ viscosity(k,cell2) = visc2
+ ! write(0,*) visc2, config_leith_parameter, config_leith_dx, config_use_leith_del2
+
tend(k,iEdge) = tend(k,iEdge) + edgeMask(k, iEdge) * u_diffusion
end do
</font>
</pre>