SICOPOLIS V3.1
 All Classes Files Functions Variables Macros
topography3_nc.F90
Go to the documentation of this file.
1 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 !
3 ! Subroutine : t o p o g r a p h y 3 _ n c
4 !
5 !> @file
6 !!
7 !! Definition of the initial surface and bedrock topography
8 !! (including gradients) and of the horizontal grid spacings dxi, deta.
9 !! For initial topography from previous simulation (read in NetCDF format).
10 !!
11 !! @section Copyright
12 !!
13 !! Copyright 2009-2013 Ralf Greve
14 !!
15 !! @section License
16 !!
17 !! This file is part of SICOPOLIS.
18 !!
19 !! SICOPOLIS is free software: you can redistribute it and/or modify
20 !! it under the terms of the GNU General Public License as published by
21 !! the Free Software Foundation, either version 3 of the License, or
22 !! (at your option) any later version.
23 !!
24 !! SICOPOLIS is distributed in the hope that it will be useful,
25 !! but WITHOUT ANY WARRANTY; without even the implied warranty of
26 !! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
27 !! GNU General Public License for more details.
28 !!
29 !! You should have received a copy of the GNU General Public License
30 !! along with SICOPOLIS. If not, see <http://www.gnu.org/licenses/>.
31 !<
32 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
33 
34 !-------------------------------------------------------------------------------
35 !> Definition of the initial surface and bedrock topography
36 !! (including gradients) and of the horizontal grid spacings dxi, deta.
37 !! For initial topography from previous simulation (read in NetCDF format).
38 !<------------------------------------------------------------------------------
39 subroutine topography3_nc(dxi, deta, z_sl, anfdatname)
40 
41 use sico_types
43 use netcdf
44 use nc_check
45 
46 implicit none
47 
48 character(len=100), intent(in) :: anfdatname
49 real(dp), intent(out) :: dxi, deta, z_sl
50 
51 integer(i4b) :: i, j, kc, kt, kr, n
52 integer(i4b) :: ios
53 character :: ch_dummy
54 
55 integer(i4b) :: ncid, ncv
56 ! ncid: ID of the output file
57 ! ncv: Variable ID
58 
59 integer(i2b), dimension(0:IMAX,0:JMAX) :: maske_conv, n_cts_conv
60 real(sp) :: time_conv, dummy_conv, z_sl_conv, &
61  v_tot_conv, a_grounded_conv, a_floating_conv, &
62  h_r_conv, &
63  xi_conv(0:imax), eta_conv(0:jmax), &
64  sigma_level_c_conv(0:kcmax), sigma_level_t_conv(0:ktmax), &
65  sigma_level_r_conv(0:krmax)
66 real(sp), dimension(0:IMAX,0:JMAX) :: lambda_conv, phi_conv, &
67  lon_conv, lat_conv, &
68  temp_s_conv, as_perp_conv, &
69  zs_conv, zm_conv, zb_conv, zl_conv, h_c_conv, h_t_conv, h_conv, &
70  q_bm_conv, q_tld_conv, &
71  am_perp_conv, &
72  qx_conv, qy_conv, &
73  dzs_dtau_conv, dzm_dtau_conv, dzb_dtau_conv, dzl_dtau_conv, &
74  dh_c_dtau_conv, dh_t_dtau_conv, dh_dtau_conv, &
75  vx_b_g_conv, vy_b_g_conv, vz_b_conv, vh_b_conv, &
76  vx_s_g_conv, vy_s_g_conv, vz_s_conv, vh_s_conv, &
77  temp_b_conv, temph_b_conv, &
78  p_b_w_conv, h_w_conv, q_gl_g_conv
79 real(sp), dimension(0:IMAX,0:JMAX,0:KCMAX) :: vx_c_conv, vy_c_conv, vz_c_conv, &
80  temp_c_conv, age_c_conv
81 real(sp), dimension(0:IMAX,0:JMAX,0:KTMAX) :: vx_t_conv, vy_t_conv, vz_t_conv, &
82  omega_t_conv, age_t_conv
83 real(sp), dimension(0:IMAX,0:JMAX,0:KRMAX) :: temp_r_conv
84 
85 !-------- Read data from time-slice file of previous simulation --------
86 
87 call check( nf90_open(anfdatpath//'/'//trim(anfdatname), nf90_nowrite, ncid) )
88 
89 call check( nf90_inq_varid(ncid, 'time', ncv) )
90 call check( nf90_get_var(ncid, ncv, time_conv) )
91 
92 if ( nf90_inq_varid(ncid, 'delta_ts', ncv) == nf90_noerr ) then
93  call check( nf90_get_var(ncid, ncv, dummy_conv) )
94 else if ( nf90_inq_varid(ncid, 'glac_index', ncv) == nf90_noerr ) then
95  call check( nf90_get_var(ncid, ncv, dummy_conv) )
96 else
97  continue
98 end if
99 
100 call check( nf90_inq_varid(ncid, 'z_sl', ncv) )
101 call check( nf90_get_var(ncid, ncv, z_sl_conv) )
102 
103 call check( nf90_inq_varid(ncid, 'xi', ncv) )
104 call check( nf90_get_var(ncid, ncv, xi_conv) )
105 
106 call check( nf90_inq_varid(ncid, 'eta', ncv) )
107 call check( nf90_get_var(ncid, ncv, eta_conv) )
108 
109 call check( nf90_inq_varid(ncid, 'sigma_level_c', ncv) )
110 call check( nf90_get_var(ncid, ncv, sigma_level_c_conv) )
111 
112 call check( nf90_inq_varid(ncid, 'sigma_level_t', ncv) )
113 call check( nf90_get_var(ncid, ncv, sigma_level_t_conv) )
114 
115 call check( nf90_inq_varid(ncid, 'sigma_level_r', ncv) )
116 call check( nf90_get_var(ncid, ncv, sigma_level_r_conv) )
117 
118 call check( nf90_inq_varid(ncid, 'lon', ncv) )
119 call check( nf90_get_var(ncid, ncv, lon_conv) )
120 
121 call check( nf90_inq_varid(ncid, 'lat', ncv) )
122 call check( nf90_get_var(ncid, ncv, lat_conv) )
123 
124 call check( nf90_inq_varid(ncid, 'lambda', ncv) )
125 call check( nf90_get_var(ncid, ncv, lambda_conv) )
126 
127 call check( nf90_inq_varid(ncid, 'phi', ncv) )
128 call check( nf90_get_var(ncid, ncv, phi_conv) )
129 
130 call check( nf90_inq_varid(ncid, 'temp_s', ncv) )
131 call check( nf90_get_var(ncid, ncv, temp_s_conv) )
132 
133 call check( nf90_inq_varid(ncid, 'as_perp', ncv) )
134 call check( nf90_get_var(ncid, ncv, as_perp_conv) )
135 
136 call check( nf90_inq_varid(ncid, 'maske', ncv) )
137 call check( nf90_get_var(ncid, ncv, maske_conv) )
138 
139 call check( nf90_inq_varid(ncid, 'n_cts', ncv) )
140 call check( nf90_get_var(ncid, ncv, n_cts_conv) )
141 
142 call check( nf90_inq_varid(ncid, 'zs', ncv) )
143 call check( nf90_get_var(ncid, ncv, zs_conv) )
144 
145 call check( nf90_inq_varid(ncid, 'zm', ncv) )
146 call check( nf90_get_var(ncid, ncv, zm_conv) )
147 
148 call check( nf90_inq_varid(ncid, 'zb', ncv) )
149 call check( nf90_get_var(ncid, ncv, zb_conv) )
150 
151 call check( nf90_inq_varid(ncid, 'zl', ncv) )
152 call check( nf90_get_var(ncid, ncv, zl_conv) )
153 
154 call check( nf90_inq_varid(ncid, 'H_c', ncv) )
155 call check( nf90_get_var(ncid, ncv, h_c_conv) )
156 
157 call check( nf90_inq_varid(ncid, 'H_t', ncv) )
158 call check( nf90_get_var(ncid, ncv, h_t_conv) )
159 
160 call check( nf90_inq_varid(ncid, 'H', ncv) )
161 call check( nf90_get_var(ncid, ncv, h_conv) )
162 
163 call check( nf90_inq_varid(ncid, 'H_R', ncv) )
164 call check( nf90_get_var(ncid, ncv, h_r_conv) )
165 
166 call check( nf90_inq_varid(ncid, 'vx_c', ncv) )
167 call check( nf90_get_var(ncid, ncv, vx_c_conv) )
168 
169 call check( nf90_inq_varid(ncid, 'vy_c', ncv) )
170 call check( nf90_get_var(ncid, ncv, vy_c_conv) )
171 
172 call check( nf90_inq_varid(ncid, 'vz_c', ncv) )
173 call check( nf90_get_var(ncid, ncv, vz_c_conv) )
174 
175 call check( nf90_inq_varid(ncid, 'vx_t', ncv) )
176 call check( nf90_get_var(ncid, ncv, vx_t_conv) )
177 
178 call check( nf90_inq_varid(ncid, 'vy_t', ncv) )
179 call check( nf90_get_var(ncid, ncv, vy_t_conv) )
180 
181 call check( nf90_inq_varid(ncid, 'vz_t', ncv) )
182 call check( nf90_get_var(ncid, ncv, vz_t_conv) )
183 
184 call check( nf90_inq_varid(ncid, 'temp_c', ncv) )
185 call check( nf90_get_var(ncid, ncv, temp_c_conv) )
186 
187 call check( nf90_inq_varid(ncid, 'omega_t', ncv) )
188 call check( nf90_get_var(ncid, ncv, omega_t_conv) )
189 
190 call check( nf90_inq_varid(ncid, 'temp_r', ncv) )
191 call check( nf90_get_var(ncid, ncv, temp_r_conv) )
192 
193 call check( nf90_inq_varid(ncid, 'Q_bm', ncv) )
194 call check( nf90_get_var(ncid, ncv, q_bm_conv) )
195 
196 call check( nf90_inq_varid(ncid, 'Q_tld', ncv) )
197 call check( nf90_get_var(ncid, ncv, q_tld_conv) )
198 
199 call check( nf90_inq_varid(ncid, 'am_perp', ncv) )
200 call check( nf90_get_var(ncid, ncv, am_perp_conv) )
201 
202 call check( nf90_inq_varid(ncid, 'qx', ncv) )
203 call check( nf90_get_var(ncid, ncv, qx_conv) )
204 
205 call check( nf90_inq_varid(ncid, 'qy', ncv) )
206 call check( nf90_get_var(ncid, ncv, qy_conv) )
207 
208 call check( nf90_inq_varid(ncid, 'age_c', ncv) )
209 call check( nf90_get_var(ncid, ncv, age_c_conv) )
210 
211 call check( nf90_inq_varid(ncid, 'age_t', ncv) )
212 call check( nf90_get_var(ncid, ncv, age_t_conv) )
213 
214 call check( nf90_inq_varid(ncid, 'dzs_dtau', ncv) )
215 call check( nf90_get_var(ncid, ncv, dzs_dtau_conv) )
216 
217 call check( nf90_inq_varid(ncid, 'dzm_dtau', ncv) )
218 call check( nf90_get_var(ncid, ncv, dzm_dtau_conv) )
219 
220 call check( nf90_inq_varid(ncid, 'dzb_dtau', ncv) )
221 call check( nf90_get_var(ncid, ncv, dzb_dtau_conv) )
222 
223 call check( nf90_inq_varid(ncid, 'dzl_dtau', ncv) )
224 call check( nf90_get_var(ncid, ncv, dzl_dtau_conv) )
225 
226 call check( nf90_inq_varid(ncid, 'dH_c_dtau', ncv) )
227 call check( nf90_get_var(ncid, ncv, dh_c_dtau_conv) )
228 
229 call check( nf90_inq_varid(ncid, 'dH_t_dtau', ncv) )
230 call check( nf90_get_var(ncid, ncv, dh_t_dtau_conv) )
231 
232 call check( nf90_inq_varid(ncid, 'dH_dtau', ncv) )
233 call check( nf90_get_var(ncid, ncv, dh_dtau_conv) )
234 
235 if ( nf90_inq_varid(ncid, 'vx_b_g', ncv) == nf90_noerr ) then
236  call check( nf90_get_var(ncid, ncv, vx_b_g_conv) )
237 else
238  go to 110
239 end if
240 
241 if ( nf90_inq_varid(ncid, 'vy_b_g', ncv) == nf90_noerr ) then
242  call check( nf90_get_var(ncid, ncv, vy_b_g_conv) )
243 else
244  go to 110
245 end if
246 
247 if ( nf90_inq_varid(ncid, 'vz_b', ncv) == nf90_noerr ) then
248  call check( nf90_get_var(ncid, ncv, vz_b_conv) )
249 else
250  go to 110
251 end if
252 
253 if ( nf90_inq_varid(ncid, 'vh_b', ncv) == nf90_noerr ) then
254  call check( nf90_get_var(ncid, ncv, vh_b_conv) )
255 else
256  go to 110
257 end if
258 
259 if ( nf90_inq_varid(ncid, 'vx_s_g', ncv) == nf90_noerr ) then
260  call check( nf90_get_var(ncid, ncv, vx_s_g_conv) )
261 else
262  go to 110
263 end if
264 
265 if ( nf90_inq_varid(ncid, 'vy_s_g', ncv) == nf90_noerr ) then
266  call check( nf90_get_var(ncid, ncv, vy_s_g_conv) )
267 else
268  go to 110
269 end if
270 
271 if ( nf90_inq_varid(ncid, 'vz_s', ncv) == nf90_noerr ) then
272  call check( nf90_get_var(ncid, ncv, vz_s_conv) )
273 else
274  go to 110
275 end if
276 
277 if ( nf90_inq_varid(ncid, 'vh_s', ncv) == nf90_noerr ) then
278  call check( nf90_get_var(ncid, ncv, vh_s_conv) )
279 else
280  go to 110
281 end if
282 
283 if ( nf90_inq_varid(ncid, 'temp_b', ncv) == nf90_noerr ) then
284  call check( nf90_get_var(ncid, ncv, temp_b_conv) )
285 else
286  go to 110
287 end if
288 
289 if ( nf90_inq_varid(ncid, 'temph_b', ncv) == nf90_noerr ) then
290  call check( nf90_get_var(ncid, ncv, temph_b_conv) )
291 else
292  go to 110
293 end if
294 
295 if ( nf90_inq_varid(ncid, 'p_b_w', ncv) == nf90_noerr ) then
296  call check( nf90_get_var(ncid, ncv, p_b_w_conv) )
297 else
298  go to 110
299 end if
300 
301 if ( nf90_inq_varid(ncid, 'H_w', ncv) == nf90_noerr ) then
302  call check( nf90_get_var(ncid, ncv, h_w_conv) )
303 else
304  go to 110
305 end if
306 
307 if ( nf90_inq_varid(ncid, 'q_gl_g', ncv) == nf90_noerr ) then
308  call check( nf90_get_var(ncid, ncv, q_gl_g_conv) )
309 else
310  go to 110
311 end if
312 
313 go to 120
314 
315 110 continue
316 write(6,'(/1x,a)') 'topography3_nc: End-of-file condition while reading *.nc.'
317 write(6, '(1x,a)') ' Some variables will be undefined.'
318 
319 120 continue
320 
321 call check( nf90_close(ncid) )
322 
323 !-------- Convert data to real*8 and years to seconds --------
324 
325 z_sl = real(z_sl_conv,dp)
326 
327 h_r = real(h_r_conv,dp)
328 
329 do i=0, imax
330  xi(i) = real(xi_conv(i),dp)
331 end do
332 
333 do j=0, jmax
334  eta(j) = real(eta_conv(j),dp)
335 end do
336 
337 do i=0, imax
338 do j=0, jmax
339 
340  maske(j,i) = maske_conv(i,j)
341  n_cts(j,i) = n_cts_conv(i,j)
342  zs(j,i) = real(zs_conv(i,j),dp)
343  zm(j,i) = real(zm_conv(i,j),dp)
344  zb(j,i) = real(zb_conv(i,j),dp)
345  zl(j,i) = real(zl_conv(i,j),dp)
346  h_c(j,i) = real(H_c_conv(i,j),dp)
347  h_t(j,i) = real(H_t_conv(i,j),dp)
348  q_bm(j,i) = real(Q_bm_conv(i,j),dp)/year_sec
349  q_tld(j,i) = real(Q_tld_conv(i,j),dp)/year_sec
350  am_perp(j,i) = real(am_perp_conv(i,j),dp)/year_sec
351  qx(j,i) = real(qx_conv(i,j),dp)/year_sec
352  qy(j,i) = real(qy_conv(i,j),dp)/year_sec
353  dzs_dtau(j,i) = real(dzs_dtau_conv(i,j),dp)/year_sec
354  dzm_dtau(j,i) = real(dzm_dtau_conv(i,j),dp)/year_sec
355  dzb_dtau(j,i) = real(dzb_dtau_conv(i,j),dp)/year_sec
356  dzl_dtau(j,i) = real(dzl_dtau_conv(i,j),dp)/year_sec
357  dh_c_dtau(j,i) = real(dH_c_dtau_conv(i,j),dp)/year_sec
358  dh_t_dtau(j,i) = real(dH_t_dtau_conv(i,j),dp)/year_sec
359  vx_b_g(j,i) = real(vx_b_g_conv(i,j),dp)/year_sec
360  vy_b_g(j,i) = real(vy_b_g_conv(i,j),dp)/year_sec
361  vz_b(j,i) = real(vz_b_conv(i,j),dp)/year_sec
362  vx_s_g(j,i) = real(vx_s_g_conv(i,j),dp)/year_sec
363  vy_s_g(j,i) = real(vy_s_g_conv(i,j),dp)/year_sec
364  vz_s(j,i) = real(vz_s_conv(i,j),dp)/year_sec
365  temp_b(j,i) = real(temp_b_conv(i,j),dp)
366  temph_b(j,i) = real(temph_b_conv(i,j),dp)
367  p_b_w(j,i) = real(p_b_w_conv(i,j),dp)
368  h_w(j,i) = real(H_w_conv(i,j),dp)
369  q_gl_g(j,i) = real(q_gl_g_conv(i,j),dp)/year_sec
370 
371  do kr=0, krmax
372  temp_r(kr,j,i) = real(temp_r_conv(i,j,kr),dp)
373  end do
374 
375  do kt=0, ktmax
376  vx_t(kt,j,i) = real(vx_t_conv(i,j,kt),dp)/year_sec
377  vy_t(kt,j,i) = real(vy_t_conv(i,j,kt),dp)/year_sec
378  vz_t(kt,j,i) = real(vz_t_conv(i,j,kt),dp)/year_sec
379  omega_t(kt,j,i) = real(omega_t_conv(i,j,kt),dp)
380  age_t(kt,j,i) = real(age_t_conv(i,j,kt),dp)*year_sec
381  end do
382 
383  do kc=0, kcmax
384  vx_c(kc,j,i) = real(vx_c_conv(i,j,kc),dp)/year_sec
385  vy_c(kc,j,i) = real(vy_c_conv(i,j,kc),dp)/year_sec
386  vz_c(kc,j,i) = real(vz_c_conv(i,j,kc),dp)/year_sec
387  temp_c(kc,j,i) = real(temp_c_conv(i,j,kc),dp)
388  age_c(kc,j,i) = real(age_c_conv(i,j,kc),dp)*year_sec
389  end do
390 
391 end do
392 end do
393 
394 !-------- Read topography of the relaxed bedrock --------
395 
396 open(23, iostat=ios, &
397  file=inpath//'/tibet/'//zl0_file, &
398  recl=16384, status='old')
399 
400 if (ios /= 0) stop ' topography3_nc: Error when opening the zl0 file!'
401 
402 do n=1, 6; read(23,'(a)') ch_dummy; end do
403 
404 do j=jmax, 0, -1
405  read(23,*) (zl0(j,i), i=0,imax)
406 end do
407 
408 close(23, status='keep')
409 
410 !-------- Further stuff --------
411 
412 dxi = dlambda *pi_180 ! deg -> rad
413 deta = dphi *pi_180 ! deg -> rad
414 
415 do i=0, imax
416 do j=0, jmax
417  call geo_coord(phi(j,i), lambda(j,i), xi(i), eta(j))
418 end do
419 end do
420 
421 !-------- Metric tensor, gradients of the topography --------
422 
423 call metric()
424 
425 #if TOPOGRAD==0
426 call topograd_1(dxi, deta, 1)
427 #elif TOPOGRAD==1
428 call topograd_2(dxi, deta, 1)
429 #endif
430 
431 !-------- Corresponding area of grid points --------
432 
433 do i=0, imax
434 do j=0, jmax
435  area(j,i) = sq_g11_g(j,i)*sq_g22_g(j,i)*dxi*deta
436 end do
437 end do
438 
439 end subroutine topography3_nc
440 !