<p><b>duda</b> 2012-04-09 18:56:10 -0600 (Mon, 09 Apr 2012)</p><p>BRANCH COMMIT<br>
<br>
- Add a new MPAS_seekStream() routine to the mpas_io_stream module that will return the <br>
  time frame in a file matching a specified date; available matching criteria are: exact<br>
  match, latest time before the specified time, earliest time after the specified time, and<br>
  the closest time to the specified time. This routine is useful particularly when restarting<br>
  the model.<br>
<br>
- Fix a bug in the reading of 0d character fields at the stream level.<br>
<br>
- Add code in the stream IO module to handle the situation where we are reading a super-array,<br>
  but only some of the constituent arrays are available in the file; this may happen, e.g.,<br>
  in the case of moist scalars, where we provided only qc, qv, and qr in the initial condition<br>
  file, but reserve space for qs and qi.<br>
<br>
<br>
M    src/framework/Makefile<br>
M    src/framework/mpas_io_streams.F<br>
</p><hr noshade><pre><font color="gray">Modified: branches/omp_blocks/io/src/framework/Makefile
===================================================================
--- branches/omp_blocks/io/src/framework/Makefile        2012-04-09 23:31:33 UTC (rev 1764)
+++ branches/omp_blocks/io/src/framework/Makefile        2012-04-10 00:56:10 UTC (rev 1765)
@@ -53,7 +53,7 @@
 
 mpas_io.o: mpas_dmpar_types.o
 
-mpas_io_streams.o: mpas_attlist.o mpas_io.o
+mpas_io_streams.o: mpas_attlist.o mpas_grid_types.o mpas_timekeeping.o mpas_io.o
 
 mpas_io_input.o: mpas_grid_types.o mpas_dmpar.o mpas_block_decomp.o mpas_sort.o mpas_configure.o mpas_timekeeping.o mpas_io_streams.o $(ZOLTANOBJ)
 

Modified: branches/omp_blocks/io/src/framework/mpas_io_streams.F
===================================================================
--- branches/omp_blocks/io/src/framework/mpas_io_streams.F        2012-04-09 23:31:33 UTC (rev 1764)
+++ branches/omp_blocks/io/src/framework/mpas_io_streams.F        2012-04-10 00:56:10 UTC (rev 1765)
@@ -2,12 +2,15 @@
 
    use mpas_attlist
    use mpas_grid_types
+   use mpas_timekeeping
    use mpas_io
 
    type field_list_type
       integer :: field_type
       logical :: isDecomposed
       integer :: totalDimSize     ! Total size of outer dimension across all blocks for decomposed fields
+      logical, dimension(:), pointer :: isAvailable =&gt; null() ! Used for reading super-arrays where one or more 
+                                                              !   constitutent arrays may not be present in the input file
       type (field0dInteger), pointer :: int0dField =&gt; null()
       type (field1dInteger), pointer :: int1dField =&gt; null()
       type (field2dInteger), pointer :: int2dField =&gt; null()
@@ -60,7 +63,12 @@
       module procedure MPAS_writeStreamAtt_text
    end interface MPAS_writeStreamAtt
 
+   integer, parameter :: MPAS_STREAM_EXACT_TIME     = 100, &amp;
+                         MPAS_STREAM_NEAREST        = 101, &amp;
+                         MPAS_STREAM_LATEST_BEFORE  = 102, &amp;
+                         MPAS_STREAM_EARLIEST_AFTER = 103
 
+
    ! Error codes
    integer, parameter :: MPAS_STREAM_NOERR              =  0, &amp;
                          MPAS_STREAM_NOT_INITIALIZED    = -1, &amp;
@@ -114,6 +122,115 @@
    end subroutine MPAS_createStream
 
 
