3 subroutine setupds(gs_handle,nx,ny,nz,nel,melg,vertex,glo_num)
4 use kinds, only : dp, i8
12 integer(i8) :: glo_num(1),ngv
15 integer :: nx, nel, ntot, ny, nz, melg
16 integer,
parameter :: num_ds_handles = 64
17 integer,
save :: handles(num_ds_handles) = -1
20 if (nx <= num_ds_handles .and. handles(nx) >= 0)
then
21 gs_handle = handles(nx)
30 call
set_vert(glo_num,ngv,nx,nel,vertex, .false. )
34 call
gs_setup(gs_handle,glo_num,ntot,nekcomm,mp)
40 write(6,1) t1,gs_handle,nx,ngv,melg
41 1
format(
' setupds time',1pe11.4,
' seconds ',2i3,2i12)
45 if (nx <= num_ds_handles) handles(nx) = gs_handle
57 use ctimer, only : tdsum, ndsum
58 use input, only : ifldmhd
60 use tstep, only : ifield
63 real(DP),
intent(inout) :: u(*)
70 if (ifldt == ifldmhd) ifldt = 1
98 call
gs_op(gsh_fld(ifldt),u,1,1,0)
110 tdsmx=max(timee,tdsmx)
111 tdsmn=min(timee,tdsmn)
137 use input, only : ifldmhd
139 use tstep, only : ifield
148 if (ifldt == ifldmhd) ifldt = 1
152 if (op ==
'+ ') call
gs_op(gsh_fld(ifldt),u,1,1,0)
153 if (op ==
'sum') call
gs_op(gsh_fld(ifldt),u,1,1,0)
154 if (op ==
'SUM') call
gs_op(gsh_fld(ifldt),u,1,1,0)
156 if (op ==
'* ') call
gs_op(gsh_fld(ifldt),u,1,2,0)
157 if (op ==
'mul') call
gs_op(gsh_fld(ifldt),u,1,2,0)
158 if (op ==
'MUL') call
gs_op(gsh_fld(ifldt),u,1,2,0)
160 if (op ==
'm ') call
gs_op(gsh_fld(ifldt),u,1,3,0)
161 if (op ==
'min') call
gs_op(gsh_fld(ifldt),u,1,3,0)
162 if (op ==
'mna') call
gs_op(gsh_fld(ifldt),u,1,3,0)
163 if (op ==
'MIN') call
gs_op(gsh_fld(ifldt),u,1,3,0)
164 if (op ==
'MNA') call
gs_op(gsh_fld(ifldt),u,1,3,0)
166 if (op ==
'M ') call
gs_op(gsh_fld(ifldt),u,1,4,0)
167 if (op ==
'max') call
gs_op(gsh_fld(ifldt),u,1,4,0)
168 if (op ==
'mxa') call
gs_op(gsh_fld(ifldt),u,1,4,0)
169 if (op ==
'MAX') call
gs_op(gsh_fld(ifldt),u,1,4,0)
170 if (op ==
'MXA') call
gs_op(gsh_fld(ifldt),u,1,4,0)
182 use size_m
, only : ndim
183 use ctimer, only : ifsync, icalld, tvdss, tgsum, nvdss, etime1,
dnekclock
184 use ctimer, only : tdsmx, tdsmn
185 use input, only : ifldmhd
187 use tstep, only : ifield
190 REAL(DP) :: U(1),V(1),W(1)
197 if (icalld == 0) tvdss=0.0d0
198 if (icalld == 0) tgsum=0.0d0
210 if (ifldt == ifldmhd) ifldt = 1
212 call gs_op_many(gsh_fld(ifldt),u,v,w,u,u,u,ndim,1,1,0)
217 tdsmx=max(timee,tdsmx)
218 tdsmn=min(timee,tdsmn)
230 use size_m
, only : ndim
232 use input, only : ifldmhd
234 use tstep, only : ifield
237 real(DP) :: u(1),v(1),w(1)
247 if (ifldt == ifldmhd) ifldt = 1
252 if (op ==
'+ ' .OR. op ==
'sum' .OR. op ==
'SUM') &
253 call gs_op_many(gsh_fld(ifldt),u,v,w,u,u,u,ndim,1,1,0)
256 if (op ==
'* ' .OR. op ==
'mul' .OR. op ==
'MUL') &
257 call gs_op_many(gsh_fld(ifldt),u,v,w,u,u,u,ndim,1,2,0)
260 if (op ==
'm ' .OR. op ==
'min' .OR. op ==
'mna' &
261 .OR. op ==
'MIN' .OR. op ==
'MNA') &
262 call gs_op_many(gsh_fld(ifldt),u,v,w,u,u,u,ndim,1,3,0)
265 if (op ==
'M ' .OR. op ==
'max' .OR. op ==
'mxa' &
266 .OR. op ==
'MAX' .OR. op ==
'MXA') &
267 call gs_op_many(gsh_fld(ifldt),u,v,w,u,u,u,ndim,1,4,0)
subroutine dssum(u)
Direct stiffness sum.
subroutine vec_dsop(u, v, w, op)
Direct stiffness summation of the face data, for field U. Boundary condition data corresponds to comp...
real(dp) function dnekclock()
subroutine setupds(gs_handle, nx, ny, nz, nel, melg, vertex, glo_num)
setup data structures for direct stiffness operations
subroutine dsop(u, op)
generalization of dssum to other reducers.
subroutine set_vert(glo_num, ngv, nx, nel, vertex, ifcenter)
Given global array, vertex, pointing to hex vertices, set up a new array of global pointers for an nx...
subroutine vec_dssum(u, v, w)
Direct stiffness summation of the face data, for field U.