SICOPOLIS V3.0
 All Classes Files Functions Variables Macros
topograd_1.F90
Go to the documentation of this file.
1 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 !
3 ! Subroutine : t o p o g r a d _ 1
4 !
5 !> @file
6 !!
7 !! Calculation of topography gradients on the staggered grid and on the grid
8 !! points (the latter by second-order discretization).
9 !! Length rescaling with the corresponding components of the metric tensor is
10 !! already included.
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 !> Calculation of topography gradients on the staggered grid and on the grid
37 !! points (the latter by second-order discretization).
38 !! Length rescaling with the corresponding components of the metric tensor is
39 !! already included.
40 !<------------------------------------------------------------------------------
41 subroutine topograd_1(dxi, deta, n_switch)
42 
43 use sico_types
45 
46 implicit none
47 integer(i4b) :: i, j
48 integer(i4b) :: n_switch
49 real(dp) :: dxi, deta, dxi_inv, deta_inv
50 real(dp) :: zs_help(0:jmax,0:imax), zm_help(0:jmax,0:imax), &
51  zb_help(0:jmax,0:imax)
52 
53 dxi_inv = 1.0_dp/dxi
54 deta_inv = 1.0_dp/deta
55 
56 !-------- Distinguish between old and new topography data --------
57 
58 if (n_switch.eq.1) then
59 
60  do i=0, imax
61  do j=0, jmax
62  zs_help(j,i) = zs(j,i)
63  zm_help(j,i) = zm(j,i)
64  zb_help(j,i) = zb(j,i)
65  end do
66  end do
67 
68 else if (n_switch.eq.2) then
69 
70  do i=0, imax
71  do j=0, jmax
72  zs_help(j,i) = zs_neu(j,i)
73  zm_help(j,i) = zm_neu(j,i)
74  zb_help(j,i) = zb_neu(j,i)
75  end do
76  end do
77 
78 else
79  stop ' topograd_1: Wrong value for n_switch!'
80 end if
81 
82 !-------- Topography gradients on the staggered grid --------
83 
84 ! ------ x-derivatives
85 
86 do i=0, imax-1
87 do j=0, jmax
88  dzs_dxi(j,i) = (zs_help(j,i+1)-zs_help(j,i))*dxi_inv &
89  *insq_g11_sgx(j,i)
90  dzm_dxi(j,i) = (zm_help(j,i+1)-zm_help(j,i))*dxi_inv &
91  *insq_g11_sgx(j,i)
92  dzb_dxi(j,i) = (zb_help(j,i+1)-zb_help(j,i))*dxi_inv &
93  *insq_g11_sgx(j,i)
94  dh_c_dxi(j,i) = dzs_dxi(j,i)-dzm_dxi(j,i)
95  dh_t_dxi(j,i) = dzm_dxi(j,i)-dzb_dxi(j,i)
96 end do
97 end do
98 
99 ! ------ y-derivatives
100 
101 do i=0, imax
102 do j=0, jmax-1
103  dzs_deta(j,i) = (zs_help(j+1,i)-zs_help(j,i))*deta_inv &
104  *insq_g22_sgy(j,i)
105  dzm_deta(j,i) = (zm_help(j+1,i)-zm_help(j,i))*deta_inv &
106  *insq_g22_sgy(j,i)
107  dzb_deta(j,i) = (zb_help(j+1,i)-zb_help(j,i))*deta_inv &
108  *insq_g22_sgy(j,i)
109  dh_c_deta(j,i) = dzs_deta(j,i)-dzm_deta(j,i)
110  dh_t_deta(j,i) = dzm_deta(j,i)-dzb_deta(j,i)
111 end do
112 end do
113 
114 !-------- Topography gradients on the grid points --------
115 
116 ! ------ x-derivatives
117 
118 do i=1, imax-1
119 do j=0, jmax
120  dzs_dxi_g(j,i) = (zs_help(j,i+1)-zs_help(j,i-1))*0.5_dp*dxi_inv &
121  *insq_g11_g(j,i)
122  dzm_dxi_g(j,i) = (zm_help(j,i+1)-zm_help(j,i-1))*0.5_dp*dxi_inv &
123  *insq_g11_g(j,i)
124  dzb_dxi_g(j,i) = (zb_help(j,i+1)-zb_help(j,i-1))*0.5_dp*dxi_inv &
125  *insq_g11_g(j,i)
126  dh_c_dxi_g(j,i) = dzs_dxi_g(j,i)-dzm_dxi_g(j,i)
127  dh_t_dxi_g(j,i) = dzm_dxi_g(j,i)-dzb_dxi_g(j,i)
128 end do
129 end do
130 
131 do j=0, jmax
132  dzs_dxi_g(j,0) = (zs_help(j,1)-zs_help(j,0))*dxi_inv &
133  *insq_g11_g(j,0)
134  dzm_dxi_g(j,0) = (zm_help(j,1)-zm_help(j,0))*dxi_inv &
135  *insq_g11_g(j,0)
136  dzb_dxi_g(j,0) = (zb_help(j,1)-zb_help(j,0))*dxi_inv &
137  *insq_g11_g(j,0)
138  dh_c_dxi_g(j,0) = dzs_dxi_g(j,0)-dzm_dxi_g(j,0)
139  dh_t_dxi_g(j,0) = dzm_dxi_g(j,0)-dzb_dxi_g(j,0)
140  dzs_dxi_g(j,imax) = (zs_help(j,imax)-zs_help(j,imax-1)) &
141  *dxi_inv &
142  *insq_g11_g(j,imax)
143  dzm_dxi_g(j,imax) = (zm_help(j,imax)-zm_help(j,imax-1)) &
144  *dxi_inv &
145  *insq_g11_g(j,imax)
146  dzb_dxi_g(j,imax) = (zb_help(j,imax)-zb_help(j,imax-1)) &
147  *dxi_inv &
148  *insq_g11_g(j,imax)
149  dh_c_dxi_g(j,imax) = dzs_dxi_g(j,imax)-dzm_dxi_g(j,imax)
150  dh_t_dxi_g(j,imax) = dzm_dxi_g(j,imax)-dzb_dxi_g(j,imax)
151 end do
152 
153 ! ------ y-derivatives
154 
155 do i=0, imax
156 do j=1, jmax-1
157  dzs_deta_g(j,i) = (zs_help(j+1,i)-zs_help(j-1,i)) &
158  *0.5_dp*deta_inv &
159  *insq_g22_g(j,i)
160  dzm_deta_g(j,i) = (zm_help(j+1,i)-zm_help(j-1,i)) &
161  *0.5_dp*deta_inv &
162  *insq_g22_g(j,i)
163  dzb_deta_g(j,i) = (zb_help(j+1,i)-zb_help(j-1,i)) &
164  *0.5_dp*deta_inv &
165  *insq_g22_g(j,i)
166  dh_c_deta_g(j,i) = dzs_deta_g(j,i)-dzm_deta_g(j,i)
167  dh_t_deta_g(j,i) = dzm_deta_g(j,i)-dzb_deta_g(j,i)
168 end do
169 end do
170 
171 do i=0, imax
172  dzs_deta_g(0,i) = (zs_help(1,i)-zs_help(0,i))*deta_inv &
173  *insq_g22_g(0,i)
174  dzm_deta_g(0,i) = (zm_help(1,i)-zm_help(0,i))*deta_inv &
175  *insq_g22_g(0,i)
176  dzb_deta_g(0,i) = (zb_help(1,i)-zb_help(0,i))*deta_inv &
177  *insq_g22_g(0,i)
178  dh_c_deta_g(0,i) = dzs_deta_g(0,i)-dzm_deta_g(0,i)
179  dh_t_deta_g(0,i) = dzm_deta_g(0,i)-dzb_deta_g(0,i)
180  dzs_deta_g(jmax,i) = (zs_help(jmax,i)-zs_help(jmax-1,i)) &
181  *deta_inv &
182  *insq_g22_g(jmax,i)
183  dzm_deta_g(jmax,i) = (zm_help(jmax,i)-zm_help(jmax-1,i)) &
184  *deta_inv &
185  *insq_g22_g(jmax,i)
186  dzb_deta_g(jmax,i) = (zb_help(jmax,i)-zb_help(jmax-1,i)) &
187  *deta_inv &
188  *insq_g22_g(jmax,i)
189  dh_c_deta_g(jmax,i) = dzs_deta_g(jmax,i)-dzm_deta_g(jmax,i)
190  dh_t_deta_g(jmax,i) = dzm_deta_g(jmax,i)-dzb_deta_g(jmax,i)
191 end do
192 
193 end subroutine topograd_1
194 !