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