<p><b>dwj07@fsu.edu</b> 2012-01-06 11:21:56 -0700 (Fri, 06 Jan 2012)</p><p><br>
        -- TRUNK COMMIT --<br>
<br>
        Merging the time_averaging branch to the trunk. It only affects the ocean core.<br>
<br>
<br>
        Now fields can be averaged over time to give more information from a run, and allow comparisons with other ocean models.<br>
</p><hr noshade><pre><font color="gray">
Property changes on: trunk/mpas/src/core_ocean
___________________________________________________________________
Added: svn:mergeinfo
   + /branches/cam_mpas_nh/src/core_ocean:1260-1270
/branches/ocean_projects/imp_vert_mix_mrp/src/core_ocean:754-986
/branches/ocean_projects/split_explicit_mrp/src/core_ocean:1134-1138
/branches/ocean_projects/split_explicit_timestepping/src/core_ocean:1044-1097
/branches/ocean_projects/time_averaging/src/core_ocean:1271-1305
/branches/ocean_projects/vert_adv_mrp/src/core_ocean:704-745
/branches/source_renaming/src/core_ocean:1082-1113
/branches/time_manager/src/core_ocean:924-962

Modified: trunk/mpas/src/core_ocean/Makefile
===================================================================
--- trunk/mpas/src/core_ocean/Makefile        2012-01-06 18:02:59 UTC (rev 1306)
+++ trunk/mpas/src/core_ocean/Makefile        2012-01-06 18:21:56 UTC (rev 1307)
@@ -41,7 +41,8 @@
            mpas_ocn_equation_of_state.o \
            mpas_ocn_equation_of_state_jm.o \
            mpas_ocn_equation_of_state_linear.o \
-       mpas_ocn_global_diagnostics.o
+       mpas_ocn_global_diagnostics.o \
+           mpas_ocn_time_average.o
 
 
 all: core_hyd
@@ -59,10 +60,12 @@
 
 mpas_ocn_time_integration_split.o:
 
-mpas_ocn_tendency.o:
+mpas_ocn_tendency.o: mpas_ocn_time_average.o
 
 mpas_ocn_global_diagnostics.o: 
 
+mpas_ocn_time_average.o:
+
 mpas_ocn_thick_hadv.o:
 
 mpas_ocn_thick_vadv.o:
@@ -172,7 +175,8 @@
                                           mpas_ocn_equation_of_state.o \
                                           mpas_ocn_equation_of_state_jm.o \
                                           mpas_ocn_equation_of_state_linear.o \
-                                          mpas_ocn_global_diagnostics.o
+                                          mpas_ocn_global_diagnostics.o \
+                                          mpas_ocn_time_average.o
 
 clean:
         $(RM) *.o *.mod *.f90 libdycore.a

Modified: trunk/mpas/src/core_ocean/Registry
===================================================================
--- trunk/mpas/src/core_ocean/Registry        2012-01-06 18:02:59 UTC (rev 1306)
+++ trunk/mpas/src/core_ocean/Registry        2012-01-06 18:21:56 UTC (rev 1307)
@@ -258,3 +258,13 @@
 var persistent real    RiTopOfEdge ( nVertLevelsP1 nEdges Time ) 1 - RiTopOfEdge diagnostics - -
 var persistent real    vertViscTopOfEdge ( nVertLevelsP1 nEdges Time ) 1 - vertViscTopOfEdge diagnostics - -
 var persistent real    vertDiffTopOfCell ( nVertLevelsP1 nCells Time ) 1 - vertDiffTopOfCell diagnostics - -
+
+var persistent real    nAccumulate ( Time ) 1 o nAccumulate state - -
+var persistent real    acc_ssh ( nCells Time ) 1 o acc_ssh state - -
+var persistent real    acc_sshVar ( nCells Time ) 1 o acc_sshVar state - -
+var persistent real    acc_uReconstructZonal ( nVertLevels nCells Time ) 1 o acc_uReconstructZonal state - -
+var persistent real    acc_uReconstructMeridional ( nVertLevels nCells Time ) 1 o acc_uReconstructMeridional state - -
+var persistent real    acc_uReconstructZonalVar ( nVertLevels nCells Time ) 1 o acc_uReconstructZonalVar state - -
+var persistent real    acc_uReconstructMeridionalVar ( nVertLevels nCells Time ) 1 o acc_uReconstructMeridionalVar state - -
+var persistent real           acc_u ( nVertLevels nEdges Time ) 1 o acc_u state - -
+var persistent real           acc_uVar ( nVertLevels nEdges Time ) 1 o acc_uVar state - -

