1 !----------------------------------------------------------------------
2 ! Subroutine: model_res_coord
4 !----------------------------------------------------------------------
5 subroutine model_res_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,alt_id,area_id,mask_id
17 real(
kind_real),allocatable :: lon(:,:),lat(:,:),alt(:,:),area(:,:),mask(:,:)
18 character(len=1024),parameter :: subr =
'model_res_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,
'nhcells',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,
'nlev',nlev_id))
28 call mpl%ncerr(subr,nf90_inquire_dimension(ncid,nlev_id,len=model%nlev))
32 allocate(lon(model%nmg,model%nlev))
33 allocate(lat(model%nmg,model%nlev))
34 allocate(alt(model%nmg,model%nlev))
35 allocate(area(model%nmg,model%nlev))
36 allocate(mask(model%nmg,model%nlev))
38 ! Read data and close file
39 call mpl%ncerr(subr,nf90_inq_varid(ncid,
'longitude',lon_id))
40 call mpl%ncerr(subr,nf90_inq_varid(ncid,
'latitude',lat_id))
41 call mpl%ncerr(subr,nf90_inq_varid(ncid,
'altitude',alt_id))
42 call mpl%ncerr(subr,nf90_inq_varid(ncid,
'horzArea',area_id))
43 call mpl%ncerr(subr,nf90_inq_varid(ncid,
'activeCells',mask_id))
44 call mpl%ncerr(subr,nf90_get_var(ncid,lon_id,lon))
45 call mpl%ncerr(subr,nf90_get_var(ncid,lat_id,lat))
46 call mpl%ncerr(subr,nf90_get_var(ncid,alt_id,alt))
47 call mpl%ncerr(subr,nf90_get_var(ncid,area_id,area))
48 call mpl%ncerr(subr,nf90_get_var(ncid,mask_id,mask))
49 call mpl%ncerr(subr,nf90_close(ncid))
53 model%lon = (sum(lon,dim=2)/real(model%nlev,
kind_real))/
req
54 model%lat = (sum(lat,dim=2)/real(model%nlev,
kind_real))/
req
55 model%area = sum(area,dim=2)/real(model%nlev,
kind_real)/
req**2
56 model%mask = (mask(:,nam%levs(1:model%nl0))>0.5)
59 if (nam%logpres) call mpl%abort(subr,
'pressure logarithm vertical coordinate is not available for this model')
62 model%vunit(img,il0) = alt(img,nam%levs(il0))
73 end subroutine model_res_coord
75 !----------------------------------------------------------------------
76 ! Subroutine: model_res_read
78 !----------------------------------------------------------------------
79 subroutine model_res_read(model,mpl,nam,filename,fld)
84 class(model_type),intent(inout) :: model !< Model
85 type(mpl_type),intent(inout) :: mpl !< MPI data
86 type(nam_type),intent(in) :: nam !< Namelist
87 character(len=*),intent(in) :: filename !< File name
88 real(
kind_real),intent(out) :: fld(model%nmga,model%nl0,nam%nv) !< Field
92 integer :: ncid,fld_id
93 real(
kind_real) :: fld_mg(model%nmg,model%nl0)
94 character(len=1024),parameter :: subr =
'model_res_read'
97 if (mpl%
main) call mpl%ncerr(subr,nf90_open(trim(nam%datadir)
102 call mpl%ncerr(subr,nf90_inq_varid(ncid,nam%variables(iv),fld_id))
106 call mpl%ncerr(subr,nf90_get_var(ncid,fld_id,fld_mg(:,il0),(/1,nam%levs(il0)/),(/model%nmg,1/)))
109 call mpl%glb_to_loc(model%nl0,model%nmga,model%nmg,model%mga_to_mg,fld_mg,fld(:,:,iv))
113 if (mpl%
main) call mpl%ncerr(subr,nf90_close(ncid))
115 end subroutine model_res_read