SABER
model_gem.inc
Go to the documentation of this file.
1 !----------------------------------------------------------------------
2 ! Subroutine: model_gem_coord
3 !> Get GEM coordinates
4 !----------------------------------------------------------------------
5 subroutine model_gem_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=8),allocatable :: lon(:),lat(:),a(:),b(:)
18 character(len=1024),parameter :: subr = 'model_gem_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,'lon',nlon_id))
24 call mpl%ncerr(subr,nf90_inq_dimid(ncid,'lat',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,'lev',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))
36 allocate(b(model%nlev))
37 
38 ! Read data and close file
39 call mpl%ncerr(subr,nf90_inq_varid(ncid,'lon',lon_id))
40 call mpl%ncerr(subr,nf90_inq_varid(ncid,'lat',lat_id))
41 call mpl%ncerr(subr,nf90_inq_varid(ncid,'ap',a_id))
42 call mpl%ncerr(subr,nf90_inq_varid(ncid,'b',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(a(nam%levs(1:nam%nl))+b(nam%levs(1:nam%nl))*ps)
72  if (model%nl0>nam%nl) model%vunit(img,model%nl0) = log(ps)
73  else
74  model%vunit(img,:) = real(nam%levs(1:model%nl0),kind_real)
75  end if
76 end do
77 
78 ! Release memory
79 deallocate(lon)
80 deallocate(lat)
81 deallocate(a)
82 deallocate(b)
83 
84 end subroutine model_gem_coord
85 
86 !----------------------------------------------------------------------
87 ! Subroutine: model_gem_read
88 !> Read GEM field
89 !----------------------------------------------------------------------
90 subroutine model_gem_read(model,mpl,nam,filename,fld)
91 
92 implicit none
93 
94 ! Passed variables
95 class(model_type),intent(inout) :: model !< Model
96 type(mpl_type),intent(inout) :: mpl !< MPI data
97 type(nam_type),intent(in) :: nam !< Namelist
98 character(len=*),intent(in) :: filename !< File name
99 real(kind_real),intent(out) :: fld(model%nmga,model%nl0,nam%nv) !< Field
100 
101 ! Local variables
102 integer :: iv,il0,img,ilon,ilat,ndims
103 integer :: ncid,fld_id
104 integer,allocatable :: fld_tmp_int(:,:)
105 real(kind_real) :: add_offset,scale_factor
106 real(kind_real) :: fld_mg(model%nmg,model%nl0)
107 real(kind_real),allocatable :: fld_tmp(:,:,:)
108 character(len=1024),parameter :: subr = 'model_gem_read'
109 
110 if (mpl%main) then
111  ! Allocation
112  allocate(fld_tmp_int(model%nlon,model%nlat))
113  allocate(fld_tmp(model%nlon,model%nlat,model%nl0))
114 
115  ! Open file
116  call mpl%ncerr(subr,nf90_open(trim(nam%datadir)//'/'//trim(filename)//'.nc',nf90_nowrite,ncid))
117 end if
118 
119 do iv=1,nam%nv
120  if (mpl%main) then
121  ! Get variable id
122  call mpl%ncerr(subr,nf90_inq_varid(ncid,nam%variables(iv),fld_id))
123 
124  ! Check variable type and read data
125  call mpl%ncerr(subr,nf90_inquire_variable(ncid,fld_id,ndims=ndims))
126  call mpl%ncerr(subr,nf90_get_att(ncid,fld_id,'add_offset',add_offset))
127  call mpl%ncerr(subr,nf90_get_att(ncid,fld_id,'scale_factor',scale_factor))
128  select case (ndims)
129  case (2)
130  ! 2D data
131  call mpl%ncerr(subr,nf90_get_var(ncid,fld_id,fld_tmp_int,(/1,1/),(/model%nlon,model%nlat/)))
132  fld_tmp = mpl%msv%valr
133  if (trim(nam%lev2d)=='first') then
134  fld_tmp(:,:,1) = add_offset+scale_factor*real(fld_tmp_int,kind_real)
135  elseif (trim(nam%lev2d)=='last') then
136  fld_tmp(:,:,model%nl0) = add_offset+scale_factor*real(fld_tmp_int,kind_real)
137  end if
138  case (3)
139  ! 3D data
140  do il0=1,nam%nl
141  call mpl%ncerr(subr,nf90_get_var(ncid,fld_id,fld_tmp_int,(/1,1,nam%levs(il0)/),(/model%nlon,model%nlat,1/)))
142  fld_tmp(:,:,il0) = add_offset+scale_factor*real(fld_tmp_int,kind_real)
143  end do
144  case default
145  call mpl%abort(subr,'wrong number of dimensions for variable '//trim(nam%variables(iv)))
146  end select
147 
148  ! Pack
149  do il0=1,model%nl0
150  do img=1,model%nmg
151  ilon = model%mg_to_lon(img)
152  ilat = model%mg_to_lat(img)
153  fld_mg(img,il0) = fld_tmp(ilon,ilat,il0)
154  end do
155  end do
156  end if
157  call mpl%glb_to_loc(model%nl0,model%nmga,model%nmg,model%mga_to_mg,fld_mg,fld(:,:,iv))
158 end do
159 
160 if (mpl%main) then
161  ! Close file
162  call mpl%ncerr(subr,nf90_close(ncid))
163 
164  ! Release memory
165  deallocate(fld_tmp_int)
166  deallocate(fld_tmp)
167 end if
168 
169 end subroutine model_gem_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