1 !----------------------------------------------------------------------
2 ! Subroutine: model_norcpm_coord
3 !> Get NorCPM coordinates
4 !----------------------------------------------------------------------
5 subroutine model_norcpm_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,il0,ilon,ilat
16 integer :: ncid,nlon_id,nlat_id,nlev_id,lon_id,lat_id,area_id,gmask_id,dp_id
17 integer,allocatable :: gmask(:,:)
18 real(kind=8),allocatable :: lon(:,:),lat(:,:),area(:,:),dp(:,:,:)
19 character(len=1024),parameter :: subr =
'model_norcpm_coord'
21 ! Open file and get dimensions
23 call mpl%ncerr(subr,nf90_open(trim(nam%datadir)
24 call mpl%ncerr(subr,nf90_inq_dimid(ncid,
'x',nlon_id))
25 call mpl%ncerr(subr,nf90_inq_dimid(ncid,
'y',nlat_id))
26 call mpl%ncerr(subr,nf90_inquire_dimension(ncid,nlon_id,len=model%nlon))
27 call mpl%ncerr(subr,nf90_inquire_dimension(ncid,nlat_id,len=model%nlat))
28 model%nmg = model%nlon*model%nlat
29 call mpl%ncerr(subr,nf90_inq_dimid(ncid,
'kk2',nlev_id))
30 call mpl%ncerr(subr,nf90_inquire_dimension(ncid,nlev_id,len=model%nlev))
31 model%nlev = model%nlev/2
35 allocate(lon(model%nlon,model%nlat))
36 allocate(lat(model%nlon,model%nlat))
37 allocate(area(model%nlon,model%nlat))
38 allocate(gmask(model%nlon,model%nlat))
39 allocate(dp(model%nlon,model%nlat,model%nlev))
41 ! Read data and close file
42 call mpl%ncerr(subr,nf90_inq_varid(ncid,
'plon',lon_id))
43 call mpl%ncerr(subr,nf90_inq_varid(ncid,
'plat',lat_id))
44 call mpl%ncerr(subr,nf90_inq_varid(ncid,
'parea',area_id))
45 call mpl%ncerr(subr,nf90_inq_varid(ncid,
'pmask',gmask_id))
46 call mpl%ncerr(subr,nf90_inq_varid(ncid,
'dp',dp_id))
47 call mpl%ncerr(subr,nf90_get_var(ncid,lon_id,lon))
48 call mpl%ncerr(subr,nf90_get_var(ncid,lat_id,lat))
49 call mpl%ncerr(subr,nf90_get_var(ncid,area_id,area))
50 call mpl%ncerr(subr,nf90_get_var(ncid,gmask_id,gmask))
51 call mpl%ncerr(subr,nf90_get_var(ncid,dp_id,dp,(/1,1,1/),(/1,1,model%nlev/)))
52 call mpl%ncerr(subr,nf90_close(ncid))
59 model%mg_to_tile(img) = 1
60 model%mg_to_lon(img) = ilon
61 model%mg_to_lat(img) = ilat
65 model%mask(img,:) = (gmask(ilon,ilat)>0)
71 ilon = model%mg_to_lon(img)
72 ilat = model%mg_to_lat(img)
75 if (nam%levs(il0)==1) then
76 model%vunit(img,il0) =
ps+0.5*dp(ilon,ilat,1)
78 model%vunit(img,il0) =
ps+sum(dp(ilon,ilat,1:nam%levs(il0)-1))+0.5*dp(ilon,ilat,nam%levs(il0))
90 end subroutine model_norcpm_coord
92 !----------------------------------------------------------------------
93 ! Subroutine: model_norcpm_read
95 !----------------------------------------------------------------------
96 subroutine model_norcpm_read(model,mpl,nam,filename,fld)
101 class(model_type),intent(inout) :: model !< Model
102 type(mpl_type),intent(inout) :: mpl !< MPI data
103 type(nam_type),intent(in) :: nam !< Namelist
104 character(len=*),intent(in) :: filename !< File name
105 real(
kind_real),intent(out) :: fld(model%nmga,model%nl0,nam%nv) !< Field
108 integer :: iv,il0,img,ilon,ilat,ndims
109 integer :: ncid,fld_id
110 real(
kind_real) :: fld_mg(model%nmg,model%nl0)
111 real(
kind_real),allocatable :: fld_tmp(:,:,:)
112 character(len=1024),parameter :: subr =
'model_norcpm_read'
116 allocate(fld_tmp(model%nlon,model%nlat,model%nl0))
119 call mpl%ncerr(subr,nf90_open(trim(nam%datadir)
125 call mpl%ncerr(subr,nf90_inq_varid(ncid,nam%variables(iv),fld_id))
128 call mpl%ncerr(subr,nf90_inquire_variable(ncid,fld_id,ndims=ndims))
134 fld_tmp = mpl%msv%valr
135 if (trim(nam%lev2d)==
'first') then
136 call mpl%ncerr(subr,nf90_get_var(ncid,fld_id,fld_tmp(:,:,1),(/1,1,1/), &
137 & (/model%nlon,model%nlat,1/)))
138 elseif (trim(nam%lev2d)==
'last') then
139 call mpl%ncerr(subr,nf90_get_var(ncid,fld_id,fld_tmp(:,:,model%nl0),(/1,1,1/), &
140 & (/model%nlon,model%nlat,1/)))
145 call mpl%ncerr(subr,nf90_get_var(ncid,fld_id,fld_tmp(:,:,il0),(/1,1,nam%levs(il0),1/), &
146 & (/model%nlon,model%nlat,1,1/)))
149 call mpl%abort(subr,
'wrong number of dimensions for variable '
155 ilon = model%mg_to_lon(img)
156 ilat = model%mg_to_lat(img)
157 fld_mg(img,il0) = fld_tmp(ilon,ilat,il0)
161 call mpl%glb_to_loc(model%nl0,model%nmga,model%nmg,model%mga_to_mg,fld_mg,fld(:,:,iv))
166 call mpl%ncerr(subr,nf90_close(ncid))
172 end subroutine model_norcpm_read