00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039 module shr_ncread_mod
00040
00041
00042
00043 use shr_string_mod
00044 use shr_kind_mod
00045 use shr_sys_mod
00046 use shr_file_mod
00047 use shr_cal_mod
00048 use shr_log_mod, only: s_loglev => shr_log_Level
00049 use shr_log_mod, only: s_logunit => shr_log_Unit
00050 use netcdf
00051
00052 implicit none
00053
00054 private
00055
00056
00057
00058
00059
00060
00061
00062 public :: shr_ncread_varExists
00063 public :: shr_ncread_attExists
00064 public :: shr_ncread_varDimNum
00065 public :: shr_ncread_varDimSize
00066 public :: shr_ncread_varDimSizes
00067 public :: shr_ncread_dimSize
00068 public :: shr_ncread_attribute
00069 public :: shr_ncread_tCoord
00070 public :: shr_ncread_domain
00071 public :: shr_ncread_tField
00072 public :: shr_ncread_Field4dG
00073 public :: shr_ncread_print
00074 public :: shr_ncread_setAbort
00075 public :: shr_ncread_setDebug
00076 public :: shr_ncread_open
00077 public :: shr_ncread_close
00078
00079
00080
00081
00082
00083
00084
00085 interface shr_ncread_varDimSize ; module procedure &
00086 shr_ncread_varDimSizeName, &
00087 shr_ncread_varDimSizeID
00088 end interface
00089
00090 interface shr_ncread_dimSize ; module procedure &
00091 shr_ncread_dimSizeName
00092 end interface
00093
00094 interface shr_ncread_tCoord ; module procedure &
00095 shr_ncread_tCoordRC, &
00096 shr_ncread_tCoordII
00097 end interface
00098
00099 interface shr_ncread_tField ; module procedure &
00100 shr_ncread_tField2dR8, &
00101 shr_ncread_tField1dR8, &
00102 shr_ncread_tField2dIN, &
00103 shr_ncread_tField1dIN
00104 end interface
00105
00106 logical ,save :: doabort = .true.
00107 integer(SHR_KIND_IN),save :: debug = 0
00108
00109
00110 contains
00111
00112
00113
00114
00115
00116
00117
00118
00119
00120
00121
00122
00123
00124
00125
00126
00127
00128
00129 logical function shr_ncread_varExists(fileName, varName)
00130
00131 implicit none
00132
00133
00134
00135 character(*),intent(in) :: fileName
00136 character(*),intent(in) :: varName
00137
00138
00139
00140
00141 integer(SHR_KIND_IN) :: fid
00142 integer(SHR_KIND_IN) :: vid
00143 integer(SHR_KIND_IN) :: debug0
00144 integer(SHR_KIND_IN) :: rCode
00145
00146
00147 character(*),parameter :: subName = "(shr_ncread_varExists)"
00148 character(*),parameter :: F00 = "('(shr_ncread_varExists) ',4a)"
00149 character(*),parameter :: F01 = "('(shr_ncread_varExists) ',a,i6)"
00150
00151
00152
00153
00154
00155
00156 debug0 = debug
00157 call shr_ncread_setDebug(0)
00158
00159 shr_ncread_varExists = .false.
00160 call shr_ncread_open(fileName,fid,rCode)
00161 rcode = nf90_inq_varid(fid,trim(varName),vid)
00162 if (rcode == nf90_noerr) shr_ncread_varExists = .true.
00163 call shr_ncread_close(fid,rCode)
00164
00165
00166 call shr_ncread_setDebug(debug0)
00167
00168 end function shr_ncread_varExists
00169
00170
00171
00172
00173
00174
00175
00176
00177
00178
00179
00180
00181
00182
00183
00184
00185
00186 logical function shr_ncread_attExists(fileName, varName, attName)
00187
00188 implicit none
00189
00190
00191
00192 character(*),intent(in) :: fileName
00193 character(*),intent(in) :: varName
00194 character(*),intent(in) :: attName
00195
00196
00197
00198
00199 integer(SHR_KIND_IN) :: fid
00200 integer(SHR_KIND_IN) :: vid
00201 integer(SHR_KIND_IN) :: debug0
00202 integer(SHR_KIND_IN) :: rCode
00203
00204
00205 character(*),parameter :: subName = "(shr_ncread_attExists)"
00206 character(*),parameter :: F00 = "('(shr_ncread_attExists) ',4a)"
00207 character(*),parameter :: F01 = "('(shr_ncread_attExists) ',a,i6)"
00208
00209
00210
00211
00212
00213
00214 debug0 = debug
00215 call shr_ncread_setDebug(0)
00216
00217 shr_ncread_attExists = .false.
00218 call shr_ncread_open(fileName,fid,rCode)
00219 rCode = nf90_inq_varid(fid,trim(varName),vid)
00220 if (rCode == nf90_noerr) then
00221 rCode = nf90_inquire_attribute(fid, vid, trim(attName))
00222 if (rCode == nf90_noerr) shr_ncread_attExists = .true.
00223 endif
00224
00225 call shr_ncread_close(fid,rCode)
00226
00227
00228 call shr_ncread_setDebug(debug0)
00229
00230 end function shr_ncread_attExists
00231
00232
00233
00234
00235
00236
00237
00238
00239
00240
00241
00242
00243
00244
00245
00246
00247 subroutine shr_ncread_varDimNum(fileName, varName, ns, rc)
00248
00249 implicit none
00250
00251
00252
00253 character(*) ,intent(in) :: fileName
00254 character(*) ,intent(in) :: varName
00255 integer(SHR_KIND_IN),intent(out) :: ns
00256 integer(SHR_KIND_IN),intent(out),optional :: rc
00257
00258
00259
00260
00261 integer(SHR_KIND_IN) :: fid
00262 integer(SHR_KIND_IN) :: vid
00263 integer(SHR_KIND_IN) :: rCode
00264
00265
00266 character(*),parameter :: subName = "(shr_ncread_varDimNum)"
00267 character(*),parameter :: F00 = "('(shr_ncread_varDimNum) ',4a)"
00268 character(*),parameter :: F01 = "('(shr_ncread_varDimNum) ',a,i6)"
00269
00270
00271
00272
00273
00274 call shr_ncread_open(fileName,fid,rCode)
00275
00276
00277 rcode = nf90_inq_varid(fid,trim(varName),vid)
00278 call shr_ncread_handleErr(rCode, subName//" ERROR inq varid")
00279 rcode = nf90_inquire_variable(fid,vid,ndims=ns)
00280 call shr_ncread_handleErr(rCode, subName//" ERROR inq var")
00281 if (debug > 1 .and. s_loglev > 0) write(s_logunit,F01) trim(varName)//' has dims = ',ns
00282
00283 call shr_ncread_close(fid,rCode)
00284
00285 if (present(rc)) rc = rCode
00286
00287 end subroutine shr_ncread_varDimNum
00288
00289
00290
00291
00292
00293
00294
00295
00296
00297
00298
00299
00300
00301
00302
00303
00304
00305 subroutine shr_ncread_varDimSizeName(fileName, varName, dimName, ns, rc)
00306
00307 implicit none
00308
00309
00310
00311 character(*) ,intent(in) :: fileName
00312 character(*) ,intent(in) :: varName
00313 character(*) ,intent(in) :: dimName
00314 integer(SHR_KIND_IN),intent(out) :: ns
00315 integer(SHR_KIND_IN),intent(out),optional :: rc
00316
00317
00318
00319
00320 integer(SHR_KIND_IN) :: rCode
00321
00322
00323 character(*),parameter :: subName = "(shr_ncread_varDimSizeName)"
00324 character(*),parameter :: F00 = "('(shr_ncread_varDimSizeName) ',4a)"
00325 character(*),parameter :: F01 = "('(shr_ncread_varDimSizeName) ',a,i6)"
00326
00327
00328
00329
00330
00331 call shr_ncread_dimSizeName(fileName,dimName,ns,rCode)
00332
00333 if (present(rc)) rc = rCode
00334
00335 end subroutine shr_ncread_varDimSizeName
00336
00337
00338
00339
00340
00341
00342
00343
00344
00345
00346
00347
00348
00349
00350
00351
00352
00353 subroutine shr_ncread_varDimSizeID(fileName, varName, dnum, ns, rc)
00354
00355 implicit none
00356
00357
00358
00359 character(*) ,intent(in) :: fileName
00360 character(*) ,intent(in) :: varName
00361 integer(SHR_KIND_IN),intent(in) :: dnum
00362 integer(SHR_KIND_IN),intent(out) :: ns
00363 integer(SHR_KIND_IN),intent(out),optional :: rc
00364
00365
00366
00367
00368 integer(SHR_KIND_IN) :: fid
00369 integer(SHR_KIND_IN) :: vid
00370 integer(SHR_KIND_IN) :: ndims
00371 character(SHR_KIND_CS) :: dimName
00372 integer(SHR_KIND_IN),allocatable :: dids(:)
00373 integer(SHR_KIND_IN) :: rCode
00374
00375
00376 character(*),parameter :: subName = "(shr_ncread_varDimSizeID)"
00377 character(*),parameter :: F00 = "('(shr_ncread_varDimSizeID) ',4a)"
00378 character(*),parameter :: F01 = "('(shr_ncread_varDimSizeID) ',a,i6)"
00379
00380
00381
00382
00383
00384 call shr_ncread_open(fileName,fid,rCode)
00385
00386 rCode = nf90_inq_varid(fid,trim(varName),vid)
00387 call shr_ncread_handleErr(rCode,subName//' ERROR inq varid vid')
00388 rCode = nf90_inquire_variable(fid,vid,ndims=ndims)
00389 call shr_ncread_handleErr(rCode,subName//' ERROR inquire variable ndims')
00390 allocate(dids(ndims))
00391 rCode = nf90_inquire_variable(fid,vid,dimids=dids)
00392 call shr_ncread_handleErr(rCode,subName//' ERROR inquire variable dimids')
00393 rcode = nf90_inquire_dimension(fid,dids(dnum),name=dimName,len=ns)
00394 call shr_ncread_handleErr(rCode, subName//" ERROR inquire dimension")
00395 if (debug > 1 .and. s_loglev > 0) write(s_logunit,F01) trim(dimName)//' dimension has size = ',ns
00396
00397 deallocate(dids)
00398 call shr_ncread_close(fid,rCode)
00399
00400 if (present(rc)) rc = rCode
00401
00402 end subroutine shr_ncread_varDimSizeID
00403
00404
00405
00406
00407
00408
00409
00410
00411
00412
00413
00414
00415
00416
00417
00418
00419
00420 subroutine shr_ncread_varDimSizes(fileName, varName, n1, n2, n3, n4, n5, n6, rc)
00421
00422 implicit none
00423
00424
00425
00426 character(*) ,intent(in) :: fileName
00427 character(*) ,intent(in) :: varName
00428 integer(SHR_KIND_IN),intent(out),optional :: n1
00429 integer(SHR_KIND_IN),intent(out),optional :: n2
00430 integer(SHR_KIND_IN),intent(out),optional :: n3
00431 integer(SHR_KIND_IN),intent(out),optional :: n4
00432 integer(SHR_KIND_IN),intent(out),optional :: n5
00433 integer(SHR_KIND_IN),intent(out),optional :: n6
00434 Integer(SHR_KIND_IN),intent(out),optional :: rc
00435
00436
00437
00438
00439 integer(SHR_KIND_IN),parameter :: maxn = 6
00440 integer(SHR_KIND_IN) :: n
00441 integer(SHR_KIND_IN) :: fid
00442 integer(SHR_KIND_IN) :: vid
00443 integer(SHR_KIND_IN) :: ndims
00444 integer(SHR_KIND_IN),allocatable :: dids(:)
00445 integer(SHR_KIND_IN),allocatable :: ns(:)
00446 integer(SHR_KIND_IN) :: rCode
00447
00448
00449 character(*),parameter :: subName = "(shr_ncread_varDimSizes)"
00450 character(*),parameter :: F00 = "('(shr_ncread_varDimSizes) ',4a)"
00451 character(*),parameter :: F01 = "('(shr_ncread_varDimSizes) ',a,i6)"
00452
00453
00454
00455
00456
00457 call shr_ncread_open(fileName,fid,rCode)
00458
00459 rCode = nf90_inq_varid(fid,trim(varName),vid)
00460 call shr_ncread_handleErr(rCode,subName//' ERROR inq varid vid')
00461 rCode = nf90_inquire_variable(fid,vid,ndims=ndims)
00462 call shr_ncread_handleErr(rCode,subName//' ERROR inquire variable ndims')
00463 allocate(dids(ndims))
00464 allocate(ns(maxn))
00465 rCode = nf90_inquire_variable(fid,vid,dimids=dids)
00466 call shr_ncread_handleErr(rCode,subName//' ERROR inquire variable dimids')
00467
00468
00469 ns = 1
00470 do n=1,min(ndims,maxn)
00471 rcode = nf90_inquire_dimension(fid,dids(n),len=ns(n))
00472 call shr_ncread_handleErr(rCode, subName//" ERROR inquire dimension")
00473 enddo
00474
00475 call shr_ncread_close(fid,rCode)
00476
00477
00478 if (present(n1)) then
00479 n1 = ns(1)
00480 endif
00481 if (present(n2)) then
00482 n2 = ns(2)
00483 endif
00484 if (present(n3)) then
00485 n3 = ns(3)
00486 endif
00487 if (present(n4)) then
00488 n4 = ns(4)
00489 endif
00490 if (present(n5)) then
00491 n5 = ns(5)
00492 endif
00493 if (present(n6)) then
00494 n6 = ns(6)
00495 endif
00496
00497 deallocate(dids)
00498 deallocate(ns)
00499
00500 if (present(rc)) rc = rCode
00501
00502 end subroutine shr_ncread_varDimSizes
00503
00504
00505
00506
00507
00508
00509
00510
00511
00512
00513
00514
00515
00516
00517
00518
00519 subroutine shr_ncread_dimSizeName(fileName, dimName, ns, rc)
00520
00521 implicit none
00522
00523
00524
00525 character(*) ,intent(in) :: fileName
00526 character(*) ,intent(in) :: dimName
00527 integer(SHR_KIND_IN),intent(out) :: ns
00528 integer(SHR_KIND_IN),intent(out),optional :: rc
00529
00530
00531
00532
00533 integer(SHR_KIND_IN) :: fid
00534 integer(SHR_KIND_IN) :: did
00535 integer(SHR_KIND_IN) :: rCode
00536
00537
00538 character(*),parameter :: subName = "(shr_ncread_dimSizeName)"
00539 character(*),parameter :: F00 = "('(shr_ncread_dimSizeName) ',4a)"
00540 character(*),parameter :: F01 = "('(shr_ncread_dimSizeName) ',a,i6)"
00541
00542
00543
00544
00545
00546 call shr_ncread_open(fileName,fid,rCode)
00547
00548
00549 rcode = nf90_inq_dimid (fid, trim(dimName), did)
00550 call shr_ncread_handleErr(rCode, subName//" ERROR inq dimid")
00551 rcode = nf90_inquire_dimension(fid,did,len=ns)
00552 call shr_ncread_handleErr(rCode, subName//" ERROR inquire dimension")
00553 if (debug > 1 .and. s_loglev > 0) write(s_logunit,F01) trim(dimName)//' dimension has size = ',ns
00554
00555 call shr_ncread_close(fid,rCode)
00556
00557 if (present(rc)) rc = rCode
00558
00559 end subroutine shr_ncread_dimSizeName
00560
00561
00562
00563
00564
00565
00566
00567
00568
00569
00570
00571
00572
00573
00574
00575
00576
00577
00578 subroutine shr_ncread_attribute(fn,vName,aName,attrib,rc)
00579
00580 implicit none
00581
00582
00583
00584 character(*) ,intent(in) :: fn
00585 character(*) ,intent(in) :: vName
00586 character(*) ,intent(in) :: aName
00587 character(*) ,intent(out) :: attrib
00588 integer(SHR_KIND_IN),intent(out),optional :: rc
00589
00590
00591
00592
00593 integer(SHR_KIND_IN) :: n
00594 integer(SHR_KIND_IN) :: xtype
00595 integer(SHR_KIND_IN) :: len
00596 integer(SHR_KIND_IN) :: fid
00597 integer(SHR_KIND_IN) :: vid
00598 integer(SHR_KIND_IN) :: rCode
00599
00600
00601 character(*),parameter :: subName = "(shr_ncread_attribute)"
00602 character(*),parameter :: F00 = "('(shr_ncread_attribute) ',8a)"
00603 character(*),parameter :: eos = "[end-of-string]"
00604
00605
00606
00607
00608
00609
00610
00611
00612
00613
00614
00615
00616
00617 attrib = ' '
00618 rCode = 0
00619 call shr_ncread_open(fn,fid,rCode)
00620
00621 rCode = nf90_inq_varid(fid, trim(vName), vid)
00622 call shr_ncread_handleErr(rCode,subName//' nf90_inq_var')
00623
00624 rCode = nf90_inquire_attribute(fid, vid, trim(aName), xtype, len)
00625 call shr_ncread_handleErr(rCode,subName//' nf90_inq_att')
00626 if (xtype == NF90_CHAR) then
00627 rCode = nf90_get_att(fid, vid, trim(aName), attrib)
00628 n = len_trim(attrib)
00629 if (ichar(attrib(n:n)) == 0 ) then
00630 if (debug>0 .and. s_loglev > 0) then
00631 write(s_logunit,F00) 'removed null char from end of attribute...'
00632 write(s_logunit,F00) 'orig: ',trim(vName),':',trim(aName),' = ',trim(attrib),eos
00633 end if
00634 attrib(n:n) = ' '
00635 if (debug>0 .and. s_loglev > 0) then
00636 write(s_logunit,F00) 'new : ',trim(vName),':',trim(aName),' = ',trim(attrib),eos
00637 end if
00638 else
00639 if (debug>0 .and. s_loglev > 0) then
00640 write(s_logunit,F00) 'read: ',trim(vName),':',trim(aName),' = ',trim(attrib),eos
00641 end if
00642 end if
00643 call shr_ncread_handleErr(rCode,subName//' nf90_get_att attrib')
00644 else
00645 write(s_logunit,F00) 'attribute: '//trim(vName)//' '//trim(aName)//' not char'
00646 call shr_ncread_abort(subName//' attribute '//trim(aName)//' not char')
00647 endif
00648
00649 call shr_ncread_close(fid,rCode)
00650
00651 if (present(rc)) rc = rCode
00652
00653 end subroutine shr_ncread_attribute
00654
00655
00656
00657
00658
00659
00660
00661
00662
00663
00664
00665
00666
00667
00668
00669
00670
00671
00672 subroutine shr_ncread_tCoordRC(fn, tName, tvar, units, calendar, rc)
00673
00674 implicit none
00675
00676
00677
00678 character(*) ,intent(in) :: fn
00679 character(*) ,intent(in) :: tName
00680 real(SHR_KIND_R8) ,intent(out) :: tvar(:)
00681 character(*) ,intent(out) :: units
00682 character(*) ,intent(out) :: calendar
00683 integer(SHR_KIND_IN),intent(out),optional:: rc
00684
00685
00686
00687 real(SHR_KIND_R8),allocatable :: A1d(:)
00688 character(SHR_KIND_CL) :: string
00689 integer(SHR_KIND_IN) :: i
00690 integer(SHR_KIND_IN) :: ndim,nd1,pd1
00691 integer(SHR_KIND_IN) :: rCode
00692
00693
00694 character(*),parameter :: subName = "(shr_ncread_tCoordRC)"
00695 character(*),parameter :: F00 = "('(shr_ncread_tCoordRC) ',4a)"
00696 character(*),parameter :: F01 = "('(shr_ncread_tCoordRC) ',2a,3i6,2x,a)"
00697 character(*),parameter :: F02 = "('(shr_ncread_tCoordRC) ',a,i6)"
00698 character(*),parameter :: F03 = "('(shr_ncread_tCoordRC) ',a,2i6)"
00699 character(*),parameter :: F04 = "('(shr_ncread_tCoordRC) ',a,2g17.8)"
00700
00701
00702
00703
00704
00705 rCode = 0
00706
00707 if (.not.shr_ncread_varExists(fn,tName)) &
00708 call shr_ncread_abort(subName//' ERROR var does not exist '//trim(tName))
00709
00710
00711 pd1 = size(tvar,1)
00712
00713
00714 call shr_ncread_varDimNum(fn,tName,ndim)
00715 if (ndim /= 1) then
00716 write(s_logunit,F02) 'ERROR '//trim(tName)//' ndim = ',ndim
00717 call shr_ncread_abort(subName//' ERROR ndim must be 1 for '//trim(tName))
00718 endif
00719 call shr_ncread_varDimSize(fn,tName,1,nd1)
00720
00721
00722 if ( nd1 > pd1) then
00723 write(s_logunit,F03) ' nd1 pd1 error ',nd1,pd1
00724 call shr_ncread_abort(subName//' ERROR nd1 pd1 error')
00725 endif
00726
00727 call shr_ncread_attribute(fn,tName,'units',units,rc=rCode)
00728 if (shr_ncread_attExists(fn,tName,'calendar')) then
00729 call shr_ncread_attribute(fn,tName,'calendar',calendar,rc=rCode)
00730 else
00731 calendar = "gregorian"
00732 endif
00733
00734 call shr_string_leftAlign(units)
00735 call shr_string_leftAlign(calendar)
00736
00737 allocate(A1d(nd1))
00738 call shr_ncread_tfield(fn,1,tName,A1d,rc=rCode)
00739 do i=1,nd1
00740 tvar(i) = A1d(i)
00741 enddo
00742
00743 deallocate(A1d)
00744
00745 if (debug > 1 .and. s_loglev > 0) then
00746 write(s_logunit,F04) 'min/max tvar ',minval(tvar),maxval(tvar)
00747 endif
00748
00749 if (present(rc)) rc = rCode
00750
00751 end subroutine shr_ncread_tCoordRC
00752
00753
00754
00755
00756
00757
00758
00759
00760
00761
00762
00763
00764
00765
00766
00767
00768
00769
00770
00771
00772 subroutine shr_ncread_tCoordII(fn, tName, dates, secs, rc)
00773
00774 implicit none
00775
00776
00777
00778 character(*) ,intent(in) :: fn
00779 character(*) ,intent(in) :: tName
00780 integer(SHR_KIND_IN),intent(out) :: dates(:)
00781 integer(SHR_KIND_IN),intent(out) :: secs(:)
00782 integer(SHR_KIND_IN),intent(out),optional:: rc
00783
00784
00785
00786
00787 real(SHR_KIND_R8),allocatable :: tvar(:)
00788 character(SHR_KIND_CL) :: cfUnits
00789 character(SHR_KIND_CL) :: cfCalendar
00790 integer(SHR_KIND_IN) :: n
00791 integer(SHR_KIND_IN) :: nd,ns,nmax
00792 character(SHR_KIND_CS) :: units
00793 integer(SHR_KIND_IN) :: bdate
00794 real(SHR_KIND_R8) :: bsec
00795 integer(SHR_KIND_IN) :: ndate
00796 real(SHR_KIND_R8) :: nsec
00797 integer(SHR_KIND_IN) :: rCode
00798
00799
00800 character(*),parameter :: subName = "(shr_ncread_tCoordII)"
00801 character(*),parameter :: F00 = "('(shr_ncread_tCoordII) ',4a)"
00802 character(*),parameter :: F01 = "('(shr_ncread_tCoordII) ',a,i8)"
00803 character(*),parameter :: F02 = "('(shr_ncread_tCoordII) ',a,g17.7)"
00804 character(*),parameter :: eos = "[end-of-string]"
00805
00806
00807
00808
00809
00810
00811
00812
00813 rCode = 0
00814
00815 nd = size(dates)
00816 ns = size(secs)
00817 nmax = min(nd,ns)
00818 allocate(tVar(nmax))
00819
00820 call shr_ncread_tCoordRC(fn,tName,tVar,cfUnits,cfCalendar,rCode)
00821 call shr_string_parseCFtunit(cfUnits,units,bdate,bsec)
00822
00823 if (debug > 0 .and. s_loglev > 0) then
00824 write(s_logunit,F00) ' units = ',trim(units) ,eos
00825 write(s_logunit,F01) ' bdate = ',bdate
00826 write(s_logunit,F02) ' bsec = ',bsec
00827 write(s_logunit,F00) ' cfUnits = ',trim(cfUnits) ,eos
00828 write(s_logunit,F00) ' cfCalendar = ',trim(cfCalendar) ,eos
00829 endif
00830
00831 do n = 1,nMax
00832 call shr_cal_advDate(tVar(n),units,bDate,bSec,nDate,nSec,cfCalendar)
00833 dates(n) = nDate
00834 secs (n) = nSec
00835 enddo
00836
00837 deallocate(tVar)
00838
00839 if (present(rc)) rc = rCode
00840
00841 end subroutine shr_ncread_tCoordII
00842
00843
00844
00845
00846
00847
00848
00849
00850
00851
00852
00853
00854
00855
00856
00857
00858
00859
00860
00861
00862
00863
00864
00865
00866
00867
00868
00869
00870
00871
00872 subroutine shr_ncread_domain(fn, lonName, lon, latName, lat, &
00873 & maskName, mask, areaName, area, &
00874 & fracName, frac, rc)
00875
00876 implicit none
00877
00878
00879
00880 character(*) ,intent(in) :: fn
00881 character(*) ,intent(in) :: lonName
00882 real(SHR_KIND_R8) ,intent(out) :: lon(:,:)
00883 character(*) ,intent(in) :: latName
00884 real(SHR_KIND_R8) ,intent(out) :: lat(:,:)
00885 character(*) ,intent(in) ,optional:: maskName
00886 integer(SHR_KIND_IN),intent(out),optional:: mask(:,:)
00887 character(*) ,intent(in) ,optional:: areaName
00888 real(SHR_KIND_R8) ,intent(out),optional:: area(:,:)
00889 character(*) ,intent(in) ,optional:: fracName
00890 real(SHR_KIND_R8) ,intent(out),optional:: frac(:,:)
00891 integer(SHR_KIND_IN),intent(out),optional:: rc
00892
00893
00894
00895 real(SHR_KIND_R8),allocatable :: A4d(:,:,:,:)
00896 real(SHR_KIND_R8),allocatable :: P2d(:,:)
00897 character(SHR_KIND_CS) :: varName
00898 integer(SHR_KIND_IN) :: nflds
00899 integer(SHR_KIND_IN) :: n,i,j
00900 integer(SHR_KIND_IN) :: ndim,nd1,nd2
00901 integer(SHR_KIND_IN) :: pd1,pd2
00902 integer(SHR_KIND_IN) :: rCode
00903
00904
00905 character(*),parameter :: subName = "(shr_ncread_domain)"
00906 character(*),parameter :: F00 = "('(shr_ncread_domain) ',4a)"
00907 character(*),parameter :: F01 = "('(shr_ncread_domain) ',2a,3i6,2x,a)"
00908 character(*),parameter :: F02 = "('(shr_ncread_domain) ',a,i6)"
00909 character(*),parameter :: F03 = "('(shr_ncread_domain) ',a,2i6)"
00910 character(*),parameter :: F04 = "('(shr_ncread_domain) ',a,2g17.8)"
00911
00912 logical :: readmask
00913 logical :: readarea
00914 logical :: readfrac
00915
00916
00917
00918
00919 rCode = 0
00920
00921 nflds = 2
00922
00923 if (present(maskName).and.present(mask)) then
00924 nflds = nflds + 1
00925 else if (( present(maskName) .and. .not.present(mask)) .or. &
00926 (.not.present(maskName) .and. present(mask))) then
00927 write(s_logunit,F00) ' ERROR: maskName and mask must both be present or not '
00928 call shr_ncread_abort(subName//' ERROR subroutine arguments, mask')
00929 end if
00930
00931 if (present(areaName).and.present(area)) then
00932 nflds = nflds + 1
00933 else if (( present(areaName) .and. .not.present(area)) .or. &
00934 (.not.present(areaName) .and. present(area))) then
00935 write(s_logunit,F00) ' ERROR: areaName and area must both be present or not '
00936 call shr_ncread_abort(subName//' ERROR subroutine arguments, area')
00937 end if
00938
00939 if (present(fracName).and.present(frac)) then
00940 nflds = nflds + 1
00941 else if (( present(fracName) .and. .not.present(frac)) .or. &
00942 (.not.present(fracName) .and. present(frac))) then
00943 write(s_logunit,F00) ' ERROR: fracName and frac must both be present or not '
00944 call shr_ncread_abort(subName//' ERROR subroutine arguments, frac')
00945 end if
00946
00947
00948 readmask = .true.
00949 readarea = .true.
00950 readfrac = .true.
00951 do n=1,nflds
00952 if (n == 1) then
00953 varName = trim(lonName)
00954 allocate(P2d(size(lon,1),size(lon,2)))
00955 elseif (n == 2) then
00956 varName = trim(latName)
00957 allocate(P2d(size(lat,1),size(lat,2)))
00958 elseif (n > 2) then
00959 if (present(maskName) .and. readmask) then
00960 varName = trim(maskName)
00961
00962 allocate(P2d(size(mask,1),size(mask,2)))
00963 readmask = .false.
00964 else if (present(areaName) .and. readarea) then
00965 varName = trim(areaName)
00966 allocate(P2d(size(area,1),size(area,2)))
00967 readarea = .false.
00968 else if (present(fracName) .and. readfrac) then
00969 varName = trim(fracName)
00970 allocate(P2d(size(frac,1),size(frac,2)))
00971 readfrac = .false.
00972 endif
00973 end if
00974
00975 if (.not.shr_ncread_varExists(fn,varName)) &
00976 call shr_ncread_abort(subName//' ERROR var does not exist '//trim(varName))
00977
00978
00979 pd1 = size(P2d,1)
00980 pd2 = size(P2d,2)
00981
00982
00983 call shr_ncread_varDimNum(fn,varName,ndim)
00984 if (n > 2 .and. ndim /= 2) then
00985 write(s_logunit,F02) 'ERROR '//trim(varName)//' ndim = ',ndim
00986 call shr_ncread_abort(subName//' ERROR ndim must be 2 for '//trim(varName))
00987 elseif (ndim < 1 .or. ndim > 2) then
00988 write(s_logunit,F02) 'ERROR '//trim(varName)//' ndim = ',ndim
00989 call shr_ncread_abort(subName//' ERROR ndim must be 1 or 2 for '//trim(varName))
00990 endif
00991 nd1 = 1
00992 nd2 = 1
00993 if (ndim > 0) call shr_ncread_varDimSize(fn,varName,1,nd1)
00994 if (ndim > 1) call shr_ncread_varDimSize(fn,varName,2,nd2)
00995
00996
00997 if (n == 2 .and. ndim == 1) then
00998 if ( nd1 /= pd2) then
00999 write(s_logunit,F03) ' nd1 pd2 error ',nd1,pd2
01000 call shr_ncread_abort(subName//' ERROR nd1 pd2 error')
01001 endif
01002 elseif (ndim > 0 .and. nd1 /= pd1) then
01003 write(s_logunit,F03) ' nd1 pd1 error ',nd1,pd1
01004 call shr_ncread_abort(subName//' ERROR nd1 pd1 error')
01005 endif
01006 if (ndim > 1 .and. nd2 /= pd2) then
01007 write(s_logunit,F03) ' nd2 pd2 error ',nd2,pd2
01008 call shr_ncread_abort(subName//' ERROR nd2 pd2 error')
01009 endif
01010
01011
01012 allocate(A4d(nd1,nd2,1,1))
01013 A4d = 0.0_SHR_KIND_R8
01014 call shr_ncread_field4dG(fn,varName,rfld=A4d)
01015
01016
01017 do j = 1,pd2
01018 do i = 1,pd1
01019 if (n == 2 .and. ndim == 1) then
01020 P2d(i,j) = A4d(j,1,1,1)
01021 elseif (ndim == 1) then
01022 P2d(i,j) = A4d(i,1,1,1)
01023 else
01024 P2d(i,j) = A4d(i,j,1,1)
01025 endif
01026 enddo
01027 enddo
01028
01029
01030 if (n == 1) then
01031 lon(:,:) = P2d(:,:)
01032 elseif (n == 2) then
01033 lat(:,:) = P2d(:,:)
01034 elseif (n == 3) then
01035 mask(:,:) = nint(P2d(:,:))
01036 elseif (n == 4) then
01037 area(:,:) = P2d(:,:)
01038 elseif (n == 5) then
01039 frac(:,:) = P2d(:,:)
01040 endif
01041
01042
01043 deallocate(A4d)
01044 deallocate(P2d,stat=rCode)
01045
01046
01047 enddo
01048
01049 if (debug > 1 .and. s_loglev > 0) then
01050 write(s_logunit,F04) 'min/max lon ',minval(lon),maxval(lon)
01051 write(s_logunit,F04) 'min/max lat ',minval(lat),maxval(lat)
01052 write(s_logunit,F04) 'min/max mask ',minval(mask),maxval(mask)
01053 write(s_logunit,F04) 'min/max area ',minval(area),maxval(area)
01054 if (nflds >= 5 .and. s_loglev > 0) write(s_logunit,F04) 'min/max frac ',minval(frac),maxval(frac)
01055 endif
01056
01057 if (present(rc)) rc = rCode
01058
01059 end subroutine shr_ncread_domain
01060
01061
01062
01063
01064
01065
01066
01067
01068
01069
01070
01071
01072
01073
01074
01075
01076
01077
01078
01079
01080
01081
01082
01083
01084
01085
01086 subroutine shr_ncread_tField2dR8(fn, tIndex, fldName, fld, dim1, dim2, tName, fidi, rc)
01087
01088 implicit none
01089
01090
01091
01092 character(*) ,intent(in) :: fn
01093 integer(SHR_KIND_IN),intent(in) :: tIndex
01094 character(*) ,intent(in) :: fldName
01095 real(SHR_KIND_R8) ,intent(out) :: fld(:,:)
01096 character(*) ,intent(in) ,optional :: dim1
01097 character(*) ,intent(in) ,optional :: dim2
01098 character(*) ,intent(in) ,optional :: tName
01099 integer(SHR_KIND_IN),intent(in) ,optional :: fidi
01100 integer(SHR_KIND_IN),intent(out),optional :: rc
01101
01102
01103
01104
01105 real(SHR_KIND_R8),allocatable :: lfld(:,:,:,:)
01106 integer(SHR_KIND_IN) :: rCode
01107
01108
01109 character(*),parameter :: subName = "(shr_ncread_tField2dR8)"
01110 character(*),parameter :: F00 = "('(shr_ncread_tField2dR8) ',4a)"
01111
01112
01113
01114
01115 allocate(lfld(size(fld,1),size(fld,2),1,1))
01116
01117 if (present(dim1).and.present(dim2).and.present(tName)) then
01118 if (.not.present(fidi)) then
01119 call shr_ncread_field4dG(fn,fldName,rfld=lfld,dim1=dim1,dim2=dim2,dim3=tName,dim3i=tIndex,rc=rCode)
01120 else
01121 call shr_ncread_field4dG(fn,fldName,rfld=lfld,dim1=dim1,dim2=dim2,dim3=tName,dim3i=tIndex,fidi=fidi,rc=rCode)
01122 endif
01123 elseif (present(dim1).and.present(dim2)) then
01124 if (.not.present(fidi)) then
01125 call shr_ncread_field4dG(fn,fldName,rfld=lfld,dim1=dim1,dim2=dim2,dim3i=tIndex,rc=rCode)
01126 else
01127 call shr_ncread_field4dG(fn,fldName,rfld=lfld,dim1=dim1,dim2=dim2,dim3i=tIndex,fidi=fidi,rc=rCode)
01128 endif
01129 elseif (present(tName)) then
01130 if (.not.present(fidi)) then
01131 call shr_ncread_field4dG(fn,fldName,rfld=lfld,dim3=tName,dim3i=tIndex,rc=rCode)
01132 else
01133 call shr_ncread_field4dG(fn,fldName,rfld=lfld,dim3=tName,dim3i=tIndex,fidi=fidi,rc=rCode)
01134 endif
01135 elseif (.not.present(dim1).and..not.present(dim2).and..not.present(tName)) then
01136 if (.not.present(fidi)) then
01137 call shr_ncread_field4dG(fn,fldName,rfld=lfld,dim3i=tIndex,rc=rCode)
01138 else
01139 call shr_ncread_field4dG(fn,fldName,rfld=lfld,dim3i=tIndex,fidi=fidi,rc=rCode)
01140 endif
01141 else
01142 call shr_ncread_abort(subName//' ERROR argument combination not supported')
01143 endif
01144
01145 fld(:,:) = lfld(:,:,1,1)
01146 deallocate(lfld)
01147
01148 if (present(rc)) rc = rCode
01149
01150 end subroutine shr_ncread_tField2dR8
01151
01152
01153
01154
01155
01156
01157
01158
01159
01160
01161
01162
01163
01164
01165
01166
01167
01168
01169
01170
01171
01172
01173
01174
01175 subroutine shr_ncread_tField1dR8(fn, tIndex, fldName, fld, dim1, tName, fidi, rc)
01176
01177 implicit none
01178
01179
01180
01181 character(*) ,intent(in) :: fn
01182 integer(SHR_KIND_IN),intent(in) :: tIndex
01183 character(*) ,intent(in) :: fldName
01184 real(SHR_KIND_R8) ,intent(out) :: fld(:)
01185 character(*) ,intent(in) ,optional :: dim1
01186 character(*) ,intent(in) ,optional :: tName
01187 integer(SHR_KIND_IN),intent(in) ,optional :: fidi
01188 integer(SHR_KIND_IN),intent(out),optional :: rc
01189
01190
01191
01192
01193 real(SHR_KIND_R8),allocatable :: lfld(:,:,:,:)
01194 integer(SHR_KIND_IN) :: rCode
01195
01196
01197 character(*),parameter :: subName = "(shr_ncread_tField1dR8)"
01198 character(*),parameter :: F00 = "('(shr_ncread_tField1dR8) ',4a)"
01199
01200
01201
01202
01203
01204 allocate(lfld(size(fld,1),1,1,1))
01205
01206 if (present(dim1).and.present(tName)) then
01207 if (.not.present(fidi)) then
01208 call shr_ncread_field4dG(fn,fldName,rfld=lfld,dim1=dim1,dim2=tName,dim2i=tIndex,rc=rCode)
01209 else
01210 call shr_ncread_field4dG(fn,fldName,rfld=lfld,dim1=dim1,dim2=tName,dim2i=tIndex,fidi=fidi,rc=rCode)
01211 endif
01212 elseif (present(dim1)) then
01213 if (.not.present(fidi)) then
01214 call shr_ncread_field4dG(fn,fldName,rfld=lfld,dim1=dim1,dim2i=tIndex,rc=rCode)
01215 else
01216 call shr_ncread_field4dG(fn,fldName,rfld=lfld,dim1=dim1,dim2i=tIndex,fidi=fidi,rc=rCode)
01217 endif
01218 elseif (present(tName)) then
01219 if (.not.present(fidi)) then
01220 call shr_ncread_field4dG(fn,fldName,rfld=lfld,dim2=tName,dim2i=tIndex,rc=rCode)
01221 else
01222 call shr_ncread_field4dG(fn,fldName,rfld=lfld,dim2=tName,dim2i=tIndex,fidi=fidi,rc=rCode)
01223 endif
01224 elseif (.not.present(dim1).and..not.present(tName)) then
01225 if (.not.present(fidi)) then
01226 call shr_ncread_field4dG(fn,fldName,rfld=lfld,dim2i=tIndex,rc=rCode)
01227 else
01228 call shr_ncread_field4dG(fn,fldName,rfld=lfld,dim2i=tIndex,fidi=fidi,rc=rCode)
01229 endif
01230 else
01231 call shr_ncread_abort(subName//' ERROR argument combination not supported')
01232 endif
01233
01234 fld(:) = lfld(:,1,1,1)
01235 deallocate(lfld)
01236
01237 if (present(rc)) rc = rCode
01238
01239 end subroutine shr_ncread_tField1dR8
01240
01241
01242
01243
01244
01245
01246
01247
01248
01249
01250
01251
01252
01253
01254
01255
01256
01257
01258
01259
01260
01261
01262
01263
01264
01265
01266 subroutine shr_ncread_tField2dIN(fn, tIndex, fldName, fld, dim1, dim2, tName, fidi, rc)
01267
01268 implicit none
01269
01270
01271
01272 character(*) ,intent(in) :: fn
01273 integer(SHR_KIND_IN),intent(in) :: tIndex
01274 character(*) ,intent(in) :: fldName
01275 integer(SHR_KIND_IN),intent(out) :: fld(:,:)
01276 character(*) ,intent(in) ,optional :: dim1
01277 character(*) ,intent(in) ,optional :: dim2
01278 character(*) ,intent(in) ,optional :: tName
01279 integer(SHR_KIND_IN),intent(in) ,optional :: fidi
01280 integer(SHR_KIND_IN),intent(out),optional :: rc
01281
01282
01283
01284
01285 integer(SHR_KIND_IN),allocatable :: lfld(:,:,:,:)
01286 integer(SHR_KIND_IN) :: rCode
01287
01288
01289 character(*),parameter :: subName = "(shr_ncread_tField2dIN)"
01290 character(*),parameter :: F00 = "('(shr_ncread_tField2dIN) ',4a)"
01291
01292
01293
01294
01295
01296 allocate(lfld(size(fld,1),size(fld,2),1,1))
01297
01298 if (present(dim1).and.present(dim2).and.present(tName)) then
01299 if (.not.present(fidi)) then
01300 call shr_ncread_field4dG(fn,fldName,ifld=lfld,dim1=dim1,dim2=dim2,dim3=tName,dim3i=tIndex,rc=rCode)
01301 else
01302 call shr_ncread_field4dG(fn,fldName,ifld=lfld,dim1=dim1,dim2=dim2,dim3=tName,dim3i=tIndex,fidi=fidi,rc=rCode)
01303 endif
01304 elseif (present(dim1).and.present(dim2)) then
01305 if (.not.present(fidi)) then
01306 call shr_ncread_field4dG(fn,fldName,ifld=lfld,dim1=dim1,dim2=dim2,dim3i=tIndex,rc=rCode)
01307 else
01308 call shr_ncread_field4dG(fn,fldName,ifld=lfld,dim1=dim1,dim2=dim2,dim3i=tIndex,fidi=fidi,rc=rCode)
01309 endif
01310 elseif (present(tName)) then
01311 if (.not.present(fidi)) then
01312 call shr_ncread_field4dG(fn,fldName,ifld=lfld,dim3=tName,dim3i=tIndex,rc=rCode)
01313 else
01314 call shr_ncread_field4dG(fn,fldName,ifld=lfld,dim3=tName,dim3i=tIndex,fidi=fidi,rc=rCode)
01315 endif
01316 elseif (.not.present(dim1).and..not.present(dim2).and..not.present(tName)) then
01317 if (.not.present(fidi)) then
01318 call shr_ncread_field4dG(fn,fldName,ifld=lfld,dim3i=tIndex,rc=rCode)
01319 else
01320 call shr_ncread_field4dG(fn,fldName,ifld=lfld,dim3i=tIndex,fidi=fidi,rc=rCode)
01321 endif
01322 else
01323 call shr_ncread_abort(subName//' ERROR argument combination not supported')
01324 endif
01325
01326 fld(:,:) = lfld(:,:,1,1)
01327 deallocate(lfld)
01328
01329 if (present(rc)) rc = rCode
01330
01331 end subroutine shr_ncread_tField2dIN
01332
01333
01334
01335
01336
01337
01338
01339
01340
01341
01342
01343
01344
01345
01346
01347
01348
01349
01350
01351
01352
01353
01354
01355
01356 subroutine shr_ncread_tField1dIN(fn, tIndex, fldName, fld, dim1, tName, fidi, rc)
01357
01358 implicit none
01359
01360
01361
01362 character(*) ,intent(in) :: fn
01363 integer(SHR_KIND_IN),intent(in) :: tIndex
01364 character(*) ,intent(in) :: fldName
01365 integer(SHR_KIND_IN),intent(out) :: fld(:)
01366 character(*) ,intent(in) ,optional :: dim1
01367 character(*) ,intent(in) ,optional :: tName
01368 integer(SHR_KIND_IN),intent(in) ,optional :: fidi
01369 integer(SHR_KIND_IN),intent(out),optional :: rc
01370
01371
01372
01373
01374 integer(SHR_KIND_IN),allocatable :: lfld(:,:,:,:)
01375 integer(SHR_KIND_IN) :: rCode
01376
01377
01378 character(*),parameter :: subName = "(shr_ncread_tField1dIN)"
01379 character(*),parameter :: F00 = "('(shr_ncread_tField1dIN) ',4a)"
01380
01381
01382
01383
01384
01385 allocate(lfld(size(fld,1),1,1,1))
01386
01387 if (present(dim1).and.present(tName)) then
01388 if (.not.present(fidi)) then
01389 call shr_ncread_field4dG(fn,fldName,ifld=lfld,dim1=dim1,dim2=tName,dim2i=tIndex,rc=rCode)
01390 else
01391 call shr_ncread_field4dG(fn,fldName,ifld=lfld,dim1=dim1,dim2=tName,dim2i=tIndex,fidi=fidi,rc=rCode)
01392 endif
01393 elseif (present(dim1)) then
01394 if (.not.present(fidi)) then
01395 call shr_ncread_field4dG(fn,fldName,ifld=lfld,dim1=dim1,dim2i=tIndex,rc=rCode)
01396 else
01397 call shr_ncread_field4dG(fn,fldName,ifld=lfld,dim1=dim1,dim2i=tIndex,fidi=fidi,rc=rCode)
01398 endif
01399 elseif (present(tName)) then
01400 if (.not.present(fidi)) then
01401 call shr_ncread_field4dG(fn,fldName,ifld=lfld,dim2=tName,dim2i=tIndex,rc=rCode)
01402 else
01403 call shr_ncread_field4dG(fn,fldName,ifld=lfld,dim2=tName,dim2i=tIndex,fidi=fidi,rc=rCode)
01404 endif
01405 elseif (.not.present(dim1).and..not.present(tName)) then
01406 if (.not.present(fidi)) then
01407 call shr_ncread_field4dG(fn,fldName,ifld=lfld,dim2i=tIndex,rc=rCode)
01408 else
01409 call shr_ncread_field4dG(fn,fldName,ifld=lfld,dim2i=tIndex,fidi=fidi,rc=rCode)
01410 endif
01411 else
01412 call shr_ncread_abort(subName//' ERROR argument combination not supported')
01413 endif
01414
01415 fld(:) = lfld(:,1,1,1)
01416 deallocate(lfld)
01417
01418 if (present(rc)) rc = rCode
01419
01420 end subroutine shr_ncread_tField1dIN
01421
01422
01423
01424
01425
01426
01427
01428
01429
01430
01431
01432
01433
01434
01435
01436
01437
01438
01439
01440
01441
01442
01443
01444
01445
01446
01447
01448
01449
01450
01451
01452 subroutine shr_ncread_field4dG(fn, fldName, rfld, ifld, &
01453 dim1, dim1i, dim2, dim2i, dim3, dim3i, dim4, dim4i, &
01454 dim5, dim5i, dim6, dim6i, fidi, rc)
01455
01456 implicit none
01457
01458
01459
01460 character(*) ,intent(in) :: fn
01461 character(*) ,intent(in) :: fldName
01462 real(SHR_KIND_R8) ,intent(out),optional :: rfld(:,:,:,:)
01463 integer(SHR_KIND_IN),intent(out),optional :: ifld(:,:,:,:)
01464 character(*) ,intent(in) ,optional :: dim1
01465 integer(SHR_KIND_IN),intent(in) ,optional :: dim1i
01466 character(*) ,intent(in) ,optional :: dim2
01467 integer(SHR_KIND_IN),intent(in) ,optional :: dim2i
01468 character(*) ,intent(in) ,optional :: dim3
01469 integer(SHR_KIND_IN),intent(in) ,optional :: dim3i
01470 character(*) ,intent(in) ,optional :: dim4
01471 integer(SHR_KIND_IN),intent(in) ,optional :: dim4i
01472 character(*) ,intent(in) ,optional :: dim5
01473 integer(SHR_KIND_IN),intent(in) ,optional :: dim5i
01474 character(*) ,intent(in) ,optional :: dim6
01475 integer(SHR_KIND_IN),intent(in) ,optional :: dim6i
01476 integer(SHR_KIND_IN),intent(in) ,optional :: fidi
01477 integer(SHR_KIND_IN),intent(out),optional :: rc
01478
01479
01480
01481 integer(SHR_KIND_IN),parameter :: maxd = 4
01482 integer(SHR_KIND_IN) :: fid
01483 integer(SHR_KIND_IN) :: vid
01484 integer(SHR_KIND_IN) :: xtype
01485 integer(SHR_KIND_IN) :: ndims
01486 integer(SHR_KIND_IN) :: n,n1,n2,n3,n4,k
01487 integer(SHR_KIND_IN) ,allocatable :: dimid(:)
01488 integer(SHR_KIND_IN) ,allocatable :: dids(:)
01489 integer(SHR_KIND_IN) ,allocatable :: start(:)
01490 integer(SHR_KIND_IN) ,allocatable :: count(:)
01491 integer(SHR_KIND_IN) ,allocatable :: len(:)
01492 character(SHR_KIND_CS),allocatable :: name(:)
01493 real(SHR_KIND_R8) ,allocatable :: rin(:,:)
01494 integer(SHR_KIND_IN) ,allocatable :: iin(:,:)
01495 integer(SHR_KIND_IN) ,allocatable :: start2d(:)
01496 integer(SHR_KIND_IN) ,allocatable :: count2d(:)
01497 logical :: found
01498 integer(SHR_KIND_IN) :: rCode
01499
01500
01501 character(*),parameter :: subName = "(shr_ncread_field4dG) "
01502 character(*),parameter :: F00 = "('(shr_ncread_field4dG) ',4a)"
01503 character(*),parameter :: F01 = "('(shr_ncread_field4dG) ',2a,3i6,2x,a)"
01504
01505
01506
01507
01508
01509
01510 if (present(rfld).and.present(ifld)) then
01511 call shr_ncread_abort(subName//'both rfld and ifld should not be sent')
01512 endif
01513 if (.not.present(rfld).and..not.present(ifld)) then
01514 call shr_ncread_abort(subName//'either rfld or ifld must be sent')
01515 endif
01516
01517 if (.not.present(fidi)) then
01518 call shr_ncread_open(fn,fid,rCode)
01519 else
01520 fid = fidi
01521 endif
01522
01523
01524 rCode = nf90_inq_varid(fid,trim(fldName),vid)
01525 call shr_ncread_handleErr(rCode,subName//'inq varid vid: '//trim(fldName))
01526 rCode = nf90_inquire_variable(fid,vid,xtype=xtype,ndims=ndims)
01527 call shr_ncread_handleErr(rCode,subName//'inquire variable ndims: '//trim(fldName))
01528
01529
01530 n4 = max(ndims,maxd)
01531 allocate(dimid (n4)) ; dimid = 0
01532 allocate(dids (n4)) ; dids = 0
01533 allocate(name (n4)) ; name = ' '
01534 allocate(len (n4)) ; len = 1
01535 allocate(start (n4)) ; start = 1
01536 allocate(count (n4)) ; count = 1
01537 allocate(start2d(n4)) ; start2d = 1
01538 allocate(count2d(n4)) ; count2d = 1
01539
01540
01541 rCode = nf90_inquire_variable(fid,vid,dimids=dids)
01542 call shr_ncread_handleErr(rCode,subName//'inquire variable dids: '//trim(fldName))
01543 do n=1,ndims
01544 rCode = nf90_inquire_dimension(fid,dids(n),name=name(n),len=len(n))
01545 call shr_ncread_handleErr(rCode,subName//'inquire dimension len: '//trim(fldName))
01546 enddo
01547
01548
01549 if (present(dim1)) then
01550 do n=1,ndims
01551 if (trim(dim1) == trim(name(n))) dimid(1) = n
01552 enddo
01553 endif
01554 if (present(dim2)) then
01555 do n=1,ndims
01556 if (trim(dim2) == trim(name(n))) dimid(2) = n
01557 enddo
01558 endif
01559 if (present(dim3)) then
01560 do n=1,ndims
01561 if (trim(dim3) == trim(name(n))) dimid(3) = n
01562 enddo
01563 endif
01564 if (present(dim4)) then
01565 do n=1,ndims
01566 if (trim(dim4) == trim(name(n))) dimid(4) = n
01567 enddo
01568 endif
01569 if (present(dim5)) then
01570 do n=1,ndims
01571 if (trim(dim5) == trim(name(n))) dimid(5) = n
01572 enddo
01573 endif
01574 if (present(dim6)) then
01575 do n=1,ndims
01576 if (trim(dim6) == trim(name(n))) dimid(6) = n
01577 enddo
01578 endif
01579
01580
01581 do n1=1,max(maxd,ndims)
01582 k = 1
01583 do while (dimid(n1) == 0)
01584 found = .false.
01585 do n2 = 1,maxd
01586 if (dimid(n2) == k) found = .true.
01587 enddo
01588 if (found) then
01589 k = k + 1
01590 else
01591 dimid(n1) = k
01592 endif
01593 enddo
01594 enddo
01595
01596
01597 do n1=1,maxd
01598 if (dimid(n1) <= ndims) then
01599 count(dimid(n1)) = len(dimid(n1))
01600 else
01601 count(dimid(n1)) = 1
01602 endif
01603 enddo
01604
01605
01606 if (present(dim1i)) then
01607 if (dim1i < 1 .or. dim1i > len(dimid(1))) &
01608 call shr_ncread_abort(subName//'dim1i setting: '//trim(fldName))
01609 start(dimid(1)) = dim1i
01610 count(dimid(1)) = 1
01611 endif
01612 if (present(dim2i)) then
01613 if (dim2i < 1 .or. dim2i > len(dimid(2))) &
01614 call shr_ncread_abort(subName//'dim2i setting: '//trim(fldName))
01615 start(dimid(2)) = dim2i
01616 count(dimid(2)) = 1
01617 endif
01618 if (present(dim3i)) then
01619 if (dim3i < 1 .or. dim3i > len(dimid(3))) &
01620 call shr_ncread_abort(subName//'dim3i setting: '//trim(fldName))
01621 start(dimid(3)) = dim3i
01622 count(dimid(3)) = 1
01623 endif
01624 if (present(dim4i)) then
01625 if (dim4i < 1 .or. dim4i > len(dimid(4))) &
01626 call shr_ncread_abort(subName//'dim4i setting: '//trim(fldName))
01627 start(dimid(4)) = dim4i
01628 count(dimid(4)) = 1
01629 endif
01630 if (present(dim5i)) then
01631 if (dim5i < 1 .or. dim5i > len(dimid(5))) &
01632 call shr_ncread_abort(subName//'dim5i setting: '//trim(fldName))
01633 start(dimid(5)) = dim5i
01634 count(dimid(5)) = 1
01635 endif
01636 if (present(dim6i)) then
01637 if (dim6i < 1 .or. dim6i > len(dimid(6))) &
01638 call shr_ncread_abort(subName//'dim6i setting: '//trim(fldName))
01639 start(dimid(6)) = dim6i
01640 count(dimid(6)) = 1
01641 endif
01642
01643
01644 do n=1,maxd
01645 if (present(rfld)) then
01646 if (size(rfld,n) /= count(dimid(n))) then
01647 call shr_ncread_abort(subName//'fld size does not agree with count: '//trim(fldName))
01648 endif
01649 endif
01650 if (present(ifld)) then
01651 if (size(ifld,n) /= count(dimid(n))) then
01652 call shr_ncread_abort(subName//'fld size does not agree with count: '//trim(fldName))
01653 endif
01654 endif
01655 enddo
01656
01657
01658
01659 if (dimid(1) > dimid(2)) then
01660 allocate(rin(count(dimid(2)),count(dimid(1))))
01661 allocate(iin(count(dimid(2)),count(dimid(1))))
01662 else
01663 allocate(rin(count(dimid(1)),count(dimid(2))))
01664 allocate(iin(count(dimid(1)),count(dimid(2))))
01665 endif
01666 start2d = start
01667 count2d = count
01668 count2d(dimid(3)) = 1
01669 count2d(dimid(4)) = 1
01670 do n4 = 1,count(dimid(4))
01671 do n3 = 1,count(dimid(3))
01672 start2d(dimid(3)) = n3 + start(dimid(3)) - 1
01673 start2d(dimid(4)) = n4 + start(dimid(4)) - 1
01674 if (present(rfld)) then
01675 rCode = nf90_get_var(fid,vid,rin,start=start2d,count=count2d)
01676 elseif (present(ifld)) then
01677 rCode = nf90_get_var(fid,vid,iin,start=start2d,count=count2d)
01678 endif
01679 call shr_ncread_handleErr(rCode,subName//'get var: '//trim(fldName))
01680
01681
01682
01683
01684
01685
01686
01687
01688
01689
01690
01691
01692
01693
01694
01695
01696 do n2 = 1,count(dimid(2))
01697 do n1 = 1,count(dimid(1))
01698 if (dimid(1) > dimid(2)) then
01699 if (present(rfld)) then
01700 rfld(n1,n2,n3,n4) = rin(n2,n1)
01701 elseif (present(ifld)) then
01702 ifld(n1,n2,n3,n4) = iin(n2,n1)
01703 endif
01704 else
01705 if (present(rfld)) then
01706 rfld(n1,n2,n3,n4) = rin(n1,n2)
01707 elseif (present(ifld)) then
01708 ifld(n1,n2,n3,n4) = iin(n1,n2)
01709 endif
01710 endif
01711 enddo
01712 enddo
01713 enddo
01714 enddo
01715 deallocate(rin)
01716 deallocate(iin)
01717
01718 deallocate(dimid)
01719 deallocate(dids)
01720 deallocate(start)
01721 deallocate(count)
01722 deallocate(name)
01723 deallocate(len)
01724 deallocate(start2d)
01725 deallocate(count2d)
01726 if (.not.present(fidi)) then
01727 call shr_ncread_close(fid,rCode)
01728 endif
01729
01730 if (present(rc)) rc = rCode
01731
01732 end subroutine shr_ncread_field4dG
01733
01734
01735
01736
01737
01738
01739
01740
01741
01742
01743
01744
01745
01746
01747
01748
01749 subroutine shr_ncread_open(fileName,fid,rCode)
01750
01751 implicit none
01752
01753
01754
01755 character(*),intent(in) :: fileName
01756 integer(SHR_KIND_IN),intent(out) :: fid
01757 integer(SHR_KIND_IN),intent(out) :: rCode
01758
01759
01760
01761
01762 integer(SHR_KIND_IN) :: n
01763 logical :: exists
01764
01765
01766 character(*),parameter :: subName = '(shr_ncread_open) '
01767 character(*),parameter :: F00 = "('(shr_ncread_open) ',4a)"
01768
01769
01770
01771
01772
01773
01774 inquire(file=trim(fileName),exist=exists)
01775 if (.not.exists) then
01776 if (s_loglev > 0) write(s_logunit,F00) "ERROR: file does not exist: ",trim(fileName)
01777 call shr_ncread_handleErr(rCode,subName//"ERROR: file does not exist: "//trim(fileName))
01778 end if
01779
01780
01781 if (debug > 1 .and. s_loglev > 0) write(s_logunit,F00) 'open netCDF data file: ',trim(fileName)
01782 rCode = nf90_open(fileName,nf90_nowrite,fid)
01783 call shr_ncread_handleErr(rCode, subName//"ERROR opening input data file")
01784
01785 end subroutine shr_ncread_open
01786
01787
01788
01789
01790
01791
01792
01793
01794
01795
01796
01797
01798
01799
01800
01801
01802
01803 subroutine shr_ncread_attrbs(fid,vn)
01804
01805 implicit none
01806
01807
01808
01809 integer(SHR_KIND_IN),intent(in) :: fid
01810 integer(SHR_KIND_IN),intent(in) :: vn
01811
01812
01813
01814
01815 integer(SHR_KIND_IN),parameter :: strlen = SHR_KIND_CL
01816 character(len=strlen) :: cvalue
01817 integer(SHR_KIND_IN) :: natt
01818 integer(SHR_KIND_IN) :: xtype
01819 integer(SHR_KIND_IN) :: an
01820 integer(SHR_KIND_IN) :: len
01821 integer(SHR_KIND_IN) :: rCode
01822 character(SHR_KIND_CS) :: name
01823
01824
01825 character(*),parameter :: subName = "(shr_ncread_attrbs)"
01826 character(*),parameter :: F00 = "('(shr_ncread_attrbs) ',4a)"
01827 character(*),parameter :: F04 = "('(shr_ncread_attrbs) ',4x,a,i4,2a)"
01828
01829
01830
01831
01832
01833 if (vn == nf90_global) then
01834 rCode = nf90_inquire(fid,nAttributes=natt)
01835 else
01836 rCode = nf90_inquire_variable(fid,vn,nAtts=natt)
01837 endif
01838
01839 do an=1,natt
01840 rCode = nf90_inq_attname(fid,VN,an,name)
01841 call shr_ncread_handleErr(rCode,subName//' nf90_inq_attname')
01842 rCode = nf90_inquire_attribute(fid, VN, trim(name), xtype, len)
01843 call shr_ncread_handleErr(rCode,subName//' nf90_inq_att')
01844 if (xtype == NF90_CHAR) then
01845 if (len < strlen) then
01846 cvalue = ' '
01847 rCode = nf90_get_att(fid,VN,trim(name),cvalue)
01848 call shr_ncread_handleErr(rCode,subName//' nf90_get_att cvalue')
01849 if (s_loglev > 0) write(s_logunit,F04) 'attribute: ',an,' '//trim(name)//':',trim(cvalue)
01850 else
01851 if (s_loglev > 0) write(s_logunit,F04) 'attribute: ',an,' '//trim(name)//':','*** too long ***'
01852 endif
01853 else
01854 if (s_loglev > 0) write(s_logunit,F04) 'attribute: ',an,' '//trim(name)//':',' *** not char ***'
01855 endif
01856 enddo
01857
01858 end subroutine shr_ncread_attrbs
01859
01860
01861
01862
01863
01864
01865
01866
01867
01868
01869
01870
01871
01872
01873
01874
01875 subroutine shr_ncread_close(fid,rCode)
01876
01877 implicit none
01878
01879
01880
01881 integer(SHR_KIND_IN),intent(in) :: fid
01882 integer(SHR_KIND_IN),intent(out) :: rCode
01883
01884
01885
01886
01887 character(*),parameter :: subName = "(shr_ncread_close)"
01888 character(*),parameter :: F00 = "('(shr_ncread_close) ',4a)"
01889
01890
01891
01892
01893
01894
01895 if (debug > 1 .and. s_loglev > 0) write(s_logunit,F00) 'close netCDF input data file '
01896 rCode = nf90_close(fid)
01897 call shr_ncread_handleErr(rCode, subName//" ERROR closing input data file")
01898
01899 end subroutine shr_ncread_close
01900
01901
01902
01903
01904
01905
01906
01907
01908
01909
01910
01911
01912
01913
01914
01915
01916
01917 subroutine shr_ncread_print(fileName, rc)
01918
01919 implicit none
01920
01921
01922
01923 character(*) ,intent (in) :: fileName
01924 integer(SHR_KIND_IN),optional,intent (in) :: rc
01925
01926
01927
01928
01929 integer(SHR_KIND_IN) :: fid,vid,did
01930 integer(SHR_KIND_IN) :: nvar,ndim,natt
01931 integer(SHR_KIND_IN) :: vn,dn,an
01932 integer(SHR_KIND_IN) :: xtype
01933 integer(SHR_KIND_IN) :: len
01934 character(SHR_KIND_CS) :: name
01935 integer(SHR_KIND_IN) :: debug0
01936 integer(SHR_KIND_IN) :: rCode
01937
01938
01939 character(*),parameter :: subName = "(shr_ncread_print)"
01940 character(*),parameter :: F00 = "('(shr_ncread_print) ',4a)"
01941 character(*),parameter :: F01 = "('(shr_ncread_print) ',2x,a,i4,a,i4,a)"
01942 character(*),parameter :: F02 = "('(shr_ncread_print) ',2x,2a,i4,2i12,i8)"
01943 character(*),parameter :: F03 = "('(shr_ncread_print) ',2x,2a,i4,2g16.7,i8)"
01944 character(*),parameter :: F04 = "('(shr_ncread_print) ',2x,a,i4,2a)"
01945 character(*),parameter :: F05 = "('(shr_ncread_print) ',2x,a,3i6)"
01946 character(*),parameter :: F06 = "('(shr_ncread_print) ',2x,a,i4,a,i4,2a,i4)"
01947
01948
01949
01950
01951
01952 debug0 = debug
01953 call shr_ncread_setDebug(2)
01954
01955 if (s_loglev > 0) write(s_logunit,F00) fileName
01956
01957 call shr_ncread_open(fileName,fid,rCode)
01958
01959 rCode = nf90_inquire(fid,ndim,nvar,natt)
01960 call shr_ncread_handleErr(rCode,subName//' nf90_inquire')
01961
01962 if (s_loglev > 0) write(s_logunit,F05) 'ndim,nvar,natt: ',ndim,nvar,natt
01963
01964 call shr_ncread_attrbs(fid,nf90_global)
01965
01966 do dn=1,ndim
01967 rcode = nf90_inquire_dimension(fid,dn,name,len)
01968 call shr_ncread_handleErr(rCode,subName//' nf90_inquire_dim')
01969 if (s_loglev > 0) write(s_logunit,F01) 'dimension: ',dn,' '//trim(name)//'(',len,')'
01970 enddo
01971
01972 do vn=1,nvar
01973 rcode = nf90_inquire_variable(fid,vn,name,xtype,len,nAtts=natt)
01974 if (s_loglev > 0) write(s_logunit,F06) 'variable: ',vn,' '//trim(name),len,' dims',' xtype=',xtype
01975 call shr_ncread_attrbs(fid,vn)
01976 enddo
01977
01978 call shr_ncread_close(fid,rCode)
01979
01980 call shr_ncread_setDebug(debug0)
01981
01982 end subroutine shr_ncread_print
01983
01984
01985
01986
01987
01988
01989
01990
01991
01992
01993
01994
01995
01996
01997
01998
01999
02000 subroutine shr_ncread_handleErr(rCode, str)
02001
02002 implicit none
02003
02004
02005
02006 integer(SHR_KIND_IN),intent (in) :: rCode
02007 character(*) ,intent (in) :: str
02008
02009
02010
02011
02012 character(*),parameter :: F00 = "('(shr_ncread_handleErr) ',4a)"
02013
02014
02015
02016
02017
02018 if (rCode /= nf90_noerr) then
02019 write(s_logunit,F00) "netCDF error: ",trim(nf90_strerror(rCode))
02020 call shr_ncread_abort(str)
02021 end if
02022
02023 end subroutine shr_ncread_handleErr
02024
02025
02026
02027
02028
02029
02030
02031
02032
02033
02034
02035
02036
02037
02038
02039
02040 subroutine shr_ncread_setAbort(flag)
02041
02042 implicit none
02043
02044
02045
02046 logical,intent(in) :: flag
02047
02048
02049
02050
02051
02052
02053 character(*),parameter :: subName = "('shr_ncread_setAbort') "
02054 character(*),parameter :: F00 = "('(shr_ncread_setAbort) ',a) "
02055
02056
02057
02058
02059
02060 doabort = flag
02061
02062 end subroutine shr_ncread_setAbort
02063
02064
02065
02066
02067
02068
02069
02070
02071
02072
02073
02074
02075
02076
02077
02078
02079 subroutine shr_ncread_setDebug(iflag)
02080
02081 implicit none
02082
02083
02084
02085 integer,intent(in) :: iflag
02086
02087
02088
02089
02090
02091
02092 character(*),parameter :: subName = "('shr_ncread_setDebug') "
02093 character(*),parameter :: F00 = "('(shr_ncread_setDebug) ',a) "
02094
02095
02096
02097
02098
02099 debug = iflag
02100
02101 end subroutine shr_ncread_setDebug
02102
02103
02104
02105
02106
02107
02108
02109
02110
02111
02112
02113
02114
02115
02116
02117
02118 subroutine shr_ncread_abort(string)
02119
02120 implicit none
02121
02122
02123
02124 character(*),optional,intent(IN) :: string
02125
02126
02127
02128
02129 character(SHR_KIND_CL) :: lstring
02130 character(*),parameter :: subName = "(shr_ncread_abort)"
02131 character(*),parameter :: F00 = "('(shr_ncread_abort) ',a)"
02132
02133
02134
02135
02136
02137 lstring = ''
02138 if (present(string)) lstring = string
02139
02140 if (doabort) then
02141 call shr_sys_abort(lstring)
02142 else
02143 write(s_logunit,F00) ' no abort:'//trim(lstring)
02144 endif
02145
02146 end subroutine shr_ncread_abort
02147
02148
02149 end module shr_ncread_mod
02150