105 class(
schwarz_t),
target,
intent(inout) :: this
106 type(
space_t),
target,
intent(inout) :: Xh
107 type(
dofmap_t),
target,
intent(in) :: dof
108 type(
gs_t),
target,
intent(inout) :: gs_h
109 type(
mesh_t),
target,
intent(inout) :: msh
110 type(
bc_list_t),
target,
intent(inout):: bclst
114 call this%Xh_schwarz%init(
gll, xh%lx+2, xh%lx+2, xh%lx+2)
115 call this%dm_schwarz%init(msh, this%Xh_schwarz)
116 call this%gs_schwarz%init(this%dm_schwarz)
118 allocate(this%work1(this%dm_schwarz%size()))
119 allocate(this%work2(this%dm_schwarz%size()))
120 allocate(this%wt(xh%lx, xh%lx, 4, msh%gdim, msh%nelv))
122 call this%fdm%init(xh, dof, gs_h)
131 call device_map(this%work1, this%work1_d,this%dm_schwarz%size())
132 call device_map(this%work2, this%work2_d,this%dm_schwarz%size())
139 call rone(this%work1, this%dof%size())
169 integer :: enx,eny,enz, n, ie, k, ns
170 real(kind=
rp),
parameter :: zero = 0.0
171 real(kind=
rp),
parameter :: one = 1.0
172 associate(work1 => this%work1, work2 => this%work2, msh => this%msh, &
173 xh => this%Xh, xh_schwarz => this%Xh_schwarz)
180 if(.not. msh%gdim .eq. 3) enz=1
181 ns = enx*eny*enz*msh%nelv
184 call rzero(work1, ns)
188 call schwarz_extrude(work1, 0, zero, work2, 0, one , enx, eny, enz, msh%nelv)
192 call this%gs_schwarz%op(work2, ns, gs_op_add)
196 call this%gs_schwarz%op(work2, ns, gs_op_add)
198 call schwarz_extrude(work2, 0, one, work1, 0, -one, enx, eny, enz, msh%nelv)
199 call schwarz_extrude(work2, 2, one, work2, 0, one, enx, eny, enz, msh%nelv)
210 call this%gs_h%op(work1, n, gs_op_add)
214 call this%gs_h%op(work1, n, gs_op_add)
219 if (msh%gdim .eq. 2)
then
222 if (this%msh%gdim.eq. 3)
then
232 integer,
intent(in) :: n, nelv
233 real(kind=rp),
intent(inout) :: wt(n,4,2,nelv)
234 real(kind=rp),
intent(inout) :: work(n,n)
237 wt(j,1,1,ie) = 1.0_rp / work(1,j)
238 wt(j,2,1,ie) = 1.0_rp / work(2,j)
239 wt(j,3,1,ie) = 1.0_rp / work(n-1,j)
240 wt(j,4,1,ie) = 1.0_rp / work(n,j)
243 wt(i,1,2,ie) = 1.0_rp / work(i,1)
244 wt(i,2,2,ie) = 1.0_rp / work(i,2)
245 wt(i,3,2,ie) = 1.0_rp / work(i,n-1)
246 wt(i,4,2,ie) = 1.0_rp / work(i,n)
254 integer,
intent(in) ::n, nelv, ie
255 real(kind=rp),
intent(inout) :: wt(n,n,4,3,nelv)
256 real(kind=rp),
intent(inout) :: work(n,n,n)
261 wt(j,k,1,1,ie) = 1.0_rp / work(1,j,k)
262 wt(j,k,2,1,ie) = 1.0_rp / work(2,j,k)
263 wt(j,k,3,1,ie) = 1.0_rp / work(n-1,j,k)
264 wt(j,k,4,1,ie) = 1.0_rp / work(n,j,k)
270 wt(i,k,1,2,ie) = 1.0_rp / work(i,1,k)
271 wt(i,k,2,2,ie) = 1.0_rp / work(i,2,k)
272 wt(i,k,3,2,ie) = 1.0_rp / work(i,n-1,k)
273 wt(i,k,4,2,ie) = 1.0_rp / work(i,n,k)
279 wt(i,j,1,3,ie) = 1.0_rp / work(i,j,1)
280 wt(i,j,2,3,ie) = 1.0_rp / work(i,j,2)
281 wt(i,j,3,3,ie) = 1.0_rp / work(i,j,n-1)
282 wt(i,j,4,3,ie) = 1.0_rp / work(i,j,n)
327 integer,
intent(in) :: l1, l2, nx, ny, nz, nelv
328 real(kind=rp),
intent(inout) :: arr1(nx,ny,nz,nelv), arr2(nx,ny,nz,nelv)
329 real(kind=rp),
intent(in) :: f1, f2
330 integer :: i, j, k, ie, i0, i1
337 arr1(l1+1 ,j,1,ie) = f1*arr1(l1+1 ,j,1,ie) &
338 +f2*arr2(l2+1 ,j,1,ie)
339 arr1(nx-l1,j,1,ie) = f1*arr1(nx-l1,j,1,ie) &
340 +f2*arr2(nx-l2,j,1,ie)
343 arr1(i,l1+1 ,1,ie) = f1*arr1(i,l1+1 ,1,ie) &
344 +f2*arr2(i,l2+1 ,1,ie)
345 arr1(i,ny-l1,1,ie) = f1*arr1(i,ny-l1,1,ie) &
346 +f2*arr2(i,nx-l2,1,ie)
353 arr1(l1+1 ,j,k,ie) = f1*arr1(l1+1 ,j,k,ie) &
354 +f2*arr2(l2+1 ,j,k,ie)
355 arr1(nx-l1,j,k,ie) = f1*arr1(nx-l1,j,k,ie) &
356 +f2*arr2(nx-l2,j,k,ie)
361 arr1(i,l1+1 ,k,ie) = f1*arr1(i,l1+1 ,k,ie) &
362 +f2*arr2(i,l2+1 ,k,ie)
363 arr1(i,nx-l1,k,ie) = f1*arr1(i,nx-l1,k,ie) &
364 +f2*arr2(i,nx-l2,k,ie)
369 arr1(i,j,l1+1 ,ie) = f1*arr1(i,j,l1+1 ,ie) &
370 +f2*arr2(i,j,l2+1 ,ie)
371 arr1(i,j,nx-l1,ie) = f1*arr1(i,j,nx-l1,ie) &
372 +f2*arr2(i,j,nx-l2,ie)
381 real(kind=rp),
dimension(this%dof%size()),
intent(inout) :: e, r
382 integer :: n, enx, eny, enz, ns
383 real(kind=rp),
parameter :: zero = 0.0_rp
384 real(kind=rp),
parameter :: one = 1.0_rp
385 type(c_ptr) :: e_d, r_d
386 associate(work1 => this%work1, work1_d => this%work1_d,&
387 work2 => this%work2, work2_d => this%work2_d)
390 enx=this%Xh_schwarz%lx
391 eny=this%Xh_schwarz%ly
392 enz=this%Xh_schwarz%lz
393 if(.not. this%msh%gdim .eq. 3) enz=1
394 ns = enx*eny*enz*this%msh%nelv
395 if (neko_bcknd_device .eq. 1)
then
396 r_d = device_get_ptr(r)
397 e_d = device_get_ptr(e)
398 call device_event_record(this%event, glb_cmd_queue)
399 call device_stream_wait_event(aux_cmd_queue, this%event, 0)
400 call device_schwarz_toext3d(work1_d, r_d, this%Xh%lx, &
401 this%msh%nelv, aux_cmd_queue)
402 call device_schwarz_extrude(work1_d, 0, zero, work1_d, 2, one, &
403 enx,eny,enz, this%msh%nelv,aux_cmd_queue)
405 this%gs_schwarz%bcknd%gs_stream = aux_cmd_queue
406 call this%gs_schwarz%op(work1, ns, gs_op_add,this%event)
407 call device_event_sync(this%event)
408 call device_schwarz_extrude(work1_d, 0, one, work1_d, 2, -one, &
409 enx, eny, enz, this%msh%nelv, aux_cmd_queue)
411 call this%fdm%compute(work2, work1,aux_cmd_queue)
413 call device_schwarz_extrude(work1_d, 0, zero, work2_d, 0, one, &
414 enx, eny, enz, this%msh%nelv, aux_cmd_queue)
415 call this%gs_schwarz%op(work2, ns, gs_op_add,this%event)
416 call device_event_sync(this%event)
418 call device_schwarz_extrude(work2_d, 0, one, work1_d, 0, -one, &
419 enx, eny, enz, this%msh%nelv, aux_cmd_queue)
420 call device_schwarz_extrude(work2_d, 2, one, work2_d, 0, one, &
421 enx, eny, enz, this%msh%nelv, aux_cmd_queue)
422 call device_schwarz_toreg3d(e_d, work2_d, this%Xh%lx, &
423 this%msh%nelv, aux_cmd_queue)
425 call device_event_record(this%event,aux_cmd_queue)
426 call device_event_sync(this%event)
428 call this%gs_h%op(e, n, gs_op_add, this%event)
429 call this%bclst%apply_scalar(e, n)
430 call device_col2(e_d,this%wt_d, n)
431 call device_stream_wait_event(aux_cmd_queue, this%event, 0)
433 call this%bclst%apply_scalar(r, n)
438 enx, eny, enz, this%msh%nelv)
439 call this%gs_schwarz%op(work1, ns, gs_op_add)
441 enx, eny, enz, this%msh%nelv)
443 call this%fdm%compute(work2, work1)
447 enx, eny, enz, this%msh%nelv)
448 call this%gs_schwarz%op(work2, ns, gs_op_add)
450 enx, eny, enz, this%msh%nelv)
452 enx, eny, enz, this%msh%nelv)
457 call this%gs_h%op(e, n, gs_op_add)
458 call this%bclst%apply_scalar(e, n)
460 call schwarz_wt3d(e, this%wt, this%Xh%lx, this%msh%nelv)
467 integer,
intent(in) :: n, nelv
468 real(kind=rp),
intent(inout) :: e(n,n,n,nelv)
469 real(kind=rp),
intent(inout) :: wt(n,n,4,3,nelv)
470 integer :: ie, i, j, k
475 e(1 ,j,k,ie) = e(1 ,j,k,ie) * wt(j,k,1,1,ie)
476 e(2 ,j,k,ie) = e(2 ,j,k,ie) * wt(j,k,2,1,ie)
477 e(n-1,j,k,ie) = e(n-1,j,k,ie) * wt(j,k,3,1,ie)
478 e(n ,j,k,ie) = e(n ,j,k,ie) * wt(j,k,4,1,ie)
483 e(i,1 ,k,ie) = e(i,1 ,k,ie) * wt(i,k,1,2,ie)
484 e(i,2 ,k,ie) = e(i,2 ,k,ie) * wt(i,k,2,2,ie)
485 e(i,n-1,k,ie) = e(i,n-1,k,ie) * wt(i,k,3,2,ie)
486 e(i,n ,k,ie) = e(i,n ,k,ie) * wt(i,k,4,2,ie)
491 e(i,j,1 ,ie) = e(i,j,1 ,ie) * wt(i,j,1,3,ie)
492 e(i,j,2 ,ie) = e(i,j,2 ,ie) * wt(i,j,2,3,ie)
493 e(i,j,n-1,ie) = e(i,j,n-1,ie) * wt(i,j,3,3,ie)
494 e(i,j,n ,ie) = e(i,j,n ,ie) * wt(i,j,4,3,ie)