fix: improve modern diag manager performance (#1634)
[FMS.git] / test_fms / diag_manager / test_diag_manager.F90
blobdd263d2e3fb6e11a62402805ba48b201880ec856
1 !***********************************************************************
2 !*                   GNU Lesser General Public License
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* FMS is free software: you can redistribute it and/or modify it under
7 !* the terms of the GNU Lesser General Public License as published by
8 !* the Free Software Foundation, either version 3 of the License, or (at
9 !* your option) any later version.
11 !* FMS is distributed in the hope that it will be useful, but WITHOUT
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13 !* FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
14 !* for more details.
16 !* You should have received a copy of the GNU Lesser General Public
17 !* License along with FMS.  If not, see <http://www.gnu.org/licenses/>.
18 !***********************************************************************
20 ! This program runs only one of many possible tests with each execution.
21 ! Each test ends with an intentional fatal error.
22 ! diag_manager_mod is not a stateless module, and there are situations
23 ! where a fatal error leaves the module in a state that does not allow
24 ! it to function properly if used again. Therefore, the program must
25 ! be terminated after each intentional fatal error.
27 ! Each test is dependent on the diag_table, and different diag_tables
28 ! exist for each test. Depending on the test, an intentional fatal error
29 ! may be triggered upon the call to diag_manager_init, register_diag_field or send_data.
30 ! Because of this, the calls to all of those routines differ depending on the test.
32 ! The diag_table for each test is included below.
34 !--------------------------------------------------------------------------------------------------
35 ! diag_table for test 1
37 ! test_diag_manager
38 ! 1 3 1 0 0 0
39 ! #output files
40 !  "diag_test",  1, "days", 1, "days", "time",
41 ! #output variables
42 !  "test_diag_manager_mod", "dat1", "dat1", "diag_test",  "all", .false., "none", 2,
43 !--------------------------------------------------------------------------------------------------
44 ! diag_table for test 2
46 ! test_diag_manager
47 ! 1 3 1 0 0 0
48 ! #output files
49 !  "diag_test",  1, "days", 1, "days", "time",
50 ! #output variables
51 !  "test_diag_manager_mod", "dat1", "dat1", "diag_test",  "all", .false., "none", 2,
52 !--------------------------------------------------------------------------------------------------
53 ! diag_table for test 3
55 ! test_diag_manager
56 ! 1 3 1 0 0 0
57 ! #output files
58 !  "diag_test",  1, "days", 1, "days", "time",
59 ! #output variables
60 !  "test_diag_manager_mod", "dat1", "dat1", "diag_test",  "all", .false., "none", 2,
61 !--------------------------------------------------------------------------------------------------
62 ! diag_table for test 4
64 ! test_diag_manager
65 ! 1 3 1 0 0 0
66 ! #output files
67 !  "diag_test",  1, "days", 1, "days", "time",
68 !  "diag_test2", 1, "days", 1, "days", "time",
69 ! #output variables
70 !  "test_diag_manager_mod", "dat2", "dat2", "diag_test",  "all", .false., "none", 2,
71 !  "test_mod",              "dat2", "dat2", "diag_test2", "all", .false., "none", 2,
72 !--------------------------------------------------------------------------------------------------
73 ! diag_table for test 5
75 ! test_diag_manager
76 ! 1 3 1 0 0 0
77 ! #output files
78 !  "diag_test",  1, "days", 1, "days", "time",
79 !  "diag_test2", 1, "days", 1, "days", "time",
80 ! #output variables
81 !  "test_diag_manager_mod", "dat2", "dat2", "diag_test",  "all", .false., "none", 2,
82 !  "test_mod",              "dat2", "dat2", "diag_test2", "all", .false., "none", 2,
83 !--------------------------------------------------------------------------------------------------
84 ! diag_table for test 6
86 ! test_diag_manager
87 ! 1 3 1 0 0 0
88 ! #output files
89 !  "diag_test",  1, "days", 1, "days", "time",
90 !  "diag_test2", 1, "days", 1, "days", "time",
91 ! #output variables
92 !  "test_diag_manager_mod", "dat2", "dat2", "diag_test",  "all", .false., "none", 2,
93 !  "test_mod",              "dat2", "dat2", "diag_test2", "all", .false., "none", 2,
94 !--------------------------------------------------------------------------------------------------
95 ! diag_table for test 7
97 ! test_diag_manager
98 ! 1 3 1 0 0 0
99 ! #output files
100 !  "diag_test",  1, "days", 1, "days", "time",
101 ! #output variables
102 !  "test_diag_manager_mod", "dat1", "dat1", "diag_test",  "all", .false., "none", 2,
103 !--------------------------------------------------------------------------------------------------
104 ! diag_table for test 8
106 ! test_diag_manager
107 ! 1 3 1 0 0 0
108 ! #output files
109 !  "diag_test",  1, "days", 1, "days", "time",
110 !  "diag_test2", 1, "days", 1, "days", "time",
111 ! #output variables
112 !  "test_diag_manager_mod", "dat2", "dat2", "diag_test",  "all", .false., "none", 2,
113 !  "test_mod",              "dat2", "dat2", "diag_test2", "all", .false., "none", 2,
114 !--------------------------------------------------------------------------------------------------
115 ! diag_table for test 9
117 ! test_diag_manager
118 ! 1 3 1 0 0 0
119 ! #output files
120 !  "diag_test",  1, "days", 1, "days", "time",
121 ! #output variables
122 !  "test_diag_manager_mod", "bk",   "bk",   "diag_test",  "all", .false., "none", 2,
123 !--------------------------------------------------------------------------------------------------
124 ! diag_table for test 10
126 ! test_diag_manager
127 ! 1 3 1 0 0 0
128 ! #output files
129 !  "diag_test",  1, "days", 1, "days", "time",
130 ! #output variables
131 !  "test_diag_manager_mod", "bk",   "bk",   "diag_test",  "all", .false., "none", 2,
132 !--------------------------------------------------------------------------------------------------
133 ! diag_table for test 11
135 ! test_diag_manager
136 ! 1 3 1 0 0 0
137 ! #output files
138 !  "diag_test",  1, "days", 1, "days", "time",
139 ! #output variables
140 !  "test_diag_manager_mod", "dat2", "dat2", "diag_test",  "all", .false., "none", 2,
141 !--------------------------------------------------------------------------------------------------
142 ! diag_table for test 12
144 ! test_diag_manager
145 ! 1 3 1 0 0 0
146 ! #output files
147 !  "diag_test",  1, "days", 1, "days", "time",
148 ! #output variables
149 !  "test_diag_manager_mod", "dat2", "dat2", "diag_test",  "all", .false., "none", 2,
150 ! # Test of the error check that duplicate field names do not appear in same file,
151 !  "test_mod",              "dat2", "dat2", "diag_test",  "all", .false., "none", 2,
152 !--------------------------------------------------------------------------------------------------
153 ! diag_table for test 13
155 ! test_diag_manager
156 ! 1 3 1 0 0 0
157 ! #output files
158 !  "diag_test",  1, "days", 1, "days", "time",
159 !  "diag_test2", 1, "months", 1, "days", "time",
160 ! #output variables
161 !  "test_diag_manager_mod", "dat2", "dat2", "diag_test",  "all", .false., "none", 2,
162 ! # Test of WARNING message that no data is written when run length is less than output interval
163 !  "test_mod",              "dat2", "dat2", "diag_test2", "all", .false., "none", 2,
164 !--------------------------------------------------------------------------------------------------
165 ! diag_table for test 14
167 ! test_diag_manager
168 ! 1990 1 29 0 0 0
169 ! #output files
170 !  "diag_test2", 1, "months", 1, "days", "time",
171 ! #output variables
172 ! # Test of check for invalid date. (Jan 29 1990 + one month = Feb 29 1990)
173 !  "test_mod",              "dat2", "dat2", "diag_test2", "all", .false., "none", 2,
174 !--------------------------------------------------------------------------------------------------
175 ! diag_table for test 16
177 ! test_diag_manager
178 ! 1 3 1 0 0 0
179 ! #output files
180 !  "diag_test2", 1, "months", 1, "days", "time",
181 ! #output variables
182 ! # Test for output file name to be modified with appended string
183 !  "test_diag_manager_mod", "dat2", "dat2", "diag_test2", "all", .false., "none", 2,
184 !--------------------------------------------------------------------------------------------------
185 ! diag_table for test 17
187 ! test_diag_manager
188 ! 1 3 1 0 0 0
189 ! #output files
190 !  "diag_test2", 1, "days", 1, "days", "time",
191 ! #output variables
192 !  "test_diag_manager_mod", "dat2", "dat2_rms", "diag_test2", "all", "rms",  "none", 2,
193 !  "test_diag_manager_mod", "dat2", "dat2",     "diag_test2", "all", .true., "none", 2,
194 !--------------------------------------------------------------------------------------------------
195 !--------------------------------------------------------------------------------------------------
196 !> diag_table for test 23 (unstructured grid)
198 !test_diag_manager_23
199 !1990 1 1 0 0 0
200 !#output files
201 !"unstructured_diag_test", 2, "days", 2, "days", "time",
202 !#output variables
203 !"UG_unit_test", "unstructured_real_scalar_field_data", "rsf_diag_1", "unstructured_diag_test", "all", .TRUE., "none",
204 ! 1,
205 !"UG_unit_test", "unstructured_real_1D_field_data", "unstructured_real_1D_field_data", "unstructured_diag_test",
206 !"all", .TRUE., "none", 1,
207 !"UG_unit_test", "unstructured_real_2D_field_data", "unstructured_real_2D_field_data", "unstructured_diag_test",
208 !"all", .TRUE., "none", 1,
209 !"UG_unit_test", "lon", "grid_xt", "unstructured_diag_test", "all", .TRUE., "none", 1,
210 !"UG_unit_test", "lat", "grid_yt", "unstructured_diag_test", "all", .TRUE., "none", 1,
211 !--------------------------------------------------------------------------------------------------
212 PROGRAM test
213   ! This program runs only one of many possible tests with each execution.
214   ! Each test ends with an intentional fatal error.
215   ! diag_manager_mod is not a stateless module, and there are situations
216   ! where a fatal error leaves the module in a state that does not allow
217   ! it to function properly if used again. Therefore, the program must
218   ! be terminated after each intentional fatal error.
220   ! Each test is dependent on the diag_table, and different diag_tables
221   ! exist for each test. Depending on the test, an intentional fatal error
222   ! may be triggered upon the call to diag_manager_init, register_diag_field or send_data.
223   ! Because of this, the calls to all of those routines differ depending on the test.
225   USE mpp_mod, ONLY: mpp_pe, mpp_root_pe, mpp_debug, mpp_set_stack_size
226   USE mpp_domains_mod, ONLY: domain2d, mpp_define_domains, mpp_get_compute_domain
227   USE mpp_domains_mod, ONLY: mpp_define_io_domain, mpp_define_layout
228   USE mpp_domains_mod, ONLY: mpp_domains_init, mpp_domains_set_stack_size
229   USE fms_mod, ONLY: fms_init, fms_end, mpp_npes, check_nml_error
230   USE fms_mod, ONLY: error_mesg, FATAL, WARNING, NOTE, stdlog, stdout
231   USE mpp_mod, ONLY: input_nml_file
232 #ifdef use_deprecated_io
233   USE fms_io_mod, ONLY: fms_io_init, file_exist, open_file
234   USE fms_io_mod, ONLY: fms_io_exit, set_filename_appendix
235   use mpp_io_mod, only: mpp_io_init
236 #endif
237   USE constants_mod, ONLY: constants_init, PI, RAD_TO_DEG
239   USE time_manager_mod, ONLY: time_type, set_calendar_type, set_date, decrement_date, OPERATOR(+), set_time
240   USE time_manager_mod, ONLY: NOLEAP, JULIAN, GREGORIAN, THIRTY_DAY_MONTHS, OPERATOR(*), assignment(=)
241   use time_manager_mod, ONLY: OPERATOR(+), OPERATOR(-), OPERATOR(/), days_in_month
243   USE diag_manager_mod, ONLY: diag_manager_init, send_data, diag_axis_init, diag_manager_end
244   USE diag_manager_mod, ONLY: register_static_field, register_diag_field, diag_send_complete
245   USE diag_manager_mod, ONLY: diag_manager_set_time_end, diag_field_add_attribute, diag_axis_add_attribute
246   USE diag_manager_mod, ONLY: diag_field_add_cell_measures
247   USE diag_manager_mod, ONLY: get_diag_field_id, DIAG_FIELD_NOT_FOUND
248   USE diag_axis_mod, ONLY: get_axis_num
249   USE platform_mod
251   IMPLICIT NONE
253   TYPE(domain2d) :: Domain1
254   TYPE(domain2d) :: Domain2
256   REAL, ALLOCATABLE, DIMENSION(:) :: lon_global1, lonb_global1
257   REAL, ALLOCATABLE, DIMENSION(:) :: lat_global1, latb_global1
258   REAL, ALLOCATABLE, DIMENSION(:) :: lon_global2, lonb_global2
259   REAL, ALLOCATABLE, DIMENSION(:) :: lat_global2, latb_global2
260   REAL, ALLOCATABLE, DIMENSION(:) :: pfull, bk, phalf
261   REAL, ALLOCATABLE, DIMENSION(:) :: lon1, lat1, lonb1, latb1
262   REAL, ALLOCATABLE, DIMENSION(:) :: lon2, lat2, lonb2, latb2
263   REAL, ALLOCATABLE, DIMENSION(:,:,:) :: dat1, dat1h
264   REAL, ALLOCATABLE, DIMENSION(:,:,:) :: dat2, dat2h
265   REAL, ALLOCATABLE, DIMENSION(:,:) :: dat2_2d
266   REAL :: solar_constant = 1600
267   REAL :: surf_press = 1.e5
268   REAL :: dp
269   INTEGER :: id_phalf, id_pfull, id_bk
270   INTEGER :: id_lon1, id_lonb1, id_latb1, id_lat1, id_dat1
271   INTEGER :: id_lon2, id_lat2, id_dat2, id_dat2_2d, id_sol_con, id_dat2h, id_dat2h_2
272   INTEGER :: id_dat2_got, id_none_got
273   INTEGER :: i, j, k, is1, ie1, js1, je1, nml_unit, ierr, log_unit, out_unit, m
274   INTEGER :: is_in, ie_in, js_in, je_in
275   INTEGER :: is2, ie2, js2, je2, hi=1, hj=1
276   INTEGER :: nlon1, nlat1, nlon2, nlat2
277   INTEGER, DIMENSION(2) :: layout = (/1,1/)
278   INTEGER :: test_number=1
279   INTEGER :: nlon=10, nlat=10, nlev=10
280   INTEGER :: io_layout(2) = (/1,1/)
281   INTEGER :: nstep = 2
282   TYPE(time_type) :: Time, Time_step, Time_end, Time_start, Run_length
283   LOGICAL :: used, test_successful
284   CHARACTER(len=256) :: err_msg
285   INTEGER :: omp_get_num_threads
287   INTEGER :: nyc1, n, jsw, jew, isw, iew
288   INTEGER :: numthreads=1, ny_per_thread, idthread
289   INTEGER :: months=0, days=1, dt_step=1
291   ! Variables needed for test 22
292   INTEGER :: id_nv, id_nv_init
294 !!!!!! Stuff for unstrctured grid
295     integer(kind=i4_kind)              :: nx = 8                               !<Total number of grid points in the
296                                                                                !! x-dimension (longitude?)
297     integer(kind=i4_kind)              :: ny = 8           !<Total number of grid points in the y-dimension (latitude?)
298     integer(kind=i4_kind)              :: nz = 2           !<Total number of grid points in the z-dimension (height)
299     integer(kind=i4_kind)              :: nt = 2                               !<Total number of time grid points.
300     integer(kind=i4_kind)              :: io_tile_factor = 1                   !< The IO tile factor
301     integer(kind=i4_kind)              :: halo = 2                             !<Number of grid points in the halo???
302     integer(kind=i4_kind)              :: ntiles_x = 1     !<Number of tiles in the x-direction
303                                                            !! (A 2D grid of tiles is used in this test.)
304     integer(kind=i4_kind)              :: ntiles_y = 2     !<Number of tiles in the y-direction (A 2D grid of tiles is
305                                                            !! used in this test.)
306     integer(kind=i4_kind)              :: total_num_tiles  !<The total number of tiles for the run (=ntiles_x*ntiles_y)
307     integer(kind=i4_kind)              :: stackmax = 1500000 !<Default size to which the mpp stack will be set.
308     integer(kind=i4_kind)              :: stackmaxd = 500000 !<Default size to which the mpp_domains stack will be set.
309     logical(kind=l4_kind)              :: debug = .false.                      !<Flag to print debugging information.
310     character(len=64)              :: test_file = "test_unstructured_grid" !<Base filename for the unit tests.
311     character(len=64)              :: iospec = '-F cachea'                 !<Something cray related ???
312     integer(kind=i4_kind)              :: pack_size = 1    !<(Number of bits in real(kind=r8_kind))/(Number of bits
313                                                            !! in real)
314     integer(kind=i4_kind)              :: npes             !<Total number of ranks in the current pelist.
315     integer(kind=i4_kind)              :: io_status                            !<Namelist read error code.
316     real(kind=r8_kind)              :: doubledata = 0.0    !<Used to determine pack_size.  This must be kind=r8_kind.
317     real                           :: realdata = 0.0   !<Used to determine pack_size.  Do not specify a kind parameter.
318     integer(kind=i4_kind)              :: funit = 7                            !<File unit.
319     logical(kind=l4_kind)              :: fopened      !<Flag telling if a file is already open.
320     type(time_type)                :: diag_time
322     integer(kind=i4_kind)              :: output_unit=6
323 !!!!!!
327   NAMELIST /test_diag_manager_nml/ layout, test_number, nlon, nlat, nlev, io_layout, numthreads, &
328                                    dt_step, months, days
329   NAMELIST /utest_nml/nx,ny,nz,nt,ntiles_x,ntiles_y,io_tile_factor
330   ! Initialize all id* vars to be -1
331   id_nv = -1
332   id_nv_init = -1
333   id_phalf = -1
334   id_pfull = -1
335   id_bk = -1
336   id_lon1 = -1
337   id_lonb1 = -1
338   id_latb1 = -1
339   id_lat1 = -1
340   id_dat1 = -1
341   id_lon2 = -1
342   id_lat2 = -1
343   id_dat2 = -1
344   id_dat2_2d = -1
345   id_sol_con = -1
346   id_dat2h = -1
347   id_dat2h_2 = -1
348   id_dat2_got = -1
349   id_none_got = -1
351   CALL fms_init
352   log_unit = stdlog()
353   out_unit = stdout()
354   CALL constants_init
355   CALL set_calendar_type(JULIAN)
356   npes = mpp_npes()
357   READ (input_nml_file, NML=test_diag_manager_nml, IOSTAT=ierr)
358   READ (input_nml_file, NML=utest_nml, IOSTAT=i)
359   ! Check the status of reading the diag_manager_nml
360   IF ( check_nml_error(IOSTAT=ierr, NML_NAME='DIAG_MANAGER_NML') < 0 ) THEN
361      IF ( mpp_pe() == mpp_root_pe() ) THEN
362         CALL error_mesg('diag_manager_mod::diag_manager_init', &
363                & 'TEST_DIAG_MANAGER_NML not found in input.nml.  Using defaults.', WARNING)
364      END IF
365   END IF
366   WRITE (log_unit,test_diag_manager_nml)
368 SELECT CASE ( test_number ) ! Closes just before the CONTAINS block.
369   ! If the test_number == 23, then call the unstrcutured grid unit test and skip everything else.
370   CASE ( 23 )
371    !Initialize the mpp_domains module
372     if (debug) then
373         call mpp_domains_init(MPP_DEBUG)
374     else
375         call mpp_domains_init()
376     endif
378    !Initialize the mpp_io module.
379 #ifdef use_deprecated_io
380     if (debug) then
381         call mpp_io_init(MPP_DEBUG)
382     else
383         call mpp_io_init()
384     endif
385 #endif
387    !Set the mpp and mpp_domains stack sizes.
388     call mpp_set_stack_size(stackmax)
389     call mpp_domains_set_stack_size(stackmaxd)
391    !Write out test configuration parameters.
392     if (mpp_pe() .eq. mpp_root_pe()) then
393         write(output_unit,*)
394         write(output_unit,*) "Performing unstructured_io unit test with:"
395         write(output_unit,*) "Total number of ranks:                          ", &
396                              npes
397         write(output_unit,*) "Total number of grid points in the x-dimension: ", &
398                              nx
399         write(output_unit,*) "Total number of grid points in the y-dimension: ", &
400                              ny
401         write(output_unit,*) "Total number of grid points in the z-dimension: ", &
402                              nz
403         write(output_unit,*) "Total number of grid points in the t-dimension: ", &
404                              nt
405         write(output_unit,*) "Halo width (# of grid points):                  ", &
406                              halo
407         write(output_unit,*) "Using Unstructured domaintypes and calls..."
408     endif
410    !Add a suffix to the test file.
411     write(test_file,'(a,i3.3)') trim(test_file),npes
413    !Initialize the diag manager module.
414     call diag_manager_init()
416    !Set the diag_time variable to be 01/01/1990 at 00:00:00 (midnight).
417     call set_calendar_type(JULIAN)
418     time = set_date(1990,1,1,0,0,0)
419    CALL unstruct_test (nx, ny, nz, npes, ntiles_x, 1, time,io_tile_factor)
421    ! If the test_number == 12, check for the correct error and skip everything else.
422    CASE ( 12 )
423      CALL diag_manager_init(err_msg=err_msg)
424      IF ( err_msg /= '' ) THEN
425         WRITE (out_unit,'(a)') 'test12 successful: err_msg='//TRIM(err_msg)
426         CALL error_mesg('test_diag_manager','test12 successful.',NOTE)
427      ELSE
428         WRITE (out_unit,'(a)') 'test12 fails'
429         CALL error_mesg('test_diag_manager','test12 fails',FATAL)
430      END IF
432   ! If the test number is not 12 or 23, run all other tests.
433   CASE DEFAULT ! Contains all remaining code up to CONTAINS block.
434      CALL diag_manager_init
436   IF ( layout(1)*layout(2) .NE. mpp_npes() ) THEN
437      CALL mpp_define_layout((/1,nlon,1,nlat/), mpp_npes(), layout )
438   END IF
440   nlon1 = nlon
441   nlat1 = nlat
442   nlon2 = nlon * 2
443   nlat2 = nlat * 2
445   CALL mpp_define_domains((/1,nlon1,1,nlat1/), layout, Domain1, name='test_diag_manager')
446   CALL mpp_get_compute_domain(Domain1, is1, ie1, js1, je1)
447   ALLOCATE(lon_global1(nlon1), lonb_global1(nlon1+1))
448   ALLOCATE(lat_global1(nlat1), latb_global1(nlat1+1))
449   ALLOCATE(lon_global2(nlon2), lonb_global2(nlon2+1))
450   ALLOCATE(lat_global2(nlat2), latb_global2(nlat2+1))
451   ALLOCATE(pfull(nlev), bk(nlev), phalf(nlev+1))
453   ALLOCATE(lon1(is1:ie1), lat1(js1:je1), lonb1(is1:ie1+1), latb1(js1:je1+1))
454   CALL compute_grid(nlon1, nlat1, is1, ie1, js1, je1, lon_global1, lat_global1, lonb_global1, latb_global1, lon1, &
455                   & lat1, lonb1, latb1)
456   CALL mpp_define_domains((/1,nlon2,1,nlat2/), layout, Domain2, name='test_diag_manager')
457   CALL mpp_get_compute_domain(Domain2, is2, ie2, js2, je2)
458   CALL mpp_define_io_domain(Domain1, io_layout)
459   CALL mpp_define_io_domain(Domain2, io_layout)
461   ALLOCATE(lon2(is2:ie2), lat2(js2:je2), lonb2(is2:ie2+1), latb2(js2:je2+1))
462   CALL compute_grid(nlon2, nlat2, is2, ie2, js2, je2, lon_global2, lat_global2, lonb_global2, latb_global2, lon2, &
463                   & lat2, lonb2, latb2)
464   dp = surf_press/nlev
465   DO k=1, nlev+1
466      phalf(k) = dp*(k-1)
467   END DO
468   DO k=1, nlev
469      pfull(k) = .5*(phalf(k) + phalf(k+1))
470      bk(k) = pfull(k)/surf_press
471   END DO
473   ALLOCATE(dat1(is1:ie1,js1:je1,nlev))
474   ALLOCATE(dat1h(is1-hi:ie1+hi,js1-hj:je1+hj,nlev))
475   dat1h = 0.
476   dat1 = 0.
477   DO j=js1, je1
478      DO i=is1, ie1
479         dat1(i,j,1) = SIN(lon1(i))*COS(lat1(j))
480      END DO
481   END DO
482   dat1h(is1:ie1,js1:je1,1) = dat1(:,:,1)
483   dat1(:,:,2) = -dat1(:,:,1)
484   dat1h(:,:,2) = -dat1h(:,:,1)
486   ALLOCATE(dat2(is2:ie2,js2:je2,nlev))
487   ALLOCATE(dat2_2d(is2:ie2,js2:je2))
488   ALLOCATE(dat2h(is2-hi:ie2+hi,js2-hj:je2+hj,nlev))
489   dat2h = 0.
490   dat2 = 0.
491   DO j=js2, je2
492      DO i=is2, ie2
493         dat2(i,j,1) = SIN(lon2(i))*COS(lat2(j))
494      END DO
495   END DO
496   dat2h(is2:ie2,js2:je2,1) = dat2(:,:,1)
497   dat2(:,:,2) = -dat2(:,:,1)
498   dat2h(:,:,2) = -dat2h(:,:,1)
499   dat2_2d = dat2(:,:,1)
501   id_lonb1 = diag_axis_init('lonb1', RAD_TO_DEG*lonb_global1, 'degrees_E', 'x', long_name='longitude edges', &
502                           & Domain2=Domain1)
503   id_latb1 = diag_axis_init('latb1', RAD_TO_DEG*latb_global1, 'degrees_N', 'y', long_name='latitude edges',  &
504                           & Domain2=Domain1)
506   id_lon1  = diag_axis_init('lon1',  RAD_TO_DEG*lon_global1, 'degrees_E','x',long_name='longitude',Domain2=Domain1, &
507                           & edges=id_lonb1)
508   id_lat1  = diag_axis_init('lat1',  RAD_TO_DEG*lat_global1, 'degrees_N','y',long_name='latitude', Domain2=Domain1, &
509                           & edges=id_latb1)
511   id_phalf= diag_axis_init('phalf', phalf, 'Pa', 'z', long_name='half pressure level', direction=-1)
512   id_pfull= diag_axis_init('pfull', pfull, 'Pa', 'z', long_name='full pressure level', direction=-1, edges=id_phalf)
514   id_lon2 = diag_axis_init('lon2',  RAD_TO_DEG*lon_global2,  'degrees_E', 'x', long_name='longitude', Domain2=Domain2)
515   id_lat2 = diag_axis_init('lat2',  RAD_TO_DEG*lat_global2,  'degrees_N', 'y', long_name='latitude',  Domain2=Domain2)
517   IF ( test_number == 22 ) THEN
518      ! Can we get the 'nv' axis ID?
519      id_nv = get_axis_num('nv', 'nv')
520      IF ( id_nv .GT. 0 ) THEN
521         write (out_unit,'(a)') 'test22.1 Passes: id_nv has a positive value'
522      ELSE
523         write (out_unit,'(a)') 'test22.1 Failed: id_nv does not have a positive value'
524      END IF
526      ! Can I call diag_axis_init on 'nv' again, and get the same ID back?
527      id_nv_init = diag_axis_init( 'nv',(/1.,2./),'none','N','vertex number', set_name='nv')
528      IF ( id_nv_init .EQ. id_nv ) THEN
529         write (out_unit,'(a)') 'test22.2 Passes: Can call diag_axis_init on "nv" and get same ID'
530      ELSE
531         write (out_unit,'(a)') 'test22.2 Failed: Cannot call diag_axis_init on "nv" and get same ID'
532      END IF
533   END IF
535   IF ( test_number == 21 ) THEN
536      ! Testing addition of axis attributes
537      CALL diag_axis_add_attribute(id_lon1, 'real_att', 2.3)
538      CALL diag_axis_add_attribute(id_lat1, 'int_att', (/ 2, 3 /))
539      CALL diag_axis_add_attribute(id_pfull, 'char_att', 'Some string')
540   END IF
542   IF ( test_number == 14 ) THEN
543      Time = set_date(1990,1,29,0,0,0)
544   ELSE
545      Time = set_date(1990,1,1,0,0,0)
546   END IF
548   IF ( test_number == 16 ) THEN
549      ! Test 16 tests the filename appendix
550 #ifdef use_deprecated_io
551      CALL set_filename_appendix('g01')
552 #endif
553   END IF
554   id_dat1 = register_diag_field('test_diag_manager_mod', 'dat1', (/id_lon1,id_lat1,id_pfull/), Time, 'sample data','K')
555   IF ( test_number == 18 ) THEN
556      CALL diag_field_add_attribute(id_dat1, 'real_att', 2.3)
557      CALL diag_field_add_attribute(id_dat1, 'cell_methods', 'area: mean')
558      CALL diag_field_add_attribute(id_dat1, 'cell_methods', 'lon: mean')
559   END IF
560   IF ( test_number == 18 .OR. test_number == 19 ) THEN
561      id_dat2 = register_diag_field('test_diag_manager_mod', 'dat2', (/id_lon1,id_lat1,id_pfull/), Time, &
562                                  & 'sample data', 'K')
563      CALL diag_field_add_attribute(id_dat2, 'interp_method', 'none')
564      CALL diag_field_add_attribute(id_dat2, 'int_att', (/ 1, 2 /) )
565   ELSE
566      id_dat2 = register_diag_field('test_diag_manager_mod', 'dat2', (/id_lon2,id_lat2,id_pfull/), Time, &
567                                  & 'sample data', 'K')
568   END IF
569   id_sol_con = register_diag_field ('test_diag_manager_mod', 'solar_constant', Time, &
570                   'solar constant', 'watts/m2')
572   IF ( test_number == 20 ) THEN
573      id_dat2_got = get_diag_field_id('test_diag_manager_mod', 'dat2')
574      IF ( id_dat2_got == id_dat2 ) THEN
575         WRITE (out_unit,'(a)') 'test20.1 Passes, id_dat2.EQ.id_dat2_got'
576      ELSE
577         WRITE (out_unit,'(a)') 'test20.1 Failed, id_dat2.NE.id_dat2_got'
578      END IF
580      id_none_got = get_diag_field_id('no_mod', 'no_var')
581      IF ( id_none_got == DIAG_FIELD_NOT_FOUND ) THEN
582         write (out_unit,'(a)') 'test20.2 Passes, id_none_got.EQ.DIAG_FIELD_NOT_FOUND'
583      ELSE
584         write (out_unit,'(a)') 'test20.2 Failed, id_none_got.NE.DIAG_FIELD_NOT_FOUND'
585      END IF
586   END IF
588   IF ( dt_step == 0 ) CALL error_mesg ('test_diag_manager',&
589        & 'dt_step is not set', FATAL)
591   Time_step = set_time(dt_step,0)
592   Time_start = Time
593   Time_end = Time
594   DO m = 1,months
595      Time_end = Time_end + set_time(0,days_in_month(Time_end))
596   END DO
597   Time_end   = Time_end + set_time(0, days)
598   Run_length = Time_end - Time_start
599   nstep = Run_length / Time_step
601   IF ( test_number == 18 ) THEN
602      id_dat2h = register_diag_field('test_mod', 'dat2h', (/id_lon1,id_lat1,id_pfull/), Time, 'sample data', 'K',&
603           & volume=id_dat1, area=id_dat2, realm='myRealm', err_msg=err_msg)
604      IF ( err_msg /= '' .OR. id_dat2h <= 0 ) THEN
605         CALL error_mesg ('test_diag_manager',&
606              & 'Unexpected error registering dat2h '//err_msg, FATAL)
607      END IF
608      id_dat2h_2 = register_diag_field('test_mod', 'dat2h_2', (/id_lon1,id_lat1,id_pfull/), Time, 'sample data', 'K',&
609           & err_msg=err_msg)
610      CALL diag_field_add_cell_measures(id_dat2h_2, area=id_dat2, volume=id_dat1)
611   ELSE IF ( test_number == 19 ) THEN
612      id_dat2h = register_diag_field('test_mod', 'dat2h', (/id_lon1,id_lat1,id_pfull/), Time, 'sample data', 'K',&
613           & volume=id_dat1, area=id_dat1, err_msg=err_msg)
614      IF ( err_msg /= '' .OR. id_dat2h <= 0 ) THEN
615         CALL error_mesg ('test_diag_manager',&
616              & 'Expected error registering dat2h '//err_msg, NOTE)
617      END IF
618   END IF
620   IF(test_number == 16 .OR. test_number == 17 .OR. test_number == 18 .OR. test_number == 21 .OR. test_number == 22)THEN
621      is_in = 1
622      js_in = 1
623      ie_in = nlon
624      je_in = nlat
626      IF ( id_dat1 > 0 ) used = send_data(id_dat1, dat1, Time, err_msg=err_msg)
627      IF ( id_dat2 > 0 ) used = send_data(id_dat2, dat1, Time, err_msg=err_msg)
628      IF ( id_dat2h > 0 ) used = send_data(id_dat2h, dat2h, Time, is_in=is_in, js_in=js_in, ie_in=ie_in, je_in=je_in, &
629                                         & err_msg=err_msg)
630      IF ( id_dat2h_2 > 0 ) used = send_data(id_dat2h_2, dat2h, Time, is_in=is_in, js_in=js_in, ie_in=ie_in, &
631                                           & je_in=je_in, err_msg=err_msg)
632      Time = Time + set_time(0,1)
633      IF ( id_dat1 > 0 ) used = send_data(id_dat1, dat1, Time, err_msg=err_msg)
634      IF ( id_dat2 > 0 ) used = send_data(id_dat2, dat1, Time, err_msg=err_msg)
635      IF ( id_dat2h > 0 ) used = send_data(id_dat2h, dat2h, Time, is_in=is_in, js_in=js_in, ie_in=ie_in, je_in=je_in, &
636                                         & err_msg=err_msg)
637      IF ( id_dat2h_2 > 0 ) used = send_data(id_dat2h_2, dat2h, Time, is_in=is_in, js_in=js_in, ie_in=ie_in, &
638                                           & je_in=je_in, err_msg=err_msg)
639   END IF
641   !-- The following is used to test openMP
642   IF ( test_number == 15 ) THEN
643 !$      call omp_set_num_threads(numthreads)
644       nyc1 = je1 - js1 + 1
645       IF (MOD(nyc1, numthreads ) /= 0) THEN
646          CALL error_mesg ('test_diag_manager',&
647               & 'The number of OpenMP threads must be an integral multiple &
648               &of the number of rows in the compute domain', FATAL)
649      END IF
650      ny_per_thread = nyc1/numthreads
652      dat1 = 1
653      CALL diag_manager_set_time_end(Time_end)
654      DO n = 1, nstep
656         Time = Time + Time_step
657         !$OMP parallel do default(shared) private(isw, iew, jsw, jew )
659         DO jsw = js1, je1, ny_per_thread
660            jew = jsw + ny_per_thread -1
661            isw = is1
662            iew = ie1
663            if(id_dat1>0) used = send_data(id_dat1, dat1(isw:iew, jsw:jew,:), Time, &
664                                 is_in=isw-is1+1, js_in=jsw-js1+1,err_msg=err_msg)
665            if(id_sol_con>0) used = send_data(id_sol_con, solar_constant, Time, err_msg=err_msg )
666         END DO
667         !$OMP END parallel do
668         !CALL diag_send_complete(Time_step)
669      END DO
670  END IF
673   IF ( test_number == 14 ) THEN
674      id_dat2_2d = register_diag_field('test_mod', 'dat2', (/id_lon2,id_lat2/), Time, 'sample data','K',err_msg=err_msg)
675      IF ( err_msg /= '' ) THEN
676         WRITE (out_unit,'(a)') 'test14 successful. err_msg='//TRIM(err_msg)
677      ELSE
678         WRITE (out_unit,'(a)') 'test14 fails.'
679      END IF
680   ELSE
681      id_dat2_2d = register_diag_field('test_mod', 'dat2', (/id_lon2,id_lat2/), Time, 'sample data', 'K')
682   END IF
684   id_bk = register_static_field('test_diag_manager_mod', 'bk', (/id_pfull/), 'half level sigma', 'none')
686   IF ( test_number == 13 ) THEN
687      IF ( id_dat2_2d > 0 ) used=send_data(id_dat2_2d, dat2(:,:,1), Time, err_msg=err_msg)
688      IF ( err_msg == '' ) THEN
689         WRITE (out_unit,'(a)') &
690               & 'test13: successful if a WARNING message appears that refers to output interval greater than runlength'
691      ELSE
692         WRITE (out_unit,'(a)') 'test13 fails: err_msg='//TRIM(err_msg)
693      END IF
694   END IF
696   ! Note: test12 involves diag_manager_init, it does not require a call to send_data.
697   !       See call to diag_manager_init above.
699   IF ( test_number == 11 ) THEN
700      is_in = 1+hi
701      js_in = 1+hj
702      ie_in = ie2-is2+1+hi
703      je_in = je2-js2+1+hj
705      IF ( id_dat2 > 0 ) used=send_data(id_dat2, dat2h, Time, is_in=is_in, js_in=js_in, ie_in=ie_in, je_in=je_in, &
706                                      & err_msg=err_msg)
707      IF ( err_msg == '' ) THEN
708         WRITE (out_unit,'(a)') 'test11.1 successful.'
709      ELSE
710         WRITE (out_unit,'(a)') 'test11.1 fails. err_msg='//TRIM(err_msg)
711      END IF
713      ! intentional_error: je_in is missing
714      IF ( id_dat2 > 0 ) used=send_data(id_dat2, dat2h, Time, is_in=is_in, js_in=js_in, ie_in=ie_in, err_msg=err_msg)
715      IF ( err_msg == '' ) THEN
716         WRITE (out_unit,'(a)') 'test11.2 fails.'
717      ELSE
718         WRITE (out_unit,'(a)') 'test11.2 successful. err_msg='//TRIM(err_msg)
719      END IF
720   END IF
722   IF ( test_number == 10 ) THEN
723      !  1 window, no halos, static, 1 dimension, global data.
725      IF ( id_bk > 0 ) used = send_data(id_bk, bk, err_msg=err_msg)
726      IF ( err_msg == '' ) THEN
727         WRITE (out_unit,'(a)') 'test10.1 successful.'
728      ELSE
729         WRITE (out_unit,'(a)') 'test10.1 fails: err_msg='//TRIM(err_msg)
730      END IF
732      !  intentional_error: data array too large.
733      IF ( id_bk > 0 ) used = send_data(id_bk, phalf, err_msg=err_msg)
734      IF ( err_msg == '' ) THEN
735         WRITE(out_unit,'(a)') 'test10.2 fails.'
736      ELSE
737         WRITE (out_unit,'(a)') 'test10.2 successful: err_msg='//TRIM(err_msg)
738      END IF
739   END IF
741   IF ( test_number == 9 ) THEN
742      !  1 window, no halos, static, 1 dimension, global data
743      IF ( id_bk > 0 ) used = send_data(id_bk, bk, err_msg=err_msg)
744      IF ( err_msg == '' ) THEN
745         WRITE (out_unit,'(a)') 'test9.1 successful.'
746      ELSE
747         WRITE (out_unit,'(a)') 'test9.1 fails: err_msg='//TRIM(err_msg)
748      END IF
750      !  intentional_error: data array too small
751      IF ( id_bk > 0 ) used = send_data(id_bk, bk(1:nlev-1), err_msg=err_msg) ! intentional_error
752      IF ( err_msg == '' ) THEN
753         WRITE (out_unit,'(a)') 'test9.2 fails.'
754      ELSE
755         WRITE (out_unit,'(a)') 'test9.2 successful: err_msg='//TRIM(err_msg)
756      END IF
757   END IF
759   IF ( test_number == 8 ) THEN
760      !  1 window with halos
761      is_in = 1+hi
762      js_in = 1+hj
764      ie_in = ie2-is2+1+hi
765      je_in = je2-js2+1+hj
766      IF ( id_dat2 > 0 ) used=send_data(id_dat2, dat2h, Time, is_in=is_in, js_in=js_in,&
767           & ks_in=1, ie_in=ie_in, je_in=je_in, ke_in=nlev, err_msg=err_msg)
768      IF ( err_msg == '' ) THEN
769         WRITE (out_unit,'(a)') 'test8.1 successful.'
770      ELSE
771         WRITE (out_unit,'(a)') 'test8.1 fails: err_msg='//TRIM(err_msg)
772      END IF
774      !  intentional_error: data array too small in both x and y directions
775      !  Error check is done on second call to send_data. Change in value of Time triggers the check.
776      Time = Time + set_time(0,1)
777      ie_in = ie1-is1+1+hi
778      je_in = je1-js1+1+hj
779      IF ( id_dat2 > 0 ) used=send_data(id_dat2, dat1h, Time, is_in=is_in, js_in=js_in,&
780           & ks_in=1, ie_in=ie_in, je_in=je_in, ke_in=nlev, err_msg=err_msg)
781      Time = Time + set_time(0,1)
782      IF ( id_dat2 > 0 ) used=send_data(id_dat2, dat1h, Time, is_in=is_in, js_in=js_in, &
783           & ks_in=1, ie_in=ie_in, je_in=je_in, ke_in=nlev, err_msg=err_msg)
784      IF ( err_msg == '' ) THEN
785         WRITE (out_unit,'(a)') 'test8.2 fails.'
786      ELSE
787         WRITE (out_unit,'(a)') 'test8.2 successful: err_msg='//TRIM(err_msg)
788      END IF
789   END IF
791   IF ( test_number == 7 ) THEN
792      !  1 window with halos
793      is_in = 1+hi
794      js_in = 1+hj
796      ie_in = ie1-is1+1+hi
797      je_in = je1-js1+1+hj
798      IF ( id_dat1 > 0 ) used=send_data(id_dat1, dat1h, Time, is_in=is_in, js_in=js_in,&
799           & ks_in=1, ie_in=ie_in, je_in=je_in, ke_in=nlev, err_msg=err_msg)
800      IF ( err_msg == '' ) THEN
801         WRITE (out_unit,'(a)') 'test7.1 successful.'
802      ELSE
803         WRITE (out_unit,'(a)') 'test7.1 fails: err_msg='//TRIM(err_msg)
804      END IF
806      !  intentional_error: data array too large in both x and y directions
807      ie_in = ie2-is2+1+hi
808      je_in = je2-js2+1+hj
809      IF ( id_dat1 > 0 ) used=send_data(id_dat1, dat2h, Time, is_in=is_in, js_in=js_in,&
810           & ks_in=1, ie_in=ie_in, je_in=je_in, ke_in=nlev, err_msg=err_msg)
811      IF ( err_msg == '' ) THEN
812         WRITE (out_unit,'(a)') 'test7.2 fails.'
813      ELSE
814         WRITE (out_unit,'(a)') 'test7.2 successful: err_msg='//TRIM(err_msg)
815      END IF
816   END IF
818   IF ( test_number == 6 ) THEN
819      !  multiple windows, no halos
820      !  No error messages should appear at any point within either do loop for test6.1
821      test_successful = .TRUE.
822      DO i=is2, ie2
823         IF ( id_dat2 > 0 ) used = send_data(id_dat2, dat2(i:i,:,:), Time, i-is2+1, 1, err_msg=err_msg)
824         IF ( err_msg /= '' ) THEN
825            WRITE (out_unit,'(a)') 'test6.1 fails: err_msg='//TRIM(err_msg)
826            test_successful = .FALSE.
827         END IF
828      END DO
829      Time = Time + set_time(0,1)
830      DO i=is2, ie2
831         IF ( id_dat2 > 0 ) used = send_data(id_dat2, dat2(i:i,:,:), Time, i-is2+1, 1, err_msg=err_msg)
832         IF ( err_msg /= '' ) THEN
833            WRITE (out_unit,'(a)') 'test6.1 fails: err_msg='//TRIM(err_msg)
834            test_successful = .FALSE.
835         END IF
836      END DO
837      IF ( test_successful ) THEN
838         WRITE (out_unit,'(a)') 'test6.1 successful.'
839      ELSE
840         WRITE (out_unit,'(a)') 'test6.1 fails.'
841      END IF
843      !  intentional_error: data array too small in y direction
844      !  Error check is done on second call to send_data. Change in value of Time triggers the check.
845      Time = Time + set_time(0,1)
846      DO i=is2, ie2
847         IF ( id_dat2 > 0 ) used = send_data(id_dat2, dat2(i:i,js2:je2-1,:), Time, i-is2+1, 1)
848      END DO
849      Time = Time + set_time(0,1)
850      DO i=is2, ie2
851         IF ( id_dat2 > 0 ) used = send_data(id_dat2, dat2(i:i,js2:je2-1,:), Time, i-is2+1, 1, err_msg=err_msg)
852         IF ( err_msg /= '' ) EXIT ! exit immediately after error is detected. No need to continue.
853      END DO
854      IF ( err_msg == '' ) THEN
855         WRITE (out_unit,'(a)') 'test6.2 fails.'
856      ELSE
857         WRITE (out_unit,'(a)') 'test6.2 successful: err_msg='//TRIM(err_msg)
858      END IF
859   END IF
861   IF ( test_number == 5 ) THEN
862      !  multiple windows, no halos
863      !  No error messages should appear at any point within either do loop for test5.1
864      test_successful = .TRUE.
865      DO j=js2, je2
866         IF ( id_dat2 > 0 ) used = send_data(id_dat2, dat2(:,j:j,:), Time, 1, j-js2+1, err_msg=err_msg)
867         IF ( err_msg /= '' ) THEN
868            WRITE (out_unit,'(a)') 'test5.1 fails: err_msg='//TRIM(err_msg)
869            test_successful = .FALSE.
870         END IF
871      END DO
872      Time = Time + set_time(0,1)
873      DO j=js2, je2
874         IF ( id_dat2 > 0 ) used = send_data(id_dat2, dat2(:,j:j,:), Time, 1, j-js2+1, err_msg=err_msg)
875         IF ( err_msg /= '' ) THEN
876            WRITE (out_unit,'(a)') 'test5.1 fails: err_msg='//TRIM(err_msg)
877            test_successful = .FALSE.
878         END IF
879      END DO
880      IF ( test_successful ) THEN
881         WRITE (out_unit,'(a)') 'test5.1 successful.'
882      ELSE
883         WRITE (out_unit,'(a)') 'test5.1 fails.'
884      END IF
886      !  intentional_error: data array too small in x direction.
887      !  Error check is done on second call to send_data. Change in value of Time triggers the check.
888      Time = Time + set_time(0,1)
889      DO j=js2, je2
890         IF ( id_dat2 > 0 ) used = send_data(id_dat2, dat2(is2:ie2-1,j:j,:), Time, 1, j-js2+1)
891      END DO
892      Time = Time + set_time(0,1)
893      DO j=js2, je2
894         IF ( id_dat2 > 0 ) used = send_data(id_dat2, dat2(is2:ie2-1,j:j,:), Time, 1, j-js2+1, err_msg=err_msg)
895         IF ( err_msg /= '' ) EXIT ! exit immediately after error is detected. No need to continue.
896      END DO
897      IF ( err_msg == '' ) THEN
898         WRITE (out_unit,'(a)') 'test5.2 fails.'
899      ELSE
900         WRITE (out_unit,'(a)') 'test5.2 successful: err_msg='//TRIM(err_msg)
901      END IF
902   END IF
904   IF ( test_number == 4 ) THEN
905      !  1 window, no halos
906      IF ( id_dat2 > 0 ) used = send_data(id_dat2, dat2, Time, err_msg=err_msg)
907      Time = Time + set_time(0,1)
908      IF ( id_dat2 > 0 ) used = send_data(id_dat2, dat2, Time, err_msg=err_msg)
909      IF ( err_msg == '' ) THEN
910         WRITE (out_unit,'(a)') 'test4.1 successful.'
911      ELSE
912         WRITE (out_unit,'(a)') 'test4.1 fails: err_msg='//TRIM(err_msg)
913      END IF
915      !  intentional_error: data array too small in both x and y directions
916      !  Error check is done on second call to send_data. Change in value of Time triggers the check.
917      Time = Time + set_time(0,1)
918      IF ( id_dat2 > 0 ) used = send_data(id_dat2, dat1, Time, err_msg=err_msg)
919      Time = Time + set_time(0,1)
920      IF ( id_dat2 > 0 ) used = send_data(id_dat2, dat1, Time, err_msg=err_msg)
921      IF ( err_msg == '' ) THEN
922         WRITE (out_unit,'(a)') 'test4.2 fails.'
923      ELSE
924         WRITE (out_unit,'(a)') 'test4.2 successful: err_msg='//TRIM(err_msg)
925      END IF
926   END IF
928   IF ( test_number == 3 ) THEN
929      !  multiple windows, no halos
930      !  No error messages should appear at any point within do loop for test3.1
931      test_successful = .TRUE.
932      DO i=is1, ie1
933         IF ( id_dat1 > 0 ) used = send_data(id_dat1, dat1(i:i,:,:), Time, i-is1+1, 1, err_msg=err_msg)
934         IF ( err_msg /= '' ) THEN
935            WRITE (out_unit,'(a)') 'test3.1 fails: err_msg='//TRIM(err_msg)
936            test_successful = .FALSE.
937         END IF
938      END DO
939      IF ( test_successful ) THEN
940         WRITE (out_unit,'(a)') 'test3.1 successful.'
941      ELSE
942         WRITE (out_unit,'(a)') 'test3.1 fails.'
943      END IF
945      !  intentional_error: data array too large in y direction
946      DO i=is1, ie1
947         IF ( id_dat1 > 0 ) used = send_data(id_dat1, dat2(i:i,:,:), Time, i-is1+1, 1, err_msg=err_msg)
948         IF ( err_msg /= '' ) EXIT ! exit immediately after error is detected. No need to continue.
949      END DO
950      IF ( err_msg == '' ) THEN
951         WRITE (out_unit,'(a)') 'test3.2 fails.'
952      ELSE
953         WRITE (out_unit,'(a)') 'test3.2 successful: err_msg='//TRIM(err_msg)
954      END IF
955   END IF
957   IF ( test_number == 2 ) THEN
958      !  multiple windows, no halos
959      !  No error messages should appear at any point within do loop for test2.1
960      test_successful = .TRUE.
961      DO j=js1, je1
962         IF ( id_dat1 > 0 ) used = send_data(id_dat1, dat1(:,j:j,:), Time, 1, j-js1+1, err_msg=err_msg)
963         IF ( err_msg /= '' ) THEN
964            WRITE (out_unit,'(a)') 'test2.1 fails: err_msg='//TRIM(err_msg)
965            test_successful = .FALSE.
966         END IF
967      END DO
968      IF ( test_successful ) THEN
969         WRITE (out_unit,'(a)') 'test2.1 successful.'
970      ELSE
971         WRITE (out_unit,'(a)') 'test2.1 fails.'
972      END IF
974      !  intentional_error: data array too large in x direction
975      DO j=js1, je1
976         IF ( id_dat1 > 0 ) used = send_data(id_dat1, dat2(:,j:j,:), Time, 1, j-js1+1, err_msg=err_msg)
977         IF ( err_msg /= '' ) EXIT ! exit immediately after error is detected. No need to continue.
978      END DO
979      IF ( err_msg == '' ) THEN
980         WRITE (out_unit,'(a)') 'test2.2 fails.'
981      ELSE
982         WRITE (out_unit,'(a)') 'test2.2 successful: err_msg='//TRIM(err_msg)
983      END IF
984   END IF
986   IF ( test_number == 1 ) THEN
987      !  1 window, no halos
988      ! Here dat2 is too large in both x and y direction so you should get an error.
989      IF ( id_dat1 > 0 ) used = send_data(id_dat1, dat2, Time, err_msg=err_msg)
990      IF ( err_msg == '' ) THEN
991         WRITE (out_unit,'(a)') 'test1.1 fails: Intentional error not detected'
992      ELSE
993         WRITE (out_unit,'(a)') 'test1.1 successful: '//TRIM(err_msg)
994      END IF
996      ! Here dat1 has the correct shape, so no error
997      IF ( id_dat1 > 0 ) used = send_data(id_dat1, dat1, Time, err_msg=err_msg)
999      IF ( err_msg == '' ) THEN
1000         WRITE (out_unit,'(a)') 'test1.2 successful'
1001      ELSE
1002         WRITE (out_unit,'(a)') 'test1.2 fails: '//TRIM(err_msg)
1003      END IF
1004   END IF
1005   CALL diag_manager_end(Time)
1006 END SELECT ! End of case handling opened for test 12.
1008   CALL fms_end
1010 CONTAINS
1012   SUBROUTINE compute_grid(nlon, nlat, is, ie, js, je, lon_global, lat_global, lonb_global, latb_global, lon, lat, &
1013                         & lonb, latb)
1014     INTEGER, INTENT(in) :: nlon, nlat, is, ie, js, je
1015     REAL, INTENT(out), DIMENSION(:) :: lon_global, lat_global, lonb_global, latb_global, lon, lat, lonb, latb
1017     REAL :: dlon, dlat
1018     INTEGER :: i, j
1020     dlon = 2*PI/nlon
1021     dlat = PI/nlat
1023     DO i=1, nlon+1
1024        lonb_global(i) = dlon*(i-1)
1025     END DO
1026     DO j=1,nlat+1
1027        latb_global(j) = dlat*(j-1) - .5*PI
1028     END DO
1029     DO i=1,nlon
1030        lon_global(i) = .5*(lonb_global(i) + lonb_global(i+1))
1031     END DO
1032     DO j=1,nlat
1033        lat_global(j) = .5*(latb_global(j) + latb_global(j+1))
1034     END DO
1035     lon  = lon_global(is:ie)
1036     lat  = lat_global(js:je)
1037     lonb = lonb_global(is:ie+1)
1038     latb = latb_global(js:je+1)
1039   END SUBROUTINE compute_grid
1041   SUBROUTINE unstruct_test(nx, ny, nz, npes, num_domain_tiles_x, num_domain_tiles_y, diag_time,io_tile_factor)
1042         use, intrinsic :: iso_fortran_env, only: output_unit
1043         use mpp_parameter_mod,             only: FATAL
1044         use mpp_mod,                       only: mpp_error, &
1045                                                  mpp_pe, &
1046                                                  mpp_root_pe, &
1047                                                  mpp_sync, &
1048                                                  mpp_chksum
1049         use mpp_domains_mod,               only: domain2D, &
1050                                                  mpp_define_mosaic, &
1051                                                  mpp_deallocate_domain, &
1052                                                  domainUG, &
1053                                                  mpp_define_unstruct_domain, &
1054 !                                                mpp_deallocate_domainUG, &
1055                                                  mpp_get_UG_compute_domain, &
1056                                                  mpp_get_UG_domain_grid_index, &
1057                                                  mpp_get_UG_domain_ntiles
1058         use diag_axis_mod,                 only: diag_axis_init, diag_axis_add_attribute
1059         use diag_manager_mod,              only: register_diag_field, &
1060                                                  send_data
1061         use time_manager_mod,              only: time_type, &
1062                                                  set_time, &
1063                                                  operator(+), &
1064                                                  assignment(=)
1065         implicit none
1067        !Inputs/Ouputs
1068         integer(kind=i4_kind),intent(in)  :: nx                 !<The number of grid points in the x-direction.
1069         integer(kind=i4_kind),intent(in)  :: ny                 !<The number of grid points in the y-direction.
1070         integer(kind=i4_kind),intent(in)  :: nz                 !<The number of grid points in the z-direction.
1071         integer(kind=i4_kind),intent(in)  :: npes               !<The total number of ranks used in this test.
1072         integer(kind=i4_kind),intent(in)  :: num_domain_tiles_x !<The total number of domain tiles in the x-dimension
1073                                                                 !! for the 2D structured domain in this test.
1074         integer(kind=i4_kind),intent(in)  :: num_domain_tiles_y !<The total number of domain tiles in the y-dimension
1075                                                                 !! for the 2D structured domain in this test.
1076         type(time_type),intent(inout) :: diag_time          !<Time for diag_manager.
1077         integer(kind=i4_kind),intent(in)  :: io_tile_factor     !<I/O tile factor.  See below.
1079        !Local variables
1080         integer(kind=i4_kind)                              :: num_domain_tiles !<The total number of domain tiles for
1081                                                                               !! the 2D structured domain in this test.
1082         integer(kind=i4_kind)                              :: npes_per_domain_tile !<The number of ranks per domain
1083                                                                                   !! tile for the 2D structured domain.
1084         integer(kind=i4_kind)                              :: my_domain_tile_id !<The 2D structured domain tile id for
1085                                                                                 !! the current rank.
1086         logical(kind=l4_kind)                              :: is_domain_tile_root !<Flag telling if the current rank
1087                                                        !! is the root rank of its associated 2D structured domain tile.
1088         integer(kind=i4_kind),dimension(2)                 :: layout_for_full_domain !<Rank layout (2D grid) for the
1089                                                    !! full 2D structured domain.
1090                                                    !! Example: 16 ranks -> (16,1) or (8,2) or (4,4) or (2,8) or (1,16)
1091         integer(kind=i4_kind),dimension(:),allocatable     :: pe_start !<Array holding the smallest rank id assigned
1092                                                                        !! to each 2D structured domain tile.
1093         integer(kind=i4_kind),dimension(:),allocatable     :: pe_end !<Array holding the largest rank id assigned to
1094                                                                      !! each 2D structured domain tile.
1095         integer(kind=i4_kind)                              :: x_grid_points_per_domain_tile !<The number of grid
1096                                                         !! points in the x-dimension on each 2D structured domain tile.
1097         integer(kind=i4_kind)                              :: y_grid_points_per_domain_tile !<The number of grid
1098                                                         !! points in the y-dimension on each 2D structured domain tile.
1099         integer(kind=i4_kind),dimension(:,:),allocatable   :: global_indices !<Required to define the 2D structured
1100                                                                              !! domain.
1101         integer(kind=i4_kind),dimension(:,:),allocatable   :: layout2D !<Required to define the 2D structured domain.
1102         type(domain2D)                                 :: domain_2D !<A structured 2D domain.
1103         logical(kind=l4_kind),dimension(:,:,:),allocatable :: land_mask                                  !<A toy mask.
1104         integer(kind=i4_kind),dimension(:),allocatable     :: num_non_masked_grid_points_per_domain_tile !<Total number
1105                                                         !! of non-masked grid points on each 2D structured domain tile.
1106         integer(kind=i4_kind)                              :: mask_counter !<Counting variable.
1107         integer(kind=i4_kind)                              :: num_non_masked_grid_points !<Total number of non-masked
1108                                                                           !! grid points for the 2D structured domain.
1109         integer(kind=i4_kind),dimension(:),allocatable     :: num_land_tiles_per_non_masked_grid_point   !<Number of
1110                                                  !! land tiles per non-masked grid point for the 2D structured domain.
1111         integer(kind=i4_kind)                              :: num_ranks_using_unstructured_grid !<Number of ranks
1112                                                                                  !! using the unstructured domain.
1113         integer(kind=i4_kind),dimension(:),allocatable     :: unstructured_grid_point_index_map !<Array that maps
1114                                                       !! indices between the 2D structured and unstructured domains.
1115         type(domainUG)                                 :: domain_ug !<An unstructured mpp domain.
1116         integer(kind=i4_kind),dimension(:),allocatable     :: unstructured_axis_data !<Data that is registered to the
1117                                                                              !! restart file for the unstructured axis.
1118         integer(kind=i4_kind)                              :: unstructured_axis_data_size !<Size of the unstructured
1119                                                                                           !! axis data array.
1120         character(len=256)                             :: unstructured_axis_name !<Name for the unstructured axis.
1121         real,dimension(:),allocatable                  :: x_axis_data !<Data for the x-axis that is registered to the
1122                                                                       !! restart file.
1123         real,dimension(:),allocatable                  :: y_axis_data !<Data for the y-axis that is registered to the
1124                                                                       !! restart file.
1125         real,dimension(:),allocatable                  :: z_axis_data !<Data for the z-axis that is registered to the
1126                                                                       !! restart file.
1127         real                                           :: unstructured_real_scalar_field_data_ref !<Reference test
1128                                                                       !! data for an unstructured real scalar field.
1129         real,dimension(:),allocatable                  :: unstructured_real_1D_field_data_ref !<Reference
1130                                                                         !! test data for an unstructured real 1D field.
1131         real,dimension(:,:),allocatable                :: unstructured_real_2D_field_data_ref !<Reference
1132                                                                         !! test data for an unstructured real 2D field.
1133         real,dimension(:,:,:),allocatable              :: unstructured_real_3D_field_data_ref !<Reference
1134                                                                         !! test data for an unstructured real 3D field.
1135         integer                                        :: unstructured_int_scalar_field_data_ref !<Reference
1136                                                                  !! test data for an unstructured integer scalar field.
1137         integer,dimension(:),allocatable               :: unstructured_int_1D_field_data_ref !<Reference
1138                                                                      !! test data for an unstructured integer 1D field.
1139         integer,dimension(:,:),allocatable             :: unstructured_int_2D_field_data_ref !<Reference
1140                                                                      !! test data for an unstructured integer 2D field.
1141         character(len=256)                             :: unstructured_real_scalar_field_name !<Name
1142                                                                               !! for an unstructured real scalar field.
1143         real                                           :: unstructured_real_scalar_field_data !<Data
1144                                                                               !! for an unstructured real scalar field.
1145         character(len=256)                             :: unstructured_real_1D_field_name !<Name for
1146                                                                                       !! an unstructured real 1D field.
1147         real,dimension(:),allocatable                  :: unstructured_real_1D_field_data !<Data for
1148                                                                                       !! an unstructured real 1D field.
1149         character(len=256)                             :: unstructured_real_2D_field_name !<Name for
1150                                                                                       !! an unstructured real 2D field.
1151         real,dimension(:,:),allocatable                :: unstructured_real_2D_field_data !<Data for
1152                                                                                       !! an unstructured real 2D field.
1153         character(len=256)                             :: unstructured_real_3D_field_name !<Name for
1154                                                                                       !! an unstructured real 3D field.
1155         real,dimension(:,:,:),allocatable              :: unstructured_real_3D_field_data !<Data for
1156                                                                                       !! an unstructured real 3D field.
1157         character(len=256)                             :: unstructured_int_scalar_field_name !<Name for
1158                                                                                !! an unstructured integer scalar field.
1159         integer                                        :: unstructured_int_scalar_field_data !<Data for
1160                                                                                !! an unstructured integer scalar field.
1161         character(len=256)                             :: unstructured_int_1D_field_name !<Name for an
1162                                                                                       !! unstructured integer 1D field.
1163         integer,dimension(:),allocatable               :: unstructured_int_1D_field_data !<Data for an
1164                                                                                       !! unstructured integer 1D field.
1165         character(len=256)                             :: unstructured_int_2D_field_name !<Name for an
1166                                                                                       !! unstructured integer 2D field.
1167         character(len=100)                             :: unstructured_1d_alt!<Name of the unstructured 1D field if L>1
1168         integer,dimension(:,:),allocatable             :: unstructured_int_2D_field_data !<Data for an
1169                                                                                       !! unstructured integer 2D field.
1170        integer(kind=i4_kind),allocatable,dimension(:)      :: unstructured_axis_diag_id !<Id returned
1171                                                                         !! for the unstructured axis by diag_axis_init.
1172        integer(kind=i4_kind)                               :: x_axis_diag_id !<Id returned for the x-axis
1173                                                                              !! by diag_axis_init.
1174        integer(kind=i4_kind)                              :: y_axis_diag_id !<Id returned for the y-axis
1175                                                                             !! by diag_axis_init.
1176        integer(kind=i4_kind)                              :: z_axis_diag_id !<Id returned for the z-axis
1177                                                                             !! by diag_axis_init.
1178        real,allocatable,dimension(:) :: lat, lon
1179        integer(kind=i4_kind)             :: idlat
1180        integer(kind=i4_kind)                              :: idlon
1181        integer(kind=i4_kind)                              :: rsf_diag_id !<Id returned for a real scalar field
1182                                                                          !! associated with the unstructured grid by
1183                                                                          !! register_diag_field.
1184        integer(kind=i4_kind),allocatable,dimension(:)     :: rsf_diag_1d_id !<Id returned for a real 1D array field
1185                                                       !! associated with the unstructured grid by register_diag_field.
1186        integer(kind=i4_kind)                              :: rsf_diag_2d_id !<Id returned for a real 2D array field
1187                                                       !! associated with the unstructured grid by register_diag_field.
1188         integer(kind=i4_kind)                              :: num_diag_time_steps !<Number of timesteps
1189                                                                                   !! (to simulate the model running).
1190         type(time_type)                                :: diag_time_start !<Starting time for the test.
1191         type(time_type)                                :: diag_time_step !<Time step for the test.
1192         logical(kind=l4_kind)                              :: used !<Return value from send data.
1194         integer(kind=i4_kind)                              :: i !<Loop variable.
1195         integer(kind=i4_kind)                              :: j !<Loop variable.
1196         integer(kind=i4_kind)                              :: k,l=1 !<Loop variable.
1197         integer(kind=i4_kind)                              :: p !<Counting variable.
1199        !Needed to define the 2D structured domain but never used.
1200         integer(kind=i4_kind)              :: ncontacts
1201         integer(kind=i4_kind),dimension(20) :: tile1
1202         integer(kind=i4_kind),dimension(20) :: tile2
1203         integer(kind=i4_kind),dimension(20) :: istart1
1204         integer(kind=i4_kind),dimension(20) :: iend1
1205         integer(kind=i4_kind),dimension(20) :: jstart1
1206         integer(kind=i4_kind),dimension(20) :: jend1
1207         integer(kind=i4_kind),dimension(20) :: istart2
1208         integer(kind=i4_kind),dimension(20) :: iend2
1209         integer(kind=i4_kind),dimension(20) :: jstart2
1210         integer(kind=i4_kind),dimension(20) :: jend2
1212         integer(kind=i4_kind),dimension(3)  :: npes_io_group
1214        !Print out a message that the test is starting.
1215         if (mpp_pe() .eq. mpp_root_pe()) then
1216             write(output_unit,*)
1217             write(output_unit,*) "</----------------------------------------"
1218             write(output_unit,*) "Test create_unstructured_test_restart_file" &
1219                                  //" starting ..."
1220             write(output_unit,*)
1221         endif
1223        !Synchronize all ranks.
1224         call mpp_sync()
1226        !Make sure that valid inputs were passed in.
1227         if (nx .lt. 1 .or. ny .lt. 1) then
1228             call mpp_error(FATAL, &
1229                            "create_unstructured_test_restart_file:" &
1230                            //" there must be at least on grid point in the" &
1231                            //" x- and y- dimensions.")
1232         endif
1233         if (npes .gt. nx*ny) then
1234             call mpp_error(FATAL, &
1235                            "create_unstructured_test_restart_file:" &
1236                            //" the total number of ranks cannot be greater" &
1237                            //" than the total number of grid points in the" &
1238                            //" x-y plane.")
1239         endif
1240         if (num_domain_tiles_x .lt. 1 .or. num_domain_tiles_y .lt. 1) then
1241             call mpp_error(FATAL, &
1242                            "create_unstructured_test_restart_file:" &
1243                            //" there must be at least on domain tile in the" &
1244                            //" x- and y- dimensions.")
1245         endif
1246         if (mod(nx,num_domain_tiles_x) .ne. 0) then
1247             call mpp_error(FATAL, &
1248                            "create_unstructured_test_restart_file:" &
1249                            //" the total number of grid points in the" &
1250                            //" x-dimension must be evenly divisible by the" &
1251                            //" total number of domain tiles in the" &
1252                            //" x-dimension.")
1253         endif
1254         if (mod(ny,num_domain_tiles_y) .ne. 0) then
1255             call mpp_error(FATAL, &
1256                            "create_unstructured_test_restart_file:" &
1257                            //" the total number of grid points in the" &
1258                            //" y-dimension must be evenly divisible by the" &
1259                            //" total number of domain tiles in the" &
1260                            //" y-dimension.")
1261         endif
1262         if (num_domain_tiles_x*num_domain_tiles_y .gt. npes) then
1263             call mpp_error(FATAL, &
1264                            "create_unstructured_test_restart_file:" &
1265                            //" the total number of domain tiles cannot be" &
1266                            //" greater than the total number of ranks.")
1267         endif
1268         if (mod(npes,num_domain_tiles_x) .ne. 0) then
1269             call mpp_error(FATAL, &
1270                            "create_unstructured_test_restart_file:" &
1271                            //" the total number of ranks must be evenly" &
1272                            //" divisible by the total number of domain" &
1273                            //" tiles in the x-dimension.")
1274         endif
1275         if (mod(npes,num_domain_tiles_y) .ne. 0) then
1276             call mpp_error(FATAL, &
1277                            "create_unstructured_test_restart_file:" &
1278                            //" the total number of ranks must be evenly" &
1279                            //" divisible by the total number of domain" &
1280                            //" tiles in the y-dimension.")
1281         endif
1283        !Set domain tile values for the 2D structured domain.
1284         num_domain_tiles = num_domain_tiles_x*num_domain_tiles_y
1285         npes_per_domain_tile = npes/num_domain_tiles
1286         my_domain_tile_id = (mpp_pe())/npes_per_domain_tile + 1
1287         if (mpp_pe() .eq. (my_domain_tile_id-1)*npes_per_domain_tile) then
1288             is_domain_tile_root = .true.
1289         else
1290             is_domain_tile_root = .false.
1291         endif
1292         layout_for_full_domain(1) = num_domain_tiles_x
1293         layout_for_full_domain(2) = npes/layout_for_full_domain(1)
1295        !For each 2D structured domain tile, store the beginning and ending
1296        !rank ids assigned to it.  For example, if there are 8 ranks and 2
1297        !domain tiles, then tile 1 will be assigned ranks 0 - 3 and tile 2
1298        !will be assigned ranks 4 - 7.
1299         allocate(pe_start(num_domain_tiles))
1300         allocate(pe_end(num_domain_tiles))
1301         do i = 1,num_domain_tiles
1302             pe_start(i) = (i-1)*npes_per_domain_tile
1303             pe_end(i) = i*npes_per_domain_tile - 1
1304         enddo
1306        !Calculate parameters needed to construct the 2D structured domain.
1307        !All domain tiles are assumed to be the same size.
1308         x_grid_points_per_domain_tile = nx/num_domain_tiles_x
1309         y_grid_points_per_domain_tile = ny/num_domain_tiles_y
1310         allocate(global_indices(4,num_domain_tiles))
1311         do i = 1,num_domain_tiles
1312             global_indices(:,i) = (/1,x_grid_points_per_domain_tile, &
1313                                     1,y_grid_points_per_domain_tile/)
1314         enddo
1315         allocate(layout2D(2,num_domain_tiles))
1316         do i = 1,num_domain_tiles
1317             layout2D(1,i) = layout_for_full_domain(1)/num_domain_tiles_x
1318             layout2D(2,i) = layout_for_full_domain(2)/num_domain_tiles_y
1319         enddo
1321        !This test does not use the "contact" region between tiles, but
1322        !the 2D structured domain requires these inputs, so just set them
1323        !all equal to 1.
1324         ncontacts = 1
1325         tile1 = 1
1326         tile2 = 1
1327         istart1 = 1
1328         iend1 = 1
1329         jstart1 = 1
1330         jend1 = 1
1331         istart2 = 1
1332         iend2 = 1
1333         jstart2 = 1
1334         jend2 = 1
1335 !write (6,*)size(tile1)
1336        !Define the 2D structured domain.
1337         call mpp_define_mosaic(global_indices, &
1338                                layout2D, &
1339                                domain_2D, &
1340                                num_domain_tiles, &
1341                                0, &
1342                                tile1, &
1343                                tile2, &
1344                                istart1, &
1345                                iend1, &
1346                                jstart1, &
1347                                jend1, &
1348                                istart2, &
1349                                iend2, &
1350                                jstart2, &
1351                                jend2, &
1352                                pe_start, &
1353                                pe_end)
1355        !Define a toy mask to mimic what happens in the land model.
1356         allocate(land_mask(x_grid_points_per_domain_tile, &
1357                            y_grid_points_per_domain_tile, &
1358                            num_domain_tiles))
1359         allocate(num_non_masked_grid_points_per_domain_tile(num_domain_tiles))
1360         land_mask = .false.
1361         do k = 1,num_domain_tiles
1362             mask_counter = 0
1363             do j = 1,y_grid_points_per_domain_tile
1364                 do i = 1,x_grid_points_per_domain_tile
1365                     if (mod((k-1)*y_grid_points_per_domain_tile*x_grid_points_per_domain_tile + &
1366                             (j-1)*x_grid_points_per_domain_tile + &
1367                             (i-1),2) .eq. 0) then
1368                         land_mask(i,j,k) = .true.
1369                         mask_counter = mask_counter + 1
1370                     endif
1371                 enddo
1372             enddo
1373             num_non_masked_grid_points_per_domain_tile(k) = mask_counter
1374         enddo
1376        !Set the number of land tiles allowed per non-masked grid point.
1377         num_non_masked_grid_points = sum(num_non_masked_grid_points_per_domain_tile)
1378         allocate(num_land_tiles_per_non_masked_grid_point(num_non_masked_grid_points))
1379         num_land_tiles_per_non_masked_grid_point = 1
1381        !Set the number of ranks to use with the unstructured domain.  There
1382        !must be at least one grid point per rank.
1383         num_ranks_using_unstructured_grid = npes
1384         if (num_ranks_using_unstructured_grid .gt. num_non_masked_grid_points) then
1385             call mpp_error(FATAL, &
1386                            "create_unstructured_test_restart_file:" &
1387                            //" the number of ranks exceeds the number of" &
1388                            //" non-masked grid points for the unstructured" &
1389                            //" domain.")
1390         endif
1392        !Define an array used to map grid points from the "structured" 2D grid
1393        !to the "unstructured" 1D grid.  The mapping goes as follows (fortran
1394        !ording so first index is fastest):
1395        !
1396        ! 2D "structured" grid (lon,lat,tile) => 1D "unstructured" grid (p)
1397        !
1398        !where masked points are skipped.
1399         allocate(unstructured_grid_point_index_map(num_non_masked_grid_points))
1400         p = 0
1401         do k = 1,num_domain_tiles
1402             do j = 1,y_grid_points_per_domain_tile
1403                 do i = 1,x_grid_points_per_domain_tile
1404                     if (land_mask(i,j,k)) then
1405                         p = p + 1
1406                         unstructured_grid_point_index_map(p) = (j-1)*x_grid_points_per_domain_tile + i
1407                     endif
1408                 enddo
1409             enddo
1410         enddo
1411        !> Set in namelist is "I/O tile factor".  The number of ranks that
1412        !! participate in I/O for a tile is equal to:
1413        !!
1414        !! num_io_ranks_on_a_tile = num_ranks_on_the_tile / "I/O tile factor".
1415        !!
1416        !!so for:
1417        !!
1418        !! io_tile_factor = 1, all of the ranks on a tile participate in the I/O
1419        !! io_tile_factor = 2, 1/2 of the ranks on a tile participate in the I/O
1420        !! io_tile_factor = 3, 1/3 of the ranks on a tile participate in the I/O
1421        !! ...
1422        !! io_tile_factor = 0 is a special case where only one rank participates
1423        !!                  in the I/O for a tile.
1424        !! io_tile_factor = 1
1425 if (mpp_pe() == mpp_root_pe()) write(6,*) "IO_TILE_FACTOR is ",io_tile_factor
1426 allocate(unstructured_axis_diag_id(1))
1427 allocate(rsf_diag_1d_id(1))
1429        !Define the "unstructured" domain decomposition.
1430         call mpp_define_unstruct_domain(domain_ug, &
1431                                         domain_2D, &
1432                                         num_non_masked_grid_points_per_domain_tile, &
1433                                         num_land_tiles_per_non_masked_grid_point, &
1434                                         num_ranks_using_unstructured_grid, &
1435                                         io_tile_factor, &
1436                                         unstructured_grid_point_index_map)
1438 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1439 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1440 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
1441 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!      !Don't need to modify above here!      !!!!!!!!!!!!!!!!!!!!!!!!!!!!
1443        !Get the that will be registered for the unstructured axis. This should
1444        !be each rank's unstructured compute domain (I think, because a gather
1445        !is performed by the root of each I/O domain pelist.
1446         call mpp_get_UG_compute_domain(domain_ug,size=unstructured_axis_data_size)
1447         if(.not.allocated(unstructured_axis_data))allocate(unstructured_axis_data(unstructured_axis_data_size))
1448 !! THIS IS A PROBLEM !!
1449         call mpp_get_UG_domain_grid_index(domain_ug,unstructured_axis_data)
1450 !write(6,*)"ID:",mpp_pe()," DATA: ",unstructured_axis_data
1451        !Initialize the "unstructured" axis for the diagnostics.
1452         unstructured_axis_name = "ug_axis"
1454         unstructured_axis_diag_id(l) = diag_axis_init(trim(unstructured_axis_name), &
1455                                                    real(unstructured_axis_data), &
1456                                                    "none", &
1457                                                    "U", &
1458                                                    long_name="mapping indices", &
1459                                                    domainU=domain_ug)
1460    call diag_axis_add_attribute(unstructured_axis_diag_id(l),'compress','grid_xt grid_yt')
1462 !write(6,*) "ID U",unstructured_axis_diag_id
1463        !Add the x-, y-, and z-axes to the restart file.  Until a bug in
1464        !the code is resolved, I must register the unstructured axis first.
1465        !Also initialize the axes for the diagnostics.
1466         if (.not.allocated(x_axis_data)) allocate(x_axis_data(nx))
1467 !        if (.not.allocated(y_axis_data))allocate(y_axis_data(ny))
1468 !! ASSUMES 4 PEs!!!
1469 ! if (mpp_pe() > 4) call error_mesg("Diag_test_unstruct","Only 4 PEs please",fatal)
1470      do i=1,nx
1471           x_axis_data(i) = real(i)
1472      enddo
1473 !     if (mod(mpp_pe(),2).eq.0) then
1474 !        do j = 1,ny/4
1475 !            y_axis_data(j) = real(j)
1476 !        enddo
1478 !     else
1479 !        do j = 1,ny/4
1480 !            y_axis_data(j) = real(j+ny/4)
1481 !        enddo
1482 !     endif
1484        x_axis_diag_id = diag_axis_init("grid_xt", &
1485                                        x_axis_data, &
1486                                        "degrees", &
1487                                        "X", &
1488                                        long_name="longitude")
1490         if (.not.allocated(y_axis_data))allocate(y_axis_data(ny/num_domain_tiles_y))
1491         do i = 1,ny/num_domain_tiles_y
1492             y_axis_data(i) = real(i)
1493         enddo
1494        y_axis_diag_id = diag_axis_init("grid_yt", &
1495                                        y_axis_data, &
1496                                        "degrees", &
1497                                        "Y", &
1498                                        long_name="latitude")
1500         if (.not.allocated(z_axis_data))allocate(z_axis_data(nz))
1501         do i = 1,nz
1502             z_axis_data(i) = real(i*5.0)
1503         enddo
1504        z_axis_diag_id = diag_axis_init("zfull", &
1505                                        z_axis_data, &
1506                                        "km", &
1507                                        "Z", &
1508                                        long_name="dont look down")
1509 !write (6,*) z_axis_diag_id
1511        !Define some reference test data.
1513        !real scalar field.
1514         unstructured_real_scalar_field_data_ref = 1234.5678*real(l)
1516        !real 1D field.
1517         if (.not.allocated(unstructured_real_1D_field_data_ref)) &
1518               & allocate(unstructured_real_1D_field_data_ref(unstructured_axis_data_size))
1519         do i = 1,unstructured_axis_data_size
1520             unstructured_real_1D_field_data_ref(i) = real(i) *real(i)+0.1*(mpp_pe()+1)
1521         enddo
1523        !real 2D field.
1524         if (.not.allocated(unstructured_real_2D_field_data_ref)) &
1525               & allocate(unstructured_real_2D_field_data_ref(unstructured_axis_data_size,nz))
1526         do j = 1,nz
1527             do i = 1,unstructured_axis_data_size
1528                 unstructured_real_2D_field_data_ref(i,j) = real(j)+0.1*(mpp_pe()+1.0)
1529                                                            !-1.0*real((j-1)* &
1530                                                            !unstructured_axis_data_size+i) &
1531                                                            !+ 1.1111111*real(l)
1532             enddo
1533         enddo
1535        !real 3D field.
1536 !       if(.not.allocated(unstructured_real_3D_field_data_ref)
1537 !             allocate(unstructured_real_3D_field_data_ref(unstructured_axis_data_size,nz,cc_axis_size))
1538 !       do k = 1,cc_axis_size
1539 !           do j = 1,nz
1540 !               do i = 1,unstructured_axis_data_size
1541 !                   unstructured_real_3D_field_data_ref(i,j,k) = -1.0*real((k-1)*nz* &
1542 !                                                                unstructured_axis_data_size+(j-1)* &
1543 !                                                                unstructured_axis_data_size+i) &
1544 !                                                                + 2.2222222
1545 !               enddo
1546 !           enddo
1547 !       enddo
1549        !integer scalar field.
1550         unstructured_int_scalar_field_data_ref = 7654321*L
1552        !integer 1D field.
1553         if (.not.allocated(unstructured_int_1D_field_data_ref)) &
1554               & allocate(unstructured_int_1D_field_data_ref(unstructured_axis_data_size))
1555         do i = 1,unstructured_axis_data_size
1556             unstructured_int_1D_field_data_ref(i) = i - 8*l
1557         enddo
1559        !integer 2D field.
1560         if (.not.allocated(unstructured_int_2D_field_data_ref)) &
1561               & allocate(unstructured_int_2D_field_data_ref(unstructured_axis_data_size,nz))
1562         do j = 1,nz
1563             do i = 1,unstructured_axis_data_size
1564                 unstructured_int_2D_field_data_ref(i,j) = -1*((j-1)*unstructured_axis_data_size+i) + 2*L
1565             enddo
1566         enddo
1568      !> Latitude and Longitude
1569      allocate(lat(ny/num_domain_tiles_y),lon(nx))
1570      do i=1,nx
1571           lon(i) = real(i)*360.0/real(nx)
1572      enddo
1573      do j=1,ny/num_domain_tiles_y
1574           lat(j) = real(j)*180.8/real(ny)
1575      enddo
1577        !Add a real scalar field to the restart file.  Initialize it as a
1578        !diagnostic.
1579         unstructured_real_scalar_field_name = "unstructured_real_scalar_field_1"
1580         unstructured_real_scalar_field_data = unstructured_real_scalar_field_data_ref
1582        idlon = register_diag_field("UG_unit_test", &
1583                                          "lon", &
1584                                          (/x_axis_diag_id/),&
1585                                          init_time=diag_time, &
1586                                          long_name="E-W longitude", &
1587                                          units="degrees")
1588 l=SIZE(unstructured_axis_diag_id)
1590        rsf_diag_id = register_diag_field("UG_unit_test", &
1591                                          "unstructured_real_scalar_field_data", &
1592                                          init_time=diag_time, &
1593                                          long_name="rsf_diag_1", &
1594                                          units="ergs")
1595        rsf_diag_1d_id(1) = register_diag_field("UG_unit_test", &
1596                                          "unstructured_real_1D_field_data", &
1597                                          (/unstructured_axis_diag_id(1)/),&
1598                                          init_time=diag_time, &
1599                                          long_name="ONE_D_ARRAY", &
1600                                          units="ergs")
1602        rsf_diag_2d_id = register_diag_field("UG_unit_test", &
1603                                          "unstructured_real_2D_field_data", &
1604                                          (/unstructured_axis_diag_id(1), z_axis_diag_id/),&
1605                                          init_time=diag_time, &
1606                                          long_name="TWO_D_ARRAY", &
1607                                          units="ergs")
1609        idlat = register_diag_field("UG_unit_test", &
1610                                          "lat", &
1611                                          (/y_axis_diag_id/),&
1612                                          init_time=diag_time, &
1613                                          long_name="S-N latitude", &
1614                                          units="degrees")
1617 IF (l .NE. 1) THEN
1618   do l=2,3
1619    write(unstructured_1d_alt,'(a,I0)') "unstructured_real_1D",L
1620    rsf_diag_1d_id(L) = register_diag_field ("UG_unit_test", trim(unstructured_1d_alt),&
1621                                           (/unstructured_axis_diag_id(L)/),&
1622                                            init_time=diag_time, &
1623                                            long_name="OTHER"//trim(unstructured_1d_alt), &
1624                                            units="kg")
1625   enddo
1626 ENDIF !L.ne.1
1627        !Add a real 1D field to the restart file.  This field is of the form:
1628        !field = field(unstructured).
1629         unstructured_real_1D_field_name = "unstructured_real_1D_field_1"
1630         if (.not.allocated(unstructured_real_1D_field_data)) &
1631               & allocate(unstructured_real_1D_field_data(unstructured_axis_data_size))
1632         unstructured_real_1D_field_data = unstructured_real_1D_field_data_ref
1634        !Add a real 2D field to the restart file.  This field is of the form:
1635        !field = field(unstructured,z).
1636         unstructured_real_2D_field_name = "unstructured_real_2D_field_1"
1637        if (.not.allocated(unstructured_real_2D_field_data)) &
1638               & allocate(unstructured_real_2D_field_data(unstructured_axis_data_size,nz))
1639        unstructured_real_2D_field_data = unstructured_real_2D_field_data_ref
1640 !       allocate(unstructured_real_2D_field_data(unstructured_axis_data_size,nx))
1641 !       unstructured_real_2D_field_data = 1
1643        !Add a real 3D field to the restart file.  This field is of the form:
1644        !field = field(unstructured,z,cc).
1645 !       unstructured_real_3D_field_name = "unstructured_real_3D_field_1"
1646 !       if (.not.allocated(unstructured_real_3D_field_data)) &
1647 !              & allocate(unstructured_real_3D_field_data(unstructured_axis_data_size,nz,cc_axis_size))
1648 !       unstructured_real_3D_field_data = unstructured_real_3D_field_data_ref
1650        !Add an integer scalar field to the restart file.
1651         unstructured_int_scalar_field_name = "unstructured_int_scalar_field_1"
1652         unstructured_int_scalar_field_data = unstructured_int_scalar_field_data_ref
1654        !Add an integer 1D field to the restart file.  This field is of the
1655        !from: field = field(unstructured).
1656         unstructured_int_1D_field_name = "unstructured_int_1D_field_1"
1657         if (.not.allocated(unstructured_int_1D_field_data)) &
1658               & allocate(unstructured_int_1D_field_data(unstructured_axis_data_size))
1659         unstructured_int_1D_field_data = unstructured_int_1D_field_data_ref
1661        !Add an integer 2D field to the restart file.  This field is of the
1662        !form: field = field(unstructured,z).
1663         unstructured_int_2D_field_name = "unstructured_int_2D_field_1"
1664         if (.not.allocated(unstructured_int_2D_field_data)) &
1665               & allocate(unstructured_int_2D_field_data(unstructured_axis_data_size,nz))
1666         unstructured_int_2D_field_data = unstructured_int_2D_field_data_ref
1668        !Simulate the model timesteps, so that diagnostics may be written
1669        !out.
1670         num_diag_time_steps = 4
1671         diag_time_step = set_time(12*3600, 0)
1672         diag_time_start = diag_time
1673 ! used = send_data(idlat,lat,diag_time)
1674 ! used = send_data(idlon,lon,diag_time)
1675         do i = 1,num_diag_time_steps
1677            !Update the current time.
1678             diag_time = diag_time  + diag_time_step
1680            !"Evolve" the test data.
1681             unstructured_real_scalar_field_data_ref = unstructured_real_scalar_field_data_ref + &
1682                                                       real(1)
1683             unstructured_real_scalar_field_data = unstructured_real_scalar_field_data_ref
1685            !Update the data.
1686            if (rsf_diag_id .gt. 0) then
1687                used = send_data(rsf_diag_id, &
1688                                 unstructured_real_scalar_field_data, &
1689                                 diag_time)
1690            endif
1692         IF (SIZE(rsf_diag_1d_id) == 1) THEN
1693           used = send_data(rsf_diag_1d_id(1), &
1694                                 unstructured_real_1D_field_data, &
1695                                 diag_time)
1696          ELSE
1697           DO L=1,3
1698            used = send_data(rsf_diag_1d_id(L), &
1699                                 unstructured_real_1D_field_data, &
1700                                 diag_time)
1701           ENDDO
1702          ENDIF
1703           used = send_data(rsf_diag_2d_id, &
1704                                 unstructured_real_2D_field_data, &
1705                                 diag_time)
1706  used = send_data(idlat,lat,diag_time)
1707  used = send_data(idlon,lon,diag_time)
1709         enddo
1710        !Deallocate the unstructured domain.
1711         call mpp_sync()
1712 !       call mpp_deallocate_domainUG(domain_ug)
1714        !Deallocate the 2D structured domain.
1715         call mpp_deallocate_domain(domain_2D)
1717        !Deallocate local allocatables.
1718         deallocate(pe_start)
1719         deallocate(pe_end)
1720         deallocate(global_indices)
1721         deallocate(layout2D)
1722         deallocate(land_mask)
1723         deallocate(num_non_masked_grid_points_per_domain_tile)
1724         deallocate(num_land_tiles_per_non_masked_grid_point)
1725         deallocate(unstructured_grid_point_index_map)
1726         deallocate(x_axis_data)
1727         deallocate(y_axis_data)
1728         deallocate(z_axis_data)
1729         deallocate(unstructured_axis_data)
1730         deallocate(unstructured_real_1D_field_data_ref)
1731         deallocate(unstructured_real_2D_field_data_ref)
1732 !       deallocate(unstructured_real_3D_field_data_ref)
1733         deallocate(unstructured_int_1D_field_data_ref)
1734         deallocate(unstructured_int_2D_field_data_ref)
1735         deallocate(unstructured_real_1D_field_data)
1736         deallocate(unstructured_real_2D_field_data)
1737 !       deallocate(unstructured_real_3D_field_data)
1738         deallocate(unstructured_int_1D_field_data)
1739         deallocate(unstructured_int_2D_field_data)
1743        !Print out a message that the test is done.
1744         call mpp_sync()
1745         if (mpp_pe() .eq. mpp_root_pe()) then
1746             write(output_unit,*)
1747             write(output_unit,*) "Test create_unstructured_test_restart_file" &
1748                                  //" complete."
1749             write(output_unit,*) "----------------------------------------/>"
1750             write(output_unit,*)
1751         endif
1754         return
1755   END SUBROUTINE unstruct_test
1757 END PROGRAM test