SICOPOLIS V3.2
 All Classes Files Functions Variables Macros
topography3.F90
Go to the documentation of this file.
1 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 !
3 ! Subroutine : t o p o g r a p h y 3
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 binary 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 binary format).
38 !<------------------------------------------------------------------------------
39 subroutine topography3(dxi, deta, z_sl, anfdatname)
40 
41 use sico_types
43 use sico_vars
44 
45 implicit none
46 
47 character(len=100), intent(in) :: anfdatname
48 
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(len=256) :: anfdatname_with_path
54 character(len=256) :: ch_attr_title, ch_attr_institution, ch_attr_source, &
55  ch_attr_history, ch_attr_references
56 character :: ch_dummy
57 
58 integer(i2b), dimension(0:IMAX,0:JMAX) :: maske_conv, n_cts_conv, kc_cts_conv
59 real(sp) :: time_conv, dummy_conv, z_sl_conv, &
60  v_tot_conv, a_grounded_conv, a_floating_conv, &
61  h_r_conv, &
62  xi_conv(0:imax), eta_conv(0:jmax), &
63  sigma_level_c_conv(0:kcmax), sigma_level_t_conv(0:ktmax), &
64  sigma_level_r_conv(0:krmax)
65 real(sp), dimension(0:IMAX,0:JMAX) :: lambda_conv, phi_conv, &
66  lon_conv, lat_conv, &
67  temp_s_conv, as_perp_conv, &
68  zs_conv, zm_conv, zb_conv, zl_conv, &
69  h_cold_conv, h_temp_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  enth_c_conv, omega_c_conv, &
82  enh_c_conv
83 real(sp), dimension(0:IMAX,0:JMAX,0:KTMAX) :: vx_t_conv, vy_t_conv, vz_t_conv, &
84  omega_t_conv, age_t_conv, &
85  enth_t_conv, &
86  enh_t_conv
87 real(sp), dimension(0:IMAX,0:JMAX,0:KRMAX) :: temp_r_conv
88 
89 !-------- Read data from time-slice file of previous simulation --------
90 
91 anfdatname_with_path = trim(anfdatpath)//'/'//trim(anfdatname)
92 
93 open(11, iostat=ios, &
94  file=trim(anfdatname_with_path), &
95  status='old', form='unformatted')
96 if (ios /= 0) stop ' topography3: Error when opening the initial-value file!'
97 
98 read(11) ch_attr_title
99 read(11) ch_attr_institution
100 read(11) ch_attr_source
101 read(11) ch_attr_history
102 read(11) ch_attr_references
103 read(11) time_conv
104 read(11) dummy_conv
105 read(11) z_sl_conv
106 read(11) xi_conv
107 read(11) eta_conv
108 read(11) sigma_level_c_conv
109 read(11) sigma_level_t_conv
110 read(11) sigma_level_r_conv
111 read(11) lon_conv
112 read(11) lat_conv
113 read(11) lambda_conv
114 read(11) phi_conv
115 read(11) temp_s_conv
116 read(11) as_perp_conv
117 read(11) maske_conv
118 read(11) n_cts_conv
119 read(11) kc_cts_conv
120 read(11) zs_conv
121 read(11) zm_conv
122 read(11) zb_conv
123 read(11) zl_conv
124 read(11) h_cold_conv
125 read(11) h_temp_conv
126 read(11) h_conv
127 read(11) h_r_conv
128 read(11) vx_c_conv
129 read(11) vy_c_conv
130 read(11) vz_c_conv
131 read(11) vx_t_conv
132 read(11) vy_t_conv
133 read(11) vz_t_conv
134 read(11) temp_c_conv
135 read(11) omega_t_conv
136 read(11) temp_r_conv
137 read(11) enth_c_conv
138 read(11) enth_t_conv
139 read(11) omega_c_conv
140 read(11) enh_c_conv
141 read(11) enh_t_conv
142 read(11) q_bm_conv
143 read(11) q_tld_conv
144 read(11) am_perp_conv
145 read(11) qx_conv
146 read(11) qy_conv
147 read(11) age_c_conv
148 read(11) age_t_conv
149 read(11) dzs_dtau_conv
150 read(11) dzm_dtau_conv
151 read(11) dzb_dtau_conv
152 read(11) dzl_dtau_conv
153 read(11) dh_c_dtau_conv
154 read(11) dh_t_dtau_conv
155 read(11) dh_dtau_conv
156 read(11) vx_b_g_conv
157 read(11) vy_b_g_conv
158 read(11) vz_b_conv
159 read(11) vh_b_conv
160 read(11) vx_s_g_conv
161 read(11) vy_s_g_conv
162 read(11) vz_s_conv
163 read(11) vh_s_conv
164 read(11) temp_b_conv
165 read(11) temph_b_conv
166 read(11) p_b_w_conv
167 read(11) h_w_conv
168 read(11) q_gl_g_conv
169 
170 go to 120
171 
172 110 continue
173 write(6, fmt='(/1x,a)') 'topography3: End-of-file condition while reading *.erg.'
174 write(6, fmt= '(1x,a)') ' Some variables will be undefined.'
175 
176 120 continue
177 
178 close(11,status='keep')
179 
180 !-------- Convert data to real*8 and years to seconds --------
181 
182 z_sl = real(z_sl_conv,dp)
183 
184 h_r = real(h_r_conv,dp)
185 
186 do i=0, imax
187  xi(i) = real(xi_conv(i),dp)
188 end do
189 
190 do j=0, jmax
191  eta(j) = real(eta_conv(j),dp)
192 end do
193 
194 do i=0, imax
195 do j=0, jmax
196 
197  maske(j,i) = maske_conv(i,j)
198  n_cts(j,i) = n_cts_conv(i,j)
199  kc_cts(j,i) = kc_cts_conv(i,j)
200  zs(j,i) = real(zs_conv(i,j),dp)
201  zm(j,i) = real(zm_conv(i,j),dp)
202  zb(j,i) = real(zb_conv(i,j),dp)
203  zl(j,i) = real(zl_conv(i,j),dp)
204 #if (CALCMOD==1)
205  h_c(j,i) = real(H_cold_conv(i,j),dp)
206  h_t(j,i) = real(H_temp_conv(i,j),dp)
207 #elif (CALCMOD==0 || CALCMOD==2 || CALCMOD==3 || CALCMOD==-1)
208  h_c(j,i) = real(H_conv(i,j),dp)
209  h_t(j,i) = 0.0_dp
210 #endif
211  q_bm(j,i) = real(Q_bm_conv(i,j),dp)/year_sec
212  q_tld(j,i) = real(Q_tld_conv(i,j),dp)/year_sec
213  am_perp(j,i) = real(am_perp_conv(i,j),dp)/year_sec
214  qx(j,i) = real(qx_conv(i,j),dp)/year_sec
215  qy(j,i) = real(qy_conv(i,j),dp)/year_sec
216  dzs_dtau(j,i) = real(dzs_dtau_conv(i,j),dp)/year_sec
217  dzm_dtau(j,i) = real(dzm_dtau_conv(i,j),dp)/year_sec
218  dzb_dtau(j,i) = real(dzb_dtau_conv(i,j),dp)/year_sec
219  dzl_dtau(j,i) = real(dzl_dtau_conv(i,j),dp)/year_sec
220  dh_c_dtau(j,i) = real(dH_c_dtau_conv(i,j),dp)/year_sec
221  dh_t_dtau(j,i) = real(dH_t_dtau_conv(i,j),dp)/year_sec
222  vx_b_g(j,i) = real(vx_b_g_conv(i,j),dp)/year_sec
223  vy_b_g(j,i) = real(vy_b_g_conv(i,j),dp)/year_sec
224  vz_b(j,i) = real(vz_b_conv(i,j),dp)/year_sec
225  vx_s_g(j,i) = real(vx_s_g_conv(i,j),dp)/year_sec
226  vy_s_g(j,i) = real(vy_s_g_conv(i,j),dp)/year_sec
227  vz_s(j,i) = real(vz_s_conv(i,j),dp)/year_sec
228  temp_b(j,i) = real(temp_b_conv(i,j),dp)
229  temph_b(j,i) = real(temph_b_conv(i,j),dp)
230  p_b_w(j,i) = real(p_b_w_conv(i,j),dp)
231  h_w(j,i) = real(H_w_conv(i,j),dp)
232  q_gl_g(j,i) = real(q_gl_g_conv(i,j),dp)/year_sec
233 
234  do kr=0, krmax
235  temp_r(kr,j,i) = real(temp_r_conv(i,j,kr),dp)
236  end do
237 
238  do kt=0, ktmax
239  vx_t(kt,j,i) = real(vx_t_conv(i,j,kt),dp)/year_sec
240  vy_t(kt,j,i) = real(vy_t_conv(i,j,kt),dp)/year_sec
241  vz_t(kt,j,i) = real(vz_t_conv(i,j,kt),dp)/year_sec
242  omega_t(kt,j,i) = real(omega_t_conv(i,j,kt),dp)
243  age_t(kt,j,i) = real(age_t_conv(i,j,kt),dp)*year_sec
244  enth_t(kt,j,i) = real(enth_t_conv(i,j,kt),dp)
245  enh_t(kt,j,i) = real(enh_t_conv(i,j,kt),dp)
246  end do
247 
248  do kc=0, kcmax
249  vx_c(kc,j,i) = real(vx_c_conv(i,j,kc),dp)/year_sec
250  vy_c(kc,j,i) = real(vy_c_conv(i,j,kc),dp)/year_sec
251  vz_c(kc,j,i) = real(vz_c_conv(i,j,kc),dp)/year_sec
252  temp_c(kc,j,i) = real(temp_c_conv(i,j,kc),dp)
253  age_c(kc,j,i) = real(age_c_conv(i,j,kc),dp)*year_sec
254  enth_c(kc,j,i) = real(enth_c_conv(i,j,kc),dp)
255  omega_c(kc,j,i) = real(omega_c_conv(i,j,kc),dp)
256  enh_c(kc,j,i) = real(enh_c_conv(i,j,kc),dp)
257  end do
258 
259 end do
260 end do
261 
262 !-------- Read topography of the relaxed bedrock --------
263 
264 open(23, iostat=ios, &
265  file=inpath//'/'//trim(ch_domain_short)//'/'//zl0_file, &
266  recl=16384, status='old')
267 
268 if (ios /= 0) stop ' topography3: Error when opening the zl0 file!'
269 
270 do n=1, 6; read(23, fmt='(a)') ch_dummy; end do
271 
272 do j=jmax, 0, -1
273  read(23, fmt=*) (zl0(j,i), i=0,imax)
274 end do
275 
276 close(23, status='keep')
277 
278 !-------- Further stuff --------
279 
280 dxi = dx *1000.0_dp ! km -> m
281 deta = dx *1000.0_dp ! km -> m
282 
283 do i=0, imax
284 do j=0, jmax
285  zl0(j,i) = zl0(j,i)*1000.0_dp ! km --> m
286  call geo_coord(phi(j,i), lambda(j,i), xi(i), eta(j))
287 end do
288 end do
289 
290 !-------- Metric tensor, gradients of the topography --------
291 
292 call metric()
293 
294 #if TOPOGRAD==0
295 call topograd_1(dxi, deta, 1)
296 #elif TOPOGRAD==1
297 call topograd_2(dxi, deta, 1)
298 #endif
299 
300 !-------- Corresponding area of grid points --------
301 
302 do i=0, imax
303 do j=0, jmax
304  area(j,i) = sq_g11_g(j,i)*sq_g22_g(j,i)*dxi*deta
305 end do
306 end do
307 
308 end subroutine topography3
309 !
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 topography3(dxi, deta, z_sl, anfdatname)
Definition of the initial surface and bedrock topography (including gradients) and of the horizontal ...
Definition: topography3.F90:39
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
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.