SICOPOLIS V3.1
 All Classes Files Functions Variables Macros
calc_elra.F90
Go to the documentation of this file.
1 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 !
3 ! Subroutine : c a l c _ e l r a
4 !
5 !> @file
6 !!
7 !! Computation of the isostatic steady-state displacement of the lithosphere
8 !! (wss) for the ELRA model.
9 !!
10 !! @section Copyright
11 !!
12 !! Copyright 2009-2013 Ralf Greve, Sascha Knell
13 !!
14 !! @section License
15 !!
16 !! This file is part of SICOPOLIS.
17 !!
18 !! SICOPOLIS is free software: you can redistribute it and/or modify
19 !! it under the terms of the GNU General Public License as published by
20 !! the Free Software Foundation, either version 3 of the License, or
21 !! (at your option) any later version.
22 !!
23 !! SICOPOLIS is distributed in the hope that it will be useful,
24 !! but WITHOUT ANY WARRANTY; without even the implied warranty of
25 !! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
26 !! GNU General Public License for more details.
27 !!
28 !! You should have received a copy of the GNU General Public License
29 !! along with SICOPOLIS. If not, see <http://www.gnu.org/licenses/>.
30 !<
31 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
32 
33 !-------------------------------------------------------------------------------
34 !> Computation of the isostatic steady-state displacement of the lithosphere
35 !! (wss) for the ELRA model.
36 !<------------------------------------------------------------------------------
37 subroutine calc_elra(z_sl, dxi, deta)
38 
39 use sico_types
41 
42 implicit none
43 
44 real(dp), intent(in) :: z_sl, dxi, deta
45 
46 integer(i4b) :: i, j, ir, jr, il, jl, n
47 integer(i4b) :: ir_max, jr_max, min_imax_jmax
48 integer(i4b) :: il_begin, il_end, jl_begin, jl_end
49 real(dp) :: l_r, l_r_inv, rho_g, rho_sw_g, fac_wss
50 real(dp) :: ra_max
51 real(dp) :: dxi_inv, deta_inv
52 real(dp) :: kei_r_incr_inv
53 real(dp), allocatable, dimension(:,:) :: f_0
54 
55 real(dp), parameter :: r_infl = 8.0_dp
56  ! Radius of non-locality influence of the elastic
57  ! lithosphere plate (in units of l_r)
58 
59 !-------- Initialisations --------
60 
61 l_r = (flex_rig_l/(rho_a*g))**0.25_dp
62 
63 ra_max = r_infl*l_r ! Radius of non-locality influence (in m)
64 
65 l_r_inv = 1.0_dp/l_r
66 rho_g = rho*g
67 rho_sw_g = rho_sw*g
68 
69 fac_wss = l_r**2.0_dp/(2.0_dp*pi*flex_rig_l)
70 
71 dxi_inv = 1.0_dp/dxi
72 deta_inv = 1.0_dp/deta
73 
74 kei_r_incr_inv = 1.0_dp/kei_r_incr
75 
76 ir_max = floor(ra_max*dxi_inv)
77 jr_max = floor(ra_max*deta_inv)
78  ! Radius of non-locality influence (in grid points)
79 
80 min_imax_jmax = min(imax, jmax)
81 ir_max = min(ir_max, min_imax_jmax)
82 jr_max = min(jr_max, min_imax_jmax)
83  ! ir_max and jr_max constrained by size of domain
84 
85 il_begin = -ir_max
86 il_end = ir_max+imax
87 jl_begin = -jr_max
88 jl_end = jr_max+jmax
89 
90 !-------- Ice/water load --------
91 
92 allocate(f_0(jl_begin:jl_end, il_begin:il_end))
93 
94 do il=il_begin, il_end
95 do jl=jl_begin, jl_end
96 
97  i = min(max(il, 0), imax)
98  j = min(max(jl, 0), jmax)
99 
100  if (maske(j,i)==0) then
101  f_0(jl,il) = rho_g * area(j,i) * (h_c(j,i) + h_t(j,i))
102  else if (maske(j,i)==1) then
103  f_0(jl,il) = 0.0_dp
104  else ! (maske(j,i)>=2)
105  f_0(jl,il) = rho_sw_g * area(j,i) * z_sl
106  ! Water load relative to the present sea-level stand (0 m)
107  ! -> can be positive or negative
108  end if
109 
110 end do
111 end do
112 
113 !-------- Computation of the steady-state displacement of the lithosphere
114 ! (wss, positive downward) --------
115 
116 do i=0, imax
117 do j=0, jmax
118 
119  wss(j,i) = 0.0_dp
120 
121  do ir=-ir_max, ir_max
122  do jr=-jr_max, jr_max
123 
124 #if GRID==0
125  n = nint( dist_dxdy(jr,ir)*l_r_inv*kei_r_incr_inv )
126 #elif GRID==1
127  n = nint( sq_g11_g(j,i)*dist_dxdy(jr,ir)*l_r_inv*kei_r_incr_inv )
128  ! The correction with g11 only is OK because g11=g22 for this case
129 #elif GRID==2
130  n = nint( dist_dxdy(jr,ir)*l_r_inv*kei_r_incr_inv )
131  ! The distortion due to g11 and g22 is already accounted for in
132  ! the variable dist_dxdy(jr,ir)
133 #endif
134 
135  wss(j,i) = wss(j,i) - fac_wss*f_0(jr+j,ir+i)*kei(n)
136 
137  end do
138  end do
139 
140 end do
141 end do
142 
143 deallocate (f_0)
144 
145 end subroutine calc_elra
146 !