SABER
model_gfs.inc
Go to the documentation of this file.
1 !----------------------------------------------------------------------
2 ! Subroutine: model_gfs_coord
3 !> Get GFS coordinates
4 !----------------------------------------------------------------------
5 subroutine model_gfs_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,ilon,ilat
16 integer :: ncid,nlon_id,nlat_id,nlev_id,lon_id,lat_id,a_id,b_id
17 real(kind_real),allocatable :: lon(:),lat(:),a(:),b(:)
18 character(len=1024),parameter :: subr = 'model_gfs_coord'
19 
20 ! Open file and get dimensions
21 model%ntile = 1
22 call mpl%ncerr(subr,nf90_open(trim(nam%datadir)//'/grid.nc',nf90_share,ncid))
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))
30 
31 ! Allocation
32 call model%alloc
33 allocate(lon(model%nlon))
34 allocate(lat(model%nlat))
35 allocate(a(model%nlev+1))
36 allocate(b(model%nlev+1))
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,'ak',a_id))
42 call mpl%ncerr(subr,nf90_inq_varid(ncid,'bk',b_id))
43 call mpl%ncerr(subr,nf90_get_var(ncid,lon_id,lon))
44 call mpl%ncerr(subr,nf90_get_var(ncid,lat_id,lat))
45 call mpl%ncerr(subr,nf90_get_var(ncid,a_id,a))
46 call mpl%ncerr(subr,nf90_get_var(ncid,b_id,b))
47 call mpl%ncerr(subr,nf90_close(ncid))
48 
49 ! Convert to radian
50 lon = lon*deg2rad
51 lat = lat*deg2rad
52 
53 ! Model grid
54 img = 0
55 do ilon=1,model%nlon
56  do ilat=1,model%nlat
57  img = img+1
58  model%mg_to_tile(img) = 1
59  model%mg_to_lon(img) = ilon
60  model%mg_to_lat(img) = ilat
61  model%lon(img) = lon(ilon)
62  model%lat(img) = lat(ilat)
63  end do
64 end do
65 model%area = 4.0*pi/real(model%nmg,kind_real)
66 model%mask = .true.
67 
68 ! Vertical unit
69 do img=1,model%nmg
70  if (nam%logpres) then
71  model%vunit(img,1:nam%nl) = log(0.5*(a(nam%levs(1:nam%nl))+a(nam%levs(1:nam%nl)+1)) &
72  & +0.5*(b(nam%levs(1:nam%nl))+b(nam%levs(1:nam%nl)+1))*ps)
73  if (model%nl0>nam%nl) model%vunit(img,model%nl0) = log(ps)
74  else
75  model%vunit(img,:) = real(nam%levs(1:model%nl0),kind_real)
76  end if
77 end do
78 
79 ! Release memory
80 deallocate(lon)
81 deallocate(lat)
82 deallocate(a)
83 deallocate(b)
84 
85 end subroutine model_gfs_coord
86 
87 !----------------------------------------------------------------------
88 ! Subroutine: model_gfs_read
89 !> Read GFS field
90 !----------------------------------------------------------------------
91 subroutine model_gfs_read(model,mpl,nam,filename,fld)
92 
93 implicit none
94 
95 ! Passed variables
96 class(model_type),intent(inout) :: model !< Model
97 type(mpl_type),intent(inout) :: mpl !< MPI data
98 type(nam_type),intent(in) :: nam !< Namelist
99 character(len=*),intent(in) :: filename !< File name
100 real(kind_real),intent(out) :: fld(model%nmga,model%nl0,nam%nv) !< Field
101 
102 ! Local variables
103 integer :: iv,il0,img,ilon,ilat,ndims
104 integer :: ncid,fld_id
105 real(kind_real) :: fld_mg(model%nmg,model%nl0)
106 real(kind_real),allocatable :: fld_tmp(:,:,:)
107 character(len=1024),parameter :: subr = 'model_gfs_read'
108 
109 if (mpl%main) then
110  ! Allocation
111  allocate(fld_tmp(model%nlon,model%nlat,model%nl0))
112 
113  ! Open file
114  call mpl%ncerr(subr,nf90_open(trim(nam%datadir)//'/'//trim(filename)//'.nc',nf90_nowrite,ncid))
115 end if
116 
117 do iv=1,nam%nv
118  if (mpl%main) then
119  ! Get variable id
120  call mpl%ncerr(subr,nf90_inq_varid(ncid,nam%variables(iv),fld_id))
121 
122 
123  ! Read data
124  select case (ndims)
125  case (2)
126  ! 2D data
127  fld_tmp = mpl%msv%valr
128  if (trim(nam%lev2d)=='first') then
129  call mpl%ncerr(subr,nf90_get_var(ncid,fld_id,fld_tmp(:,:,1),(/1,1/),(/model%nlon,model%nlat/)))
130  elseif (trim(nam%lev2d)=='last') then
131  call mpl%ncerr(subr,nf90_get_var(ncid,fld_id,fld_tmp(:,:,model%nl0),(/1,1/),(/model%nlon,model%nlat/)))
132  end if
133  case (3)
134  ! 3D data
135  do il0=1,nam%nl
136  call mpl%ncerr(subr,nf90_get_var(ncid,fld_id,fld_tmp(:,:,il0),(/1,1,nam%levs(il0)/),(/model%nlon,model%nlat,1/)))
137  end do
138  case default
139  call mpl%abort(subr,'wrong number of dimensions for variable '//trim(nam%variables(iv)))
140  end select
141 
142  ! Pack
143  do il0=1,model%nl0
144  do img=1,model%nmg
145  ilon = model%mg_to_lon(img)
146  ilat = model%mg_to_lat(img)
147  fld_mg(img,il0) = fld_tmp(ilon,ilat,il0)
148  end do
149  end do
150  end if
151  call mpl%glb_to_loc(model%nl0,model%nmga,model%nmg,model%mga_to_mg,fld_mg,fld(:,:,iv))
152 end do
153 
154 if (mpl%main) then
155  ! Close file
156  call mpl%ncerr(subr,nf90_close(ncid))
157 
158  ! Release memory
159  deallocate(fld_tmp)
160 end if
161 
162 end subroutine model_gfs_read
main
int main(int argc, char **argv)
Definition: bump_main.cc:17
tools_const::ps
real(kind_real), parameter, public ps
Definition: tools_const.F90:19
tools_const::deg2rad
real(kind_real), parameter, public deg2rad
Definition: tools_const.F90:15
tools_const::pi
real(kind_real), parameter, public pi
Definition: tools_const.F90:14
tools_kinds::kind_real
integer, parameter, public kind_real
Definition: tools_kinds.F90:18
type_rng::a
integer(kind=int64), parameter a
Definition: type_rng.F90:18