SCALE-RM
scale_external_io.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
15  !-----------------------------------------------------------------------------
16  !
17  !++ used modules
18  !
19  use gtool_file, only: &
21  use scale_precision
22  use scale_stdio
24  use scale_process, only: &
26  !-----------------------------------------------------------------------------
27  implicit none
28  private
29  !-----------------------------------------------------------------------------
30  !
31  !++ Public procedure
32  !
33  public :: externalfilegetshape
34  public :: externalfilegetglobalattv
36  public :: externalfilevarexistence
37  public :: externalfileread
38  public :: externalfilereadoffset
39 
40  interface externalfilegetglobalattv
41  module procedure externalfilegetglobalattvinteger
42  module procedure externalfilegetglobalattvrealsp
43  module procedure externalfilegetglobalattvrealdp
44  end interface externalfilegetglobalattv
45 
46  interface externalfileread
47  module procedure externalfileread2drealsp
48  module procedure externalfileread2drealdp
49  module procedure externalfileread3drealsp
50  module procedure externalfileread3drealdp
51  module procedure externalfileread4drealsp
52  module procedure externalfileread4drealdp
53  end interface externalfileread
54 
55  interface externalfilereadoffset
56  !module procedure ExternalFileRead2DRealSP
57  !module procedure ExternalFileRead2DRealDP
58  module procedure externalfilereadoffset3drealsp
59  module procedure externalfilereadoffset3drealdp
60  module procedure externalfilereadoffset4drealsp
61  module procedure externalfilereadoffset4drealdp
62  end interface externalfilereadoffset
63 
64  !-----------------------------------------------------------------------------
65  !
66  !++ Public parameters & variables
67  !
68  integer, public, parameter :: iscale = 1 ! use gtool, coz it's not external
69  integer, public, parameter :: iwrfarw = 2
70  integer, public, parameter :: inicam = 3
71  integer, public, parameter :: igrads = 4
72 
73 
74  !-----------------------------------------------------------------------------
75  !
76  !++ Private procedure
77  !
78  private :: externalfilemakefname
79  private :: externaltakedimension
80  private :: convertarrayorder
81  private :: handle_err
82 
83  interface convertarrayorder
84  module procedure convertarrayorderwrf2dsp
85  module procedure convertarrayorderwrf2ddp
86  module procedure convertarrayorderwrf3dsp
87  module procedure convertarrayorderwrf3ddp
88  module procedure convertarrayorderwrf4dsp
89  module procedure convertarrayorderwrf4ddp
90  end interface convertarrayorder
91 
92  !-----------------------------------------------------------------------------
93  !
94  !++ Private parameters & variables
95  !
96  !-----------------------------------------------------------------------------
97 contains
98 
99  !-----------------------------------------------------------------------------
100  ! ExternalFileGetShape
101  !-----------------------------------------------------------------------------
102  subroutine externalfilegetshape( &
103  dims, & ! (out)
104  timelen, & ! (out)
105  mdlid, & ! (in)
106  basename, & ! (in)
107  myrank, & ! (in)
108  single & ! (in) optional
109  )
110  use netcdf ![external lib]
111  implicit none
112 
113  integer, intent(out) :: dims(:)
114  integer, intent(out) :: timelen
115  integer, intent( in) :: mdlid
116  character(len=*), intent( in) :: basename
117  integer, intent( in) :: myrank
118  logical, intent( in), optional :: single
119 
120  integer :: status
121  integer :: ncid, unlimid
122  integer :: dims_org(7)
123  character(len=NF90_MAX_NAME) :: tname
124  character(len=H_LONG) :: fname = ''
125  logical :: single_ = .false.
126 
127  intrinsic size
128  intrinsic shape
129  !---------------------------------------------------------------------------
130 
131  if ( present(single) ) then
132  single_ = single
133  else
134  single_ = .false.
135  endif
136 
137  call externalfilemakefname( fname,mdlid,basename,myrank,single_ )
138 
139  status = nf90_open( trim(fname), nf90_nowrite, ncid )
140  if (status /= nf90_noerr) call handle_err(status, __line__)
141 
142  status = nf90_inquire( ncid, unlimiteddimid=unlimid )
143  if (status /= nf90_noerr) call handle_err(status, __line__)
144 
145  status = nf90_inquire_dimension( ncid, unlimid, tname, timelen )
146  if (status /= nf90_noerr) call handle_err(status, __line__)
147 
148  if( trim(tname)=='time' .OR. trim(tname)=='Time' ) then
149  if( io_l ) write(io_fid_log,*) 'Time Dimension Name: '//trim(tname)
150  else
151  write(*,*) 'xxx Not appropriate time dimension is used in the external file. Check!'
152  call prc_mpistop
153  endif
154 
155  call externaltakedimension( dims_org(:),ncid,mdlid )
156 
157  ! convert dimension order for return to scale-system
158  if( mdlid == iwrfarw )then !MODEL ID: WRF-ARW
159  dims(1) = dims_org(3)
160  dims(2) = dims_org(1)
161  dims(3) = dims_org(2)
162  dims(4) = dims_org(6)
163  dims(5) = dims_org(4)
164  dims(6) = dims_org(5)
165  dims(7) = dims_org(7)
166  else
167  dims(:) = dims_org(:)
168  endif
169 
170  status = nf90_close(ncid)
171  if (status /= nf90_noerr) call handle_err(status, __line__)
172 
173  return
174  end subroutine externalfilegetshape
175 
176  !-----------------------------------------------------------------------------
177  ! ExternalFileGet Global Attribute (value, real, single precision)
178  !-----------------------------------------------------------------------------
180  var, & ! (out)
181  mdlid, & ! (in)
182  basename, & ! (in)
183  attname, & ! (in)
184  myrank, & ! (in)
185  single & ! (in) optional
186  )
187  use netcdf ![external lib]
188  implicit none
189 
190  integer, intent(out) :: var(:)
191  integer, intent( in) :: mdlid
192  character(len=*), intent( in) :: basename
193  character(len=*), intent( in) :: attname
194  integer, intent( in) :: myrank
195  logical, intent( in), optional :: single
196 
197  integer, allocatable :: work(:)
198 
199  integer :: status
200  integer :: i, ncid, length
201  character(len=H_LONG) :: fname = ''
202  logical :: single_ = .false.
203 
204  intrinsic size
205  intrinsic shape
206  !---------------------------------------------------------------------------
207 
208  if ( present(single) ) then
209  single_ = single
210  else
211  single_ = .false.
212  endif
213 
214  call externalfilemakefname( fname,mdlid,basename,myrank,single_ )
215 
216  status = nf90_open( trim(fname), nf90_nowrite, ncid )
217  if (status /= nf90_noerr) call handle_err(status, __line__)
218 
219  status = nf90_inquire_attribute(ncid, nf90_global, trim(attname), len=length)
220  if (status /= nf90_noerr) call handle_err(status, __line__)
221 
222  allocate( work(length) )
223 
224  status = nf90_get_att(ncid, nf90_global, trim(attname), work)
225  if (status /= nf90_noerr) call handle_err(status, __line__)
226 
227  do i = 1, length
228  var(i) = work(i)
229  enddo
230 
231  status = nf90_close(ncid)
232  if (status /= nf90_noerr) call handle_err(status, __line__)
233  deallocate( work )
234 
235  return
236  end subroutine externalfilegetglobalattvinteger
237 
238  !-----------------------------------------------------------------------------
239  ! ExternalFileGet Global Attribute (value, real, single precision)
240  !-----------------------------------------------------------------------------
241  subroutine externalfilegetglobalattvrealsp( &
242  var, & ! (out)
243  mdlid, & ! (in)
244  basename, & ! (in)
245  attname, & ! (in)
246  myrank, & ! (in)
247  single & ! (in) optional
248  )
249  use netcdf ![external lib]
250  implicit none
251 
252  real(SP), intent(out) :: var(:)
253  integer, intent( in) :: mdlid
254  character(len=*), intent( in) :: basename
255  character(len=*), intent( in) :: attname
256  integer, intent( in) :: myrank
257  logical, intent( in), optional :: single
258 
259  real(SP), allocatable :: work(:)
260 
261  integer :: status
262  integer :: i, ncid, length
263  character(len=H_LONG) :: fname = ''
264  logical :: single_ = .false.
265 
266  intrinsic size
267  intrinsic shape
268  !---------------------------------------------------------------------------
269 
270  if ( present(single) ) then
271  single_ = single
272  else
273  single_ = .false.
274  endif
275 
276  call externalfilemakefname( fname,mdlid,basename,myrank,single_ )
277 
278  status = nf90_open( trim(fname), nf90_nowrite, ncid )
279  if (status /= nf90_noerr) call handle_err(status, __line__)
280 
281  status = nf90_inquire_attribute(ncid, nf90_global, trim(attname), len=length)
282  if (status /= nf90_noerr) call handle_err(status, __line__)
283 
284  allocate( work(length) )
285 
286  status = nf90_get_att(ncid, nf90_global, trim(attname), work)
287  if (status /= nf90_noerr) call handle_err(status, __line__)
288 
289  do i = 1, length
290  var(i) = work(i)
291  enddo
292 
293  status = nf90_close(ncid)
294  if (status /= nf90_noerr) call handle_err(status, __line__)
295  deallocate( work )
296 
297  return
298  end subroutine externalfilegetglobalattvrealsp
299 
300  !-----------------------------------------------------------------------------
301  ! ExternalFileGet Global Attribute (value, real, single precision)
302  !-----------------------------------------------------------------------------
303  subroutine externalfilegetglobalattvrealdp( &
304  var, & ! (out)
305  mdlid, & ! (in)
306  basename, & ! (in)
307  attname, & ! (in)
308  myrank, & ! (in)
309  single & ! (in) optional
310  )
311  use netcdf ![external lib]
312  implicit none
313 
314  real(DP), intent(out) :: var(:)
315  integer, intent( in) :: mdlid
316  character(len=*), intent( in) :: basename
317  character(len=*), intent( in) :: attname
318  integer, intent( in) :: myrank
319  logical, intent( in), optional :: single
320 
321  real(DP), allocatable :: work(:)
322 
323  integer :: status
324  integer :: i, ncid, length
325  character(len=H_LONG) :: fname = ''
326  logical :: single_ = .false.
327 
328  intrinsic size
329  intrinsic shape
330  !---------------------------------------------------------------------------
331 
332  if ( present(single) ) then
333  single_ = single
334  else
335  single_ = .false.
336  endif
337 
338  call externalfilemakefname( fname,mdlid,basename,myrank,single_ )
339 
340  status = nf90_open( trim(fname), nf90_nowrite, ncid )
341  if (status /= nf90_noerr) call handle_err(status, __line__)
342 
343  status = nf90_inquire_attribute(ncid, nf90_global, trim(attname), len=length)
344  if (status /= nf90_noerr) call handle_err(status, __line__)
345 
346  allocate( work(length) )
347 
348  status = nf90_get_att(ncid, nf90_global, trim(attname), work)
349  if (status /= nf90_noerr) call handle_err(status, __line__)
350 
351  do i = 1, length
352  var(i) = work(i)
353  enddo
354 
355  status = nf90_close(ncid)
356  if (status /= nf90_noerr) call handle_err(status, __line__)
357  deallocate( work )
358 
359  return
360  end subroutine externalfilegetglobalattvrealdp
361 
362  !-----------------------------------------------------------------------------
363  ! ExternalFileGet Global Attribute (character)
364  !-----------------------------------------------------------------------------
365  subroutine externalfilegetglobalattc( &
366  chr, & ! (out)
367  mdlid, & ! (in)
368  basename, & ! (in)
369  attname, & ! (in)
370  myrank, & ! (in)
371  single & ! (in) optional
372  )
373  use netcdf ![external lib]
374  implicit none
375 
376  character(len=*), intent(out) :: chr(:)
377  integer, intent( in) :: mdlid
378  character(len=*), intent( in) :: basename
379  character(len=*), intent( in) :: attname
380  integer, intent( in) :: myrank
381  logical, intent( in), optional :: single
382 
383  integer :: status
384  integer :: ncid, length
385  character(len=H_LONG) :: fname = ''
386  character(len=80) :: work
387  logical :: single_ = .false.
388 
389  intrinsic size
390  intrinsic shape
391  !---------------------------------------------------------------------------
392 
393  if ( present(single) ) then
394  single_ = single
395  else
396  single_ = .false.
397  endif
398 
399  call externalfilemakefname( fname,mdlid,basename,myrank,single_ )
400 
401  status = nf90_open( trim(fname), nf90_nowrite, ncid )
402  if (status /= nf90_noerr) call handle_err(status, __line__)
403 
404  status = nf90_inquire_attribute(ncid, nf90_global, trim(attname), len=length)
405  if (status /= nf90_noerr) call handle_err(status, __line__)
406 
407  if( len(work) < length ) then
408  write(*,*) 'xxx Not enough space to put attribute values. [externalio/scalelib]'
409  call prc_mpistop
410  endif
411 
412  status = nf90_get_att(ncid, nf90_global, trim(attname), work)
413  if (status /= nf90_noerr) call handle_err(status, __line__)
414 
415  chr = trim(work)
416 
417  return
418  end subroutine externalfilegetglobalattc
419 
420  !-----------------------------------------------------------------------------
422  !-----------------------------------------------------------------------------
423  subroutine externalfilevarexistence( &
424  existence, & ! (out)
425  basename, & ! (in)
426  varname, & ! (in)
427  myrank, & ! (in)
428  mdlid, & ! (in)
429  single & ! (in) optional
430  )
431  use netcdf ![external lib]
432  implicit none
433 
434  logical, intent(out) :: existence
435  character(len=*), intent( in) :: basename
436  character(len=*), intent( in) :: varname
437  integer, intent( in) :: myrank
438  integer, intent( in) :: mdlid
439  logical, intent( in), optional :: single
440 
441  integer :: ncid, varid
442  integer :: status
443 
444  character(len=H_LONG) :: fname = ''
445  logical :: single_ = .false.
446 
447  intrinsic size
448  intrinsic shape
449  !---------------------------------------------------------------------------
450 
451  if ( present(single) ) then
452  single_ = single
453  else
454  single_ = .false.
455  endif
456 
457  call externalfilemakefname( fname,mdlid,basename,myrank,single_ )
458 
459  status = nf90_open( trim(fname), nf90_nowrite, ncid )
460  if (status /= nf90_noerr) call handle_err(status, __line__)
461 
462  status = nf90_inq_varid( ncid, trim(varname), varid )
463  if (status == nf90_noerr) then
464  existence = .true.
465  else
466  existence = .false.
467  if( io_l ) write(io_fid_log,*) '+++ not exist variable: ', trim(varname)
468  endif
469 
470  status = nf90_close(ncid)
471  if (status /= nf90_noerr) call handle_err(status, __line__)
472 
473  return
474  end subroutine externalfilevarexistence
475 
476  !-----------------------------------------------------------------------------
478  !-----------------------------------------------------------------------------
479  subroutine externalfileread2drealsp( &
480  var, & ! (out)
481  basename, & ! (in)
482  varname, & ! (in)
483  ts, & ! (in)
484  te, & ! (in)
485  myrank, & ! (in)
486  mdlid, & ! (in)
487  nx, & ! (in)
488  single & ! (in) optional
489  )
490  use netcdf ![external lib]
491  implicit none
492 
493  real(SP), intent(out) :: var(:,:)
494  character(len=*), intent( in) :: basename
495  character(len=*), intent( in) :: varname
496  integer, intent( in) :: ts
497  integer, intent( in) :: te
498  integer, intent( in) :: myrank
499  integer, intent( in) :: mdlid
500  integer, intent( in) :: nx
501  logical, intent( in), optional :: single
502 
503  real(SP), allocatable :: var_org(:,:)
504  integer :: ncid, varid
505  integer :: status
506  integer :: precis
507 
508  integer :: tcount
509  character(len=H_LONG) :: fname = ''
510  logical :: single_ = .false.
511 
512  intrinsic size
513  intrinsic shape
514  !---------------------------------------------------------------------------
515 
516  tcount = te - ts + 1
517 
518  if ( present(single) ) then
519  single_ = single
520  else
521  single_ = .false.
522  endif
523 
524  call externalfilemakefname( fname,mdlid,basename,myrank,single_ )
525 
526  status = nf90_open( trim(fname), nf90_nowrite, ncid )
527  if (status /= nf90_noerr) call handle_err(status, __line__)
528 
529  ! based on specified dimension size
530  allocate( var_org(nx,tcount) )
531 
532  status = nf90_inq_varid( ncid, trim(varname), varid )
533  if (status /= nf90_noerr) call handle_err(status, __line__)
534 
535  status = nf90_inquire_variable( ncid, varid, xtype=precis )
536  if(status /= nf90_noerr) call handle_err(status, __line__)
537  if(precis /= nf90_float) then
538  write(*,*) 'xxx Internal Error: [scale_external_io]/[ExternalFileRead2DSP]'
539  call prc_mpistop
540  endif
541 
542  status = nf90_get_var( ncid, varid, var_org(:,:), start = (/ 1,ts /), &
543  count = (/ nx,tcount /) )
544  if (status /= nf90_noerr) call handle_err(status, __line__)
545 
546  status = nf90_close(ncid)
547  if (status /= nf90_noerr) call handle_err(status, __line__)
548 
549  call convertarrayorder( var,var_org,tcount,nx )
550 
551  deallocate( var_org )
552 
553  return
554  end subroutine externalfileread2drealsp
555  subroutine externalfileread2drealdp( &
556  var, & ! (out)
557  basename, & ! (in)
558  varname, & ! (in)
559  ts, & ! (in)
560  te, & ! (in)
561  myrank, & ! (in)
562  mdlid, & ! (in)
563  nx, & ! (in)
564  single & ! (in) optional
565  )
566  use netcdf ![external lib]
567  implicit none
568 
569  real(DP), intent(out) :: var(:,:)
570  character(len=*), intent( in) :: basename
571  character(len=*), intent( in) :: varname
572  integer, intent( in) :: ts
573  integer, intent( in) :: te
574  integer, intent( in) :: myrank
575  integer, intent( in) :: mdlid
576  integer, intent( in) :: nx
577  logical, intent( in), optional :: single
578 
579  real(DP), allocatable :: var_org(:,:)
580  integer :: ncid, varid
581  integer :: status
582  integer :: precis
583 
584  integer :: tcount
585  character(len=H_LONG) :: fname = ''
586  logical :: single_ = .false.
587 
588  intrinsic size
589  intrinsic shape
590  !---------------------------------------------------------------------------
591 
592  tcount = te - ts + 1
593 
594  if ( present(single) ) then
595  single_ = single
596  else
597  single_ = .false.
598  endif
599 
600  call externalfilemakefname( fname,mdlid,basename,myrank,single_ )
601 
602  status = nf90_open( trim(fname), nf90_nowrite, ncid )
603  if (status /= nf90_noerr) call handle_err(status, __line__)
604 
605  ! based on specified dimension size
606  allocate( var_org(nx,tcount) )
607 
608  status = nf90_inq_varid( ncid, trim(varname), varid )
609  if (status /= nf90_noerr) call handle_err(status, __line__)
610 
611  status = nf90_inquire_variable( ncid, varid, xtype=precis )
612  if(status /= nf90_noerr) call handle_err(status, __line__)
613  if(precis /= nf90_double) then
614  write(*,*) 'xxx Internal Error: [scale_external_io]/[ExternalFileRead2DDP]'
615  call prc_mpistop
616  endif
617 
618  status = nf90_get_var( ncid, varid, var_org(:,:), start = (/ 1,ts /), &
619  count = (/ nx,tcount /) )
620  if (status /= nf90_noerr) call handle_err(status, __line__)
621 
622  status = nf90_close(ncid)
623  if (status /= nf90_noerr) call handle_err(status, __line__)
624 
625  call convertarrayorder( var,var_org,tcount,nx )
626 
627  deallocate( var_org )
628 
629  return
630  end subroutine externalfileread2drealdp
631  subroutine externalfileread3drealsp( &
632  var, & ! (out)
633  basename, & ! (in)
634  varname, & ! (in)
635  ts, & ! (in)
636  te, & ! (in)
637  myrank, & ! (in)
638  mdlid, & ! (in)
639  single, & ! (in) optional
640  xstag, & ! (in) optional
641  ystag, & ! (in) optional
642  option & ! (in) optional
643  )
644  use netcdf ![external lib]
645  implicit none
646 
647  real(SP), intent(out) :: var(:,:,:)
648  character(len=*), intent( in) :: basename
649  character(len=*), intent( in) :: varname
650  integer, intent( in) :: ts
651  integer, intent( in) :: te
652  integer, intent( in) :: myrank
653  integer, intent( in) :: mdlid
654  logical, intent( in), optional :: single
655  logical, intent( in), optional :: xstag
656  logical, intent( in), optional :: ystag
657  logical, intent( in), optional :: option
658 
659  real(SP), allocatable :: var_org(:,:,:)
660  integer :: ncid, varid
661  integer :: status
662  integer :: precis
663  integer :: nx, ny
664  integer :: dims(7)
665 
666  integer :: tcount
667  character(len=H_LONG) :: fname = ''
668  logical :: single_
669  logical :: option_
670 
671  intrinsic size
672  intrinsic shape
673  !---------------------------------------------------------------------------
674 
675  single_ = .false.
676  option_ = .false.
677 
678  tcount = te - ts + 1
679 
680  if ( present(single) ) then
681  single_ = single
682  endif
683 
684  if ( present(option) ) then
685  option_ = option
686  end if
687 
688  call externalfilemakefname( fname,mdlid,basename,myrank,single_ )
689 
690  status = nf90_open( trim(fname), nf90_nowrite, ncid )
691  if (status /= nf90_noerr) then
692  if ( option_ ) then
693  return
694  else
695  call handle_err(status, __line__)
696  end if
697  end if
698 
699  ! retrieve dimension size in data original order
700  call externaltakedimension( dims(:),ncid,mdlid )
701  nx = dims(1)
702  if ( present(xstag) ) then
703  if ( xstag ) then
704  nx = dims(4)
705  endif
706  endif
707  ny = dims(2)
708  if ( present(ystag) ) then
709  if ( ystag ) then
710  ny = dims(5)
711  endif
712  endif
713  allocate( var_org(nx,ny,tcount) )
714 
715  status = nf90_inq_varid( ncid, trim(varname), varid )
716  if (status /= nf90_noerr) call handle_err(status, __line__)
717 
718  status = nf90_inquire_variable( ncid, varid, xtype=precis )
719  if(status /= nf90_noerr) call handle_err(status, __line__)
720  if(precis /= nf90_float) then
721  write(*,*) 'xxx Internal Error: [scale_external_io]/[ExternalFileRead3DSP]'
722  call prc_mpistop
723  endif
724 
725  status = nf90_get_var( ncid, varid, var_org(:,:,:), start = (/ 1,1,ts /), &
726  count = (/ nx,ny,tcount /) )
727  if (status /= nf90_noerr) call handle_err(status, __line__)
728 
729  status = nf90_close(ncid)
730  if (status /= nf90_noerr) call handle_err(status, __line__)
731 
732  call convertarrayorder( var,var_org,tcount,nx,ny )
733 
734  deallocate( var_org )
735 
736  return
737  end subroutine externalfileread3drealsp
738  subroutine externalfileread3drealdp( &
739  var, & ! (out)
740  basename, & ! (in)
741  varname, & ! (in)
742  ts, & ! (in)
743  te, & ! (in)
744  myrank, & ! (in)
745  mdlid, & ! (in)
746  single, & ! (in) optional
747  xstag, & ! (in) optional
748  ystag & ! (in) optional
749  )
750  use netcdf ![external lib]
751  implicit none
752 
753  real(DP), intent(out) :: var(:,:,:)
754  character(len=*), intent( in) :: basename
755  character(len=*), intent( in) :: varname
756  integer, intent( in) :: ts
757  integer, intent( in) :: te
758  integer, intent( in) :: myrank
759  integer, intent( in) :: mdlid
760  logical, intent( in), optional :: single
761  logical, intent( in), optional :: xstag
762  logical, intent( in), optional :: ystag
763 
764  real(DP), allocatable :: var_org(:,:,:)
765  integer :: ncid, varid
766  integer :: status
767  integer :: precis
768  integer :: nx, ny
769  integer :: dims(7)
770 
771  integer :: tcount
772  character(len=H_LONG) :: fname = ''
773  logical :: single_ = .false.
774 
775  intrinsic size
776  intrinsic shape
777  !---------------------------------------------------------------------------
778 
779  tcount = te - ts + 1
780 
781  if ( present(single) ) then
782  single_ = single
783  else
784  single_ = .false.
785  endif
786 
787  call externalfilemakefname( fname,mdlid,basename,myrank,single )
788 
789  status = nf90_open( trim(fname), nf90_nowrite, ncid )
790  if (status /= nf90_noerr) call handle_err(status, __line__)
791 
792  ! retrieve dimension size in data original order
793  call externaltakedimension( dims(:),ncid,mdlid )
794  nx = dims(1)
795  if ( present(xstag) ) then
796  if ( xstag ) then
797  nx = dims(4)
798  endif
799  endif
800  ny = dims(2)
801  if ( present(ystag) ) then
802  if ( ystag ) then
803  ny = dims(5)
804  endif
805  endif
806  allocate( var_org(nx,ny,tcount) )
807 
808  status = nf90_inq_varid( ncid, trim(varname), varid )
809  if (status /= nf90_noerr) call handle_err(status, __line__)
810 
811  status = nf90_inquire_variable( ncid, varid, xtype=precis )
812  if(status /= nf90_noerr) call handle_err(status, __line__)
813  if(precis /= nf90_double) then
814  write(*,*) 'xxx Internal Error: [scale_external_io]/[ExternalFileRead3DDP]'
815  call prc_mpistop
816  endif
817 
818  status = nf90_get_var( ncid, varid, var_org(:,:,:), start = (/ 1,1,ts /), &
819  count = (/ nx,ny,tcount /) )
820  if (status /= nf90_noerr) call handle_err(status, __line__)
821 
822  status = nf90_close(ncid)
823  if (status /= nf90_noerr) call handle_err(status, __line__)
824 
825  call convertarrayorder( var,var_org,tcount,nx,ny )
826 
827  deallocate( var_org )
828 
829  return
830  end subroutine externalfileread3drealdp
831  subroutine externalfileread4drealsp( &
832  var, & ! (out)
833  basename, & ! (in)
834  varname, & ! (in)
835  ts, & ! (in)
836  te, & ! (in)
837  myrank, & ! (in)
838  mdlid, & ! (in)
839  single, & ! (in) optional
840  xstag, & ! (in) optional
841  ystag, & ! (in) optional
842  zstag, & ! (in) optional
843  landgrid, & ! (in) optional
844  option & ! (in) optional
845  )
846  use netcdf ![external lib]
847  implicit none
848 
849  real(SP), intent(out) :: var(:,:,:,:)
850  character(len=*), intent( in) :: basename
851  character(len=*), intent( in) :: varname
852  integer, intent( in) :: ts
853  integer, intent( in) :: te
854  integer, intent( in) :: myrank
855  integer, intent( in) :: mdlid
856  logical, intent( in), optional :: single
857  logical, intent( in), optional :: xstag
858  logical, intent( in), optional :: ystag
859  logical, intent( in), optional :: zstag
860  logical, intent( in), optional :: landgrid
861  logical, intent( in), optional :: option
862 
863  real(SP), allocatable :: var_org(:,:,:,:)
864  integer :: ncid, varid
865  integer :: status
866  integer :: precis
867  integer :: nx, ny, nz
868  integer :: dims(7)
869 
870  integer :: tcount
871  character(len=H_LONG) :: fname = ''
872  logical :: single_
873  logical :: option_
874 
875  intrinsic size
876  intrinsic shape
877  !---------------------------------------------------------------------------
878 
879  single_ = .false.
880  option_ = .false.
881 
882  tcount = te - ts + 1
883 
884  if ( present(single) ) then
885  single_ = single
886  endif
887 
888  if ( present(option) ) then
889  option_ = option
890  endif
891 
892  call externalfilemakefname( fname,mdlid,basename,myrank,single_ )
893 
894  status = nf90_open( trim(fname), nf90_nowrite, ncid )
895  if (status /= nf90_noerr) then
896  if ( option_ ) then
897  return
898  else
899  write(*,*) trim(fname)
900  call handle_err(status, __line__)
901  end if
902  end if
903 
904  ! retrieve dimension size in data original order
905  call externaltakedimension( dims(:),ncid,mdlid )
906  nx = dims(1)
907  if ( present(xstag) ) then
908  if ( xstag ) then
909  nx = dims(4)
910  endif
911  endif
912  ny = dims(2)
913  if ( present(ystag) ) then
914  if ( ystag ) then
915  ny = dims(5)
916  endif
917  endif
918  nz = dims(3)
919  if ( present(zstag) ) then
920  if ( zstag ) then
921  nz = dims(6)
922  endif
923  endif
924  if ( present(landgrid) ) then
925  if ( landgrid ) then
926  nz = dims(7)
927  endif
928  endif
929  allocate( var_org(nx,ny,nz,tcount) )
930 
931  status = nf90_inq_varid( ncid, trim(varname), varid )
932  if (status /= nf90_noerr) call handle_err(status, __line__)
933 
934  status = nf90_inquire_variable( ncid, varid, xtype=precis )
935  if(status /= nf90_noerr) call handle_err(status, __line__)
936  if(precis /= nf90_float) then
937  write(*,*) 'xxx Internal Error: [scale_external_io]/[ExternalFileRead4DSP]'
938  call prc_mpistop
939  endif
940 
941  status = nf90_get_var( ncid, varid, var_org(:,:,:,:), start = (/ 1,1,1,ts /), &
942  count = (/ nx,ny,nz,tcount /) )
943  if (status /= nf90_noerr) call handle_err(status, __line__)
944 
945  status = nf90_close(ncid)
946  if (status /= nf90_noerr) call handle_err(status, __line__)
947 
948  call convertarrayorder( var,var_org,tcount,nz,nx,ny )
949 
950  deallocate( var_org )
951 
952  return
953  end subroutine externalfileread4drealsp
954  subroutine externalfileread4drealdp( &
955  var, & ! (out)
956  basename, & ! (in)
957  varname, & ! (in)
958  ts, & ! (in)
959  te, & ! (in)
960  myrank, & ! (in)
961  mdlid, & ! (in)
962  single, & ! (in) optional
963  xstag, & ! (in) optional
964  ystag, & ! (in) optional
965  zstag, & ! (in) optional
966  landgrid & ! (in) optional
967  )
968  use netcdf ![external lib]
969  implicit none
970 
971  real(DP), intent(out) :: var(:,:,:,:)
972  character(len=*), intent( in) :: basename
973  character(len=*), intent( in) :: varname
974  integer, intent( in) :: ts
975  integer, intent( in) :: te
976  integer, intent( in) :: myrank
977  integer, intent( in) :: mdlid
978  logical, intent( in), optional :: single
979  logical, intent( in), optional :: xstag
980  logical, intent( in), optional :: ystag
981  logical, intent( in), optional :: zstag
982  logical, intent( in), optional :: landgrid
983 
984  real(DP), allocatable :: var_org(:,:,:,:)
985  integer :: ncid, varid
986  integer :: status
987  integer :: precis
988  integer :: nx, ny, nz
989  integer :: dims(7)
990 
991  integer :: tcount
992  character(len=H_LONG) :: fname = ''
993  logical :: single_ = .false.
994 
995  intrinsic size
996  intrinsic shape
997  !---------------------------------------------------------------------------
998 
999  tcount = te - ts + 1
1000 
1001  if ( present(single) ) then
1002  single_ = single
1003  else
1004  single_ = .false.
1005  endif
1006 
1007  call externalfilemakefname( fname,mdlid,basename,myrank,single )
1008 
1009  status = nf90_open( trim(fname), nf90_nowrite, ncid )
1010  if (status /= nf90_noerr) call handle_err(status, __line__)
1011 
1012  ! retrieve dimension size in data original order
1013  call externaltakedimension( dims(:),ncid,mdlid )
1014  nx = dims(1)
1015  if ( present(xstag) ) then
1016  if ( xstag ) then
1017  nx = dims(4)
1018  endif
1019  endif
1020  ny = dims(2)
1021  if ( present(ystag) ) then
1022  if ( ystag ) then
1023  ny = dims(5)
1024  endif
1025  endif
1026  nz = dims(3)
1027  if ( present(zstag) ) then
1028  if ( zstag ) then
1029  nz = dims(6)
1030  endif
1031  endif
1032  if ( present(landgrid) ) then
1033  if ( landgrid ) then
1034  nz = dims(7)
1035  endif
1036  endif
1037  allocate( var_org(nx,ny,nz,tcount) )
1038 
1039  status = nf90_inq_varid( ncid, trim(varname), varid )
1040  if (status /= nf90_noerr) call handle_err(status, __line__)
1041 
1042  status = nf90_inquire_variable( ncid, varid, xtype=precis )
1043  if(status /= nf90_noerr) call handle_err(status, __line__)
1044  if(precis /= nf90_double) then
1045  write(*,*) 'xxx Internal Error: [scale_external_io]/[ExternalFileRead4DDP]'
1046  call prc_mpistop
1047  endif
1048 
1049  status = nf90_get_var( ncid, varid, var_org(:,:,:,:), start = (/ 1,1,1,ts /), &
1050  count = (/ nx,ny,nz,tcount /) )
1051  if (status /= nf90_noerr) call handle_err(status, __line__)
1052 
1053  status = nf90_close(ncid)
1054  if (status /= nf90_noerr) call handle_err(status, __line__)
1055 
1056  call convertarrayorder( var,var_org,tcount,nz,nx,ny )
1057 
1058  deallocate( var_org )
1059 
1060  return
1061  end subroutine externalfileread4drealdp
1062 
1063  subroutine externalfilereadoffset3drealsp( &
1064  var, & ! (out)
1065  basename, & ! (in)
1066  varname, & ! (in)
1067  ts, & ! (in)
1068  te, & ! (in)
1069  myrank, & ! (in)
1070  mdlid, & ! (in)
1071  single, & ! (in) optional
1072  xstag, & ! (in) optional
1073  ystag & ! (in) optional
1074  )
1075  use netcdf ![external lib]
1076  implicit none
1077 
1078  real(SP), intent(out) :: var(:,:,:)
1079  character(len=*), intent( in) :: basename
1080  character(len=*), intent( in) :: varname
1081  integer, intent( in) :: ts
1082  integer, intent( in) :: te
1083  integer, intent( in) :: myrank
1084  integer, intent( in) :: mdlid
1085  logical, intent( in), optional :: single
1086  logical, intent( in), optional :: xstag
1087  logical, intent( in), optional :: ystag
1088 
1089  real(SP), allocatable :: var_org(:,:,:)
1090  integer(2), allocatable :: short(:,:,:)
1091 
1092  real(4) :: scale_factor, add_offset
1093 
1094  integer :: ncid, varid
1095  integer :: status
1096  integer :: precis
1097  integer :: nx, ny
1098  integer :: dims(7)
1099 
1100  integer :: tcount
1101  character(len=H_LONG) :: fname = ''
1102  logical :: single_ = .false.
1103 
1104  intrinsic size
1105  intrinsic shape
1106  !---------------------------------------------------------------------------
1107 
1108  tcount = te - ts + 1
1109 
1110  if ( present(single) ) then
1111  single_ = single
1112  else
1113  single_ = .false.
1114  endif
1115 
1116  call externalfilemakefname( fname,mdlid,basename,myrank,single_ )
1117 
1118  status = nf90_open( trim(fname), nf90_nowrite, ncid )
1119  if (status /= nf90_noerr) call handle_err(status, __line__)
1120 
1121  ! retrieve dimension size in data original order
1122  call externaltakedimension( dims(:),ncid,mdlid )
1123  nx = dims(1)
1124  if ( present(xstag) ) then
1125  if ( xstag ) then
1126  nx = dims(4)
1127  endif
1128  endif
1129  ny = dims(2)
1130  if ( present(ystag) ) then
1131  if ( ystag ) then
1132  ny = dims(5)
1133  endif
1134  endif
1135  allocate( var_org(nx,ny,tcount) )
1136  allocate( short(nx,ny,tcount) )
1137 
1138  status = nf90_inq_varid( ncid, trim(varname), varid )
1139  if (status /= nf90_noerr) call handle_err(status, __line__)
1140 
1141  status = nf90_inquire_variable( ncid, varid, xtype=precis )
1142  if(status /= nf90_noerr) call handle_err(status, __line__)
1143 
1144  if(precis /= nf90_short) then
1145  status = nf90_get_var( ncid, varid, var_org(:,:,:), start = (/ 1,1,ts /), &
1146  count = (/ nx,ny,tcount /) )
1147  if (status /= nf90_noerr) call handle_err(status, __line__)
1148  else
1149  status = nf90_get_att(ncid, varid, "scale_factor", scale_factor)
1150  if (status /= nf90_noerr) call handle_err(status, __line__)
1151 
1152  status = nf90_get_att(ncid, varid, "add_offset", add_offset)
1153  if (status /= nf90_noerr) call handle_err(status, __line__)
1154 
1155  status = nf90_get_var( ncid, varid, short(:,:,:), start = (/ 1,1,ts /), &
1156  count = (/ nx,ny,tcount /) )
1157  if (status /= nf90_noerr) call handle_err(status, __line__)
1158 
1159  var_org(:,:,:) = real( short(:,:,:),kind=sp )*scale_factor + add_offset
1160  end if
1161 
1162  status = nf90_close(ncid)
1163  if (status /= nf90_noerr) call handle_err(status, __line__)
1164 
1165  call convertarrayorder( var,var_org,tcount,nx,ny )
1166 
1167  deallocate( var_org )
1168  deallocate( short )
1169 
1170  return
1171  end subroutine externalfilereadoffset3drealsp
1172  subroutine externalfilereadoffset3drealdp( &
1173  var, & ! (out)
1174  basename, & ! (in)
1175  varname, & ! (in)
1176  ts, & ! (in)
1177  te, & ! (in)
1178  myrank, & ! (in)
1179  mdlid, & ! (in)
1180  single, & ! (in) optional
1181  xstag, & ! (in) optional
1182  ystag & ! (in) optional
1183  )
1184  use netcdf ![external lib]
1185  implicit none
1186 
1187  real(DP), intent(out) :: var(:,:,:)
1188  character(len=*), intent( in) :: basename
1189  character(len=*), intent( in) :: varname
1190  integer, intent( in) :: ts
1191  integer, intent( in) :: te
1192  integer, intent( in) :: myrank
1193  integer, intent( in) :: mdlid
1194  logical, intent( in), optional :: single
1195  logical, intent( in), optional :: xstag
1196  logical, intent( in), optional :: ystag
1197 
1198  real(DP), allocatable :: var_org(:,:,:)
1199  integer(2), allocatable :: short(:,:,:)
1200 
1201  real(4) :: scale_factor, add_offset
1202 
1203  integer :: ncid, varid
1204  integer :: status
1205  integer :: precis
1206  integer :: nx, ny
1207  integer :: dims(7)
1208 
1209  integer :: tcount
1210  character(len=H_LONG) :: fname = ''
1211  logical :: single_ = .false.
1212 
1213  intrinsic size
1214  intrinsic shape
1215  !---------------------------------------------------------------------------
1216 
1217  tcount = te - ts + 1
1218 
1219  if ( present(single) ) then
1220  single_ = single
1221  else
1222  single_ = .false.
1223  endif
1224 
1225  call externalfilemakefname( fname,mdlid,basename,myrank,single_ )
1226 
1227  status = nf90_open( trim(fname), nf90_nowrite, ncid )
1228  if (status /= nf90_noerr) call handle_err(status, __line__)
1229 
1230  ! retrieve dimension size in data original order
1231  call externaltakedimension( dims(:),ncid,mdlid )
1232  nx = dims(1)
1233  if ( present(xstag) ) then
1234  if ( xstag ) then
1235  nx = dims(4)
1236  endif
1237  endif
1238  ny = dims(2)
1239  if ( present(ystag) ) then
1240  if ( ystag ) then
1241  ny = dims(5)
1242  endif
1243  endif
1244  allocate( var_org(nx,ny,tcount) )
1245  allocate( short(nx,ny,tcount) )
1246 
1247  status = nf90_inq_varid( ncid, trim(varname), varid )
1248  if (status /= nf90_noerr) call handle_err(status, __line__)
1249 
1250  status = nf90_inquire_variable( ncid, varid, xtype=precis )
1251  if(status /= nf90_noerr) call handle_err(status, __line__)
1252 
1253  if(precis /= nf90_short) then
1254  status = nf90_get_var( ncid, varid, var_org(:,:,:), start = (/ 1,1,ts /), &
1255  count = (/ nx,ny,tcount /) )
1256  if (status /= nf90_noerr) call handle_err(status, __line__)
1257  else
1258  status = nf90_get_att(ncid, varid, "scale_factor", scale_factor)
1259  if (status /= nf90_noerr) call handle_err(status, __line__)
1260 
1261  status = nf90_get_att(ncid, varid, "add_offset", add_offset)
1262  if (status /= nf90_noerr) call handle_err(status, __line__)
1263 
1264  status = nf90_get_var( ncid, varid, short(:,:,:), start = (/ 1,1,ts /), &
1265  count = (/ nx,ny,tcount /) )
1266  if (status /= nf90_noerr) call handle_err(status, __line__)
1267 
1268  var_org(:,:,:) = real( real(short(:,:,:),kind=sp)*scale_factor + add_offset, kind=dp )
1269  endif
1270 
1271  status = nf90_close(ncid)
1272  if (status /= nf90_noerr) call handle_err(status, __line__)
1273 
1274  call convertarrayorder( var,var_org,tcount,nx,ny )
1275 
1276  deallocate( var_org )
1277  deallocate( short )
1278 
1279  return
1280  end subroutine externalfilereadoffset3drealdp
1281  subroutine externalfilereadoffset4drealsp( &
1282  var, & ! (out)
1283  basename, & ! (in)
1284  varname, & ! (in)
1285  ts, & ! (in)
1286  te, & ! (in)
1287  myrank, & ! (in)
1288  mdlid, & ! (in)
1289  single, & ! (in) optional
1290  xstag, & ! (in) optional
1291  ystag, & ! (in) optional
1292  zstag, & ! (in) optional
1293  landgrid & ! (in) optional
1294  )
1295  use netcdf ![external lib]
1296  implicit none
1297 
1298  real(SP), intent(out) :: var(:,:,:,:)
1299  character(len=*), intent( in) :: basename
1300  character(len=*), intent( in) :: varname
1301  integer, intent( in) :: ts
1302  integer, intent( in) :: te
1303  integer, intent( in) :: myrank
1304  integer, intent( in) :: mdlid
1305  logical, intent( in), optional :: single
1306  logical, intent( in), optional :: xstag
1307  logical, intent( in), optional :: ystag
1308  logical, intent( in), optional :: zstag
1309  logical, intent( in), optional :: landgrid
1310 
1311  real(SP), allocatable :: var_org(:,:,:,:)
1312  integer(2), allocatable :: short(:,:,:,:)
1313 
1314  real(4) :: scale_factor, add_offset
1315 
1316  integer :: ncid, varid
1317  integer :: status
1318  integer :: precis
1319  integer :: nx, ny, nz
1320  integer :: dims(7)
1321 
1322  integer :: tcount
1323  character(len=H_LONG) :: fname = ''
1324  logical :: single_ = .false.
1325 
1326  intrinsic size
1327  intrinsic shape
1328  !---------------------------------------------------------------------------
1329 
1330  tcount = te - ts + 1
1331 
1332  if ( present(single) ) then
1333  single_ = single
1334  else
1335  single_ = .false.
1336  endif
1337 
1338  call externalfilemakefname( fname,mdlid,basename,myrank,single_ )
1339 
1340  status = nf90_open( trim(fname), nf90_nowrite, ncid )
1341  if (status /= nf90_noerr) call handle_err(status, __line__)
1342 
1343  ! retrieve dimension size in data original order
1344  call externaltakedimension( dims(:),ncid,mdlid )
1345  nx = dims(1)
1346  if ( present(xstag) ) then
1347  if ( xstag ) then
1348  nx = dims(4)
1349  endif
1350  endif
1351  ny = dims(2)
1352  if ( present(ystag) ) then
1353  if ( ystag ) then
1354  ny = dims(5)
1355  endif
1356  endif
1357  nz = dims(3)
1358  if ( present(zstag) ) then
1359  if ( zstag ) then
1360  nz = dims(6)
1361  endif
1362  endif
1363  if ( present(landgrid) ) then
1364  if ( landgrid ) then
1365  nz = dims(7)
1366  endif
1367  endif
1368  allocate( var_org(nx,ny,nz,tcount) )
1369  allocate( short(nx,ny,nz,tcount) )
1370 
1371  status = nf90_inq_varid( ncid, trim(varname), varid )
1372  if (status /= nf90_noerr) call handle_err(status, __line__)
1373 
1374  status = nf90_inquire_variable( ncid, varid, xtype=precis )
1375  if (status /= nf90_noerr) call handle_err(status, __line__)
1376 
1377  if(precis /= nf90_short) then
1378  status = nf90_get_var( ncid, varid, var_org(:,:,:,:), start = (/ 1,1,1,ts /), &
1379  count = (/ nx,ny,nz,tcount /) )
1380  if (status /= nf90_noerr) call handle_err(status, __line__)
1381  else
1382  status = nf90_get_att(ncid, varid, "scale_factor", scale_factor)
1383  if (status /= nf90_noerr) call handle_err(status, __line__)
1384 
1385  status = nf90_get_att(ncid, varid, "add_offset", add_offset)
1386  if (status /= nf90_noerr) call handle_err(status, __line__)
1387 
1388  status = nf90_get_var( ncid, varid, short(:,:,:,:), start = (/ 1,1,1,ts /), &
1389  count = (/ nx,ny,nz,tcount /) )
1390  if (status /= nf90_noerr) call handle_err(status, __line__)
1391 
1392  var_org(:,:,:,:) = real( short(:,:,:,:),kind=sp )*scale_factor + add_offset
1393  end if
1394 
1395  status = nf90_close(ncid)
1396  if (status /= nf90_noerr) call handle_err(status, __line__)
1397 
1398  call convertarrayorder( var,var_org,tcount,nz,nx,ny )
1399 
1400  deallocate( var_org )
1401  deallocate( short )
1402 
1403  return
1404  end subroutine externalfilereadoffset4drealsp
1405  subroutine externalfilereadoffset4drealdp( &
1406  var, & ! (out)
1407  basename, & ! (in)
1408  varname, & ! (in)
1409  ts, & ! (in)
1410  te, & ! (in)
1411  myrank, & ! (in)
1412  mdlid, & ! (in)
1413  single, & ! (in) optional
1414  xstag, & ! (in) optional
1415  ystag, & ! (in) optional
1416  zstag, & ! (in) optional
1417  landgrid & ! (in) optional
1418  )
1419  use netcdf ![external lib]
1420  implicit none
1421 
1422  real(DP), intent(out) :: var(:,:,:,:)
1423  character(len=*), intent( in) :: basename
1424  character(len=*), intent( in) :: varname
1425  integer, intent( in) :: ts
1426  integer, intent( in) :: te
1427  integer, intent( in) :: myrank
1428  integer, intent( in) :: mdlid
1429  logical, intent( in), optional :: single
1430  logical, intent( in), optional :: xstag
1431  logical, intent( in), optional :: ystag
1432  logical, intent( in), optional :: zstag
1433  logical, intent( in), optional :: landgrid
1434 
1435  real(DP), allocatable :: var_org(:,:,:,:)
1436  integer(2), allocatable :: short(:,:,:,:)
1437 
1438  real(4) :: scale_factor, add_offset
1439 
1440  integer :: ncid, varid
1441  integer :: status
1442  integer :: precis
1443  integer :: nx, ny, nz
1444  integer :: dims(7)
1445 
1446  integer :: tcount
1447  character(len=H_LONG) :: fname = ''
1448  logical :: single_ = .false.
1449 
1450  intrinsic size
1451  intrinsic shape
1452  !---------------------------------------------------------------------------
1453 
1454  tcount = te - ts + 1
1455 
1456  if ( present(single) ) then
1457  single_ = single
1458  else
1459  single_ = .false.
1460  endif
1461 
1462  call externalfilemakefname( fname,mdlid,basename,myrank,single_ )
1463 
1464  status = nf90_open( trim(fname), nf90_nowrite, ncid )
1465  if (status /= nf90_noerr) call handle_err(status, __line__)
1466 
1467  ! retrieve dimension size in data original order
1468  call externaltakedimension( dims(:),ncid,mdlid )
1469  nx = dims(1)
1470  if ( present(xstag) ) then
1471  if ( xstag ) then
1472  nx = dims(4)
1473  endif
1474  endif
1475  ny = dims(2)
1476  if ( present(ystag) ) then
1477  if ( ystag ) then
1478  ny = dims(5)
1479  endif
1480  endif
1481  nz = dims(3)
1482  if ( present(zstag) ) then
1483  if ( zstag ) then
1484  nz = dims(6)
1485  endif
1486  endif
1487  if ( present(landgrid) ) then
1488  if ( landgrid ) then
1489  nz = dims(7)
1490  endif
1491  endif
1492  allocate( var_org(nx,ny,nz,tcount) )
1493  allocate( short(nx,ny,nz,tcount) )
1494 
1495  status = nf90_inq_varid( ncid, trim(varname), varid )
1496  if (status /= nf90_noerr) call handle_err(status, __line__)
1497 
1498  status = nf90_inquire_variable( ncid, varid, xtype=precis )
1499  if(status /= nf90_noerr) call handle_err(status, __line__)
1500 
1501  if(precis /= nf90_short) then
1502  status = nf90_get_var( ncid, varid, var_org(:,:,:,:), start = (/ 1,1,1,ts /), &
1503  count = (/ nx,ny,nz,tcount /) )
1504  if (status /= nf90_noerr) call handle_err(status, __line__)
1505  else
1506  status = nf90_get_att(ncid, varid, "scale_factor", scale_factor)
1507  if (status /= nf90_noerr) call handle_err(status, __line__)
1508 
1509  status = nf90_get_att(ncid, varid, "add_offset", add_offset)
1510  if (status /= nf90_noerr) call handle_err(status, __line__)
1511 
1512  status = nf90_get_var( ncid, varid, short(:,:,:,:), start = (/ 1,1,1,ts /), &
1513  count = (/ nx,ny,nz,tcount /) )
1514  if (status /= nf90_noerr) call handle_err(status, __line__)
1515 
1516  var_org(:,:,:,:) = real( real(short(:,:,:,:),kind=sp)*scale_factor + add_offset, kind=dp )
1517  end if
1518 
1519  status = nf90_close(ncid)
1520  if (status /= nf90_noerr) call handle_err(status, __line__)
1521 
1522  call convertarrayorder( var,var_org,tcount,nz,nx,ny )
1523 
1524  deallocate( var_org )
1525  deallocate( short )
1526 
1527  return
1528  end subroutine externalfilereadoffset4drealdp
1529 
1530  !-----------------------------------------------------------------------------
1531  ! ExternalMakeFName
1532  !-----------------------------------------------------------------------------
1533  subroutine externalfilemakefname( &
1534  fname, & ! (out)
1535  mdlid, & ! (in)
1536  basename, & ! (in)
1537  myrank, & ! (in)
1538  single & ! (in) optional
1539  )
1540  implicit none
1541 
1542  character(len=*), intent(out) :: fname
1543  integer, intent( in) :: mdlid
1544  character(len=*), intent( in) :: basename
1545  integer, intent( in) :: myrank
1546  logical, intent( in) :: single
1547 
1548  intrinsic size
1549  intrinsic shape
1550  !---------------------------------------------------------------------------
1551 
1552  if( mdlid == iwrfarw )then !TYPE: WRF-ARW
1553  if ( single ) then
1554  fname = trim(basename)
1555  else
1556  call filemakefname(fname,trim(basename),'_',myrank,4)
1557  endif
1558  elseif( mdlid == inicam )then !TYPE: NICAM-NETCDF
1559  if ( single ) then
1560  fname = trim(basename)//'.peall.nc'
1561  else
1562  call filemakefname(fname,trim(basename),'anl.pe',myrank,6)
1563  endif
1564  !elseif( mdlid == iGrADS )then !TYPE: GrADS
1565  ! if ( single ) then
1566  ! fname = trim(basename)//'.anl'
1567  ! else
1568  ! call FileMakeFname(fname,trim(basename),'anl.pe',myrank,6)
1569  ! endif
1570  else
1571  write(*,*) 'xxx failed, wrong filetype: [scale_external_io]/[ExternalFileMakeFName]'
1572  call prc_mpistop
1573  endif
1574 
1575  return
1576  end subroutine externalfilemakefname
1577 
1578  !-----------------------------------------------------------------------------
1579  ! ExternalMakeFName
1580  !-----------------------------------------------------------------------------
1581  subroutine externaltakedimension( &
1582  dims, & ! (out)
1583  ncid, & ! (in)
1584  mdlid & ! (in)
1585  )
1586  use netcdf ![external lib]
1587  implicit none
1588 
1589  integer, intent(out) :: dims(:)
1590  integer, intent( in) :: ncid
1591  integer, intent( in) :: mdlid
1592 
1593  integer :: dimid
1594  integer :: status
1595 
1596  intrinsic size
1597  intrinsic shape
1598  !---------------------------------------------------------------------------
1599 
1600  if( mdlid == iwrfarw )then !MODEL ID: WRF-ARW
1601  status = nf90_inq_dimid( ncid, "west_east", dimid )
1602  status = nf90_inquire_dimension( ncid, dimid,len=dims(1) )
1603  status = nf90_inq_dimid( ncid, "south_north", dimid )
1604  status = nf90_inquire_dimension( ncid, dimid,len=dims(2) )
1605  status = nf90_inq_dimid( ncid, "bottom_top", dimid )
1606  status = nf90_inquire_dimension( ncid, dimid,len=dims(3) )
1607  status = nf90_inq_dimid( ncid, "west_east_stag", dimid )
1608  status = nf90_inquire_dimension( ncid, dimid,len=dims(4) )
1609  status = nf90_inq_dimid( ncid, "south_north_stag", dimid )
1610  status = nf90_inquire_dimension( ncid, dimid,len=dims(5) )
1611  status = nf90_inq_dimid( ncid, "bottom_top_stag", dimid )
1612  status = nf90_inquire_dimension( ncid, dimid,len=dims(6) )
1613  status = nf90_inq_dimid( ncid, "soil_layers_stag", dimid )
1614  status = nf90_inquire_dimension( ncid, dimid,len=dims(7) )
1615 
1616  elseif( mdlid == inicam )then !MODEL ID: NICAM-NETCDF
1617  status = nf90_inq_dimid( ncid, "lon", dimid )
1618  status = nf90_inquire_dimension( ncid, dimid,len=dims(1) )
1619  status = nf90_inq_dimid( ncid, "lat", dimid )
1620  status = nf90_inquire_dimension( ncid, dimid,len=dims(2) )
1621  status = nf90_inq_dimid( ncid, "lev", dimid )
1622  status = nf90_inquire_dimension( ncid, dimid,len=dims(3) )
1623 
1624  else
1625  write(*,*) 'xxx This external file format is not supported, Sorry.'
1626  call prc_mpistop
1627  endif
1628 
1629  return
1630  end subroutine externaltakedimension
1631 
1632  !-----------------------------------------------------------------------------
1633  subroutine convertarrayorderwrf2dsp( &
1634  var, & ! (out)
1635  var_org, & ! (in)
1636  tcount, & ! (in)
1637  nx & ! (in)
1638  )
1639  implicit none
1640 
1641  real(SP), intent(out) :: var(:,:)
1642  real(SP), intent( in) :: var_org(:,:)
1643  integer, intent( in) :: tcount
1644  integer, intent( in) :: nx
1645  integer :: n, i
1646  intrinsic shape
1647 
1648  do n = 1, tcount
1649  do i = 1, nx
1650  var(i,n) = var_org(i,n)
1651  end do
1652  end do
1653 
1654  return
1655  end subroutine convertarrayorderwrf2dsp
1656  subroutine convertarrayorderwrf2ddp( &
1657  var, & ! (out)
1658  var_org, & ! (in)
1659  tcount, & ! (in)
1660  nx & ! (in)
1661  )
1662  implicit none
1663 
1664  real(DP), intent(out) :: var(:,:)
1665  real(DP), intent( in) :: var_org(:,:)
1666  integer, intent( in) :: tcount
1667  integer, intent( in) :: nx
1668  integer :: n, i
1669  intrinsic shape
1670 
1671  do n = 1, tcount
1672  do i = 1, nx
1673  var(i,n) = var_org(i,n)
1674  end do
1675  end do
1676 
1677  return
1678  end subroutine convertarrayorderwrf2ddp
1679  subroutine convertarrayorderwrf3dsp( &
1680  var, & ! (out)
1681  var_org, & ! (in)
1682  tcount, & ! (in)
1683  nx, & ! (in)
1684  ny & ! (in)
1685  )
1686  implicit none
1687 
1688  real(SP), intent(out) :: var(:,:,:)
1689  real(SP), intent( in) :: var_org(:,:,:)
1690  integer, intent( in) :: tcount
1691  integer, intent( in) :: nx
1692  integer, intent( in) :: ny
1693  integer :: n, i, j
1694  intrinsic shape
1695 
1696  do n = 1, tcount
1697  do j = 1, ny
1698  do i = 1, nx
1699  var(i,j,n) = var_org(i,j,n)
1700  end do
1701  end do
1702  end do
1703 
1704  return
1705  end subroutine convertarrayorderwrf3dsp
1706  subroutine convertarrayorderwrf3ddp( &
1707  var, & ! (out)
1708  var_org, & ! (in)
1709  tcount, & ! (in)
1710  nx, & ! (in)
1711  ny & ! (in)
1712  )
1713  implicit none
1714 
1715  real(DP), intent(out) :: var(:,:,:)
1716  real(DP), intent( in) :: var_org(:,:,:)
1717  integer, intent( in) :: tcount
1718  integer, intent( in) :: nx
1719  integer, intent( in) :: ny
1720  integer :: n, i, j
1721  intrinsic shape
1722 
1723  do n = 1, tcount
1724  do j = 1, ny
1725  do i = 1, nx
1726  var(i,j,n) = var_org(i,j,n)
1727  end do
1728  end do
1729  end do
1730 
1731  return
1732  end subroutine convertarrayorderwrf3ddp
1733  subroutine convertarrayorderwrf4dsp( &
1734  var, & ! (out)
1735  var_org, & ! (in)
1736  tcount, & ! (in)
1737  nz, & ! (in)
1738  nx, & ! (in)
1739  ny & ! (in)
1740  )
1741  implicit none
1742 
1743  real(SP), intent(out) :: var(:,:,:,:)
1744  real(SP), intent( in) :: var_org(:,:,:,:)
1745  integer, intent( in) :: tcount
1746  integer, intent( in) :: nz
1747  integer, intent( in) :: nx
1748  integer, intent( in) :: ny
1749  integer :: n, k, i, j
1750  intrinsic shape
1751 
1752  do n = 1, tcount
1753  do j = 1, ny
1754  do i = 1, nx
1755  do k = 1, nz
1756  var(k,i,j,n) = var_org(i,j,k,n)
1757  end do
1758  end do
1759  end do
1760  end do
1761 
1762  return
1763  end subroutine convertarrayorderwrf4dsp
1764  subroutine convertarrayorderwrf4ddp( &
1765  var, & ! (out)
1766  var_org, & ! (in)
1767  tcount, & ! (in)
1768  nz, & ! (in)
1769  nx, & ! (in)
1770  ny & ! (in)
1771  )
1772  implicit none
1773 
1774  real(DP), intent(out) :: var(:,:,:,:)
1775  real(DP), intent( in) :: var_org(:,:,:,:)
1776  integer, intent( in) :: tcount
1777  integer, intent( in) :: nz
1778  integer, intent( in) :: nx
1779  integer, intent( in) :: ny
1780  integer :: n, k, i, j
1781  intrinsic shape
1782 
1783  do n = 1, tcount
1784  do j = 1, ny
1785  do i = 1, nx
1786  do k = 1, nz
1787  var(k,i,j,n) = var_org(i,j,k,n)
1788  end do
1789  end do
1790  end do
1791  end do
1792 
1793  return
1794  end subroutine convertarrayorderwrf4ddp
1795 
1796  !-----------------------------------------------------------------------------
1797  subroutine handle_err(status, line)
1798  use netcdf ![external lib]
1799  use scale_process, only: &
1800  prc_mpistop
1801  implicit none
1802  integer, intent(in) :: status
1803  integer, intent(in) :: line
1804 
1805  write(*,*) 'xxx Error in scale_external_io.f90 at line', line
1806  write(*,*) nf90_strerror(status)
1807  call prc_mpistop
1808 
1809  return
1810  end subroutine handle_err
1811 
1812 end module scale_external_io
module GTOOL_FILE
Definition: gtool_file.f90:17
subroutine externalfileread2drealdp(var, basename, varname, ts, te, myrank, mdlid, nx, single)
subroutine, public externalfilevarexistence(existence, basename, varname, myrank, mdlid, single)
Check Existence of a Variable.
subroutine, public prc_mpistop
Abort MPI.
integer, parameter, public inicam
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:61
subroutine convertarrayorderwrf2dsp(var, var_org, tcount, nx)
module STDIO
Definition: scale_stdio.F90:12
integer, parameter, public igrads
integer, parameter, public iwrfarw
integer, parameter, public iscale
subroutine externalfileread3drealdp(var, basename, varname, ts, te, myrank, mdlid, single, xstag, ystag)
subroutine externalfileread2drealsp(var, basename, varname, ts, te, myrank, mdlid, nx, single)
File Read.
module FILE I/O (netcdf)
subroutine externalfileread4drealdp(var, basename, varname, ts, te, myrank, mdlid, single, xstag, ystag, zstag, landgrid)
module grid index
subroutine externalfilereadoffset3drealdp(var, basename, varname, ts, te, myrank, mdlid, single, xstag, ystag)
subroutine, public externalfilegetshape(dims, timelen, mdlid, basename, myrank, single)
subroutine, public filemakefname(fname, basename, prefix, myrank, len)
module PROCESS
subroutine externalfilegetglobalattvrealsp(var, mdlid, basename, attname, myrank, single)
subroutine, public externalfilegetglobalattc(chr, mdlid, basename, attname, myrank, single)
subroutine externalfileread4drealsp(var, basename, varname, ts, te, myrank, mdlid, single, xstag, ystag, zstag, landgrid, option)
subroutine externalfilereadoffset3drealsp(var, basename, varname, ts, te, myrank, mdlid, single, xstag, ystag)
subroutine externalfilereadoffset4drealsp(var, basename, varname, ts, te, myrank, mdlid, single, xstag, ystag, zstag, landgrid)
module PRECISION
integer, parameter, public sp
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
subroutine externalfileread3drealsp(var, basename, varname, ts, te, myrank, mdlid, single, xstag, ystag, option)
subroutine externalfilereadoffset4drealdp(var, basename, varname, ts, te, myrank, mdlid, single, xstag, ystag, zstag, landgrid)
subroutine externalfilegetglobalattvrealdp(var, mdlid, basename, attname, myrank, single)
integer, parameter, public dp
subroutine externalfilegetglobalattvinteger(var, mdlid, basename, attname, myrank, single)