SICOPOLIS V3.2
 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-2016 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 sico_vars
44 use netcdf
45 use nc_check
46 
47 implicit none
48 
49 character(len=100), intent(in) :: anfdatname
50 
51 real(dp), intent(out) :: dxi, deta, z_sl
52 
53 integer(i4b) :: i, j, kc, kt, kr
54 integer(i4b) :: ios
55 character(len=256) :: anfdatname_with_path
56 
57 integer(i4b) :: ncid, ncv
58 ! ncid: ID of the output file
59 ! ncv: Variable ID
60 
61 integer(i2b), dimension(0:IMAX,0:JMAX) :: maske_conv, n_cts_conv, kc_cts_conv
62 real(sp) :: time_conv, dummy_conv, z_sl_conv, &
63  v_tot_conv, a_grounded_conv, a_floating_conv, &
64  h_r_conv, &
65  xi_conv(0:imax), eta_conv(0:jmax), &
66  sigma_level_c_conv(0:kcmax), sigma_level_t_conv(0:ktmax), &
67  sigma_level_r_conv(0:krmax)
68 real(sp), dimension(0:IMAX,0:JMAX) :: lambda_conv, phi_conv, &
69  lon_conv, lat_conv, &
70  temp_s_conv, as_perp_conv, &
71  zs_conv, zm_conv, zb_conv, zl_conv, &
72  h_cold_conv, h_temp_conv, h_conv, &
73  q_bm_conv, q_tld_conv, &
74  am_perp_conv, &
75  qx_conv, qy_conv, &
76  dzs_dtau_conv, dzm_dtau_conv, dzb_dtau_conv, dzl_dtau_conv, &
77  dh_c_dtau_conv, dh_t_dtau_conv, dh_dtau_conv, &
78  vx_b_g_conv, vy_b_g_conv, vz_b_conv, vh_b_conv, &
79  vx_s_g_conv, vy_s_g_conv, vz_s_conv, vh_s_conv, &
80  temp_b_conv, temph_b_conv, &
81  p_b_w_conv, h_w_conv, q_gl_g_conv
82 real(sp), dimension(0:IMAX,0:JMAX,0:KCMAX) :: vx_c_conv, vy_c_conv, vz_c_conv, &
83  temp_c_conv, age_c_conv, &
84  enth_c_conv, omega_c_conv, &
85  enh_c_conv
86 real(sp), dimension(0:IMAX,0:JMAX,0:KTMAX) :: vx_t_conv, vy_t_conv, vz_t_conv, &
87  omega_t_conv, age_t_conv, &
88  enth_t_conv, &
89  enh_t_conv
90 real(sp), dimension(0:IMAX,0:JMAX,0:KRMAX) :: temp_r_conv
91 
92 !-------- Read data from time-slice file of previous simulation --------
93 
94 anfdatname_with_path = trim(anfdatpath)//'/'//trim(anfdatname)
95 
96 call check( nf90_open(trim(anfdatname_with_path), nf90_nowrite, ncid) )
97 
98 call check( nf90_inq_varid(ncid, 'time', ncv) )
99 call check( nf90_get_var(ncid, ncv, time_conv) )
100 
101 if ( nf90_inq_varid(ncid, 'delta_ts', ncv) == nf90_noerr ) then
102  call check( nf90_get_var(ncid, ncv, dummy_conv) )
103 else if ( nf90_inq_varid(ncid, 'glac_index', ncv) == nf90_noerr ) then
104  call check( nf90_get_var(ncid, ncv, dummy_conv) )
105 else
106  continue
107 end if
108 
109 call check( nf90_inq_varid(ncid, 'z_sl', ncv) )
110 call check( nf90_get_var(ncid, ncv, z_sl_conv) )
111 
112 call check( nf90_inq_varid(ncid, 'xi', ncv) )
113 call check( nf90_get_var(ncid, ncv, xi_conv) )
114 
115 call check( nf90_inq_varid(ncid, 'eta', ncv) )
116 call check( nf90_get_var(ncid, ncv, eta_conv) )
117 
118 call check( nf90_inq_varid(ncid, 'sigma_level_c', ncv) )
119 call check( nf90_get_var(ncid, ncv, sigma_level_c_conv) )
120 
121 call check( nf90_inq_varid(ncid, 'sigma_level_t', ncv) )
122 call check( nf90_get_var(ncid, ncv, sigma_level_t_conv) )
123 
124 call check( nf90_inq_varid(ncid, 'sigma_level_r', ncv) )
125 call check( nf90_get_var(ncid, ncv, sigma_level_r_conv) )
126 
127 call check( nf90_inq_varid(ncid, 'lon', ncv) )
128 call check( nf90_get_var(ncid, ncv, lon_conv) )
129 
130 call check( nf90_inq_varid(ncid, 'lat', ncv) )
131 call check( nf90_get_var(ncid, ncv, lat_conv) )
132 
133 call check( nf90_inq_varid(ncid, 'lambda', ncv) )
134 call check( nf90_get_var(ncid, ncv, lambda_conv) )
135 
136 call check( nf90_inq_varid(ncid, 'phi', ncv) )
137 call check( nf90_get_var(ncid, ncv, phi_conv) )
138 
139 call check( nf90_inq_varid(ncid, 'temp_s', ncv) )
140 call check( nf90_get_var(ncid, ncv, temp_s_conv) )
141 
142 call check( nf90_inq_varid(ncid, 'as_perp', ncv) )
143 call check( nf90_get_var(ncid, ncv, as_perp_conv) )
144 
145 call check( nf90_inq_varid(ncid, 'maske', ncv) )
146 call check( nf90_get_var(ncid, ncv, maske_conv) )
147 
148 call check( nf90_inq_varid(ncid, 'n_cts', ncv) )
149 call check( nf90_get_var(ncid, ncv, n_cts_conv) )
150 
151 if ( nf90_inq_varid(ncid, 'kc_cts', ncv) == nf90_noerr ) then
152  call check( nf90_get_var(ncid, ncv, kc_cts_conv) )
153 else
154  go to 110
155 end if
156 
157 call check( nf90_inq_varid(ncid, 'zs', ncv) )
158 call check( nf90_get_var(ncid, ncv, zs_conv) )
159 
160 call check( nf90_inq_varid(ncid, 'zm', ncv) )
161 call check( nf90_get_var(ncid, ncv, zm_conv) )
162 
163 call check( nf90_inq_varid(ncid, 'zb', ncv) )
164 call check( nf90_get_var(ncid, ncv, zb_conv) )
165 
166 call check( nf90_inq_varid(ncid, 'zl', ncv) )
167 call check( nf90_get_var(ncid, ncv, zl_conv) )
168 
169 call check( nf90_inq_varid(ncid, 'H_cold', ncv) )
170 call check( nf90_get_var(ncid, ncv, h_cold_conv) )
171 
172 call check( nf90_inq_varid(ncid, 'H_temp', ncv) )
173 call check( nf90_get_var(ncid, ncv, h_temp_conv) )
174 
175 call check( nf90_inq_varid(ncid, 'H', ncv) )
176 call check( nf90_get_var(ncid, ncv, h_conv) )
177 
178 call check( nf90_inq_varid(ncid, 'H_R', ncv) )
179 call check( nf90_get_var(ncid, ncv, h_r_conv) )
180 
181 call check( nf90_inq_varid(ncid, 'vx_c', ncv) )
182 call check( nf90_get_var(ncid, ncv, vx_c_conv) )
183 
184 call check( nf90_inq_varid(ncid, 'vy_c', ncv) )
185 call check( nf90_get_var(ncid, ncv, vy_c_conv) )
186 
187 call check( nf90_inq_varid(ncid, 'vz_c', ncv) )
188 call check( nf90_get_var(ncid, ncv, vz_c_conv) )
189 
190 call check( nf90_inq_varid(ncid, 'vx_t', ncv) )
191 call check( nf90_get_var(ncid, ncv, vx_t_conv) )
192 
193 call check( nf90_inq_varid(ncid, 'vy_t', ncv) )
194 call check( nf90_get_var(ncid, ncv, vy_t_conv) )
195 
196 call check( nf90_inq_varid(ncid, 'vz_t', ncv) )
197 call check( nf90_get_var(ncid, ncv, vz_t_conv) )
198 
199 call check( nf90_inq_varid(ncid, 'temp_c', ncv) )
200 call check( nf90_get_var(ncid, ncv, temp_c_conv) )
201 
202 call check( nf90_inq_varid(ncid, 'omega_t', ncv) )
203 call check( nf90_get_var(ncid, ncv, omega_t_conv) )
204 
205 call check( nf90_inq_varid(ncid, 'temp_r', ncv) )
206 call check( nf90_get_var(ncid, ncv, temp_r_conv) )
207 
208 if ( nf90_inq_varid(ncid, 'enth_c', ncv) == nf90_noerr ) then
209  call check( nf90_get_var(ncid, ncv, enth_c_conv) )
210 else
211  go to 110
212 end if
213 
214 if ( nf90_inq_varid(ncid, 'enth_t', ncv) == nf90_noerr ) then
215  call check( nf90_get_var(ncid, ncv, enth_t_conv) )
216 else
217  go to 110
218 end if
219 
220 if ( nf90_inq_varid(ncid, 'omega_c', ncv) == nf90_noerr ) then
221  call check( nf90_get_var(ncid, ncv, omega_c_conv) )
222 else
223  go to 110
224 end if
225 
226 if ( nf90_inq_varid(ncid, 'enh_c', ncv) == nf90_noerr ) then
227  call check( nf90_get_var(ncid, ncv, enh_c_conv) )
228 else
229  go to 110
230 end if
231 
232 if ( nf90_inq_varid(ncid, 'enh_t', ncv) == nf90_noerr ) then
233  call check( nf90_get_var(ncid, ncv, enh_t_conv) )
234 else
235  go to 110
236 end if
237 
238 call check( nf90_inq_varid(ncid, 'Q_bm', ncv) )
239 call check( nf90_get_var(ncid, ncv, q_bm_conv) )
240 
241 call check( nf90_inq_varid(ncid, 'Q_tld', ncv) )
242 call check( nf90_get_var(ncid, ncv, q_tld_conv) )
243 
244 call check( nf90_inq_varid(ncid, 'am_perp', ncv) )
245 call check( nf90_get_var(ncid, ncv, am_perp_conv) )
246 
247 call check( nf90_inq_varid(ncid, 'qx', ncv) )
248 call check( nf90_get_var(ncid, ncv, qx_conv) )
249 
250 call check( nf90_inq_varid(ncid, 'qy', ncv) )
251 call check( nf90_get_var(ncid, ncv, qy_conv) )
252 
253 call check( nf90_inq_varid(ncid, 'age_c', ncv) )
254 call check( nf90_get_var(ncid, ncv, age_c_conv) )
255 
256 call check( nf90_inq_varid(ncid, 'age_t', ncv) )
257 call check( nf90_get_var(ncid, ncv, age_t_conv) )
258 
259 call check( nf90_inq_varid(ncid, 'dzs_dtau', ncv) )
260 call check( nf90_get_var(ncid, ncv, dzs_dtau_conv) )
261 
262 call check( nf90_inq_varid(ncid, 'dzm_dtau', ncv) )
263 call check( nf90_get_var(ncid, ncv, dzm_dtau_conv) )
264 
265 call check( nf90_inq_varid(ncid, 'dzb_dtau', ncv) )
266 call check( nf90_get_var(ncid, ncv, dzb_dtau_conv) )
267 
268 call check( nf90_inq_varid(ncid, 'dzl_dtau', ncv) )
269 call check( nf90_get_var(ncid, ncv, dzl_dtau_conv) )
270 
271 call check( nf90_inq_varid(ncid, 'dH_c_dtau', ncv) )
272 call check( nf90_get_var(ncid, ncv, dh_c_dtau_conv) )
273 
274 call check( nf90_inq_varid(ncid, 'dH_t_dtau', ncv) )
275 call check( nf90_get_var(ncid, ncv, dh_t_dtau_conv) )
276 
277 call check( nf90_inq_varid(ncid, 'dH_dtau', ncv) )
278 call check( nf90_get_var(ncid, ncv, dh_dtau_conv) )
279 
280 if ( nf90_inq_varid(ncid, 'vx_b_g', ncv) == nf90_noerr ) then
281  call check( nf90_get_var(ncid, ncv, vx_b_g_conv) )
282 else
283  go to 110
284 end if
285 
286 if ( nf90_inq_varid(ncid, 'vy_b_g', ncv) == nf90_noerr ) then
287  call check( nf90_get_var(ncid, ncv, vy_b_g_conv) )
288 else
289  go to 110
290 end if
291 
292 if ( nf90_inq_varid(ncid, 'vz_b', ncv) == nf90_noerr ) then
293  call check( nf90_get_var(ncid, ncv, vz_b_conv) )
294 else
295  go to 110
296 end if
297 
298 if ( nf90_inq_varid(ncid, 'vh_b', ncv) == nf90_noerr ) then
299  call check( nf90_get_var(ncid, ncv, vh_b_conv) )
300 else
301  go to 110
302 end if
303 
304 if ( nf90_inq_varid(ncid, 'vx_s_g', ncv) == nf90_noerr ) then
305  call check( nf90_get_var(ncid, ncv, vx_s_g_conv) )
306 else
307  go to 110
308 end if
309 
310 if ( nf90_inq_varid(ncid, 'vy_s_g', ncv) == nf90_noerr ) then
311  call check( nf90_get_var(ncid, ncv, vy_s_g_conv) )
312 else
313  go to 110
314 end if
315 
316 if ( nf90_inq_varid(ncid, 'vz_s', ncv) == nf90_noerr ) then
317  call check( nf90_get_var(ncid, ncv, vz_s_conv) )
318 else
319  go to 110
320 end if
321 
322 if ( nf90_inq_varid(ncid, 'vh_s', ncv) == nf90_noerr ) then
323  call check( nf90_get_var(ncid, ncv, vh_s_conv) )
324 else
325  go to 110
326 end if
327 
328 if ( nf90_inq_varid(ncid, 'temp_b', ncv) == nf90_noerr ) then
329  call check( nf90_get_var(ncid, ncv, temp_b_conv) )
330 else
331  go to 110
332 end if
333 
334 if ( nf90_inq_varid(ncid, 'temph_b', ncv) == nf90_noerr ) then
335  call check( nf90_get_var(ncid, ncv, temph_b_conv) )
336 else
337  go to 110
338 end if
339 
340 if ( nf90_inq_varid(ncid, 'p_b_w', ncv) == nf90_noerr ) then
341  call check( nf90_get_var(ncid, ncv, p_b_w_conv) )
342 else
343  go to 110
344 end if
345 
346 if ( nf90_inq_varid(ncid, 'H_w', ncv) == nf90_noerr ) then
347  call check( nf90_get_var(ncid, ncv, h_w_conv) )
348 else
349  go to 110
350 end if
351 
352 if ( nf90_inq_varid(ncid, 'q_gl_g', ncv) == nf90_noerr ) then
353  call check( nf90_get_var(ncid, ncv, q_gl_g_conv) )
354 else
355  go to 110
356 end if
357 
358 go to 120
359 
360 110 continue
361 write(6,'(/1x,a)') 'topography3_nc: End-of-file condition while reading *.nc.'
362 write(6, '(1x,a)') ' Some variables will be undefined.'
363 
364 120 continue
365 
366 call check( nf90_close(ncid) )
367 
368 !-------- Convert data to real*8 and years to seconds --------
369 
370 z_sl = real(z_sl_conv,dp)
371 
372 h_r = real(h_r_conv,dp)
373 
374 do i=0, imax
375  xi(i) = real(xi_conv(i),dp)
376 end do
377 
378 do j=0, jmax
379  eta(j) = real(eta_conv(j),dp)
380 end do
381 
382 do i=0, imax
383 do j=0, jmax
384 
385  maske(j,i) = maske_conv(i,j)
386  n_cts(j,i) = n_cts_conv(i,j)
387  kc_cts(j,i) = kc_cts_conv(i,j)
388  zs(j,i) = real(zs_conv(i,j),dp)
389  zm(j,i) = real(zm_conv(i,j),dp)
390  zb(j,i) = real(zb_conv(i,j),dp)
391  zl(j,i) = real(zl_conv(i,j),dp)
392 #if (CALCMOD==1)
393  h_c(j,i) = real(H_cold_conv(i,j),dp)
394  h_t(j,i) = real(H_temp_conv(i,j),dp)
395 #elif (CALCMOD==0 || CALCMOD==2 || CALCMOD==3 || CALCMOD==-1)
396  h_c(j,i) = real(H_conv(i,j),dp)
397  h_t(j,i) = 0.0_dp
398 #endif
399  q_bm(j,i) = real(Q_bm_conv(i,j),dp)/year_sec
400  q_tld(j,i) = real(Q_tld_conv(i,j),dp)/year_sec
401  am_perp(j,i) = real(am_perp_conv(i,j),dp)/year_sec
402  qx(j,i) = real(qx_conv(i,j),dp)/year_sec
403  qy(j,i) = real(qy_conv(i,j),dp)/year_sec
404  dzs_dtau(j,i) = real(dzs_dtau_conv(i,j),dp)/year_sec
405  dzm_dtau(j,i) = real(dzm_dtau_conv(i,j),dp)/year_sec
406  dzb_dtau(j,i) = real(dzb_dtau_conv(i,j),dp)/year_sec
407  dzl_dtau(j,i) = real(dzl_dtau_conv(i,j),dp)/year_sec
408  dh_c_dtau(j,i) = real(dH_c_dtau_conv(i,j),dp)/year_sec
409  dh_t_dtau(j,i) = real(dH_t_dtau_conv(i,j),dp)/year_sec
410  vx_b_g(j,i) = real(vx_b_g_conv(i,j),dp)/year_sec
411  vy_b_g(j,i) = real(vy_b_g_conv(i,j),dp)/year_sec
412  vz_b(j,i) = real(vz_b_conv(i,j),dp)/year_sec
413  vx_s_g(j,i) = real(vx_s_g_conv(i,j),dp)/year_sec
414  vy_s_g(j,i) = real(vy_s_g_conv(i,j),dp)/year_sec
415  vz_s(j,i) = real(vz_s_conv(i,j),dp)/year_sec
416  temp_b(j,i) = real(temp_b_conv(i,j),dp)
417  temph_b(j,i) = real(temph_b_conv(i,j),dp)
418  p_b_w(j,i) = real(p_b_w_conv(i,j),dp)
419  h_w(j,i) = real(H_w_conv(i,j),dp)
420  q_gl_g(j,i) = real(q_gl_g_conv(i,j),dp)/year_sec
421 
422  do kr=0, krmax
423  temp_r(kr,j,i) = real(temp_r_conv(i,j,kr),dp)
424  end do
425 
426  do kt=0, ktmax
427  vx_t(kt,j,i) = real(vx_t_conv(i,j,kt),dp)/year_sec
428  vy_t(kt,j,i) = real(vy_t_conv(i,j,kt),dp)/year_sec
429  vz_t(kt,j,i) = real(vz_t_conv(i,j,kt),dp)/year_sec
430  omega_t(kt,j,i) = real(omega_t_conv(i,j,kt),dp)
431  age_t(kt,j,i) = real(age_t_conv(i,j,kt),dp)*year_sec
432  enth_t(kt,j,i) = real(enth_t_conv(i,j,kt),dp)
433  enh_t(kt,j,i) = real(enh_t_conv(i,j,kt),dp)
434  end do
435 
436  do kc=0, kcmax
437  vx_c(kc,j,i) = real(vx_c_conv(i,j,kc),dp)/year_sec
438  vy_c(kc,j,i) = real(vy_c_conv(i,j,kc),dp)/year_sec
439  vz_c(kc,j,i) = real(vz_c_conv(i,j,kc),dp)/year_sec
440  temp_c(kc,j,i) = real(temp_c_conv(i,j,kc),dp)
441  age_c(kc,j,i) = real(age_c_conv(i,j,kc),dp)*year_sec
442  enth_c(kc,j,i) = real(enth_c_conv(i,j,kc),dp)
443  omega_c(kc,j,i) = real(omega_c_conv(i,j,kc),dp)
444  enh_c(kc,j,i) = real(enh_c_conv(i,j,kc),dp)
445  end do
446 
447 end do
448 end do
449 
450 !-------- Set topography of the relaxed bedrock --------
451 
452 zl0 = 0.0_dp
453 
454 !-------- Further stuff --------
455 
456 dxi = dx *1000.0_dp ! km -> m
457 deta = dx *1000.0_dp ! km -> m
458 
459 do i=0, imax
460 do j=0, jmax
461  call geo_coord(phi(j,i), lambda(j,i), xi(i), eta(j))
462 end do
463 end do
464 
465 !-------- Metric tensor, gradients of the topography --------
466 
467 call metric()
468 
469 #if TOPOGRAD==0
470 call topograd_1(dxi, deta, 1)
471 #elif TOPOGRAD==1
472 call topograd_2(dxi, deta, 1)
473 #endif
474 
475 !-------- Corresponding area of grid points --------
476 
477 do i=0, imax
478 do j=0, jmax
479  area(j,i) = sq_g11_g(j,i)*sq_g22_g(j,i)*dxi*deta
480 end do
481 end do
482 
483 end subroutine topography3_nc
484 !
NetCDF error capturing.
Definition: nc_check.F90:35
subroutine topograd_2(dxi, deta, n_switch)
Calculation of topography gradients on the staggered grid and on the grid points (the latter by fourt...
Definition: topograd_2.F90:41
Declarations of kind types for SICOPOLIS.
Definition: sico_types.F90:35
subroutine topograd_1(dxi, deta, n_switch)
Calculation of topography gradients on the staggered grid and on the grid points (the latter by secon...
Definition: topograd_1.F90:41
subroutine topography3_nc(dxi, deta, z_sl, anfdatname)
Definition of the initial surface and bedrock topography (including gradients) and of the horizontal ...
subroutine check(status, ch_calling_routine)
NetCDF error capturing.
Definition: nc_check.F90:45
Declarations of global variables for SICOPOLIS (for the ANT domain).
Definition: sico_vars.F90:35
subroutine geo_coord(phi_val, lambda_val, x_val, y_val)
Computation of longitude lambda and latitude phi for position (x,y) in the numerical domain...
Definition: geo_coord.F90:37
subroutine metric()
Definition of the components g11 and g22 of the metric tensor of the applied coordinates.
Definition: metric.F90:37
Declarations of global variables for SICOPOLIS.