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