1 | |
---|
2 | !> Provides reusable IO routines for OASIS |
---|
3 | |
---|
4 | MODULE mod_oasis_io |
---|
5 | |
---|
6 | USE mod_oasis_kinds |
---|
7 | USE mod_oasis_data |
---|
8 | USE mod_oasis_parameters |
---|
9 | USE mod_oasis_mpi |
---|
10 | USE mod_oasis_sys |
---|
11 | USE mod_oasis_ioshr |
---|
12 | USE mct_mod |
---|
13 | USE netcdf |
---|
14 | |
---|
15 | implicit none |
---|
16 | |
---|
17 | private |
---|
18 | |
---|
19 | !--- interfaces --- |
---|
20 | public :: oasis_io_varexists |
---|
21 | public :: oasis_io_read_avfld |
---|
22 | public :: oasis_io_read_avfile |
---|
23 | public :: oasis_io_write_avfile |
---|
24 | public :: oasis_io_read_array |
---|
25 | public :: oasis_io_write_array |
---|
26 | public :: oasis_io_read_avfbf |
---|
27 | public :: oasis_io_write_avfbf |
---|
28 | public :: oasis_io_write_2dgridint_fromroot |
---|
29 | public :: oasis_io_write_2dgridfld_fromroot |
---|
30 | public :: oasis_io_write_3dgridfld_fromroot |
---|
31 | public :: oasis_io_write_fldattr |
---|
32 | public :: oasis_io_read_field_fromroot |
---|
33 | |
---|
34 | integer(kind=ip_i4_p) :: oasis_io_64bit_offset = nf90_64bit_offset |
---|
35 | #ifdef CDF_64BIT_DATA |
---|
36 | integer(kind=ip_i4_p) :: oasis_io_64bit_data = nf90_64bit_data |
---|
37 | #else |
---|
38 | integer(kind=ip_i4_p) :: oasis_io_64bit_data = nf90_64bit_offset ! set to offset by default |
---|
39 | #endif |
---|
40 | |
---|
41 | !=========================================================================== |
---|
42 | CONTAINS |
---|
43 | !=========================================================================== |
---|
44 | |
---|
45 | !=============================================================================== |
---|
46 | |
---|
47 | !> Checks whether the var fldname is in the file |
---|
48 | |
---|
49 | logical function oasis_io_varexists(filename,fldname) |
---|
50 | |
---|
51 | implicit none |
---|
52 | |
---|
53 | character(len=*), intent(in) :: filename !< filename |
---|
54 | character(len=*), intent(in) :: fldname !< fldname |
---|
55 | |
---|
56 | !--- local --- |
---|
57 | logical :: exists ! file exist flag |
---|
58 | integer(ip_i4_p) :: ncid ! netcf file id |
---|
59 | integer(ip_i4_p) :: varid ! netcf variable id |
---|
60 | integer(ip_i4_p) :: status ! error code |
---|
61 | character(len=*),parameter :: subname = '(oasis_io_varexists)' |
---|
62 | |
---|
63 | !------------------------------------------------------------------------------- |
---|
64 | call oasis_debug_enter(subname) |
---|
65 | |
---|
66 | oasis_io_varexists = .false. |
---|
67 | inquire(file=trim(filename),exist=exists) |
---|
68 | if (exists) then |
---|
69 | status = nf90_open(filename,NF90_NOWRITE,ncid) |
---|
70 | call check_status(status,subname,__FILE__,__LINE__) |
---|
71 | status = nf90_inq_varid(ncid,trim(fldname),varid) |
---|
72 | if (status == NF90_NOERR) then |
---|
73 | oasis_io_varexists = .true. |
---|
74 | endif |
---|
75 | status = nf90_close(ncid) |
---|
76 | call check_status(status,subname,__FILE__,__LINE__) |
---|
77 | endif |
---|
78 | |
---|
79 | call oasis_debug_exit(subname) |
---|
80 | |
---|
81 | end function oasis_io_varexists |
---|
82 | |
---|
83 | !=============================================================================== |
---|
84 | |
---|
85 | !> Reads single field from a file into an attribute Vector |
---|
86 | |
---|
87 | subroutine oasis_io_read_avfld(filename,av,gsmap,mpicom,avfld,filefld,fldtype) |
---|
88 | |
---|
89 | implicit none |
---|
90 | |
---|
91 | character(len=*), intent(in) :: filename !< file name |
---|
92 | type(mct_aVect) , intent(inout) :: av !< avect |
---|
93 | type(mct_gsmap) , intent(in) :: gsmap !< gsmap decomp |
---|
94 | integer(ip_i4_p), intent(in) :: mpicom !< mpi comm |
---|
95 | character(len=*), intent(in) :: avfld !< av field name |
---|
96 | character(len=*), intent(in) :: filefld !< file field name |
---|
97 | character(len=*), intent(in),optional :: fldtype !< kind |
---|
98 | |
---|
99 | !--- local --- |
---|
100 | integer(ip_i4_p) :: n,n1,i,j,fk,fk1 ! index |
---|
101 | integer(ip_i4_p) :: nx ! 2d global size nx |
---|
102 | integer(ip_i4_p) :: ny ! 2d global size ny |
---|
103 | type(mct_aVect) :: av_g ! avect global data |
---|
104 | integer(ip_i4_p) :: master_task,iam,ierr ! mpi info |
---|
105 | integer(ip_i4_p) :: ncid,dimid,dimid2(2),varid ! netcdf info |
---|
106 | integer(ip_i4_p) :: dlen ! dimension length |
---|
107 | integer(ip_i4_p) :: status ! error code |
---|
108 | logical :: exists ! file existance |
---|
109 | real(ip_double_p),allocatable :: array2(:,:) |
---|
110 | integer(ip_i4_p) ,allocatable :: array2i(:,:) |
---|
111 | integer(ip_i4_p) :: ifldtype ! field type int (1) or real (2) |
---|
112 | |
---|
113 | character(len=*),parameter :: subname = '(oasis_io_read_avfld)' |
---|
114 | |
---|
115 | !------------------------------------------------------------------------------- |
---|
116 | ! |
---|
117 | !------------------------------------------------------------------------------- |
---|
118 | |
---|
119 | if (mpicom == MPI_COMM_NULL) return |
---|
120 | |
---|
121 | ! empty filename, just return |
---|
122 | |
---|
123 | call oasis_debug_enter(subname) |
---|
124 | if (len_trim(filename) < 1) then |
---|
125 | call oasis_debug_exit(subname) |
---|
126 | return |
---|
127 | endif |
---|
128 | |
---|
129 | master_task = 0 |
---|
130 | call MPI_COMM_RANK(mpicom,iam,ierr) |
---|
131 | |
---|
132 | ifldtype = 2 ! real default |
---|
133 | if (present(fldtype)) then |
---|
134 | ifldtype = 0 |
---|
135 | if (trim(fldtype) == 'int') ifldtype = 1 |
---|
136 | if (trim(fldtype) == 'real') ifldtype = 2 |
---|
137 | if (ifldtype == 0) then |
---|
138 | write(nulprt,*) subname,estr,'in fldtype argument' |
---|
139 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
140 | endif |
---|
141 | endif |
---|
142 | |
---|
143 | call mct_aVect_gather(av,av_g,gsmap,master_task,mpicom) |
---|
144 | |
---|
145 | if (iam == master_task) then |
---|
146 | |
---|
147 | inquire(file=trim(filename),exist=exists) |
---|
148 | if (exists) then |
---|
149 | status = nf90_open(trim(filename),NF90_NOWRITE,ncid) |
---|
150 | call check_status(status,subname,__FILE__,__LINE__) |
---|
151 | else |
---|
152 | write(nulprt,*) subname,estr,'file missing ',trim(filename) |
---|
153 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
154 | endif |
---|
155 | |
---|
156 | status = nf90_inq_varid(ncid,trim(filefld),varid) |
---|
157 | call check_status(status,subname,__FILE__,__LINE__) |
---|
158 | status = nf90_inquire_variable(ncid,varid,ndims=dlen,dimids=dimid2) |
---|
159 | call check_status(status,subname,__FILE__,__LINE__) |
---|
160 | if (dlen > 2) then |
---|
161 | write(nulprt,*) subname,estr,'variable ndims gt 2 ',trim(filefld),dlen |
---|
162 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
163 | endif |
---|
164 | status = nf90_inquire_dimension(ncid,dimid2(1),len=nx) |
---|
165 | call check_status(status,subname,__FILE__,__LINE__) |
---|
166 | ny = 1 |
---|
167 | if (dlen == 2) then |
---|
168 | status = nf90_inquire_dimension(ncid,dimid2(2),len=ny) |
---|
169 | call check_status(status,subname,__FILE__,__LINE__) |
---|
170 | endif |
---|
171 | |
---|
172 | if (size(av_g%rAttr,dim=2) /= nx*ny) then |
---|
173 | write(nulprt,*) subname,estr,'av gsize nx ny mismatch in file :',& |
---|
174 | trim(filename),SIZE(av_g%rAttr,dim=2),nx,ny |
---|
175 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
176 | endif |
---|
177 | |
---|
178 | if (ifldtype == 1) then |
---|
179 | allocate(array2i(nx,ny)) |
---|
180 | status = nf90_get_var(ncid,varid,array2i) |
---|
181 | call check_status(status,subname,__FILE__,__LINE__) |
---|
182 | |
---|
183 | n = mct_avect_indexIA(av_g,trim(avfld)) |
---|
184 | n1 = 0 |
---|
185 | do j = 1,ny |
---|
186 | do i = 1,nx |
---|
187 | n1 = n1 + 1 |
---|
188 | av_g%iAttr(n,n1) = array2i(i,j) |
---|
189 | enddo |
---|
190 | enddo |
---|
191 | deallocate(array2i) |
---|
192 | else |
---|
193 | allocate(array2(nx,ny)) |
---|
194 | status = nf90_get_var(ncid,varid,array2) |
---|
195 | call check_status(status,subname,__FILE__,__LINE__) |
---|
196 | |
---|
197 | n = mct_avect_indexRA(av_g,trim(avfld)) |
---|
198 | n1 = 0 |
---|
199 | do j = 1,ny |
---|
200 | do i = 1,nx |
---|
201 | n1 = n1 + 1 |
---|
202 | av_g%rAttr(n,n1) = array2(i,j) |
---|
203 | enddo |
---|
204 | enddo |
---|
205 | deallocate(array2) |
---|
206 | endif |
---|
207 | |
---|
208 | status = nf90_close(ncid) |
---|
209 | call check_status(status,subname,__FILE__,__LINE__) |
---|
210 | |
---|
211 | endif |
---|
212 | |
---|
213 | call mct_aVect_scatter(av_g,av,gsmap,master_task,mpicom) |
---|
214 | if (iam == master_task) then |
---|
215 | call mct_aVect_clean(av_g) |
---|
216 | endif |
---|
217 | |
---|
218 | call oasis_debug_exit(subname) |
---|
219 | |
---|
220 | end subroutine oasis_io_read_avfld |
---|
221 | |
---|
222 | !=============================================================================== |
---|
223 | |
---|
224 | !> Write an attribute keyword value pair to a field |
---|
225 | |
---|
226 | subroutine oasis_io_write_fldattr(filename,fldname,keyword,value) |
---|
227 | |
---|
228 | ! --------------------------------------- |
---|
229 | ! Designed to work with oasis3 write_grid . |
---|
230 | ! --------------------------------------- |
---|
231 | |
---|
232 | implicit none |
---|
233 | |
---|
234 | character(len=*), intent(in) :: filename !< file name |
---|
235 | character(len=*), intent(in) :: fldname !< field name |
---|
236 | character(len=*), intent(in) :: keyword !< attribute name |
---|
237 | character(len=*), intent(in) :: value !< attribute value |
---|
238 | |
---|
239 | !--- local --- |
---|
240 | integer(ip_i4_p) :: ncid,dimid,dimid2(2),varid ! cdf info |
---|
241 | integer(ip_i4_p) :: status ! error code |
---|
242 | integer(ip_i4_p) :: ind ! string index |
---|
243 | logical :: exists ! file existance |
---|
244 | character(len=ic_med) :: gridname ! grid name derived from fldname |
---|
245 | |
---|
246 | character(len=*),parameter :: subname = '(oasis_io_write_fldattr)' |
---|
247 | |
---|
248 | !------------------------------------------------------------------------------- |
---|
249 | ! |
---|
250 | !------------------------------------------------------------------------------- |
---|
251 | |
---|
252 | call oasis_debug_enter(subname) |
---|
253 | |
---|
254 | ! expects to run only on 1 pe. |
---|
255 | ! if (iam == master_task) then |
---|
256 | |
---|
257 | if (OASIS_debug >= 5) then |
---|
258 | write(nulprt,*) subname,' write ',trim(filename),' ',trim(fldname) |
---|
259 | write(nulprt,*) subname,' write ',trim(keyword),' ',trim(value) |
---|
260 | endif |
---|
261 | |
---|
262 | ! open file, must already exist |
---|
263 | inquire(file=trim(filename),exist=exists) |
---|
264 | if (exists) then |
---|
265 | status = nf90_open(filename,NF90_WRITE,ncid) |
---|
266 | call check_status(status,subname,__FILE__,__LINE__) |
---|
267 | status = nf90_redef(ncid) |
---|
268 | call check_status(status,subname,__FILE__,__LINE__) |
---|
269 | else |
---|
270 | write(nulprt,*) subname,estr,'file does not exist' |
---|
271 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
272 | endif |
---|
273 | |
---|
274 | ! get varid, must already exist |
---|
275 | status = nf90_inq_varid(ncid,trim(fldname),varid) |
---|
276 | if (status == nf90_noerr) then |
---|
277 | status = nf90_put_att(ncid,varid,trim(keyword),trim(value)) |
---|
278 | call check_status(status,subname,__FILE__,__LINE__) |
---|
279 | else |
---|
280 | write(nulprt,*) subname,estr,'variable does not exist ',trim(fldname) |
---|
281 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
282 | endif |
---|
283 | |
---|
284 | status = nf90_close(ncid) |
---|
285 | call check_status(status,subname,__FILE__,__LINE__) |
---|
286 | |
---|
287 | ! endif |
---|
288 | |
---|
289 | call oasis_debug_exit(subname) |
---|
290 | |
---|
291 | end subroutine oasis_io_write_fldattr |
---|
292 | |
---|
293 | !=============================================================================== |
---|
294 | |
---|
295 | !> Writes all fields from an attribute vector to a file |
---|
296 | |
---|
297 | subroutine oasis_io_write_avfile(rstfile,av,gsmap,mpicom,nx,ny,nampre) |
---|
298 | |
---|
299 | implicit none |
---|
300 | |
---|
301 | character(len=*), intent(in) :: rstfile !< filename |
---|
302 | type(mct_aVect) , intent(in) :: av !< avect |
---|
303 | type(mct_gsmap) , intent(in) :: gsmap !< gsmap decomp |
---|
304 | integer(ip_i4_p), intent(in) :: mpicom !< mpi comm |
---|
305 | integer(ip_i4_p), intent(in) :: nx !< 2d global nx size |
---|
306 | integer(ip_i4_p), intent(in) :: ny !< 2d global ny size |
---|
307 | character(len=*), intent(in),optional :: nampre !< field name prepend string on file |
---|
308 | |
---|
309 | !--- local --- |
---|
310 | integer(ip_i4_p) :: n,n1,i,j,fk,fk1 ! index |
---|
311 | integer(ip_i4_p) :: nxf,nyf ! field size on file |
---|
312 | type(mct_aVect) :: av_g ! avect global data |
---|
313 | type(mct_string) :: mstring ! mct char type |
---|
314 | character(ic_med) :: itemc ! string converted to char |
---|
315 | character(ic_med) :: lnampre ! local nampre |
---|
316 | character(ic_med) :: lstring ! local filename |
---|
317 | integer(ip_i4_p) :: master_task,iam,ierr ! mpi info |
---|
318 | integer(ip_i4_p) :: ncid,dimid,dimid2(2),varid ! netcdf info |
---|
319 | integer(ip_i4_p) :: dlen ! dimension length |
---|
320 | integer(ip_i4_p) :: iflag ! netcdf control |
---|
321 | integer(ip_i4_p) :: status ! error code |
---|
322 | logical :: exists ! file existance |
---|
323 | real(ip_double_p),allocatable :: array2(:,:) |
---|
324 | |
---|
325 | character(len=*),parameter :: subname = '(oasis_io_write_avfile)' |
---|
326 | |
---|
327 | !------------------------------------------------------------------------------- |
---|
328 | ! |
---|
329 | !------------------------------------------------------------------------------- |
---|
330 | |
---|
331 | if (mpicom == MPI_COMM_NULL) return |
---|
332 | |
---|
333 | call oasis_debug_enter(subname) |
---|
334 | |
---|
335 | ! empty filename, just return |
---|
336 | |
---|
337 | if (len_trim(rstfile) < 1) then |
---|
338 | call oasis_debug_exit(subname) |
---|
339 | return |
---|
340 | endif |
---|
341 | |
---|
342 | lnampre = "" |
---|
343 | if (present(nampre)) then |
---|
344 | lnampre = trim(nampre) |
---|
345 | endif |
---|
346 | |
---|
347 | master_task = 0 |
---|
348 | call MPI_COMM_RANK(mpicom,iam,ierr) |
---|
349 | |
---|
350 | call mct_aVect_gather(av,av_g,gsmap,master_task,mpicom) |
---|
351 | |
---|
352 | if (iam == master_task) then |
---|
353 | if (size(av_g%rAttr,dim=2) /= nx*ny) then |
---|
354 | write(nulprt,*) subname,estr,'av gsize nx ny mismatch in file :',& |
---|
355 | trim(lstring),SIZE(av_g%rAttr,dim=2),nx,ny |
---|
356 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
357 | endif |
---|
358 | |
---|
359 | inquire(file=trim(rstfile),exist=exists) |
---|
360 | if (exists) then |
---|
361 | status = nf90_open(trim(rstfile),NF90_WRITE,ncid) |
---|
362 | call check_status(status,subname,__FILE__,__LINE__) |
---|
363 | status = nf90_redef(ncid) |
---|
364 | call check_status(status,subname,__FILE__,__LINE__) |
---|
365 | else |
---|
366 | iflag = NF90_CLOBBER |
---|
367 | if (cdf_filetype == 'cdf2') iflag = ior(iflag,oasis_io_64bit_offset) |
---|
368 | if (cdf_filetype == 'cdf5') iflag = ior(iflag,oasis_io_64bit_data) |
---|
369 | status = nf90_create(trim(rstfile),iflag,ncid) |
---|
370 | call check_status(status,subname,__FILE__,__LINE__) |
---|
371 | endif |
---|
372 | |
---|
373 | do n = 1,mct_aVect_nRAttr(av_g) |
---|
374 | call mct_aVect_getRList(mstring,n,av_g) |
---|
375 | itemc = mct_string_toChar(mstring) |
---|
376 | itemc = trim(lnampre)//trim(itemc) |
---|
377 | call mct_string_clean(mstring) |
---|
378 | |
---|
379 | status = nf90_inq_dimid(ncid,trim(itemc)//'_nx',dimid2(1)) |
---|
380 | if (status /= NF90_NOERR) then |
---|
381 | status = nf90_def_dim(ncid,trim(itemc)//'_nx',nx,dimid2(1)) |
---|
382 | call check_status(status,subname,__FILE__,__LINE__) |
---|
383 | endif |
---|
384 | |
---|
385 | status = nf90_inq_dimid(ncid,trim(itemc)//'_ny',dimid2(2)) |
---|
386 | if (status /= NF90_NOERR) then |
---|
387 | status = nf90_def_dim(ncid,trim(itemc)//'_ny',ny,dimid2(2)) |
---|
388 | call check_status(status,subname,__FILE__,__LINE__) |
---|
389 | endif |
---|
390 | |
---|
391 | status = nf90_inquire_dimension(ncid,dimid2(1),len=dlen) |
---|
392 | call check_status(status,subname,__FILE__,__LINE__) |
---|
393 | if (dlen /= nx) then |
---|
394 | write(nulprt,*) subname,wstr,'dlen ne nx ',dlen,nx |
---|
395 | call oasis_flush(nulprt) |
---|
396 | ! call oasis_abort(file=__FILE__,line=__LINE__) |
---|
397 | endif |
---|
398 | |
---|
399 | status = nf90_inquire_dimension(ncid,dimid2(2),len=dlen) |
---|
400 | call check_status(status,subname,__FILE__,__LINE__) |
---|
401 | if (dlen /= ny) then |
---|
402 | write(nulprt,*) subname,wstr,'dlen ne ny ',dlen,ny |
---|
403 | call oasis_flush(nulprt) |
---|
404 | ! call oasis_abort(file=__FILE__,line=__LINE__) |
---|
405 | endif |
---|
406 | |
---|
407 | status = nf90_inq_varid(ncid,trim(itemc),varid) |
---|
408 | if (status /= NF90_NOERR) then |
---|
409 | status = nf90_def_var(ncid,trim(itemc),NF90_DOUBLE,dimid2,varid) |
---|
410 | call check_status(status,subname,__FILE__,__LINE__) |
---|
411 | endif |
---|
412 | |
---|
413 | enddo |
---|
414 | |
---|
415 | status = nf90_enddef(ncid) |
---|
416 | call check_status(status,subname,__FILE__,__LINE__) |
---|
417 | |
---|
418 | nxf = 0 |
---|
419 | nyf = 0 |
---|
420 | do n = 1,mct_aVect_nRAttr(av_g) |
---|
421 | call mct_aVect_getRList(mstring,n,av_g) |
---|
422 | itemc = mct_string_toChar(mstring) |
---|
423 | itemc = trim(lnampre)//trim(itemc) |
---|
424 | call mct_string_clean(mstring) |
---|
425 | status = nf90_inq_varid(ncid,trim(itemc),varid) |
---|
426 | call check_status(status,subname,__FILE__,__LINE__) |
---|
427 | if (n == 1) then |
---|
428 | status = nf90_inquire_variable(ncid,varid,ndims=dlen,dimids=dimid2) |
---|
429 | call check_status(status,subname,__FILE__,__LINE__) |
---|
430 | status = nf90_inquire_dimension(ncid,dimid2(1),len=nxf) |
---|
431 | call check_status(status,subname,__FILE__,__LINE__) |
---|
432 | status = nf90_inquire_dimension(ncid,dimid2(2),len=nyf) |
---|
433 | call check_status(status,subname,__FILE__,__LINE__) |
---|
434 | if (dlen /= 2 .or. nx*ny /= nxf*nyf) then |
---|
435 | write(nulprt,*) subname,estr,'ndims and size does not match on file' |
---|
436 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
437 | endif |
---|
438 | allocate(array2(nxf,nyf)) |
---|
439 | endif |
---|
440 | |
---|
441 | n1 = 0 |
---|
442 | do j = 1,nyf |
---|
443 | do i = 1,nxf |
---|
444 | n1 = n1 + 1 |
---|
445 | array2(i,j) = av_g%rAttr(n,n1) |
---|
446 | enddo |
---|
447 | enddo |
---|
448 | |
---|
449 | status = nf90_put_var(ncid,varid,array2) |
---|
450 | call check_status(status,subname,__FILE__,__LINE__) |
---|
451 | enddo |
---|
452 | deallocate(array2) |
---|
453 | call mct_aVect_clean(av_g) |
---|
454 | |
---|
455 | status = nf90_close(ncid) |
---|
456 | call check_status(status,subname,__FILE__,__LINE__) |
---|
457 | |
---|
458 | endif |
---|
459 | |
---|
460 | call oasis_debug_exit(subname) |
---|
461 | |
---|
462 | end subroutine oasis_io_write_avfile |
---|
463 | |
---|
464 | !=============================================================================== |
---|
465 | |
---|
466 | !> Reads all fields for an attribute vector from a file |
---|
467 | |
---|
468 | subroutine oasis_io_read_avfile(rstfile,av,gsmap,mpicom,abort,nampre,didread) |
---|
469 | |
---|
470 | implicit none |
---|
471 | |
---|
472 | character(len=*), intent(in) :: rstfile !< filename |
---|
473 | type(mct_aVect) , intent(inout) :: av !< avect |
---|
474 | type(mct_gsmap) , intent(in) :: gsmap !< gsmap decomp |
---|
475 | integer(ip_i4_p), intent(in) :: mpicom !< mpi comm |
---|
476 | logical , intent(in) ,optional :: abort !< abort on error flag, default is true |
---|
477 | character(len=*), intent(in) ,optional :: nampre !< name prepend string for fields on file |
---|
478 | logical , intent(out),optional :: didread !< flag indicating that read was successful |
---|
479 | |
---|
480 | !--- local --- |
---|
481 | integer(ip_i4_p) :: n,n1,i,j,fk,fk1 ! index |
---|
482 | integer(ip_i4_p) :: nx ! 2d global size nx |
---|
483 | integer(ip_i4_p) :: ny ! 2d global size ny |
---|
484 | type(mct_aVect) :: av_g ! avect global data |
---|
485 | type(mct_string) :: mstring ! mct char type |
---|
486 | character(ic_med) :: itemc ! string converted to char |
---|
487 | character(ic_med) :: lnampre ! local nampre |
---|
488 | character(ic_med) :: lstring ! local filename |
---|
489 | integer(ip_i4_p) :: master_task,iam,ierr ! mpi info |
---|
490 | integer(ip_i4_p) :: ncid,dimid,dimid2(2),varid ! netcdf info |
---|
491 | integer(ip_i4_p) :: dlen ! dimension length |
---|
492 | integer(ip_i4_p) :: status ! error code |
---|
493 | logical :: exists ! file existance |
---|
494 | logical :: labort ! local abort flag |
---|
495 | real(ip_double_p),allocatable :: array2(:,:) |
---|
496 | |
---|
497 | character(len=*),parameter :: subname = '(oasis_io_read_avfile)' |
---|
498 | |
---|
499 | !------------------------------------------------------------------------------- |
---|
500 | ! |
---|
501 | !------------------------------------------------------------------------------- |
---|
502 | |
---|
503 | if (mpicom == MPI_COMM_NULL) return |
---|
504 | |
---|
505 | call oasis_debug_enter(subname) |
---|
506 | |
---|
507 | if (present(didread)) didread = .false. |
---|
508 | |
---|
509 | ! empty filename, just return |
---|
510 | |
---|
511 | if (len_trim(rstfile) < 1) then |
---|
512 | call oasis_debug_exit(subname) |
---|
513 | return |
---|
514 | endif |
---|
515 | |
---|
516 | labort = .true. |
---|
517 | if (present(abort)) then |
---|
518 | labort = abort |
---|
519 | endif |
---|
520 | |
---|
521 | lnampre = "" |
---|
522 | if (present(nampre)) then |
---|
523 | lnampre = trim(nampre) |
---|
524 | endif |
---|
525 | |
---|
526 | master_task = 0 |
---|
527 | call MPI_COMM_RANK(mpicom,iam,ierr) |
---|
528 | |
---|
529 | call mct_aVect_gather(av,av_g,gsmap,master_task,mpicom) |
---|
530 | |
---|
531 | if (iam == master_task) then |
---|
532 | |
---|
533 | inquire(file=trim(rstfile),exist=exists) |
---|
534 | if (.not.exists) then |
---|
535 | if (labort) THEN |
---|
536 | write(nulprt,*) subname,estr,'file missing ',trim(rstfile) |
---|
537 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
538 | else |
---|
539 | write(nulprt,*) subname,wstr,'file missing ',trim(rstfile) |
---|
540 | call oasis_flush(nulprt) |
---|
541 | endif |
---|
542 | else |
---|
543 | status = nf90_open(trim(rstfile),NF90_NOWRITE,ncid) |
---|
544 | call check_status(status,subname,__FILE__,__LINE__) |
---|
545 | |
---|
546 | do n = 1,mct_aVect_nRAttr(av_g) |
---|
547 | call mct_aVect_getRList(mstring,n,av_g) |
---|
548 | itemc = mct_string_toChar(mstring) |
---|
549 | itemc = trim(lnampre)//trim(itemc) |
---|
550 | call mct_string_clean(mstring) |
---|
551 | |
---|
552 | status = nf90_inq_varid(ncid,trim(itemc),varid) |
---|
553 | if (status /= NF90_NOERR) then |
---|
554 | if (labort) THEN |
---|
555 | write(nulprt,*) subname,estr,'var missing on file = ',trim(itemc),':',trim(nf90_strerror(status)) |
---|
556 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
557 | ! else |
---|
558 | ! write(nulprt,*) subname,wstr,'var missing on file = ',trim(itemc),':',trim(nf90_strerror(status)) |
---|
559 | ! call oasis_flush(nulprt) |
---|
560 | endif |
---|
561 | |
---|
562 | else |
---|
563 | status = nf90_inquire_variable(ncid,varid,ndims=dlen,dimids=dimid2) |
---|
564 | call check_status(status,subname,__FILE__,__LINE__) |
---|
565 | if (dlen > 2) then |
---|
566 | write(nulprt,*) subname,estr,'variable ndims gt 2 on file ',trim(itemc),dlen |
---|
567 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
568 | endif |
---|
569 | status = nf90_inquire_dimension(ncid,dimid2(1),len=nx) |
---|
570 | call check_status(status,subname,__FILE__,__LINE__) |
---|
571 | ny = 1 |
---|
572 | if (dlen == 2) then |
---|
573 | status = nf90_inquire_dimension(ncid,dimid2(2),len=ny) |
---|
574 | call check_status(status,subname,__FILE__,__LINE__) |
---|
575 | endif |
---|
576 | |
---|
577 | if (size(av_g%rAttr,dim=2) /= nx*ny) then |
---|
578 | write(nulprt,*) subname,estr,'av gsize nx ny mismatch in file = ',& |
---|
579 | trim(rstfile),SIZE(av_g%rAttr,dim=2),nx,ny |
---|
580 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
581 | endif |
---|
582 | |
---|
583 | allocate(array2(nx,ny)) |
---|
584 | |
---|
585 | status = nf90_get_var(ncid,varid,array2) |
---|
586 | call check_status(status,subname,__FILE__,__LINE__) |
---|
587 | |
---|
588 | n1 = 0 |
---|
589 | do j = 1,ny |
---|
590 | do i = 1,nx |
---|
591 | n1 = n1 + 1 |
---|
592 | av_g%rAttr(n,n1) = array2(i,j) |
---|
593 | enddo |
---|
594 | enddo |
---|
595 | if (present(didread)) didread = .true. |
---|
596 | |
---|
597 | deallocate(array2) |
---|
598 | endif ! varid valid |
---|
599 | enddo |
---|
600 | |
---|
601 | status = nf90_close(ncid) |
---|
602 | call check_status(status,subname,__FILE__,__LINE__) |
---|
603 | |
---|
604 | endif ! file exists |
---|
605 | endif |
---|
606 | |
---|
607 | call mct_aVect_scatter(av_g,av,gsmap,master_task,mpicom) |
---|
608 | if (iam == master_task) then |
---|
609 | call mct_aVect_clean(av_g) |
---|
610 | endif |
---|
611 | |
---|
612 | call oasis_debug_exit(subname) |
---|
613 | |
---|
614 | end subroutine oasis_io_read_avfile |
---|
615 | |
---|
616 | !=============================================================================== |
---|
617 | |
---|
618 | !> Reads an integer or real field from a file into an array |
---|
619 | |
---|
620 | subroutine oasis_io_read_array(rstfile,mpicom,iarray,ivarname,rarray,rvarname,abort) |
---|
621 | |
---|
622 | implicit none |
---|
623 | |
---|
624 | character(len=*), intent(in) :: rstfile !< filename |
---|
625 | integer(ip_i4_p), intent(in) :: mpicom !< mpi comm |
---|
626 | integer(ip_i4_p), intent(inout),optional :: iarray(:) !< integer data on root |
---|
627 | character(len=*), intent(in),optional :: ivarname !< integer variable name on file |
---|
628 | real(ip_double_p),intent(inout),optional :: rarray(:) !< real data on root |
---|
629 | character(len=*), intent(in),optional :: rvarname !< real variable name on file |
---|
630 | logical , intent(in),optional :: abort !< abort on error flag, default is true |
---|
631 | |
---|
632 | !--- local --- |
---|
633 | integer(ip_i4_p) :: ncnt |
---|
634 | integer(ip_i4_p) :: master_task,iam,ierr ! mpi info |
---|
635 | integer(ip_i4_p) :: ncid,dimid,dimid1(1),varid ! netcdf info |
---|
636 | integer(ip_i4_p) :: dlen ! dimension length |
---|
637 | integer(ip_i4_p) :: status ! error code |
---|
638 | logical :: exists ! file existance |
---|
639 | logical :: labort ! local abort flag |
---|
640 | |
---|
641 | character(len=*),parameter :: subname = '(oasis_io_read_array)' |
---|
642 | |
---|
643 | !------------------------------------------------------------------------------- |
---|
644 | ! |
---|
645 | !------------------------------------------------------------------------------- |
---|
646 | |
---|
647 | if (mpicom == MPI_COMM_NULL) return |
---|
648 | |
---|
649 | call oasis_debug_enter(subname) |
---|
650 | |
---|
651 | ! empty filename, just return |
---|
652 | |
---|
653 | if (len_trim(rstfile) < 1) then |
---|
654 | call oasis_debug_exit(subname) |
---|
655 | return |
---|
656 | endif |
---|
657 | |
---|
658 | labort = .true. |
---|
659 | if (present(abort)) then |
---|
660 | labort = abort |
---|
661 | endif |
---|
662 | |
---|
663 | master_task = 0 |
---|
664 | call MPI_COMM_RANK(mpicom,iam,ierr) |
---|
665 | |
---|
666 | if (iam == master_task) then |
---|
667 | |
---|
668 | inquire(file=trim(rstfile),exist=exists) |
---|
669 | if (.not.exists) then |
---|
670 | if (labort) THEN |
---|
671 | write(nulprt,*) subname,estr,'file missing ',trim(rstfile) |
---|
672 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
673 | else |
---|
674 | write(nulprt,*) subname,wstr,'file missing ',trim(rstfile) |
---|
675 | call oasis_flush(nulprt) |
---|
676 | endif |
---|
677 | else |
---|
678 | status = nf90_open(trim(rstfile),NF90_NOWRITE,ncid) |
---|
679 | call check_status(status,subname,__FILE__,__LINE__) |
---|
680 | |
---|
681 | if (present(iarray)) then |
---|
682 | if (.not. present(ivarname)) then |
---|
683 | write(nulprt,*) subname,estr,'iarray must have ivarname set' |
---|
684 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
685 | endif |
---|
686 | |
---|
687 | ncnt = size(iarray) |
---|
688 | |
---|
689 | status = nf90_inq_varid(ncid,trim(ivarname),varid) |
---|
690 | if (status /= NF90_NOERR) then |
---|
691 | if (labort) THEN |
---|
692 | write(nulprt,*) subname,estr,'var missing on file = ',trim(ivarname),':',trim(nf90_strerror(status)) |
---|
693 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
694 | ! else |
---|
695 | ! write(nulprt,*) subname,wstr,'var missing on file = ',trim(ivarname),':',trim(nf90_strerror(status)) |
---|
696 | ! call oasis_flush(nulprt) |
---|
697 | endif |
---|
698 | else |
---|
699 | status = nf90_inquire_variable(ncid,varid,ndims=dlen,dimids=dimid1) |
---|
700 | call check_status(status,subname,__FILE__,__LINE__) |
---|
701 | if (dlen /= 1) then |
---|
702 | write(nulprt,*) subname,estr,'variable ndims ne 1 ',trim(ivarname),dlen |
---|
703 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
704 | endif |
---|
705 | status = nf90_inquire_dimension(ncid,dimid1(1),len=dlen) |
---|
706 | call check_status(status,subname,__FILE__,__LINE__) |
---|
707 | |
---|
708 | if (ncnt /= dlen) then |
---|
709 | write(nulprt,*) subname,estr,'iarray ncnt dlen mismatch ',ncnt,dlen |
---|
710 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
711 | endif |
---|
712 | |
---|
713 | status = nf90_get_var(ncid,varid,iarray) |
---|
714 | call check_status(status,subname,__FILE__,__LINE__) |
---|
715 | endif |
---|
716 | endif |
---|
717 | |
---|
718 | if (present(rarray)) then |
---|
719 | if (.not. present(rvarname)) then |
---|
720 | write(nulprt,*) subname,estr,'rarray must have rvarname set' |
---|
721 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
722 | endif |
---|
723 | |
---|
724 | ncnt = size(rarray) |
---|
725 | |
---|
726 | status = nf90_inq_varid(ncid,trim(rvarname),varid) |
---|
727 | if (status /= NF90_NOERR) then |
---|
728 | if (labort) THEN |
---|
729 | write(nulprt,*) subname,estr,'var missing on file = ',trim(rvarname),':',trim(nf90_strerror(status)) |
---|
730 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
731 | ! else |
---|
732 | ! write(nulprt,*) subname,wstr,'var missing on file = ',trim(rvarname),':',trim(nf90_strerror(status)) |
---|
733 | ! call oasis_flush(nulprt) |
---|
734 | endif |
---|
735 | else |
---|
736 | status = nf90_inquire_variable(ncid,varid,ndims=dlen,dimids=dimid1) |
---|
737 | call check_status(status,subname,__FILE__,__LINE__) |
---|
738 | if (dlen /= 1) then |
---|
739 | write(nulprt,*) subname,estr,'variable ndims ne 1 ',trim(rvarname),dlen |
---|
740 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
741 | endif |
---|
742 | status = nf90_inquire_dimension(ncid,dimid1(1),len=dlen) |
---|
743 | call check_status(status,subname,__FILE__,__LINE__) |
---|
744 | |
---|
745 | if (ncnt /= dlen) then |
---|
746 | write(nulprt,*) subname,estr,'rarray ncnt dlen mismatch ',ncnt,dlen |
---|
747 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
748 | endif |
---|
749 | |
---|
750 | status = nf90_get_var(ncid,varid,rarray) |
---|
751 | call check_status(status,subname,__FILE__,__LINE__) |
---|
752 | endif |
---|
753 | endif |
---|
754 | |
---|
755 | status = nf90_close(ncid) |
---|
756 | call check_status(status,subname,__FILE__,__LINE__) |
---|
757 | |
---|
758 | endif |
---|
759 | endif |
---|
760 | |
---|
761 | if (present(iarray)) then |
---|
762 | call oasis_mpi_bcast(iarray,mpicom,trim(subname)//':iarray',master_task) |
---|
763 | endif |
---|
764 | |
---|
765 | if (present(rarray)) then |
---|
766 | call oasis_mpi_bcast(rarray,mpicom,trim(subname)//':rarray',master_task) |
---|
767 | endif |
---|
768 | |
---|
769 | call oasis_debug_exit(subname) |
---|
770 | |
---|
771 | end subroutine oasis_io_read_array |
---|
772 | |
---|
773 | !=============================================================================== |
---|
774 | |
---|
775 | !> Writes a real or integer array to a file |
---|
776 | |
---|
777 | subroutine oasis_io_write_array(rstfile,mpicom,iarray,ivarname,rarray,rvarname) |
---|
778 | |
---|
779 | implicit none |
---|
780 | |
---|
781 | character(len=*), intent(in) :: rstfile !< filename |
---|
782 | integer(ip_i4_p), intent(in) :: mpicom !< mpi comm |
---|
783 | integer(ip_i4_p), intent(in),optional :: iarray(:) !< integer data on root |
---|
784 | character(len=*), intent(in),optional :: ivarname !< integer variable name on file |
---|
785 | real(ip_double_p),intent(in),optional :: rarray(:) !< real data on root |
---|
786 | character(len=*), intent(in),optional :: rvarname !< real variable name on file |
---|
787 | |
---|
788 | !--- local --- |
---|
789 | integer(ip_i4_p) :: ncnt |
---|
790 | integer(ip_i4_p) :: master_task,iam,ierr ! mpi info |
---|
791 | integer(ip_i4_p) :: ncid,dimid,dimid1(1),varid ! netcdf info |
---|
792 | integer(ip_i4_p) :: dlen ! dimension length |
---|
793 | integer(ip_i4_p) :: iflag ! netcdf control |
---|
794 | integer(ip_i4_p) :: status ! error code |
---|
795 | logical :: exists ! file existance |
---|
796 | |
---|
797 | character(len=*),parameter :: subname = '(oasis_io_write_array)' |
---|
798 | |
---|
799 | !------------------------------------------------------------------------------- |
---|
800 | ! |
---|
801 | !------------------------------------------------------------------------------- |
---|
802 | |
---|
803 | if (mpicom == MPI_COMM_NULL) return |
---|
804 | |
---|
805 | call oasis_debug_enter(subname) |
---|
806 | |
---|
807 | ! empty filename, just return |
---|
808 | |
---|
809 | if (len_trim(rstfile) < 1) then |
---|
810 | call oasis_debug_exit(subname) |
---|
811 | return |
---|
812 | endif |
---|
813 | |
---|
814 | master_task = 0 |
---|
815 | call MPI_COMM_RANK(mpicom,iam,ierr) |
---|
816 | |
---|
817 | if (iam == master_task) then |
---|
818 | |
---|
819 | inquire(file=trim(rstfile),exist=exists) |
---|
820 | if (exists) then |
---|
821 | status = nf90_open(trim(rstfile),NF90_WRITE,ncid) |
---|
822 | call check_status(status,subname,__FILE__,__LINE__) |
---|
823 | status = nf90_redef(ncid) |
---|
824 | call check_status(status,subname,__FILE__,__LINE__) |
---|
825 | else |
---|
826 | iflag = NF90_CLOBBER |
---|
827 | if (cdf_filetype == 'cdf2') iflag = ior(iflag,oasis_io_64bit_offset) |
---|
828 | if (cdf_filetype == 'cdf5') iflag = ior(iflag,oasis_io_64bit_data) |
---|
829 | status = nf90_create(trim(rstfile),iflag,ncid) |
---|
830 | call check_status(status,subname,__FILE__,__LINE__) |
---|
831 | endif |
---|
832 | |
---|
833 | if (present(iarray)) then |
---|
834 | if (.not. present(ivarname)) then |
---|
835 | write(nulprt,*) subname,estr,'iarray must have ivarname set' |
---|
836 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
837 | endif |
---|
838 | |
---|
839 | ncnt = size(iarray) |
---|
840 | |
---|
841 | status = nf90_inq_dimid(ncid,trim(ivarname)//'_ncnt',dimid1(1)) |
---|
842 | if (status /= NF90_NOERR) then |
---|
843 | status = nf90_def_dim(ncid,trim(ivarname)//'_ncnt',ncnt,dimid1(1)) |
---|
844 | call check_status(status,subname,__FILE__,__LINE__) |
---|
845 | endif |
---|
846 | |
---|
847 | status = nf90_inquire_dimension(ncid,dimid1(1),len=dlen) |
---|
848 | call check_status(status,subname,__FILE__,__LINE__) |
---|
849 | if (dlen /= ncnt) then |
---|
850 | write(nulprt,*) subname,estr,'iarray dlen ne ncnt ',dlen,ncnt |
---|
851 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
852 | endif |
---|
853 | |
---|
854 | status = nf90_inq_varid(ncid,trim(ivarname),varid) |
---|
855 | if (status /= NF90_NOERR) then |
---|
856 | status = nf90_def_var(ncid,trim(ivarname),NF90_INT,dimid1,varid) |
---|
857 | call check_status(status,subname,__FILE__,__LINE__) |
---|
858 | endif |
---|
859 | endif |
---|
860 | |
---|
861 | if (present(rarray)) then |
---|
862 | if (.not. present(rvarname)) then |
---|
863 | write(nulprt,*) subname,estr,'rarray must have rvarname set' |
---|
864 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
865 | endif |
---|
866 | |
---|
867 | ncnt = size(rarray) |
---|
868 | |
---|
869 | status = nf90_inq_dimid(ncid,trim(rvarname)//'_ncnt',dimid1(1)) |
---|
870 | if (status /= NF90_NOERR) then |
---|
871 | status = nf90_def_dim(ncid,trim(rvarname)//'_ncnt',ncnt,dimid1(1)) |
---|
872 | call check_status(status,subname,__FILE__,__LINE__) |
---|
873 | endif |
---|
874 | |
---|
875 | status = nf90_inquire_dimension(ncid,dimid1(1),len=dlen) |
---|
876 | call check_status(status,subname,__FILE__,__LINE__) |
---|
877 | if (dlen /= ncnt) then |
---|
878 | write(nulprt,*) subname,estr,'rarray dlen ne ncnt ',dlen,ncnt |
---|
879 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
880 | endif |
---|
881 | |
---|
882 | status = nf90_inq_varid(ncid,trim(rvarname),varid) |
---|
883 | if (status /= NF90_NOERR) then |
---|
884 | status = nf90_def_var(ncid,trim(rvarname),NF90_DOUBLE,dimid1,varid) |
---|
885 | call check_status(status,subname,__FILE__,__LINE__) |
---|
886 | endif |
---|
887 | endif |
---|
888 | |
---|
889 | status = nf90_enddef(ncid) |
---|
890 | call check_status(status,subname,__FILE__,__LINE__) |
---|
891 | |
---|
892 | if (present(iarray)) then |
---|
893 | status = nf90_inq_varid(ncid,trim(ivarname),varid) |
---|
894 | call check_status(status,subname,__FILE__,__LINE__) |
---|
895 | status = nf90_put_var(ncid,varid,iarray) |
---|
896 | call check_status(status,subname,__FILE__,__LINE__) |
---|
897 | endif |
---|
898 | |
---|
899 | if (present(rarray)) then |
---|
900 | status = nf90_inq_varid(ncid,trim(rvarname),varid) |
---|
901 | call check_status(status,subname,__FILE__,__LINE__) |
---|
902 | status = nf90_put_var(ncid,varid,rarray) |
---|
903 | call check_status(status,subname,__FILE__,__LINE__) |
---|
904 | endif |
---|
905 | |
---|
906 | status = nf90_close(ncid) |
---|
907 | call check_status(status,subname,__FILE__,__LINE__) |
---|
908 | |
---|
909 | endif |
---|
910 | |
---|
911 | call oasis_debug_exit(subname) |
---|
912 | |
---|
913 | end subroutine oasis_io_write_array |
---|
914 | |
---|
915 | !=============================================================================== |
---|
916 | |
---|
917 | !> Write each field in an attribute vector to an individual files |
---|
918 | |
---|
919 | subroutine oasis_io_write_avfbf(av,gsmap,mpicom,nx,ny,msec,f_string,filename) |
---|
920 | |
---|
921 | ! --------------------------------------- |
---|
922 | ! This works only for a single av to a file |
---|
923 | ! Optionally can specify time info, and filename info |
---|
924 | ! --------------------------------------- |
---|
925 | |
---|
926 | implicit none |
---|
927 | |
---|
928 | type(mct_aVect) , intent(in) :: av !< avect |
---|
929 | type(mct_gsmap) , intent(in) :: gsmap !< gsmap decomp |
---|
930 | integer(ip_i4_p), intent(in) :: mpicom !< mpi comm |
---|
931 | integer(ip_i4_p), intent(in) :: nx !< 2d global nx size |
---|
932 | integer(ip_i4_p), intent(in) :: ny !< 2d global ny size |
---|
933 | integer(ip_i4_p), intent(in),optional :: msec !< optional time info in seconds |
---|
934 | character(len=*), intent(in),optional :: f_string !< optional f_string to append to filename |
---|
935 | character(len=*), intent(in),optional :: filename !< optional output filename |
---|
936 | |
---|
937 | !--- local --- |
---|
938 | integer(ip_i4_p) :: n,n1,i,j,fk,fk1 ! index |
---|
939 | type(mct_aVect) :: av_g ! avect global data |
---|
940 | type(mct_string) :: mstring ! mct char type |
---|
941 | character(ic_med) :: itemc ! f_string converted to char |
---|
942 | character(ic_med) :: lfn ! local filename |
---|
943 | character(ic_med) :: lstring ! local filename |
---|
944 | integer(ip_i4_p) :: master_task,iam,ierr ! mpi info |
---|
945 | integer(ip_i4_p) :: ncid,dimid,dimid3(3),varid ! netcdf info |
---|
946 | integer(ip_i4_p) :: start3(3),count3(3) ! netcdf info |
---|
947 | integer(ip_i4_p) :: start1(1),count1(1) ! netcdf info |
---|
948 | integer(ip_i4_p) :: lmsec(1) ! local msec value |
---|
949 | integer(ip_i4_p) :: dlen ! dimension length |
---|
950 | integer(ip_i4_p) :: iflag ! netcdf control |
---|
951 | integer(ip_i4_p) :: status ! error code |
---|
952 | logical :: exists ! file existance |
---|
953 | logical :: newfile ! to create a new history file |
---|
954 | logical :: whead,wdata ! for writing restart/history cdf files |
---|
955 | integer(ip_i4_p) ,allocatable :: time(:) |
---|
956 | real(ip_double_p),allocatable :: array3(:,:,:) |
---|
957 | real(ip_double_p) :: tbnds(2) |
---|
958 | |
---|
959 | character(len=*),parameter :: subname = '(oasis_io_write_avfbf)' |
---|
960 | |
---|
961 | !------------------------------------------------------------------------------- |
---|
962 | ! |
---|
963 | !------------------------------------------------------------------------------- |
---|
964 | |
---|
965 | if (mpicom == MPI_COMM_NULL) return |
---|
966 | |
---|
967 | call oasis_debug_enter(subname) |
---|
968 | |
---|
969 | lmsec = 0 |
---|
970 | if (present(msec)) then |
---|
971 | lmsec = msec |
---|
972 | endif |
---|
973 | |
---|
974 | lstring = " " |
---|
975 | if (present(f_string)) then |
---|
976 | lstring = trim(f_string) |
---|
977 | endif |
---|
978 | |
---|
979 | master_task = 0 |
---|
980 | call MPI_COMM_RANK(mpicom,iam,ierr) |
---|
981 | |
---|
982 | #if (PIO_DEFINED) |
---|
983 | ! tcraig, not working as of Oct 2011 |
---|
984 | call oasis_ioshr_wopen(lfn,clobber=.true.,cdf64=.true.) |
---|
985 | |
---|
986 | do fk = fk1,2 |
---|
987 | if (fk == 1) then |
---|
988 | whead = .true. |
---|
989 | wdata = .false. |
---|
990 | elseif (fk == 2) then |
---|
991 | whead = .false. |
---|
992 | wdata = .true. |
---|
993 | else |
---|
994 | write(nulprt,*) subname,estr,'fk illegal' |
---|
995 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
996 | end if |
---|
997 | |
---|
998 | call oasis_ioshr_write(lfn,& |
---|
999 | time_units='seconds',time_cal='seconds',time_val=real(lmsec,ip_double_p),& |
---|
1000 | nt=1,whead=whead,wdata=wdata) |
---|
1001 | |
---|
1002 | call oasis_ioshr_write(lfn, gsmap, av, 'pout', & |
---|
1003 | whead=whead,wdata=wdata,nx=nx,ny=ny,nt=1, & |
---|
1004 | use_float=.false.) |
---|
1005 | |
---|
1006 | if (fk == 1) call oasis_ioshr_enddef(lfn) |
---|
1007 | enddo |
---|
1008 | |
---|
1009 | call oasis_ioshr_close(lfn) |
---|
1010 | #else |
---|
1011 | |
---|
1012 | call mct_aVect_gather(av,av_g,gsmap,master_task,mpicom) |
---|
1013 | if (iam == master_task) then |
---|
1014 | if (size(av_g%rAttr,dim=2) /= nx*ny) then |
---|
1015 | write(nulprt,*) subname,estr,'av gsize nx ny mismatch in file :',& |
---|
1016 | trim(filename),SIZE(av_g%rAttr,dim=2),nx,ny |
---|
1017 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
1018 | endif |
---|
1019 | |
---|
1020 | allocate(array3(nx,ny,1)) |
---|
1021 | do n = 1,mct_aVect_nRAttr(av_g) |
---|
1022 | call mct_aVect_getRList(mstring,n,av_g) |
---|
1023 | itemc = mct_string_toChar(mstring) |
---|
1024 | call mct_string_clean(mstring) |
---|
1025 | if (present(filename)) then |
---|
1026 | lfn = trim(filename) |
---|
1027 | else |
---|
1028 | lfn = trim(itemc)//trim(lstring)//'.nc' |
---|
1029 | endif |
---|
1030 | n1 = 0 |
---|
1031 | do j = 1,ny |
---|
1032 | do i = 1,nx |
---|
1033 | n1 = n1 + 1 |
---|
1034 | array3(i,j,1) = av_g%rAttr(n,n1) |
---|
1035 | enddo |
---|
1036 | enddo |
---|
1037 | |
---|
1038 | start1 = 1 |
---|
1039 | count1 = 1 |
---|
1040 | start3 = 1 |
---|
1041 | count3(1) = nx |
---|
1042 | count3(2) = ny |
---|
1043 | count3(3) = 1 |
---|
1044 | |
---|
1045 | newfile = .true. |
---|
1046 | inquire(file=trim(lfn),exist=exists) |
---|
1047 | if (exists) then |
---|
1048 | newfile = .false. |
---|
1049 | status = nf90_open(lfn,NF90_WRITE,ncid) |
---|
1050 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1051 | status = nf90_inq_dimid(ncid,'time',dimid) |
---|
1052 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1053 | status = nf90_inquire_dimension(ncid,dimid,len=dlen) |
---|
1054 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1055 | allocate(time(dlen)) |
---|
1056 | status = nf90_inq_varid(ncid,'time',varid) |
---|
1057 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1058 | status = nf90_get_var(ncid,varid,time) |
---|
1059 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1060 | |
---|
1061 | !--- check whether the current time is less than the current file time axis |
---|
1062 | !--- if so, then assume this is an old file and clobber it and start new |
---|
1063 | do i = 1,dlen |
---|
1064 | if (time(i) >= lmsec(1)) newfile=.true. |
---|
1065 | enddo |
---|
1066 | deallocate(time) |
---|
1067 | if (.not.newfile) then |
---|
1068 | start1(1) = dlen + 1 |
---|
1069 | start3(3) = start1(1) |
---|
1070 | endif |
---|
1071 | endif |
---|
1072 | |
---|
1073 | if (newfile) then |
---|
1074 | iflag = NF90_CLOBBER |
---|
1075 | if (cdf_filetype == 'cdf2') iflag = ior(iflag,oasis_io_64bit_offset) |
---|
1076 | if (cdf_filetype == 'cdf5') iflag = ior(iflag,oasis_io_64bit_data) |
---|
1077 | status = nf90_create(lfn,iflag,ncid) |
---|
1078 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1079 | status = nf90_def_dim(ncid,'nx',nx,dimid3(1)) |
---|
1080 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1081 | status = nf90_def_dim(ncid,'ny',ny,dimid3(2)) |
---|
1082 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1083 | status = nf90_def_dim(ncid,'time',NF90_UNLIMITED,dimid) |
---|
1084 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1085 | dimid3(3) = dimid |
---|
1086 | status = nf90_def_var(ncid,'time',NF90_INT,dimid,varid) |
---|
1087 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1088 | status = nf90_def_var(ncid,trim(itemc),NF90_DOUBLE,dimid3,varid) |
---|
1089 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1090 | status = nf90_enddef(ncid) |
---|
1091 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1092 | endif |
---|
1093 | |
---|
1094 | status = nf90_inq_varid(ncid,'time',varid) |
---|
1095 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1096 | status = nf90_put_var(ncid,varid,lmsec,start1,count1) |
---|
1097 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1098 | status = nf90_inq_varid(ncid,trim(itemc),varid) |
---|
1099 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1100 | status = nf90_put_var(ncid,varid,array3,start3,count3) |
---|
1101 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1102 | status = nf90_close(ncid) |
---|
1103 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1104 | enddo |
---|
1105 | deallocate(array3) |
---|
1106 | call mct_aVect_clean(av_g) |
---|
1107 | endif |
---|
1108 | |
---|
1109 | |
---|
1110 | #endif |
---|
1111 | |
---|
1112 | call oasis_debug_exit(subname) |
---|
1113 | |
---|
1114 | end subroutine oasis_io_write_avfbf |
---|
1115 | |
---|
1116 | !=============================================================================== |
---|
1117 | |
---|
1118 | !> Read each field in an attribute vector from individual files |
---|
1119 | |
---|
1120 | subroutine oasis_io_read_avfbf(av,gsmap,mpicom,msec,f_string,filename) |
---|
1121 | |
---|
1122 | ! --------------------------------------- |
---|
1123 | ! This works only for a single av from a file |
---|
1124 | ! Optionally can specify time info, and filename info |
---|
1125 | ! --------------------------------------- |
---|
1126 | |
---|
1127 | implicit none |
---|
1128 | |
---|
1129 | type(mct_aVect) , intent(inout) :: av !< avect |
---|
1130 | type(mct_gsmap) , intent(in) :: gsmap !< gsmap decomp |
---|
1131 | integer(ip_i4_p), intent(in) :: mpicom !< mpi comm |
---|
1132 | integer(ip_i4_p), intent(in),optional :: msec !< optional time info in seconds |
---|
1133 | character(len=*), intent(in),optional :: f_string !< optional f_string to append to filename |
---|
1134 | character(len=*), intent(in),optional :: filename !< optional input filename |
---|
1135 | |
---|
1136 | !--- local --- |
---|
1137 | integer(ip_i4_p) :: n,n1,i,j,fk,fk1 ! index |
---|
1138 | integer(ip_i4_p) :: nx,ny ! grid size from file |
---|
1139 | type(mct_aVect) :: av_g ! avect global data |
---|
1140 | type(mct_string) :: mstring ! mct char type |
---|
1141 | character(ic_med) :: itemc ! f_string converted to char |
---|
1142 | character(ic_med) :: lfn ! local filename |
---|
1143 | character(ic_med) :: lstring ! local filename |
---|
1144 | integer(ip_i4_p) :: master_task,iam,ierr ! mpi info |
---|
1145 | integer(ip_i4_p) :: ncid,dimid,dimid3(3),varid ! netcdf info |
---|
1146 | integer(ip_i4_p) :: start3(3),count3(3) ! netcdf info |
---|
1147 | integer(ip_i4_p) :: lmsec(1) ! local msec value |
---|
1148 | integer(ip_i4_p) :: dlen ! dimension length |
---|
1149 | integer(ip_i4_p) :: status ! error code |
---|
1150 | logical :: exists ! file existance |
---|
1151 | logical :: whead,wdata ! for writing restart/history cdf files |
---|
1152 | real(ip_double_p),allocatable :: array3(:,:,:) |
---|
1153 | integer(ip_i4_p) ,allocatable :: time(:) |
---|
1154 | real(ip_double_p) :: tbnds(2) |
---|
1155 | |
---|
1156 | character(len=*),parameter :: subname = '(oasis_io_read_avfbf)' |
---|
1157 | |
---|
1158 | !------------------------------------------------------------------------------- |
---|
1159 | ! |
---|
1160 | !------------------------------------------------------------------------------- |
---|
1161 | |
---|
1162 | if (mpicom == MPI_COMM_NULL) return |
---|
1163 | |
---|
1164 | call oasis_debug_enter(subname) |
---|
1165 | |
---|
1166 | lmsec = 0 |
---|
1167 | if (present(msec)) then |
---|
1168 | lmsec = msec |
---|
1169 | endif |
---|
1170 | |
---|
1171 | lstring = " " |
---|
1172 | if (present(f_string)) then |
---|
1173 | lstring = trim(f_string) |
---|
1174 | endif |
---|
1175 | |
---|
1176 | master_task = 0 |
---|
1177 | call MPI_COMM_RANK(mpicom,iam,ierr) |
---|
1178 | |
---|
1179 | call mct_aVect_gather(av,av_g,gsmap,master_task,mpicom) |
---|
1180 | if (iam == master_task) then |
---|
1181 | do n = 1,mct_aVect_nRAttr(av_g) |
---|
1182 | call mct_aVect_getRList(mstring,n,av_g) |
---|
1183 | itemc = mct_string_toChar(mstring) |
---|
1184 | call mct_string_clean(mstring) |
---|
1185 | if (present(filename)) then |
---|
1186 | lfn = trim(filename) |
---|
1187 | else |
---|
1188 | lfn = trim(itemc)//trim(lstring)//'.nc' |
---|
1189 | endif |
---|
1190 | |
---|
1191 | inquire(file=trim(lfn),exist=exists) |
---|
1192 | if (.not.exists) then |
---|
1193 | write(nulprt,*) subname,estr,'file not found ',trim(lfn) |
---|
1194 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
1195 | endif |
---|
1196 | |
---|
1197 | status = nf90_open(lfn,NF90_NOWRITE,ncid) |
---|
1198 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1199 | |
---|
1200 | status = nf90_inq_dimid(ncid,'time',dimid) |
---|
1201 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1202 | status = nf90_inquire_dimension(ncid,dimid,len=dlen) |
---|
1203 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1204 | allocate(time(dlen)) |
---|
1205 | status = nf90_inq_varid(ncid,'time',varid) |
---|
1206 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1207 | status = nf90_get_var(ncid,varid,time) |
---|
1208 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1209 | n1 = 0 |
---|
1210 | do j = 1,dlen |
---|
1211 | if (time(j) == lmsec(1)) n1 = j |
---|
1212 | enddo |
---|
1213 | deallocate(time) |
---|
1214 | if (n1 < 1) then |
---|
1215 | write(nulprt,*) subname,estr,'time not found on file ',trim(lfn),lmsec |
---|
1216 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
1217 | endif |
---|
1218 | |
---|
1219 | status = nf90_inq_varid(ncid,trim(itemc),varid) |
---|
1220 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1221 | status = nf90_inquire_variable(ncid,varid,dimids=dimid3) |
---|
1222 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1223 | status = nf90_inquire_dimension(ncid,dimid3(1),len=nx) |
---|
1224 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1225 | status = nf90_inquire_dimension(ncid,dimid3(2),len=ny) |
---|
1226 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1227 | |
---|
1228 | if (size(av_g%rAttr,dim=2) /= nx*ny) then |
---|
1229 | write(nulprt,*) subname,estr,'av gsize nx ny mismatch in file :',& |
---|
1230 | trim(filename),SIZE(av_g%rAttr,dim=2),nx,ny |
---|
1231 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
1232 | endif |
---|
1233 | |
---|
1234 | start3 = 1 |
---|
1235 | count3(1) = nx |
---|
1236 | count3(2) = ny |
---|
1237 | count3(3) = 1 |
---|
1238 | start3(3) = n1 |
---|
1239 | allocate(array3(nx,ny,1)) |
---|
1240 | |
---|
1241 | status = nf90_get_var(ncid,varid,array3,start3,count3) |
---|
1242 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1243 | status = nf90_close(ncid) |
---|
1244 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1245 | |
---|
1246 | n1 = 0 |
---|
1247 | do j = 1,ny |
---|
1248 | do i = 1,nx |
---|
1249 | n1 = n1 + 1 |
---|
1250 | av_g%rAttr(n,n1) = array3(i,j,1) |
---|
1251 | enddo |
---|
1252 | enddo |
---|
1253 | |
---|
1254 | deallocate(array3) |
---|
1255 | |
---|
1256 | enddo |
---|
1257 | endif |
---|
1258 | |
---|
1259 | call mct_aVect_scatter(av_g,av,gsmap,master_task,mpicom) |
---|
1260 | if (iam == master_task) then |
---|
1261 | call mct_aVect_clean(av_g) |
---|
1262 | endif |
---|
1263 | |
---|
1264 | call oasis_debug_exit(subname) |
---|
1265 | |
---|
1266 | end subroutine oasis_io_read_avfbf |
---|
1267 | |
---|
1268 | !=============================================================================== |
---|
1269 | |
---|
1270 | !> Read a field on the root task from a file into an array. |
---|
1271 | |
---|
1272 | subroutine oasis_io_read_field_fromroot(filename,fldname,ifld2,fld2,fld3,nx,ny,nz) |
---|
1273 | |
---|
1274 | ! --------------------------------------- |
---|
1275 | ! Supports 2d integer or real and 3d real arrays. |
---|
1276 | ! --------------------------------------- |
---|
1277 | |
---|
1278 | implicit none |
---|
1279 | |
---|
1280 | character(len=*) , intent(in) :: filename !< file name |
---|
1281 | character(len=*) , intent(in) :: fldname !< field name |
---|
1282 | integer(ip_i4_p) , intent(inout),optional :: ifld2(:,:) !< 2d integer array |
---|
1283 | real(ip_realwp_p), intent(inout),optional :: fld2(:,:) !< 2d real array |
---|
1284 | real(ip_realwp_p), intent(inout),optional :: fld3(:,:,:) !< 3d real array |
---|
1285 | integer(ip_i4_p) , intent(inout),optional :: nx !< global nx size |
---|
1286 | integer(ip_i4_p) , intent(inout),optional :: ny !< global ny size |
---|
1287 | integer(ip_i4_p) , intent(inout),optional :: nz !< global nz size |
---|
1288 | |
---|
1289 | !--- local --- |
---|
1290 | integer(ip_i4_p) :: ncid,varid ! cdf info |
---|
1291 | integer(ip_i4_p) :: n,ndims,xtype |
---|
1292 | integer(ip_i4_p),allocatable :: dimid(:),nd(:) |
---|
1293 | integer(ip_i4_p) :: status ! error code |
---|
1294 | integer(ip_i4_p) :: ind ! string index |
---|
1295 | logical :: exists ! file existance |
---|
1296 | character(len=ic_med) :: gridname ! grid name derived from fldname |
---|
1297 | |
---|
1298 | character(len=*),parameter :: subname = '(oasis_io_read_field_fromroot)' |
---|
1299 | |
---|
1300 | !------------------------------------------------------------------------------- |
---|
1301 | ! |
---|
1302 | !------------------------------------------------------------------------------- |
---|
1303 | |
---|
1304 | call oasis_debug_enter(subname) |
---|
1305 | |
---|
1306 | ! expects to run only on 1 pe. |
---|
1307 | ! if (iam == master_task) then |
---|
1308 | |
---|
1309 | ! EM Since parallel version, must be run only on 1 pe. |
---|
1310 | if ( mpi_rank_local == 0 ) then |
---|
1311 | |
---|
1312 | if (OASIS_debug >= 5) then |
---|
1313 | write(nulprt,*) subname,' read ',trim(filename),' ',trim(fldname) |
---|
1314 | endif |
---|
1315 | |
---|
1316 | inquire(file=trim(filename),exist=exists) |
---|
1317 | if (exists) then |
---|
1318 | status = nf90_open(filename,NF90_NOWRITE,ncid) |
---|
1319 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1320 | else |
---|
1321 | write(nulprt,*) subname,estr,'in filename ',trim(filename) |
---|
1322 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
1323 | endif |
---|
1324 | |
---|
1325 | status = nf90_inq_varid(ncid,trim(fldname),varid) |
---|
1326 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1327 | |
---|
1328 | status = nf90_inquire_variable(ncid,varid,ndims=ndims,xtype=xtype) |
---|
1329 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1330 | |
---|
1331 | allocate(dimid(ndims),nd(ndims)) |
---|
1332 | |
---|
1333 | status = nf90_inquire_variable(ncid,varid,dimids=dimid) |
---|
1334 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1335 | do n = 1,ndims |
---|
1336 | status = nf90_inquire_dimension(ncid,dimid(n),len=nd(n)) |
---|
1337 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1338 | enddo |
---|
1339 | |
---|
1340 | if (present(ifld2) .or. present(fld2) .or. present(fld3)) then |
---|
1341 | if (xtype == NF90_INT .and. ndims == 2 .and. present(ifld2)) then |
---|
1342 | status = nf90_get_var(ncid,varid,ifld2) |
---|
1343 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1344 | elseif (xtype /= NF90_INT .and. ndims == 2 .and. present(fld2)) then |
---|
1345 | status = nf90_get_var(ncid,varid,fld2) |
---|
1346 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1347 | elseif (xtype /= NF90_INT .and. ndims == 3 .and. present(fld3)) then |
---|
1348 | status = nf90_get_var(ncid,varid,fld3) |
---|
1349 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1350 | else |
---|
1351 | write(nulprt,*) subname,estr,'mismatch in field and data' |
---|
1352 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
1353 | endif |
---|
1354 | endif |
---|
1355 | |
---|
1356 | status = nf90_close(ncid) |
---|
1357 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1358 | |
---|
1359 | if (present(nx)) then |
---|
1360 | nx = nd(1) |
---|
1361 | endif |
---|
1362 | |
---|
1363 | if (present(ny)) then |
---|
1364 | ny = nd(2) |
---|
1365 | endif |
---|
1366 | |
---|
1367 | if (present(nz)) then |
---|
1368 | nz = nd(3) |
---|
1369 | endif |
---|
1370 | |
---|
1371 | deallocate(dimid,nd) |
---|
1372 | |
---|
1373 | endif |
---|
1374 | |
---|
1375 | ! Scatter info |
---|
1376 | ! EM Since parallel version, must be broadcasted on a subset of mpi ranks |
---|
1377 | call oasis_mpi_bcast(xtype,mpi_comm_map,subname//' xtype') |
---|
1378 | call oasis_mpi_bcast(ndims,mpi_comm_map,subname//' ndims') |
---|
1379 | |
---|
1380 | if (present(ifld2) .or. present(fld2) .or. present(fld3)) then |
---|
1381 | if (xtype == NF90_INT .and. ndims == 2 .and. present(ifld2)) then |
---|
1382 | call oasis_mpi_bcast(ifld2,mpi_comm_map,subname//' mask') |
---|
1383 | elseif (xtype /= NF90_INT .and. ndims == 2 .and. present(fld2)) then |
---|
1384 | call oasis_mpi_bcast(fld2,mpi_comm_map,subname//' others') |
---|
1385 | elseif (xtype /= NF90_INT .and. ndims == 3 .and. present(fld3)) then |
---|
1386 | call oasis_mpi_bcast(fld3,mpi_comm_map,subname//' corners') |
---|
1387 | endif |
---|
1388 | endif |
---|
1389 | |
---|
1390 | if (present(nx)) call oasis_mpi_bcast(nx,mpi_comm_map,subname//' nx') |
---|
1391 | if (present(ny)) call oasis_mpi_bcast(ny,mpi_comm_map,subname//' nx') |
---|
1392 | if (present(nz)) call oasis_mpi_bcast(nz,mpi_comm_map,subname//' nx') |
---|
1393 | |
---|
1394 | ! endif |
---|
1395 | |
---|
1396 | call oasis_debug_exit(subname) |
---|
1397 | |
---|
1398 | end subroutine oasis_io_read_field_fromroot |
---|
1399 | |
---|
1400 | !=============================================================================== |
---|
1401 | |
---|
1402 | !> Write a real array named field from the root task to a file. |
---|
1403 | |
---|
1404 | subroutine oasis_io_write_2dgridfld_fromroot(filename,fldname,fld,nx,ny) |
---|
1405 | |
---|
1406 | ! --------------------------------------- |
---|
1407 | ! Designed to work with oasis3 write_grid . |
---|
1408 | ! --------------------------------------- |
---|
1409 | |
---|
1410 | implicit none |
---|
1411 | |
---|
1412 | character(len=*), intent(in) :: filename !< file name |
---|
1413 | character(len=*), intent(in) :: fldname !< field name |
---|
1414 | real(ip_realwp_p),intent(in) :: fld(:,:) !< 2d real field |
---|
1415 | integer(ip_i4_p), intent(in) :: nx !< 2d global nx size |
---|
1416 | integer(ip_i4_p), intent(in) :: ny !< 2d global ny size |
---|
1417 | |
---|
1418 | !--- local --- |
---|
1419 | integer(ip_i4_p) :: ncid,dimid,dimid2(2),varid ! cdf info |
---|
1420 | integer(ip_i4_p) :: iflag ! netcdf control |
---|
1421 | integer(ip_i4_p) :: status ! error code |
---|
1422 | integer(ip_i4_p) :: ind ! string index |
---|
1423 | logical :: exists ! file existance |
---|
1424 | character(len=ic_med) :: gridname ! grid name derived from fldname |
---|
1425 | |
---|
1426 | character(len=*),parameter :: subname = '(oasis_io_write_2dgridfld_fromroot)' |
---|
1427 | |
---|
1428 | !------------------------------------------------------------------------------- |
---|
1429 | ! |
---|
1430 | !------------------------------------------------------------------------------- |
---|
1431 | |
---|
1432 | call oasis_debug_enter(subname) |
---|
1433 | |
---|
1434 | ! expects to run only on 1 pe. |
---|
1435 | ! if (iam == master_task) then |
---|
1436 | |
---|
1437 | if (OASIS_debug >= 5) then |
---|
1438 | write(nulprt,*) subname,' write ',trim(filename),' ',trim(fldname),nx,ny |
---|
1439 | endif |
---|
1440 | |
---|
1441 | ind = index(trim(fldname),'.') |
---|
1442 | if (ind < 2) then |
---|
1443 | write(nulprt,*) subname,estr,'in fldname ',trim(fldname) |
---|
1444 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
1445 | endif |
---|
1446 | gridname = fldname(1:ind-1) |
---|
1447 | |
---|
1448 | inquire(file=trim(filename),exist=exists) |
---|
1449 | if (exists) then |
---|
1450 | status = nf90_open(filename,NF90_WRITE,ncid) |
---|
1451 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1452 | status = nf90_redef(ncid) |
---|
1453 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1454 | else |
---|
1455 | iflag = NF90_CLOBBER |
---|
1456 | if (cdf_filetype == 'cdf2') iflag = ior(iflag,oasis_io_64bit_offset) |
---|
1457 | if (cdf_filetype == 'cdf5') iflag = ior(iflag,oasis_io_64bit_data) |
---|
1458 | status = nf90_create(filename,iflag,ncid) |
---|
1459 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1460 | endif |
---|
1461 | |
---|
1462 | ! define x dimension if it doesn't exist |
---|
1463 | status = nf90_inq_dimid(ncid,'x_'//trim(gridname),dimid2(1)) |
---|
1464 | if (status /= NF90_NOERR) then |
---|
1465 | status = nf90_def_dim(ncid,'x_'//trim(gridname),nx,dimid2(1)) |
---|
1466 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1467 | endif |
---|
1468 | |
---|
1469 | ! define y dimension if it doesn't exist |
---|
1470 | status = nf90_inq_dimid(ncid,'y_'//trim(gridname),dimid2(2)) |
---|
1471 | if (status /= NF90_NOERR) then |
---|
1472 | status = nf90_def_dim(ncid,'y_'//trim(gridname),ny,dimid2(2)) |
---|
1473 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1474 | endif |
---|
1475 | |
---|
1476 | ! define var if it doesn't exist |
---|
1477 | status = nf90_inq_varid(ncid,trim(fldname),varid) |
---|
1478 | if (status /= NF90_NOERR) then |
---|
1479 | status = nf90_def_var(ncid,trim(fldname),NF90_DOUBLE,dimid2,varid) |
---|
1480 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1481 | status = nf90_enddef(ncid) |
---|
1482 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1483 | status = nf90_put_var(ncid,varid,fld) |
---|
1484 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1485 | else |
---|
1486 | status = nf90_enddef(ncid) |
---|
1487 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1488 | endif |
---|
1489 | |
---|
1490 | status = nf90_close(ncid) |
---|
1491 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1492 | |
---|
1493 | ! endif |
---|
1494 | |
---|
1495 | call oasis_debug_exit(subname) |
---|
1496 | |
---|
1497 | end subroutine oasis_io_write_2dgridfld_fromroot |
---|
1498 | |
---|
1499 | !=============================================================================== |
---|
1500 | |
---|
1501 | !> Write an integer array named field from the root task to a file. |
---|
1502 | |
---|
1503 | subroutine oasis_io_write_2dgridint_fromroot(filename,fldname,fld,nx,ny) |
---|
1504 | |
---|
1505 | ! --------------------------------------- |
---|
1506 | ! Designed to work with oasis3 write_grid . |
---|
1507 | ! --------------------------------------- |
---|
1508 | |
---|
1509 | implicit none |
---|
1510 | |
---|
1511 | character(len=*), intent(in) :: filename !< file name |
---|
1512 | character(len=*), intent(in) :: fldname !< field name |
---|
1513 | integer(ip_i4_p), intent(in) :: fld(:,:) !< integer field |
---|
1514 | integer(ip_i4_p), intent(in) :: nx !< 2d global nx size |
---|
1515 | integer(ip_i4_p), intent(in) :: ny !< 2d global ny size |
---|
1516 | |
---|
1517 | !--- local --- |
---|
1518 | integer(ip_i4_p) :: ncid,dimid,dimid2(2),varid ! cdf info |
---|
1519 | integer(ip_i4_p) :: iflag ! netcdf control |
---|
1520 | integer(ip_i4_p) :: status ! error code |
---|
1521 | integer(ip_i4_p) :: ind ! string index |
---|
1522 | logical :: exists ! file existance |
---|
1523 | character(len=ic_med) :: gridname ! grid name derived from fldname |
---|
1524 | |
---|
1525 | character(len=*),parameter :: subname = '(oasis_io_write_2dgridint_fromroot)' |
---|
1526 | |
---|
1527 | !------------------------------------------------------------------------------- |
---|
1528 | ! |
---|
1529 | !------------------------------------------------------------------------------- |
---|
1530 | |
---|
1531 | call oasis_debug_enter(subname) |
---|
1532 | |
---|
1533 | ! expects to run only on 1 pe. |
---|
1534 | ! if (iam == master_task) then |
---|
1535 | |
---|
1536 | if (OASIS_debug >= 5) then |
---|
1537 | write(nulprt,*) subname,' write ',trim(filename),' ',trim(fldname),nx,ny |
---|
1538 | endif |
---|
1539 | |
---|
1540 | ind = index(trim(fldname),'.') |
---|
1541 | if (ind < 2) then |
---|
1542 | write(nulprt,*) subname,estr,'in fldname ',trim(fldname) |
---|
1543 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
1544 | endif |
---|
1545 | gridname = fldname(1:ind-1) |
---|
1546 | |
---|
1547 | inquire(file=trim(filename),exist=exists) |
---|
1548 | if (exists) then |
---|
1549 | status = nf90_open(filename,NF90_WRITE,ncid) |
---|
1550 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1551 | status = nf90_redef(ncid) |
---|
1552 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1553 | else |
---|
1554 | iflag = NF90_CLOBBER |
---|
1555 | if (cdf_filetype == 'cdf2') iflag = ior(iflag,oasis_io_64bit_offset) |
---|
1556 | if (cdf_filetype == 'cdf5') iflag = ior(iflag,oasis_io_64bit_data) |
---|
1557 | status = nf90_create(filename,iflag,ncid) |
---|
1558 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1559 | endif |
---|
1560 | |
---|
1561 | ! define x dimension if it doesn't exist |
---|
1562 | status = nf90_inq_dimid(ncid,'x_'//trim(gridname),dimid2(1)) |
---|
1563 | if (status /= NF90_NOERR) then |
---|
1564 | status = nf90_def_dim(ncid,'x_'//trim(gridname),nx,dimid2(1)) |
---|
1565 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1566 | endif |
---|
1567 | |
---|
1568 | ! define y dimension if it doesn't exist |
---|
1569 | status = nf90_inq_dimid(ncid,'y_'//trim(gridname),dimid2(2)) |
---|
1570 | if (status /= NF90_NOERR) then |
---|
1571 | status = nf90_def_dim(ncid,'y_'//trim(gridname),ny,dimid2(2)) |
---|
1572 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1573 | endif |
---|
1574 | |
---|
1575 | ! define var if it doesn't exist |
---|
1576 | status = nf90_inq_varid(ncid,trim(fldname),varid) |
---|
1577 | if (status /= NF90_NOERR) then |
---|
1578 | status = nf90_def_var(ncid,trim(fldname),NF90_INT,dimid2,varid) |
---|
1579 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1580 | status = nf90_enddef(ncid) |
---|
1581 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1582 | status = nf90_put_var(ncid,varid,fld) |
---|
1583 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1584 | else |
---|
1585 | status = nf90_enddef(ncid) |
---|
1586 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1587 | endif |
---|
1588 | |
---|
1589 | status = nf90_close(ncid) |
---|
1590 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1591 | |
---|
1592 | ! endif |
---|
1593 | |
---|
1594 | call oasis_debug_exit(subname) |
---|
1595 | |
---|
1596 | end subroutine oasis_io_write_2dgridint_fromroot |
---|
1597 | |
---|
1598 | !=============================================================================== |
---|
1599 | |
---|
1600 | !> Write a 3d real array named field from the root task to a file. |
---|
1601 | |
---|
1602 | subroutine oasis_io_write_3dgridfld_fromroot(filename,fldname,fld,nx,ny,nc) |
---|
1603 | |
---|
1604 | ! --------------------------------------- |
---|
1605 | ! Designed to work with oasis3 write_grid. |
---|
1606 | ! --------------------------------------- |
---|
1607 | |
---|
1608 | implicit none |
---|
1609 | |
---|
1610 | character(len=*), intent(in) :: filename !< file name |
---|
1611 | character(len=*), intent(in) :: fldname !< field name |
---|
1612 | real(ip_realwp_p), intent(in) :: fld(:,:,:)!< 3d real array |
---|
1613 | integer(ip_i4_p), intent(in) :: nx !< 3d global nx size |
---|
1614 | integer(ip_i4_p), intent(in) :: ny !< 3d global ny size |
---|
1615 | integer(ip_i4_p), intent(in) :: nc !< 3d global nz size or nc size for corners |
---|
1616 | |
---|
1617 | !--- local --- |
---|
1618 | integer(ip_i4_p) :: ncid,dimid,dimid3(3),varid ! cdf info |
---|
1619 | integer(ip_i4_p) :: iflag ! netcdf control |
---|
1620 | integer(ip_i4_p) :: status ! error code |
---|
1621 | integer(ip_i4_p) :: ind ! string index |
---|
1622 | logical :: exists ! file existance |
---|
1623 | character(len=ic_med) :: gridname ! grid name derived from fldname |
---|
1624 | |
---|
1625 | character(len=*),parameter :: subname = '(oasis_io_write_3dgridfld_fromroot)' |
---|
1626 | |
---|
1627 | !------------------------------------------------------------------------------- |
---|
1628 | ! |
---|
1629 | !------------------------------------------------------------------------------- |
---|
1630 | |
---|
1631 | call oasis_debug_enter(subname) |
---|
1632 | |
---|
1633 | ! expects to run only on 1 pe. |
---|
1634 | ! if (iam == master_task) then |
---|
1635 | |
---|
1636 | if (OASIS_debug >= 5) then |
---|
1637 | write(nulprt,*) subname,' write ',trim(filename),' ',trim(fldname),nx,ny,nc |
---|
1638 | endif |
---|
1639 | |
---|
1640 | ind = index(trim(fldname),'.') |
---|
1641 | if (ind < 2) then |
---|
1642 | write(nulprt,*) subname,estr,'in fldname ',trim(fldname) |
---|
1643 | call oasis_abort(file=__FILE__,line=__LINE__) |
---|
1644 | endif |
---|
1645 | gridname = fldname(1:ind-1) |
---|
1646 | |
---|
1647 | inquire(file=trim(filename),exist=exists) |
---|
1648 | if (exists) then |
---|
1649 | status = nf90_open(filename,NF90_WRITE,ncid) |
---|
1650 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1651 | status = nf90_redef(ncid) |
---|
1652 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1653 | else |
---|
1654 | iflag = NF90_CLOBBER |
---|
1655 | if (cdf_filetype == 'cdf2') iflag = ior(iflag,oasis_io_64bit_offset) |
---|
1656 | if (cdf_filetype == 'cdf5') iflag = ior(iflag,oasis_io_64bit_data) |
---|
1657 | status = nf90_create(filename,iflag,ncid) |
---|
1658 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1659 | endif |
---|
1660 | |
---|
1661 | ! define x dimension if it doesn't exist |
---|
1662 | status = nf90_inq_dimid(ncid,'x_'//trim(gridname),dimid3(1)) |
---|
1663 | if (status /= NF90_NOERR) then |
---|
1664 | status = nf90_def_dim(ncid,'x_'//trim(gridname),nx,dimid3(1)) |
---|
1665 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1666 | endif |
---|
1667 | |
---|
1668 | ! define y dimension if it doesn't exist |
---|
1669 | status = nf90_inq_dimid(ncid,'y_'//trim(gridname),dimid3(2)) |
---|
1670 | if (status /= NF90_NOERR) then |
---|
1671 | status = nf90_def_dim(ncid,'y_'//trim(gridname),ny,dimid3(2)) |
---|
1672 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1673 | endif |
---|
1674 | |
---|
1675 | ! define crn dimension if it doesn't exist |
---|
1676 | status = nf90_inq_dimid(ncid,'crn_'//trim(gridname),dimid3(3)) |
---|
1677 | if (status /= NF90_NOERR) then |
---|
1678 | status = nf90_def_dim(ncid,'crn_'//trim(gridname),nc,dimid3(3)) |
---|
1679 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1680 | endif |
---|
1681 | |
---|
1682 | ! define var if it doesn't exist |
---|
1683 | status = nf90_inq_varid(ncid,trim(fldname),varid) |
---|
1684 | if (status /= NF90_NOERR) then |
---|
1685 | status = nf90_def_var(ncid,trim(fldname),NF90_DOUBLE,dimid3,varid) |
---|
1686 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1687 | status = nf90_enddef(ncid) |
---|
1688 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1689 | status = nf90_put_var(ncid,varid,fld) |
---|
1690 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1691 | else |
---|
1692 | status = nf90_enddef(ncid) |
---|
1693 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1694 | endif |
---|
1695 | |
---|
1696 | status = nf90_close(ncid) |
---|
1697 | call check_status(status,subname,__FILE__,__LINE__) |
---|
1698 | |
---|
1699 | ! endif |
---|
1700 | |
---|
1701 | call oasis_debug_exit(subname) |
---|
1702 | |
---|
1703 | end subroutine oasis_io_write_3dgridfld_fromroot |
---|
1704 | |
---|
1705 | !------------------------------------------------------------------- |
---|
1706 | |
---|
1707 | subroutine check_status(stat,subn,file,line) |
---|
1708 | |
---|
1709 | integer(ip_i4_p), intent(in) :: stat |
---|
1710 | character(len=*), intent(in) :: subn |
---|
1711 | character(len=*), intent(in) :: file |
---|
1712 | integer(ip_i4_p), intent(in) :: line |
---|
1713 | |
---|
1714 | character(len=*),parameter :: subname = '(check_status)' |
---|
1715 | |
---|
1716 | if (stat /= NF90_NOERR) then |
---|
1717 | write(nulprt,*) subname,estr,' netcdf status' |
---|
1718 | write(nulprt,*) trim(subn),' nf90_strerror = ',trim(nf90_strerror(stat)) |
---|
1719 | call oasis_abort(file=file,line=line) |
---|
1720 | endif |
---|
1721 | |
---|
1722 | end subroutine check_status |
---|
1723 | |
---|
1724 | !------------------------------------------------------------------- |
---|
1725 | |
---|
1726 | END MODULE mod_oasis_io |
---|