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