<p><b>duda</b> 2010-05-21 11:20:37 -0600 (Fri, 21 May 2010)</p><p>Add code in module_io_input.F and module_io_output.F to read and write <br>
two new global attributes, on_a_sphere and sphere_radius. If these attributes <br>
are not present in an input file, on_a_sphere will default to "YES" and <br>
sphere_radius will default to 1.0.<br>
<br>
The values of these attributes are accessible in the grid data structure;<br>
that is, they can be accessed in the code as grid % on_a_sphere<br>
(or mesh % on_a_sphere in some places) and grid % sphere_radius.<br>
<br>
Note: In the grid data structure, on_a_sphere is a logical, rather than<br>
a character string, as it is stored in the netCDF files.<br>
<br>
<br>
M module_grid_types.F<br>
M module_io_input.F<br>
M module_io_output.F<br>
</p><hr noshade><pre><font color="gray">Modified: trunk/mpas/src/framework/module_grid_types.F
===================================================================
--- trunk/mpas/src/framework/module_grid_types.F        2010-05-20 23:43:36 UTC (rev 297)
+++ trunk/mpas/src/framework/module_grid_types.F        2010-05-21 17:20:37 UTC (rev 298)
@@ -70,6 +70,9 @@
#include "field_dimensions.inc"
+ logical :: on_a_sphere
+ real (kind=RKIND) :: sphere_radius
+
#include "time_invariant_fields.inc"
end type grid_meta
Modified: trunk/mpas/src/framework/module_io_input.F
===================================================================
--- trunk/mpas/src/framework/module_io_input.F        2010-05-20 23:43:36 UTC (rev 297)
+++ trunk/mpas/src/framework/module_io_input.F        2010-05-21 17:20:37 UTC (rev 298)
@@ -49,6 +49,9 @@
integer :: i, j, k
type (io_input_object) :: input_obj
#include "dim_decls.inc"
+
+ character (len=16) :: c_on_a_sphere
+ real (kind=RKIND) :: r_sphere_radius
integer :: readCellStart, readCellEnd, nReadCells
integer :: readEdgeStart, readEdgeEnd, nReadEdges
@@ -729,6 +732,18 @@
#include "dim_dummy_args.inc"
)
+ !
+ ! Read attributes
+ !
+ call io_input_get_att_text(input_obj, 'on_a_sphere', c_on_a_sphere)
+ call io_input_get_att_real(input_obj, 'sphere_radius', r_sphere_radius)
+ if (index(c_on_a_sphere, 'YES') /= 0) then
+ domain % blocklist % mesh % on_a_sphere = .true.
+ else
+ domain % blocklist % mesh % on_a_sphere = .false.
+ end if
+ domain % blocklist % mesh % sphere_radius = r_sphere_radius
+
if (.not. config_do_restart) then
input_obj % time = 1
else
@@ -1100,7 +1115,59 @@
end subroutine io_input_get_dimension
+
+ subroutine io_input_get_att_real(input_obj, attname, attvalue)
+
+ implicit none
+ type (io_input_object), intent(in) :: input_obj
+ character (len=*), intent(in) :: attname
+ real (kind=RKIND), intent(out) :: attvalue
+
+ include 'netcdf.inc'
+
+ integer :: nferr
+
+ if (RKIND == 8) then
+ nferr = nf_get_att_double(input_obj % rd_ncid, NF_GLOBAL, attname, attvalue)
+ else
+ nferr = nf_get_att_real(input_obj % rd_ncid, NF_GLOBAL, attname, attvalue)
+ end if
+ if (nferr /= NF_NOERR) then
+ write(0,*) 'Warning: Attribute '//trim(attname)//' not found in '//trim(input_obj % filename)
+ if (index(attname, 'sphere_radius') /= 0) then
+ write(0,*) ' Setting '//trim(attname)//' to 1.0'
+ attvalue = 1.0
+ end if
+ end if
+
+ end subroutine io_input_get_att_real
+
+
+ subroutine io_input_get_att_text(input_obj, attname, attvalue)
+
+ implicit none
+
+ type (io_input_object), intent(in) :: input_obj
+ character (len=*), intent(in) :: attname
+ character (len=*), intent(out) :: attvalue
+
+ include 'netcdf.inc'
+
+ integer :: nferr
+
+ nferr = nf_get_att_text(input_obj % rd_ncid, NF_GLOBAL, attname, attvalue)
+ if (nferr /= NF_NOERR) then
+ write(0,*) 'Warning: Attribute '//trim(attname)//' not found in '//trim(input_obj % filename)
+ if (index(attname, 'on_a_sphere') /= 0) then
+ write(0,*) ' Setting '//trim(attname)//' to ''YES'''
+ attvalue = 'YES'
+ end if
+ end if
+
+ end subroutine io_input_get_att_text
+
+
subroutine io_input_field0dReal(input_obj, field)
implicit none
Modified: trunk/mpas/src/framework/module_io_output.F
===================================================================
--- trunk/mpas/src/framework/module_io_output.F        2010-05-20 23:43:36 UTC (rev 297)
+++ trunk/mpas/src/framework/module_io_output.F        2010-05-21 17:20:37 UTC (rev 298)
@@ -90,6 +90,7 @@
! although in future, work needs to be done to write model state
! from many distributed blocks
call io_output_init(output_obj, domain % dminfo, &
+ block_ptr % mesh, &
#include "output_dim_actual_args.inc"
)
@@ -326,6 +327,7 @@
subroutine io_output_init( output_obj, &
dminfo, &
+ mesh, &
#include "dim_dummy_args.inc"
)
@@ -335,6 +337,7 @@
type (io_output_object), intent(inout) :: output_obj
type (dm_info), intent(in) :: dminfo
+ type (grid_meta), intent(in) :: mesh
#include "dim_dummy_decls.inc"
integer :: nferr
@@ -348,6 +351,17 @@
#endif
#include "netcdf_def_dims_vars.inc"
+
+ if (mesh % on_a_sphere) then
+ nferr = nf_put_att_text(output_obj % wr_ncid, NF_GLOBAL, 'on_a_sphere', 16, 'YES ')
+ else
+ nferr = nf_put_att_text(output_obj % wr_ncid, NF_GLOBAL, 'on_a_sphere', 16, 'NO ')
+ end if
+ if (RKIND == 8) then
+ nferr = nf_put_att_double(output_obj % wr_ncid, NF_GLOBAL, 'sphere_radius', NF_DOUBLE, 1, mesh % sphere_radius)
+ else
+ nferr = nf_put_att_real(output_obj % wr_ncid, NF_GLOBAL, 'sphere_radius', NF_FLOAT, 1, mesh % sphere_radius)
+ end if
nferr = nf_enddef(output_obj % wr_ncid)
end if
</font>
</pre>