<p><b>dwj07@fsu.edu</b> 2011-09-16 12:29:23 -0600 (Fri, 16 Sep 2011)</p><p><br>
        -- BRANCH COMMIT --<br>
<br>
        Adding modules for forcings in ocean core.<br>
<br>
        Wind stress and bottom drag.<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 16:39:49 UTC (rev 1004)
+++ branches/ocean_projects/performance/src/core_ocean/Makefile        2011-09-16 18:29:23 UTC (rev 1005)
@@ -8,6 +8,9 @@
         module_OcnVelHmix.o \
         module_OcnVelHmixDel2.o \
         module_OcnVelHmixDel4.o \
+         module_OcnVelForcing.o \
+         module_OcnVelForcingWindStress.o \
+         module_OcnVelForcingBottomDrag.o \
         module_OcnVelPressureGrad.o \
module_time_integration.o \
module_global_diagnostics.o
@@ -35,6 +38,12 @@
module_OcnVelHmixDel4.o:
+module_OcnVelForcing.o: module_OcnVelForcingWindStress.o module_OcnVelForcingBottomDrag.o
+
+module_OcnVelForcingWindStress.o:
+
+module_OcnVelForcingBottomDrag.o:
+
module_OcnVelCoriolis.o:
module_mpas_core.o: module_advection.o \
@@ -45,6 +54,9 @@
                                        module_OcnVelHmix.o \
                                        module_OcnVelHmixDel2.o \
                                        module_OcnVelHmixDel4.o \
+                                        module_OcnVelForcing.o \
+                                        module_OcnVelForcingWindStress.o \
+                                        module_OcnVelForcingBottomDrag.o \
                                        module_OcnVelPressureGrad.o \
                                        module_time_integration.o
