<p><b>dwj07@fsu.edu</b> 2011-12-21 16:07:11 -0700 (Wed, 21 Dec 2011)</p><p><br>
        -- BRANCH COMMIT --<br>
<br>
        Adding a time averaging module.<br>
<br>
        Right now it is always zero for some reason. That needs to be fixed.<br>
</p><hr noshade><pre><font color="gray">Modified: branches/ocean_projects/time_averaging/src/core_ocean/Makefile
===================================================================
--- branches/ocean_projects/time_averaging/src/core_ocean/Makefile        2011-12-21 22:40:19 UTC (rev 1272)
+++ branches/ocean_projects/time_averaging/src/core_ocean/Makefile        2011-12-21 23:07:11 UTC (rev 1273)
@@ -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: branches/ocean_projects/time_averaging/src/core_ocean/Registry
===================================================================
--- branches/ocean_projects/time_averaging/src/core_ocean/Registry        2011-12-21 22:40:19 UTC (rev 1272)
+++ branches/ocean_projects/time_averaging/src/core_ocean/Registry        2011-12-21 23:07:11 UTC (rev 1273)
@@ -258,3 +258,11 @@
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 - -
Modified: branches/ocean_projects/time_averaging/src/core_ocean/mpas_ocn_mpas_core.F
===================================================================
--- branches/ocean_projects/time_averaging/src/core_ocean/mpas_ocn_mpas_core.F        2011-12-21 22:40:19 UTC (rev 1272)
+++ branches/ocean_projects/time_averaging/src/core_ocean/mpas_ocn_mpas_core.F        2011-12-21 23:07:11 UTC (rev 1273)
@@ -23,6 +23,8 @@
use ocn_vmix
+ use ocn_time_average
+
type (io_output_object) :: restart_obj
integer :: current_outfile_frames
@@ -207,6 +209,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)
@@ -299,7 +302,13 @@
write(0,*) 'Initial time ', timeStamp
call write_output_frame(output_obj, output_frame, domain)
+ block_ptr => domain % blocklist
+ do while(associated(block_ptr))
+ call ocn_time_average_init(block_ptr % state % time_levs(1) % state)
+ block_ptr => 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
@@ -326,7 +335,18 @@
call mpas_output_state_finalize(output_obj, domain % dminfo)
call mpas_output_state_init(output_obj, domain, "OUTPUT", trim(timeStamp))
end if
+
+ block_ptr => domain % blocklist
+ do while (associated(block_ptr))
+ call ocn_time_average_normalize(block_ptr % state % time_levs(1) % state)
+ block_ptr => block_ptr % next
+ end do
call write_output_frame(output_obj, output_frame, domain)
+ block_ptr => domain % blocklist
+ do while (associated(block_ptr))
+ call ocn_time_average_init(block_ptr % state % time_levs(1) % state)
+ block_ptr => block_ptr % next
+ end do
end if
if (mpas_is_alarm_ringing(clock, restartAlarmID, ierr=ierr)) then
Modified: branches/ocean_projects/time_averaging/src/core_ocean/mpas_ocn_tendency.F
===================================================================
--- branches/ocean_projects/time_averaging/src/core_ocean/mpas_ocn_tendency.F        2011-12-21 22:40:19 UTC (rev 1272)
+++ branches/ocean_projects/time_averaging/src/core_ocean/mpas_ocn_tendency.F        2011-12-21 23:07:11 UTC (rev 1273)
@@ -38,6 +38,8 @@
use ocn_equation_of_state
use ocn_vmix
+ use ocn_time_average
+
implicit none
private
save
@@ -1107,6 +1109,8 @@
call ocn_wtop(s,grid)
+ call ocn_time_average_accumulate(s)
+
call mpas_timer_stop("ocn_diagnostic_solve")
end subroutine ocn_diagnostic_solve!}}}
Added: branches/ocean_projects/time_averaging/src/core_ocean/mpas_ocn_time_average.F
===================================================================
--- branches/ocean_projects/time_averaging/src/core_ocean/mpas_ocn_time_average.F         (rev 0)
+++ branches/ocean_projects/time_averaging/src/core_ocean/mpas_ocn_time_average.F        2011-12-21 23:07:11 UTC (rev 1273)
@@ -0,0 +1,98 @@
+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
+
+ nAccumulate => state % nAccumulate % scalar
+
+ acc_ssh => state % acc_ssh % array
+ acc_sshVar => state % acc_sshVar % array
+ acc_uReconstructZonal => state % acc_uReconstructZonal % array
+ acc_uReconstructMeridional => state % acc_uReconstructMeridional % array
+ acc_uReconstructZonalVar => state % acc_uReconstructZonalVar % array
+ acc_uReconstructMeridionalVar => state % acc_uReconstructMeridionalVar % 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
+
+ end subroutine ocn_time_average_init!}}}
+
+ subroutine ocn_time_average_accumulate(state)!{{{
+ type (state_type), intent(inout) :: state
+
+ real, pointer :: nAccumulate
+
+ real (kind=RKIND), dimension(:), pointer :: ssh
+ real (kind=RKIND), dimension(:,:), pointer :: uReconstructZonal, uReconstructMeridional
+
+ real (kind=RKIND), dimension(:,:), pointer :: acc_uReconstructZonal, acc_uReconstructMeridional, acc_uReconstructZonalVar, acc_uReconstructMeridionalVar
+ real (kind=RKIND), dimension(:), pointer :: acc_ssh, acc_sshVar
+
+ nAccumulate => state % nAccumulate % scalar
+
+ ssh => state % ssh % array
+ uReconstructZonal => state % uReconstructZonal % array
+ uReconstructMeridional => state % uReconstructMeridional % array
+
+ acc_ssh => state % acc_ssh % array
+ acc_sshVar => state % acc_sshVar % array
+ acc_uReconstructZonal => state % acc_uReconstructZonal % array
+ acc_uReconstructMeridional => state % acc_uReconstructMeridional % array
+ acc_uReconstructZonalVar => state % acc_uReconstructZonalVar % array
+ acc_uReconstructMeridionalVar => state % acc_uReconstructMeridionalVar % array
+
+ acc_ssh = acc_ssh + ssh
+ acc_sshVar = acc_sshVar + ssh**2
+ acc_uReconstructZonal = acc_uReconstructZonal + uReconstructZonal
+ acc_uReconstructMeridional = acc_uReconstructMeridional + uReconstructMeridional
+ acc_uReconstructZonalVar = acc_uReconstructZonalVar + uReconstructZonal**2
+ acc_uReconstructMeridionalVar = acc_uReconstructMeridionalVar + uReconstructMeridional**2
+
+ nAccumulate = 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
+
+ nAccumulate => state % nAccumulate % scalar
+
+ acc_ssh => state % acc_ssh % array
+ acc_sshVar => state % acc_sshVar % array
+ acc_uReconstructZonal => state % acc_uReconstructZonal % array
+ acc_uReconstructMeridional => state % acc_uReconstructMeridional % array
+ acc_uReconstructZonalVar => state % acc_uReconstructZonalVar % array
+ acc_uReconstructMeridionalVar => state % acc_uReconstructMeridionalVar % 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
+ end subroutine ocn_time_average_normalize!}}}
+
+end module ocn_time_average
</font>
</pre>