SICOPOLIS V3.2
 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-2016 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 use sico_vars
46 
47 implicit none
48 
49 integer(i4b), intent(in) :: n_switch
50 real(dp), intent(in) :: dxi, deta
51 
52 integer(i4b) :: i, j
53 real(dp) :: dxi_inv, deta_inv
54 real(dp) :: zs_help(0:jmax,0:imax), zm_help(0:jmax,0:imax), &
55  zb_help(0:jmax,0:imax)
56 
57 dxi_inv = 1.0_dp/dxi
58 deta_inv = 1.0_dp/deta
59 
60 !-------- Distinguish between old and new topography data --------
61 
62 if (n_switch == 1) then
63 
64  do i=0, imax
65  do j=0, jmax
66  zs_help(j,i) = zs(j,i)
67  zm_help(j,i) = zm(j,i)
68  zb_help(j,i) = zb(j,i)
69  end do
70  end do
71 
72 else if (n_switch == 2) then
73 
74  do i=0, imax
75  do j=0, jmax
76  zs_help(j,i) = zs_neu(j,i)
77  zm_help(j,i) = zm_neu(j,i)
78  zb_help(j,i) = zb_neu(j,i)
79  end do
80  end do
81 
82 else
83  stop ' topograd_1: Wrong value for n_switch!'
84 end if
85 
86 !-------- Topography gradients on the staggered grid --------
87 
88 ! ------ x-derivatives
89 
90 do i=0, imax-1
91 do j=0, jmax
92  dzs_dxi(j,i) = (zs_help(j,i+1)-zs_help(j,i))*dxi_inv &
93  *insq_g11_sgx(j,i)
94  dzm_dxi(j,i) = (zm_help(j,i+1)-zm_help(j,i))*dxi_inv &
95  *insq_g11_sgx(j,i)
96  dzb_dxi(j,i) = (zb_help(j,i+1)-zb_help(j,i))*dxi_inv &
97  *insq_g11_sgx(j,i)
98  dh_c_dxi(j,i) = dzs_dxi(j,i)-dzm_dxi(j,i)
99  dh_t_dxi(j,i) = dzm_dxi(j,i)-dzb_dxi(j,i)
100 end do
101 end do
102 
103 ! ------ y-derivatives
104 
105 do i=0, imax
106 do j=0, jmax-1
107  dzs_deta(j,i) = (zs_help(j+1,i)-zs_help(j,i))*deta_inv &
108  *insq_g22_sgy(j,i)
109  dzm_deta(j,i) = (zm_help(j+1,i)-zm_help(j,i))*deta_inv &
110  *insq_g22_sgy(j,i)
111  dzb_deta(j,i) = (zb_help(j+1,i)-zb_help(j,i))*deta_inv &
112  *insq_g22_sgy(j,i)
113  dh_c_deta(j,i) = dzs_deta(j,i)-dzm_deta(j,i)
114  dh_t_deta(j,i) = dzm_deta(j,i)-dzb_deta(j,i)
115 end do
116 end do
117 
118 !-------- Topography gradients on the grid points --------
119 
120 ! ------ x-derivatives
121 
122 do i=1, imax-1
123 do j=0, jmax
124  dzs_dxi_g(j,i) = (zs_help(j,i+1)-zs_help(j,i-1))*0.5_dp*dxi_inv &
125  *insq_g11_g(j,i)
126  dzm_dxi_g(j,i) = (zm_help(j,i+1)-zm_help(j,i-1))*0.5_dp*dxi_inv &
127  *insq_g11_g(j,i)
128  dzb_dxi_g(j,i) = (zb_help(j,i+1)-zb_help(j,i-1))*0.5_dp*dxi_inv &
129  *insq_g11_g(j,i)
130  dh_c_dxi_g(j,i) = dzs_dxi_g(j,i)-dzm_dxi_g(j,i)
131  dh_t_dxi_g(j,i) = dzm_dxi_g(j,i)-dzb_dxi_g(j,i)
132 end do
133 end do
134 
135 do j=0, jmax
136  dzs_dxi_g(j,0) = (zs_help(j,1)-zs_help(j,0))*dxi_inv &
137  *insq_g11_g(j,0)
138  dzm_dxi_g(j,0) = (zm_help(j,1)-zm_help(j,0))*dxi_inv &
139  *insq_g11_g(j,0)
140  dzb_dxi_g(j,0) = (zb_help(j,1)-zb_help(j,0))*dxi_inv &
141  *insq_g11_g(j,0)
142  dh_c_dxi_g(j,0) = dzs_dxi_g(j,0)-dzm_dxi_g(j,0)
143  dh_t_dxi_g(j,0) = dzm_dxi_g(j,0)-dzb_dxi_g(j,0)
144  dzs_dxi_g(j,imax) = (zs_help(j,imax)-zs_help(j,imax-1)) &
145  *dxi_inv &
146  *insq_g11_g(j,imax)
147  dzm_dxi_g(j,imax) = (zm_help(j,imax)-zm_help(j,imax-1)) &
148  *dxi_inv &
149  *insq_g11_g(j,imax)
150  dzb_dxi_g(j,imax) = (zb_help(j,imax)-zb_help(j,imax-1)) &
151  *dxi_inv &
152  *insq_g11_g(j,imax)
153  dh_c_dxi_g(j,imax) = dzs_dxi_g(j,imax)-dzm_dxi_g(j,imax)
154  dh_t_dxi_g(j,imax) = dzm_dxi_g(j,imax)-dzb_dxi_g(j,imax)
155 end do
156 
157 ! ------ y-derivatives
158 
159 do i=0, imax
160 do j=1, jmax-1
161  dzs_deta_g(j,i) = (zs_help(j+1,i)-zs_help(j-1,i)) &
162  *0.5_dp*deta_inv &
163  *insq_g22_g(j,i)
164  dzm_deta_g(j,i) = (zm_help(j+1,i)-zm_help(j-1,i)) &
165  *0.5_dp*deta_inv &
166  *insq_g22_g(j,i)
167  dzb_deta_g(j,i) = (zb_help(j+1,i)-zb_help(j-1,i)) &
168  *0.5_dp*deta_inv &
169  *insq_g22_g(j,i)
170  dh_c_deta_g(j,i) = dzs_deta_g(j,i)-dzm_deta_g(j,i)
171  dh_t_deta_g(j,i) = dzm_deta_g(j,i)-dzb_deta_g(j,i)
172 end do
173 end do
174 
175 do i=0, imax
176  dzs_deta_g(0,i) = (zs_help(1,i)-zs_help(0,i))*deta_inv &
177  *insq_g22_g(0,i)
178  dzm_deta_g(0,i) = (zm_help(1,i)-zm_help(0,i))*deta_inv &
179  *insq_g22_g(0,i)
180  dzb_deta_g(0,i) = (zb_help(1,i)-zb_help(0,i))*deta_inv &
181  *insq_g22_g(0,i)
182  dh_c_deta_g(0,i) = dzs_deta_g(0,i)-dzm_deta_g(0,i)
183  dh_t_deta_g(0,i) = dzm_deta_g(0,i)-dzb_deta_g(0,i)
184  dzs_deta_g(jmax,i) = (zs_help(jmax,i)-zs_help(jmax-1,i)) &
185  *deta_inv &
186  *insq_g22_g(jmax,i)
187  dzm_deta_g(jmax,i) = (zm_help(jmax,i)-zm_help(jmax-1,i)) &
188  *deta_inv &
189  *insq_g22_g(jmax,i)
190  dzb_deta_g(jmax,i) = (zb_help(jmax,i)-zb_help(jmax-1,i)) &
191  *deta_inv &
192  *insq_g22_g(jmax,i)
193  dh_c_deta_g(jmax,i) = dzs_deta_g(jmax,i)-dzm_deta_g(jmax,i)
194  dh_t_deta_g(jmax,i) = dzm_deta_g(jmax,i)-dzb_deta_g(jmax,i)
195 end do
196 
197 end subroutine topograd_1
198 !
Declarations of kind types for SICOPOLIS.
Definition: sico_types.F90:35
subroutine topograd_1(dxi, deta, n_switch)
Calculation of topography gradients on the staggered grid and on the grid points (the latter by secon...
Definition: topograd_1.F90:41
Declarations of global variables for SICOPOLIS (for the ANT domain).
Definition: sico_vars.F90:35
Declarations of global variables for SICOPOLIS.