<p><b>duda</b> 2011-10-11 15:03:15 -0600 (Tue, 11 Oct 2011)</p><p>Pass uReconstruct{X,Y,Z,Zonal,Meridional} and u (normal velocity) as separate<br>
array arguments to the reconstruct() routine. This allows individual cores to<br>
store these fields in different derived types, making the reconstruction code<br>
effectively core-independent.<br>
<br>
<br>
M    src/core_hyd_atmos/module_mpas_core.F<br>
M    src/core_hyd_atmos/Registry<br>
M    src/core_hyd_atmos/module_time_integration.F<br>
M    src/core_sw/module_mpas_core.F<br>
M    src/core_sw/module_time_integration.F<br>
M    src/core_ocean/mpas_ocn_time_integration_rk4.F<br>
M    src/core_ocean/mpas_ocn_time_integration_split.F<br>
M    src/core_ocean/mpas_ocn_mpas_core.F<br>
M    src/operators/module_vector_reconstruction.F<br>
</p><hr noshade><pre><font color="gray">Modified: trunk/mpas/src/core_hyd_atmos/Registry
===================================================================
--- trunk/mpas/src/core_hyd_atmos/Registry        2011-10-11 21:01:57 UTC (rev 1063)
+++ trunk/mpas/src/core_hyd_atmos/Registry        2011-10-11 21:03:15 UTC (rev 1064)
@@ -137,11 +137,11 @@
 var persistent real    ke ( nVertLevels nCells Time ) 2 o ke state - -
 var persistent real    pv_vertex ( nVertLevels nVertices Time ) 2 o pv_vertex state - -
 var persistent real    pv_cell ( nVertLevels nCells Time ) 2 o pv_cell state - -
-var persistent real    uReconstructX ( nVertLevels nCells Time ) 2 o uReconstructX state - -
-var persistent real    uReconstructY ( nVertLevels nCells Time ) 2 o uReconstructY state - -
-var persistent real    uReconstructZ ( nVertLevels nCells Time ) 2 o uReconstructZ state - -
-var persistent real    uReconstructZonal ( nVertLevels nCells Time ) 2 o uReconstructZonal state - -
-var persistent real    uReconstructMeridional ( nVertLevels nCells Time ) 2 o uReconstructMeridional state - -
+var persistent real    uReconstructX ( nVertLevels nCells Time ) 1 o uReconstructX diag - -
+var persistent real    uReconstructY ( nVertLevels nCells Time ) 1 o uReconstructY diag - -
+var persistent real    uReconstructZ ( nVertLevels nCells Time ) 1 o uReconstructZ diag - -
+var persistent real    uReconstructZonal ( nVertLevels nCells Time ) 1 o uReconstructZonal diag - -
+var persistent real    uReconstructMeridional ( nVertLevels nCells Time ) 1 o uReconstructMeridional diag - -
 
 # Tendency variables
 var persistent real    tend_h ( nVertLevels nCells Time ) 1 - h tend - -

Modified: trunk/mpas/src/core_hyd_atmos/module_mpas_core.F
===================================================================
--- trunk/mpas/src/core_hyd_atmos/module_mpas_core.F        2011-10-11 21:01:57 UTC (rev 1063)
+++ trunk/mpas/src/core_hyd_atmos/module_mpas_core.F        2011-10-11 21:03:15 UTC (rev 1064)
@@ -124,7 +124,14 @@
       call initialize_advection_rk(mesh)
       call rbfInterp_initialize(mesh)
       call init_reconstruct(mesh)
-      call reconstruct(block % state % time_levs(1) % state, mesh)
+      call reconstruct(mesh, block % state % time_levs(1) % state % u % array, &amp;
+                       block % diag % uReconstructX % array,                   &amp;
+                       block % diag % uReconstructY % array,                   &amp;
+                       block % diag % uReconstructZ % array,                   &amp;
+                       block % diag % uReconstructZonal % array,               &amp;
+                       block % diag % uReconstructMeridional % array           &amp;
+                      )
+
   
    end subroutine mpas_init_block
    

Modified: trunk/mpas/src/core_hyd_atmos/module_time_integration.F
===================================================================
--- trunk/mpas/src/core_hyd_atmos/module_time_integration.F        2011-10-11 21:01:57 UTC (rev 1063)
+++ trunk/mpas/src/core_hyd_atmos/module_time_integration.F        2011-10-11 21:03:15 UTC (rev 1064)
@@ -321,7 +321,14 @@
       !
       block =&gt; domain % blocklist
       do while (associated(block))
