SICOPOLIS V3.1
 All Classes Files Functions Variables Macros
topography1.F90
Go to the documentation of this file.
1 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 !
3 ! Subroutine : t o p o g r a p h y 1
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 an initial topography with a thin ice layer (thickness 2*H_MIN)
10 !! everywhere on the land area.
11 !!
12 !! @section Copyright
13 !!
14 !! Copyright 2009-2013 Ralf Greve
15 !!
16 !! @section License
17 !!
18 !! This file is part of SICOPOLIS.
19 !!
20 !! SICOPOLIS is free software: you can redistribute it and/or modify
21 !! it under the terms of the GNU General Public License as published by
22 !! the Free Software Foundation, either version 3 of the License, or
23 !! (at your option) any later version.
24 !!
25 !! SICOPOLIS is distributed in the hope that it will be useful,
26 !! but WITHOUT ANY WARRANTY; without even the implied warranty of
27 !! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
28 !! GNU General Public License for more details.
29 !!
30 !! You should have received a copy of the GNU General Public License
31 !! along with SICOPOLIS. If not, see <http://www.gnu.org/licenses/>.
32 !<
33 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
34 
35 !-------------------------------------------------------------------------------
36 !> Definition of the initial surface and bedrock topography
37 !! (including gradients) and of the horizontal grid spacings dxi, deta.
38 !! For an initial topography with a thin ice layer (thickness 2*H_MIN)
39 !! everywhere on the land area.
40 !<------------------------------------------------------------------------------
41 subroutine topography1(dxi, deta)
42 
43 use sico_types
45 
46 implicit none
47 
48 real(dp), intent(out) :: dxi, deta
49 
50 integer(i4b) :: i, j, n
51 integer(i4b) :: ios, n_dummy
52 real(dp) :: d_dummy
53 real(dp) :: xi0, eta0
54 character :: ch_dummy
55 
56 character(len= 8) :: ch_imax
57 character(len=128) :: fmt4
58 
59 write(ch_imax, fmt='(i8)') imax
60 write(fmt4, fmt='(a)') '('//trim(adjustl(ch_imax))//'(i1),i1)'
61 
62 !-------- Set topography --------
63 
64 zl0 = 0.0_dp
65 
66 open(24, iostat=ios, &
67  file=inpath//'/heino/'//mask_present_file, &
68  recl=1024, status='old')
69 
70 if (ios /= 0) stop ' topography1: Error when opening the mask file!'
71 
72 do n=1, 6; read(24, fmt='(a)') ch_dummy; end do
73 
74 do j=jmax, 0, -1
75  read(24, fmt=trim(fmt4)) (maske(j,i), i=0,imax)
76 end do
77 
78 close(24, status='keep')
79 
80 !-------- Further stuff --------
81 
82 #if (GRID==0 || GRID==1)
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 #elif GRID==2
91 
92 stop ' topography1: GRID==2 not allowed for this application!'
93 
94 #endif
95 
96 do i=0, imax
97 do j=0, jmax
98 
99  zs(j,i) = zl0(j,i)
100  zb(j,i) = zl0(j,i)
101  zl(j,i) = zl0(j,i)
102 
103  if (maske(j,i) <= 1) then
104  maske(j,i) = 0
105  zs(j,i) = zs(j,i) + 2.0_dp*h_min
106  end if
107 
108  xi(i) = xi0 + real(i,dp)*dxi
109  eta(j) = eta0 + real(j,dp)*deta
110 
111  call geo_coord(phi(j,i), lambda(j,i), xi(i), eta(j))
112 
113  zm(j,i) = zb(j,i)
114  n_cts(j,i) = -1
115 
116  h_c(j,i) = zs(j,i)-zm(j,i)
117  h_t(j,i) = 0.0_dp
118 
119  dzs_dtau(j,i) = 0.0_dp
120  dzm_dtau(j,i) = 0.0_dp
121  dzb_dtau(j,i) = 0.0_dp
122  dzl_dtau(j,i) = 0.0_dp
123  dh_c_dtau(j,i) = 0.0_dp
124  dh_t_dtau(j,i) = 0.0_dp
125 
126 end do
127 end do
128 
129 !-------- Metric tensor, gradients of the topography --------
130 
131 call metric()
132 
133 #if TOPOGRAD==0
134 call topograd_1(dxi, deta, 1)
135 #elif TOPOGRAD==1
136 call topograd_2(dxi, deta, 1)
137 #endif
138 
139 !-------- Corresponding area of grid points --------
140 
141 do i=0, imax
142 do j=0, jmax
143  area(j,i) = sq_g11_g(j,i)*sq_g22_g(j,i)*dxi*deta
144 end do
145 end do
146 
147 end subroutine topography1
148 !