SICOPOLIS V3.2
 All Classes Files Functions Variables Macros
calc_vz_ssa.F90
Go to the documentation of this file.
1 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 !
3 ! Subroutine : c a l c _ v z _ s s a
4 !
5 !> @file
6 !!
7 !! Computation of the vertical velocity vz in the shallow shelf approximation.
8 !!
9 !! @section Copyright
10 !!
11 !! Copyright 2009-2016 Ralf Greve, Tatsuru Sato
12 !!
13 !! @section License
14 !!
15 !! This file is part of SICOPOLIS.
16 !!
17 !! SICOPOLIS is free software: you can redistribute it and/or modify
18 !! it under the terms of the GNU General Public License as published by
19 !! the Free Software Foundation, either version 3 of the License, or
20 !! (at your option) any later version.
21 !!
22 !! SICOPOLIS is distributed in the hope that it will be useful,
23 !! but WITHOUT ANY WARRANTY; without even the implied warranty of
24 !! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
25 !! GNU General Public License for more details.
26 !!
27 !! You should have received a copy of the GNU General Public License
28 !! along with SICOPOLIS. If not, see <http://www.gnu.org/licenses/>.
29 !<
30 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
31 
32 !-------------------------------------------------------------------------------
33 !> Computation of the vertical velocity vz in the shallow shelf approximation.
34 !<------------------------------------------------------------------------------
35 subroutine calc_vz_ssa(z_sl, dxi, deta, dzeta_c)
36 
37 use sico_types
39 use sico_vars
40 
41 implicit none
42 
43 real(dp), intent(in) :: z_sl
44 real(dp), intent(in) :: dxi, deta, dzeta_c
45 
46 integer(i4b) :: i, j, kt, kc
47 real(dp), dimension(0:JMAX,0:IMAX) :: vz_sl
48 real(dp), dimension(0:KCMAX) :: zeta_c_sgz, eaz_c_sgz, eaz_c_quotient_sgz
49 real(dp) :: dxi_inv, deta_inv
50 real(dp) :: dvx_dxi, dvy_deta
51 
52 dxi_inv = 1.0_dp/dxi
53 deta_inv = 1.0_dp/deta
54 
55 !-------- Abbreviations --------
56 
57 do kc=0, kcmax-1
58 
59  zeta_c_sgz(kc) = (kc+0.5_dp)*dzeta_c
60 
61  if (flag_aa_nonzero) then
62  eaz_c_sgz(kc) = exp(aa*zeta_c_sgz(kc))
63  eaz_c_quotient_sgz(kc) = (eaz_c_sgz(kc)-1.0_dp)/(ea-1.0_dp)
64  else
65  eaz_c_sgz(kc) = 1.0_dp
66  eaz_c_quotient_sgz(kc) = zeta_c_sgz(kc)
67  end if
68 
69 end do
70 
71 !-------- Computation of vz --------
72 
73 do i=1, imax-1
74 do j=1, jmax-1
75 
76  if (maske(j,i)==3_i2b) then ! floating ice
77 
78 ! ------ Derivatives of the horizontal velocity
79 
80  dvx_dxi = insq_g11_g(j,i)*insq_g22_g(j,i) &
81  *(vx_m(j,i)*sq_g22_sgx(j,i)-vx_m(j,i-1)*sq_g22_sgx(j,i-1)) &
82  *dxi_inv
83 
84  dvy_deta = insq_g11_g(j,i)*insq_g22_g(j,i) &
85  *(vy_m(j,i)*sq_g11_sgy(j,i)-vy_m(j-1,i)*sq_g11_sgy(j-1,i)) &
86  *deta_inv
87 
88 ! ------ Basal velocity vz_b
89 
90  vz_b(j,i) = 0.5_dp*(vx_m(j,i)+vx_m(j,i-1))*dzb_dxi_g(j,i) &
91  +0.5_dp*(vy_m(j,i)+vy_m(j-1,i))*dzb_deta_g(j,i) &
92  +dzb_dtau(j,i)-q_b_tot(j,i)
93  ! kinematic boundary condition at the ice base
94 
95 ! ------ Velocity at sea level vz_sl
96 
97  vz_sl(j,i) = vz_b(j,i) - (z_sl-zb(j,i))*(dvx_dxi+dvy_deta)
98 
99 ! ------ Surface velocity vz_s
100 
101  vz_s(j,i) = vz_sl(j,i) - (zs(j,i)-z_sl)*(dvx_dxi+dvy_deta)
102 
103 ! ------ Velocity vz_m at the interface between
104 ! the upper (kc) and the lower (kt) domain
105 
106  if ((n_cts(j,i) == -1).or.(n_cts(j,i) == 0)) then
107  ! cold ice base, temperate ice base
108 
109  vz_m(j,i) = vz_b(j,i)
110 
111  else ! n_cts(j,i) == 1, temperate ice layer
112 
113  vz_m(j,i) = vz_b(j,i) - h_t(j,i)*(dvx_dxi+dvy_deta)
114 
115  end if
116 
117 ! ------ 3-D velocity vz_c and vz_t
118 
119  do kc=0, kcmax-1
120  vz_c(kc,j,i) = vz_sl(j,i) &
121  -(zm(j,i)+eaz_c_quotient_sgz(kc)*h_c(j,i)-z_sl) &
122  *(dvx_dxi+dvy_deta)
123  end do
124 
125  if ((n_cts(j,i) == -1).or.(n_cts(j,i) == 0)) then
126  ! cold ice base, temperate ice base
127 
128  do kt=0, ktmax-1
129  vz_t(kt,j,i) = vz_b(j,i)
130  end do
131 
132  else ! n_cts(j,i) == 1, temperate ice layer
133 
134  do kt=0, ktmax-1
135  vz_t(kt,j,i) = vz_sl(j,i) &
136  -(zb(j,i) &
137  +0.5_dp*(zeta_t(kt)+zeta_t(kt+1))*h_t(j,i) &
138  -z_sl) &
139  *(dvx_dxi+dvy_deta)
140  end do
141 
142  end if
143 
144  end if
145 
146 end do
147 end do
148 
149 end subroutine calc_vz_ssa
150 !
subroutine calc_vz_ssa(z_sl, dxi, deta, dzeta_c)
Computation of the vertical velocity vz in the shallow shelf approximation.
Definition: calc_vz_ssa.F90:35
Declarations of kind types for SICOPOLIS.
Definition: sico_types.F90:35
Declarations of global variables for SICOPOLIS (for the ANT domain).
Definition: sico_vars.F90:35
Declarations of global variables for SICOPOLIS.