00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018 module shr_scam_mod
00019
00020
00021
00022 use shr_kind_mod
00023 use shr_sys_mod
00024 use shr_file_mod
00025 use shr_kind_mod, only : R8=>SHR_KIND_R8,IN=>SHR_KIND_IN,CL=>SHR_KIND_CL
00026 use shr_log_mod, only : s_loglev => shr_log_Level
00027 use shr_log_mod, only : s_logunit => shr_log_Unit
00028 use netcdf
00029
00030 implicit none
00031
00032 private
00033
00034
00035
00036
00037
00038
00039
00040 public :: shr_scam_getCloseLatLon
00041 public :: shr_scam_checkSurface
00042
00043
00044
00045
00046
00047
00048
00049 CONTAINS
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070 subroutine shr_scam_getCloseLatLon(ncid, targetLat, targetLon, closeLat, closeLon, &
00071 closeLatIdx, closeLonIdx)
00072 implicit none
00073
00074
00075
00076 integer(IN),intent(in) :: ncid
00077 real (R8),intent(in) :: targetLat
00078 real (R8),intent(in) :: targetLon
00079 real (R8),intent(out) :: closeLat
00080 real (R8),intent(out) :: closeLon
00081 integer(IN),intent(out) :: closeLatIdx
00082 integer(IN),intent(out) :: closeLonIdx
00083
00084
00085
00086
00087 real (R8),allocatable :: lats(:),lons(:),poslons(:)
00088 real (R8) :: postargetlon
00089 integer(IN) :: rcode
00090 integer(IN) :: i
00091 integer(IN) :: len
00092 integer(IN) :: ndims
00093 integer(IN) :: nvars
00094 integer(IN) :: nvarid
00095 integer(IN) :: ndimid
00096 integer(IN) :: strt(nf90_max_var_dims),count(nf90_max_var_dims)
00097 integer(IN) :: nlon,nlat
00098 integer(IN), dimension(nf90_max_var_dims) :: dimids
00099 character(len=80) :: name,var_name
00100 character(*),parameter :: subname = "(shr_scam_getCloseLatLon) "
00101
00102
00103
00104
00105
00106
00107
00108 rcode = nf90_inquire(ncid, nVariables=nvars)
00109 if (rcode /= nf90_noerr) then
00110 call shr_sys_abort(subname//"ERROR from nf90_inquire ")
00111 endif
00112
00113
00114
00115 nlat=0
00116 nlon=0
00117 nvarid=0
00118
00119
00120
00121 do while (nvarid < nvars .and.(nlon.eq.0 .or. nlat.eq.0))
00122 nvarid=nvarid+1
00123 rcode = nf90_inquire_variable(ncid, nvarid, var_name, ndims=ndims,dimids = dimids)
00124 if (rcode /= nf90_noerr) then
00125 call shr_sys_abort(subname//"ERROR inquiring about variable "//trim(var_name))
00126 endif
00127
00128
00129
00130 if ( var_name .eq. 'lat'.or. var_name .eq. 'latixy'.or. var_name .eq. 'yc'.or.var_name.eq.'lsmlat'.or.&
00131 var_name .eq. 'LAT'.or. var_name .eq. 'LATIXY'.or. var_name .eq. 'YC'.or.var_name.eq.'LSMLAT' ) then
00132
00133
00134
00135 do ndimid = 1,ndims
00136 rcode = nf90_inquire_dimension(ncid, dimids(ndimid), name, len)
00137 if (rcode /= nf90_noerr) then
00138 call shr_sys_abort(subname//"ERROR: Cant read netcdf latitude variable dimension")
00139 endif
00140 if ( name .eq. 'lat'.or. name .eq. 'latixy'.or. name .eq. 'nj'.or. name .eq. 'lsmlat' .or. &
00141 name .eq. 'LAT'.or. name .eq. 'LATIXY'.or. name .eq. 'NJ'.or. name .eq. 'LSMLAT' ) then
00142 strt(ndimid) = 1
00143 count(ndimid) = len
00144 nlat=len
00145 else
00146 strt(ndimid) = 1
00147 count(ndimid) = 1
00148 endif
00149 end do
00150 if (nlat.eq.0) then
00151 call shr_sys_abort( subname//"ERROR: Cant find a useable latitude dimension (lat,nj,latixy, or lsmlat")
00152 endif
00153 allocate(lats(nlat))
00154 rcode= nf90_get_var(ncid, nvarid ,lats, start = strt, count = count)
00155 if (rcode /= nf90_noerr) then
00156 call shr_sys_abort( subname//"ERROR: Cant read netcdf latitude variable dimension")
00157 endif
00158 end if
00159
00160
00161
00162 if ( var_name .eq. 'lon'.or. var_name .eq. 'longxy'.or. var_name .eq. 'xc'.or.var_name.eq.'lsmlon'.or.&
00163 var_name .eq. 'LON'.or. var_name .eq. 'LONGXY'.or. var_name .eq. 'XC'.or.var_name.eq.'LSMLON' ) then
00164 do ndimid = 1,ndims
00165 rcode = nf90_inquire_dimension(ncid, dimids(ndimid), name, len)
00166 if (rcode /= nf90_noerr) then
00167 call shr_sys_abort( subname//"ERROR: Cant read netcdf latitude variable dimension")
00168 endif
00169 if ( name .eq. 'lon'.or. name .eq. 'longxy'.or. name .eq. 'ni'.or. name .eq. 'lsmlon' .or. &
00170 name .eq. 'LON'.or. name .eq. 'LONGXY'.or. name .eq. 'NI'.or. name .eq. 'LSMLON' ) then
00171 strt(ndimid) = 1
00172 count(ndimid) = len
00173 nlon=len
00174 else
00175 strt(ndimid) = 1
00176 count(ndimid) = 1
00177 endif
00178 end do
00179 if (nlon.eq.0) then
00180 call shr_sys_abort( subname//"ERROR: Cant find a useable longitude dimension (lon,ni,longxy, or lsmlon")
00181 endif
00182 allocate(lons(nlon))
00183 allocate(poslons(nlon))
00184 rcode= nf90_get_var(ncid, nvarid ,lons, start = strt, count = count)
00185 if (rcode /= nf90_noerr) then
00186 call shr_sys_abort( subname//"ERROR: Cant read netcdf latitude variable dimension")
00187 endif
00188 end if
00189 end do
00190
00191
00192
00193 if (nlon.eq.0) then
00194 call shr_sys_abort( subname//"ERROR: Couldnt find a longitude coordinate variable")
00195 end if
00196 if (nlat.eq.0) then
00197 call shr_sys_abort( subname//"ERROR: Couldnt find a latitude coordinate variable")
00198 end if
00199
00200
00201
00202 poslons=mod(lons+360._r8,360._r8)
00203 postargetlon=mod(targetlon+360._r8,360._r8)
00204
00205
00206
00207 closelonidx=(MINLOC(abs(poslons-postargetlon),dim=1))
00208 closelatidx=(MINLOC(abs(lats-targetlat),dim=1))
00209 closelon=lons(closelonidx)
00210 closelat=lats(closelatidx)
00211
00212
00213
00214 deallocate(lats)
00215 deallocate(lons)
00216 deallocate(poslons)
00217
00218 return
00219
00220 end subroutine shr_scam_getCloseLatLon
00221
00222
00223
00224
00225
00226
00227
00228
00229
00230
00231
00232
00233
00234
00235
00236
00237 subroutine shr_scam_checkSurface(scmlon, scmlat, ocn_compid, ocn_mpicom, lnd_present, ocn_present, ice_present)
00238
00239
00240 use shr_strdata_mod
00241 use shr_dmodel_mod
00242 use mct_mod
00243 implicit none
00244
00245
00246
00247 real(R8), intent(in) :: scmlon,scmlat
00248 integer(IN), intent(in) :: ocn_compid
00249 integer(IN), intent(in) :: ocn_mpicom
00250 logical, optional, intent(out) :: lnd_present
00251 logical, optional, intent(out) :: ice_present
00252 logical, optional, intent(out) :: ocn_present
00253
00254
00255
00256
00257 integer(IN) :: rcode
00258 integer(IN) :: ncid_ocn
00259 integer(IN) :: fracid
00260 integer(IN) :: closeLatIdx
00261 integer(IN) :: closeLonIdx
00262 integer(IN) :: unitn
00263 real (R8) :: ocn_frac(1,1)
00264 real (R8) :: closeLat
00265 real (R8) :: closeLon
00266 character(len=CL) :: nrevsn = ' '
00267 character(len=CL) :: rest_pfile = './rpointer.dom'
00268 character(len=CL) :: bndtvs
00269 character(len=CL) :: focndomain
00270 logical :: sstcyc
00271 logical :: docn_exists
00272 logical :: ocn_exists
00273 logical :: exists
00274 logical :: aqua_planet
00275 logical :: single_column
00276
00277
00278 character(*),parameter :: subname = "(shr_scam_checkSurface) "
00279 character(*),parameter :: F00 = "('(shr_scam_checkSurface) ',8a)"
00280 type(shr_strdata_type) :: SDAT
00281 character(len=CL) :: decomp = '1d'
00282 character(len=CL) :: restfilm = 'unset'
00283 character(len=CL) :: restfils = 'unset'
00284 character(len=CL) :: ocn_in = 'unset'
00285 integer(IN) :: nfrac
00286 namelist /dom_inparm/ sstcyc, nrevsn, rest_pfile, bndtvs, focndomain
00287 namelist / docn_nml / ocn_in, decomp, restfilm, restfils
00288
00289
00290
00291
00292
00293 inquire( file='ocn_in', exist=ocn_exists )
00294 inquire( file='docn_in', exist=docn_exists )
00295 if (ocn_exists) then
00296
00297
00298 unitn = shr_file_getUnit()
00299 open( unitn, file='ocn_in', status='old' )
00300 rcode = 1
00301 do while ( rcode /= 0 )
00302 read(unitn, dom_inparm, iostat=rcode)
00303 if (rcode < 0) then
00304 call shr_sys_abort( 'shr_scam_checkSurface encountered end-of-file on namelist read' )
00305 endif
00306 end do
00307 close( unitn )
00308 call shr_file_freeUnit(unitn)
00309
00310
00311
00312 inquire(file=trim(focndomain),exist=exists)
00313 if (.not.exists) call shr_sys_abort(subName//"ERROR: file does not exist: "//trim(focndomain))
00314 rcode = nf90_open(focndomain,nf90_nowrite,ncid_ocn)
00315 if (rCode /= nf90_noerr) call shr_sys_abort(subName//"ERROR opening data file : "//trim(focndomain))
00316 if (s_loglev > 0) write(s_logunit,F00) 'opened netCDF data file: ',trim(focndomain)
00317
00318
00319
00320 call shr_scam_getCloseLatLon(ncid_ocn,scmlat,scmlon,closelat,closelon,closelatidx,closelonidx)
00321 rcode = nf90_inq_varid(ncid_ocn, 'frac', fracid)
00322 if (rcode /= nf90_noerr) then
00323 call shr_sys_abort(subname//"ERROR getting varid from variable frac in file "//trim(focndomain))
00324 end if
00325 rcode = nf90_get_var(ncid_ocn,fracid,ocn_frac,start=(/closelonidx,closelatidx/),count=(/1,1/))
00326 if (rcode /= nf90_noerr) then
00327 call shr_sys_abort(subname//"ERROR getting ocean fraction from "//trim(focndomain))
00328 end if
00329
00330
00331
00332 if ( present(ocn_present) ) ocn_present=(ocn_frac(1,1).gt.0.)
00333 if ( present(ocn_present).and.present(ice_present)) ice_present=ocn_present
00334 if ( present(lnd_present)) lnd_present=(ocn_frac(1,1).lt.1.)
00335 else if (docn_exists) then
00336
00337
00338 unitn = shr_file_getUnit()
00339 open( unitn, file='docn_in', status='old' )
00340 rcode = 1
00341 do while ( rcode /= 0 )
00342 read (unitn,nml=docn_nml,iostat=rcode)
00343 if (rcode < 0) then
00344 call shr_sys_abort( 'shr_scam_checkSurface encountered end-of-file on namelist read' )
00345 endif
00346 end do
00347 close( unitn )
00348 call shr_file_freeUnit(unitn)
00349 call shr_strdata_readnml(SDAT,ocn_in)
00350 call shr_dmodel_readgrid(SDAT%grid,SDAT%gsmap,SDAT%nxg,SDAT%nyg, &
00351 SDAT%domainfile, ocn_compid, ocn_mpicom, '1d', readfrac=.true., &
00352 scmmode=.true.,scmlon=scmlon,scmlat=scmlat)
00353 nfrac = mct_aVect_indexRA(SDAT%grid%data,'frac')
00354 if ( present(ocn_present) ) ocn_present=(SDAT%grid%data%rAttr(nfrac,1).gt.0.)
00355 if ( present(ocn_present).and.present(ice_present)) ice_present=ocn_present
00356 if ( present(lnd_present)) lnd_present=(SDAT%grid%data%rAttr(nfrac,1).lt.1.)
00357 else
00358
00359 if ( present(ocn_present) ) ocn_present=.false.
00360 if ( present(ice_present) ) ice_present=.false.
00361 if ( present(lnd_present) ) lnd_present=.true.
00362 end if
00363
00364 end subroutine shr_scam_checkSurface
00365
00366
00367
00368 end module shr_scam_mod
00369