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