Modified: trunk/mpas/src/core_ocean/mpas_ocn_mpas_core.F
===================================================================
--- trunk/mpas/src/core_ocean/mpas_ocn_mpas_core.F        2012-01-06 18:02:59 UTC (rev 1306)
+++ trunk/mpas/src/core_ocean/mpas_ocn_mpas_core.F        2012-01-06 18:21:56 UTC (rev 1307)
@@ -23,6 +23,8 @@
 
    use ocn_vmix
 
+   use ocn_time_average
+
    type (io_output_object) :: restart_obj
 
    integer :: current_outfile_frames
@@ -209,6 +211,7 @@
       real (kind=RKIND), intent(in) :: dt
       integer :: i, iEdge, iCell, k
    
+      call ocn_time_average_init(block % state % time_levs(1) % state)
    
       call ocn_diagnostic_solve(dt, block % state % time_levs(1) % state, mesh)
 
@@ -301,7 +304,13 @@
       write(0,*) 'Initial time ', timeStamp
 
       call write_output_frame(output_obj, output_frame, domain)
+      block_ptr =&gt; domain % blocklist
 
+      do while(associated(block_ptr))
+        call ocn_time_average_init(block_ptr % state % time_levs(1) % state)
+        block_ptr =&gt; block_ptr % next
+      end do
+
       ! During integration, time level 1 stores the model state at the beginning of the
       !   time step, and time level 2 stores the state advanced dt in time by timestep(...)
       itimestep = 0
@@ -328,7 +337,18 @@
                call mpas_output_state_finalize(output_obj, domain % dminfo)
                call mpas_output_state_init(output_obj, domain, &quot;OUTPUT&quot;, trim(timeStamp))
             end if
+
+            block_ptr =&gt; domain % blocklist
+            do while (associated(block_ptr))
+                call ocn_time_average_normalize(block_ptr % state % time_levs(1) % state)
+                block_ptr =&gt; block_ptr % next
+            end do
             call write_output_frame(output_obj, output_frame, domain)
+            block_ptr =&gt; domain % blocklist
+            do while (associated(block_ptr))
+                call ocn_time_average_init(block_ptr % state % time_levs(1) % state)
+                block_ptr =&gt; block_ptr % next
+            end do
          end if
 
          if (mpas_is_alarm_ringing(clock, restartAlarmID, ierr=ierr)) then

Modified: trunk/mpas/src/core_ocean/mpas_ocn_tendency.F
===================================================================
--- trunk/mpas/src/core_ocean/mpas_ocn_tendency.F        2012-01-06 18:02:59 UTC (rev 1306)
+++ trunk/mpas/src/core_ocean/mpas_ocn_tendency.F        2012-01-06 18:21:56 UTC (rev 1307)
@@ -38,6 +38,8 @@
    use ocn_equation_of_state
    use ocn_vmix
 
+   use ocn_time_average
+
    implicit none
    private
    save

