<p><b>dwj07@fsu.edu</b> 2011-09-16 13:27:37 -0600 (Fri, 16 Sep 2011)</p><p><br>
        -- BRANCH COMMIT --<br>
<br>
        Adding modules for horizontal and vertical advection of thickness.<br>
</p><hr noshade><pre><font color="gray">Modified: branches/ocean_projects/performance/src/core_ocean/Makefile
===================================================================
--- branches/ocean_projects/performance/src/core_ocean/Makefile        2011-09-16 18:29:23 UTC (rev 1005)
+++ branches/ocean_projects/performance/src/core_ocean/Makefile        2011-09-16 19:27:37 UTC (rev 1006)
@@ -3,6 +3,8 @@
OBJS = module_mpas_core.o \
module_test_cases.o \
module_advection.o \
+         module_OcnThickHadv.o \
+         module_OcnThickVadv.o \
         module_OcnVelCoriolis.o \
         module_OcnVelVadv.o \
         module_OcnVelHmix.o \
@@ -28,6 +30,10 @@
module_global_diagnostics.o:
+module_OcnThickHadv.o:
+
+module_OcnThickVadv.o:
+
module_OcnVelPressureGrad.o:
module_OcnVelVadv.o:
@@ -47,6 +53,8 @@
module_OcnVelCoriolis.o:
module_mpas_core.o: module_advection.o \
+                                        module_OcnThickHadv.o \
+                                        module_OcnThickVadv.o \
                                        module_global_diagnostics.o \
                                        module_test_cases.o \
                                        module_OcnVelCoriolis.o \
