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 module shr_stream_mod
00032
00033 use shr_sys_mod
00034 use shr_kind_mod
00035 use shr_const_mod
00036 use shr_string_mod
00037 use shr_ncread_mod
00038 use shr_mpi_mod
00039 use shr_file_mod
00040 use shr_cal_mod
00041
00042 use shr_log_mod, only : s_loglev => shr_log_Level
00043 use shr_log_mod, only : s_logunit => shr_log_Unit
00044 use perf_mod
00045
00046 implicit none
00047
00048 private
00049
00050
00051
00052 public :: shr_stream_streamType
00053 public :: shr_stream_fileType
00054
00055
00056
00057 public :: shr_stream_init
00058 public :: shr_stream_set
00059 public :: shr_stream_parseInput
00060 public :: shr_stream_findBounds
00061 public :: shr_stream_getFileFieldList
00062 public :: shr_stream_getModelFieldList
00063 public :: shr_stream_getFileFieldName
00064 public :: shr_stream_getModelFieldName
00065 public :: shr_stream_getFirstFileName
00066 public :: shr_stream_getNextFileName
00067 public :: shr_stream_getPrevFileName
00068 public :: shr_stream_getFilePath
00069
00070 public :: shr_stream_getDataSource
00071 public :: shr_stream_getDomainInfo
00072 public :: shr_stream_getFile
00073 public :: shr_stream_getNFiles
00074 public :: shr_stream_dataDump
00075 public :: shr_stream_restWrite
00076 public :: shr_stream_restRead
00077 public :: shr_stream_setDebug
00078 public :: shr_stream_setAbort
00079 public :: shr_stream_getDebug
00080
00081
00082
00083
00084
00085
00086
00087
00088 character(SHR_KIND_CS),parameter,public :: shr_stream_taxis_cycle = 'cycle'
00089 character(SHR_KIND_CS),parameter,public :: shr_stream_taxis_extend = 'extend'
00090 character(SHR_KIND_CS),parameter,public :: shr_stream_taxis_limit = 'limit'
00091 character(SHR_KIND_CS),parameter,public :: shr_stream_file_null = 'not_set'
00092
00093
00094 type shr_stream_fileType
00095 character(SHR_KIND_CL) :: name = ''
00096 logical :: haveData = .false.
00097 integer (SHR_KIND_IN) :: nt = 0
00098 integer (SHR_KIND_IN),pointer :: date(:)
00099 integer (SHR_KIND_IN),pointer :: secs(:)
00100 end type shr_stream_fileType
00101
00102
00103 integer(SHR_KIND_IN),parameter :: nFileMax = 1000
00104
00105 type shr_stream_streamType
00106
00107
00108 logical :: init = .false.
00109 integer (SHR_KIND_IN) :: nFiles = 0
00110 character(SHR_KIND_CS) :: dataSource = 'undefined'
00111 character(SHR_KIND_CL) :: filePath = ''
00112 type(shr_stream_fileType) :: file(nFileMax)
00113
00114
00115 integer(SHR_KIND_IN) :: yearFirst = 0
00116 integer(SHR_KIND_IN) :: yearLast = 0
00117 integer(SHR_KIND_IN) :: yearAlign = 0
00118 integer(SHR_KIND_IN) :: offset = 0
00119 character(SHR_KIND_CS) :: taxMode = trim(shr_stream_taxis_cycle)
00120
00121
00122 integer(SHR_KIND_IN) :: k_lvd,n_lvd
00123 logical :: found_lvd = .false.
00124 integer(SHR_KIND_IN) :: k_gvd,n_gvd
00125 logical :: found_gvd = .false.
00126
00127
00128 character(SHR_KIND_CX) :: fldListFile = ''
00129 character(SHR_KIND_CX) :: fldListModel = ''
00130 character(SHR_KIND_CL) :: domFilePath = ''
00131 character(SHR_KIND_CL) :: domFileName = ''
00132 character(SHR_KIND_CS) :: domTvarName = ''
00133 character(SHR_KIND_CS) :: domXvarName = ''
00134 character(SHR_KIND_CS) :: domYvarName = ''
00135 character(SHR_KIND_CS) :: domAreaName = ''
00136 character(SHR_KIND_CS) :: domMaskName = ''
00137
00138 character(SHR_KIND_CS) :: tInterpAlgo = 'unused'
00139 end type shr_stream_streamType
00140
00141
00142 real(SHR_KIND_R8) ,parameter :: spd = SHR_CONST_CDAY
00143
00144 integer(SHR_KIND_IN),save :: debug = 0
00145 logical ,save :: doabort = .true.
00146
00147
00148 contains
00149
00150
00151
00152
00153
00154
00155
00156
00157
00158
00159
00160
00161
00162
00163
00164
00165
00166 subroutine shr_stream_init(strm,infoFile,yearFirst,yearLast,yearAlign,taxMode,rc)
00167
00168
00169
00170 type(shr_stream_streamType) ,intent(out) :: strm
00171 character(*) ,intent(in) :: infoFile
00172 integer (SHR_KIND_IN) ,intent(in) :: yearFirst
00173 integer (SHR_KIND_IN) ,intent(in) :: yearLast
00174 integer (SHR_KIND_IN) ,intent(in) :: yearAlign
00175 character(*) ,optional,intent(in) :: taxMode
00176 integer (SHR_KIND_IN),optional,intent(out) :: rc
00177
00178
00179
00180
00181 integer (SHR_KIND_IN) :: n
00182 character(SHR_KIND_CL) :: str
00183 integer (SHR_KIND_IN) :: int
00184 character(SHR_KIND_CL) :: subStr
00185 integer (SHR_KIND_IN) :: nUnit
00186 character(SHR_KIND_CS) :: startTag
00187 character(SHR_KIND_CS) :: endTag
00188 character(SHR_KIND_CS) :: fldNameFile
00189 character(SHR_KIND_CS) :: fldNameModel
00190 character(SHR_KIND_CX) :: fldListFile
00191 character(SHR_KIND_CX) :: fldListModel
00192 integer (SHR_KIND_IN) :: rCode, rCode2
00193
00194
00195 character(*),parameter :: subName = '(shr_stream_init) '
00196 character(*),parameter :: F00 = "('(shr_stream_init) ',8a)"
00197
00198
00199
00200
00201
00202
00203
00204
00205
00206
00207
00208 rCode = 0
00209 write(s_logunit,F00) 'Reading file ',trim(infoFile)
00210
00211
00212
00213
00214 strm%nFiles = 0
00215 do n=1,nFileMax
00216 strm%file%name = trim(shr_stream_file_null)
00217 strm%file%nt = 0
00218 strm%file%haveData = .false.
00219
00220
00221 end do
00222 strm%yearFirst = yearFirst
00223 strm%yearLast = yearLast
00224 strm%yearAlign = yearAlign
00225 if (present(taxMode)) then
00226 strm%taxMode = trim(taxMode)
00227 else
00228 strm%taxMode = trim(shr_stream_taxis_cycle)
00229 endif
00230 strm%fldListFile = ' '
00231 strm%fldListModel = ' '
00232 strm%filePath = ' '
00233 strm%domFilePath = ' '
00234 strm%k_lvd = -1
00235 strm%n_lvd = -1
00236 strm%found_lvd = .false.
00237 strm%k_gvd = -1
00238 strm%n_gvd = -1
00239 strm%found_gvd = .false.
00240
00241
00242 if (debug>0 .and. s_loglev>0) write(s_logunit,F00) ' reading data source'
00243
00244
00245 nUnit = shr_file_getUnit()
00246
00247
00248 startTag = "<dataSource>"
00249 endTag = "</dataSource>"
00250 open(nUnit,file=infoFile,STATUS='OLD',FORM='FORMATTED',ACTION='READ')
00251 call shr_stream_readUpToTag(nUnit,startTag,rc=rCode)
00252 if (rCode /= 0) goto 999
00253
00254
00255 read(nUnit,'(a)',END=999) str
00256 call shr_string_leftAlign(str)
00257 strm%dataSource = str
00258 if (debug>0 .and. s_loglev>0) write(s_logunit,F00) ' * format = ', trim(strm%dataSource)
00259
00260 close(nUnit)
00261 call shr_file_freeUnit(nUnit)
00262
00263
00264 if (debug>0 .and. s_loglev>0) write(s_logunit,F00) ' reading field data variable names'
00265
00266
00267
00268 open(nUnit,file=infoFile,STATUS='OLD',FORM='FORMATTED',ACTION='READ')
00269 startTag = "<fieldInfo>"
00270 endTag = "</fieldInfo>"
00271 call shr_stream_readUpToTag(nUnit,startTag,rc=rCode)
00272 if (rCode /= 0) goto 999
00273 startTag = "<variableNames>"
00274 endTag = "</variableNames>"
00275 call shr_stream_readUpToTag(nUnit,startTag,rc=rCode)
00276 if (rCode /= 0) goto 999
00277
00278
00279 n=0
00280 do while (.true.)
00281 read(nUnit,'(a)',END=999) str
00282 call shr_string_leftAlign(str)
00283 n=n+1
00284 if (str(1:len_trim(endTag)) == trim(endTag)) exit
00285 fldNameFile = ""
00286 fldNameModel = ""
00287 read(str,*,iostat=rCode) fldNameFile,fldNameModel
00288 if (len_trim(fldNameFile)==0 .or. len_trim(fldNameModel)==0 ) then
00289 rCode = 1
00290 write(s_logunit,F00) "ERROR: reading field names"
00291 write(s_logunit,F00) '* fldNameFile = ',trim(fldNameFile)
00292 write(s_logunit,F00) '* fldNameModel = ',trim(fldNameModel)
00293 call shr_stream_abort(subName//"ERROR: reading field names")
00294 end if
00295 if (n==1) then
00296 strm%fldListFile = trim(fldNameFile )
00297 strm%fldListModel = trim(fldNameModel)
00298 else
00299 strm%fldListFile = trim(strm%fldListFile ) // ":" // trim(fldNameFile )
00300 strm%fldListModel = trim(strm%fldListModel) // ":" // trim(fldNameModel)
00301 end if
00302 end do
00303 if (debug>0 .and. s_loglev>0) write(s_logunit,F00) ' * file field list = ',trim(strm%fldListFile )
00304 if (debug>0 .and. s_loglev>0) write(s_logunit,F00) ' * model field list = ',trim(strm%fldListModel)
00305 if (n==0) then
00306 rCode = 1
00307 write(s_logunit,F00) "ERROR: no input field names"
00308 call shr_stream_abort(subName//"ERROR: no input field names")
00309 end if
00310
00311 close(nUnit)
00312
00313
00314 if (debug>0 .and. s_loglev>0) write(s_logunit,F00) ' reading time-interpolation alogrithm '
00315
00316
00317
00318 open(nUnit,file=infoFile,STATUS='OLD',FORM='FORMATTED',ACTION='READ')
00319 startTag = "<fieldInfo>"
00320 endTag = "</fieldInfo>"
00321 call shr_stream_readUpToTag(nUnit,startTag,rc=rCode2)
00322 if (rCode2 /= 0)then
00323 rCode = rCode2
00324 goto 999
00325 end if
00326
00327
00328
00329
00330
00331
00332
00333
00334
00335
00336
00337
00338
00339 close(nUnit)
00340
00341
00342 if (debug>0 .and. s_loglev>0) write(s_logunit,F00) ' reading offset'
00343
00344
00345
00346 open(nUnit,file=infoFile,STATUS='OLD',FORM='FORMATTED',ACTION='READ')
00347 startTag = "<fieldInfo>"
00348 endTag = "</fieldInfo>"
00349 call shr_stream_readUpToTag(nUnit,startTag,rc=rCode2)
00350 if (rCode2 /= 0)then
00351 rCode = rCode2
00352 goto 999
00353 end if
00354 startTag = "<offset>"
00355 endTag = "</offset>"
00356 call shr_stream_readUpToTag(nUnit,startTag,optionalTag=.true.,rc=rCode2)
00357 if (rCode2 == 0) then
00358
00359 read(nUnit,*,END=999) int
00360 strm%offset = int
00361 else
00362 strm%offset = 0
00363 end if
00364 if (debug>0 .and. s_loglev>0) write(s_logunit,F00) ' * offset ',strm%offset
00365
00366 close(nUnit)
00367
00368
00369 if (debug>0 .and. s_loglev>0) write(s_logunit,F00) ' reading data file path'
00370
00371
00372
00373 open(nUnit,file=infoFile,STATUS='OLD',FORM='FORMATTED',ACTION='READ')
00374 startTag = "<fieldInfo>"
00375 endTag = "</fieldInfo>"
00376 call shr_stream_readUpToTag(nUnit,startTag,rc=rCode2)
00377 if (rCode2 /= 0)then
00378 rCode = rCode2
00379 goto 999
00380 end if
00381 startTag = "<filePath>"
00382 endTag = "</filePath>"
00383 call shr_stream_readUpToTag(nUnit,startTag,rc=rCode2)
00384 if (rCode2 /= 0)then
00385 rCode = rCode2
00386 goto 999
00387 end if
00388
00389
00390 read(nUnit,'(a)',END=999) str
00391 call shr_string_leftAlign(str)
00392 n = len_trim(str)
00393 if (n>0 .and. str(n:n) /= '/') str(n+1:n+2) = "/ "
00394 if (n==0) str = "./ "
00395 strm%FilePath = str
00396 if (debug>0 .and. s_loglev>0) write(s_logunit,F00) ' * data file path = ', trim(strm%FilePath)
00397
00398 close(nUnit)
00399
00400
00401 if (debug>0 .and. s_loglev>0) write(s_logunit,F00) ' reading field data file names'
00402
00403
00404
00405 open(nUnit,file=infoFile,STATUS='OLD',FORM='FORMATTED',ACTION='READ')
00406 startTag = "<fieldInfo>"
00407 endTag = "</fieldInfo>"
00408 call shr_stream_readUpToTag(nUnit,startTag,rc=rCode2)
00409 if (rCode2 /= 0)then
00410 rCode = rCode2
00411 goto 999
00412 end if
00413 startTag = "<fileNames>"
00414 endTag = "</fileNames>"
00415 call shr_stream_readUpToTag(nUnit,startTag,rc=rCode2)
00416 if (rCode2 /= 0)then
00417 rCode = rCode2
00418 goto 999
00419 end if
00420
00421
00422 n=0
00423 do while (.true.)
00424 read(nUnit,'(a)',END=999) str
00425 call shr_string_leftAlign(str)
00426 if (str(1:len_trim(endTag)) == trim(endTag)) exit
00427 n=n+1
00428 if (n > nFileMax) then
00429 rCode = 1
00430 write(s_logunit,F00) "ERROR: exceeded max number of files"
00431 call shr_stream_abort(subName//"ERROR: exceeded max number of files")
00432 if ( present(rc) ) rc = rCode
00433 close(nUnit)
00434 call shr_file_freeUnit(nUnit)
00435 return
00436 end if
00437 strm%file(n)%name = str
00438 if (debug>0 .and. s_loglev>0) write(s_logunit,F00) ' * ',trim(strm%file(n)%name)
00439 end do
00440 strm%nFiles = n
00441 if (n==0) then
00442 rCode = 1
00443 write(s_logunit,F00) "ERROR: no input file names"
00444 call shr_stream_abort(subName//"ERROR: no input file names")
00445 end if
00446
00447 close(nUnit)
00448
00449
00450 if (debug>0 .and. s_loglev>0) write(s_logunit,F00) ' reading domain data variable names'
00451
00452
00453
00454 open(nUnit,file=infoFile,STATUS='OLD',FORM='FORMATTED',ACTION='READ')
00455 startTag = "<domainInfo>"
00456 endTag = "</domainInfo>"
00457 call shr_stream_readUpToTag(nUnit,startTag,rc=rCode2)
00458 if (rCode2 /= 0)then
00459 rCode = rCode2
00460 goto 999
00461 end if
00462 startTag = "<variableNames>"
00463 endTag = "</variableNames>"
00464 call shr_stream_readUpToTag(nUnit,startTag,rc=rCode2)
00465 if (rCode2 /= 0)then
00466 rCode = rCode2
00467 goto 999
00468 end if
00469
00470
00471 n=0
00472 do while (.true.)
00473 read(nUnit,'(a)',END=999) str
00474 call shr_string_leftAlign(str)
00475 n=n+1
00476 if (str(1:len_trim(endTag)) == trim(endTag)) exit
00477 fldNameFile = ""
00478 fldNameModel = ""
00479 read(str,*,iostat=rCode2) fldNameFile,fldNameModel
00480 if (len_trim(fldNameFile)==0 .or. len_trim(fldNameModel)==0 ) then
00481 rCode = 1
00482 write(s_logunit,F00) "ERROR: reading field names"
00483 write(s_logunit,F00) '* fldNameFile = ',trim(fldNameFile)
00484 write(s_logunit,F00) '* fldNameModel = ',trim(fldNameModel)
00485 call shr_stream_abort(subName//"ERROR: reading field names")
00486 end if
00487 if (n==1) then
00488 fldListFile = trim(fldNameFile )
00489 fldListModel = trim(fldNameModel)
00490 else
00491 fldListFile = trim(fldListFile ) // ":" // trim(fldNameFile )
00492 fldListModel = trim(fldListModel) // ":" // trim(fldNameModel)
00493 end if
00494 end do
00495 if (debug>0 .and. s_loglev>0) write(s_logunit,F00) ' * file field list = ',trim(fldListFile )
00496 if (debug>0 .and. s_loglev>0) write(s_logunit,F00) ' * model field list = ',trim(fldListModel)
00497
00498 if (n==0) then
00499 rCode = 1
00500 write(s_logunit,F00) "ERROR: no input field names"
00501 call shr_stream_abort(subName//"ERROR: no input field names")
00502 else
00503
00504 n = shr_string_listGetIndexF(fldListModel,"time")
00505 if (n==0) then
00506 rCode = 1
00507 write(s_logunit,F00) "ERROR: no input field names"
00508 call shr_stream_abort(subName//"ERROR: no time variable name")
00509 else
00510 call shr_string_listGetName (fldListFile,n,substr,rc)
00511 strm%domTvarName = subStr
00512 endif
00513
00514 n = shr_string_listGetIndexF(fldListModel,"lon")
00515 if (n==0) then
00516 rCode = 1
00517 write(s_logunit,F00) "ERROR: no input field names"
00518 call shr_stream_abort(subName//"ERROR: no lon variable name")
00519 else
00520 call shr_string_listGetName (fldListFile,n,substr,rc)
00521 strm%domXvarName = subStr
00522 endif
00523
00524 n = shr_string_listGetIndexF(fldListModel,"lat")
00525 if (n==0) then
00526 rCode = 1
00527 write(s_logunit,F00) "ERROR: no input field names"
00528 call shr_stream_abort(subName//"ERROR: no time variable name")
00529 else
00530 call shr_string_listGetName (fldListFile,n,substr,rc)
00531 strm%domYvarName = subStr
00532 endif
00533
00534 n = shr_string_listGetIndexF(fldListModel,"area")
00535 if (n==0) then
00536 rCode = 1
00537 write(s_logunit,F00) "ERROR: no input field names"
00538 call shr_stream_abort(subName//"ERROR: no time variable name")
00539 else
00540 call shr_string_listGetName (fldListFile,n,substr,rc)
00541 strm%domAreaName = subStr
00542 endif
00543
00544 n = shr_string_listGetIndexF(fldListModel,"mask")
00545 if (n==0) then
00546 rCode = 1
00547 write(s_logunit,F00) "ERROR: no input field names"
00548 call shr_stream_abort(subName//"ERROR: no time variable name")
00549 else
00550 call shr_string_listGetName (fldListFile,n,substr,rc)
00551 strm%domMaskName = subStr
00552 endif
00553 end if
00554
00555 close(nUnit)
00556
00557
00558 if (debug>0 .and. s_loglev>0) write(s_logunit,F00) ' reading domain data file path'
00559
00560
00561
00562 startTag = "<domainInfo>"
00563 endTag = "</domainInfo>"
00564 open(nUnit,file=infoFile,STATUS='OLD',FORM='FORMATTED',ACTION='READ')
00565 call shr_stream_readUpToTag(nUnit,startTag,rc=rCode2)
00566 if (rCode2 /= 0)then
00567 rCode = rCode2
00568 goto 999
00569 end if
00570 startTag = "<filePath>"
00571 endTag = "</filePath>"
00572 call shr_stream_readUpToTag(nUnit,startTag,rc=rCode2)
00573 if (rCode2 /= 0)then
00574 rCode = rCode2
00575 goto 999
00576 end if
00577
00578
00579 read(nUnit,'(a)',END=999) str
00580 call shr_string_leftAlign(str)
00581 n = len_trim(str)
00582 if (n>0 .and. str(n:n) /= '/') str(n+1:n+2) = "/ "
00583 if (n==0) str = "./ "
00584 strm%domFilePath = str
00585 if (debug>0 .and. s_loglev>0) write(s_logunit,F00) ' * data file path = ', trim(strm%domFilePath)
00586
00587 close(nUnit)
00588
00589 if (debug>0 .and. s_loglev>0) write(s_logunit,F00) ' reading domain data file name'
00590
00591
00592
00593 open(nUnit,file=infoFile,STATUS='OLD',FORM='FORMATTED',ACTION='READ')
00594 startTag = "<domainInfo>"
00595 endTag = "</domainInfo>"
00596 call shr_stream_readUpToTag(nUnit,startTag,rc=rCode2)
00597 if (rCode2 /= 0)then
00598 rCode = rCode2
00599 goto 999
00600 end if
00601 startTag = "<fileNames>"
00602 endTag = "</fileNames>"
00603 call shr_stream_readUpToTag(nUnit,startTag,rc=rCode2)
00604 if (rCode2 /= 0)then
00605 rCode = rCode2
00606 goto 999
00607 end if
00608
00609
00610 read(nUnit,'(a)',END=999) str
00611 call shr_string_leftAlign(str)
00612 strm%domFileName = str
00613 if (debug>0 .and. s_loglev>0) write(s_logunit,F00) ' * ',trim(strm%domFileName)
00614
00615 close(nUnit)
00616
00617
00618
00619
00620 strm%init = .true.
00621 if ( present(rc) ) rc = rCode
00622 call shr_file_freeUnit(nUnit)
00623 return
00624
00625 999 continue
00626 write(s_logunit,F00) "ERROR: unexpected end-of-file while reading ",trim(startTag)
00627 write(s_logunit,F00) " error code = ", rCode
00628 call shr_stream_abort(subName//"ERROR: unexpected end-of-file")
00629 close(nUnit)
00630 if ( present(rc) ) rc = rCode
00631 call shr_file_freeUnit(nUnit)
00632
00633 end subroutine shr_stream_init
00634
00635
00636
00637
00638
00639
00640
00641
00642
00643
00644
00645
00646
00647
00648
00649
00650 subroutine shr_stream_set(strm,yearFirst,yearLast,yearAlign,offset,taxMode, &
00651 fldListFile,fldListModel,domFilePath,domFileName, &
00652 domTvarName,domXvarName,domYvarName,domAreaName,domMaskName, &
00653 filePath,filename,dataSource,rc)
00654
00655
00656
00657 type(shr_stream_streamType) ,intent(inout) :: strm
00658 integer (SHR_KIND_IN),optional,intent(in) :: yearFirst
00659 integer (SHR_KIND_IN),optional,intent(in) :: yearLast
00660 integer (SHR_KIND_IN),optional,intent(in) :: yearAlign
00661 integer (SHR_KIND_IN),optional,intent(in) :: offset
00662 character(*) ,optional,intent(in) :: taxMode
00663 character(*) ,optional,intent(in) :: fldListFile
00664 character(*) ,optional,intent(in) :: fldListModel
00665 character(*) ,optional,intent(in) :: domFilePath
00666 character(*) ,optional,intent(in) :: domFileName
00667 character(*) ,optional,intent(in) :: domTvarName
00668 character(*) ,optional,intent(in) :: domXvarName
00669 character(*) ,optional,intent(in) :: domYvarName
00670 character(*) ,optional,intent(in) :: domAreaName
00671 character(*) ,optional,intent(in) :: domMaskName
00672 character(*) ,optional,intent(in) :: filePath
00673 character(*) ,optional,intent(in) :: filename(:)
00674 character(*) ,optional,intent(in) :: dataSource
00675 integer (SHR_KIND_IN),optional,intent(out) :: rc
00676
00677
00678
00679
00680 integer(SHR_KIND_IN) :: n
00681
00682
00683 character(*),parameter :: subName = '(shr_stream_set) '
00684 character(*),parameter :: F00 = "('(shr_stream_set) ',8a)"
00685 character(*),parameter :: F01 = "('(shr_stream_set) ',1a,i6)"
00686
00687
00688
00689
00690
00691
00692 strm%nFiles = 0
00693 do n=1,nFileMax
00694 strm%file%name = trim(shr_stream_file_null)
00695 strm%file%nt = 0
00696 strm%file%haveData = .false.
00697
00698
00699 end do
00700 strm%yearFirst = yearFirst
00701 strm%yearLast = yearLast
00702 strm%yearAlign = yearAlign
00703 if (present(taxMode)) then
00704 strm%taxMode = trim(taxMode)
00705 else
00706 strm%taxMode = trim(shr_stream_taxis_cycle)
00707 endif
00708 strm%fldListFile = ' '
00709 strm%fldListModel = ' '
00710 strm%filePath = ' '
00711 strm%domFilePath = ' '
00712 strm%k_lvd = -1
00713 strm%n_lvd = -1
00714 strm%found_lvd = .false.
00715 strm%k_gvd = -1
00716 strm%n_gvd = -1
00717 strm%found_gvd = .false.
00718 strm%init = .true.
00719
00720 if ( present(rc) ) rc = 0
00721
00722 if (present(yearFirst)) then
00723 strm%yearFirst = yearFirst
00724 endif
00725 if (present(yearLast)) then
00726 strm%yearLast = yearLast
00727 endif
00728 if (present(yearAlign)) then
00729 strm%yearAlign = yearAlign
00730 endif
00731 if (present(offset)) then
00732 strm%offset = offset
00733 endif
00734 if (present(taxMode)) then
00735 strm%taxMode = trim(taxMode)
00736 endif
00737 if (present(fldListFile)) then
00738 strm%fldListFile = trim(fldListFile)
00739 endif
00740 if (present(fldListModel)) then
00741 strm%fldListModel = trim(fldListModel)
00742 endif
00743 if (present(domFilePath)) then
00744 strm%domFilePath = trim(domFilePath)
00745 endif
00746 if (present(domFileName)) then
00747 strm%domFileName = trim(domFileName)
00748 endif
00749 if (present(domTvarName)) then
00750 strm%domTvarName = trim(domTvarName)
00751 endif
00752 if (present(domXvarName)) then
00753 strm%domXvarName = trim(domXvarName)
00754 endif
00755 if (present(domYvarName)) then
00756 strm%domYvarName = trim(domYvarName)
00757 endif
00758 if (present(domAreaName)) then
00759 strm%domAreaName = trim(domAreaName)
00760 endif
00761 if (present(domMaskName)) then
00762 strm%domMaskName = trim(domMaskName)
00763 endif
00764 if (present(filePath)) then
00765 strm%filePath = trim(filePath)
00766 endif
00767 if (present(filename)) then
00768 write(s_logunit,F01) "size of filename = ",size(filename)
00769 write(s_logunit,F00) "filename = ",filename
00770 do n = 1,size(filename)
00771 if (trim(filename(n)) == trim(shr_stream_file_null)) then
00772
00773 else
00774 if (n > nFileMax) then
00775 write(s_logunit,F00) "ERROR: exceeded max number of files"
00776 call shr_stream_abort(subName//"ERROR: exceeded max number of files")
00777 if ( present(rc) ) rc = 1
00778 return
00779 endif
00780 strm%nFiles = n
00781 strm%file(n)%name = trim(filename(n))
00782 endif
00783 enddo
00784 endif
00785
00786 end subroutine shr_stream_set
00787
00788
00789 subroutine shr_stream_readUpToTag(nUnit,tag,optionalTag,rc)
00790
00791
00792 integer(SHR_KIND_IN),intent(in ) :: nUnit
00793 character(*) ,intent(in ) :: tag
00794 logical, optional ,intent(in ) :: optionalTag
00795 integer(SHR_KIND_IN),intent(out) :: rc
00796
00797
00798 character(SHR_KIND_CL) :: str
00799 logical :: localOptionalTag
00800
00801
00802 character(*),parameter :: subName = '(shr_stream_readUpToTag) '
00803 character(*),parameter :: F00 = "('(shr_stream_readUpToTag) ',8a)"
00804
00805
00806
00807
00808
00809 rc = 1
00810 localOptionalTag = .false.
00811 if (present(optionalTag)) localOptionalTag = optionalTag
00812 do while (.true.)
00813 read(nUnit,'(a)',END=999) str
00814 str = adjustL(str)
00815 if (str(1:len_trim(adjustL(tag))) == trim(adjustL(tag))) then
00816 rc = 0
00817 exit
00818 end if
00819 end do
00820
00821 999 continue
00822
00823 if (rc /= 0 .and. .not. localOptionalTag ) then
00824 write(s_logunit,F00) "ERROR: tag not found: ",trim(tag)
00825 call shr_stream_abort(subName//"ERROR: tag not found")
00826 end if
00827
00828 end subroutine shr_stream_readUpToTag
00829
00830
00831
00832
00833
00834
00835
00836
00837
00838
00839
00840
00841
00842
00843
00844
00845
00846 subroutine shr_stream_parseInput(str,fileName,yearAlign,yearFirst,yearLast,rc)
00847
00848
00849
00850 character(*) ,intent(in) :: str
00851 character(*) ,intent(out) :: fileName
00852 integer (SHR_KIND_IN) ,intent(out) :: yearFirst
00853 integer (SHR_KIND_IN) ,intent(out) :: yearLast
00854 integer (SHR_KIND_IN) ,intent(out) :: yearAlign
00855 integer (SHR_KIND_IN),optional,intent(out) :: rc
00856
00857
00858
00859
00860 integer (SHR_KIND_IN) :: n
00861 character(SHR_KIND_CL) :: str2
00862
00863
00864 character(*),parameter :: subName = '(shr_stream_parseInput) '
00865 character(*),parameter :: F00 = "('(shr_stream_parseInput) ',8a)"
00866 character(*),parameter :: F01 = "('(shr_stream_parseInput) ',a,3i10)"
00867
00868
00869
00870
00871
00872
00873
00874
00875 if (debug>1 .and. s_loglev > 0) write(s_logunit,F00) "str = ",trim(str)
00876
00877 str2 = adjustL(str)
00878 n = index(str2," ")
00879 fileName = str2(:n)
00880 read(str2(n:),*) yearAlign,yearFirst,yearLast
00881
00882 if (debug>1 .and. s_loglev > 0) then
00883 write(s_logunit,F00) "fileName = ",trim(fileName)
00884 write(s_logunit,F01) "yearAlign = ",yearAlign
00885 write(s_logunit,F01) "yearFirst = ",yearFirst
00886 write(s_logunit,F01) "yearLast = ",yearLast
00887 end if
00888
00889 if (present(rc)) rc = 0
00890
00891 end subroutine shr_stream_parseInput
00892
00893
00894
00895
00896
00897
00898
00899
00900
00901
00902
00903
00904
00905
00906
00907
00908
00909
00910 subroutine shr_stream_findBounds(strm,mDateIn, secIn, &
00911 & mDateLB,dDateLB,secLB,n_lb,fileLB, &
00912 & mDateUB,dDateUB,secUB,n_ub,fileUB )
00913
00914
00915
00916 type(shr_stream_streamType),intent(inout):: strm
00917 integer(SHR_KIND_IN) ,intent(in) :: mDateIn
00918 integer(SHR_KIND_IN) ,intent(in) :: secIn
00919 integer(SHR_KIND_IN) ,intent(out) :: mDateLB
00920 integer(SHR_KIND_IN) ,intent(out) :: dDateLB
00921 integer(SHR_KIND_IN) ,intent(out) :: secLB
00922 integer(SHR_KIND_IN) ,intent(out) :: n_lb
00923 character(*) ,intent(out) :: fileLB
00924 integer(SHR_KIND_IN) ,intent(out) :: mDateUB
00925 integer(SHR_KIND_IN) ,intent(out) :: dDateUB
00926 integer(SHR_KIND_IN) ,intent(out) :: secUB
00927 integer(SHR_KIND_IN) ,intent(out) :: n_ub
00928 character(*) ,intent(out) :: fileUB
00929
00930
00931
00932
00933 character(SHR_KIND_CL) :: fileName
00934 integer (SHR_KIND_IN) :: nt
00935 integer (SHR_KIND_IN) :: dDateIn
00936 integer (SHR_KIND_IN) :: dDateF
00937 integer (SHR_KIND_IN) :: dDateL
00938 integer (SHR_KIND_IN) :: n,nf
00939 integer (SHR_KIND_IN) :: k,kf
00940 integer (SHR_KIND_IN) :: k_ub,k_lb
00941 integer (SHR_KIND_IN) :: rCode
00942
00943 integer (SHR_KIND_IN) :: mYear
00944 integer (SHR_KIND_IN) :: yrFirst
00945 integer (SHR_KIND_IN) :: yrLast
00946 integer (SHR_KIND_IN) :: yrAlign
00947 integer (SHR_KIND_IN) :: nYears
00948 integer (SHR_KIND_IN) :: dYear
00949 integer (SHR_KIND_IN) :: yy,mm,dd
00950 real (SHR_KIND_R8) :: rDateIn
00951 real (SHR_KIND_R8) :: rDate1
00952 real (SHR_KIND_R8) :: rDate2
00953 real (SHR_KIND_R8) :: rDatelvd
00954 real (SHR_KIND_R8) :: rDategvd
00955 logical :: cycle
00956 logical :: limit
00957
00958
00959 character(*),parameter :: subName = '(shr_stream_findBounds) '
00960 character(*),parameter :: F00 = "('(shr_stream_findBounds) ',8a)"
00961 character(*),parameter :: F01 = "('(shr_stream_findBounds) ',a,i9.8,a)"
00962 character(*),parameter :: F02 = "('(shr_stream_findBounds) ',a,2i9.8,i6,i5,1x,a)"
00963 character(*),parameter :: F03 = "('(shr_stream_findBounds) ',a,i4)"
00964 character(*),parameter :: F04 = "('(shr_stream_findBounds) ',2a,i4)"
00965
00966
00967
00968
00969
00970
00971
00972
00973 if (debug>0 .and. s_loglev > 0) write(s_logunit,F02) "DEBUG: ---------- enter ------------------"
00974
00975 rCode = 0
00976 if ( .not. strm%init ) then
00977 rCode = 1
00978 call shr_stream_abort(trim(subName)//" ERROR: trying to find bounds of uninitialized stream")
00979 return
00980 end if
00981
00982 if (trim(strm%taxMode) == trim(shr_stream_taxis_cycle)) then
00983 cycle = .true.
00984 limit = .false.
00985 elseif (trim(strm%taxMode) == trim(shr_stream_taxis_extend)) then
00986 cycle = .false.
00987 limit = .false.
00988 elseif (trim(strm%taxMode) == trim(shr_stream_taxis_limit)) then
00989 cycle = .false.
00990 limit = .true.
00991 else
00992 write(s_logunit,*) trim(subName),' ERROR: illegal taxMode = ',trim(strm%taxMode)
00993 call shr_stream_abort(trim(subName)//' ERROR: illegal taxMode = '//trim(strm%taxMode))
00994 endif
00995
00996
00997
00998
00999
01000 mYear = mDateIn/10000
01001 yrFirst = strm%yearFirst
01002 yrLast = strm%yearLast
01003 yrAlign = strm%yearAlign
01004 nYears = yrLast - yrFirst + 1
01005 dDateF = yrFirst * 10000 + 101
01006 dDateL = (yrLast+1) * 10000 + 101
01007
01008 if (cycle) then
01009 dYear = yrFirst + modulo(mYear-yrAlign+(2*nYears),nYears)
01010 else
01011 dYear = yrFirst + mYear - yrAlign
01012 endif
01013
01014 if (dYear < 0) then
01015 write(s_logunit,*) trim(subName),' ERROR: dyear lt zero = ',dYear
01016 call shr_stream_abort(trim(subName)//' ERROR: dyear lt one')
01017 endif
01018
01019 dDateIn = dYear*10000 + modulo(mDateIn,10000)
01020 rDateIn = dDateIn + secIn/spd
01021
01022
01023
01024
01025
01026
01027
01028
01029
01030 if (.not. strm%found_lvd) then
01031 A: do k=1,strm%nFiles
01032 if (.not. strm%file(k)%haveData) then
01033 call shr_stream_readtCoord(strm, k, rCode)
01034 if ( rCode /= 0 )then
01035 call shr_stream_abort(trim(subName)//" ERROR: readtCoord1")
01036 return
01037 end if
01038 end if
01039 do n=1,strm%file(k)%nt
01040 if ( dDateF <= strm%file(k)%date(n) ) then
01041
01042 strm%k_lvd = k
01043 strm%n_lvd = n
01044 strm%found_lvd = .true.
01045 exit A
01046 end if
01047 end do
01048 end do A
01049 if (.not. strm%found_lvd) then
01050 rCode = 1
01051 write(s_logunit,F00) "ERROR: LVD not found, all data is before yearFirst"
01052 call shr_stream_abort(trim(subName)//" ERROR: LVD not found, all data is before yearFirst")
01053 else
01054
01055 if ( dDateL <= strm%file(strm%k_lvd)%date(strm%n_lvd) ) then
01056 rCode = 1
01057 write(s_logunit,F00) "ERROR: LVD not found, all data is after yearLast"
01058 call shr_stream_abort(trim(subName)//" ERROR: LVD not found, all data is after yearLast")
01059 end if
01060 end if
01061 if (debug>1 .and. s_loglev > 0) then
01062 if (strm%found_lvd) write(s_logunit,F01) "DEBUG: found LVD = ",strm%file(k)%date(n)
01063 end if
01064 end if
01065
01066 if (strm%found_lvd) then
01067 k = strm%k_lvd
01068 n = strm%n_lvd
01069 rDatelvd = strm%file(k)%date(n) + strm%file(k)%secs(n)/spd
01070 else
01071 write(s_logunit,F00) "ERROR: LVD not found yet"
01072 call shr_stream_abort(trim(subName)//" ERROR: LVD not found yet")
01073 endif
01074
01075 if (strm%found_gvd) then
01076 k = strm%k_gvd
01077 n = strm%n_gvd
01078 rDategvd = strm%file(k)%date(n) + strm%file(k)%secs(n)/spd
01079 else
01080 rDategvd = 99991231.0
01081 endif
01082
01083
01084
01085
01086
01087
01088
01089
01090
01091
01092
01093 if (rDateIn < rDatelvd) then
01094 if (limit) then
01095 write(s_logunit,*) trim(subName)," ERROR: limit on and rDateIn lt rDatelvd",rDateIn,rDatelvd
01096 call shr_stream_abort(trim(subName)//" ERROR: rDateIn lt rDatelvd limit true")
01097 return
01098 endif
01099
01100 if (.not.cycle) then
01101 k_lb = strm%k_lvd
01102 n_lb = strm%n_lvd
01103 dDateLB = 00000101
01104 mDateLB = 00000101
01105 secLB = 0
01106 fileLB = strm%file(k_lb)%name
01107
01108 k_ub = strm%k_lvd
01109 n_ub = strm%n_lvd
01110 dDateUB = strm%file(k_ub)%date(n_ub)
01111 call shr_cal_date2ymd(dDateUB,yy,mm,dd)
01112 yy = yy + (mYear-dYear)
01113 call shr_cal_ymd2date(yy,mm,dd,mDateUB)
01114 secUB = strm%file(k_ub)%secs(n_ub)
01115 fileUB = strm%file(k_ub)%name
01116
01117
01118 return
01119 endif
01120
01121 if (cycle) then
01122
01123 if (.not. strm%found_gvd) then
01124
01125 B: do k=strm%nFiles,1,-1
01126
01127 if (.not. strm%file(k)%haveData) then
01128 call shr_stream_readtCoord(strm, k, rCode)
01129 if ( rCode /= 0 )then
01130 call shr_stream_abort(trim(subName)//" ERROR: readtCoord2")
01131 return
01132 end if
01133 end if
01134
01135 do n=strm%file(k)%nt,1,-1
01136 if ( strm%file(k)%date(n) < dDateL ) then
01137 strm%k_gvd = k
01138 strm%n_gvd = n
01139 strm%found_gvd = .true.
01140 rDategvd = strm%file(k)%date(n) + strm%file(k)%secs(n)/spd
01141 if (debug>1 .and. s_loglev > 0) write(s_logunit,F01) "DEBUG: found GVD ",strm%file(k)%date(n)
01142 exit B
01143 end if
01144 end do
01145 end do B
01146 end if
01147
01148 if (.not. strm%found_gvd) then
01149 write(s_logunit,F00) "ERROR: GVD not found1"
01150 call shr_stream_abort(trim(subName)//" ERROR: GVD not found1")
01151 endif
01152
01153 k_lb = strm%k_gvd
01154 n_lb = strm%n_gvd
01155 dDateLB = strm%file(k_lb)%date(n_lb)
01156 call shr_cal_date2ymd(dDateLB,yy,mm,dd)
01157 yy = yy + (mYear-dYear-nYears)
01158 call shr_cal_ymd2date(yy,mm,dd,mDateLB)
01159 secLB = strm%file(k_lb)%secs(n_lb)
01160 fileLB = strm%file(k_lb)%name
01161
01162 k_ub = strm%k_lvd
01163 n_ub = strm%n_lvd
01164 dDateUB = strm%file(k_ub)%date(n_ub)
01165 call shr_cal_date2ymd(dDateUB,yy,mm,dd)
01166 yy = yy + (mYear-dYear)
01167 call shr_cal_ymd2date(yy,mm,dd,mDateUB)
01168 secUB = strm%file(k_ub)%secs(n_ub)
01169 fileUB = strm%file(k_ub)%name
01170
01171
01172 return
01173 endif
01174
01175
01176
01177
01178
01179
01180
01181
01182 else if (strm%found_gvd .and. rDateIn >= rDategvd) then
01183 if (limit) then
01184 write(s_logunit,*) trim(subName)," ERROR: limit on and rDateIn gt rDategvd",rDateIn,rDategvd
01185 call shr_stream_abort(trim(subName)//" ERROR: rDateIn gt rDategvd limit true")
01186 return
01187 endif
01188
01189 if (.not.cycle) then
01190 k_lb = strm%k_gvd
01191 n_lb = strm%n_gvd
01192 dDateLB = strm%file(k_lb)%date(n_lb)
01193 call shr_cal_date2ymd(dDateLB,yy,mm,dd)
01194 yy = yy + (mYear-dYear)
01195 call shr_cal_ymd2date(yy,mm,dd,mDateLB)
01196 secLB = strm%file(k_lb)%secs(n_lb)
01197 fileLB = strm%file(k_lb)%name
01198
01199 k_ub = strm%k_gvd
01200 n_ub = strm%n_gvd
01201 dDateUB = 99991231
01202 mDateUB = 99991231
01203 secUB = 0
01204 fileUB = strm%file(k_ub)%name
01205
01206
01207 return
01208 endif
01209
01210 if (cycle) then
01211 k_lb = strm%k_gvd
01212 n_lb = strm%n_gvd
01213 dDateLB = strm%file(k_lb)%date(n_lb)
01214 call shr_cal_date2ymd(dDateLB,yy,mm,dd)
01215 yy = yy + (mYear-dYear)
01216 call shr_cal_ymd2date(yy,mm,dd,mDateLB)
01217 secLB = strm%file(k_lb)%secs(n_lb)
01218 fileLB = strm%file(k_lb)%name
01219
01220 k_ub = strm%k_lvd
01221 n_ub = strm%n_lvd
01222 dDateUB = strm%file(k_ub)%date(n_ub)
01223 call shr_cal_date2ymd(dDateUB,yy,mm,dd)
01224 yy = yy + (mYear-dYear+nYears)
01225 call shr_cal_ymd2date(yy,mm,dd,mDateUB)
01226 secUB = strm%file(k_ub)%secs(n_ub)
01227 fileUB = strm%file(k_ub)%name
01228
01229
01230 return
01231 endif
01232
01233 else
01234
01235
01236
01237
01238 k_lb = strm%k_lvd
01239 n_lb = strm%n_lvd
01240 C: do k=strm%k_lvd,strm%nFiles
01241
01242 if (.not. strm%file(k)%haveData) then
01243 call shr_stream_readtCoord(strm, k, rCode)
01244 if ( rCode /= 0 )then
01245 call shr_stream_abort(trim(subName)//" ERROR: readtCoord3")
01246 return
01247 end if
01248 end if
01249
01250 n = strm%file(k)%nt
01251 rDate1 = strm%file(k)%date(n) + strm%file(k)%secs(n)/spd
01252
01253 if (.not. strm%found_gvd) then
01254 n = strm%file(k)%nt
01255 if (dDateL <= strm%file(k)%date(n)) then
01256
01257 if (k > 1) then
01258 strm%k_gvd = k-1
01259 strm%n_gvd = strm%file(k-1)%nt
01260 strm%found_gvd = .true.
01261 endif
01262 do n=1,strm%file(k)%nt
01263 if ( strm%file(k)%date(n) < dDateL ) then
01264 strm%k_gvd = k
01265 strm%n_gvd = n
01266 strm%found_gvd = .true.
01267 endif
01268 enddo
01269 elseif (k == strm%nFiles) then
01270 strm%k_gvd = k
01271 strm%n_gvd = strm%file(k)%nt
01272 strm%found_gvd = .true.
01273 end if
01274 if (strm%found_gvd) then
01275 kf = strm%k_gvd
01276 nf = strm%n_gvd
01277 rDategvd = strm%file(kf)%date(nf) + strm%file(kf)%secs(nf)/spd
01278 endif
01279 end if
01280
01281
01282
01283
01284
01285
01286
01287
01288 if (strm%found_gvd .and. rDateIn >= rDategvd) then
01289 if (limit) then
01290 write(s_logunit,*) trim(subName)," ERROR: limit on and rDateIn gt rDategvd",rDateIn,rDategvd
01291 call shr_stream_abort(trim(subName)//" ERROR: rDateIn gt rDategvd limit true")
01292 return
01293 endif
01294
01295 if (.not.cycle) then
01296 k_lb = strm%k_gvd
01297 n_lb = strm%n_gvd
01298 dDateLB = strm%file(k_lb)%date(n_lb)
01299 call shr_cal_date2ymd(dDateLB,yy,mm,dd)
01300 yy = yy + (mYear-dYear)
01301 call shr_cal_ymd2date(yy,mm,dd,mDateLB)
01302 secLB = strm%file(k_lb)%secs(n_lb)
01303 fileLB = strm%file(k_lb)%name
01304
01305 k_ub = strm%k_gvd
01306 n_ub = strm%n_gvd
01307 dDateUB = 99991231
01308 mDateUB = 99991231
01309 secUB = 0
01310 fileUB = strm%file(k_ub)%name
01311
01312
01313 return
01314 endif
01315
01316 if (cycle) then
01317 k_lb = strm%k_gvd
01318 n_lb = strm%n_gvd
01319 dDateLB = strm%file(k_lb)%date(n_lb)
01320 call shr_cal_date2ymd(dDateLB,yy,mm,dd)
01321 yy = yy + (mYear-dYear)
01322 call shr_cal_ymd2date(yy,mm,dd,mDateLB)
01323 secLB = strm%file(k_lb)%secs(n_lb)
01324 fileLB = strm%file(k_lb)%name
01325
01326 k_ub = strm%k_lvd
01327 n_ub = strm%n_lvd
01328 dDateUB = strm%file(k_ub)%date(n_ub)
01329 call shr_cal_date2ymd(dDateUB,yy,mm,dd)
01330 yy = yy + (mYear-dYear+nYears)
01331 call shr_cal_ymd2date(yy,mm,dd,mDateUB)
01332 secUB = strm%file(k_ub)%secs(n_ub)
01333 fileUB = strm%file(k_ub)%name
01334
01335
01336 return
01337 endif
01338
01339 endif
01340
01341 if ( rDate1 < rDateIn ) then
01342
01343 k_lb = k
01344 n_lb = strm%file(k)%nt
01345 else
01346
01347 do n=1,strm%file(k)%nt
01348 rDate2 = strm%file(k)%date(n) + strm%file(k)%secs(n)/spd
01349 if ( rDate2 <= rDateIn ) then
01350
01351 k_lb = k
01352 n_lb = n
01353 else
01354
01355 k_ub = k
01356 n_ub = n
01357
01358 dDateLB = strm%file(k_lb)%date(n_lb)
01359 call shr_cal_date2ymd(dDateLB,yy,mm,dd)
01360 yy = yy + (mYear-dYear)
01361 call shr_cal_ymd2date(yy,mm,dd,mDateLB)
01362 secLB = strm%file(k_lb)%secs(n_lb)
01363 fileLB = strm%file(k_lb)%name
01364
01365 dDateUB = strm%file(k_ub)%date(n_ub)
01366 call shr_cal_date2ymd(dDateUB,yy,mm,dd)
01367 yy = yy + (mYear-dYear)
01368 call shr_cal_ymd2date(yy,mm,dd,mDateUB)
01369 secUB = strm%file(k_ub)%secs(n_ub)
01370 fileUB = strm%file(k_ub)%name
01371
01372
01373 return
01374 endif
01375 enddo
01376 endif
01377 end do C
01378 endif
01379
01380 call shr_stream_abort(trim(subName)//' ERROR: findBounds failed')
01381 return
01382
01383 end subroutine shr_stream_findBounds
01384
01385
01386
01387
01388
01389
01390
01391
01392
01393
01394
01395
01396
01397
01398
01399 subroutine shr_stream_readTCoord(strm,k,rc)
01400
01401 use shr_cal_mod, only: shr_cal_advDate
01402 use netcdf
01403
01404 implicit none
01405
01406
01407
01408 type(shr_stream_streamType) ,intent(inout) :: strm
01409 integer(SHR_KIND_IN) ,intent(in) :: k
01410 integer(SHR_KIND_IN),optional,intent(out) :: rc
01411
01412
01413
01414
01415 character(SHR_KIND_CL) :: fileName
01416 integer(SHR_KIND_IN) :: nt
01417 integer(SHR_KIND_IN) :: num,n
01418 integer(SHR_KIND_IN) :: din,dout
01419 real(SHR_KIND_R8) :: sin,sout,offin
01420 integer(SHR_KIND_IN) :: lrc
01421 integer(SHR_KIND_IN) :: fid,vid,ndims,rcode
01422 integer(SHR_KIND_IN),allocatable :: dids(:)
01423 character(SHR_KIND_CL) :: units,calendar
01424 character(SHR_KIND_CS) :: bunits
01425 integer(SHR_KIND_IN) :: bdate
01426 real(SHR_KIND_R8) :: bsec
01427 integer(SHR_KIND_IN) :: ndate
01428 real(SHR_KIND_R8) :: nsec
01429 real(SHR_KIND_R8),allocatable :: tvar(:)
01430
01431 character(*),parameter :: subname = '(shr_stream_readTCoord) '
01432 character(*),parameter :: F01 = "('(shr_stream_readTCoord) ',a,2i7)"
01433
01434
01435
01436 lrc = 0
01437
01438
01439 call shr_stream_getFile(strm%filePath,strm%file(k)%name,fileName)
01440 rCode = nf90_open(fileName,nf90_nowrite,fid)
01441 if (rcode /= nf90_noerr) call shr_sys_abort(subname//' ERROR: nf90_open')
01442 rCode = nf90_inq_varid(fid,trim(strm%domTvarName),vid)
01443 if (rcode /= nf90_noerr) call shr_sys_abort(subname//' ERROR: nf90_inq_varid')
01444 rCode = nf90_inquire_variable(fid,vid,ndims=ndims)
01445 if (rcode /= nf90_noerr) call shr_sys_abort(subname//' ERROR: nf90_inquire_variable1')
01446 allocate(dids(ndims))
01447 rCode = nf90_inquire_variable(fid,vid,dimids=dids)
01448 if (rcode /= nf90_noerr) call shr_sys_abort(subname//' ERROR: nf90_inquire_variable2')
01449 rCode = nf90_inquire_dimension(fid,dids(1),len=nt)
01450 if (rcode /= nf90_noerr) call shr_sys_abort(subname//' ERROR: nf90_inquire_dimension')
01451 deallocate(dids)
01452
01453 allocate(strm%file(k)%date(nt),strm%file(k)%secs(nt))
01454 strm%file(k)%nt = nt
01455
01456 units = ' '
01457 calendar = ' '
01458 rCode = nf90_get_att(fid, vid, 'units', units)
01459 if (rcode /= nf90_noerr) call shr_sys_abort(subname//' ERROR: nf90_get_att units')
01460 rCode = nf90_inquire_attribute(fid, vid, 'calendar')
01461 if (rCode == nf90_noerr) then
01462 rCode = nf90_get_att(fid, vid, 'calendar', calendar)
01463 if (rcode /= nf90_noerr) call shr_sys_abort(subname//' ERROR: nf90_get_att calendar')
01464 else
01465 calendar = 'gregorian'
01466 endif
01467 n = len_trim(units)
01468 if (ichar(units(n:n)) == 0 ) units(n:n) = ' '
01469 n = len_trim(calendar)
01470 if (ichar(calendar(n:n)) == 0 ) calendar(n:n) = ' '
01471 call shr_string_leftalign(units)
01472 call shr_string_leftalign(calendar)
01473 call shr_string_parseCFtunit(units,bunits,bdate,bsec)
01474
01475 allocate(tvar(nt))
01476 rcode = nf90_get_var(fid,vid,tvar)
01477 if (rcode /= nf90_noerr) call shr_sys_abort(subname//' ERROR: nf90_get_var')
01478 rCode = nf90_close(fid)
01479 if (rcode /= nf90_noerr) call shr_sys_abort(subname//' ERROR: nf90_close')
01480 do n = 1,nt
01481 call shr_cal_advDate(tvar(n),bunits,bdate,bsec,ndate,nsec,calendar)
01482 strm%file(k)%date(n) = ndate
01483 strm%file(k)%secs(n) = nint(nsec)
01484 enddo
01485 deallocate(tvar)
01486
01487 if (strm%offset /= 0) then
01488 if (size(strm%file(k)%date) /= size(strm%file(k)%secs)) then
01489
01490 write(s_logunit,F01) "Incompatable date and secs sizes",size(strm%file(k)%date),size(strm%file(k)%secs)
01491 call shr_sys_abort()
01492 endif
01493 num = size(strm%file(k)%date)
01494 offin = strm%offset * 1.0_SHR_KIND_R8
01495 do n = 1,num
01496 din = strm%file(k)%date(n)
01497 sin = strm%file(k)%secs(n) * 1.0_SHR_KIND_R8
01498 call shr_cal_advDate(offin,'seconds',din,sin,dout,sout)
01499
01500 strm%file(k)%date(n) = dout
01501 strm%file(k)%secs(n) = sout
01502 enddo
01503 endif
01504
01505 strm%file(k)%haveData = .true.
01506 call shr_stream_verifyTCoord(strm,k,lrc)
01507
01508 if (present(rc)) then
01509 rc = lrc
01510 endif
01511
01512 end subroutine shr_stream_readTCoord
01513
01514
01515
01516
01517
01518
01519
01520
01521
01522
01523
01524
01525
01526
01527 subroutine shr_stream_verifyTCoord(strm,k,rc)
01528
01529
01530
01531 type(shr_stream_streamType),intent(in) :: strm
01532 integer(SHR_KIND_IN) :: k
01533 integer(SHR_KIND_IN) :: rc
01534
01535
01536
01537
01538 integer(SHR_KIND_IN) :: n
01539 integer(SHR_KIND_IN) :: nt
01540 integer(SHR_KIND_IN) :: date1,secs1
01541 integer(SHR_KIND_IN) :: date2,secs2
01542 logical :: checkIt
01543
01544
01545 character(*),parameter :: subName = '(shr_stream_verifyTCoord) '
01546 character(*),parameter :: F00 = "('(shr_stream_verifyTCoord) ',8a)"
01547 character(*),parameter :: F01 = "('(shr_stream_verifyTCoord) ',a,2i7)"
01548 character(*),parameter :: F02 = "('(shr_stream_verifyTCoord) ',a,2i9.8)"
01549
01550
01551
01552
01553
01554
01555
01556
01557
01558
01559
01560 rc = 0
01561 if (debug>1 .and. s_loglev > 0) write(s_logunit,F01) "checking t-coordinate data for file k =",k
01562
01563 if ( .not. strm%file(k)%haveData) then
01564 rc = 1
01565 write(s_logunit,F01) "Don't have data for file ",k
01566 call shr_stream_abort(subName//"ERROR: can't check -- file not read.")
01567 return
01568 end if
01569
01570 do n=1,strm%file(k)%nt+1
01571 checkIt = .false.
01572
01573
01574 if (n==1) then
01575
01576 if (k>1) then
01577 if ( strm%file(k-1)%haveData ) then
01578 nt = strm%file(k-1)%nt
01579 date1 = strm%file(k-1)%date(nt)
01580 secs1 = strm%file(k-1)%secs(nt)
01581 date2 = strm%file(k )%date(n)
01582 secs2 = strm%file(k )%secs(n)
01583 checkIt = .true.
01584 if (debug>1 .and. s_loglev > 0) write(s_logunit,F01) "comparing with previous file for file k =",k
01585 end if
01586 end if
01587 else if (n==strm%file(k)%nt+1) then
01588
01589 if (k<strm%nFiles) then
01590 if ( strm%file(k+1)%haveData ) then
01591 nt = strm%file(k )%nt
01592 date1 = strm%file(k )%date(nt)
01593 secs1 = strm%file(k )%secs(nt)
01594 date2 = strm%file(k+1)%date(1)
01595 secs2 = strm%file(k+1)%secs(1)
01596 checkIt = .true.
01597 if (debug>1 .and. s_loglev > 0) write(s_logunit,F01) "comparing with next file for file k =",k
01598 end if
01599 end if
01600 else
01601
01602 date1 = strm%file(k)%date(n-1)
01603 secs1 = strm%file(k)%secs(n-1)
01604 date2 = strm%file(k)%date(n )
01605 secs2 = strm%file(k)%secs(n )
01606 checkIt = .true.
01607 end if
01608
01609
01610 if (checkIt) then
01611 if ( date1 > date2 ) then
01612 rc = 1
01613 write(s_logunit,F01) "ERROR: calendar dates must be increasing"
01614 write(s_logunit,F02) "date(n), date(n+1) = ",date1,date2
01615 call shr_stream_abort(subName//"ERROR: calendar dates must be increasing")
01616 return
01617 else if ( date1 == date2 ) then
01618 if ( secs1 >= secs2 ) then
01619 rc = 1
01620 write(s_logunit,F01) "ERROR: elapsed seconds on a date must be strickly increasing"
01621 write(s_logunit,F02) "secs(n), secs(n+1) = ",secs1,secs2
01622 call shr_stream_abort(subName//"ERROR: elapsed seconds must be increasing")
01623 return
01624 end if
01625 end if
01626 if ( secs1 < 0 .or. spd < secs1 ) then
01627 rc = 1
01628 write(s_logunit,F01) "ERROR: elapsed seconds out of valid range [0,spd]"
01629 write(s_logunit,F02) "secs(n) = ",secs1
01630 call shr_stream_abort(subName//"ERROR: elapsed seconds out of range")
01631 return
01632 end if
01633 end if
01634 end do
01635
01636 if (debug>0 .and. s_loglev > 0) write(s_logunit,F01) "data is OK (non-decreasing) for file k =",k
01637
01638 end subroutine shr_stream_verifyTCoord
01639
01640
01641
01642
01643
01644
01645
01646
01647
01648
01649
01650
01651
01652
01653
01654
01655 subroutine shr_stream_getFileFieldList(stream,list,rc)
01656
01657 implicit none
01658
01659
01660
01661 type(shr_stream_streamType) ,intent(in) :: stream
01662 character(*) ,intent(out) :: list
01663 integer(SHR_KIND_IN),optional,intent(out) :: rc
01664
01665
01666
01667
01668 integer(SHR_KIND_IN) :: rCode
01669
01670
01671 character(*),parameter :: subName = '(shr_stream_getFileFieldList) '
01672 character(*),parameter :: F00 = "('(shr_stream_getFileFieldList) ',4a)"
01673
01674
01675
01676
01677
01678 rCode = 0
01679
01680 list = stream%fldListFile
01681
01682 if (present(rc)) rc = rCode
01683
01684 end subroutine shr_stream_getFileFieldList
01685
01686
01687
01688
01689
01690
01691
01692
01693
01694
01695
01696
01697
01698
01699
01700
01701 subroutine shr_stream_getModelFieldList(stream,list,rc)
01702
01703 implicit none
01704
01705
01706
01707 type(shr_stream_streamType) ,intent(in) :: stream
01708 character(*) ,intent(out) :: list
01709 integer(SHR_KIND_IN),optional,intent(out) :: rc
01710
01711
01712
01713
01714 integer(SHR_KIND_IN) :: rCode
01715
01716
01717 character(*),parameter :: subName = '(shr_stream_getModelFieldList) '
01718 character(*),parameter :: F00 = "('(shr_stream_getModelFieldList) ',4a)"
01719
01720
01721
01722
01723
01724 rCode = 0
01725
01726 list = stream%fldListModel
01727
01728 if (present(rc)) rc = rCode
01729
01730 end subroutine shr_stream_getModelFieldList
01731
01732
01733
01734
01735
01736
01737
01738
01739
01740
01741
01742
01743
01744
01745
01746
01747 subroutine shr_stream_getFileFieldName(stream,k,name,rc)
01748
01749 implicit none
01750
01751
01752
01753 type(shr_stream_streamType) ,intent(in) :: stream
01754 integer(SHR_KIND_IN) ,intent(in) :: k
01755 character(*) ,intent(out) :: name
01756 integer(SHR_KIND_IN),optional,intent(out) :: rc
01757
01758
01759
01760
01761 integer(SHR_KIND_IN) :: rCode
01762
01763
01764 character(*),parameter :: subName = '(shr_stream_getFileFieldName) '
01765 character(*),parameter :: F00 = "('(shr_stream_getFileFieldName) ',4a)"
01766
01767
01768
01769
01770
01771 rCode = 0
01772
01773 call shr_string_listGetName(stream%fldListFile,k,name,rCode)
01774
01775 if (present(rc)) rc = rCode
01776
01777 end subroutine shr_stream_getFileFieldName
01778
01779
01780
01781
01782
01783
01784
01785
01786
01787
01788
01789
01790
01791
01792
01793
01794 subroutine shr_stream_getModelFieldName(stream,k,name,rc)
01795
01796 implicit none
01797
01798
01799
01800 type(shr_stream_streamType) ,intent(in) :: stream
01801 integer(SHR_KIND_IN) ,intent(in) :: k
01802 character(*) ,intent(out) :: name
01803 integer(SHR_KIND_IN),optional,intent(out) :: rc
01804
01805
01806
01807
01808 integer(SHR_KIND_IN) :: rCode
01809
01810
01811 character(*),parameter :: subName = '(shr_stream_getModelFieldName) '
01812 character(*),parameter :: F00 = "('(shr_stream_getModelFieldName) ',4a)"
01813
01814
01815
01816
01817
01818 rCode = 0
01819
01820 call shr_string_listGetName(stream%fldListModel,k,name,rCode)
01821
01822 if (present(rc)) rc = rCode
01823
01824 end subroutine shr_stream_getModelFieldName
01825
01826
01827
01828
01829
01830
01831
01832
01833
01834
01835
01836
01837
01838
01839 subroutine shr_stream_getFilepath(strm,path)
01840
01841
01842
01843 type(shr_stream_streamType),intent(in) :: strm
01844 character(*) ,intent(out) :: path
01845
01846
01847
01848
01849
01850
01851
01852 path = strm%filePath
01853
01854 end subroutine shr_stream_getFilePath
01855
01856
01857
01858
01859
01860
01861
01862
01863
01864
01865
01866
01867
01868
01869 subroutine shr_stream_getTInterpAlgo(strm,algo)
01870
01871
01872
01873 type(shr_stream_streamType),intent(in) :: strm
01874 character(*) ,intent(out) :: algo
01875
01876
01877
01878
01879
01880
01881
01882 algo = strm%tInterpAlgo
01883
01884 end subroutine shr_stream_getTInterpAlgo
01885
01886
01887
01888
01889
01890
01891
01892
01893
01894
01895
01896
01897
01898
01899 subroutine shr_stream_getDataSource(strm,str)
01900
01901
01902
01903 type(shr_stream_streamType),intent(in) :: strm
01904 character(*) ,intent(out) :: str
01905
01906
01907
01908
01909
01910
01911
01912 str = strm%dataSource
01913
01914 end subroutine shr_stream_getDataSource
01915
01916
01917
01918
01919
01920
01921
01922
01923
01924
01925
01926
01927
01928
01929 subroutine shr_stream_getDomainInfo(strm,filePath,fileName,timeName,lonName,latName,maskName,areaName)
01930
01931
01932
01933 type(shr_stream_streamType),intent(in) :: strm
01934 character(*) ,intent(out) :: filePath
01935 character(*) ,intent(out) :: fileName
01936 character(*) ,intent(out) :: timeName
01937 character(*) ,intent(out) :: lonName
01938 character(*) ,intent(out) :: latName
01939 character(*) ,intent(out) :: maskName
01940 character(*) ,intent(out) :: areaName
01941
01942
01943
01944
01945
01946
01947
01948 filePath = strm%domFilePath
01949 fileName = strm%domFileName
01950 timeName = strm%domTvarName
01951 lonName = strm%domXvarName
01952 latName = strm%domYvarName
01953 maskName = strm%domMaskName
01954 areaName = strm%domAreaName
01955
01956 end subroutine shr_stream_getDomainInfo
01957
01958
01959
01960
01961
01962
01963
01964
01965
01966
01967
01968
01969
01970
01971
01972
01973 subroutine shr_stream_getFile(filePath,fileName,localFile,rc)
01974
01975 use shr_file_mod, only: shr_file_queryPrefix, shr_file_noPrefix
01976
01977 implicit none
01978
01979
01980
01981 character(*) ,intent(in) :: filePath
01982 character(*) ,intent(inout) :: fileName
01983 character(*) ,optional,intent(out) :: localFile
01984 integer(SHR_KIND_IN),optional,intent(out) :: rc
01985
01986
01987
01988
01989 character(SHR_KIND_CL) :: localFn
01990 integer (SHR_KIND_IN) :: rCode
01991
01992
01993 character(*),parameter :: subName = '(shr_stream_getFile) '
01994 character(*),parameter :: F00 = "('(shr_stream_getFile) ',4a)"
01995
01996
01997
01998
01999
02000
02001
02002
02003
02004
02005
02006
02007
02008
02009
02010
02011
02012 rCode = 0
02013
02014 if ( shr_file_queryPrefix(filePath) /= shr_file_noPrefix ) then
02015 localFn = fileName
02016 call shr_file_get(rCode,localFn, trim(filePath)//fileName)
02017 else
02018 localFn = trim(filePath)//fileName
02019 end if
02020
02021 if (debug>0 .and. s_loglev > 0) then
02022 write(s_logunit,F00) "DEBUG: remote file : ",trim(filePath)//trim(fileName)
02023 write(s_logunit,F00) "DEBUG: local file : ",trim(localFn)
02024 end if
02025
02026 if (.not. present(localFile)) fileName = localFn
02027 if ( present(localFile)) localFile = localFn
02028
02029 if (present(rc)) rc = rCode
02030
02031 end subroutine shr_stream_getFile
02032
02033
02034
02035
02036
02037
02038
02039
02040
02041
02042
02043
02044
02045
02046 subroutine shr_stream_getFirstFileName(strm,file,path)
02047
02048
02049
02050 type(shr_stream_streamType),intent(in) :: strm
02051 character(*) ,intent(out) :: file
02052 character(*),optional ,intent(out) :: path
02053
02054
02055
02056
02057
02058
02059
02060 if (present(path)) path = strm%filePath
02061 file = strm%file(1)%name
02062
02063 end subroutine shr_stream_getFirstFileName
02064
02065
02066
02067
02068
02069
02070
02071
02072
02073
02074
02075
02076
02077
02078 subroutine shr_stream_getNextFileName(strm,fn,fnNext,path,rc)
02079
02080
02081
02082 type(shr_stream_streamType),intent(in) :: strm
02083 character(*) ,intent(in) :: fn
02084 character(*) ,intent(out) :: fnNext
02085 character(*),optional ,intent(out) :: path
02086 integer(SHR_KIND_IN),optional,intent(out) :: rc
02087
02088
02089
02090
02091 integer (SHR_KIND_IN) :: rCode
02092 integer(SHR_KIND_IN) :: n
02093 logical :: found
02094
02095
02096 character(*),parameter :: subName = '(shr_stream_getNextFileName) '
02097 character(*),parameter :: F00 = "('(shr_stream_getNextFileName) ',8a)"
02098
02099
02100
02101
02102
02103
02104 rCode = 0
02105 if (present(path)) path = strm%filePath
02106
02107
02108 found = .false.
02109 do n = 1,strm%nFiles
02110 if ( trim(fn) == trim(strm%file(n)%name)) then
02111 found = .true.
02112 exit
02113 end if
02114 end do
02115 if (.not. found) then
02116 rCode = 1
02117 write(s_logunit,F00) "ERROR: input file name is not in stream: ",trim(fn)
02118 call shr_stream_abort(subName//"ERROR: file name not in stream: "//trim(fn))
02119 end if
02120
02121
02122 n = n+1
02123 if (strm%found_lvd .and. strm%found_gvd) then
02124 if (n > strm%k_gvd) n = strm%k_lvd
02125 else if (strm%found_lvd ) then
02126 if (n > strm%nFiles) n = strm%k_lvd
02127 else if (n > strm%nFiles ) then
02128 n = 1
02129 end if
02130
02131 fnNext = trim(strm%file(n)%name)
02132 if ( present(rc) ) rc = rCode
02133
02134 end subroutine shr_stream_getNextFileName
02135
02136
02137
02138
02139
02140
02141
02142
02143
02144
02145
02146
02147
02148
02149 subroutine shr_stream_getPrevFileName(strm,fn,fnPrev,path,rc)
02150
02151
02152
02153 type(shr_stream_streamType) ,intent(in) :: strm
02154 character(*) ,intent(in) :: fn
02155 character(*) ,intent(out) :: fnPrev
02156 character(*),optional ,intent(out) :: path
02157 integer(SHR_KIND_IN),optional,intent(out) :: rc
02158
02159
02160
02161
02162 integer (SHR_KIND_IN) :: rCode
02163 integer(SHR_KIND_IN) :: n
02164 logical :: found
02165
02166
02167 character(*),parameter :: subName = '(shr_stream_getPrevFileName) '
02168 character(*),parameter :: F00 = "('(shr_stream_getPrevFileName) ',8a)"
02169
02170
02171
02172
02173
02174
02175 rCode = 0
02176 if (present(path)) path = strm%filePath
02177
02178
02179 found = .false.
02180 do n = 1,strm%nFiles
02181 if ( trim(fn) == trim(strm%file(n)%name)) then
02182 found = .true.
02183 exit
02184 end if
02185 end do
02186 if (.not. found) then
02187 rCode = 1
02188 write(s_logunit,F00) "ERROR: input file name is not in stream: ",trim(fn)
02189 call shr_stream_abort(subName//"ERROR: file name not in stream: "//trim(fn))
02190 end if
02191
02192
02193 n = n-1
02194 if (strm%found_lvd .and. strm%found_gvd) then
02195 if ( n < strm%k_lvd) n = strm%k_gvd
02196 end if
02197 if (n>0) then
02198 fnPrev = trim(strm%file(n)%name)
02199 else
02200 fnPrev = "unknown "
02201 end if
02202 if ( present(rc) ) rc = rCode
02203
02204 end subroutine shr_stream_getPrevFileName
02205
02206
02207
02208
02209
02210
02211
02212
02213
02214
02215
02216
02217
02218
02219 subroutine shr_stream_getNFiles(strm,nfiles)
02220
02221
02222
02223 type(shr_stream_streamType),intent(in) :: strm
02224 integer(SHR_KIND_IN) ,intent(out) :: nfiles
02225
02226
02227
02228
02229
02230
02231
02232 nfiles = strm%nfiles
02233
02234 end subroutine shr_stream_getNFiles
02235
02236
02237
02238
02239
02240
02241
02242
02243
02244
02245
02246
02247
02248
02249 subroutine shr_stream_restWrite(strm,fileName,caseName,caseDesc,nstrms,rc)
02250
02251 implicit none
02252
02253
02254
02255 type(shr_stream_streamType) ,intent(in) :: strm(:)
02256 character(*) ,intent(in) :: fileName
02257 character(*) ,intent(in) :: caseName
02258 character(*) ,intent(in) :: caseDesc
02259 integer(SHR_KIND_IN),optional,intent(in) :: nstrms
02260 integer(SHR_KIND_IN),optional,intent(out) :: rc
02261
02262
02263
02264
02265 integer (SHR_KIND_IN) :: rCode
02266 integer(SHR_KIND_IN) :: nStreams
02267 integer(SHR_KIND_IN) :: k,n
02268 character( 8) :: dStr
02269 character(10) :: tStr
02270 character(SHR_KIND_CS) :: str
02271 integer(SHR_KIND_IN) :: nUnit
02272 integer(SHR_KIND_IN) :: nt
02273
02274
02275 character(*),parameter :: subName = '(shr_stream_restWrite) '
02276 character(*),parameter :: F00 = "('(shr_stream_restWrite) ',16a) "
02277 character(*),parameter :: F01 = "('(shr_stream_restWrite) ',a,i5,a,5a) "
02278 character(*),parameter :: F02 = "('(shr_stream_restWrite) ',a,i5,a,5i8) "
02279 character(*),parameter :: F03 = "('(shr_stream_restWrite) ',a,i5,a,5l3) "
02280
02281
02282
02283
02284
02285 rCode = 0
02286
02287 if (present(nstrms)) then
02288 if (size(strm) < nstrms) then
02289 write(s_logunit,F02) "ERROR: nstrms too large for strm",size(strm),nstrms
02290 call shr_stream_abort(subname//": ERROR: nstrms too large for strm")
02291 endif
02292 nStreams = nstrms
02293 else
02294 nStreams = size(strm)
02295 endif
02296 call date_and_time(dStr,tStr)
02297
02298
02299 if (s_loglev > 0) then
02300 write(s_logunit,F00) "case name : ",trim(caseName)
02301 write(s_logunit,F00) "case description : ",trim(caseDesc)
02302 write(s_logunit,F00) "File created : ",dStr(1:4)//'-'//dStr(5:6)//'-'//dStr(7:8)//' ' &
02303 & //tStr(1:2)//':'//tStr(3:4)//':'//tStr(5:6)
02304 write(s_logunit,F01) "Number of streams ",nStreams
02305 endif
02306
02307
02308
02309
02310
02311 nUnit = shr_file_getUnit()
02312 open(nUnit,file=trim(fileName),form="unformatted",action="write")
02313
02314 str = "case name : "//caseName
02315 write(nUnit) str
02316 str = "case description : "//caseDesc
02317 write(nUnit) str
02318 str = 'File created : '//dStr(1:4)//'-'//dStr(5:6)//'-'//dStr(7:8)//' ' &
02319 & //tStr(1:2)//':'//tStr(3:4)//':'//tStr(5:6)
02320 write(nUnit) str
02321
02322 write(nUnit) nStreams
02323 do k = 1,nStreams
02324 if (.not. strm(k)%init) then
02325 rCode = 1
02326 write(s_logunit,F01) "ERROR: can't write uninitialized stream to a restart file, k = ",k
02327 call shr_stream_abort(subName//": ERROR: given uninitialized stream")
02328 end if
02329
02330 write(nUnit) strm(k)%init
02331 write(nUnit) strm(k)%nFiles
02332 write(nUnit) strm(k)%dataSource
02333 write(nUnit) strm(k)%filePath
02334
02335 if (s_loglev > 0) write(s_logunit,F01) "* stream ",k," first file name = ",trim(strm(k)%file(1)%name)
02336 if (s_loglev > 0) write(s_logunit,F03) "* stream ",k," first have data = ",strm(k)%file(1)%haveData
02337 if (s_loglev > 0) write(s_logunit,F02) "* stream ",k," first nt = ",strm(k)%file(1)%nt
02338 nt = strm(k)%file(1)%nt
02339 if (strm(k)%file(1)%haveData) then
02340 if (s_loglev > 0) write(s_logunit,F02) "* stream ",k," first date secs = ", &
02341 strm(k)%file(1)%date(1),strm(k)%file(1)%secs(1)
02342 if (s_loglev > 0) write(s_logunit,F02) "* stream ",k," last date secs = ", &
02343 strm(k)%file(1)%date(nt),strm(k)%file(1)%secs(nt)
02344 endif
02345 do n=1,strm(k)%nFiles
02346 write(nUnit) strm(k)%file(n)%name
02347 write(nUnit) strm(k)%file(n)%haveData
02348 write(nUnit) strm(k)%file(n)%nt
02349 if (strm(k)%file(n)%haveData) then
02350 write(nUnit) strm(k)%file(n)%date(:)
02351 write(nUnit) strm(k)%file(n)%secs(:)
02352 end if
02353 end do
02354
02355 write(nUnit) strm(k)%yearFirst
02356 write(nUnit) strm(k)%yearLast
02357 write(nUnit) strm(k)%yearAlign
02358 write(nUnit) strm(k)%offset
02359
02360
02361 write(nUnit) strm(k)%k_lvd
02362 write(nUnit) strm(k)%n_lvd
02363 write(nUnit) strm(k)%found_lvd
02364 write(nUnit) strm(k)%k_gvd
02365 write(nUnit) strm(k)%n_gvd
02366 write(nUnit) strm(k)%found_gvd
02367
02368 write(nUnit) strm(k)%fldListFile
02369 write(nUnit) strm(k)%fldListModel
02370 write(nUnit) strm(k)%tInterpAlgo
02371 write(nUnit) strm(k)%domFileName
02372 write(nUnit) strm(k)%domFilePath
02373 write(nUnit) strm(k)%domTvarName
02374 write(nUnit) strm(k)%domXvarName
02375 write(nUnit) strm(k)%domYvarName
02376 write(nUnit) strm(k)%domAreaName
02377 write(nUnit) strm(k)%domMaskName
02378
02379 end do
02380
02381 close(nUnit)
02382 call shr_file_freeUnit(nUnit)
02383 if ( present(rc) ) rc = rCode
02384
02385 end subroutine shr_stream_restWrite
02386
02387
02388
02389
02390
02391
02392
02393
02394
02395
02396
02397
02398
02399
02400
02401
02402 subroutine shr_stream_restRead(strm,fileName,nstrms,rc)
02403
02404 implicit none
02405
02406
02407
02408 type(shr_stream_streamType) ,intent(inout) :: strm(:)
02409 character(*) ,intent(in) :: fileName
02410 integer(SHR_KIND_IN),optional,intent(in) :: nstrms
02411 integer(SHR_KIND_IN),optional,intent(out) :: rc
02412
02413
02414
02415
02416 integer (SHR_KIND_IN) :: rCode
02417 integer(SHR_KIND_IN) :: nStreams
02418 integer(SHR_KIND_IN) :: k,n
02419 integer(SHR_KIND_IN) :: nt
02420 character(SHR_KIND_CS) :: str
02421 integer(SHR_KIND_IN) :: nUnit
02422
02423
02424 character(*),parameter :: subName = '(shr_stream_restRead) '
02425 character(*),parameter :: F00 = "('(shr_stream_restRead) ',16a) "
02426 character(*),parameter :: F01 = "('(shr_stream_restRead) ',a,i5,a,5a) "
02427 character(*),parameter :: F02 = "('(shr_stream_restRead) ',a,i5,a,5i8) "
02428 character(*),parameter :: F03 = "('(shr_stream_restRead) ',a,i5,a,5l3) "
02429
02430
02431
02432
02433
02434 rCode = 0
02435
02436
02437
02438
02439 nUnit = shr_file_getUnit()
02440 open(nUnit,file=trim(fileName),form="unformatted",status="old",action="read", iostat=rCode)
02441 if ( rCode /= 0 )then
02442 call shr_file_freeUnit(nUnit)
02443 call shr_stream_abort(subName//": ERROR: error opening file: "//trim(fileName) )
02444 if ( present(rc) ) rc = rCode
02445 return
02446 end if
02447
02448 read(nUnit) str
02449 if (s_loglev > 0) write(s_logunit,F00) trim(str)
02450 read(nUnit) str
02451 if (s_loglev > 0) write(s_logunit,F00) trim(str)
02452 read(nUnit) str
02453 if (s_loglev > 0) write(s_logunit,F00) trim(str)
02454
02455 read(nUnit) nStreams
02456 if (present(nstrms)) then
02457 if (nstrms /= nStreams) then
02458 write(s_logunit,F02) "ERROR: nstrms ne nStreams on restart",nstrms,nStreams
02459 call shr_stream_abort(subname//": ERROR: nstrms ne nStreams on restart")
02460 endif
02461 nStreams = nstrms
02462 endif
02463 if (s_loglev > 0) write(s_logunit,F01) "Number of streams ",nStreams
02464
02465 do k = 1,nStreams
02466 read(nUnit) strm(k)%init
02467 if (.not. strm(k)%init) then
02468 rCode = 1
02469 write(s_logunit,F01) "ERROR: uninitialized stream in restart file, k = ",k
02470 call shr_stream_abort(subName//": ERROR: reading uninitialized stream")
02471 end if
02472
02473 read(nUnit) strm(k)%nFiles
02474 read(nUnit) strm(k)%dataSource
02475 read(nUnit) strm(k)%filePath
02476
02477 do n=1,strm(k)%nFiles
02478 read(nUnit) strm(k)%file(n)%name
02479 read(nUnit) strm(k)%file(n)%haveData
02480 read(nUnit) strm(k)%file(n)%nt
02481 nt = strm(k)%file(n)%nt
02482 if (strm(k)%file(n)%haveData) then
02483 allocate(strm(k)%file(n)%date(nt))
02484 allocate(strm(k)%file(n)%secs(nt))
02485 read(nUnit) strm(k)%file(n)%date(:)
02486 read(nUnit) strm(k)%file(n)%secs(:)
02487 end if
02488 end do
02489
02490 if (s_loglev > 0) write(s_logunit,F01) "* stream ",k," first file name = ",trim(strm(k)%file(1)%name)
02491 if (s_loglev > 0) write(s_logunit,F03) "* stream ",k," first have data = ",strm(k)%file(1)%haveData
02492 if (s_loglev > 0) write(s_logunit,F02) "* stream ",k," first nt = ",strm(k)%file(1)%nt
02493 if (strm(k)%file(1)%haveData) then
02494 nt = strm(k)%file(1)%nt
02495 if (s_loglev > 0) write(s_logunit,F02) "* stream ",k," first date secs = ", &
02496 strm(k)%file(1)%date(1),strm(k)%file(1)%secs(1)
02497 if (s_loglev > 0) write(s_logunit,F02) "* stream ",k," last date secs = ", &
02498 strm(k)%file(1)%date(nt),strm(k)%file(1)%secs(nt)
02499 endif
02500
02501 read(nUnit) strm(k)%yearFirst
02502 read(nUnit) strm(k)%yearLast
02503 read(nUnit) strm(k)%yearAlign
02504 read(nUnit) strm(k)%offset
02505
02506
02507 read(nUnit) strm(k)%k_lvd
02508 read(nUnit) strm(k)%n_lvd
02509 read(nUnit) strm(k)%found_lvd
02510 read(nUnit) strm(k)%k_gvd
02511 read(nUnit) strm(k)%n_gvd
02512 read(nUnit) strm(k)%found_gvd
02513
02514 read(nUnit) strm(k)%fldListFile
02515 read(nUnit) strm(k)%fldListModel
02516 read(nUnit) strm(k)%tInterpAlgo
02517 read(nUnit) strm(k)%domFileName
02518 read(nUnit) strm(k)%domFilePath
02519 read(nUnit) strm(k)%domTvarName
02520 read(nUnit) strm(k)%domXvarName
02521 read(nUnit) strm(k)%domYvarName
02522 read(nUnit) strm(k)%domAreaName
02523 read(nUnit) strm(k)%domMaskName
02524
02525 end do
02526
02527 close(nUnit)
02528 call shr_file_freeUnit(nUnit)
02529 if ( present(rc) ) rc = rCode
02530
02531 end subroutine shr_stream_restRead
02532
02533
02534
02535
02536
02537
02538
02539
02540
02541
02542
02543
02544
02545
02546 subroutine shr_stream_dataDump(strm)
02547
02548
02549
02550 type(shr_stream_streamType),intent(in) :: strm
02551
02552
02553
02554
02555 integer(SHR_KIND_IN) :: k
02556
02557
02558 character(*),parameter :: subName = '(shr_stream_dataDump) '
02559 character(*),parameter :: F00 = "('(shr_stream_dataDump) ',8a)"
02560 character(*),parameter :: F01 = "('(shr_stream_dataDump) ',a,3i5)"
02561 character(*),parameter :: F02 = "('(shr_stream_dataDump) ',a,365i9.8)"
02562 character(*),parameter :: F03 = "('(shr_stream_dataDump) ',a,365i6)"
02563
02564
02565
02566
02567
02568 if (s_loglev <= 0) return
02569
02570 write(s_logunit,F00) "dump internal data for debugging..."
02571
02572
02573
02574
02575 write(s_logunit,F01) "nFiles = ", strm%nFiles
02576 write(s_logunit,F00) "filePath = ", trim(strm%filePath)
02577 do k=1,strm%nFiles
02578 write(s_logunit,F01) "data for file k = ",k
02579 write(s_logunit,F00) "* file(k)%name = ", trim(strm%file(k)%name)
02580 if ( strm%file(k)%haveData ) then
02581 write(s_logunit,F01) "* file(k)%nt = ", strm%file(k)%nt
02582 write(s_logunit,F02) "* file(k)%date(:) = ", strm%file(k)%date(:)
02583 write(s_logunit,F03) "* file(k)%Secs(:) = ", strm%file(k)%secs(:)
02584 else
02585 write(s_logunit,F00) "* time coord data not read in yet for this file"
02586 end if
02587 end do
02588 write(s_logunit,F01) "yearF/L/A = ", strm%yearFirst,strm%yearLast,strm%yearAlign
02589 write(s_logunit,F01) "offset = ", strm%offset
02590 write(s_logunit,F00) "taxMode = ", trim(strm%taxMode)
02591
02592 write(s_logunit,F00) "fldListFile = ", trim(strm%fldListFile)
02593 write(s_logunit,F00) "fldListModel = ", trim(strm%fldListModel)
02594 write(s_logunit,F00) "domFileName = ", trim(strm%domFileName)
02595 write(s_logunit,F00) "domFilePath = ", trim(strm%domFilePath)
02596 write(s_logunit,F00) "domTvarName = ", trim(strm%domTvarName)
02597 write(s_logunit,F00) "domXvarName = ", trim(strm%domXvarName)
02598 write(s_logunit,F00) "domYvarName = ", trim(strm%domYvarName)
02599 write(s_logunit,F00) "domAreaName = ", trim(strm%domAreaName)
02600 write(s_logunit,F00) "domMaskName = ", trim(strm%domMaskName)
02601
02602 end subroutine shr_stream_dataDump
02603
02604
02605
02606
02607
02608
02609
02610
02611
02612
02613
02614
02615
02616
02617
02618
02619 subroutine shr_stream_setDebug(level)
02620
02621 implicit none
02622
02623
02624
02625 integer,intent(in) :: level
02626
02627
02628
02629
02630 character(*),parameter :: subName = '(shr_stream_setDebug) '
02631 character(*),parameter :: F00 = "('(shr_stream_setDebug) ',a) "
02632 character(*),parameter :: F01 = "('(shr_stream_setDebug) ',a,i4) "
02633
02634
02635
02636
02637
02638 debug = level
02639 if (s_loglev > 0) write(s_logunit,F01) "debug level reset to ",level
02640
02641 end subroutine shr_stream_setDebug
02642
02643
02644
02645
02646
02647
02648
02649
02650
02651
02652
02653
02654
02655
02656
02657
02658 subroutine shr_stream_getDebug(level)
02659
02660 implicit none
02661
02662
02663
02664 integer,intent(out) :: level
02665
02666
02667
02668
02669 character(*),parameter :: subName = '(shr_stream_getDebug) '
02670 character(*),parameter :: F00 = "('(shr_stream_getDebug) ',a) "
02671
02672
02673
02674
02675
02676 level = debug
02677
02678 end subroutine shr_stream_getDebug
02679
02680
02681
02682
02683
02684
02685
02686
02687
02688
02689
02690
02691
02692
02693
02694
02695 subroutine shr_stream_setAbort(flag)
02696
02697 implicit none
02698
02699
02700
02701 logical,intent(in) :: flag
02702
02703
02704
02705
02706 character(*),parameter :: subName = '(shr_stream_setAbort) '
02707
02708
02709
02710
02711
02712 doabort = flag
02713
02714 end subroutine shr_stream_setAbort
02715
02716
02717
02718
02719
02720
02721
02722
02723
02724
02725
02726
02727
02728
02729
02730 subroutine shr_stream_abort( msg )
02731
02732 implicit none
02733
02734
02735 character(len=*), optional, intent(IN) :: msg
02736
02737
02738
02739 character(SHR_KIND_CL) :: lmsg
02740
02741
02742 character(*),parameter :: subName = '(shr_stream_abort) '
02743 character(*),parameter :: F00 = "('(shr_stream_abort) ',a) "
02744
02745
02746
02747
02748
02749 lmsg = ' '
02750 if (present(msg)) lmsg= msg
02751
02752 if (doabort) then
02753 call shr_sys_abort(lmsg)
02754 else
02755 write(s_logunit,F00) trim(lmsg)
02756 endif
02757
02758 end subroutine shr_stream_abort
02759
02760
02761
02762
02763
02764
02765
02766
02767
02768
02769
02770
02771
02772
02773
02774
02775
02776 subroutine shr_stream_bcast(stream,comm,rc)
02777
02778 implicit none
02779
02780
02781
02782 type(shr_stream_streamType), intent(inout) :: stream
02783 integer(SHR_KIND_IN), intent(in) :: comm
02784 integer(SHR_KIND_IN),optional,intent(out) :: rc
02785
02786
02787
02788
02789 integer :: n,nt
02790 integer :: pid
02791
02792
02793 character(*),parameter :: subName = '(shr_stream_bcast) '
02794 character(*),parameter :: F00 = "('(shr_stream_bcast) ',a) "
02795
02796
02797
02798
02799
02800 if ( present(rc) ) rc = 0
02801 call shr_mpi_commRank(comm,pid,subName)
02802
02803 call shr_mpi_bcast(stream%init ,comm,subName)
02804 call shr_mpi_bcast(stream%nFiles ,comm,subName)
02805 call shr_mpi_bcast(stream%dataSource ,comm,subName)
02806 call shr_mpi_bcast(stream%filePath ,comm,subName)
02807 call shr_mpi_bcast(stream%yearFirst ,comm,subName)
02808 call shr_mpi_bcast(stream%yearLast ,comm,subName)
02809 call shr_mpi_bcast(stream%yearAlign ,comm,subName)
02810 call shr_mpi_bcast(stream%offset ,comm,subName)
02811 call shr_mpi_bcast(stream%taxMode ,comm,subName)
02812 call shr_mpi_bcast(stream%k_lvd ,comm,subName)
02813 call shr_mpi_bcast(stream%n_lvd ,comm,subName)
02814 call shr_mpi_bcast(stream%found_lvd ,comm,subName)
02815 call shr_mpi_bcast(stream%k_gvd ,comm,subName)
02816 call shr_mpi_bcast(stream%n_gvd ,comm,subName)
02817 call shr_mpi_bcast(stream%found_gvd ,comm,subName)
02818 call shr_mpi_bcast(stream%fldListFile ,comm,subName)
02819 call shr_mpi_bcast(stream%fldListModel,comm,subName)
02820
02821 call shr_mpi_bcast(stream%domFileName ,comm,subName)
02822 call shr_mpi_bcast(stream%domFilePath ,comm,subName)
02823 call shr_mpi_bcast(stream%domTvarName ,comm,subName)
02824 call shr_mpi_bcast(stream%domXvarName ,comm,subName)
02825 call shr_mpi_bcast(stream%domYvarName ,comm,subName)
02826 call shr_mpi_bcast(stream%domMaskName ,comm,subName)
02827
02828 do n = 1,stream%nFiles
02829 call shr_mpi_bcast(stream%file(n)%name ,comm,subName)
02830 call shr_mpi_bcast(stream%file(n)%haveData,comm,subName)
02831 call shr_mpi_bcast(stream%file(n)%nt ,comm,subName)
02832 nt = stream%file(n)%nt
02833 if (pid /= 0) allocate(stream%file(n)%date(nt),stream%file(n)%secs(nt))
02834 call shr_mpi_bcast(stream%file(n)%date ,comm,subName)
02835 call shr_mpi_bcast(stream%file(n)%secs ,comm,subName)
02836 enddo
02837
02838 end subroutine shr_stream_bcast
02839
02840
02841 end module shr_stream_mod
02842
02843