<p><b>ringler@lanl.gov</b> 2010-03-30 12:06:43 -0600 (Tue, 30 Mar 2010)</p><p><br>
changes to support new data fields needed for lateral boundaries<br>
<br>
Makefile: turn off expand_levels to run in single-layer shallow-water mode<br>
mpas_interface: turn off vector reconstruction (will have to made consistent with new formulation)<br>
Registry: added boundaryEdge and boundaryVertex data fields<br>
module_time_integration: change &quot;uBC&quot; to &quot;boundaryEdge&quot;, turn off vector reconstruction<br>
</p><hr noshade><pre><font color="gray">Modified: branches/lateral_boundary_conditions/Makefile
===================================================================
--- branches/lateral_boundary_conditions/Makefile        2010-03-29 17:50:58 UTC (rev 169)
+++ branches/lateral_boundary_conditions/Makefile        2010-03-30 18:06:43 UTC (rev 170)
@@ -1,7 +1,7 @@
 #MODEL_FORMULATION = -DNCAR_FORMULATION
 MODEL_FORMULATION = -DLANL_FORMULATION
 
-EXPAND_LEVELS = -DEXPAND_LEVELS=26
+#EXPAND_LEVELS = -DEXPAND_LEVELS=26
 #FILE_OFFSET = -DOFFSET64BIT
 
 #########################

Modified: branches/lateral_boundary_conditions/src/core_sw/Registry
===================================================================
--- branches/lateral_boundary_conditions/src/core_sw/Registry        2010-03-29 17:50:58 UTC (rev 169)
+++ branches/lateral_boundary_conditions/src/core_sw/Registry        2010-03-30 18:06:43 UTC (rev 170)
@@ -83,7 +83,8 @@
 var real    rbf_value ( maxEdges nCells ) - rbf_value - -
 
 # Boundary conditions: read from input, saved in restart and written to output
-var integer uBC ( nVertLevels nEdges ) iro uBC - -
+var integer boundaryEdge ( nVertLevels nEdges ) iro boundaryEdge - -
+var integer boundaryVertex ( nVertLevels nVertices ) iro boundaryVertex - -
 
 # Prognostic variables: read from input, saved in restart, and written to output
 var real    u ( nVertLevels nEdges Time ) iro u - -

Modified: branches/lateral_boundary_conditions/src/core_sw/module_time_integration.F
===================================================================
--- branches/lateral_boundary_conditions/src/core_sw/module_time_integration.F        2010-03-29 17:50:58 UTC (rev 169)
+++ branches/lateral_boundary_conditions/src/core_sw/module_time_integration.F        2010-03-30 18:06:43 UTC (rev 170)
@@ -1,6 +1,6 @@
 module time_integration
 
-   use vector_reconstruction
+!  use vector_reconstruction
    use grid_types
    use configure
    use constants
@@ -125,7 +125,7 @@
         do while (associated(block))
            call compute_tend(block % intermediate_step(TEND), block % intermediate_step(PROVIS), block % mesh)
            call compute_scalar_tend(block % intermediate_step(TEND), block % intermediate_step(PROVIS), block % mesh)
-           call enforce_uBC(block % intermediate_step(TEND), block % mesh)
+           call enforce_boundaryEdge(block % intermediate_step(TEND), block % mesh)
            block =&gt; block % next
         end do
 
@@ -214,7 +214,7 @@
 
          call compute_solve_diagnostics(dt, block % time_levs(2) % state, block % mesh)
 
-         call reconstruct(block % time_levs(2) % state, block % mesh)
+!        call reconstruct(block % time_levs(2) % state, block % mesh)
 
          block =&gt; block % next
       end do
@@ -450,7 +450,7 @@
       real (kind=RKIND), dimension(:), pointer :: h_s, fVertex, fEdge, dvEdge, dcEdge, areaCell, areaTriangle
       real (kind=RKIND), dimension(:,:), pointer :: vh, weightsOnEdge, kiteAreasOnVertex, h_edge, h, u, v, tend_h, tend_u, &amp;
                                                     circulation, vorticity, ke, pv_edge, pv_vertex, pv_cell, gradPVn, gradPVt, divergence
