<p><b>mpetersen@lanl.gov</b> 2012-02-16 09:41:42 -0700 (Thu, 16 Feb 2012)</p><p>Add tendency subroutines for high frequency thickness and low frequency divergence. Zlevel using RK4 and split KE results remains identical to the trunk (p75p-u). RK4 can now run in ztilde mode.<br>
</p><hr noshade><pre><font color="gray">Modified: branches/ocean_projects/ztilde/src/core_ocean/Registry
===================================================================
--- branches/ocean_projects/ztilde/src/core_ocean/Registry        2012-02-14 22:01:52 UTC (rev 1509)
+++ branches/ocean_projects/ztilde/src/core_ocean/Registry        2012-02-16 16:41:42 UTC (rev 1510)
@@ -23,7 +23,7 @@
namelist character grid config_pressure_type pressure
namelist real grid config_rho0 1028
namelist real grid config_div_low_freq_time 86400.0
-namelist real grid config_h_restore_time 86400.0
+namelist real grid config_hhf_restore_time 86400.0
namelist real grid config_h_high_freq_diff 1.0e-5
namelist integer split_explicit_ts config_n_ts_iter 2
namelist integer split_explicit_ts config_n_bcl_iter_beg 2
Modified: branches/ocean_projects/ztilde/src/core_ocean/mpas_ocn_tendency.F
===================================================================
--- branches/ocean_projects/ztilde/src/core_ocean/mpas_ocn_tendency.F        2012-02-14 22:01:52 UTC (rev 1509)
+++ branches/ocean_projects/ztilde/src/core_ocean/mpas_ocn_tendency.F        2012-02-16 16:41:42 UTC (rev 1510)
@@ -41,6 +41,9 @@
use ocn_time_average
+ ! mrp 120215 added temporarily for hHighFreq tendency
+ use ocn_tracer_hmix_del2
+
implicit none
private
save
@@ -64,6 +67,8 @@
public :: ocn_tend_h, &
ocn_tend_u, &
+ ocn_tend_hHighFreq, &
+ ocn_tend_divergenceLowFreq, &
ocn_tend_scalar, &
ocn_diagnostic_solve, &
ocn_div_hu, &
@@ -409,7 +414,7 @@
!
!-----------------------------------------------------------------------
- subroutine ocn_tend_hHighFreq(tend, s, d, grid)!{{{
+ subroutine ocn_tend_hHighFreq(tend, s, grid)!{{{
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! Compute high frequency thickness tendency for z-tilde vertical grid
!
@@ -423,41 +428,53 @@
type (tend_type), intent(inout) :: tend
type (state_type), intent(inout) :: s
- type (diagnostics_type), intent(in) :: d
type (mesh_type), intent(in) :: grid
- real (kind=RKIND), dimension(:,:), pointer :: tend_hHighFreq
-
integer :: err
+ integer :: nCells, nVertLevels, iCell, k
+ integer, dimension(:), pointer :: maxLevelCell
+ real (kind=RKIND), dimension(:), pointer :: div_hu_btr
+ real (kind=RKIND), dimension(:,:), pointer :: &
+ tend_hHighFreq, hHighFreq, div_hu, divergenceLowFreq
+
call mpas_timer_start("ocn_tend_hHighFreq")
+ err = 0
- tend_hHighFreq => tend % hHighFreq % array
-
+ nCells = grid % nCells
+ nVertLevels = grid % nVertLevels
+
+ divergenceLowFreq => s % divergenceLowFreq % array
+ hHighFreq => s % hHighFreq % array
+ div_hu => s % div_hu % array
+ div_hu_btr => s % div_hu_btr % array
+
+ maxLevelCell => grid % maxLevelCell % array
+
+ tend_hHighFreq => tend % hHighFreq % array
+
!
- ! height tendency: start accumulating tendency terms
+ ! Low Frequency Divergence Tendency
!
tend_hHighFreq = 0.0
+ do iCell=1,nCells
+ do k=1,maxLevelCell(iCell)
+ tend_hHighFreq(k,iCell) = &
+ - div_hu(k,iCell) + div_hu_btr(iCell) + divergenceLowFreq(k,iCell) &
+ -2.0*pii/config_hhf_restore_time * hHighFreq(k,iCell)
+ end do
+ end do
!
- ! height tendency: horizontal advection term -</font>
<font color="blue">abla\cdot ( hu)
+ ! tracer tendency: del2 horizontal hhf diffusion, div(\kappa_{hf} </font>
<font color="gray">abla h^{hf})
!
- ! See Ringler et al. (2010) jcp paper, eqn 19, 21, and fig. 3.
- ! for explanation of divergence operator.
- !
-! call mpas_timer_start("hadv", .false., thickHadvTimer)
-! call ocn_thick_hadv_tend(grid, uTransport, h_edge, tend_h, err)
-! call mpas_timer_stop("hadv", thickHadvTimer)
+ call mpas_timer_start("hmix", .false., tracerHmixTimer)
+ call ocn_hhf_hmix_del2_tend(grid, hHighFreq, tend_hHighFreq, err)
+ call mpas_timer_stop("hmix", tracerHmixTimer)
- !
- ! height tendency: vertical advection term -d/dz(hw)
- !
-! call mpas_timer_start("vadv", .false., thickVadvTimer)
-! call ocn_thick_vadv_tend(grid, wtop, tend_h, err)
-! call mpas_timer_stop("vadv", thickVadvTimer)
call mpas_timer_stop("ocn_tend_hHighFreq")
-
+
end subroutine ocn_tend_hHighFreq!}}}
!***********************************************************************
@@ -473,7 +490,7 @@
!
!-----------------------------------------------------------------------
- subroutine ocn_tend_divergenceLowFreq(tend, s, d, grid)!{{{
+ subroutine ocn_tend_divergenceLowFreq(tend, s, grid)!{{{
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! Compute height and normal wind tendencies, as well as diagnostic variables
!
@@ -487,22 +504,42 @@
type (tend_type), intent(inout) :: tend
type (state_type), intent(inout) :: s
- type (diagnostics_type), intent(in) :: d
type (mesh_type), intent(in) :: grid
- real (kind=RKIND), dimension(:,:), pointer :: tend_divergenceLowFreq
-
integer :: err
+ integer :: nCells, nVertLevels, iCell, k
+ integer, dimension(:), pointer :: maxLevelCell
+ real (kind=RKIND), dimension(:), pointer :: div_hu_btr
+ real (kind=RKIND), dimension(:,:), pointer :: tend_divergenceLowFreq, div_hu, divergenceLowFreq
+
call mpas_timer_start("ocn_tend_divergenceLowFreq")
+ err = 0
- tend_divergenceLowFreq => tend % divergenceLowFreq % array
+ nCells = grid % nCells
+ nVertLevels = grid % nVertLevels
+ divergenceLowFreq => s % divergenceLowFreq % array
+ div_hu => s % div_hu % array
+ div_hu_btr => s % div_hu_btr % array
+
+ maxLevelCell => grid % maxLevelCell % array
+
+ tend_divergenceLowFreq => tend % divergenceLowFreq % array
+
!
- ! height tendency: start accumulating tendency terms
+ ! Low Frequency Divergence Tendency
!
tend_divergenceLowFreq = 0.0
+ do iCell=1,nCells
+ do k=1,maxLevelCell(iCell)
+ tend_divergenceLowFreq(k,iCell) = &
+ -2.0*pii/config_div_low_freq_time &
+ *(divergenceLowFreq(k,iCell) - div_hu(k,iCell) + div_hu_btr(iCell))
+ end do
+ end do
+
call mpas_timer_stop("ocn_tend_divergenceLowFreq")
end subroutine ocn_tend_divergenceLowFreq!}}}
Modified: branches/ocean_projects/ztilde/src/core_ocean/mpas_ocn_time_integration_rk4.F
===================================================================
--- branches/ocean_projects/ztilde/src/core_ocean/mpas_ocn_time_integration_rk4.F        2012-02-14 22:01:52 UTC (rev 1509)
+++ branches/ocean_projects/ztilde/src/core_ocean/mpas_ocn_time_integration_rk4.F        2012-02-16 16:41:42 UTC (rev 1510)
@@ -112,7 +112,7 @@
block % state % time_levs(2) % state % u % array(:,:) = block % state % time_levs(1) % state % u % array(:,:)
block % state % time_levs(2) % state % h % array(:,:) = block % state % time_levs(1) % state % h % array(:,:)
if (config_vert_grid_type.eq.'ztilde') then
- block % state % time_levs(2) % state % hHighFreq% array(:,:) = block % state % time_levs(1) % state % hHighFreq% array(:,:)
+ block % state % time_levs(2) % state % hHighFreq % array(:,:) = block % state % time_levs(1) % state % hHighFreq % array(:,:)
block % state % time_levs(2) % state % divergenceLowFreq % array(:,:) = block % state % time_levs(1) % state % divergenceLowFreq % array(:,:)
endif
do iCell=1,block % mesh % nCells ! couple tracers to h
@@ -173,8 +173,8 @@
do while (associated(block))
if (config_vert_grid_type.eq.'ztilde') then
- !call ocn_tend_hHighFreq(block % tend, provis, block % diagnostics, block % mesh)
- !call ocn_tend_divergenceLowFreq(block % tend, provis, block % diagnostics, block % mesh)
+ call ocn_tend_hHighFreq(block % tend, provis, block % mesh)
+ call ocn_tend_divergenceLowFreq(block % tend, provis, block % mesh)
endif
call ocn_div_hu(provis, provis, block % mesh)
Modified: branches/ocean_projects/ztilde/src/core_ocean/mpas_ocn_tracer_hmix_del2.F
===================================================================
--- branches/ocean_projects/ztilde/src/core_ocean/mpas_ocn_tracer_hmix_del2.F        2012-02-14 22:01:52 UTC (rev 1509)
+++ branches/ocean_projects/ztilde/src/core_ocean/mpas_ocn_tracer_hmix_del2.F        2012-02-16 16:41:42 UTC (rev 1510)
@@ -36,6 +36,7 @@
!--------------------------------------------------------------------
public :: ocn_tracer_hmix_del2_tend, &
+ ocn_hhf_hmix_del2_tend, &
ocn_tracer_hmix_del2_init
!--------------------------------------------------------------------
@@ -158,7 +159,7 @@
do k=1,maxLevelEdgeTop(iEdge)
do iTracer=1,num_tracers
- ! \kappa_2 </font>
<font color="blue">abla \phi on edge
+ ! </font>
<font color="black">abla \phi on edge
tracer_turb_flux = tracers(iTracer,k,cell2) - tracers(iTracer,k,cell1)
! div(h \kappa_2 </font>
<font color="gray">abla \phi) at cell center
@@ -176,6 +177,121 @@
!***********************************************************************
!
+! routine ocn_hhf_hmix_del2_tend
+!
+!> \brief Computes laplacian tendency term for horizontal hhf mixing
+!> \author Doug Jacobsen
+!> \date 19 September 2011
+!> \version SVN:$Id$
+!> \details
+!> This routine computes the horizontal mixing tendency for
+!> high frequency thickness
+!> based on current state using a laplacian parameterization.
+!
+!-----------------------------------------------------------------------
+
+ subroutine ocn_hhf_hmix_del2_tend(grid, hHighFreq, tend, err)!{{{
+
+ !-----------------------------------------------------------------
+ !
+ ! input variables
+ !
+ !-----------------------------------------------------------------
+
+ type (mesh_type), intent(in) :: &
+ grid !< Input: grid information
+
+ real (kind=RKIND), dimension(:,:), intent(in) :: &
+ hHighFreq !< Input: high frequency thickness
+
+ !-----------------------------------------------------------------
+ !
+ ! input/output variables
+ !
+ !-----------------------------------------------------------------
+
+ real (kind=RKIND), dimension(:,:), intent(inout) :: &
+ tend !< Input/Output: high freq thickness tendency
+
+ !-----------------------------------------------------------------
+ !
+ ! output variables
+ !
+ !-----------------------------------------------------------------
+
+ integer, intent(out) :: err !< Output: error flag
+
+ !-----------------------------------------------------------------
+ !
+ ! local variables
+ !
+ !-----------------------------------------------------------------
+
+ integer :: iEdge, nEdges, nVertLevels, cell1, cell2, k
+
+ integer, dimension(:,:), allocatable :: boundaryMask
+
+ integer, dimension(:), pointer :: maxLevelEdgeTop
+ integer, dimension(:,:), pointer :: cellsOnEdge, edgeMask
+
+ real (kind=RKIND) :: invAreaCell1, invAreaCell2
+ real (kind=RKIND) :: hhf_turb_flux, flux, r_tmp
+
+ real (kind=RKIND), dimension(:), pointer :: areaCell, dvEdge, dcEdge
+ real (kind=RKIND), dimension(:), pointer :: meshScalingDel2
+
+ !-----------------------------------------------------------------
+ !
+ ! call relevant routines for computing tendencies
+ ! note that the user can choose multiple options and the
+ ! tendencies will be added together
+ !
+ !-----------------------------------------------------------------
+
+ err = 0
+
+ if (.not.del2On) return
+
+ nEdges = grid % nEdges
+ nVertLevels = grid % nVertLevels
+
+ maxLevelEdgeTop => grid % maxLevelEdgeTop % array
+ cellsOnEdge => grid % cellsOnEdge % array
+ edgeMask => grid % edgeMask % array
+ areaCell => grid % areaCell % array
+ dvEdge => grid % dvEdge % array
+ dcEdge => grid % dcEdge % array
+ meshScalingDel2 => grid % meshScalingDel2 % array
+
+ !
+ ! compute a boundary mask to enforce insulating boundary conditions in the horizontal
+ !
+ do iEdge=1,nEdges
+ cell1 = cellsOnEdge(1,iEdge)
+ cell2 = cellsOnEdge(2,iEdge)
+ invAreaCell1 = 1.0/areaCell(cell1)
+ invAreaCell2 = 1.0/areaCell(cell2)
+
+ r_tmp = meshScalingDel2(iEdge) * config_h_high_freq_diff * dvEdge(iEdge) / dcEdge(iEdge)
+
+ do k=1,maxLevelEdgeTop(iEdge)
+ ! </font>
<font color="blue">abla h^{hf} on edge
+ hhf_turb_flux = hHighFreq(k,cell2) - hHighFreq(k,cell1)
+
+ ! div(\kappa_{hf} </font>
<font color="blue">abla h^{hf}) at cell center
+ flux = hhf_turb_flux * edgeMask(k, iEdge) * r_tmp
+
+ tend(k,cell1) = tend(k,cell1) + flux * invAreaCell1
+ tend(k,cell2) = tend(k,cell2) - flux * invAreaCell2
+ end do
+
+ end do
+ !--------------------------------------------------------------------
+
+ end subroutine ocn_hhf_hmix_del2_tend!}}}
+
+!***********************************************************************
+!
! routine ocn_tracer_hmix_del2_init
!
!> \brief Initializes ocean tracer horizontal mixing quantities
</font>
</pre>