00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017 module shr_strdata_mod
00018
00019 use shr_const_mod, only: SHR_CONST_PI
00020 use shr_kind_mod, only: IN=>SHR_KIND_IN, R8=>SHR_KIND_R8, &
00021 CS=>SHR_KIND_CS, CL=>SHR_KIND_CL, &
00022 CX=>SHR_KIND_CX
00023 use shr_sys_mod, only : shr_sys_abort, shr_sys_flush
00024 use shr_mpi_mod, only : shr_mpi_bcast
00025 use shr_file_mod, only : shr_file_getunit, shr_file_freeunit
00026 use shr_log_mod, only: loglev => shr_log_Level
00027 use shr_log_mod, only: logunit => shr_log_Unit
00028 use shr_stream_mod
00029 use shr_string_mod
00030 use shr_map_mod
00031 use shr_cal_mod, only: shr_cal_date2eday
00032 use shr_orb_mod, only: shr_orb_decl, shr_orb_cosz, shr_orb_undef_real
00033 use shr_tinterp_mod
00034
00035 use shr_dmodel_mod
00036
00037 use shr_mct_mod
00038 use mct_mod
00039 use perf_mod
00040 use pio
00041
00042 implicit none
00043
00044 private
00045
00046
00047
00048 public shr_strdata_type
00049
00050
00051
00052 public shr_strdata_readnml
00053 public shr_strdata_bcastnml
00054 public shr_strdata_restRead
00055 public shr_strdata_restWrite
00056 public shr_strdata_setOrbs
00057 public shr_strdata_print
00058 public shr_strdata_init
00059 public shr_strdata_create
00060 public shr_strdata_advance
00061 public shr_strdata_clean
00062 public shr_strdata_setlogunit
00063
00064
00065
00066
00067
00068
00069 integer(IN),parameter :: nStrMax = 30
00070 integer(IN),parameter :: nVecMax = 30
00071 character(len=*),public,parameter :: shr_strdata_nullstr = 'null'
00072 character(len=*),parameter :: shr_strdata_unset = 'NOT_SET'
00073 real(R8),parameter,private :: dtlimit_default = 1.5_R8
00074
00075 type shr_strdata_type
00076
00077 character(CL) :: dataMode
00078 character(CL) :: domainFile
00079 character(CL) :: streams (nStrMax)
00080 character(CL) :: taxMode (nStrMax) = trim(shr_stream_taxis_cycle)
00081 real(R8) :: dtlimit (nStrMax) = dtlimit_default
00082 character(CL) :: vectors (nVecMax)
00083 character(CL) :: fillalgo(nStrMax)
00084 character(CL) :: fillmask(nStrMax)
00085 character(CL) :: fillread(nStrMax)
00086 character(CL) :: fillwrit(nStrMax)
00087 character(CL) :: mapalgo (nStrMax)
00088 character(CL) :: mapmask (nStrMax)
00089 character(CL) :: mapread(nStrMax)
00090 character(CL) :: mapwrit(nStrMax)
00091 character(CL) :: tintalgo(nStrMax)
00092 integer(IN) :: io_type
00093 integer(IN) :: num_iotasks
00094 integer(IN) :: io_root
00095 integer(IN) :: io_stride
00096 integer(IN) :: num_agg
00097
00098 real(R8) :: eccen
00099 real(R8) :: mvelpp
00100 real(R8) :: lambm0
00101 real(R8) :: obliqr
00102 integer(IN) :: modeldt
00103
00104 integer(IN) :: nxg = 0
00105 integer(IN) :: nyg = 0
00106 integer(IN) :: lsize
00107 type(mct_gsmap) :: gsmap
00108 type(mct_ggrid) :: grid
00109 type(mct_avect) :: avs(nStrMax)
00110
00111 type(shr_stream_streamType) :: stream(nStrMax)
00112 type(iosystem_desc_t) :: pio_subsystem
00113 type(io_desc_t) :: pio_iodesc(nStrMax)
00114 integer(IN) :: nstreams = 0
00115 integer(IN) :: strnxg(nStrMax)
00116 integer(IN) :: strnyg(nStrMax)
00117 logical :: dofill(nStrMax)
00118 logical :: domaps(nStrMax)
00119 integer(IN) :: lsizeR(nStrMax)
00120 type(mct_gsmap) :: gsmapR(nStrMax)
00121 type(mct_rearr) :: rearrR(nStrMax)
00122 type(mct_ggrid) :: gridR(nStrMax)
00123 type(mct_avect) :: avRLB(nStrMax)
00124 type(mct_avect) :: avRUB(nStrMax)
00125 type(mct_avect) :: avFUB(nStrMax)
00126 type(mct_avect) :: avFLB(nStrMax)
00127 type(mct_avect) :: avCoszen(nStrMax)
00128 type(mct_sMatP) :: sMatPf(nStrMax)
00129 type(mct_sMatP) :: sMatPs(nStrMax)
00130 integer(IN) :: ymdLB(nStrMax),todLB(nStrMax)
00131 integer(IN) :: ymdUB(nStrMax),todUB(nStrMax)
00132 real(R8) :: dtmin(nStrMax) = 1.0e30
00133 real(R8) :: dtmax(nStrMax) = 0.0
00134 integer(IN) :: ymd ,tod
00135 integer(IN) :: nvectors
00136 integer(IN) :: ustrm (nVecMax)
00137 integer(IN) :: vstrm (nVecMax)
00138 character(CL) :: allocstring
00139 end type shr_strdata_type
00140
00141 real(R8),parameter,private :: deg2rad = SHR_CONST_PI/180.0_R8
00142 character(len=*),parameter :: allocstring_value = 'strdata_allocated'
00143
00144
00145 contains
00146
00147
00148 subroutine shr_strdata_init(SDAT,mpicom,compid,name,scmmode,scmlon,scmlat, &
00149 gsmap,ggrid,nxg,nyg)
00150
00151 implicit none
00152
00153 type(shr_strdata_type),intent(inout) :: SDAT
00154 integer(IN) ,intent(in) :: mpicom
00155 integer(IN) ,intent(in) :: compid
00156 character(len=*) ,intent(in),optional :: name
00157 logical ,intent(in),optional :: scmmode
00158 real(R8) ,intent(in),optional :: scmlon
00159 real(R8) ,intent(in),optional :: scmlat
00160 type(mct_gsmap) ,intent(in),optional :: gsmap
00161 type(mct_ggrid) ,intent(in),optional :: ggrid
00162 integer(IN) ,intent(in),optional :: nxg
00163 integer(IN) ,intent(in),optional :: nyg
00164
00165 integer(IN) :: n,m,k
00166 integer(IN) :: nu,nv
00167 integer(IN) :: my_task,npes
00168 integer(IN),parameter :: master_task = 0
00169 character(CS) :: lname
00170 character(CL) :: filePath
00171 character(CL) :: fileName
00172 character(CS) :: timeName
00173 character(CS) :: lonName
00174 character(CS) :: latName
00175 character(CS) :: maskName
00176 character(CS) :: areaName
00177 character(CS) :: uname
00178 character(CS) :: vname
00179 character(CX) :: fldList
00180 integer(IN) :: lsize
00181 integer(IN) :: nfiles
00182 integer(IN) :: ierr
00183 integer(IN) :: method
00184 integer(IN), pointer :: dof(:)
00185 type(mct_sMat):: sMati
00186 logical :: lscmmode
00187
00188 character(len=*),parameter :: subname = "(shr_strdata_init) "
00189 character(*),parameter :: F00 = "('(shr_strdata_init) ',8a)"
00190
00191
00192
00193
00194
00195 call MPI_COMM_SIZE(mpicom,npes,ierr)
00196 call MPI_COMM_RANK(mpicom,my_task,ierr)
00197
00198
00199 if (my_task == master_task) then
00200 do n=1,nStrMax
00201
00202 if (trim(SDAT%streams(n)) /= trim(shr_strdata_nullstr)) SDAT%nstreams = max(SDAT%nstreams,n)
00203
00204 call shr_stream_getNFiles(SDAT%stream(n),nfiles)
00205 if (nfiles > 0) SDAT%nstreams = max(SDAT%nstreams,n)
00206
00207 if (trim(SDAT%taxMode(n)) == trim(shr_stream_taxis_extend)) SDAT%dtlimit(n) = 1.0e30
00208 end do
00209 SDAT%nvectors = 0
00210 do n=1,nVecMax
00211 if (trim(SDAT%vectors(n)) /= trim(shr_strdata_nullstr)) SDAT%nvectors = n
00212 end do
00213 endif
00214
00215 call shr_mpi_bcast(SDAT%nstreams ,mpicom,'nstreams')
00216 call shr_mpi_bcast(SDAT%nvectors ,mpicom,'nvectors')
00217 call shr_mpi_bcast(SDAT%dtlimit ,mpicom,'dtlimit')
00218
00219 if (SDAT%nstreams < 1) return
00220
00221 n = 0
00222 if (present(gsmap)) then
00223 n = n + 1
00224 endif
00225 if (present(ggrid)) then
00226 n = n + 1
00227 endif
00228 if (present(nxg)) then
00229 n = n + 1
00230 endif
00231 if (present(nyg)) then
00232 n = n + 1
00233 endif
00234
00235 if ( n == 0 .or. n == 4) then
00236
00237 else
00238 write(logunit,*) subname,' ERROR: gsmap, ggrid, nxg, and nyg must be specified together'
00239 call shr_sys_abort(subname//' ERROR: gsmap, ggrid, nxg, nyg set together')
00240 endif
00241
00242 lscmmode = .false.
00243 if (present(scmmode)) then
00244 lscmmode = scmmode
00245 if (lscmmode) then
00246 if (.not.present(scmlon) .or. .not.present(scmlat)) then
00247 write(logunit,*) subname,' ERROR: scmmode requires scmlon and scmlat'
00248 call shr_sys_abort(subname//' ERROR: scmmode1 lon lat')
00249 endif
00250 endif
00251 endif
00252
00253 lname = ""
00254 if (present(name)) then
00255 lname = "_"//trim(name)
00256 endif
00257
00258
00259
00260 call pio_init(my_task, mpicom, SDAT%num_iotasks, SDAT%num_agg, SDAT%io_stride, &
00261 PIO_REARR_BOX, SDAT%pio_subsystem, base=SDAT%io_root)
00262
00263 do n = 1,SDAT%nstreams
00264 if (my_task == master_task) then
00265 call shr_stream_getDomainInfo(SDAT%stream(n),filePath,fileName,timeName,lonName, &
00266 latName,maskName,areaName)
00267 call shr_stream_getFile(filePath,fileName)
00268 endif
00269 call shr_mpi_bcast(fileName,mpicom)
00270 call shr_mpi_bcast(lonName,mpicom)
00271 call shr_mpi_bcast(latName,mpicom)
00272 call shr_mpi_bcast(maskName,mpicom)
00273 call shr_mpi_bcast(areaName,mpicom)
00274 call shr_dmodel_readgrid(SDAT%gridR(n),SDAT%gsmapR(n),SDAT%strnxg(n),SDAT%strnyg(n), &
00275 fileName, compid, mpicom, '1d', lonName, latName, maskName, areaName)
00276 SDAT%lsizeR(n) = mct_gsmap_lsize(SDAT%gsmapR(n),mpicom)
00277
00278 call mct_gsmap_OrderedPoints(SDAT%gsmapR(n), my_task, dof)
00279 call pio_initdecomp(SDAT%pio_subsystem, pio_double, &
00280 (/SDAT%strnxg(n),SDAT%strnyg(n)/), dof, SDAT%pio_iodesc(n))
00281 deallocate(dof)
00282 enddo
00283
00284
00285
00286 if (present(gsmap)) then
00287 SDAT%nxg = nxg
00288 SDAT%nyg = nyg
00289 lsize = mct_gsmap_lsize(gsmap,mpicom)
00290 call mct_gsmap_Copy(gsmap,SDAT%gsmap)
00291 call mct_ggrid_init(SDAT%grid, ggrid, lsize)
00292 call mct_aVect_copy(ggrid%data, SDAT%grid%data)
00293 else
00294 if (trim(SDAT%domainfile) == trim(shr_strdata_nullstr)) then
00295 if (SDAT%nstreams > 0) then
00296 if (my_task == master_task) then
00297 call shr_stream_getDomainInfo(SDAT%stream(1),filePath,fileName,timeName,lonName, &
00298 latName,maskName,areaName)
00299 call shr_stream_getFile(filePath,fileName)
00300 endif
00301 call shr_mpi_bcast(fileName,mpicom)
00302 call shr_mpi_bcast(lonName,mpicom)
00303 call shr_mpi_bcast(latName,mpicom)
00304 call shr_mpi_bcast(maskName,mpicom)
00305 call shr_mpi_bcast(areaName,mpicom)
00306 if (lscmmode) then
00307 call shr_dmodel_readgrid(SDAT%grid,SDAT%gsmap,SDAT%nxg,SDAT%nyg, &
00308 fileName, compid, mpicom, '1d', lonName, latName, maskName, areaName, &
00309 scmmode=lscmmode,scmlon=scmlon,scmlat=scmlat)
00310 else
00311 call shr_dmodel_readgrid(SDAT%grid,SDAT%gsmap,SDAT%nxg,SDAT%nyg, &
00312 fileName, compid, mpicom, '1d', lonName, latName, maskName, areaName)
00313 endif
00314 endif
00315 else
00316 if (lscmmode) then
00317 call shr_dmodel_readgrid(SDAT%grid,SDAT%gsmap,SDAT%nxg,SDAT%nyg, &
00318 SDAT%domainfile, compid, mpicom, '1d', readfrac=.true., &
00319 scmmode=lscmmode,scmlon=scmlon,scmlat=scmlat)
00320 else
00321 call shr_dmodel_readgrid(SDAT%grid,SDAT%gsmap,SDAT%nxg,SDAT%nyg, &
00322 SDAT%domainfile, compid, mpicom, '1d', readfrac=.true.)
00323 endif
00324 endif
00325 endif
00326 SDAT%lsize = mct_gsmap_lsize(SDAT%gsmap,mpicom)
00327
00328
00329
00330 do n = 1,SDAT%nstreams
00331 if (shr_dmodel_gGridCompare(SDAT%gridR(n),SDAT%gsmapR(n),SDAT%grid,SDAT%gsmap, &
00332 shr_dmodel_gGridCompareMaskSubset,mpicom) .or. trim(SDAT%fillalgo(n))=='none') then
00333 SDAT%dofill(n) = .false.
00334 else
00335 SDAT%dofill(n) = .true.
00336 endif
00337 if (trim(SDAT%mapmask(n)) == 'dstmask') then
00338 method = shr_dmodel_gGridCompareXYabsMask
00339 else
00340 method = shr_dmodel_gGridCompareXYabs
00341 endif
00342 if (shr_dmodel_gGridCompare(SDAT%gridR(n),SDAT%gsmapR(n),SDAT%grid,SDAT%gsmap, &
00343 method,mpicom,0.01_r8) .or. trim(SDAT%mapalgo(n))=='none') then
00344 SDAT%domaps(n) = .false.
00345 else
00346 SDAT%domaps(n) = .true.
00347 endif
00348
00349 if (SDAT%dofill(n)) then
00350 if (trim(SDAT%fillread(n)) == trim(shr_strdata_unset)) then
00351 if (my_task == master_task) write(logunit,F00) ' calling shr_dmodel_mapSet for fill'
00352 call shr_dmodel_mapSet(SDAT%sMatPf(n), &
00353 SDAT%gridR(n),SDAT%gsmapR(n),SDAT%strnxg(n),SDAT%strnyg(n), &
00354 SDAT%gridR(n),SDAT%gsmapR(n),SDAT%strnxg(n),SDAT%strnyg(n), &
00355 name='mapFill', type='cfill', &
00356 algo=trim(SDAT%fillalgo(n)),mask=trim(SDAT%fillmask(n)),vect='scalar', &
00357 compid=compid,mpicom=mpicom)
00358 if (trim(SDAT%fillwrit(n)) /= trim(shr_strdata_unset)) then
00359 if (my_task == master_task) write(logunit,F00) ' writing ',trim(SDAT%fillwrit(n))
00360 call shr_mct_sMatWritednc(SDAT%sMatPf(n)%Matrix,SDAT%fillwrit(n),compid,mpicom)
00361 endif
00362 else
00363 if (my_task == master_task) write(logunit,F00) ' reading ',trim(SDAT%fillread(n))
00364 call shr_mct_sMatReaddnc(sMati,SDAT%gsmapR(n),SDAT%gsmapR(n),'src', &
00365 filename=trim(SDAT%fillread(n)),mytask=my_task,mpicom=mpicom)
00366 call mct_sMatP_Init(SDAT%sMatPf(n),sMati,SDAT%gsMapR(n),SDAT%gsmapR(n),0, mpicom, compid)
00367 call mct_sMat_Clean(sMati)
00368 endif
00369 endif
00370 if (SDAT%domaps(n)) then
00371 if (trim(SDAT%mapread(n)) == trim(shr_strdata_unset)) then
00372 if (my_task == master_task) write(logunit,F00) ' calling shr_dmodel_mapSet for remap'
00373 call shr_dmodel_mapSet(SDAT%sMatPs(n), &
00374 SDAT%gridR(n),SDAT%gsmapR(n),SDAT%strnxg(n),SDAT%strnyg(n), &
00375 SDAT%grid ,SDAT%gsmap ,SDAT%nxg ,SDAT%nyg, &
00376 name='mapScalar', type='remap', &
00377 algo=trim(SDAT%mapalgo(n)),mask=trim(SDAT%mapmask(n)), vect='scalar', &
00378 compid=compid,mpicom=mpicom)
00379 if (trim(SDAT%mapwrit(n)) /= trim(shr_strdata_unset)) then
00380 if (my_task == master_task) write(logunit,F00) ' writing ',trim(SDAT%mapwrit(n))
00381 call shr_mct_sMatWritednc(SDAT%sMatPs(n)%Matrix,SDAT%mapwrit(n),compid,mpicom)
00382 endif
00383 else
00384 if (my_task == master_task) write(logunit,F00) ' reading ',trim(SDAT%mapread(n))
00385 call shr_mct_sMatReaddnc(sMati,SDAT%gsmapR(n),SDAT%gsmap,'src', &
00386 filename=trim(SDAT%mapread(n)),mytask=my_task,mpicom=mpicom)
00387 call mct_sMatP_Init(SDAT%sMatPs(n),sMati,SDAT%gsMapR(n),SDAT%gsmap,0, mpicom, compid)
00388 call mct_sMat_Clean(sMati)
00389 endif
00390 else
00391 call mct_rearr_init(SDAT%gsmapR(n), SDAT%gsmap, mpicom, SDAT%rearrR(n))
00392 endif
00393 enddo
00394
00395
00396
00397 do n = 1,SDAT%nstreams
00398 if (my_task == master_task) then
00399 call shr_stream_getModelFieldList(SDAT%stream(n),fldList)
00400 endif
00401 call shr_mpi_bcast(fldList,mpicom)
00402 call mct_aVect_init(SDAT%avs(n) ,rlist=fldList,lsize=SDAT%lsize)
00403 call mct_aVect_init(SDAT%avFLB(n),rlist=fldList,lsize=SDAT%lsize)
00404 call mct_aVect_init(SDAT%avFUB(n),rlist=fldList,lsize=SDAT%lsize)
00405 call mct_aVect_init(SDAT%avRLB(n),rlist=fldList,lsize=SDAT%lsizeR(n))
00406 call mct_aVect_init(SDAT%avRUB(n),rlist=fldList,lsize=SDAT%lsizeR(n))
00407 if (trim(SDAT%tintalgo(n)) == 'coszen') then
00408 call mct_aVect_init(SDAT%avCoszen(n),rlist="tavCosz",lsize=SDAT%lsize)
00409 endif
00410 enddo
00411
00412
00413
00414 do m = 1,SDAT%nvectors
00415 if (.not. shr_string_listIsValid(SDAT%vectors(m))) then
00416 write(logunit,*) trim(subname),' vec fldlist invalid m=',m,trim(SDAT%vectors(m))
00417 call shr_sys_abort(subname//': vec fldlist invalid:'//trim(SDAT%vectors(m)))
00418 endif
00419 if (shr_string_listGetNum(SDAT%vectors(m)) /= 2) then
00420 write(logunit,*) trim(subname),' vec fldlist ne 2 m=',m,trim(SDAT%vectors(m))
00421 call shr_sys_abort(subname//': vec fldlist ne 2:'//trim(SDAT%vectors(m)))
00422 endif
00423 call shr_string_listGetName(SDAT%vectors(m),1,uname)
00424 call shr_string_listGetName(SDAT%vectors(m),2,vname)
00425 nu = 0
00426 nv = 0
00427 do n = 1,SDAT%nstreams
00428 k = mct_aVect_indexRA(SDAT%avRLB(n),trim(uname),perrWith='quiet')
00429 if (k > 0) nu = n
00430 k = mct_aVect_indexRA(SDAT%avRLB(n),trim(vname),perrWith='quiet')
00431 if (k > 0) nv = n
00432 enddo
00433 if (nu == 0 .or. nv == 0) then
00434 write(logunit,*) trim(subname),' vec flds not found m=',m,trim(SDAT%vectors(m))
00435 call shr_sys_abort(subname//': vec flds not found:'//trim(SDAT%vectors(m)))
00436 endif
00437 if (nu /= nv) then
00438 if ((.not. shr_dmodel_gGridCompare(SDAT%gridR(nu),SDAT%gsmapR(nu), &
00439 SDAT%gridR(nv),SDAT%gsmapR(nv), &
00440 shr_dmodel_gGridCompareXYabs,mpicom,0.01_r8)) .or. &
00441 (.not. shr_dmodel_gGridCompare(SDAT%gridR(nu),SDAT%gsmapR(nu), &
00442 SDAT%gridR(nv),SDAT%gsmapR(nv), &
00443 shr_dmodel_gGridCompareMaskZeros,mpicom))) then
00444 write(logunit,*) trim(subname),' vec fld doms not same m=',m,trim(SDAT%vectors(m))
00445 call shr_sys_abort(subname//': vec fld doms not same:'//trim(SDAT%vectors(m)))
00446 endif
00447 endif
00448 SDAT%ustrm(m) = nu
00449 SDAT%vstrm(m) = nv
00450 enddo
00451
00452 end subroutine shr_strdata_init
00453
00454
00455
00456 subroutine shr_strdata_advance(SDAT,ymd,tod,mpicom,istr,timers)
00457
00458 implicit none
00459
00460 type(shr_strdata_type),intent(inout) :: SDAT
00461 integer(IN),intent(in) :: ymd
00462 integer(IN),intent(in) :: tod
00463 integer(IN),intent(in) :: mpicom
00464 character(len=*),intent(in),optional :: istr
00465 logical ,intent(in),optional :: timers
00466
00467 integer(IN) :: n,m,i,k,l,kf
00468 integer(IN) :: my_task,npes
00469 integer(IN),parameter :: master_task = 0
00470 logical :: mssrmlf
00471 logical,allocatable :: newData(:)
00472 integer(IN) :: ierr
00473 integer(IN) :: nu,nv
00474 integer(IN) :: lsize,lsizeR,lsizeF
00475 integer(IN) :: yy,mm,dd
00476 type(mct_avect) :: avRtmp
00477 type(mct_avect) :: avRV,avFV
00478 character(len=32) :: lstr
00479 logical :: ltimers
00480 real(R8) :: flb,fub
00481
00482
00483 real(R8) :: calday
00484 real(R8) :: declin
00485 real(R8) :: eccf
00486 real(R8),pointer :: lonr(:)
00487 real(R8),pointer :: latr(:)
00488 real(R8),pointer :: cosz(:)
00489 real(R8),pointer :: tavCosz(:)
00490 real(R8),pointer :: xlon(:),ylon(:)
00491 real(R8),parameter :: solZenMin = 0.001_R8
00492
00493 integer(IN) :: edayLB,edayUB
00494 real(R8) :: dtime
00495 integer(IN) :: uvar,vvar
00496 logical :: someNewData
00497 character(CS) :: uname
00498 character(CS) :: vname
00499 character(len=*),parameter :: timname = "_strd_adv"
00500 integer(IN),parameter :: tadj = 2
00501
00502
00503 character(*),parameter :: subname = "(shr_strdata_advance) "
00504
00505
00506
00507
00508
00509 if (SDAT%nstreams < 1) return
00510
00511 lstr = ''
00512 if (present(istr)) then
00513 lstr = trim(istr)
00514 endif
00515
00516 ltimers = .true.
00517 if (present(timers)) then
00518 ltimers = timers
00519 endif
00520
00521 if (.not.ltimers) call t_adj_detailf(tadj)
00522
00523 call MPI_COMM_SIZE(mpicom,npes,ierr)
00524 call MPI_COMM_RANK(mpicom,my_task,ierr)
00525
00526 mssrmlf = .false.
00527
00528 SDAT%ymd = ymd
00529 SDAT%tod = tod
00530
00531 if (SDAT%nstreams > 0) then
00532 allocate(newData(SDAT%nstreams))
00533
00534 do n = 1,SDAT%nstreams
00535 call t_barrierf(trim(lstr)//trim(timname)//'_readLBUB_BARRIER',mpicom)
00536 call t_startf(trim(lstr)//trim(timname)//'_readLBUB')
00537 call shr_dmodel_readLBUB(SDAT%stream(n),SDAT%pio_subsystem,SDAT%io_type,SDAT%pio_iodesc(n), &
00538 ymd,tod,mpicom,SDAT%gsmapR(n),&
00539 SDAT%avRLB(n),SDAT%ymdLB(n),SDAT%todLB(n), &
00540 SDAT%avRUB(n),SDAT%ymdUB(n),SDAT%todUB(n),newData(n), &
00541 istr=trim(lstr)//'_readLBUB')
00542 if (newData(n)) then
00543 call shr_cal_date2eday(SDAT%ymdLB(n),edayLB)
00544 call shr_cal_date2eday(SDAT%ymdUB(n),edayUB)
00545 dtime = abs(real(edayUB-edayLB,R8) + real(SDAT%todUB(n)-SDAT%todLB(n),R8)/shr_const_cDay)
00546 SDAT%dtmin(n) = min(SDAT%dtmin(n),dtime)
00547 SDAT%dtmax(n) = max(SDAT%dtmax(n),dtime)
00548 if ((SDAT%dtmax(n)/SDAT%dtmin(n)) > SDAT%dtlimit(n)) then
00549 write(logunit,*) trim(subName),' ERROR: dt limit1 ',SDAT%dtmax(n),SDAT%dtmin(n),SDAT%dtlimit(n)
00550 write(logunit,*) trim(subName),' ERROR: dt limit2 ',SDAT%ymdLB(n),SDAT%todLB(n), &
00551 SDAT%ymdUB(n),SDAT%todUB(n)
00552 call shr_sys_abort(trim(subName)//' ERROR dt limit')
00553 endif
00554 endif
00555 call t_stopf(trim(lstr)//trim(timname)//'_readLBUB')
00556 enddo
00557
00558 do n = 1,SDAT%nstreams
00559 if (newData(n)) then
00560
00561 if (SDAT%doFill(n)) then
00562 call t_startf(trim(lstr)//trim(timname)//'_fill')
00563 lsize = mct_aVect_lsize(SDAT%avRLB(n))
00564 call mct_aVect_init(avRtmp,SDAT%avRLB(n),lsize)
00565 call mct_aVect_copy(SDAT%avRLB(n),avRtmp)
00566 call mct_sMat_avMult(avRtmp,SDAT%sMatPf(n),SDAT%avRLB(n))
00567 call mct_aVect_copy(SDAT%avRUB(n),avRtmp)
00568 call mct_sMat_avMult(avRtmp,SDAT%sMatPf(n),SDAT%avRUB(n))
00569 call mct_aVect_clean(avRtmp)
00570 call t_stopf(trim(lstr)//trim(timname)//'_fill')
00571 endif
00572
00573 if (SDAT%domaps(n)) then
00574 call t_startf(trim(lstr)//trim(timname)//'_map')
00575 call mct_sMat_avMult(SDAT%avRLB(n),SDAT%sMatPs(n),SDAT%avFLB(n))
00576 call mct_sMat_avMult(SDAT%avRUB(n),SDAT%sMatPs(n),SDAT%avFUB(n))
00577 call t_stopf(trim(lstr)//trim(timname)//'_map')
00578 else
00579 call t_startf(trim(lstr)//trim(timname)//'_rearr')
00580 call mct_rearr_rearrange(SDAT%avRLB(n),SDAT%avFLB(n),SDAT%rearrR(n))
00581 call mct_rearr_rearrange(SDAT%avRUB(n),SDAT%avFUB(n),SDAT%rearrR(n))
00582 call t_stopf(trim(lstr)//trim(timname)//'_rearr')
00583 endif
00584
00585 endif
00586 enddo
00587
00588 do m = 1,SDAT%nvectors
00589 nu = SDAT%ustrm(m)
00590 nv = SDAT%vstrm(m)
00591 if ((SDAT%domaps(nu) .or. SDAT%domaps(nv)) .and. &
00592 (newdata(nu) .or. newdata(nv))) then
00593
00594 call t_startf(trim(lstr)//trim(timname)//'_vect')
00595 call shr_string_listGetName(SDAT%vectors(m),1,uname)
00596 call shr_string_listGetName(SDAT%vectors(m),2,vname)
00597 lsizeR = mct_aVect_lsize(SDAT%avRLB(nu))
00598 lsizeF = mct_aVect_lsize(SDAT%avFLB(nu))
00599 call mct_aVect_init(avRV,rlist=SDAT%vectors(m),lsize=lsizeR)
00600 call mct_aVect_init(avFV,rlist=SDAT%vectors(m),lsize=lsizeF)
00601 call mct_aVect_exportRattr(SDAT%gridR(nu)%data,'lon',xlon)
00602 call mct_aVect_exportRattr(SDAT%grid %data,'lon',ylon)
00603 xlon = xlon * deg2rad
00604 ylon = ylon * deg2rad
00605
00606
00607
00608 uvar = mct_aVect_indexRA(SDAT%avRLB(nu),trim(uname))
00609 vvar = mct_aVect_indexRA(SDAT%avRLB(nv),trim(vname))
00610 do i = 1,lsizeR
00611 avRV%rAttr(1,i) = SDAT%avRLB(nu)%rAttr(uvar,i) * cos(xlon(i)) &
00612 -SDAT%avRLB(nv)%rAttr(vvar,i) * sin(xlon(i))
00613 avRV%rAttr(2,i) = SDAT%avRLB(nu)%rAttr(uvar,i) * sin(xlon(i)) &
00614 +SDAT%avRLB(nv)%rAttr(vvar,i) * cos(xlon(i))
00615 enddo
00616 call mct_sMat_avMult(avRV,SDAT%sMatPs(nu),avFV)
00617
00618
00619
00620 do i = 1,lsizeF
00621 SDAT%avFLB(nu)%rAttr(uvar,i) = avFV%rAttr(1,i) * cos(ylon(i)) &
00622 +avFV%rAttr(2,i) * sin(ylon(i))
00623 SDAT%avFLB(nv)%rAttr(vvar,i) = -avFV%rAttr(1,i) * sin(ylon(i)) &
00624 +avFV%rAttr(2,i) * cos(ylon(i))
00625 enddo
00626
00627
00628
00629 uvar = mct_aVect_indexRA(SDAT%avRUB(nu),trim(uname))
00630 vvar = mct_aVect_indexRA(SDAT%avRUB(nv),trim(vname))
00631 do i = 1,lsizeR
00632 avRV%rAttr(1,i) = SDAT%avRUB(nu)%rAttr(uvar,i) * cos(xlon(i)) &
00633 -SDAT%avRUB(nv)%rAttr(vvar,i) * sin(xlon(i))
00634 avRV%rAttr(2,i) = SDAT%avRUB(nu)%rAttr(uvar,i) * sin(xlon(i)) &
00635 +SDAT%avRUB(nv)%rAttr(vvar,i) * cos(xlon(i))
00636 enddo
00637 call mct_sMat_avMult(avRV,SDAT%sMatPs(nu),avFV)
00638
00639
00640
00641 do i = 1,lsizeF
00642 SDAT%avFUB(nu)%rAttr(uvar,i) = avFV%rAttr(1,i) * cos(ylon(i)) &
00643 +avFV%rAttr(2,i) * sin(ylon(i))
00644 SDAT%avFUB(nv)%rAttr(vvar,i) = -avFV%rAttr(1,i) * sin(ylon(i)) &
00645 +avFV%rAttr(2,i) * cos(ylon(i))
00646 enddo
00647
00648 call mct_aVect_clean(avRV)
00649 call mct_aVect_clean(avFV)
00650 deallocate(xlon,ylon)
00651
00652 call t_stopf(trim(lstr)//trim(timname)//'_vect')
00653 endif
00654 enddo
00655
00656 do n = 1,SDAT%nstreams
00657
00658
00659 if (trim(SDAT%tintalgo(n)) == 'coszen') then
00660 call t_startf(trim(lstr)//trim(timname)//'_coszen')
00661
00662
00663 if (SDAT%eccen == SHR_ORB_UNDEF_REAL) then
00664 call shr_sys_abort(subname//' ERROR in orb params for coszen tinterp')
00665 else if (SDAT%modeldt < 1) then
00666 call shr_sys_abort(subname//' ERROR: model dt < 1 for coszen tinterp')
00667 endif
00668
00669
00670 lsizeF = mct_aVect_lsize(SDAT%avFLB(n))
00671 allocate(tavCosz(lsizeF),cosz(lsizeF),lonr(lsizeF),latr(lsizeF))
00672
00673
00674 kf = mct_aVect_indexRA(SDAT%grid%data,'lat')
00675 latr(:) = SDAT%grid%data%rAttr(kf,:) * deg2rad
00676 kf = mct_aVect_indexRA(SDAT%grid%data,'lon')
00677 lonr(:) = SDAT%grid%data%rAttr(kf,:) * deg2rad
00678
00679 call t_startf(trim(lstr)//trim(timname)//'_coszenC')
00680 cosz = 0.0_r8
00681 call shr_tInterp_getCosz(cosz,lonr,latr,ymd,tod, &
00682 SDAT%eccen,SDAT%mvelpp,SDAT%lambm0,SDAT%obliqr)
00683 call t_stopf(trim(lstr)//trim(timname)//'_coszenC')
00684
00685 if (newdata(n)) then
00686
00687 call t_startf(trim(lstr)//trim(timname)//'_coszenN')
00688 call shr_tInterp_getAvgCosz(tavCosz,lonr,latr, &
00689 SDAT%ymdLB(n),SDAT%todLB(n), SDAT%ymdUB(n),SDAT%todUB(n), &
00690 SDAT%eccen,SDAT%mvelpp,SDAT%lambm0,SDAT%obliqr,SDAT%modeldt)
00691 call mct_avect_importRAttr(SDAT%avCoszen(n),'tavCosz',tavCosz,lsizeF)
00692 call t_stopf(trim(lstr)//trim(timname)//'_coszenN')
00693 else
00694
00695 call mct_avect_exportRAttr(SDAT%avCoszen(n),'tavCosz',tavCosz)
00696 endif
00697
00698
00699 do i = 1,lsizeF
00700 if (cosz(i) > solZenMin) then
00701 SDAT%avs(n)%rAttr(:,i) = SDAT%avFLB(n)%rAttr(:,i)*cosz(i)/tavCosz(i)
00702 else
00703 SDAT%avs(n)%rAttr(:,i) = 0._r8
00704 endif
00705 enddo
00706 deallocate(tavCosz,cosz,lonr,latr)
00707 call t_stopf(trim(lstr)//trim(timname)//'_coszen')
00708
00709
00710 elseif (trim(SDAT%tintalgo(n)) /= trim(shr_strdata_nullstr)) then
00711
00712 call t_startf(trim(lstr)//trim(timname)//'_tint')
00713 call shr_tInterp_getFactors(SDAT%ymdlb(n),SDAT%todlb(n),SDAT%ymdub(n),SDAT%todub(n),SDAT%ymd,SDAT%tod, &
00714 flb,fub,algo=trim(SDAT%tintalgo(n)))
00715 SDAT%avs(n)%rAttr(:,:) = SDAT%avFLB(n)%rAttr(:,:)*flb + SDAT%avFUB(n)%rAttr(:,:)*fub
00716 call t_stopf(trim(lstr)//trim(timname)//'_tint')
00717
00718 else
00719 call t_startf(trim(lstr)//trim(timname)//'_zero')
00720 call mct_avect_zero(SDAT%avs(n))
00721 call t_stopf(trim(lstr)//trim(timname)//'_zero')
00722 endif
00723 enddo
00724
00725 deallocate(newData)
00726
00727 endif
00728
00729 if (.not.ltimers) call t_adj_detailf(-tadj)
00730
00731 end subroutine shr_strdata_advance
00732
00733
00734 subroutine shr_strdata_clean(SDAT)
00735
00736 implicit none
00737
00738 type(shr_strdata_type),intent(inout) :: SDAT
00739
00740 integer(IN) :: n
00741
00742
00743 character(len=*),parameter :: subname = "(shr_strdata_clean) "
00744
00745
00746
00747
00748
00749 if (SDAT%nxg * SDAT%nyg == 0) then
00750 return
00751 endif
00752
00753 SDAT%nxg = 0
00754 SDAT%nyg = 0
00755 SDAT%strnxg = 0
00756 SDAT%strnyg = 0
00757
00758 SDAT%nstreams = 0
00759 SDAT%nvectors = 0
00760 SDAT%ustrm = 0
00761 SDAT%vstrm = 0
00762
00763 SDAT%dofill = .false.
00764 SDAT%domaps = .false.
00765
00766 call mct_ggrid_clean(SDAT%grid)
00767 call mct_gsmap_clean(SDAT%gsmap)
00768 do n = 1,nStrMax
00769
00770 call mct_avect_clean(SDAT%avs(n))
00771 call mct_avect_clean(SDAT%avRLB(n))
00772 call mct_avect_clean(SDAT%avRUB(n))
00773 call mct_avect_clean(SDAT%avFLB(n))
00774 call mct_avect_clean(SDAT%avFUB(n))
00775 call mct_ggrid_clean(SDAT%gridR(n))
00776 call mct_sMatP_clean(SDAT%sMatPf(n))
00777 call mct_sMatP_clean(SDAT%sMatPs(n))
00778 call mct_gsmap_clean(SDAT%gsmapR(n))
00779 enddo
00780
00781 end subroutine shr_strdata_clean
00782
00783
00784
00785 subroutine shr_strdata_restWrite(filename,SDAT,mpicom,str1,str2)
00786
00787 implicit none
00788
00789 character(len=*) ,intent(in) :: filename
00790 type(shr_strdata_type),intent(inout) :: SDAT
00791 integer(IN) ,intent(in) :: mpicom
00792 character(len=*) ,intent(in) :: str1
00793 character(len=*) ,intent(in) :: str2
00794
00795
00796 type(shr_stream_streamtype),pointer :: streams(:)
00797 integer(IN) :: n,my_task,ier
00798
00799
00800 character(len=*),parameter :: subname = "(shr_strdata_restWrite) "
00801
00802
00803
00804
00805
00806 call MPI_COMM_RANK(mpicom,my_task,ier)
00807
00808 if (my_task == 0) then
00809 call shr_stream_restWrite(SDAT%stream,trim(filename),trim(str1),trim(str2),SDAT%nstreams)
00810 endif
00811
00812 end subroutine shr_strdata_restWrite
00813
00814
00815
00816 subroutine shr_strdata_restRead(filename,SDAT,mpicom)
00817
00818 implicit none
00819
00820 character(len=*) ,intent(in) :: filename
00821 type(shr_strdata_type),intent(inout) :: SDAT
00822 integer(IN) ,intent(in) :: mpicom
00823
00824
00825 type(shr_stream_streamtype),pointer :: streams(:)
00826 integer(IN) :: n,my_task,ier
00827
00828
00829 character(len=*),parameter :: subname = "(shr_strdata_restRead) "
00830
00831
00832
00833
00834
00835 call MPI_COMM_RANK(mpicom,my_task,ier)
00836
00837 if (my_task == 0) then
00838 call shr_stream_restRead(SDAT%stream,trim(filename),SDAT%nstreams)
00839 endif
00840
00841 end subroutine shr_strdata_restRead
00842
00843
00844
00845 subroutine shr_strdata_setOrbs(SDAT,eccen,mvelpp,lambm0,obliqr,modeldt)
00846
00847 implicit none
00848
00849 type(shr_strdata_type),intent(inout) :: SDAT
00850 real(R8),intent(in) :: eccen
00851 real(R8),intent(in) :: mvelpp
00852 real(R8),intent(in) :: lambm0
00853 real(R8),intent(in) :: obliqr
00854 integer(IN),intent(in) :: modeldt
00855
00856
00857 character(len=*),parameter :: subname = "(shr_strdata_setOrbs) "
00858
00859
00860
00861
00862
00863 SDAT%eccen = eccen
00864 SDAT%mvelpp = mvelpp
00865 SDAT%lambm0 = lambm0
00866 SDAT%obliqr = obliqr
00867 SDAT%modeldt = modeldt
00868
00869 end subroutine shr_strdata_setOrbs
00870
00871
00872
00873
00874
00875
00876
00877
00878
00879
00880
00881
00882
00883
00884
00885 subroutine shr_strdata_readnml(SDAT,file,rc,mpicom)
00886
00887 implicit none
00888
00889
00890
00891 type(shr_strdata_type),intent(inout):: SDAT
00892 character(*),optional ,intent(in) :: file
00893 integer(IN),optional ,intent(out) :: rc
00894 integer(IN),optional ,intent(in) :: mpicom
00895
00896
00897
00898 integer(IN) :: rCode
00899 integer(IN) :: nUnit
00900 integer(IN) :: n
00901 integer(IN) :: my_task
00902 integer(IN) :: master_task
00903 integer(IN) :: ntasks
00904
00905
00906 character(CL) :: dataMode
00907 character(CL) :: domainFile
00908 character(CL) :: streams(nStrMax)
00909 character(CL) :: taxMode(nStrMax)
00910 real(R8) :: dtlimit(nStrMax)
00911 character(CL) :: vectors(nVecMax)
00912 character(CL) :: fillalgo(nStrMax)
00913 character(CL) :: fillmask(nStrMax)
00914 character(CL) :: fillread(nStrMax)
00915 character(CL) :: fillwrite(nStrMax)
00916 character(CL) :: mapalgo(nStrMax)
00917 character(CL) :: mapmask(nStrMax)
00918 character(CL) :: mapread(nStrMax)
00919 character(CL) :: mapwrite(nStrMax)
00920 character(CL) :: tintalgo(nStrMax)
00921 character(CL) :: io_type
00922 integer(IN) :: num_iotasks
00923 integer(IN) :: io_root
00924 integer(IN) :: io_stride
00925 integer(IN) :: num_agg
00926 character(CL) :: fileName
00927 integer(IN) :: yearFirst
00928 integer(IN) :: yearLast
00929 integer(IN) :: yearAlign
00930
00931
00932 namelist / shr_strdata_nml / &
00933 dataMode &
00934 , domainFile &
00935 , streams &
00936 , taxMode &
00937 , dtlimit &
00938 , vectors &
00939 , fillalgo &
00940 , fillmask &
00941 , fillread &
00942 , fillwrite &
00943 , mapalgo &
00944 , mapmask &
00945 , mapread &
00946 , mapwrite &
00947 , tintalgo &
00948 , io_type &
00949 , num_iotasks &
00950 , io_root &
00951 , io_stride &
00952 , num_agg
00953
00954
00955 character(*),parameter :: subName = "(shr_strdata_readnml) "
00956 character(*),parameter :: F00 = "('(shr_strdata_readnml) ',8a)"
00957 character(*),parameter :: F01 = "('(shr_strdata_readnml) ',a,i6,a)"
00958 character(*),parameter :: F02 = "('(shr_strdata_readnml) ',a,es13.6)"
00959 character(*),parameter :: F03 = "('(shr_strdata_readnml) ',a,l6)"
00960 character(*),parameter :: F04 = "('(shr_strdata_readnml) ',a,i2,a,a)"
00961 character(*),parameter :: F20 = "('(shr_strdata_readnml) ',a,i6,a)"
00962 character(*),parameter :: F90 = "('(shr_strdata_readnml) ',58('-'))"
00963
00964
00965
00966
00967
00968 if (present(rc)) rc = 0
00969
00970 my_task = 0
00971 master_task = 0
00972 ntasks = 1
00973 if (present(mpicom)) then
00974 call MPI_COMM_RANK(mpicom,my_task,rCode)
00975 call MPI_COMM_SIZE(mpicom,ntasks,rCode)
00976 endif
00977
00978
00979 if (my_task == master_task) then
00980
00981
00982
00983
00984 dataMode = 'NULL'
00985 domainFile = trim(shr_strdata_nullstr)
00986 streams(:) = trim(shr_strdata_nullstr)
00987 taxMode(:) = trim(shr_stream_taxis_cycle)
00988 dtlimit(:) = dtlimit_default
00989 vectors(:) = trim(shr_strdata_nullstr)
00990 fillalgo(:) = 'nn'
00991 fillmask(:) = 'nomask'
00992 fillread(:) = trim(shr_strdata_unset)
00993 fillwrite(:)= trim(shr_strdata_unset)
00994 mapalgo(:) = 'bilinear'
00995 mapmask(:) = 'dstmask'
00996 mapread(:) = trim(shr_strdata_unset)
00997 mapwrite(:) = trim(shr_strdata_unset)
00998 tintalgo(:) = 'linear'
00999 io_type = 'pio_netcdf'
01000
01001 num_iotasks = 1
01002 io_root = 0
01003 io_stride = 4
01004 num_agg = 0
01005
01006
01007
01008
01009 if (present(file)) then
01010 write(logunit,F00) 'reading input namelist file: ',trim(file)
01011 nUnit = shr_file_getUnit()
01012 open (nUnit,file=trim(file),status="old",action="read")
01013 read (nUnit,nml=shr_strdata_nml,iostat=rCode)
01014 close(nUnit)
01015 call shr_file_freeUnit(nUnit)
01016 if (rCode > 0) then
01017 write(logunit,F01) 'ERROR: reading input namelist, '//trim(file)//' iostat=',rCode
01018 call shr_sys_abort(subName//": namelist read error "//trim(file))
01019 end if
01020 endif
01021
01022
01023
01024
01025 SDAT%dataMode = dataMode
01026 SDAT%domainFile = domainFile
01027 SDAT%streams(:) = streams(:)
01028 SDAT%taxMode(:) = taxMode(:)
01029 SDAT%dtlimit(:) = dtlimit(:)
01030 SDAT%vectors(:) = vectors(:)
01031 SDAT%fillalgo(:) = fillalgo(:)
01032 SDAT%fillmask(:) = fillmask(:)
01033 SDAT%fillread(:) = fillread(:)
01034 SDAT%fillwrit(:) = fillwrite(:)
01035 SDAT%mapalgo(:) = mapalgo(:)
01036 SDAT%mapmask(:) = mapmask(:)
01037 SDAT%mapread(:) = mapread(:)
01038 SDAT%mapwrit(:) = mapwrite(:)
01039 SDAT%tintalgo(:) = tintalgo(:)
01040 do n=1,nStrMax
01041 if (trim(streams(n)) /= trim(shr_strdata_nullstr)) SDAT%nstreams = max(SDAT%nstreams,n)
01042 if (trim(SDAT%taxMode(n)) == trim(shr_stream_taxis_extend)) SDAT%dtlimit(n) = 1.0e30
01043 end do
01044 SDAT%nvectors = 0
01045 do n=1,nVecMax
01046 if (trim(vectors(n)) /= trim(shr_strdata_nullstr)) SDAT%nvectors = n
01047 end do
01048
01049 do n = 1,SDAT%nstreams
01050 call shr_stream_parseInput(SDAT%streams(n),fileName,yearAlign,yearFirst,yearLast)
01051 call shr_stream_init(SDAT%stream(n),fileName,yearFirst,yearLast,yearAlign, &
01052 trim(SDAT%taxMode(n)))
01053 enddo
01054
01055 if (trim(io_type) == 'netcdf') then
01056 SDAT%io_type = iotype_std_netcdf
01057 elseif (trim(io_type) == 'pio_netcdf') then
01058 SDAT%io_type = iotype_netcdf
01059 elseif (trim(io_type) == 'pio_pnetcdf') then
01060 SDAT%io_type = iotype_pnetcdf
01061 else
01062 write(logunit,F00) 'ERROR: unknown io_type, '//trim(io_type)
01063 call shr_sys_abort(subName//": unknown io_type, "//trim(io_type))
01064 endif
01065
01066 SDAT%num_iotasks = num_iotasks
01067 SDAT%io_root = io_root
01068 SDAT%io_stride = io_stride
01069 SDAT%num_agg = num_agg
01070
01071
01072
01073 endif
01074
01075
01076 if (present(mpicom)) then
01077 call shr_strdata_bcastnml(SDAT,mpicom)
01078 endif
01079
01080 SDAT%ymdLB = -1
01081 SDAT%todLB = -1
01082 SDAT%ymdUB = -1
01083 SDAT%todUB = -1
01084 SDAT%dtmin = 1.0e30
01085 SDAT%dtmax = 0.0
01086 SDAT%eccen = SHR_ORB_UNDEF_REAL
01087 SDAT%mvelpp = SHR_ORB_UNDEF_REAL
01088 SDAT%lambm0 = SHR_ORB_UNDEF_REAL
01089 SDAT%obliqr = SHR_ORB_UNDEF_REAL
01090 SDAT%modeldt = 0
01091
01092 end subroutine shr_strdata_readnml
01093
01094
01095
01096
01097
01098
01099
01100
01101
01102
01103
01104
01105
01106
01107
01108 subroutine shr_strdata_create(SDAT, name, mpicom, compid, gsmap, ggrid, nxg, nyg, &
01109
01110 yearFirst, yearLast, yearAlign, offset, &
01111 domFilePath, domFileName, &
01112 domTvarName, domXvarName, domYvarName, domAreaName, domMaskName, &
01113 filePath, filename, fldListFile, fldListModel, &
01114
01115 taxMode, dtlimit, tintalgo, &
01116 fillalgo, fillmask, fillread, fillwrite, &
01117 mapalgo, mapmask, mapread, mapwrite, &
01118 io_type, num_iotasks, io_root, io_stride, num_agg)
01119
01120 implicit none
01121
01122
01123
01124 type(shr_strdata_type),intent(inout):: SDAT
01125 character(*) ,intent(in) :: name
01126 integer(IN) ,intent(in) :: mpicom
01127 integer(IN) ,intent(in) :: compid
01128 type(mct_gsmap) ,intent(in) :: gsmap
01129 type(mct_ggrid) ,intent(in) :: ggrid
01130 integer(IN) ,intent(in) :: nxg
01131 integer(IN) ,intent(in) :: nyg
01132
01133 integer(IN) ,intent(in) :: yearFirst
01134 integer(IN) ,intent(in) :: yearLast
01135 integer(IN) ,intent(in) :: yearAlign
01136 integer(IN) ,intent(in) :: offset
01137 character(*) ,intent(in) :: domFilePath
01138 character(*) ,intent(in) :: domFileName
01139 character(*) ,intent(in) :: domTvarName
01140 character(*) ,intent(in) :: domXvarName
01141 character(*) ,intent(in) :: domYvarName
01142 character(*) ,intent(in) :: domAreaName
01143 character(*) ,intent(in) :: domMaskName
01144 character(*) ,intent(in) :: filePath
01145 character(*) ,intent(in) :: filename(:)
01146 character(*) ,intent(in) :: fldListFile
01147 character(*) ,intent(in) :: fldListModel
01148
01149 character(*),optional ,intent(in) :: taxMode
01150 real(R8) ,optional ,intent(in) :: dtlimit
01151 character(*),optional ,intent(in) :: fillalgo
01152 character(*),optional ,intent(in) :: fillmask
01153 character(*),optional ,intent(in) :: fillread
01154 character(*),optional ,intent(in) :: fillwrite
01155 character(*),optional ,intent(in) :: mapalgo
01156 character(*),optional ,intent(in) :: mapmask
01157 character(*),optional ,intent(in) :: mapread
01158 character(*),optional ,intent(in) :: mapwrite
01159 character(*),optional ,intent(in) :: tintalgo
01160 character(*),optional ,intent(in) :: io_type
01161 integer(IN) ,optional ,intent(in) :: num_iotasks
01162 integer(IN) ,optional ,intent(in) :: io_root
01163 integer(IN) ,optional ,intent(in) :: io_stride
01164 integer(IN) ,optional ,intent(in) :: num_agg
01165
01166
01167
01168
01169
01170 character(*),parameter :: subName = "(shr_strdata_create) "
01171 character(*),parameter :: F00 = "('(shr_strdata_create) ',8a)"
01172
01173
01174
01175
01176
01177 call shr_strdata_readnml(SDAT,mpicom=mpicom)
01178
01179 SDAT%nstreams = 1
01180
01181 if (present(taxMode)) then
01182 SDAT%taxMode(1) = taxMode
01183 if (trim(SDAT%taxMode(1)) == trim(shr_stream_taxis_extend)) SDAT%dtlimit(1) = 1.0e30
01184 endif
01185 if (present(dtlimit)) then
01186 SDAT%dtlimit(1) = dtlimit
01187 endif
01188 if (present(fillalgo)) then
01189 SDAT%fillalgo(1) = fillalgo
01190 endif
01191 if (present(fillmask)) then
01192 SDAT%fillmask(1) = fillmask
01193 endif
01194 if (present(fillread)) then
01195 SDAT%fillread(1) = fillread
01196 endif
01197 if (present(fillwrite)) then
01198 SDAT%fillwrit(1) = fillwrite
01199 endif
01200 if (present(mapalgo)) then
01201 SDAT%mapalgo(1) = mapalgo
01202 endif
01203 if (present(mapmask)) then
01204 SDAT%mapmask(1) = mapmask
01205 endif
01206 if (present(mapread)) then
01207 SDAT%mapread(1) = mapread
01208 endif
01209 if (present(mapwrite)) then
01210 SDAT%mapwrit(1) = mapwrite
01211 endif
01212 if (present(tintalgo)) then
01213 SDAT%tintalgo(1) = tintalgo
01214 endif
01215 if (present(mapmask)) then
01216 SDAT%mapmask(1) = mapmask
01217 endif
01218
01219 if (present(io_type)) then
01220 if (trim(io_type) == 'netcdf') then
01221 SDAT%io_type = iotype_std_netcdf
01222 elseif (trim(io_type) == 'pio_netcdf') then
01223 SDAT%io_type = iotype_netcdf
01224 elseif (trim(io_type) == 'pio_pnetcdf') then
01225 SDAT%io_type = iotype_pnetcdf
01226 else
01227 write(logunit,F00) 'ERROR: unknown io_type, '//trim(io_type)
01228 call shr_sys_abort(subName//": unknown io_type, "//trim(io_type))
01229 endif
01230 endif
01231 if (present(num_iotasks)) then
01232 SDAT%num_iotasks = num_iotasks
01233 endif
01234 if (present(io_root)) then
01235 SDAT%io_root = io_root
01236 endif
01237 if (present(io_stride)) then
01238 SDAT%io_stride = io_stride
01239 endif
01240 if (present(num_agg)) then
01241 SDAT%num_agg = num_agg
01242 endif
01243
01244 call shr_stream_set(SDAT%stream(1),yearFirst,yearLast,yearAlign,offset,taxMode, &
01245 fldListFile,fldListModel,domFilePath,domFileName, &
01246 domTvarName,domXvarName,domYvarName,domAreaName,domMaskName, &
01247 filePath,filename,trim(name))
01248
01249 call shr_strdata_init(SDAT, mpicom, compid, &
01250 gsmap=gsmap, ggrid=ggrid, nxg=nxg, nyg=nyg)
01251
01252 end subroutine shr_strdata_create
01253
01254
01255
01256
01257
01258
01259
01260
01261
01262
01263
01264
01265
01266
01267
01268 subroutine shr_strdata_print(SDAT,name)
01269
01270 implicit none
01271
01272
01273
01274 type(shr_strdata_type) ,intent(in) :: SDAT
01275 character(len=*),optional,intent(in) :: name
01276
01277
01278
01279 integer(IN) :: n
01280 character(CL) :: lname
01281
01282
01283 character(*),parameter :: subName = "(shr_strdata_print) "
01284 character(*),parameter :: F00 = "('(shr_strdata_print) ',8a)"
01285 character(*),parameter :: F01 = "('(shr_strdata_print) ',a,i6,a)"
01286 character(*),parameter :: F02 = "('(shr_strdata_print) ',a,es13.6)"
01287 character(*),parameter :: F03 = "('(shr_strdata_print) ',a,l6)"
01288 character(*),parameter :: F04 = "('(shr_strdata_print) ',a,i2,a,a)"
01289 character(*),parameter :: F05 = "('(shr_strdata_print) ',a,i2,a,i6)"
01290 character(*),parameter :: F06 = "('(shr_strdata_print) ',a,i2,a,l2)"
01291 character(*),parameter :: F07 = "('(shr_strdata_print) ',a,i2,a,es13.6)"
01292 character(*),parameter :: F20 = "('(shr_strdata_print) ',a,i6,a)"
01293 character(*),parameter :: F90 = "('(shr_strdata_print) ',58('-'))"
01294
01295
01296
01297
01298
01299 lname = 'unknown'
01300 if (present(name)) then
01301 lname = trim(name)
01302 endif
01303
01304
01305
01306 write(logunit,F90)
01307 write(logunit,F00) "name = ",trim(lname)
01308 write(logunit,F00) "dataMode = ",trim(SDAT%dataMode)
01309 write(logunit,F00) "domainFile = ",trim(SDAT%domainFile)
01310 write(logunit,F01) "nxg = ",SDAT%nxg
01311 write(logunit,F01) "nyg = ",SDAT%nyg
01312 write(logunit,F01) "io_type = ",SDAT%io_type
01313 write(logunit,F01) "num_iotasks = ",SDAT%num_iotasks
01314 write(logunit,F01) "io_root = ",SDAT%io_root
01315 write(logunit,F01) "io_stride = ",SDAT%io_stride
01316 write(logunit,F01) "num_agg = ",SDAT%num_agg
01317 write(logunit,F02) "eccen = ",SDAT%eccen
01318 write(logunit,F02) "mvelpp = ",SDAT%mvelpp
01319 write(logunit,F02) "lambm0 = ",SDAT%lambm0
01320 write(logunit,F02) "obliqr = ",SDAT%obliqr
01321 write(logunit,F01) "nstreams = ",SDAT%nstreams
01322 do n=1, SDAT%nstreams
01323 write(logunit,F04) " streams (",n,") = ",trim(SDAT%streams(n))
01324 write(logunit,F04) " taxMode (",n,") = ",trim(SDAT%taxMode(n))
01325 write(logunit,F07) " dtlimit (",n,") = ",SDAT%dtlimit(n)
01326 write(logunit,F05) " strnxg (",n,") = ",SDAT%strnxg(n)
01327 write(logunit,F05) " strnyg (",n,") = ",SDAT%strnyg(n)
01328 write(logunit,F06) " dofill (",n,") = ",SDAT%dofill(n)
01329 write(logunit,F04) " fillalgo(",n,") = ",trim(SDAT%fillalgo(n))
01330 write(logunit,F04) " fillmask(",n,") = ",trim(SDAT%fillmask(n))
01331 write(logunit,F04) " fillread(",n,") = ",trim(SDAT%fillread(n))
01332 write(logunit,F04) " fillwrit(",n,") = ",trim(SDAT%fillwrit(n))
01333 write(logunit,F06) " domaps (",n,") = ",SDAT%domaps(n)
01334 write(logunit,F04) " mapalgo (",n,") = ",trim(SDAT%mapalgo(n))
01335 write(logunit,F04) " mapmask (",n,") = ",trim(SDAT%mapmask(n))
01336 write(logunit,F04) " mapread (",n,") = ",trim(SDAT%mapread(n))
01337 write(logunit,F04) " mapwrit (",n,") = ",trim(SDAT%mapwrit(n))
01338 write(logunit,F04) " tintalgo(",n,") = ",trim(SDAT%tintalgo(n))
01339 write(logunit,F01) " "
01340 end do
01341 write(logunit,F01) "nvectors = ",SDAT%nvectors
01342 do n=1, SDAT%nvectors
01343 write(logunit,F04) " vectors (",n,") = ",trim(SDAT%vectors(n))
01344 end do
01345 write(logunit,F90)
01346 call shr_sys_flush(logunit)
01347
01348 end subroutine shr_strdata_print
01349
01350
01351
01352
01353
01354
01355
01356
01357
01358
01359
01360
01361
01362
01363 subroutine shr_strdata_bcastnml(SDAT,mpicom,rc)
01364
01365 use shr_mpi_mod, only : shr_mpi_bcast
01366
01367 implicit none
01368
01369
01370
01371 type(shr_strdata_type),intent(inout) :: SDAT
01372 integer(IN) ,intent(in) :: mpicom
01373 integer(IN),optional ,intent(out) :: rc
01374
01375
01376
01377
01378 integer(IN) :: lrc
01379
01380
01381 character(*),parameter :: subName = "(shr_strdata_bcastnml) "
01382
01383
01384
01385
01386
01387 lrc = 0
01388
01389 call shr_mpi_bcast(SDAT%dataMode ,mpicom,'dataMode')
01390 call shr_mpi_bcast(SDAT%domainFile,mpicom,'domainFile')
01391 call shr_mpi_bcast(SDAT%nstreams ,mpicom,'nstreams')
01392 call shr_mpi_bcast(SDAT%nvectors ,mpicom,'nvectors')
01393 call shr_mpi_bcast(SDAT%streams ,mpicom,'streams')
01394 call shr_mpi_bcast(SDAT%taxMode ,mpicom,'taxMode')
01395 call shr_mpi_bcast(SDAT%dtlimit ,mpicom,'dtlimit')
01396 call shr_mpi_bcast(SDAT%vectors ,mpicom,'vectors')
01397 call shr_mpi_bcast(SDAT%fillalgo ,mpicom,'fillalgo')
01398 call shr_mpi_bcast(SDAT%fillmask ,mpicom,'fillmask')
01399 call shr_mpi_bcast(SDAT%fillread ,mpicom,'fillread')
01400 call shr_mpi_bcast(SDAT%fillwrit ,mpicom,'fillwrit')
01401 call shr_mpi_bcast(SDAT%mapalgo ,mpicom,'mapalgo')
01402 call shr_mpi_bcast(SDAT%mapmask ,mpicom,'mapmask')
01403 call shr_mpi_bcast(SDAT%mapread ,mpicom,'mapread')
01404 call shr_mpi_bcast(SDAT%mapwrit ,mpicom,'mapwrit')
01405 call shr_mpi_bcast(SDAT%tintalgo ,mpicom,'tintalgo')
01406
01407 call shr_mpi_bcast(SDAT%io_type ,mpicom,'io_type')
01408 call shr_mpi_bcast(SDAT%num_iotasks,mpicom,'num_iotasks')
01409 call shr_mpi_bcast(SDAT%io_root ,mpicom,'io_root')
01410 call shr_mpi_bcast(SDAT%io_stride ,mpicom,'io_stride')
01411 call shr_mpi_bcast(SDAT%num_agg ,mpicom,'num_agg')
01412
01413 if (present(rc)) then
01414 rc = lrc
01415 endif
01416
01417 end subroutine shr_strdata_bcastnml
01418
01419
01420
01421 subroutine shr_strdata_setlogunit(nu)
01422
01423 integer(IN),intent(in) :: nu
01424 character(len=*),parameter :: subname = "(shr_strdata_setlogunit) "
01425
01426
01427
01428 end subroutine shr_strdata_setlogunit
01429
01430
01431
01432
01433 end module shr_strdata_mod
01434