<p><b>qchen3@fsu.edu</b> 2013-02-27 15:25:28 -0700 (Wed, 27 Feb 2013)</p><p>BRANCH COMMIT<br>
<br>
GM implementation stage 1 completed:<br>
The GM module is now de-linked from the diagnostics module, <br>
and both are now called side by side in the integrators <br>
(RK4 or split-explicit). This change should make the program<br>
better organized and the gm funcationality easier to use or <br>
to modify in the future.<br>
</p><hr noshade><pre><font color="gray">Modified: branches/ocean_projects/gm_implement/namelist.input.ocean
===================================================================
--- branches/ocean_projects/gm_implement/namelist.input.ocean        2013-02-27 21:40:19 UTC (rev 2513)
+++ branches/ocean_projects/gm_implement/namelist.input.ocean        2013-02-27 22:25:28 UTC (rev 2514)
@@ -2,16 +2,16 @@
         config_do_restart = .false.
         config_start_time = '0000-01-01_00:00:00'
         config_stop_time = 'none'
-        config_run_duration = '0_01:00:00'
+        config_run_duration = '1_00:00:00'
         config_calendar_type = '360day'
 /
 &amp;io
         config_input_name = 'grid.nc'
         config_output_name = 'output.nc'
         config_restart_name = 'restart.nc'
-        config_restart_interval = '0_06:00:00'
-        config_output_interval = '0_06:00:00'
-        config_stats_interval = '0_01:00:00'
+        config_restart_interval = '10_00:00:00'
+        config_output_interval = '1_00:00:00'
+        config_stats_interval = '1_00:00:00'
         config_write_stats_on_startup = .true.
         config_write_output_on_startup = .true.
         config_frames_per_outfile = 1000
@@ -19,7 +19,7 @@
         config_pio_stride = 1
 /
 &amp;time_integration
-        config_dt = 20.0
+        config_dt = 500.0
         config_time_integrator = 'split_explicit'
 /
 &amp;grid
@@ -61,7 +61,7 @@
         config_Leith_visc2_max = 2.5e3
 /
 &amp;standard_GM
-        config_h_kappa = 0.0
+        config_h_kappa = 10.0
         config_h_kappa_q = 0.0
 /
 &amp;Rayleigh_damping

Modified: branches/ocean_projects/gm_implement/src/core_ocean/mpas_ocn_diagnostics.F
===================================================================
--- branches/ocean_projects/gm_implement/src/core_ocean/mpas_ocn_diagnostics.F        2013-02-27 21:40:19 UTC (rev 2513)
+++ branches/ocean_projects/gm_implement/src/core_ocean/mpas_ocn_diagnostics.F        2013-02-27 22:25:28 UTC (rev 2514)
@@ -19,7 +19,6 @@
    use mpas_constants
    use mpas_timer
 
-   use ocn_gm
    use ocn_equation_of_state
 
    implicit none
@@ -102,15 +101,13 @@
         weightsOnEdge, kiteAreasOnVertex, h_edge, h, u, v, pressure,&amp;
         circulation, vorticity, ke, ke_edge, MontPot, wTop, zMid, &amp;
         Vor_edge, Vor_vertex, Vor_cell, gradVor_n, gradVor_t, divergence, &amp;
-        rho, rhoDisplaced, temperature, salinity, kev, kevc, uBolusGM, uTransport, &amp;
+        rho, rhoDisplaced, temperature, salinity, kev, kevc, &amp;
         vertVelocityTop, BruntVaisalaFreqTop
       real (kind=RKIND), dimension(:,:,:), pointer :: tracers, deriv_two
       character :: c1*6
 
       h           =&gt; s % h % array
       u           =&gt; s % u % array
-      uTransport  =&gt; s % uTransport % array
-      uBolusGM    =&gt; s % uBolusGM % array
       v           =&gt; s % v % array
       h_edge      =&gt; s % h_edge % array
       circulation =&gt; s % circulation % array
@@ -475,16 +472,6 @@
 
       end do
 
-      !
-      ! Apply the GM closure as a bolus velocity
-      !
-      if (config_h_kappa .GE. epsilon(0D0)) then
-         call ocn_gm_compute_uBolus(s,grid)
-      else
-         ! mrp efficiency note: if uBolusGM is guaranteed to be zero, this can be removed.
-         uBolusGM = 0.0
-      end if
-
    end subroutine ocn_diagnostic_solve!}}}
 
 !***********************************************************************

Modified: branches/ocean_projects/gm_implement/src/core_ocean/mpas_ocn_gm.F
===================================================================
--- branches/ocean_projects/gm_implement/src/core_ocean/mpas_ocn_gm.F        2013-02-27 21:40:19 UTC (rev 2513)
+++ branches/ocean_projects/gm_implement/src/core_ocean/mpas_ocn_gm.F        2013-02-27 22:25:28 UTC (rev 2514)
@@ -61,7 +61,7 @@
       else
 
          ! Nothing for now for all other grid types (zlevel, zstar, ztilde)
