SABER
model_res.inc
Go to the documentation of this file.
1 !----------------------------------------------------------------------
2 ! Subroutine: model_res_coord
3 !> Get RES coordinates
4 !----------------------------------------------------------------------
5 subroutine model_res_coord(model,mpl,nam)
6 
7 implicit none
8 
9 ! Passed variables
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
13 
14 ! Local variables
15 integer :: img,il0
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'
19 
20 ! Open file and get dimensions
21 model%ntile = 1
22 model%nlon = mpl%msv%vali
23 model%nlat = mpl%msv%vali
24 call mpl%ncerr(subr,nf90_open(trim(nam%datadir)//'/grid.nc',nf90_share,ncid))
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))
29 
30 ! Allocation
31 call model%alloc
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))
37 
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))
50 
51 ! Model grid
52 model%mg_to_tile = 1
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)
57 
58 ! Vertical unit
59 if (nam%logpres) call mpl%abort(subr,'pressure logarithm vertical coordinate is not available for this model')
60 do il0=1,model%nl0
61  do img=1,model%nmg
62  model%vunit(img,il0) = alt(img,nam%levs(il0))
63  end do
64 end do
65 
66 ! Release memory
67 deallocate(lon)
68 deallocate(lat)
69 deallocate(alt)
70 deallocate(area)
71 deallocate(mask)
72 
73 end subroutine model_res_coord
74 
75 !----------------------------------------------------------------------
76 ! Subroutine: model_res_read
77 !> Read RES field
78 !----------------------------------------------------------------------
79 subroutine model_res_read(model,mpl,nam,filename,fld)
80 
81 implicit none
82 
83 ! Passed variables
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
89 
90 ! Local variables
91 integer :: iv,il0
92 integer :: ncid,fld_id
93 real(kind_real) :: fld_mg(model%nmg,model%nl0)
94 character(len=1024),parameter :: subr = 'model_res_read'
95 
96 ! Open file
97 if (mpl%main) call mpl%ncerr(subr,nf90_open(trim(nam%datadir)//'/'//trim(filename)//'.nc',nf90_nowrite,ncid))
98 
99 do iv=1,nam%nv
100  if (mpl%main) then
101  ! Get variable id
102  call mpl%ncerr(subr,nf90_inq_varid(ncid,nam%variables(iv),fld_id))
103 
104  ! Read data
105  do il0=1,nam%nl
106  call mpl%ncerr(subr,nf90_get_var(ncid,fld_id,fld_mg(:,il0),(/1,nam%levs(il0)/),(/model%nmg,1/)))
107  end do
108  end if
109  call mpl%glb_to_loc(model%nl0,model%nmga,model%nmg,model%mga_to_mg,fld_mg,fld(:,:,iv))
110 end do
111 
112 ! Close file
113 if (mpl%main) call mpl%ncerr(subr,nf90_close(ncid))
114 
115 end subroutine model_res_read
main
int main(int argc, char **argv)
Definition: bump_main.cc:17
tools_const::req
real(kind_real), parameter, public req
Definition: tools_const.F90:17
tools_kinds::kind_real
integer, parameter, public kind_real
Definition: tools_kinds.F90:18