Added: branches/ocean_projects/performance/src/core_ocean/module_OcnVelForcing.F
===================================================================
--- branches/ocean_projects/performance/src/core_ocean/module_OcnVelForcing.F         (rev 0)
+++ branches/ocean_projects/performance/src/core_ocean/module_OcnVelForcing.F        2011-09-16 18:29:23 UTC (rev 1005)
@@ -0,0 +1,180 @@
+!|||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||
+!
+! OcnVelForcing
+!
+!> \brief MPAS ocean forcing driver
+!> \author Doug Jacobsen
+!> \date 16 September 2011
+!> \version SVN:$Id:$
+!> \details
+!> This module contains the main driver routine for computing
+!> tendencies from forcings.
+!
+!-----------------------------------------------------------------------
+
+module OcnVelForcing
+
+ use grid_types
+ use configure
+
+ use OcnVelForcingWindStress
+ use OcnVelForcingBottomDrag
+
+ implicit none
+ private
+ save
+
+ !--------------------------------------------------------------------
+ !
+ ! Public parameters
+ !
+ !--------------------------------------------------------------------
+
+ !--------------------------------------------------------------------
+ !
+ ! Public member functions
+ !
+ !--------------------------------------------------------------------
+
+ public :: OcnVelForcingTend, &
+ OcnVelForcingInit
+
+ !--------------------------------------------------------------------
+ !
+ ! Private module variables
+ !
+ !--------------------------------------------------------------------
+
+
+!***********************************************************************
+
+contains
+
+!***********************************************************************
+!
+! routine OcnVelForcingTend
+!
+!> \brief Computes tendency term from forcings
+!> \author Doug Jacobsen
+!> \date 15 September 2011
+!> \version SVN:$Id$
+!> \details
+!> This routine computes the forcing tendency for momentum
+!> based on current state and user choices of forcings.
+!> Multiple forcings may be chosen and added together. These
+!> tendencies are generally computed by calling the specific routine
+!> for the chosen forcing, so this routine is primarily a
+!> driver for managing these choices.
+!
+!-----------------------------------------------------------------------
+
+ subroutine OcnVelForcingTend(grid, u, u_src, ke_edge, h_edge, tend, err)
+
+ !-----------------------------------------------------------------
+ !
+ ! input variables
+ !
+ !-----------------------------------------------------------------
+
+ real (kind=RKIND), dimension(:,:), intent(in) :: &
+ u !< Input: velocity
+
+ real (kind=RKIND), dimension(:,:), intent(in) :: &
+ u_src !< Input: wind stress
+
+ real (kind=RKIND), dimension(:,:), intent(in) :: &
+ ke_edge !< Input: kinetic energy at edge
+
+ 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 :: err1, err2
+
+ !-----------------------------------------------------------------
+ !
+ ! call relevant routines for computing tendencies
+ ! note that the user can choose multiple options and the
+ ! tendencies will be added together
+ !
+ !-----------------------------------------------------------------
+
+ call OcnVelForcingWindStressTend(grid, u_src, h_edge, tend, err1)
+ call OcnVelForcingBottomDragTend(grid, u, ke_edge, h_edge, tend, err2)
+
+ err = err1 .or. err2
+
+ !--------------------------------------------------------------------
+
+ end subroutine OcnVelForcingTend
+
+!***********************************************************************
+!
+! routine OcnVelForcingInit
+!
+!> \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 OcnVelForcingInit(err)
+
+ !--------------------------------------------------------------------
+
+ !-----------------------------------------------------------------
+ !
+ ! call individual init routines for each parameterization
+ !
+ !-----------------------------------------------------------------
+
+ integer, intent(out) :: err
+
+ integer :: err1, err2
+
+ call OcnVelForcingWindStressInit(err1)
+ call OcnVelForcingBottomDragInit(err2)
+
+ err = err1 .or. err2
+
+ !--------------------------------------------------------------------
+
+ end subroutine OcnVelForcingInit
+
+!***********************************************************************
+
+end module OcnVelForcing
+
+!|||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||
Added: branches/ocean_projects/performance/src/core_ocean/module_OcnVelForcingBottomDrag.F
===================================================================
--- branches/ocean_projects/performance/src/core_ocean/module_OcnVelForcingBottomDrag.F         (rev 0)
+++ branches/ocean_projects/performance/src/core_ocean/module_OcnVelForcingBottomDrag.F        2011-09-16 18:29:23 UTC (rev 1005)
@@ -0,0 +1,201 @@
+!|||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||
+!
+! OcnVelForcingBottomDrag
+!
+!> \brief MPAS ocean bottom drag
+!> \author Doug Jacobsen
+!> \date 16 September 2011
+!> \version SVN:$Id:$
+!> \details
+!> This module contains the routine for computing
+!> tendencies from bottom drag.
+!
+!-----------------------------------------------------------------------
+
+module OcnVelForcingBottomDrag
+
+ use grid_types
+ use configure
+
+ implicit none
+ private
+ save
+
+ !--------------------------------------------------------------------
+ !
+ ! Public parameters
+ !
+ !--------------------------------------------------------------------
+
+ !--------------------------------------------------------------------
+ !
+ ! Public member functions
+ !
+ !--------------------------------------------------------------------
+
+ public :: OcnVelForcingBottomDragTend, &
+ OcnVelForcingBottomDragInit
+
+ !--------------------------------------------------------------------
+ !
+ ! Private module variables
+ !
+ !--------------------------------------------------------------------
+
+ logical :: bottomDragOn
+ real (kind=RKIND) :: bottomDragCoef
+
+
+!***********************************************************************
+
+contains
+
+!***********************************************************************
+!
+! routine OcnVelForcingBottomDragTend
+!
+!> \brief Computes tendency term from bottom drag
+!> \author Doug Jacobsen
+!> \date 15 September 2011
+!> \version SVN:$Id$
+!> \details
+!> This routine computes the bottom drag tendency for momentum
+!> based on current state.
+!
+!-----------------------------------------------------------------------
+
+ subroutine OcnVelForcingBottomDragTend(grid, u, ke_edge, h_edge, tend, err)
+
+ !-----------------------------------------------------------------
+ !
+ ! input variables
+ !
+ !-----------------------------------------------------------------
+
+ real (kind=RKIND), dimension(:,:), intent(in) :: &
+ u !< Input: velocity
+
+ real (kind=RKIND), dimension(:,:), intent(in) :: &
+ ke_edge !< Input: kinetic energy at edge
+
+ 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, nEdgesSolve, k
+ integer, dimension(:), pointer :: maxLevelEdgeTop
+
+ !-----------------------------------------------------------------
+ !
+ ! 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.bottomDragOn) return
+
+ nEdgesSolve = grid % nEdgesSolve
+ maxLevelEdgeTop => grid % maxLevelEdgeTop % array
+
+ do iEdge=1,grid % nEdgesSolve
+
+ k = maxLevelEdgeTop(iEdge)
+
+ ! efficiency note: it would be nice to avoid this
+ ! if within a do. This could be done with
+ ! k = max(maxLevelEdgeTop(iEdge),1)
+ ! and then tend_u(1,iEdge) is just not used for land cells.
+
+ if (k>0) then
+ ! bottom drag is the same as POP:
+ ! -c |u| u where c is unitless and 1.0e-3.
+ ! see POP Reference guide, section 3.4.4.
+
+ tend(k,iEdge) = tend(k,iEdge) &
+ -bottomDragCoef*u(k,iEdge) &
+ *sqrt(2.0*ke_edge(k,iEdge))/h_edge(k,iEdge)
+
+ endif
+
+ enddo
+
+
+
+ !--------------------------------------------------------------------
+
+ end subroutine OcnVelForcingBottomDragTend
+
+!***********************************************************************
+!
+! routine OcnVelForcingBottomDragInit
+!
+!> \brief Initializes ocean bottom drag
+!> \author Doug Jacobsen
+!> \date 16 September 2011
+!> \version SVN:$Id$
+!> \details
+!> This routine initializes quantities related to bottom drag
+!> in the ocean.
+!
+!-----------------------------------------------------------------------
+
+
+ subroutine OcnVelForcingBottomDragInit(err)
+
+ !--------------------------------------------------------------------
+
+ !-----------------------------------------------------------------
+ !
+ ! call individual init routines for each parameterization
+ !
+ !-----------------------------------------------------------------
+
+ integer, intent(out) :: err
+
+
+ err = 0
+
+ bottomDragOn = .false.
+
+ if (.not.config_implicit_vertical_mix) then
+ bottomDragOn = .true.
+ bottomDragCoef = config_bottom_drag_coeff
+ endif
+
+ !--------------------------------------------------------------------
+
+ end subroutine OcnVelForcingBottomDragInit
+
+!***********************************************************************
+
+end module OcnVelForcingBottomDrag
+
+!|||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||
Added: branches/ocean_projects/performance/src/core_ocean/module_OcnVelForcingWindStress.F
===================================================================
--- branches/ocean_projects/performance/src/core_ocean/module_OcnVelForcingWindStress.F         (rev 0)
+++ branches/ocean_projects/performance/src/core_ocean/module_OcnVelForcingWindStress.F        2011-09-16 18:29:23 UTC (rev 1005)
@@ -0,0 +1,190 @@
+!|||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||
+!
+! OcnVelForcingWindStress
+!
+!> \brief MPAS ocean wind stress
+!> \author Doug Jacobsen
+!> \date 16 September 2011
+!> \version SVN:$Id:$
+!> \details
+!> This module contains the routine for computing
+!> tendencies from wind stress.
+!
+!-----------------------------------------------------------------------
+
+module OcnVelForcingWindStress
+
+ use grid_types
+ use configure
+
+ implicit none
+ private
+ save
+
+ !--------------------------------------------------------------------
+ !
+ ! Public parameters
+ !
+ !--------------------------------------------------------------------
+
+ !--------------------------------------------------------------------
+ !
+ ! Public member functions
+ !
+ !--------------------------------------------------------------------
+
+ public :: OcnVelForcingWindStressTend, &
+ OcnVelForcingWindStressInit
+
+ !--------------------------------------------------------------------
+ !
+ ! Private module variables
+ !
+ !--------------------------------------------------------------------
+
+ logical :: windStressOn
+ real (kind=RKIND) :: rho_ref
+
+
+!***********************************************************************
+
+contains
+
+!***********************************************************************
+!
+! routine OcnVelForcingWindStressTend
+!
+!> \brief Computes tendency term from wind stress
+!> \author Doug Jacobsen
+!> \date 15 September 2011
+!> \version SVN:$Id$
+!> \details
+!> This routine computes the wind stress tendency for momentum
+!> based on current state.
+!
+!-----------------------------------------------------------------------
+
+ subroutine OcnVelForcingWindStressTend(grid, u_src, h_edge, tend, err)
+
+ !-----------------------------------------------------------------
+ !
+ ! input variables
+ !
+ !-----------------------------------------------------------------
+
+ real (kind=RKIND), dimension(:,:), intent(in) :: &
+ u_src !< Input: wind stress
+
+ 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, nEdgesSolve, k
+ integer, dimension(:), pointer :: maxLevelEdgeTop
+
+
+ !-----------------------------------------------------------------
+ !
+ ! 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.windStressOn) return
+
+ nEdgesSolve = grid % nEdgesSolve
+ maxLevelEdgeTop => grid % maxLevelEdgeTop % array
+
+ do iEdge=1,nEdgesSolve
+
+ k = maxLevelEdgeTop(iEdge)
+
+ ! efficiency note: it would be nice to avoid this
+ ! if within a do. This could be done with
+ ! k = max(maxLevelEdgeTop(iEdge),1)
+ ! and then tend_u(1,iEdge) is just not used for land cells.
+
+ if (k>0) then
+ ! forcing in top layer only
+ tend(1,iEdge) = tend(1,iEdge) &
+ + u_src(1,iEdge)/rho_ref/h_edge(1,iEdge)
+ endif
+
+ enddo
+
+
+ !--------------------------------------------------------------------
+
+ end subroutine OcnVelForcingWindStressTend
+
+!***********************************************************************
+!
+! routine OcnVelForcingWindStressInit
+!
+!> \brief Initializes ocean wind stress forcing
+!> \author Doug Jacobsen
+!> \date 16 September 2011
+!> \version SVN:$Id$
+!> \details
+!> This routine initializes quantities related to wind stress
+!> in the ocean.
+!
+!-----------------------------------------------------------------------
+
+
+ subroutine OcnVelForcingWindStressInit(err)
+
+ !--------------------------------------------------------------------
+
+ !-----------------------------------------------------------------
+ !
+ ! call individual init routines for each parameterization
+ !
+ !-----------------------------------------------------------------
+
+ integer, intent(out) :: err
+
+
+ windStressOn = .true.
+ rho_ref = 1000.0
+
+ err = 0
+
+
+ !--------------------------------------------------------------------
+
+ end subroutine OcnVelForcingWindStressInit
+
+!***********************************************************************
+
+end module OcnVelForcingWindStress
+
+!|||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||||
Modified: branches/ocean_projects/performance/src/core_ocean/module_mpas_core.F
===================================================================
--- branches/ocean_projects/performance/src/core_ocean/module_mpas_core.F        2011-09-16 16:39:49 UTC (rev 1004)
+++ branches/ocean_projects/performance/src/core_ocean/module_mpas_core.F        2011-09-16 18:29:23 UTC (rev 1005)
@@ -8,6 +8,7 @@
use OcnVelPressureGrad
use OcnVelVadv
use OcnVelHmix
+ use OcnVelForcing
type (io_output_object) :: restart_obj
integer :: restart_frame
@@ -78,6 +79,7 @@
call OcnVelPressureGradInit(err)
call OcnVelVadvInit(err)
call OcnVelHmixInit(err)
+ call OcnVelForcingInit(err)
! mrp 100316 In order for this to work, we need to pass domain % dminfo as an
! input arguement into mpas_init. Ask about that later. For now, there will be
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 16:39:49 UTC (rev 1004)
+++ branches/ocean_projects/performance/src/core_ocean/module_time_integration.F        2011-09-16 18:29:23 UTC (rev 1005)
@@ -12,6 +12,7 @@
use OcnVelPressureGrad
use OcnVelVadv
use OcnVelHmix
+ use OcnVelForcing
contains
@@ -1870,34 +1871,9 @@
! know the bottom edge with nonzero velocity and place the drag there.
call timer_start("compute_tend_u-forcings")
- do iEdge=1,grid % nEdgesSolve
- k = maxLevelEdgeTop(iEdge)
+ call OcnVelForcingTend(grid, u, u_src, ke_edge, h_edge, tend_u, err)
- ! efficiency note: it would be nice to avoid this
- ! if within a do. This could be done with
- ! k = max(maxLevelEdgeTop(iEdge),1)
- ! and then tend_u(1,iEdge) is just not used for land cells.
-
- if (k>0) then
-
- ! forcing in top layer only
- tend_u(1,iEdge) = tend_u(1,iEdge) &
- + u_src(1,iEdge)/rho_ref/h_edge(1,iEdge)
-
- ! bottom drag is the same as POP:
- ! -c |u| u where c is unitless and 1.0e-3.
- ! see POP Reference guide, section 3.4.4.
-
- if (.not.config_implicit_vertical_mix) then
- tend_u(k,iEdge) = tend_u(k,iEdge) &
- - config_bottom_drag_coeff*u(k,iEdge) &
- *sqrt(2.0*ke_edge(k,iEdge))/h_edge(k,iEdge)
- end if
-
- endif
-
- enddo
call timer_stop("compute_tend_u-forcings")
!
</font>
</pre>