<p><b>jgraham@lanl.gov</b> 2012-11-15 12:28:25 -0700 (Thu, 15 Nov 2012)</p><p>Save tendencies to output files for calculation of energy transfer functions.<br>
Now all that remains is the weighted velocities at edges.<br>
</p><hr noshade><pre><font color="gray">Modified: branches/ocean_projects/turbulent_transfer/src/core_ocean/Registry
===================================================================
--- branches/ocean_projects/turbulent_transfer/src/core_ocean/Registry        2012-11-15 19:19:35 UTC (rev 2306)
+++ branches/ocean_projects/turbulent_transfer/src/core_ocean/Registry        2012-11-15 19:28:25 UTC (rev 2307)
@@ -287,6 +287,15 @@
 var persistent real    rhoDisplaced ( nVertLevels nCells Time ) 2 - rhoDisplaced state - -
 var persistent real    viscosity ( nVertLevels nEdges Time ) 2 o viscosity state - -
 
+% jpg 121101 diagnostic variables for energy transfer analysis
+var persistent real    uWeighted ( nVertLevels nEdges Time ) 2 o uWeighted state - -
+var persistent real    uCoriolisTend ( nVertLevels nEdges Time ) 2 o uCoriolisTend state - -
+var persistent real    uPressureGradTend ( nVertLevels nEdges Time ) 2 o uPressureGradTend state - -
+var persistent real    uHMixTend ( nVertLevels nEdges Time ) 2 o uHMixTend state - -
+var persistent real    uForcingTend ( nVertLevels nEdges Time ) 2 o uForcingTend state - -
+var persistent real    uVAdvTend ( nVertLevels nEdges Time ) 2 o uVAdvTend state - -
+var persistent real    uVMixTend ( nVertLevels nEdges Time ) 2 o uVMixTend state - -
+
 % Other diagnostic variables: neither read nor written to any files
 var persistent real    vh ( nVertLevels nEdges Time ) 2 - vh state - -
 var persistent real    circulation ( nVertLevels nVertices Time ) 2 - circulation state - -

Modified: branches/ocean_projects/turbulent_transfer/src/core_ocean/mpas_ocn_tendency.F
===================================================================
--- branches/ocean_projects/turbulent_transfer/src/core_ocean/mpas_ocn_tendency.F        2012-11-15 19:19:35 UTC (rev 2306)
+++ branches/ocean_projects/turbulent_transfer/src/core_ocean/mpas_ocn_tendency.F        2012-11-15 19:28:25 UTC (rev 2307)
@@ -263,6 +263,118 @@
 
    end subroutine ocn_tend_u!}}}
 