-      integer, dimension(:,:), pointer :: cellsOnEdge, cellsOnVertex, verticesOnEdge, edgesOnCell, edgesOnEdge, edgesOnVertex, uBC
+      integer, dimension(:,:), pointer :: cellsOnEdge, cellsOnVertex, verticesOnEdge, edgesOnCell, edgesOnEdge, edgesOnVertex, boundaryEdge
       integer, dimension(:), pointer :: nEdgesOnCell, nEdgesOnEdge
       real (kind=RKIND) :: r, h1, h2
 
@@ -495,7 +495,7 @@
       nVertices   = grid % nVertices
       nVertLevels = grid % nVertLevels
 
-      uBC =&gt; grid % uBC % array
+      boundaryEdge =&gt; grid % boundaryEdge % array
 
       !
       ! Compute height on cell edges at velocity locations
@@ -712,12 +712,12 @@
       !
       ! set pv_edge = fEdge / h_edge at boundary points
       !
-      if (maxval(uBC).gt.0) then
+      if (maxval(boundaryEdge).gt.0) then
       do iEdge = 1,nEdges
          cell1 = cellsOnEdge(1,iEdge)
          cell2 = cellsOnEdge(2,iEdge)
          do k = 1,nVertLevels
-           if(uBC(k,iEdge).eq.1) then
+           if(boundaryEdge(k,iEdge).eq.1) then
              if(cell1.gt.0) h1 = h(k,cell1)
              if(cell2.gt.0) h2 = h(k,cell2)
              pv_edge(k,iEdge) = fEdge(iEdge) / ( max(h1,h2) )
@@ -731,13 +731,13 @@
    end subroutine compute_solve_diagnostics
 
 
-   subroutine enforce_uBC(tend, grid)
+   subroutine enforce_boundaryEdge(tend, grid)
    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
    ! Enforce any boundary conditions on the normal velocity at each edge
    !
    ! Input: grid - grid metadata
    !
-   ! Output: tend_u set to zero at uBC == 1 locations
+   ! Output: tend_u set to zero at boundaryEdge == 1 locations
    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
 
 
@@ -746,7 +746,7 @@
       type (grid_state), intent(inout) :: tend
       type (grid_meta), intent(in) :: grid
 
-      integer, dimension(:,:), pointer :: uBC
+      integer, dimension(:,:), pointer :: boundaryEdge
       real (kind=RKIND), dimension(:,:), pointer :: tend_u
       integer :: nCells, nEdges, nVertices, nVertLevels
       integer :: iEdge, k
@@ -756,22 +756,22 @@
       nVertices   = grid % nVertices
       nVertLevels = grid % nVertLevels
 
-      uBC         =&gt; grid % uBC % array
-      tend_u      =&gt; tend % u % array
+      boundaryEdge         =&gt; grid % boundaryEdge % array
+      tend_u               =&gt; tend % u % array
 
-      if(maxval(uBC).le.0) return
+      if(maxval(boundaryEdge).le.0) return
 
       do iEdge = 1,nEdges
         do k = 1,nVertLevels
 
-          if(uBC(k,iEdge).eq.1) then
+          if(boundaryEdge(k,iEdge).eq.1) then
              tend_u(k,iEdge) = 0.0
           endif
 
         enddo
        enddo
 
-   end subroutine enforce_uBC
+   end subroutine enforce_boundaryEdge
 
 
 end module time_integration

Modified: branches/lateral_boundary_conditions/src/core_sw/mpas_interface.F
===================================================================
--- branches/lateral_boundary_conditions/src/core_sw/mpas_interface.F        2010-03-29 17:50:58 UTC (rev 169)
+++ branches/lateral_boundary_conditions/src/core_sw/mpas_interface.F        2010-03-30 18:06:43 UTC (rev 170)
@@ -2,7 +2,7 @@
 
    use grid_types
    use test_cases
-   use vector_reconstruction
+!  use vector_reconstruction
 
    implicit none
 
@@ -25,7 +25,7 @@
    real (kind=RKIND), intent(in) :: dt
 
    call compute_solve_diagnostics(dt, block % time_levs(1) % state, mesh)
-   call init_reconstruct(mesh)
+!  call init_reconstruct(mesh)
 
 end subroutine mpas_init
 

</font>
</pre>