1 !----------------------------------------------------------------------
2 ! Subroutine: model_ifs_coord
4 !----------------------------------------------------------------------
5 subroutine model_ifs_coord(model,mpl,nam)
10 class(model_type),intent(inout) :: model !< Model
11 type(mpl_type),intent(inout) :: mpl !< MPI data
12 type(nam_type),intent(in) :: nam !< Namelist
15 integer :: img,ilon,ilat
16 integer :: ncid,nlon_id,nlat_id,nlev_id,lon_id,lat_id,pres_id
17 real(
kind_real),allocatable :: lon(:),lat(:),pres(:)
18 character(len=1024),parameter :: subr =
'model_ifs_coord'
20 ! Open file and get dimensions
22 call mpl%ncerr(subr,nf90_open(trim(nam%datadir)
23 call mpl%ncerr(subr,nf90_inq_dimid(ncid,
'longitude',nlon_id))
24 call mpl%ncerr(subr,nf90_inq_dimid(ncid,
'latitude',nlat_id))
25 call mpl%ncerr(subr,nf90_inquire_dimension(ncid,nlon_id,len=model%nlon))
26 call mpl%ncerr(subr,nf90_inquire_dimension(ncid,nlat_id,len=model%nlat))
27 model%nmg = model%nlon*model%nlat
28 call mpl%ncerr(subr,nf90_inq_dimid(ncid,
'level',nlev_id))
29 call mpl%ncerr(subr,nf90_inquire_dimension(ncid,nlev_id,len=model%nlev))
33 allocate(lon(model%nlon))
34 allocate(lat(model%nlat))
35 allocate(pres(model%nlev))
37 ! Read data and close file
38 call mpl%ncerr(subr,nf90_inq_varid(ncid,
'longitude',lon_id))
39 call mpl%ncerr(subr,nf90_inq_varid(ncid,
'latitude',lat_id))
40 call mpl%ncerr(subr,nf90_inq_varid(ncid,
'pf',pres_id))
41 call mpl%ncerr(subr,nf90_get_var(ncid,lon_id,lon))
42 call mpl%ncerr(subr,nf90_get_var(ncid,lat_id,lat))
43 call mpl%ncerr(subr,nf90_get_var(ncid,pres_id,pres))
44 call mpl%ncerr(subr,nf90_close(ncid))
55 model%mg_to_tile(img) = 1
56 model%mg_to_lon(img) = ilon
57 model%mg_to_lat(img) = ilat
58 model%lon(img) = lon(ilon)
59 model%lat(img) = lat(ilat)
68 model%vunit(img,1:nam%nl) = log(pres(nam%levs(1:nam%nl)))
69 if (model%nl0>nam%nl) model%vunit(img,model%nl0) = log(
ps)
71 model%vunit(img,:) = real(nam%levs(1:model%nl0),
kind_real)
78 end subroutine model_ifs_coord
80 !----------------------------------------------------------------------
81 ! Subroutine: model_ifs_read
83 !----------------------------------------------------------------------
84 subroutine model_ifs_read(model,mpl,nam,filename,fld)
89 class(model_type),intent(inout) :: model !< Model
90 type(mpl_type),intent(inout) :: mpl !< MPI data
91 type(nam_type),intent(in) :: nam !< Namelist
92 character(len=*),intent(in) :: filename !< File name
93 real(
kind_real),intent(out) :: fld(model%nmga,model%nl0,nam%nv) !< Field
96 integer :: iv,il0,img,ilon,ilat
97 integer :: ncid,fld_id
98 real(
kind_real) :: fld_mg(model%nmg,model%nl0)
99 real(
kind_real),allocatable :: fld_tmp(:,:,:)
100 character(len=1024),parameter :: subr =
'model_ifs_read'
104 allocate(fld_tmp(model%nlon,model%nlat,model%nl0))
107 call mpl%ncerr(subr,nf90_open(trim(nam%datadir)
113 call mpl%ncerr(subr,nf90_inq_varid(ncid,nam%variables(iv),fld_id))
117 call mpl%ncerr(subr,nf90_get_var(ncid,fld_id,fld_tmp(:,:,il0),(/1,1,nam%levs(il0),1/),(/model%nlon,model%nlat,1,1/)))
123 ilon = model%mg_to_lon(img)
124 ilat = model%mg_to_lat(img)
125 fld_mg(img,il0) = fld_tmp(ilon,ilat,il0)
129 call mpl%glb_to_loc(model%nl0,model%nmga,model%nmg,model%mga_to_mg,fld_mg,fld(:,:,iv))
134 call mpl%ncerr(subr,nf90_close(ncid))
140 end subroutine model_ifs_read