+! jpg 121101 save tendencies for post-processing transfer analysis
+   subroutine ocn_tend_u_diagnostic(s, d, grid)!{{{
+      implicit none
+
+      type (state_type), intent(inout) :: s !&lt; Input: State information
+      type (diagnostics_type), intent(in) :: d !&lt; Input: Diagnostic information
+      type (mesh_type), intent(in) :: grid !&lt; Input: Grid information
+
+      real (kind=RKIND), dimension(:,:), pointer :: &amp;
+        h_edge, h, u, rho, zMid, pressure, &amp;
+        circulation, vorticity, viscosity, ke, ke_edge, Vor_edge, &amp;
+        MontPot, wTop, divergence, vertViscTopOfEdge, tend_u_cor, &amp;
+        tend_u_press, tend_u_hmix, tend_u_force, tend_u_vadv, tend_u_vmix
+
+      real (kind=RKIND), dimension(:,:), pointer :: u_src
+
+      integer :: err
+
+      call mpas_timer_start(&quot;ocn_tend_u_diagnostic&quot;)
+
+      u           =&gt; s % u % array
+      rho         =&gt; s % rho % array
+      wTop        =&gt; s % wTop % array
+      zMid        =&gt; s % zMid % array
+      h_edge      =&gt; s % h_edge % array
+      viscosity   =&gt; s % viscosity % array
+      vorticity   =&gt; s % vorticity % array
+      divergence  =&gt; s % divergence % array
+      ke          =&gt; s % ke % array
+      ke_edge     =&gt; s % ke_edge % array
+      Vor_edge     =&gt; s % Vor_edge % array
+      MontPot     =&gt; s % MontPot % array
+      pressure    =&gt; s % pressure % array
+      vertViscTopOfEdge =&gt; d % vertViscTopOfEdge % array
+
+      tend_u_cor      =&gt; s % uCoriolisTend % array
+      tend_u_press      =&gt; s % uPressureGradTend % array
+      tend_u_hmix      =&gt; s % uHMixTend % array
+      tend_u_force      =&gt; s % uForcingTend % array
+      tend_u_vadv      =&gt; s % uVAdvTend % array
+      tend_u_vmix      =&gt; s % uVMixTend % array
+                  
+      u_src =&gt; grid % u_src % array
+
+      !
+      ! velocity tendency terms
+      !
+      ! jpg 121101 efficiency: I'm guessing these memory addresses need to be zeroed each call, but I'm not sure how time state 2 is setup
+      tend_u_cor(:,:) = 0.0
+      tend_u_press(:,:) = 0.0
+      tend_u_hmix(:,:) = 0.0
+      tend_u_force(:,:) = 0.0
+      tend_u_vadv(:,:) = 0.0
+      tend_u_vmix(:,:) = 0.0
+
+      !
+      ! velocity tendency: nonlinear Coriolis term and grad of kinetic energy
+      !
+
+      call mpas_timer_start(&quot;coriolis&quot;, .false., velCorTimer)
+      call ocn_vel_coriolis_tend(grid, Vor_edge, h_edge, u, ke, tend_u_cor, err)
+      call mpas_timer_stop(&quot;coriolis&quot;, velCorTimer)
+
+      !
+      ! velocity tendency: vertical advection term -w du/dz
+      !
+      call mpas_timer_start(&quot;vadv&quot;, .false., velVadvTimer)
+      call ocn_vel_vadv_tend(grid, u, h_edge, wtop, tend_u_vadv, err)
+      call mpas_timer_stop(&quot;vadv&quot;, velVadvTimer)
+
+      !
+      ! velocity tendency: pressure gradient
+      !
+      call mpas_timer_start(&quot;pressure grad&quot;, .false., velPgradTimer)
+      if (config_pressure_type.eq.'MontgomeryPotential') then
+          call ocn_vel_pressure_grad_tend(grid, MontPot,  zMid, rho, tend_u_press, err)
+      else
+          call ocn_vel_pressure_grad_tend(grid, pressure, zMid, rho, tend_u_press, err)
+      end if
+      call mpas_timer_stop(&quot;pressure grad&quot;, velPgradTimer)
+
+      !
+      ! velocity tendency: del2 dissipation, </font>
<font color="black">u_2 </font>
<font color="blue">abla^2 u
+      !   computed as </font>
<font color="black">u( </font>
<font color="black">abla divergence + k \times </font>
<font color="gray">abla vorticity )
+      !   strictly only valid for config_h_mom_eddy_visc2 == constant
+      !
+      call mpas_timer_start(&quot;hmix&quot;, .false., velHmixTimer)
+      call ocn_vel_hmix_tend(grid, divergence, vorticity, viscosity, tend_u_hmix, err)
+      call mpas_timer_stop(&quot;hmix&quot;, velHmixTimer)
+
+      !
+      ! velocity tendency: forcing and bottom drag
+      !
+      ! mrp 101115 note: in order to include flux boundary conditions, we will need to 
+      ! know the bottom edge with nonzero velocity and place the drag there.
+
+      call mpas_timer_start(&quot;forcings&quot;, .false., velForceTimer)
+      call ocn_vel_forcing_tend(grid, u, u_src, ke_edge, h_edge, tend_u_force, err)
+      call mpas_timer_stop(&quot;forcings&quot;, velForceTimer)
+
+      !
+      ! velocity tendency: vertical mixing d/dz( nu_v du/dz))
+      !
+      if (.not.config_implicit_vertical_mix) then
+          call mpas_timer_start(&quot;explicit vmix&quot;, .false., velExpVmixTimer)
+          call ocn_vel_vmix_tend_explicit(grid, u, h_edge, vertvisctopofedge, tend_u_vmix, err)
+          call mpas_timer_stop(&quot;explicit vmix&quot;, velExpVmixTimer)
+      endif
+      call mpas_timer_stop(&quot;ocn_tend_u_diagnostic&quot;)
+
+      end subroutine ocn_tend_u_diagnostic !}}}
+
 !***********************************************************************
 !
 !  routine ocn_tendSalar

Modified: branches/ocean_projects/turbulent_transfer/src/core_ocean/mpas_ocn_time_integration_rk4.F
===================================================================
--- branches/ocean_projects/turbulent_transfer/src/core_ocean/mpas_ocn_time_integration_rk4.F        2012-11-15 19:19:35 UTC (rev 2306)
+++ branches/ocean_projects/turbulent_transfer/src/core_ocean/mpas_ocn_time_integration_rk4.F        2012-11-15 19:28:25 UTC (rev 2307)
@@ -167,6 +167,12 @@
            endif
 
            call ocn_tend_scalar(block % tend, block % provis, block % diagnostics, block % mesh, dt)
+
+! jpg 121101 save tendencies for post-processing transfer analysis
+           if (rk_step == 1) then
+              call ocn_tend_u_diagnostic(block % state % time_levs(2) % state, block % diagnostics, block % mesh)
+           endif
+
            block =&gt; block % next
         end do
         call mpas_timer_stop(&quot;RK4-tendency computations&quot;)

</font>
</pre>