SPEED
MAKE_EXPL_SOURCE.f90 File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine make_expl_source (nn, ct, ww, dd, sxx, syy, szz, syz, szx, sxy, nl_expl, check_ne, check_dist_ne, length_cne, ielem, facsexpl, func_type, func_indx, func_data, nfdata, nf, t
 Computes the explosive moment tensor.
 

Function/Subroutine Documentation

◆ make_expl_source()

subroutine make_expl_source ( integer*4  nn,
real*8, dimension(nn)  ct,
real*8, dimension(nn)  ww,
real*8, dimension(nn,nn)  dd,
real*8, dimension(nn,nn,nn)  sxx,
real*8, dimension(nn,nn,nn)  syy,
real*8, dimension(nn,nn,nn)  szz,
real*8, dimension(nn,nn,nn)  syz,
real*8, dimension(nn,nn,nn)  szx,
real*8, dimension(nn,nn,nn)  sxy,
integer*4  nl_expl,
integer*4, dimension(length_cne,4)  check_ne,
real*8, dimension(length_cne,1)  check_dist_ne,
integer*4  length_cne,
integer*4  ielem,
real*8, dimension(nl_expl,6)  facsexpl,
integer*4, dimension(nf)  func_type,
integer*4, dimension(nf+1)  func_indx,
real*8, dimension(nfdata)  func_data,
integer*4  nfdata,
integer*4  nf,
  t 
)

Computes the explosive moment tensor.

Author
Ilario Mazzieri
Date
September, 2013
Version
1.0
Parameters
[in]nnnumber of 1-D GLL nodes
[in]xqGLL nodes
[in]wqGLL weights
[in]ddspectral derivatives matrix
[in]nl_explnumber of explosive loads
[in]length_cnelength of check_ns
[in]check_neinfo about sismic load
[in]check_dist_nedistance from epicenter
[in]ielemelement index
[in]facsexplseismic moment factor
[in]numberof fuctions
[in]func_typefunction types
[in]func_indxfunction indixes for function data
[in]func_datafunction data
[in]t_stresstime
[in]cs_nnzlength of cs
[in]csspectral connectivity vector
[in]tag_funcfunctin label
[in]nelem_locnumber of local elements
[in]local_el_numlocal element numeration
[in]nn_locnumber of local nodes
[in]local_n_numlocal node numeration
[in,out]sxxnodal values for the stress tensor
[in,out]syynodal values for the stress tensor
[in,out]szznodal values for the stress tensor
[in,out]syznodal values for the stress tensor
[in,out]szxnodal values for the stress tensor
[in,out]sxynodal values for the stress tensor

Definition at line 53 of file MAKE_EXPL_SOURCE.f90.

59 cs_nnz,cs,tag_func, &
60 nelem_loc, local_el_num,&
61 nn_loc, local_n_num)
62
63 implicit none
64
65 integer*4 :: nn, nelem_loc, nn_loc
66 integer*4 :: i,j,k,p,q,r
67 integer*4 :: length_cne,ielem,nl_expl,nfdata
68 integer*4 :: nf
69 integer*4 :: is,in,fn, ie
70 integer*4 :: cs_nnz
71
72 integer*4, dimension(nelem_loc) :: local_el_num
73 integer*4, dimension(nn_loc) :: local_n_num
74 integer*4, dimension(nf) :: tag_func
75 integer*4, dimension(nf) :: func_type
76 integer*4, dimension(nf+1) :: func_indx
77
78 integer*4, dimension(0:cs_nnz) :: cs
79 integer*4, dimension(length_cne,4) :: check_ne
80
81 real*8 :: get_func_value_sism
82 real*8 :: t_stress
83
84 real*8, dimension(nfdata) :: func_data
85 real*8, dimension(nn) :: ct,ww
86
87 real*8, dimension(nn,nn) :: dd
88 real*8, dimension(length_cne,1) :: check_dist_ne
89 real*8, dimension(nl_expl,6) :: facsexpl
90
91 real*8, dimension(nn,nn,nn) :: sxx,syy,szz,syz,szx,sxy
92
93! STRESS CALCULATION
94
95
96
97 if ((ielem.ge.check_ne(1,4)).and. (ielem.le.check_ne(length_cne,4))) then
98
99 do fn = 1,nf
100
101 do i = 1,length_cne
102 if (ielem.eq.check_ne(i,4)) then
103
104
105 do r = 1,nn
106 do q = 1,nn
107 do p = 1,nn
108
109
110 in = local_n_num(cs(cs(ielem -1) +is))
111
112 if ((tag_func(fn).eq.(check_ne(i,2))) .and.(check_ne(i,1).eq.in)) then
113
114
115 !STRESS RE-CALCULATION
116
117
118 sxx(p,q,r) = sxx(p,q,r) &
119 - get_func_value_sism(nf,func_type,func_indx,func_data, nfdata,&
120 fn,t_stress,check_dist_ne(i,1),0.d0) &
121 * facsexpl(check_ne(i,3),1)
122
123 syy(p,q,r) = syy(p,q,r) &
124 - get_func_value_sism(nf,func_type,func_indx,func_data, nfdata,&
125 fn,t_stress,check_dist_ne(i,1),0.d0) &
126 * facsexpl(check_ne(i,3),2)
127
128 szz(p,q,r) = szz(p,q,r) &
129 - get_func_value_sism(nf,func_type,func_indx,func_data, nfdata,&
130 fn,t_stress,check_dist_ne(i,1),0.d0) &
131 * facsexpl(check_ne(i,3),3)
132
133 syz(p,q,r) = syz(p,q,r) &
134 - get_func_value_sism(nf,func_type,func_indx,func_data, nfdata,&
135 fn,t_stress,check_dist_ne(i,1),0.d0) &
136 * facsexpl(check_ne(i,3),4)
137
138 szx(p,q,r) = szx(p,q,r) &
139 - get_func_value_sism(nf,func_type,func_indx,func_data, nfdata,&
140 fn,t_stress,check_dist_ne(i,1),0.d0) &
141 * facsexpl(check_ne(i,3),5)
142
143 sxy(p,q,r) = sxy(p,q,r) &
144 - get_func_value_sism(nf,func_type,func_indx,func_data, nfdata,&
145 fn,t_stress,check_dist_ne(i,1),0.d0) &
146 * facsexpl(check_ne(i,3),6)
147
148 endif ! if ((tag_func(fn).eq.(check_node_expl(i,2))) &...
149
150 enddo
151 enddo
152 enddo
153
154 endif !ielem
155 enddo !i
156 enddo !fn
157
158 endif !ielem
159
160
161 return
162
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.

References get_func_value_sism().

Here is the call graph for this function: