<p><b>dwj07@fsu.edu</b> 2012-10-24 14:22:05 -0600 (Wed, 24 Oct 2012)</p><p><br>
        -- TRUNK COMMIT --<br>
<br>
        Addition of scratch variables.<br>
<br>
        Now within registry you can specify scratch as the persistence of variables in addition to persistent.<br>
        This causes the % array pointer of the field to not be allocated upon start-up.<br>
<br>
        Routines are added to allocate and deallocate the % array pointer of a field. Either for all blocks in blocklist or for a single block.<br>
        These routines are:<br>
                mpas_allocate_scratch_field<br>
                mpas_deallocate_scratch_field<br>
<br>
        Scratch variables cannot be in a super array.<br>
<br>
<br>
        NOTE: If you want to create an additional variable group (state, mesh, etc) for scratch variables <br>
              it needs to have more than 1 variable or else registry will fail.<br>
</p><hr noshade><pre><font color="gray">Modified: trunk/mpas/src/framework/mpas_grid_types.F
===================================================================
--- trunk/mpas/src/framework/mpas_grid_types.F        2012-10-24 17:07:59 UTC (rev 2264)
+++ trunk/mpas/src/framework/mpas_grid_types.F        2012-10-24 20:22:05 UTC (rev 2265)
@@ -363,6 +363,26 @@
       type (dm_info), pointer :: dminfo
    end type domain_type
 
+   interface mpas_allocate_scratch_field
+     module procedure mpas_allocate_scratch_field1d_integer
+     module procedure mpas_allocate_scratch_field2d_integer
+     module procedure mpas_allocate_scratch_field3d_integer
+     module procedure mpas_allocate_scratch_field1d_real
+     module procedure mpas_allocate_scratch_field2d_real
+     module procedure mpas_allocate_scratch_field3d_real
+     module procedure mpas_allocate_scratch_field1d_char
+   end interface
+
+   interface mpas_deallocate_scratch_field
+     module procedure mpas_deallocate_scratch_field1d_integer
+     module procedure mpas_deallocate_scratch_field2d_integer
+     module procedure mpas_deallocate_scratch_field3d_integer
+     module procedure mpas_deallocate_scratch_field1d_real
+     module procedure mpas_deallocate_scratch_field2d_real
+     module procedure mpas_deallocate_scratch_field3d_real
+     module procedure mpas_deallocate_scratch_field1d_char
+   end interface
+
    interface mpas_deallocate_field
      module procedure mpas_deallocate_field0d_integer
      module procedure mpas_deallocate_field1d_integer
@@ -444,6 +464,406 @@
 
    end subroutine mpas_deallocate_domain!}}}
 
