1 !----------------------------------------------------------------------
2 ! Subroutine: model_nemo_coord
3 !> Get NEMO coordinates
4 !----------------------------------------------------------------------
5 subroutine model_nemo_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
15 integer :: il0,img,ilat,ilon
16 integer :: ncid,nlon_id,nlat_id,nlev_id,lon_id,lat_id,tmask_id,e1t_id,e2t_id,e3t_id
17 integer(kind=1),allocatable :: tmask(:,:,:)
18 real(
kind_real),allocatable :: lon(:,:),lat(:,:),e1t(:,:),e2t(:,:),e3t(:,:,:)
19 character(len=1024),parameter :: subr =
'model_nemo_coord'
21 ! Open file and get dimensions
23 call mpl%ncerr(subr,nf90_open(trim(nam%datadir)
24 call mpl%ncerr(subr,nf90_inq_dimid(ncid,
'x',nlon_id))
25 call mpl%ncerr(subr,nf90_inq_dimid(ncid,
'y',nlat_id))
26 call mpl%ncerr(subr,nf90_inquire_dimension(ncid,nlon_id,len=model%nlon))
27 call mpl%ncerr(subr,nf90_inquire_dimension(ncid,nlat_id,len=model%nlat))
28 model%nmg = model%nlon*model%nlat
29 call mpl%ncerr(subr,nf90_inq_dimid(ncid,
'z',nlev_id))
30 call mpl%ncerr(subr,nf90_inquire_dimension(ncid,nlev_id,len=model%nlev))
34 allocate(lon(model%nlon,model%nlat))
35 allocate(lat(model%nlon,model%nlat))
36 allocate(tmask(model%nlon,model%nlat,model%nl0))
37 allocate(e1t(model%nlon,model%nlat))
38 allocate(e2t(model%nlon,model%nlat))
39 allocate(e3t(model%nlon,model%nlat,model%nlev))
41 ! Read data and close file
42 call mpl%ncerr(subr,nf90_inq_varid(ncid,
'nav_lon',lon_id))
43 call mpl%ncerr(subr,nf90_inq_varid(ncid,
'nav_lat',lat_id))
44 call mpl%ncerr(subr,nf90_inq_varid(ncid,
'tmask',tmask_id))
45 call mpl%ncerr(subr,nf90_inq_varid(ncid,
'e1t',e1t_id))
46 call mpl%ncerr(subr,nf90_inq_varid(ncid,
'e2t',e2t_id))
47 call mpl%ncerr(subr,nf90_inq_varid(ncid,
'e3t',e3t_id))
48 call mpl%ncerr(subr,nf90_get_var(ncid,lon_id,lon))
49 call mpl%ncerr(subr,nf90_get_var(ncid,lat_id,lat))
51 call mpl%ncerr(subr,nf90_get_var(ncid,tmask_id,tmask(:,:,il0),(/1,1,nam%levs(il0),1/),(/model%nlon,model%nlat,1,1/)))
53 call mpl%ncerr(subr,nf90_get_var(ncid,e1t_id,e1t,(/1,1,1/),(/model%nlon,model%nlat,1/)))
54 call mpl%ncerr(subr,nf90_get_var(ncid,e2t_id,e2t,(/1,1,1/),(/model%nlon,model%nlat,1/)))
56 call mpl%ncerr(subr,nf90_get_var(ncid,e3t_id,e3t(:,:,il0),(/1,1,il0,1/),(/model%nlon,model%nlat,1,1/)))
58 call mpl%ncerr(subr,nf90_close(ncid))
69 model%mg_to_tile(img) = 1
70 model%mg_to_lon(img) = ilon
71 model%mg_to_lat(img) = ilat
72 model%lon(img) = lon(ilon,ilat)
73 model%lat(img) = lat(ilon,ilat)
74 model%area(img) = e1t(ilon,ilat)*e2t(ilon,ilat)/
req**2
76 model%mask(img,il0) = (tmask(ilon,ilat,il0)>0)
84 ilon = model%mg_to_lon(img)
85 ilat = model%mg_to_lat(img)
87 if (nam%levs(il0)==1) then
88 model%vunit(img,il0) = -0.5*e3t(ilon,ilat,1)
90 model%vunit(img,il0) = -sum(e3t(ilon,ilat,1:nam%levs(il0)-1))-0.5*e3t(ilon,ilat,nam%levs(il0))
93 if (model%nl0>nam%nl) model%vunit(img,model%nl0) = 0.0
95 model%vunit(img,:) = real(nam%levs(1:model%nl0),
kind_real)
104 end subroutine model_nemo_coord
106 !----------------------------------------------------------------------
107 ! Subroutine: model_nemo_read
109 !----------------------------------------------------------------------
110 subroutine model_nemo_read(model,mpl,nam,filename,fld)
115 class(model_type),intent(inout) :: model !< Model
116 type(mpl_type),intent(inout) :: mpl !< MPI data
117 type(nam_type),intent(in) :: nam !< Namelist
118 character(len=*),intent(in) :: filename !< File name
119 real(
kind_real),intent(out) :: fld(model%nmga,model%nl0,nam%nv) !< Field
122 integer :: iv,il0,img,ilon,ilat,ndims
123 integer :: ncid,fld_id
124 real(
kind_real) :: fld_tmp2,fld_mg(model%nmg,model%nl0)
125 real(
kind_real),allocatable :: fld_tmp(:,:,:)
126 character(len=1024),parameter :: subr =
'model_nemo_read'
130 allocate(fld_tmp(model%nlon,model%nlat,model%nl0))
133 call mpl%ncerr(subr,nf90_open(trim(nam%datadir)
139 call mpl%ncerr(subr,nf90_inq_varid(ncid,nam%variables(iv),fld_id))
142 call mpl%ncerr(subr,nf90_inquire_variable(ncid,fld_id,ndims=ndims))
148 fld_tmp = mpl%msv%valr
149 if (trim(nam%lev2d)==
'first') then
150 call mpl%ncerr(subr,nf90_get_var(ncid,fld_id,fld_tmp(:,:,1),(/1,1,1/),(/model%nlon,model%nlat,1/)))
151 elseif (trim(nam%lev2d)==
'last') then
152 call mpl%ncerr(subr,nf90_get_var(ncid,fld_id,fld_tmp(:,:,model%nl0),(/1,1,1/),(/model%nlon,model%nlat,1/)))
157 call mpl%ncerr(subr,nf90_get_var(ncid,fld_id,fld_tmp(:,:,il0),(/1,1,nam%levs(il0),1/),(/model%nlon,model%nlat,1,1/)))
158 select case (trim(nam%variables(iv)))
163 call mpl%ncerr(subr,nf90_get_var(ncid,fld_id,fld_tmp2,(/model%nlon,ilat,nam%levs(il0),1/)))
165 call mpl%ncerr(subr,nf90_get_var(ncid,fld_id,fld_tmp2,(/ilon-1,ilat,nam%levs(il0),1/)))
167 fld_tmp(ilon,ilat,il0) = 0.5*(fld_tmp(ilon,ilat,il0)+fld_tmp2)
174 call mpl%ncerr(subr,nf90_get_var(ncid,fld_id,fld_tmp2,(/ilon,model%nlat,nam%levs(il0),1/)))
176 call mpl%ncerr(subr,nf90_get_var(ncid,fld_id,fld_tmp2,(/ilon,ilat-1,nam%levs(il0),1/)))
178 fld_tmp(ilon,ilat,il0) = 0.5*(fld_tmp(ilon,ilat,il0)+fld_tmp2)
184 call mpl%abort(subr,
'wrong number of dimensions for variable '
190 ilon = model%mg_to_lon(img)
191 ilat = model%mg_to_lat(img)
192 fld_mg(img,il0) = fld_tmp(ilon,ilat,il0)
196 call mpl%glb_to_loc(model%nl0,model%nmga,model%nmg,model%mga_to_mg,fld_mg,fld(:,:,iv))
201 call mpl%ncerr(subr,nf90_close(ncid))
207 end subroutine model_nemo_read