SPEED
CHECK_EXPL.f90 File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine check_expl (cs_nnz, cs, nm, tm, sd, nl_expl, num_ne, max_num_ne, sour_ne, dist_sour_ne, check_ne, check_dist_ne, length_cne, fun_expl, nf, tag_func, val_expl, nn_loc, local_n_num)
 Fills array check_ns for explosive force.
 

Function/Subroutine Documentation

◆ check_expl()

subroutine check_expl ( integer*4  cs_nnz,
integer*4, dimension(0:cs_nnz)  cs,
integer*4  nm,
integer*4, dimension(nm)  tm,
integer*4, dimension(nm)  sd,
integer*4  nl_expl,
integer*4, dimension(nl_expl)  num_ne,
integer*4  max_num_ne,
integer*4, dimension(max_num_ne,nl_expl)  sour_ne,
real*8, dimension(max_num_ne,nl_expl)  dist_sour_ne,
integer*4, dimension(length_cne,4)  check_ne,
real*8, dimension(length_cne,1)  check_dist_ne,
integer*4  length_cne,
integer*4, dimension(nl_expl)  fun_expl,
integer*4  nf,
integer*4, dimension(nf)  tag_func,
real*8, dimension(nl_expl,20)  val_expl,
integer*4  nn_loc,
integer*4, dimension(nn_loc)  local_n_num 
)

Fills array check_ns for explosive force.

Author
Ilario Mazzieri
Date
September, 2013
Version
1.0
Parameters
[in]cs_nnzlength of cs
[in]csspectral connectivity vector
[in]nmnumber of materials
[in]tmlabel for materias
[in]sdpolynomial degree vector
[in]nl_explnumber of explosive nodes
[in]num_nenumber of explosive nodes
[in]max_num_nemax number of explosive nodes
[in]sour_neexplosive source
[in]dist_sour_nedistance of the node from the explosive source
[in]length_cnelength check explosive nodes (useless)
[in]fun_explfunction associeted with the explosive load
[in]nfnumber of functions
[in]tag_funclabel for explosive functions
[in]val_explvalues for the explosive loads
[in]nn_locnumber of local nodes
[in]local_n_numlocal numeration vector
[out]check_necheck_ne(i,1) explosive source node, check_ne(i,2) explosive function for the node i, check_ne(i,3) i for explosive load, check_ne(i,4) number of the element
[out]check_dist_necheck_dist_ne(i,1) = dist_sour_ne / val_expl

Definition at line 47 of file CHECK_EXPL.f90.

56
57 implicit none
58
59 integer*4 :: cs_nnz,nm,ne,nl_expl,nf, nn_loc
60 integer*4 :: max_num_ne,im,ie,iexpl,nn,is,in,ip,i,j,k,h,fn
61 integer*4 :: length_cne
62
63 integer*4, dimension(0:cs_nnz) :: cs
64 integer*4, dimension(nm) :: tm,sd
65 integer*4, dimension(nl_expl) :: num_ne
66 integer*4, dimension(nl_expl) :: fun_expl
67 integer*4, dimension(nn_loc) :: local_n_num
68 integer*4, dimension(nf) :: tag_func
69
70 integer*4, dimension(length_cne,4) :: check_ne
71 integer*4, dimension(max_num_ne,nl_expl) :: sour_ne
72
73 real*8 :: vel_prop
74
75 real*8, dimension(:), allocatable :: ct,ww
76
77 real*8, dimension(:,:), allocatable :: dd
78 real*8, dimension(nl_expl,20) :: val_expl
79 real*8, dimension(max_num_ne,nl_expl) :: dist_sour_ne
80 real*8, dimension(length_cne,1) :: check_dist_ne
81
82
83
84 h = 0
85 nn = 2
86 ne = cs(0) -1
87
88 allocate(ct(nn),ww(nn),dd(nn,nn))
89 call make_lgl_nw(nn,ct,ww,dd)
90
91 do im = 1,nm
92 if ((sd(im) +1).ne.nn) then
93 deallocate(ct,ww,dd)
94 nn = sd(im) +1
95 allocate(ct(nn),ww(nn),dd(nn,nn))
96 call make_lgl_nw(nn,ct,ww,dd)
97 endif
98
99 do ie = 1,ne
100 if (cs(cs(ie -1) +0).eq.tm(im)) then
101 do k = 1,nn
102 do j = 1,nn
103 do i = 1,nn
104 is = nn*nn*(k -1) +nn*(j -1) +i
105 in = cs(cs(ie -1) +is)
106
107 do iexpl = 1,nl_expl
108 do ip = 1,num_ne(iexpl)
109
110 if (local_n_num(in) .eq. sour_ne(ip,iexpl)) then
111
112 h = h + 1
113 check_ne(h,1) = sour_ne(ip,iexpl)
114 check_ne(h,2) = fun_expl(iexpl)
115 check_ne(h,3) = iexpl
116 check_ne(h,4) = ie
117 check_dist_ne(h,1) = dist_sour_ne(ip,iexpl) / val_expl(iexpl,19)
118
119 endif
120 enddo !ip
121 enddo !iexpl
122 enddo !i
123 enddo !j
124 enddo !k
125
126
127 endif !if (cs....)
128
129 enddo !ie = 1,ne
130
131 enddo !im = 1,nm
132
133 return
134
subroutine make_lgl_nw(nb_pnt, xq, wq, dd)
Makes Gauss-Legendre-Lobatto nodes, weigths and spectral derivatives.

References make_lgl_nw().

Here is the call graph for this function: