SICOPOLIS V3.0
 All Classes Files Functions Variables Macros
topography2.F90
Go to the documentation of this file.
1 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 !
3 ! Subroutine : t o p o g r a p h y 2
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 ice-free initial topography with relaxed lithosphere.
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 ice-free initial topography with relaxed lithosphere.
38 !<------------------------------------------------------------------------------
39 subroutine topography2(dxi, deta)
40 
41 use sico_types
43 
44 implicit none
45 integer(i4b) :: i, j, n
46 integer(i4b) :: ios
47 real(dp) :: dxi, deta
48 real(dp) :: xi0, eta0
49 character :: ch_dummy
50 
51 !-------- Read topography --------
52 
53 open(23, iostat=ios, &
54  file=inpath//'/scand/'//zl0_file, &
55  recl=8192, status='old')
56 
57 if (ios.ne.0) stop ' topography2: Error when opening the zl0 file!'
58 
59 open(24, iostat=ios, &
60  file=inpath//'/scand/'//mask_present_file, &
61  recl=1024, status='old')
62 
63 if (ios.ne.0) stop ' topography2: Error when opening the mask file!'
64 
65 do n=1, 6; read(23,'(a)') ch_dummy; end do
66 
67 do j=jmax, 0, -1
68  read(23,*) (zl0(j,i), i=0,imax)
69 end do
70 
71 do n=1, 6; read(24,'(a)') ch_dummy; end do
72 
73 do j=jmax, 0, -1
74  read(24,2300) (maske(j,i), i=0,imax)
75 end do
76 
77 close(23, status='keep')
78 close(24, status='keep')
79 
80 2300 format(imax(i1),i1)
81 
82 !-------- Further stuff --------
83 
84 dxi = dx *1000.0_dp ! km -> m
85 deta = dx *1000.0_dp ! km -> m
86 
87 xi0 = x0 *1000.0_dp ! km -> m
88 eta0 = y0 *1000.0_dp ! km -> m
89 
90 do i=0, imax
91 do j=0, jmax
92 
93  if (maske(j,i) <= 1 ) then
94  maske(j,i) = 1
95  zs(j,i) = zl0(j,i)
96  zb(j,i) = zl0(j,i)
97  zl(j,i) = zl0(j,i)
98  else ! (maske(j,i) >= 2 )
99  maske(j,i) = 2
100 #if ( MARGIN==1 || MARGIN==2 )
101  zs(j,i) = zl0(j,i)
102  zb(j,i) = zl0(j,i)
103 #elif ( MARGIN==3 )
104  zs(j,i) = 0.0_dp ! present-day
105  zb(j,i) = 0.0_dp ! sea level
106 #endif
107  zl(j,i) = zl0(j,i)
108  end if
109 
110  xi(i) = xi0 + real(i,dp)*dxi
111  eta(j) = eta0 + real(j,dp)*deta
112 
113  call geo_coord(phi(j,i), lambda(j,i), xi(i), eta(j))
114 
115  zm(j,i) = zb(j,i)
116  n_cts(j,i) = -1
117 
118  h_c(j,i) = 0.0_dp
119  h_t(j,i) = 0.0_dp
120 
121  dzs_dtau(j,i) = 0.0_dp
122  dzm_dtau(j,i) = 0.0_dp
123  dzb_dtau(j,i) = 0.0_dp
124  dzl_dtau(j,i) = 0.0_dp
125  dh_c_dtau(j,i) = 0.0_dp
126  dh_t_dtau(j,i) = 0.0_dp
127 
128 end do
129 end do
130 
131 !-------- Metric tensor, gradients of the topography --------
132 
133 call metric
134 
135 #if TOPOGRAD==0
136 call topograd_1(dxi, deta, 1)
137 #elif TOPOGRAD==1
138 call topograd_2(dxi, deta, 1)
139 #endif
140 
141 !-------- Corresponding area of grid points --------
142 
143 do i=0, imax
144 do j=0, jmax
145  area(j,i) = sq_g11_g(j,i)*sq_g22_g(j,i)*dxi*deta
146 end do
147 end do
148 
149 end subroutine topography2
150 !