+   integer function MPAS_seekStream(stream, seekTime, seekPosition, actualTime, ierr)
+
+      implicit none
+
+      type (MPAS_Stream_type), intent(inout) :: stream
+      character (len=*), intent(in) :: seekTime
+      integer, intent(in) :: seekPosition
+      character (len=*), intent(out) :: actualTime
+      integer, intent(out), optional :: ierr
+      
+      integer :: io_err
+      integer :: i
+      integer :: timeDim
+      character (len=64), dimension(:), pointer :: xtimes
+      character (len=32) :: strTemp
+      type (MPAS_Time_type) :: sliceTime, startTime
+      type (MPAS_TimeInterval_type) :: timeDiff, minTimeDiff
+
+      write(0,*) 'Called MPAS_seekStream'
+
+      !
+      ! Initialize output arguments
+      !
+      if (present(ierr)) ierr = MPAS_STREAM_NOERR
+      MPAS_seekStream = 0
+
+
+      !
+      ! Sanity checks
+      !
+      if (.not. stream % isInitialized) then
+         if (present(ierr)) ierr = MPAS_IO_ERR
+         return
+      end if
+
+      call MPAS_io_inq_dim(stream % fileHandle, 'Time', timeDim, io_err)
+      if (timeDim &lt;= 0 .or. io_err /= MPAS_IO_NOERR) then
+         if (present(ierr)) ierr = MPAS_IO_ERR
+         return
+      end if
+
+write(0,*) 'Found ', timeDim, ' times in file' 
+
+      allocate(xtimes(timeDim))
+
+      do i=1,timeDim
+         call MPAS_io_set_frame(stream % fileHandle, i, io_err)
+         call MPAS_io_get_var(stream % fileHandle, 'xtime', xtimes(i), io_err)
+      end do
+
+      if (len(seekTime) &gt; 32) then
+         write(strTemp, '(a)') seekTime(1:32)
+      else
+         write(strTemp, '(a)') trim(seekTime)
+      end if
+      call mpas_set_timeInterval(interval=minTimeDiff, DD=10000)
+      call mpas_set_time(curr_time=startTime, dateTimeString=strTemp)
+
+      do i=1,timeDim
+         write(strTemp, '(a)') trim(xtimes(i)(1:32))
+         call mpas_set_time(curr_time=sliceTime, dateTimeString=strTemp)
+         if (seekPosition == MPAS_STREAM_EXACT_TIME) then
+            if (sliceTime == startTime) then
+               minTimeDiff = timeDiff
+               MPAS_seekStream = i
+            end if
+         else if (seekPosition == MPAS_STREAM_NEAREST) then
+            timeDiff = abs(sliceTime - startTime)
+            if (timeDiff &lt; minTimeDiff) then
+               minTimeDiff = timeDiff
+               MPAS_seekStream = i
+            end if
+         else if (seekPosition == MPAS_STREAM_LATEST_BEFORE) then
+            if (sliceTime &lt;= startTime) then
+               timeDiff = abs(sliceTime - startTime)
+               if (timeDiff &lt; minTimeDiff) then
+                  minTimeDiff = timeDiff
+                  MPAS_seekStream = i
+               end if
+            end if
+         else if (seekPosition == MPAS_STREAM_EARLIEST_AFTER) then
+            if (sliceTime &gt;= startTime) then
+               timeDiff = abs(sliceTime - startTime)
+               if (timeDiff &lt; minTimeDiff) then
+                  minTimeDiff = timeDiff
+                  MPAS_seekStream = i
+               end if
+            end if
+         else
+            write(0,*) 'Error in MPAS_seekStream: unrecognized seekPosition'
+            deallocate(xtimes)
+            if (present(ierr)) ierr = MPAS_IO_ERR
+            return
+         end if
+      end do
+
+      if (MPAS_seekStream == 0) then
+         deallocate(xtimes)
+         if (present(ierr)) ierr = MPAS_IO_ERR
+         return
+      end if
+
+      write(actualTime, '(a)') trim(xtimes(MPAS_seekStream)(1:32))
+
+      deallocate(xtimes)
+
+   end function MPAS_seekStream
+
+
    subroutine MPAS_streamAddField_0dInteger(stream, field, ierr)
 
       implicit none
@@ -214,6 +331,8 @@
       integer, dimension(:), pointer :: indices
       type (field_list_type), pointer :: field_list_cursor
       type (field_list_type), pointer :: new_field_list_node