+   subroutine mpas_allocate_scratch_field1d_integer(f, single_block_in)!{{{
+       type (field1dInteger), pointer :: f
+       logical, intent(in), optional :: single_block_in
+       logical :: single_block
+       type (field1dInteger), pointer :: f_cursor
+
+       if(present(single_block_in)) then
+          single_block = single_block_in
+       else
+          single_block = .false.
+       end if
+
+       if(.not. single_block) then
+          f_cursor =&gt; f
+          do while(associated(f_cursor))
+             if(.not.associated(f_cursor % array)) then
+                allocate(f_cursor % array(f_cursor % dimSizes(1)))
+             end if
+             f_cursor =&gt; f_cursor % next
+          end do
+       else
+          if(.not.associated(f % array)) then
+            allocate(f % array(f % dimSizes(1)))
+          end if
+       end if
+
+   end subroutine mpas_allocate_scratch_field1d_integer!}}}
+
+   subroutine mpas_allocate_scratch_field2d_integer(f, single_block_in)!{{{
+       type (field2dInteger), pointer :: f
+       logical, intent(in), optional :: single_block_in
+       logical :: single_block
+       type (field2dInteger), pointer :: f_cursor
+
+       if(present(single_block_in)) then
+          single_block = single_block_in
+       else
+          single_block = .false.
+       end if
+
+       if(.not. single_block) then
+          f_cursor =&gt; f
+          do while(associated(f_cursor))
+             if(.not.associated(f_cursor % array)) then
+                allocate(f_cursor % array(f_cursor % dimSizes(1), f_cursor % dimSizes(2)))
+             end if
+             f_cursor =&gt; f_cursor % next
+          end do
+       else
+          if(.not.associated(f % array)) then
+            allocate(f % array(f % dimSizes(1), f % dimSizes(2)))
+          end if
+       end if
+
+   end subroutine mpas_allocate_scratch_field2d_integer!}}}
+
+   subroutine mpas_allocate_scratch_field3d_integer(f, single_block_in)!{{{
+       type (field3dInteger), pointer :: f
+       logical, intent(in), optional :: single_block_in
+       logical :: single_block
+       type (field3dInteger), pointer :: f_cursor
+
+       if(present(single_block_in)) then
+          single_block = single_block_in
+       else
+          single_block = .false.
+       end if
+
+       if(.not. single_block) then
+          f_cursor =&gt; f
+          do while(associated(f_cursor))
+             if(.not.associated(f_cursor % array)) then
+                allocate(f_cursor % array(f_cursor % dimSizes(1), f_cursor % dimSizes(2), f_cursor % dimSizes(3)))
+             end if
+             f_cursor =&gt; f_cursor % next
+          end do
+       else
+          if(.not.associated(f % array)) then
+            allocate(f % array(f % dimSizes(1), f % dimSizes(2), f % dimSizes(3)))
+          end if
+       end if
+
+   end subroutine mpas_allocate_scratch_field3d_integer!}}}
+
+   subroutine mpas_allocate_scratch_field1d_real(f, single_block_in)!{{{
+       type (field1dReal), pointer :: f
+       logical, intent(in), optional :: single_block_in
+       logical :: single_block
+       type (field1dReal), pointer :: f_cursor
+
+       if(present(single_block_in)) then
+          single_block = single_block_in
+       else
+          single_block = .false.
+       end if
+
+       if(.not. single_block) then
+          f_cursor =&gt; f
+          do while(associated(f_cursor))
+             if(.not.associated(f_cursor % array)) then
+                allocate(f_cursor % array(f_cursor % dimSizes(1)))
+             end if
+             f_cursor =&gt; f_cursor % next
+          end do
+       else
+          if(.not.associated(f % array)) then
+            allocate(f % array(f % dimSizes(1)))
+          end if
+       end if
+
+   end subroutine mpas_allocate_scratch_field1d_real!}}}
+
+   subroutine mpas_allocate_scratch_field2d_real(f, single_block_in)!{{{
+       type (field2dReal), pointer :: f
+       logical, intent(in), optional :: single_block_in
+       logical :: single_block
+       type (field2dReal), pointer :: f_cursor
+
+       if(present(single_block_in)) then
+          single_block = single_block_in
+       else
+          single_block = .false.
+       end if
+
+       if(.not. single_block) then
+          f_cursor =&gt; f
+          do while(associated(f_cursor))
+             if(.not.associated(f_cursor % array)) then
+                allocate(f_cursor % array(f_cursor % dimSizes(1), f_cursor % dimSizes(2)))
+             end if
+             f_cursor =&gt; f_cursor % next
+          end do
+       else
+          if(.not.associated(f % array)) then
+            allocate(f % array(f % dimSizes(1), f % dimSizes(2)))
+          end if
+       end if
+
+   end subroutine mpas_allocate_scratch_field2d_real!}}}
+
+   subroutine mpas_allocate_scratch_field3d_real(f, single_block_in)!{{{
+       type (field3dReal), pointer :: f
+       logical, intent(in), optional :: single_block_in
+       logical :: single_block
+       type (field3dReal), pointer :: f_cursor
+
+       if(present(single_block_in)) then
+          single_block = single_block_in
+       else
+          single_block = .false.
+       end if
+
+       if(.not. single_block) then
+          f_cursor =&gt; f
+          do while(associated(f_cursor))
+             if(.not.associated(f_cursor % array)) then
+                allocate(f_cursor % array(f_cursor % dimSizes(1), f_cursor % dimSizes(2), f_cursor % dimSizes(3)))
+             end if
+             f_cursor =&gt; f_cursor % next
+          end do
+       else
+          if(.not.associated(f % array)) then
+            allocate(f % array(f % dimSizes(1), f % dimSizes(2), f % dimSizes(3)))
+          end if
+       end if
+
+   end subroutine mpas_allocate_scratch_field3d_real!}}}
+
+   subroutine mpas_allocate_scratch_field1d_char(f, single_block_in)!{{{
+       type (field1dChar), pointer :: f
+       logical, intent(in), optional :: single_block_in
+       logical :: single_block
+       type (field1dChar), pointer :: f_cursor
+
+       if(present(single_block_in)) then
+          single_block = single_block_in
+       else
+          single_block = .false.
+       end if
+
+       if(.not. single_block) then
+          f_cursor =&gt; f
+          do while(associated(f_cursor))
+             if(.not.associated(f_cursor % array)) then
+                allocate(f_cursor % array(f_cursor % dimSizes(1)))
+             end if
+             f_cursor =&gt; f_cursor % next
+          end do
+       else
+          if(.not.associated(f % array)) then
+            allocate(f % array(f % dimSizes(1)))
+          end if
+       end if
+
+   end subroutine mpas_allocate_scratch_field1d_char!}}}
+
+   subroutine mpas_deallocate_scratch_field1d_integer(f, single_block_in)!{{{
+       type (field1dInteger), pointer :: f
+       logical, intent(in), optional :: single_block_in
+       logical :: single_block
+       type (field1dInteger), pointer :: f_cursor
+
+       if(present(single_block_in)) then
+          single_block = single_block_in
+       else
+          single_block = .false.
+       end if
+
+       if(.not.single_block) then
+          f_cursor =&gt; f
+          do while(associated(f_cursor))
+            if(associated(f_cursor % array)) then
+              deallocate(f_cursor % array)
+            end if
+   
+            f_cursor =&gt; f_cursor % next
+          end do
+       else
+          if(associated(f % array)) then
+             deallocate(f % array)
+          end if
+       end if
+
+   end subroutine mpas_deallocate_scratch_field1d_integer!}}}
+
+   subroutine mpas_deallocate_scratch_field2d_integer(f, single_block_in)!{{{
+       type (field2dInteger), pointer :: f
+       logical, intent(in), optional :: single_block_in
+       logical :: single_block
+       type (field2dInteger), pointer :: f_cursor
+
+       if(present(single_block_in)) then
+          single_block = single_block_in
+       else
+          single_block = .false.
+       end if
+
+       if(.not.single_block) then
+          f_cursor =&gt; f
+          do while(associated(f_cursor))
+            if(associated(f_cursor % array)) then
+              deallocate(f_cursor % array)
+            end if
+   
+            f_cursor =&gt; f_cursor % next
+          end do
+       else
+          if(associated(f % array)) then
+             deallocate(f % array)
+          end if
+       end if
+
+   end subroutine mpas_deallocate_scratch_field2d_integer!}}}
+
+   subroutine mpas_deallocate_scratch_field3d_integer(f, single_block_in)!{{{
+       type (field3dInteger), pointer :: f
+       logical, intent(in), optional :: single_block_in
+       logical :: single_block
+       type (field3dInteger), pointer :: f_cursor
+
+       if(present(single_block_in)) then
+          single_block = single_block_in
+       else
+          single_block = .false.
+       end if
+
+       if(.not.single_block) then
+          f_cursor =&gt; f
+          do while(associated(f_cursor))
+            if(associated(f_cursor % array)) then
+              deallocate(f_cursor % array)
+            end if
+   
+            f_cursor =&gt; f_cursor % next
+          end do
+       else
+          if(associated(f % array)) then
+             deallocate(f % array)
+          end if
+       end if
+
+   end subroutine mpas_deallocate_scratch_field3d_integer!}}}
+
+   subroutine mpas_deallocate_scratch_field1d_real(f, single_block_in)!{{{
+       type (field1dReal), pointer :: f
+       logical, intent(in), optional :: single_block_in
+       logical :: single_block
+       type (field1dReal), pointer :: f_cursor
+
+       if(present(single_block_in)) then
+          single_block = single_block_in
+       else
+          single_block = .false.
+       end if
+
+       if(.not.single_block) then
+          f_cursor =&gt; f
+          do while(associated(f_cursor))
+            if(associated(f_cursor % array)) then
+              deallocate(f_cursor % array)
+            end if
+   
+            f_cursor =&gt; f_cursor % next
+          end do
+       else
+          if(associated(f % array)) then
+             deallocate(f % array)
+          end if
+       end if
+
+   end subroutine mpas_deallocate_scratch_field1d_real!}}}
+
+   subroutine mpas_deallocate_scratch_field2d_real(f, single_block_in)!{{{
+       type (field2dReal), pointer :: f
+       logical, intent(in), optional :: single_block_in
+       logical :: single_block
+       type (field2dReal), pointer :: f_cursor
+
+       if(present(single_block_in)) then
+          single_block = single_block_in
+       else
+          single_block = .false.
+       end if
+
+       if(.not.single_block) then
+          f_cursor =&gt; f
+          do while(associated(f_cursor))
+            if(associated(f_cursor % array)) then
+              deallocate(f_cursor % array)
+            end if
+   
+            f_cursor =&gt; f_cursor % next
+          end do
+       else
+          if(associated(f % array)) then
+             deallocate(f % array)
+          end if
+       end if
+
+   end subroutine mpas_deallocate_scratch_field2d_real!}}}
+
+   subroutine mpas_deallocate_scratch_field3d_real(f, single_block_in)!{{{
+       type (field3dReal), pointer :: f
+       logical, intent(in), optional :: single_block_in
+       logical :: single_block
+       type (field3dReal), pointer :: f_cursor
+
+       if(present(single_block_in)) then
+          single_block = single_block_in
+       else
+          single_block = .false.
+       end if
+
+       if(.not.single_block) then
+          f_cursor =&gt; f
+          do while(associated(f_cursor))
+            if(associated(f_cursor % array)) then
+              deallocate(f_cursor % array)
+            end if
+   
+            f_cursor =&gt; f_cursor % next
+          end do
+       else
+          if(associated(f % array)) then
+             deallocate(f % array)
+          end if
+       end if
+
+   end subroutine mpas_deallocate_scratch_field3d_real!}}}
+
+   subroutine mpas_deallocate_scratch_field1d_char(f, single_block_in)!{{{
+       type (field1dChar), pointer :: f
+       logical, intent(in), optional :: single_block_in
+       logical :: single_block
+       type (field1dChar), pointer :: f_cursor
+
+       if(present(single_block_in)) then
+          single_block = single_block_in
+       else
+          single_block = .false.
+       end if
+
+       if(.not.single_block) then
+          f_cursor =&gt; f
+          do while(associated(f_cursor))
+            if(associated(f_cursor % array)) then
+              deallocate(f_cursor % array)
+            end if
+   
+            f_cursor =&gt; f_cursor % next
+          end do
+       else
+          if(associated(f % array)) then
+             deallocate(f % array)
+          end if
+       end if
+
+   end subroutine mpas_deallocate_scratch_field1d_char!}}}
+
+
    subroutine mpas_deallocate_field0d_integer(f)!{{{
        type (field0dInteger), pointer :: f
        type (field0dInteger), pointer :: f_cursor

Modified: trunk/mpas/src/registry/gen_inc.c
===================================================================
--- trunk/mpas/src/registry/gen_inc.c        2012-10-24 17:07:59 UTC (rev 2264)
+++ trunk/mpas/src/registry/gen_inc.c        2012-10-24 20:22:05 UTC (rev 2265)
@@ -813,6 +813,10 @@
             fortprintf(fd, &quot;      %s %% %s %% fieldName = \'%s\'</font>
<font color="black">&quot;, group_ptr-&gt;name, var_ptr-&gt;name_in_code, var_ptr-&gt;name_in_file);
             fortprintf(fd, &quot;      %s %% %s %% isSuperArray = .false.</font>
<font color="blue">&quot;, group_ptr-&gt;name, var_ptr-&gt;name_in_code);
             if (var_ptr-&gt;ndims &gt; 0) {
+                            if(var_ptr-&gt;persistence == SCRATCH){
+                                  fortprintf(fd, &quot;      ! SCRATCH VARIABLE</font>
<font color="blue">&quot;);
+                                  fortprintf(fd, &quot;      nullify(%s %% %s %% array)</font>
<font color="gray">&quot;, group_ptr-&gt;name, var_ptr-&gt;name_in_code); 
+                          } else if(var_ptr-&gt;persistence == PERSISTENT){
                fortprintf(fd, &quot;      allocate(%s %% %s %% array(&quot;, group_ptr-&gt;name, var_ptr-&gt;name_in_code);
                dimlist_ptr = var_ptr-&gt;dimlist;
                if (!strncmp(dimlist_ptr-&gt;dim-&gt;name_in_file, &quot;nCells&quot;, 1024) ||
@@ -843,6 +847,7 @@
                else if (var_ptr-&gt;vtype == CHARACTER)
                   fortprintf(fd, &quot;      %s %% %s %% array = \'\'</font>
<font color="gray">&quot;, group_ptr-&gt;name, var_ptr-&gt;name_in_code ); /* initialize field to zero */
 
+                          }
                dimlist_ptr = var_ptr-&gt;dimlist;
                i = 1;
                while (dimlist_ptr) {
@@ -869,7 +874,7 @@
                   i++;
                   dimlist_ptr = dimlist_ptr-&gt;next;
                }
-            }
+                        }
 
             if (var_ptr-&gt;timedim) fortprintf(fd, &quot;      %s %% %s %% hasTimeDimension = .true.</font>
<font color="black">&quot;, group_ptr-&gt;name, var_ptr-&gt;name_in_code);
             else fortprintf(fd, &quot;      %s %% %s %% hasTimeDimension = .false.</font>
<font color="gray">&quot;, group_ptr-&gt;name, var_ptr-&gt;name_in_code);
@@ -934,14 +939,18 @@
                var_list_ptr2 = var_list_ptr;
                var_list_ptr = var_list_ptr-&gt;next;
             }
-            fortprintf(fd, &quot;      deallocate(%s %% %s %% array)</font>
<font color="blue">&quot;, group_ptr-&gt;name, var_list_ptr2-&gt;var-&gt;super_array);
+            fortprintf(fd, &quot;      if(associated(%s %% %s %% array)) then</font>
<font color="blue">&quot;, group_ptr-&gt;name, var_list_ptr2-&gt;var-&gt;super_array);
+            fortprintf(fd, &quot;         deallocate(%s %% %s %% array)</font>
<font color="blue">&quot;, group_ptr-&gt;name, var_list_ptr2-&gt;var-&gt;super_array);
+            fortprintf(fd, &quot;      end if</font>
<font color="black">&quot;);
             fortprintf(fd, &quot;      deallocate(%s %% %s %% ioinfo)</font>
<font color="black">&quot;, group_ptr-&gt;name, var_list_ptr2-&gt;var-&gt;super_array);
             fortprintf(fd, &quot;      call mpas_deallocate_attlist(%s %% %s %% attList)</font>
<font color="black">&quot;, group_ptr-&gt;name, var_list_ptr2-&gt;var-&gt;super_array);
             fortprintf(fd, &quot;      deallocate(%s %% %s)</font>
<font color="black"></font>
<font color="red">&quot;, group_ptr-&gt;name, var_list_ptr2-&gt;var-&gt;super_array);
          }
          else {
             if (var_ptr-&gt;ndims &gt; 0) {
-               fortprintf(fd, &quot;      deallocate(%s %% %s %% array)</font>
<font color="blue">&quot;, group_ptr-&gt;name, var_ptr-&gt;name_in_code);
+               fortprintf(fd, &quot;      if(associated(%s %% %s %% array)) then</font>
<font color="blue">&quot;, group_ptr-&gt;name, var_ptr-&gt;name_in_code);
+               fortprintf(fd, &quot;         deallocate(%s %% %s %% array)</font>
<font color="blue">&quot;, group_ptr-&gt;name, var_ptr-&gt;name_in_code);
+               fortprintf(fd, &quot;      end if</font>
<font color="black">&quot;);
                fortprintf(fd, &quot;      deallocate(%s %% %s %% ioinfo)</font>
<font color="black">&quot;, group_ptr-&gt;name, var_ptr-&gt;name_in_code);
                fortprintf(fd, &quot;      call mpas_deallocate_attlist(%s %% %s %% attList)</font>
<font color="black">&quot;, group_ptr-&gt;name, var_ptr-&gt;name_in_code);
                fortprintf(fd, &quot;      deallocate(%s %% %s)</font>
<font color="black"></font>
<font color="gray">&quot;, group_ptr-&gt;name, var_ptr-&gt;name_in_code);
@@ -1111,7 +1120,6 @@
      var_list_ptr = group_ptr-&gt;vlist;
      var_list_ptr = var_list_ptr-&gt;next;
      var_ptr = var_list_ptr-&gt;var;
-
      
      int ntime_levs = 1;
      
@@ -2126,6 +2134,7 @@
 
          dimlist_ptr = var_ptr-&gt;dimlist;
          i = 1;
+                 if(var_ptr-&gt;persistence == PERSISTENT){
          while (dimlist_ptr) {
             if (i == var_ptr-&gt;ndims) { 
 
@@ -2172,6 +2181,7 @@
             i++;
             dimlist_ptr = dimlist_ptr -&gt; next;
          }
+                 }
 
          if (var_list_ptr) var_list_ptr = var_list_ptr-&gt;next;
       }

</font>
</pre>