-         uBolusGM(:,:) = 0.0
+         uBolusGM(:,:) = 0.000000001     ! A tiny number just for testing for the moment
 
       end if
 

Modified: branches/ocean_projects/gm_implement/src/core_ocean/mpas_ocn_time_integration_rk4.F
===================================================================
--- branches/ocean_projects/gm_implement/src/core_ocean/mpas_ocn_time_integration_rk4.F        2013-02-27 21:40:19 UTC (rev 2513)
+++ branches/ocean_projects/gm_implement/src/core_ocean/mpas_ocn_time_integration_rk4.F        2013-02-27 22:25:28 UTC (rev 2514)
@@ -24,6 +24,7 @@
    use ocn_tendency
    use ocn_diagnostics
 
+   use ocn_gm
    use ocn_equation_of_state
    use ocn_vmix
    use ocn_time_average
@@ -215,7 +216,13 @@
 
               call ocn_diagnostic_solve(dt, block % provis, block % mesh)
 
-              ! Compute velocity transport, used in advection terms of h and tracer tendency
+              ! Compute uBolusGM and then tracer transport velocity 
+              if (config_h_kappa .GE. epsilon(0D0)) then
+                  call ocn_gm_compute_uBolus(block % provis, block % mesh)
+              else
+                  block % provis % uBolusGM % array(:,:) = 0.0
+              end if

               block % provis % uTransport % array(:,:) &amp;
                     = block % provis % u          % array(:,:) &amp;
                     + block % provis % uBolusGM   % array(:,:)
@@ -316,8 +323,14 @@
          end if
 
          call ocn_diagnostic_solve(dt, block % state % time_levs(2) % state, block % mesh)
+   
+         ! Compute uBolusGM and the tracer transport velocity
+         if (config_h_kappa .GE. epsilon(0D0)) then
+             call ocn_gm_compute_uBolus(block % state % time_levs(2) % state, block % mesh)
+         else
+             block % state % time_levs(2) % state % uBolusGM % array(:,:) = 0.0
+         end if
 
-         ! Compute velocity transport, used in advection terms of h and tracer tendency
             block % state % time_levs(2) % state % uTransport % array(:,:) &amp;
           = block % state % time_levs(2) % state % u % array(:,:) &amp;
           + block % state % time_levs(2) % state % uBolusGM % array(:,:)

Modified: branches/ocean_projects/gm_implement/src/core_ocean/mpas_ocn_time_integration_split.F
===================================================================
--- branches/ocean_projects/gm_implement/src/core_ocean/mpas_ocn_time_integration_split.F        2013-02-27 21:40:19 UTC (rev 2513)
+++ branches/ocean_projects/gm_implement/src/core_ocean/mpas_ocn_time_integration_split.F        2013-02-27 22:25:28 UTC (rev 2514)
@@ -25,6 +25,7 @@
 
    use ocn_tendency
    use ocn_diagnostics
+   use ocn_gm
 
    use ocn_equation_of_state
    use ocn_vmix
@@ -880,6 +881,15 @@
                ! I can par this down later.
                call ocn_diagnostic_solve(dt, block % state % time_levs(2) % state, block % mesh)
 
+              ! Compute uBolusGM
+              ! QC note: Not sure if uBolus needs to be computed here. If not, then the next few lines can be deleted.
+              !if (config_h_kappa .GE. epsilon(0D0)) then
+              !    call ocn_gm_compute_uBolus(block % state % time_levs(2) % state, block % mesh)
+              !else
+              !    block % state % time_levs(2) % state % uBolusGM % array(:,:) = 0.0
+              !end if
+
+
             !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
             !
             !  If large iteration complete, compute all variables at time n+1
@@ -953,6 +963,14 @@
           ! implicit vmix routine that follows. mrp 121023.
           call ocn_diagnostic_solve(dt, block % state % time_levs(2) % state, block % mesh)
 
+          ! Compute uBolusGM
+          ! QC note: probably uBolusGM needs not to be computed here.
+          !if (config_h_kappa .GE. epsilon(0D0)) then
+          !    call ocn_gm_compute_uBolus(block % state % time_levs(2) % state, block % mesh)
+          !else
+          !    block % state % time_levs(2) % state % uBolusGM % array(:,:) = 0.0
+          !end if
+
           call ocn_vmix_implicit(dt, block % mesh, block % diagnostics, block % state % time_levs(2) % state, err)
 
           block =&gt; block % next
@@ -987,6 +1005,13 @@
 
          call ocn_diagnostic_solve(dt, block % state % time_levs(2) % state, block % mesh)
 
+         ! Compute uBolusGM 
+         if (config_h_kappa .GE. epsilon(0D0)) then
+             call ocn_gm_compute_uBolus(block % state % time_levs(2) % state, block % mesh)
+         else
+             block % state % time_levs(2) % state % uBolusGM % array(:,:) = 0.0
+         end if
+
          ! Compute velocity transport, used in advection terms of h and tracer tendency
             block % state % time_levs(2) % state % uTransport % array(:,:) &amp;
           = block % state % time_levs(2) % state % u % array(:,:) &amp;

</font>
</pre>