Added: branches/ocean_projects/performance/src/core_ocean/module_OcnThickHadv.F
===================================================================
--- branches/ocean_projects/performance/src/core_ocean/module_OcnThickHadv.F         (rev 0)
+++ branches/ocean_projects/performance/src/core_ocean/module_OcnThickHadv.F        2011-09-16 19:27:37 UTC (rev 1006)
@@ -0,0 +1,210 @@
+!|||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||
+!
+! OcnThickHadv
+!
+!> \brief MPAS ocean horizontal advection for thickness
+!> \author Doug Jacobsen
+!> \date 16 September 2011
+!> \version SVN:$Id:$
+!> \details
+!> This module contains the routine for computing
+!> tendencies for thickness from horizontal advection
+!
+!-----------------------------------------------------------------------
+
+module OcnThickHadv
+
+ use grid_types
+ use configure
+
+ implicit none
+ private
+ save
+
+ !--------------------------------------------------------------------
+ !
+ ! Public parameters
+ !
+ !--------------------------------------------------------------------
+
+ !--------------------------------------------------------------------
+ !
+ ! Public member functions
+ !
+ !--------------------------------------------------------------------
+
+ public :: OcnThickHadvTend, &
+ OcnThickHadvInit
+
+ !--------------------------------------------------------------------
+ !
+ ! Private module variables
+ !
+ !--------------------------------------------------------------------
+
+
+!***********************************************************************
+
+contains
+
+!***********************************************************************
+!
+! routine OcnThickHadvTend
+!
+!> \brief Computes tendency term from horizontal advection of thickness
+!> \author Doug Jacobsen
+!> \date 15 September 2011
+!> \version SVN:$Id$
+!> \details
+!> This routine computes the horizontal advection tendency for
+!> thicknes based on current state and user choices of forcings.
+!
+!-----------------------------------------------------------------------
+
+ subroutine OcnThickHadvTend(grid, u, h_edge, tend, err)
+
+ !-----------------------------------------------------------------
+ !
+ ! input variables
+ !
+ !-----------------------------------------------------------------
+
+ real (kind=RKIND), dimension(:,:), intent(in) :: &
+ u !< Input: velocity
+
+ real (kind=RKIND), dimension(:,:), intent(in) :: &
+ h_edge !< Input: thickness at edge
+
+ type (mesh_type), intent(in) :: &
+ grid !< Input: grid information
+
+ !-----------------------------------------------------------------
+ !
+ ! input/output variables
+ !
+ !-----------------------------------------------------------------
+
+ real (kind=RKIND), dimension(:,:), intent(inout) :: &
+ tend !< Input/Output: velocity tendency
+
+ !-----------------------------------------------------------------
+ !
+ ! output variables
+ !
+ !-----------------------------------------------------------------
+
+ integer, intent(out) :: err
+
+ !-----------------------------------------------------------------
+ !
+ ! local variables
+ !
+ !-----------------------------------------------------------------
+
+ integer :: iEdge, nEdges, cell1, cell2, nVertLevels, k
+ integer :: iCell, nCells
+
+ integer, dimension(:), pointer :: maxLevelEdgeTop
+ integer, dimension(:,:), pointer :: cellsOnEdge
+
+ real (kind=RKIND) :: flux
+ real (kind=RKIND), dimension(:), pointer :: dvEdge, areaCell
+
+ !-----------------------------------------------------------------
+ !
+ ! call relevant routines for computing tendencies
+ ! note that the user can choose multiple options and the
+ ! tendencies will be added together
+ !
+ !-----------------------------------------------------------------
+
+ err = 0
+
+ nEdges = grid % nEdges
+ nCells = grid % nCells
+ nVertLevels = grid % nVertLevels
+
+ maxLevelEdgeTop => grid % maxLevelEdgeTop % array
+ cellsOnEdge => grid % cellsOnEdge % array
+ dvEdge => grid % dvEdge % array
+ areaCell => grid % areaCell % array
+
+ if (config_vert_grid_type.eq.'isopycnal') then
+
+ do iEdge=1,nEdges
+ cell1 = cellsOnEdge(1,iEdge)
+ cell2 = cellsOnEdge(2,iEdge)
+ do k=1,nVertLevels
+ flux = u(k,iEdge) * dvEdge(iEdge) * h_edge(k,iEdge)
+ tend(k,cell1) = tend(k,cell1) - flux
+ tend(k,cell2) = tend(k,cell2) + flux
+ end do
+ end do
+ do iCell=1,nCells
+ do k=1,nVertLevels
+ tend(k,iCell) = tend(k,iCell) / areaCell(iCell)
+ end do
+ end do
+
+ elseif (config_vert_grid_type.eq.'zlevel') then
+
+ do iEdge=1,nEdges
+ cell1 = cellsOnEdge(1,iEdge)
+ cell2 = cellsOnEdge(2,iEdge)
+ do k=1,min(1,maxLevelEdgeTop(iEdge))
+ flux = u(k,iEdge) * dvEdge(iEdge) * h_edge(k,iEdge)
+ tend(k,cell1) = tend(k,cell1) - flux
+ tend(k,cell2) = tend(k,cell2) + flux
+ end do
+ end do
+ do iCell=1,nCells
+ tend(1,iCell) = tend(1,iCell) / areaCell(iCell)
+ end do
+
+ endif ! config_vert_grid_type
+
+
+ !--------------------------------------------------------------------
+
+ end subroutine OcnThickHadvTend
+
+!***********************************************************************
+!
+! routine OcnThickHadvInit
+!
+!> \brief Initializes ocean forcings
+!> \author Doug Jacobsen
+!> \date 16 September 2011
+!> \version SVN:$Id$
+!> \details
+!> This routine initializes quantities related to forcings
+!> in the ocean. Since a multiple forcings are available,
+!> this routine primarily calls the
+!> individual init routines for each forcing.
+!
+!-----------------------------------------------------------------------
+
+
+ subroutine OcnThickHadvInit(err)
+
+ !--------------------------------------------------------------------
+
+ !-----------------------------------------------------------------
+ !
+ ! call individual init routines for each parameterization
+ !
+ !-----------------------------------------------------------------
+
+ integer, intent(out) :: err
+
+ err = 0
+
+ !--------------------------------------------------------------------
+
+ end subroutine OcnThickHadvInit
+
+!***********************************************************************
+
+end module OcnThickHadv
+
+!|||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||
Added: branches/ocean_projects/performance/src/core_ocean/module_OcnThickVadv.F
===================================================================
--- branches/ocean_projects/performance/src/core_ocean/module_OcnThickVadv.F         (rev 0)
+++ branches/ocean_projects/performance/src/core_ocean/module_OcnThickVadv.F        2011-09-16 19:27:37 UTC (rev 1006)
@@ -0,0 +1,165 @@
+!|||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||
+!
+! OcnThickVadv
+!
+!> \brief MPAS ocean vertical advection for thickness
+!> \author Doug Jacobsen
+!> \date 16 September 2011
+!> \version SVN:$Id:$
+!> \details
+!> This module contains the routine for computing
+!> tendencies for thickness from vertical advection
+!
+!-----------------------------------------------------------------------
+
+module OcnThickVadv
+
+ use grid_types
+ use configure
+
+ implicit none
+ private
+ save
+
+ !--------------------------------------------------------------------
+ !
+ ! Public parameters
+ !
+ !--------------------------------------------------------------------
+
+ !--------------------------------------------------------------------
+ !
+ ! Public member functions
+ !
+ !--------------------------------------------------------------------
+
+ public :: OcnThickVadvTend, &
+ OcnThickVadvInit
+
+ !--------------------------------------------------------------------
+ !
+ ! Private module variables
+ !
+ !--------------------------------------------------------------------
+
+
+!***********************************************************************
+
+contains
+
+!***********************************************************************
+!
+! routine OcnThickVadvTend
+!
+!> \brief Computes tendency term from vertical advection of thickness
+!> \author Doug Jacobsen
+!> \date 15 September 2011
+!> \version SVN:$Id$
+!> \details
+!> This routine computes the vertical advection tendency for
+!> thicknes based on current state and user choices of forcings.
+!
+!-----------------------------------------------------------------------
+
+ subroutine OcnThickVadvTend(grid, wTop, tend, err)
+
+ !-----------------------------------------------------------------
+ !
+ ! input variables
+ !
+ !-----------------------------------------------------------------
+
+ real (kind=RKIND), dimension(:,:), intent(in) :: &
+ wTop !< Input: vertical velocity on top layer
+
+ type (mesh_type), intent(in) :: &
+ grid !< Input: grid information
+
+ !-----------------------------------------------------------------
+ !
+ ! input/output variables
+ !
+ !-----------------------------------------------------------------
+
+ real (kind=RKIND), dimension(:,:), intent(inout) :: &
+ tend !< Input/Output: velocity tendency
+
+ !-----------------------------------------------------------------
+ !
+ ! output variables
+ !
+ !-----------------------------------------------------------------
+
+ integer, intent(out) :: err
+
+ !-----------------------------------------------------------------
+ !
+ ! local variables
+ !
+ !-----------------------------------------------------------------
+
+ integer :: iCell, nCells
+
+ !-----------------------------------------------------------------
+ !
+ ! call relevant routines for computing tendencies
+ ! note that the user can choose multiple options and the
+ ! tendencies will be added together
+ !
+ !-----------------------------------------------------------------
+
+ err = 0
+
+ nCells = grid % nCells
+
+ if (config_vert_grid_type.eq.'zlevel') then
+ do iCell=1,nCells
+ tend(1,iCell) = tend(1,iCell) + wTop(2,iCell)
+ end do
+ endif ! coordinate type
+
+
+ !--------------------------------------------------------------------
+
+ end subroutine OcnThickVadvTend
+
+!***********************************************************************
+!
+! routine OcnThickVadvInit
+!
+!> \brief Initializes ocean forcings
+!> \author Doug Jacobsen
+!> \date 16 September 2011
+!> \version SVN:$Id$
+!> \details
+!> This routine initializes quantities related to forcings
+!> in the ocean. Since a multiple forcings are available,
+!> this routine primarily calls the
+!> individual init routines for each forcing.
+!
+!-----------------------------------------------------------------------
+
+
+ subroutine OcnThickVadvInit(err)
+
+ !--------------------------------------------------------------------
+
+ !-----------------------------------------------------------------
+ !
+ ! call individual init routines for each parameterization
+ !
+ !-----------------------------------------------------------------
+
+ integer, intent(out) :: err
+
+ err = 0
+
+ !--------------------------------------------------------------------
+
+ end subroutine OcnThickVadvInit
+
+!***********************************************************************
+
+end module OcnThickVadv
+
+!|||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||
Modified: branches/ocean_projects/performance/src/core_ocean/module_time_integration.F
===================================================================
--- branches/ocean_projects/performance/src/core_ocean/module_time_integration.F        2011-09-16 18:29:23 UTC (rev 1005)
+++ branches/ocean_projects/performance/src/core_ocean/module_time_integration.F        2011-09-16 19:27:37 UTC (rev 1006)
@@ -8,6 +8,9 @@
use spline_interpolation
use timer
+ use OcnThickHadv
+ use OcnThickVadv
+
use OcnVelCoriolis
use OcnVelPressureGrad
use OcnVelVadv
@@ -1571,7 +1574,7 @@
type (mesh_type), intent(in) :: grid
integer :: iEdge, iCell, iVertex, k, cell1, cell2, &
- vertex1, vertex2, eoe, i, j
+ vertex1, vertex2, eoe, i, j, err
! mrp 110512 I just split compute_tend into compute_tend_u and compute_tend_h.
! Most of these variables can be removed, but at a later time.
@@ -1661,39 +1664,9 @@
! for z-level, only compute height tendency for top layer.
call timer_start("compute_tend_h-horiz adv")
- if (config_vert_grid_type.eq.'isopycnal') then
-
- do iEdge=1,nEdges
- cell1 = cellsOnEdge(1,iEdge)
- cell2 = cellsOnEdge(2,iEdge)
- do k=1,nVertLevels
- flux = u(k,iEdge) * dvEdge(iEdge) * h_edge(k,iEdge)
- tend_h(k,cell1) = tend_h(k,cell1) - flux
- tend_h(k,cell2) = tend_h(k,cell2) + flux
- end do
- end do
- do iCell=1,nCells
- do k=1,nVertLevels
- tend_h(k,iCell) = tend_h(k,iCell) / areaCell(iCell)
- end do
- end do
- elseif (config_vert_grid_type.eq.'zlevel') then
+ call OcnThickHadvTend(grid, u, h_edge, tend_h, err)
- do iEdge=1,nEdges
- cell1 = cellsOnEdge(1,iEdge)
- cell2 = cellsOnEdge(2,iEdge)
- do k=1,min(1,maxLevelEdgeTop(iEdge))
- flux = u(k,iEdge) * dvEdge(iEdge) * h_edge(k,iEdge)
- tend_h(k,cell1) = tend_h(k,cell1) - flux
- tend_h(k,cell2) = tend_h(k,cell2) + flux
- end do
- end do
- do iCell=1,nCells
- tend_h(1,iCell) = tend_h(1,iCell) / areaCell(iCell)
- end do
-
- endif ! config_vert_grid_type
call timer_stop("compute_tend_h-horiz adv")
!
@@ -1701,11 +1674,9 @@
!
! Vertical advection computed for top layer of a z grid only.
call timer_start("compute_tend_h-vert adv")
- if (config_vert_grid_type.eq.'zlevel') then
- do iCell=1,nCells
- tend_h(1,iCell) = tend_h(1,iCell) + wTop(2,iCell)
- end do
- endif ! coordinate type
+
+ call OcnThickVadvTend(grid, wTop, tend_h, err)
+
call timer_stop("compute_tend_h-vert adv")
call timer_stop("compute_tend_h")
</font>
</pre>