+      logical :: any_success
+      logical, dimension(:), pointer :: isAvailable
 
       if (present(ierr)) ierr = MPAS_STREAM_NOERR
 
@@ -238,19 +357,30 @@
 #include &quot;add_field_indices.inc&quot;
 
       
+      any_success = .false.
       if (field % isSuperArray) then
+         allocate(isAvailable(size(field % constituentNames)))
+         isAvailable(:) = .false.
          do i=1,size(field % constituentNames)
             call MPAS_streamAddField_generic(stream, trim(field % constituentNames(i)), MPAS_IO_INT, dimNames0, &amp;
                                              dimSizes0, field % hasTimeDimension, isDecomposed, totalDimSize, globalDimSize, &amp;
                                              indices, io_err)
+            if (io_err == MPAS_STREAM_NOERR) then
+               isAvailable(i) = .true.
+               any_success = .true.
+            end if
          end do
       else
+         nullify(isAvailable)
          call MPAS_streamAddField_generic(stream, trim(field % fieldName), MPAS_IO_INT, field % dimNames, field % dimSizes, &amp;
                                           field % hasTimeDimension, isDecomposed, totalDimSize, globalDimSize, indices, io_err)
+         if (io_err == MPAS_STREAM_NOERR) then
+            any_success = .true.
+         end if
       end if
 
       deallocate(indices)
-      if (io_err /= MPAS_STREAM_NOERR) then
+      if (.not. any_success) then
           if (present(ierr)) ierr = MPAS_IO_ERR
           return
       end if
@@ -273,6 +403,7 @@
       end do
       new_field_list_node % field_type = FIELD_1D_INT
       new_field_list_node % int1dField =&gt; field
+      new_field_list_node % isAvailable =&gt; isAvailable
 
 write(0,*) '... done adding field'
 
@@ -299,6 +430,8 @@
       integer, dimension(:), pointer :: indices
       type (field_list_type), pointer :: field_list_cursor
       type (field_list_type), pointer :: new_field_list_node
+      logical :: any_success
+      logical, dimension(:), pointer :: isAvailable
 
       if (present(ierr)) ierr = MPAS_STREAM_NOERR
 
@@ -323,19 +456,30 @@
 #include &quot;add_field_indices.inc&quot;
 
       
+      any_success = .false.
       if (field % isSuperArray) then
+         allocate(isAvailable(size(field % constituentNames)))
+         isAvailable(:) = .false.
          do i=1,size(field % constituentNames)
             call MPAS_streamAddField_generic(stream, trim(field % constituentNames(i)), MPAS_IO_INT, field % dimNames(2:ndims), &amp;
                                              field % dimSizes(2:ndims), field % hasTimeDimension, isDecomposed, totalDimSize, globalDimSize, &amp;
                                              indices, io_err)
+            if (io_err == MPAS_STREAM_NOERR) then
+               isAvailable(i) = .true.
+               any_success = .true.
+            end if
          end do
       else
+         nullify(isAvailable)
          call MPAS_streamAddField_generic(stream, trim(field % fieldName), MPAS_IO_INT, field % dimNames, field % dimSizes, &amp;
                                           field % hasTimeDimension, isDecomposed, totalDimSize, globalDimSize, indices, io_err)
+         if (io_err == MPAS_STREAM_NOERR) then
+            any_success = .true.
+         end if
       end if
 
       deallocate(indices)
-      if (io_err /= MPAS_STREAM_NOERR) then
+      if (.not. any_success) then
           if (present(ierr)) ierr = MPAS_IO_ERR
           return
       end if
@@ -358,6 +502,7 @@
       end do
       new_field_list_node % field_type = FIELD_2D_INT
       new_field_list_node % int2dField =&gt; field
+      new_field_list_node % isAvailable =&gt; isAvailable
 
 write(0,*) '... done adding field'
 
@@ -384,6 +529,8 @@
       integer, dimension(:), pointer :: indices
       type (field_list_type), pointer :: field_list_cursor
       type (field_list_type), pointer :: new_field_list_node
