SPEED
UPDATE_DAMPING_TENSOR.f90
Go to the documentation of this file.
1! Copyright (C) 2012 The SPEED FOUNDATION
2! Author: Ilario Mazzieri
3!
4! This file is part of SPEED.
5!
6! SPEED is free software; you can redistribute it and/or modify it
7! under the terms of the GNU Affero General Public License as
8! published by the Free Software Foundation, either version 3 of the
9! License, or (at your option) any later version.
10!
11! SPEED is distributed in the hope that it will be useful, but
12! WITHOUT ANY WARRANTY; without even the implied warranty of
13! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14! Affero General Public License for more details.
15!
16! You should have received a copy of the GNU Affero General Public License
17! along with SPEED. If not, see <http://www.gnu.org/licenses/>.
18
29
30
31 subroutine update_damping_tensor(nnod_loc, N_SLS, frequency_range, &
32 strain_visc, strain, deltat)
33
34
35 implicit none
36
37
38 integer*4 :: nnod_loc, N_SLS, i
39
40
41 real*8 :: deltat, frequency_range(n_sls), strain_visc(6*nnod_loc,n_sls), &
42 strain(6*nnod_loc), temp(6*nnod_loc, n_sls)
43
44
45
46 do i = 1, n_sls
47 !EXPLICIT EULER
48 !temp(:,i) = strain_visc(:,i) + deltat*frequency_range(i)*(strain - strain_visc(:,i))
49 !EXPLICIT HEUN
50 temp(:,i) = strain_visc(:,i) + &
51 (deltat*frequency_range(i) - 0.5*(deltat*frequency_range(i))**2)*(strain - strain_visc(:,i))
52
53
54 enddo
55
56 strain_visc = temp;
57
58
59 return
60
61 end subroutine update_damping_tensor
subroutine update_damping_tensor(nnod_loc, n_sls, frequency_range,
Update the damping tensor with the Heun method.