Copied: trunk/mpas/src/core_ocean/mpas_ocn_time_average.F (from rev 1305, branches/ocean_projects/time_averaging/src/core_ocean/mpas_ocn_time_average.F)
===================================================================
--- trunk/mpas/src/core_ocean/mpas_ocn_time_average.F                                (rev 0)
+++ trunk/mpas/src/core_ocean/mpas_ocn_time_average.F        2012-01-06 18:21:56 UTC (rev 1307)
@@ -0,0 +1,129 @@
+module ocn_time_average
+
+    use mpas_grid_types
+
+    implicit none
+    save
+    public
+
+    contains 
+
+    subroutine ocn_time_average_init(state)!{{{
+        type (state_type), intent(inout) :: state
+
+        real, pointer :: nAccumulate
+
+        real (kind=RKIND), dimension(:), pointer :: acc_ssh, acc_sshVar
+        real (kind=RKIND), dimension(:,:), pointer :: acc_uReconstructZonal, acc_uReconstructMeridional, acc_uReconstructZonalVar, acc_uReconstructMeridionalVar
+        real (kind=RKIND), dimension(:,:), pointer :: acc_u, acc_uVar
+
+        nAccumulate =&gt; state % nAccumulate % scalar
+
+        acc_ssh =&gt; state % acc_ssh % array
+        acc_sshVar =&gt; state % acc_sshVar % array
+        acc_uReconstructZonal =&gt; state % acc_uReconstructZonal % array
+        acc_uReconstructMeridional =&gt; state % acc_uReconstructMeridional % array
+        acc_uReconstructZonalVar =&gt; state % acc_uReconstructZonalVar % array
+        acc_uReconstructMeridionalVar =&gt; state % acc_uReconstructMeridionalVar % array
+        acc_u =&gt; state % acc_u % array
+        acc_uVar =&gt; state % acc_uVar % array
+
+        nAccumulate = 0
+
+        acc_ssh = 0.0
+        acc_sshVar = 0.0
+        acc_uReconstructZonal = 0.0
+        acc_uReconstructMeridional = 0.0
+        acc_uReconstructZonalVar = 0.0
+        acc_uReconstructMeridionalVar = 0.0
+        acc_u = 0.0
+        acc_uVar = 0.0
+
+    end subroutine ocn_time_average_init!}}}
+
+    subroutine ocn_time_average_accumulate(state, old_state)!{{{
+        type (state_type), intent(inout) :: state
+        type (state_type), intent(in) :: old_state
+
+        real, pointer :: nAccumulate, old_nAccumulate
+
+        real (kind=RKIND), dimension(:), pointer :: ssh
+        real (kind=RKIND), dimension(:,:), pointer :: uReconstructZonal, uReconstructMeridional, u
+
+        real (kind=RKIND), dimension(:,:), pointer :: acc_u, acc_uVar
+        real (kind=RKIND), dimension(:,:), pointer :: acc_uReconstructZonal, acc_uReconstructMeridional, acc_uReconstructZonalVar, acc_uReconstructMeridionalVar
+        real (kind=RKIND), dimension(:), pointer :: acc_ssh, acc_sshVar
+
+        real (kind=RKIND), dimension(:,:), pointer :: old_acc_u, old_acc_uVar
+        real (kind=RKIND), dimension(:,:), pointer :: old_acc_uReconstructZonal, old_acc_uReconstructMeridional, old_acc_uReconstructZonalVar, old_acc_uReconstructMeridionalVar
+        real (kind=RKIND), dimension(:), pointer :: old_acc_ssh, old_acc_sshVar
+
+        old_nAccumulate =&gt; old_state % nAccumulate  % scalar
+        nAccumulate =&gt; state % nAccumulate  % scalar
+
+        ssh =&gt; state % ssh % array
+        uReconstructZonal =&gt; state % uReconstructZonal % array
+        uReconstructMeridional =&gt; state % uReconstructMeridional % array
+        u =&gt; state % u % array
+
+        acc_ssh =&gt; state % acc_ssh % array
+        acc_sshVar =&gt; state % acc_sshVar % array
+        acc_uReconstructZonal =&gt; state % acc_uReconstructZonal % array
+        acc_uReconstructMeridional =&gt; state % acc_uReconstructMeridional % array
+        acc_uReconstructZonalVar =&gt; state % acc_uReconstructZonalVar % array
+        acc_uReconstructMeridionalVar =&gt; state % acc_uReconstructMeridionalVar % array
+        acc_u =&gt; state % acc_u % array
+        acc_uVar =&gt; state % acc_uVar % array
+
+        old_acc_ssh =&gt; old_state % acc_ssh % array
+        old_acc_sshVar =&gt; old_state % acc_sshVar % array
+        old_acc_uReconstructZonal =&gt; old_state % acc_uReconstructZonal % array
+        old_acc_uReconstructMeridional =&gt; old_state % acc_uReconstructMeridional % array
+        old_acc_uReconstructZonalVar =&gt; old_state % acc_uReconstructZonalVar % array
+        old_acc_uReconstructMeridionalVar =&gt; old_state % acc_uReconstructMeridionalVar % array
+        old_acc_u =&gt; old_state % acc_u % array
+        old_acc_uVar =&gt; old_state % acc_uVar % array
+
+        acc_ssh = old_acc_ssh + ssh
+        acc_sshVar = old_acc_sshVar + ssh**2
+        acc_uReconstructZonal = old_acc_uReconstructZonal + uReconstructZonal
+        acc_uReconstructMeridional = old_acc_uReconstructMeridional + uReconstructMeridional
+        acc_uReconstructZonalVar = old_acc_uReconstructZonalVar + uReconstructZonal**2
+        acc_uReconstructMeridionalVar = old_acc_uReconstructMeridionalVar + uReconstructMeridional**2
+        acc_u = old_acc_u + u
+        acc_uVar = old_acc_uVar + u**2
+
+        nAccumulate = old_nAccumulate + 1
+    end subroutine ocn_time_average_accumulate!}}}
+
+    subroutine ocn_time_average_normalize(state)!{{{
+        type (state_type), intent(inout) :: state
+
+        real, pointer :: nAccumulate
+
+        real (kind=RKIND), dimension(:), pointer :: acc_ssh, acc_sshVar
+        real (kind=RKIND), dimension(:,:), pointer :: acc_uReconstructZonal, acc_uReconstructMeridional, acc_uReconstructZonalVar, acc_uReconstructMeridionalVar
+        real (kind=RKIND), dimension(:,:), pointer :: acc_u, acc_uVar
+
+        nAccumulate =&gt; state % nAccumulate  % scalar
+
+        acc_ssh =&gt; state % acc_ssh % array
+        acc_sshVar =&gt; state % acc_sshVar % array
+        acc_uReconstructZonal =&gt; state % acc_uReconstructZonal % array
+        acc_uReconstructMeridional =&gt; state % acc_uReconstructMeridional % array
+        acc_uReconstructZonalVar =&gt; state % acc_uReconstructZonalVar % array
+        acc_uReconstructMeridionalVar =&gt; state % acc_uReconstructMeridionalVar % array
+        acc_u =&gt; state % acc_u % array
+        acc_uVar =&gt; state % acc_uVar % array
+
+        acc_ssh = acc_ssh / nAccumulate
+        acc_sshVar = acc_sshVar / nAccumulate
+        acc_uReconstructZonal = acc_uReconstructZonal / nAccumulate
+        acc_uReconstructMeridional = acc_uReconstructMeridional / nAccumulate
+        acc_uReconstructZonalVar = acc_uReconstructZonalVar / nAccumulate
+        acc_uReconstructMeridionalVar = acc_uReconstructMeridionalVar / nAccumulate
+        acc_u = acc_u / nAccumulate
+        acc_uVar = acc_uVar / nAccumulate
+    end subroutine ocn_time_average_normalize!}}}
+
+end module ocn_time_average