+      logical :: any_success
+      logical, dimension(:), pointer :: isAvailable
 
       if (present(ierr)) ierr = MPAS_STREAM_NOERR
 
@@ -407,19 +554,30 @@
       ! 
 #include &quot;add_field_indices.inc&quot;
       
+      any_success = .false.
       if (field % isSuperArray) then
+         allocate(isAvailable(size(field % constituentNames)))
+         isAvailable(:) = .false.
          do i=1,size(field % constituentNames)
             call MPAS_streamAddField_generic(stream, trim(field % constituentNames(i)), MPAS_IO_INT, field % dimNames(2:ndims), &amp;
                                              field % dimSizes(2:ndims), field % hasTimeDimension, isDecomposed, totalDimSize, globalDimSize, &amp;
                                              indices, io_err)
+            if (io_err == MPAS_STREAM_NOERR) then
+               isAvailable(i) = .true.
+               any_success = .true.
+            end if
          end do
       else
+         nullify(isAvailable)
          call MPAS_streamAddField_generic(stream, trim(field % fieldName), MPAS_IO_INT, field % dimNames, field % dimSizes, &amp;
                                           field % hasTimeDimension, isDecomposed, totalDimSize, globalDimSize, indices, io_err)
+         if (io_err == MPAS_STREAM_NOERR) then
+            any_success = .true.
+         end if
       end if
 
       deallocate(indices)
-      if (io_err /= MPAS_STREAM_NOERR) then
+      if (.not. any_success) then
           if (present(ierr)) ierr = MPAS_IO_ERR
           return
       end if
@@ -442,6 +600,7 @@
       end do
       new_field_list_node % field_type = FIELD_3D_INT
       new_field_list_node % int3dField =&gt; field
+      new_field_list_node % isAvailable =&gt; isAvailable
 
 write(0,*) '... done adding field'
 
@@ -548,6 +707,8 @@
       integer, dimension(:), pointer :: indices
       type (field_list_type), pointer :: field_list_cursor
       type (field_list_type), pointer :: new_field_list_node
+      logical :: any_success
+      logical, dimension(:), pointer :: isAvailable
 
       if (present(ierr)) ierr = MPAS_STREAM_NOERR
 
@@ -572,19 +733,30 @@
 #include &quot;add_field_indices.inc&quot;
 
       
+      any_success = .false.
       if (field % isSuperArray) then
+         allocate(isAvailable(size(field % constituentNames)))
+         isAvailable(:) = .false.
          do i=1,size(field % constituentNames)
             call MPAS_streamAddField_generic(stream, trim(field % constituentNames(i)), MPAS_IO_DOUBLE, dimNames0, &amp;
                                              dimSizes0, field % hasTimeDimension, isDecomposed, totalDimSize, globalDimSize, &amp;
                                              indices, io_err)
+            if (io_err == MPAS_STREAM_NOERR) then
+               isAvailable(i) = .true.
+               any_success = .true.
+            end if
          end do
       else
+         nullify(isAvailable)
          call MPAS_streamAddField_generic(stream, trim(field % fieldName), MPAS_IO_DOUBLE, field % dimNames, field % dimSizes, &amp;
                                           field % hasTimeDimension, isDecomposed, totalDimSize, globalDimSize, indices, io_err)
+         if (io_err == MPAS_STREAM_NOERR) then
+            any_success = .true.
+         end if
       end if
 
       deallocate(indices)
-      if (io_err /= MPAS_STREAM_NOERR) then
+      if (.not. any_success) then
           if (present(ierr)) ierr = MPAS_IO_ERR
           return
       end if
@@ -607,6 +779,7 @@
       end do
       new_field_list_node % field_type = FIELD_1D_REAL
       new_field_list_node % real1dField =&gt; field
+      new_field_list_node % isAvailable =&gt; isAvailable
 
 write(0,*) '... done adding field'
 
@@ -633,6 +806,8 @@
       integer, dimension(:), pointer :: indices
       type (field_list_type), pointer :: field_list_cursor
       type (field_list_type), pointer :: new_field_list_node
