1 !----------------------------------------------------------------------
2 ! Subroutine: model_mpas_coord
3 !> Get MPAS coordinates
4 !----------------------------------------------------------------------
5 subroutine model_mpas_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
16 integer :: ncid,nmg_id,nlev_id,lon_id,lat_id,pres_id
17 real(kind=4),allocatable :: lon(:),lat(:),pres(:)
18 character(len=1024),parameter :: subr =
'model_mpas_coord'
20 ! Open file and get dimensions
22 model%nlon = mpl%msv%vali
23 model%nlat = mpl%msv%vali
24 call mpl%ncerr(subr,nf90_open(trim(nam%datadir)
25 call mpl%ncerr(subr,nf90_inq_dimid(ncid,
'nCells',nmg_id))
26 call mpl%ncerr(subr,nf90_inquire_dimension(ncid,nmg_id,len=model%nmg))
27 call mpl%ncerr(subr,nf90_inq_dimid(ncid,
'nVertLevels',nlev_id))
28 call mpl%ncerr(subr,nf90_inquire_dimension(ncid,nlev_id,len=model%nlev))
32 allocate(lon(model%nmg))
33 allocate(lat(model%nmg))
34 allocate(pres(model%nlev))
36 ! Read data and close file
37 call mpl%ncerr(subr,nf90_inq_varid(ncid,
'lonCell',lon_id))
38 call mpl%ncerr(subr,nf90_inq_varid(ncid,
'latCell',lat_id))
39 call mpl%ncerr(subr,nf90_inq_varid(ncid,
'pressure_base',pres_id))
40 call mpl%ncerr(subr,nf90_get_var(ncid,lon_id,lon))
41 call mpl%ncerr(subr,nf90_get_var(ncid,lat_id,lat))
42 call mpl%ncerr(subr,nf90_get_var(ncid,pres_id,pres))
43 call mpl%ncerr(subr,nf90_close(ncid))
55 model%vunit(img,1:nam%nl) = log(pres(nam%levs(1:nam%nl)))
56 if (model%nl0>nam%nl) model%vunit(img,model%nl0) = log(
ps)
58 model%vunit(img,:) = real(nam%levs(1:model%nl0),
kind_real)
67 end subroutine model_mpas_coord
69 !----------------------------------------------------------------------
70 ! Subroutine: model_mpas_read
72 !----------------------------------------------------------------------
73 subroutine model_mpas_read(model,mpl,nam,filename,fld)
78 class(model_type),intent(inout) :: model !< Model
79 type(mpl_type),intent(inout) :: mpl !< MPI data
80 type(nam_type),intent(in) :: nam !< Namelist
81 character(len=*),intent(in) :: filename !< File name
82 real(
kind_real),intent(out) :: fld(model%nmga,model%nl0,nam%nv) !< Field
85 integer :: iv,il0,ndims
86 integer :: ncid,fld_id
87 real(
kind_real) :: fld_mg(model%nmg,model%nl0)
88 character(len=1024),parameter :: subr =
'model_mpas_read'
91 if (mpl%
main) call mpl%ncerr(subr,nf90_open(trim(nam%datadir)
96 call mpl%ncerr(subr,nf90_inq_varid(ncid,nam%variables(iv),fld_id))
99 call mpl%ncerr(subr,nf90_inquire_variable(ncid,fld_id,ndims=ndims))
105 fld_mg = mpl%msv%valr
106 if (trim(nam%lev2d)==
'first') then
107 call mpl%ncerr(subr,nf90_get_var(ncid,fld_id,fld_mg(:,1),(/1,1/),(/model%nmg,1/)))
108 elseif (trim(nam%lev2d)==
'last') then
109 call mpl%ncerr(subr,nf90_get_var(ncid,fld_id,fld_mg(:,model%nl0),(/1,1/),(/model%nmg,1/)))
114 call mpl%ncerr(subr,nf90_get_var(ncid,fld_id,fld_mg(:,il0),(/nam%levs(il0),1,1/),(/1,model%nmg,1/)))
117 call mpl%abort(subr,
'wrong number of dimensions for variable '
120 call mpl%glb_to_loc(model%nl0,model%nmga,model%nmg,model%mga_to_mg,fld_mg,fld(:,:,iv))
124 if (mpl%
main) call mpl%ncerr(subr,nf90_close(ncid))
126 end subroutine model_mpas_read