-         call reconstruct(block % state % time_levs(2) % state, block % mesh)
+         call reconstruct(block % mesh, block % state % time_levs(2) % state % u % array, &amp;
+                          block % diag % uReconstructX % array,                           &amp;
+                          block % diag % uReconstructY % array,                           &amp;
+                          block % diag % uReconstructZ % array,                           &amp;
+                          block % diag % uReconstructZonal % array,                       &amp;
+                          block % diag % uReconstructMeridional % array                   &amp;
+                         )
+
          call compute_w(block % state % time_levs(2) % state, block % state % time_levs(1) % state, block % mesh, dt)
          block =&gt; block % next
       end do

Modified: trunk/mpas/src/core_ocean/mpas_ocn_mpas_core.F
===================================================================
--- trunk/mpas/src/core_ocean/mpas_ocn_mpas_core.F        2011-10-11 21:01:57 UTC (rev 1063)
+++ trunk/mpas/src/core_ocean/mpas_ocn_mpas_core.F        2011-10-11 21:03:15 UTC (rev 1064)
@@ -216,7 +216,13 @@
  
       call rbfInterp_initialize(mesh)
       call init_reconstruct(mesh)
-      call reconstruct(block % state % time_levs(1) % state, mesh)
+      call reconstruct(mesh, block % state % time_levs(1) % state % u % array,                  &amp;
+                       block % state % time_levs(1) % state % uReconstructX % array,            &amp;
+                       block % state % time_levs(1) % state % uReconstructY % array,            &amp;
+                       block % state % time_levs(1) % state % uReconstructZ % array,            &amp;
+                       block % state % time_levs(1) % state % uReconstructZonal % array,        &amp;
+                       block % state % time_levs(1) % state % uReconstructMeridional % array    &amp;
+                      )
 
       ! initialize velocities and tracers on land to be -1e34
       ! The reconstructed velocity on land will have values not exactly

Modified: trunk/mpas/src/core_ocean/mpas_ocn_time_integration_rk4.F
===================================================================
--- trunk/mpas/src/core_ocean/mpas_ocn_time_integration_rk4.F        2011-10-11 21:01:57 UTC (rev 1063)
+++ trunk/mpas/src/core_ocean/mpas_ocn_time_integration_rk4.F        2011-10-11 21:03:15 UTC (rev 1064)
@@ -347,7 +347,13 @@
 
          call ocn_diagnostic_solve(dt, block % state % time_levs(2) % state, block % mesh)
 
-         call reconstruct(block % state % time_levs(2) % state, block % mesh)
+         call reconstruct(block % mesh, block % state % time_levs(2) % state % u % array,          &amp;
+                          block % state % time_levs(2) % state % uReconstructX % array,            &amp;
+                          block % state % time_levs(2) % state % uReconstructY % array,            &amp;
+                          block % state % time_levs(2) % state % uReconstructZ % array,            &amp;
+                          block % state % time_levs(2) % state % uReconstructZonal % array,        &amp;
+                          block % state % time_levs(2) % state % uReconstructMeridional % array    &amp;
+                         )
 
          block =&gt; block % next
       end do

Modified: trunk/mpas/src/core_ocean/mpas_ocn_time_integration_split.F
===================================================================
--- trunk/mpas/src/core_ocean/mpas_ocn_time_integration_split.F        2011-10-11 21:01:57 UTC (rev 1063)
+++ trunk/mpas/src/core_ocean/mpas_ocn_time_integration_split.F        2011-10-11 21:03:15 UTC (rev 1064)
@@ -1137,7 +1137,13 @@
 
          call ocn_diagnostic_solve(dt, block % state % time_levs(2) % state, block % mesh)
 
-         call reconstruct(block % state % time_levs(2) % state, block % mesh)
+         call reconstruct(block % mesh, block % state % time_levs(2) % state % u % array,          &amp;
+                          block % state % time_levs(2) % state % uReconstructX % array,            &amp;
+                          block % state % time_levs(2) % state % uReconstructY % array,            &amp;
+                          block % state % time_levs(2) % state % uReconstructZ % array,            &amp;
+                          block % state % time_levs(2) % state % uReconstructZonal % array,        &amp;
+                          block % state % time_levs(2) % state % uReconstructMeridional % array    &amp;
+                         )
 
          block =&gt; block % next
       end do

Modified: trunk/mpas/src/core_sw/module_mpas_core.F
===================================================================
--- trunk/mpas/src/core_sw/module_mpas_core.F        2011-10-11 21:01:57 UTC (rev 1063)
+++ trunk/mpas/src/core_sw/module_mpas_core.F        2011-10-11 21:03:15 UTC (rev 1064)
@@ -130,7 +130,14 @@
 
       call rbfInterp_initialize(mesh)
       call init_reconstruct(mesh)