+      logical :: any_success
+      logical, dimension(:), pointer :: isAvailable
 
       if (present(ierr)) ierr = MPAS_STREAM_NOERR
 
@@ -657,19 +832,30 @@
 #include &quot;add_field_indices.inc&quot;
 
       
+      any_success = .false.
       if (field % isSuperArray) then
+         allocate(isAvailable(size(field % constituentNames)))
+         isAvailable(:) = .false.
          do i=1,size(field % constituentNames)
             call MPAS_streamAddField_generic(stream, trim(field % constituentNames(i)), MPAS_IO_DOUBLE, field % dimNames(2:ndims), &amp;
                                              field % dimSizes(2:ndims), field % hasTimeDimension, isDecomposed, totalDimSize, globalDimSize, &amp;
                                              indices, io_err)
+            if (io_err == MPAS_STREAM_NOERR) then
+               isAvailable(i) = .true.
+               any_success = .true.
+            end if
          end do
       else
+         nullify(isAvailable)
          call MPAS_streamAddField_generic(stream, trim(field % fieldName), MPAS_IO_DOUBLE, field % dimNames, field % dimSizes, &amp;
                                           field % hasTimeDimension, isDecomposed, totalDimSize, globalDimSize, indices, io_err)
+         if (io_err == MPAS_STREAM_NOERR) then
+            any_success = .true.
+         end if
       end if
 
       deallocate(indices)
-      if (io_err /= MPAS_STREAM_NOERR) then
+      if (.not. any_success) then
           if (present(ierr)) ierr = MPAS_IO_ERR
           return
       end if
@@ -692,6 +878,7 @@
       end do
       new_field_list_node % field_type = FIELD_2D_REAL
       new_field_list_node % real2dField =&gt; field
+      new_field_list_node % isAvailable =&gt; isAvailable
 
 write(0,*) '... done adding field'
 
@@ -718,6 +905,8 @@
       integer, dimension(:), pointer :: indices
       type (field_list_type), pointer :: field_list_cursor
       type (field_list_type), pointer :: new_field_list_node
+      logical :: any_success
+      logical, dimension(:), pointer :: isAvailable
 
       if (present(ierr)) ierr = MPAS_STREAM_NOERR
 
@@ -742,20 +931,31 @@
 #include &quot;add_field_indices.inc&quot;
 
       
+      any_success = .false.
       if (field % isSuperArray) then
 write(0,*) '^^^^^^^^^^^^^^^^^^^^^^^^^^^ we are adding a super-array'
+         allocate(isAvailable(size(field % constituentNames)))
+         isAvailable(:) = .false.
          do i=1,size(field % constituentNames)
             call MPAS_streamAddField_generic(stream, trim(field % constituentNames(i)), MPAS_IO_DOUBLE, field % dimNames(2:ndims), &amp;
                                              field % dimSizes(2:ndims), field % hasTimeDimension, isDecomposed, totalDimSize, globalDimSize, &amp;
                                              indices, io_err)
+            if (io_err == MPAS_STREAM_NOERR) then
+               isAvailable(i) = .true.
+               any_success = .true.
+            end if
          end do
       else
+         nullify(isAvailable)
          call MPAS_streamAddField_generic(stream, trim(field % fieldName), MPAS_IO_DOUBLE, field % dimNames, field % dimSizes, &amp;
                                           field % hasTimeDimension, isDecomposed, totalDimSize, globalDimSize, indices, io_err)
+         if (io_err == MPAS_STREAM_NOERR) then
+            any_success = .true.
+         end if
       end if
 
       deallocate(indices)
-      if (io_err /= MPAS_STREAM_NOERR) then
+      if (.not. any_success) then
           if (present(ierr)) ierr = MPAS_IO_ERR
           return
       end if
@@ -778,8 +978,10 @@
       end do
       new_field_list_node % field_type = FIELD_3D_REAL
       new_field_list_node % real3dField =&gt; field
+      new_field_list_node % isAvailable =&gt; isAvailable
 
 write(0,*) '... done adding field'
