<p><b>duda</b> 2010-04-23 14:28:13 -0600 (Fri, 23 Apr 2010)</p><p>BRANCH COMMIT <br>
<br>
Add implementation of CAM_RESTART_TO_MPAS(...) assuming the<br>
addition of vertical grid information to the argument list.<br>
<br>
<br>
M    module_mpas_cam_interface.F<br>
</p><hr noshade><pre><font color="gray">Modified: branches/mpas_cam_coupling/src/driver_cam_interface/module_mpas_cam_interface.F
===================================================================
--- branches/mpas_cam_coupling/src/driver_cam_interface/module_mpas_cam_interface.F        2010-04-22 23:14:24 UTC (rev 206)
+++ branches/mpas_cam_coupling/src/driver_cam_interface/module_mpas_cam_interface.F        2010-04-23 20:28:13 UTC (rev 207)
@@ -13,17 +13,7 @@
    real (kind=RKIND) :: dt
    integer :: itimestep
 
-   ! TEMPORARY: Once we are initializing MPAS fields from CAM-supplied
-   ! arrays, we won't need to call mpas_setup_test_case to get initial
-   ! fields
-   interface mpas_setup_test_case
-      subroutine mpas_setup_test_case(domain)
-         use grid_types
-         type (domain_type), intent(inout) :: domain
-      end subroutine mpas_setup_test_case
-   end interface mpas_setup_test_case
 
-
    contains
 
 
@@ -310,7 +300,8 @@
    end subroutine cam_inidat_to_mpas
    
    
-   subroutine cam_restart_to_mpas(Numcols, MaxEdges, Plev, Pcnst, Psd, Phis, T, U, W, Tracer)
+   subroutine cam_restart_to_mpas(Numcols, MaxEdges, Plev, Pcnst, Pref, Ptop, &amp;
+                                  Ae, Be, Ac, Bc, Psd, Phis, T, U, W, Tracer)
    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
    ! SUBROUTINE MPAS_RESTART_TO_MPAS
    !
@@ -318,6 +309,10 @@
    !        MaxEdges - maximum number of edges per cell
    !        Plev - number of vertical levels
    !        Pcnst - number of tracers
+   !        Pref - reference surface pressure (Pa)
+   !        Ptop - pressure at top of model (Pa)
+   !        Ae, Be - A and B arrays at vertical edges
+   !        Ac, Bc - A and B arrays at vertical centers
    !        Psd - dry surface pressure (Pa)
    !        Phis - surface geopotential (m^2/s^2)
    !        T - temperature (K)
@@ -332,6 +327,10 @@
       integer, intent(in) :: MaxEdges
       integer, intent(in) :: Plev
       integer, intent(in) :: Pcnst
+      real (kind=RKIND), intent(in) :: Pref
+      real (kind=RKIND), intent(in) :: Ptop
+      real (kind=RKIND), dimension(Plev+1), intent(in) :: Ae, Be
+      real (kind=RKIND), dimension(Plev), intent(in) :: Ac, Bc
       real (kind=RKIND), dimension(Numcols), intent(in) :: Psd
       real (kind=RKIND), dimension(Numcols), intent(in) :: Phis
       real (kind=RKIND), dimension(Numcols,Plev), intent(in) :: T
@@ -339,6 +338,7 @@
       real (kind=RKIND), dimension(Numcols,Plev+1), intent(in) :: W
       real (kind=RKIND), dimension(Numcols,Plev,Pcnst), intent(in) :: Tracer
 
+      integer :: iCell, k, i
       type (block_type), pointer :: block
    
       write(0,*) 'Called CAM_RESTART_TO_MPAS'
@@ -350,48 +350,42 @@
       !
       if (Numcols /= block % mesh % nCellsSolve) then
          write(0,*) 'Error: mismatch between Numcols and nCellsSolve: ', Numcols, block % mesh % nCellsSolve
+         return
       end if
       if (Plev /= block % mesh % nVertLevels) then
          write(0,*) 'Error: mismatch between Plev and nVertLevels: ', Plev, block % mesh % nVertLevels
+         return
       end if
       if (Pcnst /= num_scalars) then
          write(0,*) 'Error: mismatch between Pcnst and num_scalars: ', Pcnst, num_scalars
+         return
       end if
 
 
       !
       ! Fill in MPAS fields in block % time_levs(1) % state arrays from arguments 
+      ! Initialize most variables using cam_inidat_to_mpas
       !
+      call cam_inidat_to_mpas(Numcols, MaxEdges, Plev, Pcnst, Pref, Ptop, &amp;
+                              Ae, Be, Ac, Bc, Psd, Phis, T, U, Tracer)
 
-! TEMPORARY: until we are actually given fields, just call one of the test case
-! setup routines to initialize the MPAS prognostic fields
-      call mpas_setup_test_case(domain)
+      do iCell=1,Numcols
+         do k=1,Plev+1
+            block % time_levs(1) % state % ww % array(k,iCell) = W(iCell,k)
+         end do
+      end do
 
 
       !
       ! Do a ghost-cell update
       !
-      call dmpar_exch_halo_field1dReal(domain % dminfo, block % time_levs(1) % state % surface_pressure % array(:), &amp;
-                                       block % mesh % nCells, &amp;
-                                       block % parinfo % cellsToSend, block % parinfo % cellsToRecv)
-      call dmpar_exch_halo_field2dReal(domain % dminfo, block % time_levs(1) % state % theta % array(:,:), &amp;
-                                       block % mesh % nVertLevels, block % mesh % nCells, &amp;
-                                       block % parinfo % cellsToSend, block % parinfo % cellsToRecv)
-      call dmpar_exch_halo_field2dReal(domain % dminfo, block % time_levs(1) % state % u % array(:,:), &amp;
-                                       block % mesh % nVertLevels, block % mesh % nEdges, &amp;
-                                       block % parinfo % edgesToSend, block % parinfo % edgesToRecv)
-      call dmpar_exch_halo_field2dReal(domain % dminfo, block % time_levs(1) % state % pressure % array(:,:), &amp;
-                                       block % mesh % nVertLevels+1, block % mesh % nCells, &amp;
-                                       block % parinfo % cellsToSend, block % parinfo % cellsToRecv)
-      call dmpar_exch_halo_field2dReal(domain % dminfo, block % time_levs(1) % state % geopotential % array(:,:), &amp;
-                                       block % mesh % nVertLevels+1, block % mesh % nCells, &amp;
-                                       block % parinfo % cellsToSend, block % parinfo % cellsToRecv)
       call dmpar_exch_halo_field2dReal(domain % dminfo, block % time_levs(1) % state % ww % array(:,:), &amp;
                                        block % mesh % nVertLevels+1, block % mesh % nCells, &amp;
                                        block % parinfo % cellsToSend, block % parinfo % cellsToRecv)
-      call dmpar_exch_halo_field3dReal(domain % dminfo, block % time_levs(1) % state % scalars % array(:,:,:), &amp;
-                                       num_scalars, block % mesh % nVertLevels, block % mesh % nCells, &amp;
-                                       block % parinfo % cellsToSend, block % parinfo % cellsToRecv)
+
+      do i=2,nTimeLevs
+         block % time_levs(i) % state % ww % array(:,:) = block % time_levs(1) % state % ww % array(:,:)
+      end do
    
    end subroutine cam_restart_to_mpas
    

</font>
</pre>