-      call reconstruct(block % state % time_levs(1) % state, mesh)
+      call reconstruct(mesh, block % state % time_levs(1) % state % u % array,                  &amp;
+                       block % state % time_levs(1) % state % uReconstructX % array,            &amp;
+                       block % state % time_levs(1) % state % uReconstructY % array,            &amp;
+                       block % state % time_levs(1) % state % uReconstructZ % array,            &amp;
+                       block % state % time_levs(1) % state % uReconstructZonal % array,        &amp;
+                       block % state % time_levs(1) % state % uReconstructMeridional % array    &amp;
+                      )
+
    
    end subroutine mpas_init_block
    

Modified: trunk/mpas/src/core_sw/module_time_integration.F
===================================================================
--- trunk/mpas/src/core_sw/module_time_integration.F        2011-10-11 21:01:57 UTC (rev 1063)
+++ trunk/mpas/src/core_sw/module_time_integration.F        2011-10-11 21:03:15 UTC (rev 1064)
@@ -229,7 +229,13 @@
 
          call compute_solve_diagnostics(dt, block % state % time_levs(2) % state, block % mesh)
 
-         call reconstruct(block % state % time_levs(2) % state, block % mesh)
+         call reconstruct(block % mesh, block % state % time_levs(2) % state % u % array,          &amp;
+                          block % state % time_levs(2) % state % uReconstructX % array,            &amp;
+                          block % state % time_levs(2) % state % uReconstructY % array,            &amp;
+                          block % state % time_levs(2) % state % uReconstructZ % array,            &amp;
+                          block % state % time_levs(2) % state % uReconstructZonal % array,        &amp;
+                          block % state % time_levs(2) % state % uReconstructMeridional % array    &amp;
+                         )
 
          block =&gt; block % next
       end do

Modified: trunk/mpas/src/operators/module_vector_reconstruction.F
===================================================================
--- trunk/mpas/src/operators/module_vector_reconstruction.F        2011-10-11 21:01:57 UTC (rev 1063)
+++ trunk/mpas/src/operators/module_vector_reconstruction.F        2011-10-11 21:03:15 UTC (rev 1064)
@@ -112,7 +112,7 @@
 
   end subroutine init_reconstruct
 
-  subroutine reconstruct(state, grid)
+  subroutine reconstruct(grid, u, uReconstructX, uReconstructY, uReconstructZ, uReconstructZonal, uReconstructMeridional)
   !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
   ! Purpose: reconstruct vector field at cell centers based on radial basis functions
   !
@@ -123,18 +123,17 @@
 
     implicit none
 
-    type (state_type), intent(inout) :: state 
     type (mesh_type), intent(in) :: grid
+    real (kind=RKIND), dimension(:,:), intent(in) :: u
+    real (kind=RKIND), dimension(:,:), intent(out) :: uReconstructX, uReconstructY, uReconstructZ
+    real (kind=RKIND), dimension(:,:), intent(out) :: uReconstructZonal, uReconstructMeridional
 
     !   temporary arrays needed in the compute procedure
     integer :: nCellsSolve
     integer, dimension(:,:), pointer :: edgesOnCell
     integer, dimension(:), pointer :: nEdgesOnCell
     integer :: iCell,iEdge, i
-    real (kind=RKIND), dimension(:,:), pointer :: u
     real(kind=RKIND), dimension(:), pointer :: latCell, lonCell
-    real (kind=RKIND), dimension(:,:), pointer :: uReconstructX, uReconstructY, uReconstructZ
-    real (kind=RKIND), dimension(:,:), pointer :: uReconstructZonal, uReconstructMeridional
 
     real (kind=RKIND), dimension(:,:,:), pointer :: coeffs_reconstruct
 
@@ -150,15 +149,9 @@
     edgesOnCell =&gt; grid % edgesOnCell % array
     nEdgesOnCell=&gt; grid % nEdgesOnCell % array
     nCellsSolve = grid % nCellsSolve
-    u =&gt; state % u % array
-    uReconstructX =&gt; state % uReconstructX % array
-    uReconstructY =&gt; state % uReconstructY % array
-    uReconstructZ =&gt; state % uReconstructZ % array
 
     latCell       =&gt; grid % latCell % array
     lonCell       =&gt; grid % lonCell % array
-    uReconstructZonal =&gt; state % uReconstructZonal % array
-    uReconstructMeridional =&gt; state % uReconstructMeridional % array
     on_a_sphere = grid % on_a_sphere
 
     ! init the intent(out)

</font>
</pre>