SPEED
MAKE_SEISMIC_MOMENT.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
52
53 subroutine make_seismic_moment(nn,ct,ww,dd,&
54 sxx,syy,szz,syz,szx,sxy,&
55 nl_sism,&
56 check_ns,check_dist_ns,&
57 length_cns,ielem,facsmom,&
58 tausmom,&
59 func_type,func_indx,func_data,nfdata,nf,t_stress,&
60 cs_nnz,cs,tag_func, &
61 nelem_loc, local_el_num,&
62 nn_loc, local_n_num)
63
64 implicit none
65
66 integer*4 :: nn, nelem_loc, nn_loc
67 integer*4 :: i,j,k,p,q,r
68 integer*4 :: cs_nnz
69 integer*4 :: is,in,fn,ie
70 integer*4 :: length_cns,ielem,nl_sism
71 integer*4 :: nf,nfdata
72
73 integer*4, dimension(nelem_loc) :: local_el_num
74 integer*4, dimension(nn_loc) :: local_n_num
75 integer*4, dimension(nf) :: func_type
76 integer*4, dimension(nf +1) :: func_indx
77 integer*4, dimension(nf) :: tag_func
78 integer*4, dimension(0:cs_nnz) :: cs
79
80 integer*4, dimension(length_cns,4) :: check_ns
81
82 real*8 :: t_stress
83 real*8 :: get_func_value_sism, val_sism
84
85 real*8, dimension(nn) :: ct,ww
86 real*8, dimension(nfdata) :: func_data
87
88 real*8, dimension(nn,nn) :: dd
89 real*8, dimension(length_cns,1) :: check_dist_ns
90 real*8, dimension(nl_sism,6) :: facsmom
91 real*8, dimension(nl_sism,1) :: tausmom
92
93 real*8, dimension(nn,nn,nn) :: sxx,syy,szz,syz,szx,sxy
94 logical :: exist
95
96 if ((ielem .ge. check_ns(1,4)) .and. (ielem .le. check_ns(length_cns,4))) then
97
98 do i = 1,length_cns
99
100 if (ielem .eq. check_ns(i,4)) then
101
102
103 do r = 1,nn
104 do q = 1,nn
105 do p = 1,nn
106
107 is = nn*nn*(r -1) +nn*(q -1) +p
108
109
110 in = local_n_num(cs(cs(ielem -1) +is))
111
112 do fn = 1,nf
113
114 if ((tag_func(fn) .eq. (check_ns(i,2))) .and. (check_ns(i,1) .eq. in)) then
115
116
117 !STRESS RE-CALCULATION
118
119 !t_stress = current time,
120 !check_dist_ns(i,1) = rupture time
121 !tausmom = rise time
122 !facsmom = seismic moment
123 val_sism = get_func_value_sism(nf,func_type,func_indx,func_data, nfdata, &
124 fn,t_stress,check_dist_ns(i,1), &
125 tausmom(check_ns(i,3),1))
126
127 !write(*,*) VAL_SISM
128 !read(*,*)
129
130! sxx(p,q,r) = sxx(p,q,r) - GET_FUNC_VALUE_SISM(nf,func_type,func_indx,func_data, nfdata, &
131! fn,t_stress,check_dist_ns(i,1), &
132! tausmom(check_ns(i,3),1)) * facsmom(check_ns(i,3),1)
133
134! syy(p,q,r) = syy(p,q,r) - GET_FUNC_VALUE_SISM(nf,func_type,func_indx,func_data, nfdata,&
135! fn,t_stress,check_dist_ns(i,1), &
136! tausmom(check_ns(i,3),1)) * facsmom(check_ns(i,3),2)
137
138! szz(p,q,r) = szz(p,q,r) - GET_FUNC_VALUE_SISM(nf,func_type,func_indx,func_data, nfdata,&
139! fn,t_stress,check_dist_ns(i,1), &
140! tausmom(check_ns(i,3),1)) * facsmom(check_ns(i,3),3)
141
142! syz(p,q,r) = syz(p,q,r) - GET_FUNC_VALUE_SISM(nf,func_type,func_indx,func_data, nfdata,&
143! fn,t_stress,check_dist_ns(i,1), &
144! tausmom(check_ns(i,3),1)) * facsmom(check_ns(i,3),4)
145
146! szx(p,q,r) = szx(p,q,r) - GET_FUNC_VALUE_SISM(nf,func_type,func_indx,func_data, nfdata,&
147! fn,t_stress,check_dist_ns(i,1), &
148! tausmom(check_ns(i,3),1)) * facsmom(check_ns(i,3),5)
149
150! sxy(p,q,r) = sxy(p,q,r) - GET_FUNC_VALUE_SISM(nf,func_type,func_indx,func_data, nfdata,&
151! fn,t_stress,check_dist_ns(i,1), &
152! tausmom(check_ns(i,3),1)) * facsmom(check_ns(i,3),6)
153
154 sxx(p,q,r) = sxx(p,q,r) - val_sism * facsmom(check_ns(i,3),1)
155
156 syy(p,q,r) = syy(p,q,r) - val_sism * facsmom(check_ns(i,3),2)
157
158 szz(p,q,r) = szz(p,q,r) - val_sism * facsmom(check_ns(i,3),3)
159
160 syz(p,q,r) = syz(p,q,r) - val_sism * facsmom(check_ns(i,3),4)
161
162 szx(p,q,r) = szx(p,q,r) - val_sism * facsmom(check_ns(i,3),5)
163
164 sxy(p,q,r) = sxy(p,q,r) - val_sism * facsmom(check_ns(i,3),6)
165
166
167 ! inquire(file="sft_th_debug.txt", exist=exist)
168 ! if (exist) then
169 ! open(82, file="sft_th_debug.txt", status="old", position="append", action="write")
170 ! else
171 ! open(82, file="sft_th_debug.txt", status="new", action="write")
172 ! end if
173 ! write(82,*) 'time', t_stress, 'funct_value', GET_FUNC_VALUE_SISM(nf,func_type,func_indx,&
174 ! func_data, nfdata, fn,t_stress,check_dist_ns(i,1), &
175 ! tausmom(check_ns(i,3),1))
176 ! close(82)
177
178
179 endif ! if ((tag_func(fn).eq.(check_node_sism(i,2))) &...
180
181 enddo !fn = 1,nf
182
183 enddo !p
184 enddo !q
185 enddo !r
186
187 endif !ielem
188
189 enddo !i
190
191
192 endif !ielem
193
194 return
195
196 end subroutine make_seismic_moment
197
real *8 function get_func_value_sism(nb_fnc, type_fnc, ind_fnc, data_fnc, nb_data_fnc, id_fnc,
Computes time evolution function for seismic source.
subroutine make_seismic_moment(nn, ct, ww, dd, sxx, syy, szz, syz, szx, sxy, nl_sism, check_ns, check_dist_ns, length_cns, ielem, facsmom, tausmom, func_type, func_indx, func_data, nfdata, nf, t
Creates the seismic moment tensor.