Nek5000
SEM for Incompressible NS
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
dssum.F90
Go to the documentation of this file.
1 !-----------------------------------------------------------------------
3 subroutine setupds(gs_handle,nx,ny,nz,nel,melg,vertex,glo_num)
4  use kinds, only : dp, i8
5  use size_m, only : nid
6  use ctimer, only : dnekclock
7  use parallel, only : mp=>np, nekcomm
8  implicit none
9 
10  integer :: gs_handle
11  integer :: vertex(1)
12  integer(i8) :: glo_num(1),ngv
13 
14  real(DP) :: t0, t1
15  integer :: nx, nel, ntot, ny, nz, melg
16  integer, parameter :: num_ds_handles = 64
17  integer, save :: handles(num_ds_handles) = -1
18 
19  ! If we've already set up it, return an existing handle
20  if (nx <= num_ds_handles .and. handles(nx) >= 0) then
21  gs_handle = handles(nx)
22  return
23 
24  ! Make a new handle
25  else
26 
27  t0 = dnekclock()
28 
29  ! Global-to-local mapping for gs
30  call set_vert(glo_num,ngv,nx,nel,vertex, .false. )
31 
32  ! Initialize gather-scatter code
33  ntot = nx*ny*nz*nel
34  call gs_setup(gs_handle,glo_num,ntot,nekcomm,mp)
35 
36  ! call gs_chkr(glo_num)
37 
38  t1 = dnekclock() - t0
39  if (nid == 0) then
40  write(6,1) t1,gs_handle,nx,ngv,melg
41  1 format(' setupds time',1pe11.4,' seconds ',2i3,2i12)
42  endif
43 
44  ! save the handle
45  if (nx <= num_ds_handles) handles(nx) = gs_handle
46 
47  endif
48 
49  return
50 end subroutine setupds
51 
52 !-----------------------------------------------------------------------
54 subroutine dssum(u)
55  use kinds, only : dp
56  use ctimer, only : ifsync, icalld, tdsmx, tdsmn, etime1, dnekclock
57  use ctimer, only : tdsum, ndsum
58  use input, only : ifldmhd
59  use parallel, only : gsh_fld
60  use tstep, only : ifield
61  implicit none
62 
63  real(DP), intent(inout) :: u(*)
64 
65  integer :: ifldt
66  real(DP) :: timee
67 
68  ifldt = ifield
69 ! if (ifldt.eq.0) ifldt = 1
70  if (ifldt == ifldmhd) ifldt = 1
71 ! write(6,*) ifldt,ifield,gsh_fld(ifldt),imesh,' ifldt'
72 
73  if (ifsync) call nekgsync()
74 
75 #ifndef NOTIMER
76  if (icalld == 0) then
77  tdsmx=0.
78  tdsmn=0.
79  endif
80  icalld=icalld+1
81  etime1=dnekclock()
82 #endif
83 
84 
85 ! T ~ ~T T
86 ! Implement QQ := J Q Q J
87 
88 
89 ! T
90 ! This is the J part, translating child data
91 ! call apply_Jt(u,nx,ny,nz,nel)
92 
93 
94 
95 ! ~ ~T
96 ! This is the Q Q part
97 
98  call gs_op(gsh_fld(ifldt),u,1,1,0) ! 1 ==> +
99 
100 
101 
102 ! This is the J part, interpolating parent solution onto child
103 ! call apply_J(u,nx,ny,nz,nel)
104 
105 
106 #ifndef NOTIMER
107  timee=(dnekclock()-etime1)
108  tdsum=tdsum+timee
109  ndsum=icalld
110  tdsmx=max(timee,tdsmx)
111  tdsmn=min(timee,tdsmn)
112 #endif
113 
114  return
115 end subroutine dssum
116 
117 !-----------------------------------------------------------------------
134 subroutine dsop(u,op)
135  use kinds, only : dp
136  use ctimer, only : ifsync
137  use input, only : ifldmhd
138  use parallel, only : gsh_fld
139  use tstep, only : ifield
140  implicit none
141 
142  real(DP) :: u(1)
143  character(3) :: op
144  integer :: ifldt
145 
146  ifldt = ifield
147 ! if (ifldt.eq.0) ifldt = 1
148  if (ifldt == ifldmhd) ifldt = 1
149 
150  if (ifsync) call nekgsync()
151 
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)
155 
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)
159 
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)
165 
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)
171 
172  return
173 end subroutine dsop
174 
175 !-----------------------------------------------------------------------
180 subroutine vec_dssum(u,v,w)
181  use kinds, only : dp
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
186  use parallel, only : gsh_fld
187  use tstep, only : ifield
188  implicit none
189 
190  REAL(DP) :: U(1),V(1),W(1)
191  integer :: ifldt
192  real(DP) :: timee
193 
194  if(ifsync) call nekgsync()
195 
196 #ifndef NOTIMER
197  if (icalld == 0) tvdss=0.0d0
198  if (icalld == 0) tgsum=0.0d0
199  icalld=icalld+1
200  nvdss=icalld
201  etime1=dnekclock()
202 #endif
203 
204 !============================================================================
205 ! execution phase
206 !============================================================================
207 
208  ifldt = ifield
209 ! if (ifldt.eq.0) ifldt = 1
210  if (ifldt == ifldmhd) ifldt = 1
211 
212  call gs_op_many(gsh_fld(ifldt),u,v,w,u,u,u,ndim,1,1,0)
213 
214 #ifndef NOTIMER
215  timee=(dnekclock()-etime1)
216  tvdss=tvdss+timee
217  tdsmx=max(timee,tdsmx)
218  tdsmn=min(timee,tdsmn)
219 #endif
220 
221  return
222 end subroutine vec_dssum
223 
224 !-----------------------------------------------------------------------
228 subroutine vec_dsop(u,v,w,op)
229  use kinds, only : dp
230  use size_m, only : ndim
231  use ctimer, only : ifsync
232  use input, only : ifldmhd
233  use parallel, only : gsh_fld
234  use tstep, only : ifield
235  implicit none
236 
237  real(DP) :: u(1),v(1),w(1)
238  character(3) :: op
239  integer :: ifldt
240 
241 !============================================================================
242 ! execution phase
243 !============================================================================
244 
245  ifldt = ifield
246 ! if (ifldt.eq.0) ifldt = 1
247  if (ifldt == ifldmhd) ifldt = 1
248 
249 ! write(6,*) 'opdsop: ',op,ifldt,ifield
250  if(ifsync) call nekgsync()
251 
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)
254 
255 
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)
258 
259 
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)
263 
264 
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)
268 
269 
270  return
271 end subroutine vec_dsop
272 
273 !-----------------------------------------------------------------------
subroutine dssum(u)
Direct stiffness sum.
Definition: dssum.F90:54
cleaned
Definition: tstep_mod.F90:2
Input parameters from preprocessors.
Definition: input_mod.F90:11
subroutine vec_dsop(u, v, w, op)
Direct stiffness summation of the face data, for field U. Boundary condition data corresponds to comp...
Definition: dssum.F90:228
#define gs_setup
Definition: gs.c:29
real(dp) function dnekclock()
Definition: ctimer_mod.F90:103
subroutine setupds(gs_handle, nx, ny, nz, nel, melg, vertex, glo_num)
setup data structures for direct stiffness operations
Definition: dssum.F90:3
subroutine dsop(u, op)
generalization of dssum to other reducers.
Definition: dssum.F90:134
cleaned
Definition: parallel_mod.F90:2
#define gs_op
Definition: gs.c:15
subroutine nekgsync()
Definition: comm_mpi.F90:318
subroutine vec_dssum(u, v, w)
Direct stiffness summation of the face data, for field U.
Definition: dssum.F90:180
static uint np
Definition: findpts_test.c:63