+write(0,*) 'DEBUGGING : Finished adding 3d real field '//trim(field % fieldName)
 
    end subroutine MPAS_streamAddField_3dReal
 
@@ -1108,6 +1310,7 @@
       do while (associated(field_cursor))
          if (field_cursor % field_type == FIELD_0D_INT) then
 
+write(0,*) 'DEBUGGING : *************** '//trim(field_cursor % int0dField % fieldName)
 write(0,*) 'Reading in field '//trim(field_cursor % int0dField % fieldName)
 write(0,*) '   &gt; is the field decomposed? ', field_cursor % isDecomposed
 write(0,*) '   &gt; outer dimension size ', field_cursor % totalDimSize
@@ -1130,6 +1333,7 @@
             end do
 
          else if (field_cursor % field_type == FIELD_1D_INT) then
+write(0,*) 'DEBUGGING : *************** '//trim(field_cursor % int1dField % fieldName)
 
             if (field_cursor % int1dField % isSuperArray) then
                ncons = size(field_cursor % int1dField % constituentNames)
@@ -1139,6 +1343,7 @@
             end if
 
             do j=1,ncons
+               if ((field_cursor % int1dField % isSuperArray) .and. (.not. field_cursor % isAvailable(j))) cycle
                if (field_cursor % int1dField % isSuperArray) then
                   call MPAS_io_get_var(stream % fileHandle, field_cursor % int1dField % constituentNames(j), int0d_temp, io_err)
                else
@@ -1202,6 +1407,7 @@
 
          else if (field_cursor % field_type == FIELD_2D_INT) then
 
+write(0,*) 'DEBUGGING : *************** '//trim(field_cursor % int2dField % fieldName)
             if (field_cursor % int2dField % isSuperArray) then
                ncons = size(field_cursor % int2dField % constituentNames)
                allocate(int1d_temp(field_cursor % totalDimSize))
@@ -1212,6 +1418,7 @@
             end if
 
             do j=1,ncons
+               if ((field_cursor % int2dField % isSuperArray) .and. (.not. field_cursor % isAvailable(j))) cycle
                if (field_cursor % int2dField % isSuperArray) then
                   call MPAS_io_get_var(stream % fileHandle, field_cursor % int2dField % constituentNames(j), int1d_temp, io_err)
                else
@@ -1279,6 +1486,7 @@
 
          else if (field_cursor % field_type == FIELD_3D_INT) then
 