Modified: trunk/mpas/src/core_ocean/mpas_ocn_time_integration_rk4.F
===================================================================
--- trunk/mpas/src/core_ocean/mpas_ocn_time_integration_rk4.F        2012-01-06 18:02:59 UTC (rev 1306)
+++ trunk/mpas/src/core_ocean/mpas_ocn_time_integration_rk4.F        2012-01-06 18:21:56 UTC (rev 1307)
@@ -24,7 +24,8 @@
    use ocn_tendency
 
    use ocn_equation_of_state
-   use ocn_Vmix
+   use ocn_vmix
+   use ocn_time_average
 
    implicit none
    private
@@ -356,6 +357,8 @@
                           block % state % time_levs(2) % state % uReconstructMeridional % array    &amp;
                          )
 
+         call ocn_time_average_accumulate(block % state % time_levs(2) % state, block % state % time_levs(1) % state)
+
          block =&gt; block % next
       end do
       call mpas_timer_stop(&quot;RK4-cleaup phase&quot;)

Modified: trunk/mpas/src/core_ocean/mpas_ocn_time_integration_split.F
===================================================================
--- trunk/mpas/src/core_ocean/mpas_ocn_time_integration_split.F        2012-01-06 18:02:59 UTC (rev 1306)
+++ trunk/mpas/src/core_ocean/mpas_ocn_time_integration_split.F        2012-01-06 18:21:56 UTC (rev 1307)
@@ -27,8 +27,8 @@
 
    use ocn_equation_of_state
    use ocn_vmix
+   use ocn_time_average
 
-
    implicit none
    private
    save
@@ -1251,6 +1251,8 @@
                           block % state % time_levs(2) % state % uReconstructMeridional % array    &amp;
                          )
 
+         call ocn_time_average_accumulate(block % state % time_levs(2) % state, block % state % time_levs(1) % state)
+
          block =&gt; block % next
       end do
       call mpas_timer_stop(&quot;split_explicit_timestep&quot;)

</font>
</pre>