Octopus
nl_operator.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
24 use accel_oct_m
25 use batch_oct_m
27 use debug_oct_m
28 use global_oct_m
29 use index_oct_m
30 use iso_c_binding
31 use math_oct_m
32 use mesh_oct_m
34 use mpi_oct_m
39 use parser_oct_m
41 use space_oct_m
43 use types_oct_m
45
46 implicit none
47
48 private
49 public :: &
71
77 private
78 integer :: nri = 0
79 integer, allocatable :: imin(:)
80 integer, allocatable :: imax(:)
81 integer, allocatable :: ri(:, :)
82 end type nl_operator_index_t
83
85 type nl_operator_t
86 private
87 type(stencil_t), public :: stencil
88 type(mesh_t), pointer :: mesh => null()
89 integer, allocatable :: nn(:)
90 integer, public :: np = 0
91 ! When running in parallel mode, the next three arrays are unique on each node.
92 real(real64), allocatable, public :: w(:,:)
93
94 logical, public :: const_w = .true.
95
96 type(accel_mem_t), public :: buff_weights
97 type(accel_mem_t), public :: buff_half_weights
98
99 character(len=40) :: label
100
102 integer, public :: nri = 0
103 integer, allocatable, public :: ri(:,:)
104 integer, allocatable, public :: rimap(:)
105 integer, allocatable, public :: rimap_inv(:)
106
107 integer :: ninner = 0
108 integer :: nouter = 0
109
110 type(nl_operator_index_t) :: inner
111 type(nl_operator_index_t) :: outer
112
113 type(accel_kernel_t) :: kernel
114 type(accel_mem_t) :: buff_imin
115 type(accel_mem_t) :: buff_imax
116 type(accel_mem_t) :: buff_ri
117 type(accel_mem_t) :: buff_map
118 type(accel_mem_t) :: buff_all
119 type(accel_mem_t) :: buff_inner
120 type(accel_mem_t) :: buff_outer
121 type(accel_mem_t) :: buff_stencil
122 type(accel_mem_t) :: buff_ip_to_xyz
123 type(accel_mem_t) :: buff_xyz_to_ip
124
125 ! For multigrid solvers
126 type(nl_operator_t), public, pointer :: coarser => null()
127
128 end type nl_operator_t
129
130 integer, parameter :: &
131 OP_FORTRAN = 0, &
132 op_vec = 1, &
133 op_min = op_fortran, &
134 op_max = op_vec
135
136 integer, parameter :: &
137 OP_INVMAP = 1, &
138 op_map = 2, &
139 op_nomap = 3
140
141 integer, public, parameter :: OP_ALL = 3, op_inner = 1, op_outer = 2
142
143 logical :: compact_boundaries
144
145 interface
146 integer function op_is_available(opid, type)
147 implicit none
148 integer, intent(in) :: opid, type
149 end function op_is_available
150 end interface
151
152 integer :: dfunction_global = -1
153 integer :: zfunction_global = -1
154 integer :: function_opencl
155
156contains
157
158 ! ---------------------------------------------------------
161 subroutine nl_operator_global_init(namespace)
162 type(namespace_t), intent(in) :: namespace
163
164 integer :: default
165
167
168 !%Variable OperateDouble
169 !%Type integer
170 !%Section Execution::Optimization
171 !%Default optimized
172 !%Description
173 !% This variable selects the subroutine used to apply non-local
174 !% operators over the grid for real functions.
175 !%Option fortran 0
176 !% The standard Fortran function.
177 !%Option optimized 1
178 !% This version is optimized using vector primitives (if available).
179 !%End
181 !%Variable OperateComplex
182 !%Type integer
183 !%Section Execution::Optimization
184 !%Default optimized
185 !%Description
186 !% This variable selects the subroutine used to apply non-local
187 !% operators over the grid for complex functions.
188 !%Option fortran 0
189 !% The standard Fortran function.
190 !%Option optimized 1
191 !% This version is optimized using vector primitives (if available).
192 !%End
193
194 default = op_vec
196 call parse_variable(namespace, 'OperateDouble', default, dfunction_global)
197 if (.not. varinfo_valid_option('OperateDouble', dfunction_global)) call messages_input_error(namespace, 'OperateDouble')
199 call parse_variable(namespace, 'OperateComplex', default, zfunction_global)
200 if (.not. varinfo_valid_option('OperateComplex', zfunction_global)) call messages_input_error(namespace, 'OperateComplex')
202
205 !%Variable OperateAccel
206 !%Type integer
207 !%Default map
208 !%Section Execution::Optimization
209 !%Description
210 !% This variable selects the subroutine used to apply non-local
211 !% operators over the grid when an accelerator device is used.
212 !%Option invmap 1
213 !% The standard implementation ported to OpenCL.
214 !%Option map 2
215 !% A different version, more suitable for GPUs.
216 !%End
217 call parse_variable(namespace, 'OperateAccel', op_map, function_opencl)
218
219 call messages_obsolete_variable(namespace, 'OperateOpenCL', 'OperateAccel')
220
221 end if
222
223 !%Variable NLOperatorCompactBoundaries
224 !%Type logical
225 !%Default no
226 !%Section Execution::Optimization
227 !%Description
228 !% (Experimental) When set to yes, for finite systems Octopus will
229 !% map boundary points for finite-differences operators to a few
230 !% memory locations. This increases performance, however it is
231 !% experimental and has not been thoroughly tested.
232 !%End
233
234 call parse_variable(namespace, 'NLOperatorCompactBoundaries', .false., compact_boundaries)
235
236 if (compact_boundaries) then
237 call messages_experimental('NLOperatorCompactBoundaries')
238 end if
241 end subroutine nl_operator_global_init
242
243 ! ---------------------------------------------------------
244
249 end subroutine nl_operator_global_end
250
251 ! ---------------------------------------------------------
255 subroutine nl_operator_init(op, label)
256 type(nl_operator_t), intent(inout) :: op
257 character(len=*), intent(in) :: label
258
259 push_sub(nl_operator_init)
260
261 op%label = label
262
263 pop_sub(nl_operator_init)
264 end subroutine nl_operator_init
265
266
267 ! ---------------------------------------------------------
268 subroutine nl_operator_copy(opo, opi)
269 type(nl_operator_t), intent(inout) :: opo
270 type(nl_operator_t), target, intent(in) :: opi
271
272 push_sub(nl_operator_copy)
273
274 ! We cannot currently copy the GPU kernel for the nl_operator
275 assert(.not. accel_is_enabled())
276
277 call nl_operator_end(opo)
278 call nl_operator_init(opo, opi%label)
279
280 call stencil_copy(opi%stencil, opo%stencil)
281
282 opo%np = opi%np
283 opo%mesh => opi%mesh
284
285 safe_allocate_source_a(opo%nn, opi%nn)
286 safe_allocate_source_a(opo%w, opi%w)
287
288 opo%const_w = opi%const_w
289
290 opo%nri = opi%nri
291 assert(allocated(opi%ri))
292
293 safe_allocate_source_a(opo%ri, opi%ri)
294 safe_allocate_source_a(opo%rimap, opi%rimap)
295 safe_allocate_source_a(opo%rimap_inv, opi%rimap_inv)
296
297 if (opi%mesh%parallel_in_domains) then
298 opo%inner%nri = opi%inner%nri
299 safe_allocate_source_a(opo%inner%imin, opi%inner%imin)
300 safe_allocate_source_a(opo%inner%imax, opi%inner%imax)
301 safe_allocate_source_a(opo%inner%ri, opi%inner%ri)
302
303 opo%outer%nri = opi%outer%nri
304 safe_allocate_source_a(opo%outer%imin, opi%outer%imin)
305 safe_allocate_source_a(opo%outer%imax, opi%outer%imax)
306 safe_allocate_source_a(opo%outer%ri, opi%outer%ri)
307 end if
308
309 ! We create the corresponding GPU buffers
310 if (accel_is_enabled() .and. opo%const_w) then
311 call accel_create_buffer(opo%buff_weights, accel_mem_read_only, type_float, opo%stencil%size)
312 call accel_write_buffer(opo%buff_weights, opo%stencil%size, opo%w(:, 1))
313 call accel_create_buffer(opo%buff_half_weights, accel_mem_read_only, type_float, opo%stencil%size)
314 call accel_write_buffer(opo%buff_half_weights, opo%stencil%size, -m_half*opo%w(:, 1))
315 end if
316
317
318 pop_sub(nl_operator_copy)
319 end subroutine nl_operator_copy
320
321
322 ! ---------------------------------------------------------
323 subroutine nl_operator_build(space, mesh, op, np, const_w, regenerate)
324 class(space_t), intent(in) :: space
325 type(mesh_t), target, intent(in) :: mesh
326 type(nl_operator_t), intent(inout) :: op
327 integer, intent(in) :: np
328 logical, optional, intent(in) :: const_w
329 logical, optional, intent(in) :: regenerate
330
331 integer :: ii, jj, p1(space%dim), time, current
332 integer, allocatable :: st1(:), st2(:), st1r(:)
333 integer :: nn
334 integer :: ir, maxp, iinner, iouter
335 logical :: change, force_change
336 character(len=200) :: flags
337 integer, allocatable :: inner_points(:), outer_points(:), all_points(:)
339 push_sub(nl_operator_build)
340
341 op%const_w = optional_default(const_w, .false.)
342 if (mesh%parallel_in_domains .and. .not. op%const_w) then
343 call messages_experimental('Domain parallelization with curvilinear coordinates')
344 end if
345
346 assert(np > 0)
347
348 ! store values in structure
349 op%np = np
350 op%mesh => mesh
351
352 ! allocate weights op%w
353 if (op%const_w) then
354 safe_allocate(op%w(1:op%stencil%size, 1))
355 message(1) = 'Debug: nl_operator_build: working with constant weights.'
356 call messages_info(1, debug_only=.true.)
357 else
358 safe_allocate(op%w(1:op%stencil%size, 1:op%np))
359 message(1) = 'Debug: nl_operator_build: working with non-constant weights.'
360 call messages_info(1, debug_only=.true.)
361 end if
362
363 ! set initially to zero
364 op%w = m_zero
365
366 ! Build lookup table
367 safe_allocate(st1(1:op%stencil%size))
368 safe_allocate(st1r(1:op%stencil%size))
369 safe_allocate(st2(1:op%stencil%size))
370
371 op%nri = 0
372 do time = 1, 2
373 st2 = 0
374 do ii = 1, np
375 p1 = 0
376 call mesh_local_index_to_coords(mesh, ii, p1)
377
378 do jj = 1, op%stencil%size
379 ! Get local index of p1 plus current stencil point.
380 st1(jj) = mesh_local_index_from_coords(mesh, p1 + op%stencil%points(:, jj))
381
382 ! if boundary conditions are zero, we can remap boundary
383 ! points to reduce memory accesses. We cannot do this for the
384 ! first point, since it is used to build the weights, so it
385 ! has to have the positions right
386 if (ii > 1 .and. compact_boundaries .and. mesh_compact_boundaries(mesh) .and. .not. space%is_periodic()) then
387 st1(jj) = min(st1(jj), mesh%np + 1)
388 end if
389 assert(st1(jj) > 0)
390 end do
391
392 st1(1:op%stencil%size) = st1(1:op%stencil%size) - ii
393
394 change = any(st1 /= st2)
395
396 !the next is to detect when we move from a point that does not
397 !have boundary points as neighbours to one that has
398 force_change = any(st1 + ii > mesh%np) .and. all(st2 + ii - 1 <= mesh%np)
399
400 if (change .and. compact_boundaries .and. mesh_compact_boundaries(mesh) .and. .not. space%is_periodic()) then
401 !try to repair it by changing the boundary points
402 do jj = 1, op%stencil%size
403 if (st1(jj) + ii > mesh%np .and. st2(jj) + ii - 1 > mesh%np .and. st2(jj) + ii <= mesh%np_part) then
404 st1r(jj) = st2(jj)
405 else
406 st1r(jj) = st1(jj)
407 end if
408 end do
409
410 change = any(st1r /= st2)
411
412 if (.not. change) st1(:) = st1r(:)
413 end if
414
415 ! if the stencil changes
416 if (change .or. force_change) then
417 !store it
418 st2(:) = st1(:)
419
420 !first time, just count
421 if (time == 1) op%nri = op%nri + 1
422
423 !second time, store
424 if (time == 2) then
425 current = current + 1
426 op%ri(1:op%stencil%size, current) = st1(1:op%stencil%size)
427 end if
428 end if
429
430 if (time == 2) op%rimap(ii) = current
431
432 end do
433
434 !after counting, allocate
435 if (time == 1) then
436 safe_deallocate_a(op%ri)
437 safe_deallocate_a(op%rimap)
438 safe_deallocate_a(op%rimap_inv)
439
440 safe_allocate(op%ri(1:op%stencil%size, 1:op%nri))
441 safe_allocate(op%rimap(1:op%np))
442 safe_allocate(op%rimap_inv(1:op%nri + 1))
443 op%ri = 0
444 op%rimap = 0
445 op%rimap_inv = 0
446 current = 0
447
448 ! the sizes
449 if (mesh%use_curvilinear) then
450 safe_allocate(op%nn(1:op%nri))
451 ! for the moment all the sizes are the same
452 op%nn = op%stencil%size
453 end if
454 end if
455
456 end do
457
458 !the inverse mapping
459 op%rimap_inv(1) = 0
460 do jj = 1, op%np
461 op%rimap_inv(op%rimap(jj) + 1) = jj
462 end do
463 op%rimap_inv(op%nri + 1) = op%np
464
465 safe_deallocate_a(st1)
466 safe_deallocate_a(st1r)
467 safe_deallocate_a(st2)
468
469 do jj = 1, op%nri
470 nn = op%rimap_inv(jj + 1) - op%rimap_inv(jj)
471 end do
472
473 if (op%mesh%parallel_in_domains) then
474 !now build the arrays required to apply the nl_operator by parts
475
476 !count points
477 op%inner%nri = 0
478 op%outer%nri = 0
479 do ir = 1, op%nri
480 maxp = op%rimap_inv(ir + 1) + maxval(op%ri(1:op%stencil%size, ir))
481 if (maxp <= np) then
482 !inner point
483 op%inner%nri = op%inner%nri + 1
484 assert(op%inner%nri <= op%nri)
485 else
486 !outer point
487 op%outer%nri = op%outer%nri + 1
488 assert(op%outer%nri <= op%nri)
489 end if
490 end do
491
492 assert(op%inner%nri + op%outer%nri == op%nri)
493
494 if (optional_default(regenerate, .false.)) then
495 safe_deallocate_a(op%inner%imin)
496 safe_deallocate_a(op%inner%imax)
497 safe_deallocate_a(op%inner%ri)
498 safe_deallocate_a(op%outer%imin)
499 safe_deallocate_a(op%outer%imax)
500 safe_deallocate_a(op%outer%ri)
501 end if
502 safe_allocate(op%inner%imin(1:op%inner%nri + 1))
503 safe_allocate(op%inner%imax(1:op%inner%nri))
504 safe_allocate(op%inner%ri(1:op%stencil%size, 1:op%inner%nri))
505
506 safe_allocate(op%outer%imin(1:op%outer%nri + 1))
507 safe_allocate(op%outer%imax(1:op%outer%nri))
508 safe_allocate(op%outer%ri(1:op%stencil%size, 1:op%outer%nri))
509
510 !now populate the arrays
511 iinner = 0
512 iouter = 0
513 do ir = 1, op%nri
514 maxp = op%rimap_inv(ir + 1) + maxval(op%ri(1:op%stencil%size, ir))
515 if (maxp <= np) then
516 !inner point
517 iinner = iinner + 1
518 op%inner%imin(iinner) = op%rimap_inv(ir)
519 op%inner%imax(iinner) = op%rimap_inv(ir + 1)
520 op%inner%ri(1:op%stencil%size, iinner) = op%ri(1:op%stencil%size, ir)
521 else
522 !outer point
523 iouter = iouter + 1
524 op%outer%imin(iouter) = op%rimap_inv(ir)
525 op%outer%imax(iouter) = op%rimap_inv(ir + 1)
526 op%outer%ri(1:op%stencil%size, iouter) = op%ri(1:op%stencil%size, ir)
527 end if
528 end do
529
530 !verify that all points in the inner operator are actually inner
531 do ir = 1, op%inner%nri
532 do ii = op%inner%imin(ir) + 1, op%inner%imax(ir)
533 assert(all(ii + op%inner%ri(1:op%stencil%size, ir) <= mesh%np))
534 end do
535 end do
536
537 end if
538
539 if (accel_is_enabled() .and. op%const_w) then
540
541 write(flags, '(i5)') op%stencil%size
542 flags='-DNDIM=3 -DSTENCIL_SIZE='//trim(adjustl(flags))
543
544 if (op%mesh%parallel_in_domains) flags = '-DINDIRECT '//trim(flags)
545
546 select case (function_opencl)
547 case (op_invmap)
548 call accel_kernel_build(op%kernel, 'operate.cl', 'operate', flags)
549 case (op_map)
550 call accel_kernel_build(op%kernel, 'operate.cl', 'operate_map', flags)
551 end select
552
553 ! conversion to i8 needed to avoid integer overflow
554 call accel_create_buffer(op%buff_ri, accel_mem_read_only, type_integer, int(op%nri, int64)*op%stencil%size)
555 call accel_write_buffer(op%buff_ri, int(op%nri, int64)*op%stencil%size, op%ri)
556
557 select case (function_opencl)
558 case (op_invmap)
559 call accel_create_buffer(op%buff_imin, accel_mem_read_only, type_integer, op%nri)
560 call accel_write_buffer(op%buff_imin, op%nri, op%rimap_inv(1:))
561 call accel_create_buffer(op%buff_imax, accel_mem_read_only, type_integer, op%nri)
562 call accel_write_buffer(op%buff_imax, op%nri, op%rimap_inv(2:))
563
564 case (op_map)
565
567 call accel_write_buffer(op%buff_map, op%mesh%np, (op%rimap - 1)*op%stencil%size)
568
569 if (op%mesh%parallel_in_domains) then
570
571 safe_allocate(inner_points(1:op%mesh%np))
572 safe_allocate(outer_points(1:op%mesh%np))
573 safe_allocate(all_points(1:op%mesh%np))
574
575 op%ninner = 0
576 op%nouter = 0
577
578 do ii = 1, op%mesh%np
579 all_points(ii) = ii - 1
580 maxp = ii + maxval(op%ri(1:op%stencil%size, op%rimap(ii)))
581 if (maxp <= op%mesh%np) then
582 op%ninner = op%ninner + 1
583 inner_points(op%ninner) = ii - 1
584 else
585 op%nouter = op%nouter + 1
586 outer_points(op%nouter) = ii - 1
587 end if
588 end do
589
591 call accel_write_buffer(op%buff_all, op%mesh%np, all_points)
592
594 call accel_write_buffer(op%buff_inner, op%ninner, inner_points)
595
597 call accel_write_buffer(op%buff_outer, op%nouter, outer_points)
598
599 safe_deallocate_a(inner_points)
600 safe_deallocate_a(outer_points)
601 safe_deallocate_a(all_points)
602
603 end if
604 end select
605 end if
606
607 pop_sub(nl_operator_build)
608
609 end subroutine nl_operator_build
610
611 ! ---------------------------------------------------------
612 subroutine nl_operator_output_weights(this)
613 type(nl_operator_t), intent(inout) :: this
614
615 integer :: istencil, idir
616
618
619 write(message(1), '(3a)') 'Debug info: Finite difference weights for ', trim(this%label), '.'
620 write(message(2), '(a)') ' Spacing:'
621 do idir = 1, this%mesh%box%dim
622 write(message(2), '(a,f16.8)') trim(message(2)), this%mesh%spacing(idir)
623 end do
624 call messages_info(2, debug_only=.true.)
625
626 do istencil = 1, this%stencil%size
627 select case(this%mesh%box%dim)
628 case(1)
629 write(message(1), '(a,i3,1i4,f25.10)') ' ', istencil, this%stencil%points(1:1, istencil), this%w(istencil, 1)
630 case(2)
631 write(message(1), '(a,i3,2i4,f25.10)') ' ', istencil, this%stencil%points(1:2, istencil), this%w(istencil, 1)
632 case(3)
633 write(message(1), '(a,i3,3i4,f25.10)') ' ', istencil, this%stencil%points(1:3, istencil), this%w(istencil, 1)
634 end select
635 call messages_info(1, debug_only=.true.)
636 end do
637
639
640 end subroutine nl_operator_output_weights
641
642 ! ---------------------------------------------------------
644 subroutine nl_operator_adjoint(op, opt, mesh, self_adjoint)
645 type(nl_operator_t), target, intent(in) :: op
646 type(nl_operator_t), target, intent(out) :: opt
647 type(mesh_t), target, intent(in) :: mesh
648 logical, intent(in) :: self_adjoint
649
650 integer :: ip, jp, kp, lp, index, ii, p1(1:mesh%box%dim)
651 real(real64), pointer, contiguous :: vol_pp(:), weights(:, :), tmp(:)
652 real(real64) :: factor
653
654 push_sub(nl_operator_adjoint)
655
656 if (self_adjoint) then
657 ! make operator self-adjoint
658 factor = m_one
659 else
660 ! make operator skew-adjoint
661 factor = -m_one
662 end if
663
664 call nl_operator_copy(opt, op)
665
666 if (mesh%parallel_in_domains) then
667 ! get ghost point values
668 safe_allocate(vol_pp(1:mesh%np_part))
669 vol_pp(1:mesh%np) = mesh%vol_pp(1:mesh%np)
670 call dpar_vec_ghost_update(mesh%pv, vol_pp)
671
672 if (.not.op%const_w) then
673 safe_allocate(weights(1:op%stencil%size, 1:mesh%np_part))
674 safe_allocate(tmp(1:mesh%np_part))
675 do ii = 1, op%stencil%size
676 tmp(1:mesh%np) = op%w(ii, 1:mesh%np)
677 call dpar_vec_ghost_update(mesh%pv, tmp)
678 weights(ii, 1:mesh%np_part) = tmp(1:mesh%np_part)
679 end do
680 safe_deallocate_p(tmp)
681 end if
682 else
683 vol_pp => mesh%vol_pp
684 weights => op%w
685 end if
686
687 if (.not.op%const_w) then
688 opt%w = m_zero
689 do ip = 1, mesh%np
690 do jp = 1, op%stencil%size
691 index = nl_operator_get_index(op, jp, ip)
692 if (index <= mesh%np) then
693 ! point is in the mesh
694 do lp = 1, op%stencil%size
695 kp = nl_operator_get_index(op, lp, index)
696 if (kp == ip) then
697 opt%w(jp, ip) = m_half*weights(jp, ip) + factor*m_half*(vol_pp(index)/vol_pp(ip))*weights(lp, index)
698 end if
699 end do
700 else if (index <= mesh%np + mesh%pv%np_ghost) then
701 ! point is in the ghost layer
702 call mesh_local_index_to_coords(mesh, index, p1)
703 do lp = 1, op%stencil%size
704 kp = mesh_local_index_from_coords(mesh, p1(1:mesh%box%dim) + &
705 op%stencil%points(1:mesh%box%dim, lp))
706 if (kp == ip) then
707 opt%w(jp, ip) = m_half*weights(jp, ip) + factor*m_half*(vol_pp(index)/vol_pp(ip))*weights(lp, index)
708 end if
709 end do
710 end if
711 end do
712 end do
713 else
714 opt%w(1:op%stencil%size, 1) = op%w(1:op%stencil%size, 1)
715 end if
716
718
719 if (mesh%parallel_in_domains) then
720 safe_deallocate_p(vol_pp)
721 safe_deallocate_p(weights)
722 end if
723
724 pop_sub(nl_operator_adjoint)
725 end subroutine nl_operator_adjoint
726
727 ! ---------------------------------------------------------
728 subroutine nl_operator_end(op)
729 type(nl_operator_t), intent(inout) :: op
730
731 push_sub(nl_operator_end)
732
733 if (accel_is_enabled() .and. op%const_w) then
734
735 call accel_release_buffer(op%buff_ri)
736 select case (function_opencl)
737 case (op_invmap)
738 call accel_release_buffer(op%buff_imin)
739 call accel_release_buffer(op%buff_imax)
740
741 case (op_map)
742 call accel_release_buffer(op%buff_map)
743 if (op%mesh%parallel_in_domains) then
744 call accel_release_buffer(op%buff_all)
745 call accel_release_buffer(op%buff_inner)
746 call accel_release_buffer(op%buff_outer)
747 end if
748
749 case (op_nomap)
750 call accel_release_buffer(op%buff_map)
751 call accel_release_buffer(op%buff_stencil)
752 call accel_release_buffer(op%buff_xyz_to_ip)
753 call accel_release_buffer(op%buff_ip_to_xyz)
754 end select
755
756 call accel_release_buffer(op%buff_weights)
757 call accel_release_buffer(op%buff_half_weights)
758 end if
759
760 safe_deallocate_a(op%inner%imin)
761 safe_deallocate_a(op%inner%imax)
762 safe_deallocate_a(op%inner%ri)
763 safe_deallocate_a(op%outer%imin)
764 safe_deallocate_a(op%outer%imax)
765 safe_deallocate_a(op%outer%ri)
766
767 safe_deallocate_a(op%w)
768
769 safe_deallocate_a(op%ri)
770 safe_deallocate_a(op%rimap)
771 safe_deallocate_a(op%rimap_inv)
772 safe_deallocate_a(op%nn)
773
774 call stencil_end(op%stencil)
775
776 pop_sub(nl_operator_end)
777 end subroutine nl_operator_end
778
779
780 ! ---------------------------------------------------------
781 integer pure function nl_operator_get_index(op, is, ip) result(res)
782 type(nl_operator_t), intent(in) :: op
783 integer, intent(in) :: is
784 integer, intent(in) :: ip
785
786 res = ip + op%ri(is, op%rimap(ip))
787 end function nl_operator_get_index
788
789 ! ---------------------------------------------------------
790
792 type(nl_operator_t), intent(inout) :: op
793
795
796 ! Update the GPU weights
797 if (accel_is_enabled() .and. op%const_w) then
798 call accel_create_buffer(op%buff_weights, accel_mem_read_only, type_float, op%stencil%size)
799 call accel_create_buffer(op%buff_half_weights, accel_mem_read_only, type_float, op%stencil%size)
800 end if
801
804
805
806 ! ---------------------------------------------------------
807
808 subroutine nl_operator_update_gpu_buffers(op)
809 type(nl_operator_t), intent(inout) :: op
810
812
813 ! Update the GPU weights
814 if (accel_is_enabled() .and. op%const_w) then
815 call accel_write_buffer(op%buff_weights, op%stencil%size, op%w(:, 1))
816 call accel_write_buffer(op%buff_half_weights, op%stencil%size, -m_half*op%w(:, 1))
817 end if
818
820 end subroutine nl_operator_update_gpu_buffers
822 ! ---------------------------------------------------------
823
824 integer pure function nl_operator_np_zero_bc(op) result(np_bc)
825 type(nl_operator_t), intent(in) :: op
826
827 integer :: jj, ii
828
829 np_bc = 0
830 do jj = 1, op%nri
831 ii = op%rimap_inv(jj + 1) + maxval(op%ri(1:op%stencil%size, jj))
832 np_bc = max(np_bc, ii)
833 end do
834
835 end function nl_operator_np_zero_bc
836
837 ! ------------------------------------------------------
838
839 logical pure function nl_operator_compact_boundaries()
840
843
844
845#include "undef.F90"
846#include "real.F90"
847#include "nl_operator_inc.F90"
848
849#include "undef.F90"
850#include "complex.F90"
851#include "nl_operator_inc.F90"
852
853end module nl_operator_oct_m
854
855!! Local Variables:
856!! mode: f90
857!! coding: utf-8
858!! End:
subroutine, public accel_release_buffer(this)
Definition: accel.F90:1248
subroutine, public accel_kernel_build(this, file_name, kernel_name, flags)
Definition: accel.F90:2053
pure logical function, public accel_is_enabled()
Definition: accel.F90:401
integer, parameter, public accel_mem_read_only
Definition: accel.F90:183
integer pure function, public accel_max_workgroup_size()
Definition: accel.F90:1474
This module implements batches of mesh functions.
Definition: batch.F90:133
Module implementing boundary conditions in Octopus.
Definition: boundaries.F90:122
subroutine, public dpar_vec_ghost_update(pv, v_local)
Updates ghost points of every node.
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public m_half
Definition: global.F90:194
real(real64), parameter, public m_one
Definition: global.F90:189
This module implements the index, used for the mesh points.
Definition: index.F90:122
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
subroutine, public weights(N, M, cc, side)
Compute the weights for finite-difference calculations:
Definition: math.F90:526
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
integer function, public mesh_local_index_from_coords(mesh, ix)
This function returns the local index of the point for a given vector of integer coordinates.
Definition: mesh.F90:939
subroutine, public mesh_local_index_to_coords(mesh, ip, ix)
Given a local point index, this function returns the set of integer coordinates of the point.
Definition: mesh.F90:951
logical pure function, public mesh_compact_boundaries(mesh)
Definition: mesh.F90:826
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1046
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:714
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1086
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:617
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:145
This module defines non-local operators.
subroutine, public dnl_operator_operate_diag(op, fo)
subroutine, public nl_operator_init(op, label)
initialize an instance of a non-local operator by setting the label
integer, parameter op_map
integer, parameter op_max
integer, parameter op_vec
subroutine, public dnl_operator_operate_batch(op, fi, fo, ghost_update, profile, points, factor, async)
subroutine, public dnl_operator_operate(op, fi, fo, ghost_update, profile, points)
subroutine, public nl_operator_update_gpu_buffers(op)
subroutine, public nl_operator_global_init(namespace)
initialize global settings for non-local operators
subroutine, public nl_operator_copy(opo, opi)
subroutine, public nl_operator_output_weights(this)
subroutine, public nl_operator_end(op)
integer, parameter, public op_inner
subroutine, public znl_operator_operate(op, fi, fo, ghost_update, profile, points)
logical pure function, public nl_operator_compact_boundaries()
subroutine, public nl_operator_global_end()
integer, parameter, public op_outer
subroutine, public nl_operator_build(space, mesh, op, np, const_w, regenerate)
logical compact_boundaries
subroutine, public znl_operator_operate_batch(op, fi, fo, ghost_update, profile, points, factor, async)
integer, parameter op_nomap
subroutine, public nl_operator_adjoint(op, opt, mesh, self_adjoint)
opt has to be initialised and built.
integer pure function, public nl_operator_np_zero_bc(op)
subroutine, public znl_operator_operate_diag(op, fo)
integer pure function, public nl_operator_get_index(op, is, ip)
subroutine, public nl_operator_allocate_gpu_buffers(op)
integer, parameter op_min
This module contains interfaces for routines in operate.c.
Definition: operate_f.F90:117
Some general things and nomenclature:
Definition: par_vec.F90:171
This module defines stencils used in Octopus.
Definition: stencil.F90:135
subroutine, public stencil_end(this)
Definition: stencil.F90:214
subroutine, public stencil_copy(input, output)
Definition: stencil.F90:196
type(type_t), parameter, public type_integer
Definition: types.F90:135
type(type_t), parameter, public type_float
Definition: types.F90:133
Describes mesh distribution to nodes.
Definition: mesh.F90:186
index type for non-local operators
data type for non local operators
int true(void)