+write(0,*) 'DEBUGGING : *************** '//trim(field_cursor % int3dField % fieldName)
             if (field_cursor % int3dField % isSuperArray) then
                ncons = size(field_cursor % int3dField % constituentNames)
                allocate(int2d_temp(field_cursor % int3dField % dimSizes(2), &amp;
@@ -1291,6 +1499,7 @@
             end if
 
             do j=1,ncons
+               if ((field_cursor % int3dField % isSuperArray) .and. (.not. field_cursor % isAvailable(j))) cycle
                if (field_cursor % int3dField % isSuperArray) then
                   call MPAS_io_get_var(stream % fileHandle, field_cursor % int3dField % constituentNames(j), int2d_temp, io_err)
                else
@@ -1358,6 +1567,7 @@
 
          else if (field_cursor % field_type == FIELD_0D_REAL) then
 
+write(0,*) 'DEBUGGING : *************** '//trim(field_cursor % real0dField % fieldName)
 write(0,*) 'Reading in field '//trim(field_cursor % real0dField % fieldName)
 write(0,*) '   &gt; is the field decomposed? ', field_cursor % isDecomposed
 write(0,*) '   &gt; outer dimension size ', field_cursor % totalDimSize
@@ -1381,6 +1591,7 @@
 
          else if (field_cursor % field_type == FIELD_1D_REAL) then
 
+write(0,*) 'DEBUGGING : *************** '//trim(field_cursor % real1dField % fieldName)
             if (field_cursor % real1dField % isSuperArray) then
                ncons = size(field_cursor % real1dField % constituentNames)
             else
@@ -1389,6 +1600,7 @@
             end if
 
             do j=1,ncons
+               if ((field_cursor % real1dField % isSuperArray) .and. (.not. field_cursor % isAvailable(j))) cycle
                if (field_cursor % real1dField % isSuperArray) then
                   call MPAS_io_get_var(stream % fileHandle, field_cursor % real1dField % constituentNames(j), real0d_temp, io_err)
                else
@@ -1452,6 +1664,7 @@
 
          else if (field_cursor % field_type == FIELD_2D_REAL) then
 
+write(0,*) 'DEBUGGING : *************** '//trim(field_cursor % real2dField % fieldName)
             if (field_cursor % real2dField % isSuperArray) then
                ncons = size(field_cursor % real2dField % constituentNames)
                allocate(real1d_temp(field_cursor % totalDimSize))
@@ -1462,6 +1675,7 @@
             end if
 
             do j=1,ncons
+               if ((field_cursor % real2dField % isSuperArray) .and. (.not. field_cursor % isAvailable(j))) cycle
                if (field_cursor % real2dField % isSuperArray) then
                   call MPAS_io_get_var(stream % fileHandle, field_cursor % real2dField % constituentNames(j), real1d_temp, io_err)
                else
@@ -1529,7 +1743,10 @@
 
          else if (field_cursor % field_type == FIELD_3D_REAL) then
 
+write(0,*) 'DEBUGGING : *************** '//trim(field_cursor % real3dField % fieldName)
+write(0,*) 'DEBUGGING : reading a 3d real array'
             if (field_cursor % real3dField % isSuperArray) then
+write(0,*) 'DEBUGGING : reading a 3d real super-array'
                ncons = size(field_cursor % real3dField % constituentNames)
                allocate(real2d_temp(field_cursor % real3dField % dimSizes(2), &amp;
                                     field_cursor % totalDimSize))
@@ -1541,7 +1758,9 @@
             end if
 
             do j=1,ncons
+               if ((field_cursor % real3dField % isSuperArray) .and. (.not. field_cursor % isAvailable(j))) cycle
                if (field_cursor % real3dField % isSuperArray) then
+write(0,*) 'DEBUGGING : calling get_var for a constitutent'
                   call MPAS_io_get_var(stream % fileHandle, field_cursor % real3dField % constituentNames(j), real2d_temp, io_err)
                else
                   call MPAS_io_get_var(stream % fileHandle, field_cursor % real3dField % fieldName, real3d_temp, io_err)
@@ -1572,6 +1791,7 @@
                   end if
                   do while (associated(field_3dreal_ptr))
                      if (field_cursor % real3dField % isSuperArray) then
+write(0,*) 'DEBUGGING : copying the temporary array'
                         field_3dreal_ptr % array(j,:,1:ownedSize) = real2d_temp(:,i:i+ownedSize-1)
                      else
                         field_3dreal_ptr % array(:,:,1:ownedSize) = real3d_temp(:,:,i:i+ownedSize-1)
@@ -1622,7 +1842,7 @@
 
 write(0,*) 'Distributing and Copying field to other blocks'
 
-            call mpas_dmpar_bcast_char(field_cursor % int0dField % block % domain % dminfo, field_cursor % char0dField % scalar)
+            call mpas_dmpar_bcast_char(field_cursor % char0dField % block % domain % dminfo, field_cursor % char0dField % scalar)
             field_0dchar_ptr =&gt; field_cursor % char0dField
             do while (associated(field_0dchar_ptr))
                field_0dchar_ptr % scalar = field_cursor % char0dField % scalar
@@ -2400,6 +2620,10 @@
 write(0,*) 'Deallocating field list'
       field_cursor =&gt; stream % fieldList
       do while (associated(field_cursor))
+         if (associated(field_cursor % isAvailable)) then
+            deallocate(field_cursor % isAvailable)
+write(0,*) 'Deallocating isAvailable array'
+         end if
          stream % fieldList =&gt; stream % fieldList % next
          deallocate(field_cursor)
          field_cursor =&gt; stream % fieldList

</font>
</pre>