Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
schwarz.f90
Go to the documentation of this file.
1! Copyright (c) 2008-2020, UCHICAGO ARGONNE, LLC.
2!
3! The UChicago Argonne, LLC as Operator of Argonne National
4! Laboratory holds copyright in the Software. The copyright holder
5! reserves all rights except those expressly granted to licensees,
6! and U.S. Government license rights.
7!
8! Redistribution and use in source and binary forms, with or without
9! modification, are permitted provided that the following conditions
10! are met:
11!
12! 1. Redistributions of source code must retain the above copyright
13! notice, this list of conditions and the disclaimer below.
14!
15! 2. Redistributions in binary form must reproduce the above copyright
16! notice, this list of conditions and the disclaimer (as noted below)
17! in the documentation and/or other materials provided with the
18! distribution.
19!
20! 3. Neither the name of ANL nor the names of its contributors
21! may be used to endorse or promote products derived from this software
22! without specific prior written permission.
23!
24! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
25! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
26! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
27! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL
28! UCHICAGO ARGONNE, LLC, THE U.S. DEPARTMENT OF
29! ENERGY OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
30! SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED
31! TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
32! DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
33! THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
34! (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
35! OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36!
37! Additional BSD Notice
38! ---------------------
39! 1. This notice is required to be provided under our contract with
40! the U.S. Department of Energy (DOE). This work was produced at
41! Argonne National Laboratory under Contract
42! No. DE-AC02-06CH11357 with the DOE.
43!
44! 2. Neither the United States Government nor UCHICAGO ARGONNE,
45! LLC nor any of their employees, makes any warranty,
46! express or implied, or assumes any liability or responsibility for the
47! accuracy, completeness, or usefulness of any information, apparatus,
48! product, or process disclosed, or represents that its use would not
49! infringe privately-owned rights.
50!
51! 3. Also, reference herein to any specific commercial products, process,
52! or services by trade name, trademark, manufacturer or otherwise does
53! not necessarily constitute or imply its endorsement, recommendation,
54! or favoring by the United States Government or UCHICAGO ARGONNE LLC.
55! The views and opinions of authors expressed
56! herein do not necessarily state or reflect those of the United States
57! Government or UCHICAGO ARGONNE, LLC, and shall
58! not be used for advertising or product endorsement purposes.
59!
61module schwarz
62 use num_types
63 use math
64 use mesh
65 use space
66 use dofmap, only : dofmap_t
67 use bc
70 use device_math
71 use fdm, only : fdm_t
72 use device
73 use neko_config
74 use bc_list, only : bc_list_t
75 use, intrinsic :: iso_c_binding
76 implicit none
77 private
78
79 type, public :: schwarz_t
80 real(kind=rp), allocatable :: work1(:)
81 real(kind=rp), allocatable :: work2(:)
82 real(kind=rp), allocatable :: wt(:,:,:,:,:)
83 type(c_ptr) :: work1_d = c_null_ptr
84 type(c_ptr) :: work2_d = c_null_ptr
85 type(c_ptr) :: wt_d = c_null_ptr
86 type(space_t) :: xh_schwarz
87 type(gs_t) :: gs_schwarz
88 type(dofmap_t) :: dm_schwarz
89 type(fdm_t) :: fdm
90 type(space_t), pointer :: xh
91 type(bc_list_t), pointer :: bclst
92 type(dofmap_t), pointer :: dof
93 type(gs_t), pointer :: gs_h
94 type(mesh_t), pointer :: msh
95 type(c_ptr) :: event
96 contains
97 procedure, pass(this) :: init => schwarz_init
98 procedure, pass(this) :: free => schwarz_free
99 procedure, pass(this) :: compute => schwarz_compute
100 end type schwarz_t
101
102contains
103
104 subroutine schwarz_init(this, Xh, dof, gs_h, bclst, msh)
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
111
112 call this%free()
113
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)
117
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))
121
122 call this%fdm%init(xh, dof, gs_h)
123
124
125 this%msh => msh
126 this%Xh => xh
127 this%bclst => bclst
128 this%dof => dof
129 this%gs_h => gs_h
130 if (neko_bcknd_device .eq. 1) then
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())
133 end if
134
135
136 call schwarz_setup_wt(this)
137 if (neko_bcknd_device .eq. 1) then
138 call device_alloc(this%wt_d,int(this%dof%size()*rp, i8))
139 call rone(this%work1, this%dof%size())
140 call schwarz_wt3d(this%work1, this%wt, xh%lx, msh%nelv)
141 call device_memcpy(this%work1, this%wt_d, this%dof%size(), &
142 host_to_device, sync=.false.)
143 call device_event_create(this%event, 2)
144 end if
145 end subroutine schwarz_init
146
147 subroutine schwarz_free(this)
148 class(schwarz_t), intent(inout) :: this
149
150 if(allocated(this%work1)) deallocate(this%work1)
151 if(allocated(this%work2)) deallocate(this%work2)
152 if(allocated(this%wt)) deallocate(this%wt)
153
154 call this%Xh_schwarz%free()
155 call this%gs_schwarz%free()
156 !why cant I do this?
157 !call dofmap_free(this%dm_schwarz)
158 call this%fdm%free()
159
160 nullify(this%Xh)
161 nullify(this%bclst)
162 nullify(this%dof)
163 nullify(this%gs_h)
164 nullify(this%msh)
165 end subroutine schwarz_free
167 subroutine schwarz_setup_wt(this)
168 class(schwarz_t), intent(inout) :: this
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)
174
175 n = this%dof%size()
176
177 enx = xh_schwarz%lx
178 eny = xh_schwarz%ly
179 enz = xh_schwarz%lz
180 if(.not. msh%gdim .eq. 3) enz=1
181 ns = enx*eny*enz*msh%nelv
182
183 call rone(work2, ns)
184 call rzero(work1, ns)
185
186 ! Sum overlap region (border excluded)
187 ! Cred to PFF for this, very clever
188 call schwarz_extrude(work1, 0, zero, work2, 0, one , enx, eny, enz, msh%nelv)
189 if (neko_bcknd_device .eq. 1) then
190 call device_memcpy(work2, this%work2_d, ns, &
191 host_to_device, sync=.false.)
192 call this%gs_schwarz%op(work2, ns, gs_op_add)
193 call device_memcpy(work2, this%work2_d, ns, &
194 device_to_host, sync=.false.)
195 else
196 call this%gs_schwarz%op(work2, ns, gs_op_add)
197 end if
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)
200
201 ! if(.not.if3d) then ! Go back to regular size array
202 ! call hsmg_schwarz_toreg2d(mg_work,mg_work(i),mg_nh(l))
203 ! else
204 call schwarz_toreg3d(work1, work2, xh%lx, msh%nelv)
205 ! endif
206
207 if (neko_bcknd_device .eq. 1) then
208 call device_memcpy(work1, this%work1_d, n, &
209 host_to_device, sync=.false.)
210 call this%gs_h%op(work1, n, gs_op_add)
211 call device_memcpy(work1, this%work1_d, n, &
212 device_to_host, sync=.true.)
213 else
214 call this%gs_h%op(work1, n, gs_op_add)
215 end if
216
217 k = 1
218 do ie = 1,msh%nelv
219 if (msh%gdim .eq. 2) then
220 call schwarz_setup_schwarz_wt2d_2(this%wt,ie,xh%lx, work1(k), msh%nelv)
221 end if
222 if (this%msh%gdim.eq. 3) then
223 call schwarz_setup_schwarz_wt3d_2(this%wt,ie,xh%lx, work1(k), msh%nelv)
224 k = k + xh%lxyz
225 end if
226 end do
227 end associate
228 end subroutine schwarz_setup_wt
229
231 subroutine schwarz_setup_schwarz_wt2d_2(wt,ie,n,work, nelv)
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)
235 integer :: ie,i,j
236 do j = 1, 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)
241 end do
242 do i = 1, n
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)
247 end do
248
249 return
250 end subroutine schwarz_setup_schwarz_wt2d_2
251
253 subroutine schwarz_setup_schwarz_wt3d_2(wt, ie, n, work, nelv)
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)
257 integer :: i,j,k
258
259 do k = 1, n
260 do j = 1, 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)
265 end do
266 end do
267
268 do k = 1, n
269 do i = 1, n
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)
274 end do
275 end do
276
277 do j = 1, n
278 do i = 1, n
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)
283 end do
284 end do
285
286 end subroutine schwarz_setup_schwarz_wt3d_2
287
289 subroutine schwarz_toreg3d(b, a, n, nelv)
290 integer, intent(in) :: n, nelv
291 real(kind=rp), intent(inout) :: a(0:n+1, 0:n+1, 0:n+1, nelv)
292 real(kind=rp), intent(inout) :: b(n,n,n,nelv)
293 integer :: i, j, k, ie
294 do ie = 1, nelv
295 do k = 1, n
296 do j = 1, n
297 do i = 1, n
298 b(i,j,k,ie) = a(i,j,k,ie)
299 end do
300 end do
301 end do
302 end do
303 end subroutine schwarz_toreg3d
304
306 subroutine schwarz_toext3d(a, b, n, nelv)
307 integer, intent(in) :: n, nelv
308 real (kind=rp), intent(inout) :: a(0:n+1,0:n+1,0:n+1,nelv), b(n,n,n,nelv)
309 integer :: i,j,k,ie
310
311 call rzero(a, (n+2)*(n+2)*(n+2)*nelv)
312 do ie = 1, nelv
313 do k = 1, n
314 do j = 1, n
315 do i = 1, n
316 a(i,j,k,ie) = b(i,j,k,ie)
317 end do
318 end do
319 end do
320 end do
321 end subroutine schwarz_toext3d
322
326 subroutine schwarz_extrude(arr1, l1, f1, arr2, l2, f2, nx, ny, nz, nelv)
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
331 i0=2
332 i1=nx-1
333
334 if(nz .eq. 1) then
335 do ie = 1, nelv
336 do j = 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)
341 end do
342 do i = i0, i1
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)
347 end do
348 end do
349 else
350 do ie = 1, nelv
351 do k = i0, i1
352 do j = i0, i1
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)
357 end do
358 end do
359 do k = i0, i1
360 do i = i0, i1
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)
365 end do
366 end do
367 do j = i0, i1
368 do i = i0, i1
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)
373 end do
374 end do
375 end do
376 endif
377 end subroutine schwarz_extrude
378
379 subroutine schwarz_compute(this, e, r)
380 class(schwarz_t), intent(inout) :: this
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)
388
389 n = this%dof%size()
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)
404
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)
410
411 call this%fdm%compute(work2, work1,aux_cmd_queue) ! do local solves
412
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)
417
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)
424
425 call device_event_record(this%event,aux_cmd_queue)
426 call device_event_sync(this%event)
427
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)
432 else
433 call this%bclst%apply_scalar(r, n)
434 call schwarz_toext3d(work1, r, this%Xh%lx, this%msh%nelv)
435
436 ! exchange interior nodes
437 call schwarz_extrude(work1, 0, zero, work1, 2, one, &
438 enx, eny, enz, this%msh%nelv)
439 call this%gs_schwarz%op(work1, ns, gs_op_add)
440 call schwarz_extrude(work1, 0, one, work1, 2, -one, &
441 enx, eny, enz, this%msh%nelv)
442
443 call this%fdm%compute(work2, work1) ! do local solves
444
445 ! Sum overlap region (border excluded)
446 call schwarz_extrude(work1, 0, zero, work2, 0, one, &
447 enx, eny, enz, this%msh%nelv)
448 call this%gs_schwarz%op(work2, ns, gs_op_add)
449 call schwarz_extrude(work2, 0, one, work1, 0, -one, &
450 enx, eny, enz, this%msh%nelv)
451 call schwarz_extrude(work2, 2, one, work2, 0, one, &
452 enx, eny, enz, this%msh%nelv)
453
454 call schwarz_toreg3d(e, work2, this%Xh%lx, this%msh%nelv)
455
456 ! sum border nodes
457 call this%gs_h%op(e, n, gs_op_add)
458 call this%bclst%apply_scalar(e, n)
459
460 call schwarz_wt3d(e, this%wt, this%Xh%lx, this%msh%nelv)
461 end if
462 end associate
463 end subroutine schwarz_compute
464
465 !Apply schwarz weights along the boundary of each element.
466 subroutine schwarz_wt3d(e,wt,n, 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
471
472 do ie = 1, nelv
473 do k = 1, n
474 do j = 1, n
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)
479 end do
480 end do
481 do k = 1, n
482 do i = 3, n-2
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)
487 end do
488 end do
489 do j = 3, n-2
490 do i = 3, n-2
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)
495 end do
496 end do
497 end do
498 end subroutine schwarz_wt3d
499end module schwarz
Map a Fortran array to a device (allocate and associate)
Definition device.F90:57
Copy data between host and device (or device and device)
Definition device.F90:51
Defines a list of bc_t.
Definition bc_list.f90:34
Defines a boundary condition.
Definition bc.f90:34
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public host_to_device
Definition device.F90:47
integer, parameter, public device_to_host
Definition device.F90:47
subroutine, public device_alloc(x_d, s)
Allocate memory on the device.
Definition device.F90:164
subroutine, public device_event_create(event, flags)
Create a device event queue.
Definition device.F90:1164
Defines a mapping of the degrees of freedom.
Definition dofmap.f90:35
Type for the Fast Diagonalization connected with the schwarz overlapping solves.
Definition fdm.f90:61
Gather-scatter.
Definition math.f90:60
subroutine, public rone(a, n)
Set all elements to one.
Definition math.f90:227
subroutine, public rzero(a, n)
Zero a real vector.
Definition math.f90:194
Defines a mesh.
Definition mesh.f90:34
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public i8
Definition num_types.f90:7
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Overlapping schwarz solves.
Definition schwarz.f90:61
subroutine schwarz_toreg3d(b, a, n, nelv)
convert array a from extended size to regular
Definition schwarz.f90:290
subroutine schwarz_setup_wt(this)
setup weights
Definition schwarz.f90:168
subroutine schwarz_init(this, xh, dof, gs_h, bclst, msh)
Definition schwarz.f90:105
subroutine schwarz_extrude(arr1, l1, f1, arr2, l2, f2, nx, ny, nz, nelv)
Sum values along rows l1, l2 with weights f1, f2 and store along row l1. Helps us avoid complicated c...
Definition schwarz.f90:327
subroutine schwarz_setup_schwarz_wt3d_2(wt, ie, n, work, nelv)
Setup schwarz weights, 3d, second step.
Definition schwarz.f90:254
subroutine schwarz_compute(this, e, r)
Definition schwarz.f90:380
subroutine schwarz_wt3d(e, wt, n, nelv)
Definition schwarz.f90:467
subroutine schwarz_toext3d(a, b, n, nelv)
convert array a from original size to size extended array with border
Definition schwarz.f90:307
subroutine schwarz_setup_schwarz_wt2d_2(wt, ie, n, work, nelv)
Setup schwarz weights, 2d, second step.
Definition schwarz.f90:232
subroutine schwarz_free(this)
Definition schwarz.f90:148
Defines a function space.
Definition space.f90:34
integer, parameter, public gll
Definition space.f90:48
A list of allocatable `bc_t`. Follows the standard interface of lists.
Definition bc_list.f90:46
The function space for the SEM solution fields.
Definition space.f90:62