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 if (mesh%parallel_in_domains .and. .not. const_w) then
342 call messages_experimental('Domain parallelization with curvilinear coordinates')
343 end if
344
345 assert(np > 0)
346
347 ! store values in structure
348 op%np = np
349 op%mesh => mesh
350 op%const_w = .false.
351 if (present(const_w)) op%const_w = const_w
352
353 ! allocate weights op%w
354 if (op%const_w) then
355 safe_allocate(op%w(1:op%stencil%size, 1))
356 if (debug%info) then
357 message(1) = 'Info: nl_operator_build: working with constant weights.'
358 call messages_info(1)
359 end if
360 else
361 safe_allocate(op%w(1:op%stencil%size, 1:op%np))
362 if (debug%info) then
363 message(1) = 'Info: nl_operator_build: working with non-constant weights.'
364 call messages_info(1)
365 end if
366 end if
367
368 ! set initially to zero
369 op%w = m_zero
370
371 ! Build lookup table
372 safe_allocate(st1(1:op%stencil%size))
373 safe_allocate(st1r(1:op%stencil%size))
374 safe_allocate(st2(1:op%stencil%size))
375
376 op%nri = 0
377 do time = 1, 2
378 st2 = 0
379 do ii = 1, np
380 p1 = 0
381 call mesh_local_index_to_coords(mesh, ii, p1)
382
383 do jj = 1, op%stencil%size
384 ! Get local index of p1 plus current stencil point.
385 st1(jj) = mesh_local_index_from_coords(mesh, p1 + op%stencil%points(:, jj))
386
387 ! if boundary conditions are zero, we can remap boundary
388 ! points to reduce memory accesses. We cannot do this for the
389 ! first point, since it is used to build the weights, so it
390 ! has to have the positions right
391 if (ii > 1 .and. compact_boundaries .and. mesh_compact_boundaries(mesh) .and. .not. space%is_periodic()) then
392 st1(jj) = min(st1(jj), mesh%np + 1)
393 end if
394 assert(st1(jj) > 0)
395 end do
396
397 st1(1:op%stencil%size) = st1(1:op%stencil%size) - ii
398
399 change = any(st1 /= st2)
400
401 !the next is to detect when we move from a point that does not
402 !have boundary points as neighbours to one that has
403 force_change = any(st1 + ii > mesh%np) .and. all(st2 + ii - 1 <= mesh%np)
404
405 if (change .and. compact_boundaries .and. mesh_compact_boundaries(mesh) .and. .not. space%is_periodic()) then
406 !try to repair it by changing the boundary points
407 do jj = 1, op%stencil%size
408 if (st1(jj) + ii > mesh%np .and. st2(jj) + ii - 1 > mesh%np .and. st2(jj) + ii <= mesh%np_part) then
409 st1r(jj) = st2(jj)
410 else
411 st1r(jj) = st1(jj)
412 end if
413 end do
414
415 change = any(st1r /= st2)
417 if (.not. change) st1 = st1r
418 end if
419
420 ! if the stencil changes
421 if (change .or. force_change) then
422 !store it
423 st2 = st1
424
425 !first time, just count
426 if (time == 1) op%nri = op%nri + 1
427
428 !second time, store
429 if (time == 2) then
430 current = current + 1
431 op%ri(1:op%stencil%size, current) = st1(1:op%stencil%size)
432 end if
433 end if
434
435 if (time == 2) op%rimap(ii) = current
436
437 end do
438
439 !after counting, allocate
440 if (time == 1) then
441 safe_deallocate_a(op%ri)
442 safe_deallocate_a(op%rimap)
443 safe_deallocate_a(op%rimap_inv)
444
445 safe_allocate(op%ri(1:op%stencil%size, 1:op%nri))
446 safe_allocate(op%rimap(1:op%np))
447 safe_allocate(op%rimap_inv(1:op%nri + 1))
448 op%ri = 0
449 op%rimap = 0
450 op%rimap_inv = 0
451 current = 0
452
453 ! the sizes
454 if (mesh%use_curvilinear) then
455 safe_allocate(op%nn(1:op%nri))
456 ! for the moment all the sizes are the same
457 op%nn = op%stencil%size
458 end if
459 end if
460
461 end do
462
463 !the inverse mapping
464 op%rimap_inv(1) = 0
465 do jj = 1, op%np
466 op%rimap_inv(op%rimap(jj) + 1) = jj
467 end do
468 op%rimap_inv(op%nri + 1) = op%np
469
470 safe_deallocate_a(st1)
471 safe_deallocate_a(st1r)
472 safe_deallocate_a(st2)
473
474 do jj = 1, op%nri
475 nn = op%rimap_inv(jj + 1) - op%rimap_inv(jj)
476 end do
477
478 if (op%mesh%parallel_in_domains) then
479 !now build the arrays required to apply the nl_operator by parts
480
481 !count points
482 op%inner%nri = 0
483 op%outer%nri = 0
484 do ir = 1, op%nri
485 maxp = op%rimap_inv(ir + 1) + maxval(op%ri(1:op%stencil%size, ir))
486 if (maxp <= np) then
487 !inner point
488 op%inner%nri = op%inner%nri + 1
489 assert(op%inner%nri <= op%nri)
490 else
491 !outer point
492 op%outer%nri = op%outer%nri + 1
493 assert(op%outer%nri <= op%nri)
494 end if
495 end do
496
497 assert(op%inner%nri + op%outer%nri == op%nri)
498
499 if (optional_default(regenerate, .false.)) then
500 safe_deallocate_a(op%inner%imin)
501 safe_deallocate_a(op%inner%imax)
502 safe_deallocate_a(op%inner%ri)
503 safe_deallocate_a(op%outer%imin)
504 safe_deallocate_a(op%outer%imax)
505 safe_deallocate_a(op%outer%ri)
506 end if
507 safe_allocate(op%inner%imin(1:op%inner%nri + 1))
508 safe_allocate(op%inner%imax(1:op%inner%nri))
509 safe_allocate(op%inner%ri(1:op%stencil%size, 1:op%inner%nri))
510
511 safe_allocate(op%outer%imin(1:op%outer%nri + 1))
512 safe_allocate(op%outer%imax(1:op%outer%nri))
513 safe_allocate(op%outer%ri(1:op%stencil%size, 1:op%outer%nri))
514
515 !now populate the arrays
516 iinner = 0
517 iouter = 0
518 do ir = 1, op%nri
519 maxp = op%rimap_inv(ir + 1) + maxval(op%ri(1:op%stencil%size, ir))
520 if (maxp <= np) then
521 !inner point
522 iinner = iinner + 1
523 op%inner%imin(iinner) = op%rimap_inv(ir)
524 op%inner%imax(iinner) = op%rimap_inv(ir + 1)
525 op%inner%ri(1:op%stencil%size, iinner) = op%ri(1:op%stencil%size, ir)
526 else
527 !outer point
528 iouter = iouter + 1
529 op%outer%imin(iouter) = op%rimap_inv(ir)
530 op%outer%imax(iouter) = op%rimap_inv(ir + 1)
531 op%outer%ri(1:op%stencil%size, iouter) = op%ri(1:op%stencil%size, ir)
532 end if
533 end do
534
535 !verify that all points in the inner operator are actually inner
536 do ir = 1, op%inner%nri
537 do ii = op%inner%imin(ir) + 1, op%inner%imax(ir)
538 assert(all(ii + op%inner%ri(1:op%stencil%size, ir) <= mesh%np))
539 end do
540 end do
541
542 end if
543
544 if (accel_is_enabled() .and. op%const_w) then
545
546 write(flags, '(i5)') op%stencil%size
547 flags='-DNDIM=3 -DSTENCIL_SIZE='//trim(adjustl(flags))
548
549 if (op%mesh%parallel_in_domains) flags = '-DINDIRECT '//trim(flags)
550
551 select case (function_opencl)
552 case (op_invmap)
553 call accel_kernel_build(op%kernel, 'operate.cl', 'operate', flags)
554 case (op_map)
555 call accel_kernel_build(op%kernel, 'operate.cl', 'operate_map', flags)
556 end select
557
558 ! conversion to i8 needed to avoid integer overflow
559 call accel_create_buffer(op%buff_ri, accel_mem_read_only, type_integer, int(op%nri, int64)*op%stencil%size)
560 call accel_write_buffer(op%buff_ri, int(op%nri, int64)*op%stencil%size, op%ri)
561
562 select case (function_opencl)
563 case (op_invmap)
564 call accel_create_buffer(op%buff_imin, accel_mem_read_only, type_integer, op%nri)
565 call accel_write_buffer(op%buff_imin, op%nri, op%rimap_inv(1:))
566 call accel_create_buffer(op%buff_imax, accel_mem_read_only, type_integer, op%nri)
567 call accel_write_buffer(op%buff_imax, op%nri, op%rimap_inv(2:))
568
569 case (op_map)
570
572 call accel_write_buffer(op%buff_map, op%mesh%np, (op%rimap - 1)*op%stencil%size)
573
574 if (op%mesh%parallel_in_domains) then
575
576 safe_allocate(inner_points(1:op%mesh%np))
577 safe_allocate(outer_points(1:op%mesh%np))
578 safe_allocate(all_points(1:op%mesh%np))
579
580 op%ninner = 0
581 op%nouter = 0
582
583 do ii = 1, op%mesh%np
584 all_points(ii) = ii - 1
585 maxp = ii + maxval(op%ri(1:op%stencil%size, op%rimap(ii)))
586 if (maxp <= op%mesh%np) then
587 op%ninner = op%ninner + 1
588 inner_points(op%ninner) = ii - 1
589 else
590 op%nouter = op%nouter + 1
591 outer_points(op%nouter) = ii - 1
592 end if
593 end do
594
596 call accel_write_buffer(op%buff_all, op%mesh%np, all_points)
597
599 call accel_write_buffer(op%buff_inner, op%ninner, inner_points)
600
602 call accel_write_buffer(op%buff_outer, op%nouter, outer_points)
603
604 safe_deallocate_a(inner_points)
605 safe_deallocate_a(outer_points)
606 safe_deallocate_a(all_points)
607
608 end if
609 end select
610 end if
611
612 pop_sub(nl_operator_build)
613
614 end subroutine nl_operator_build
615
616 ! ---------------------------------------------------------
617 subroutine nl_operator_output_weights(this)
618 type(nl_operator_t), intent(inout) :: this
619
620 integer :: istencil, idir
621
623
624 if (debug%info) then
625
626 write(message(1), '(3a)') 'Debug info: Finite difference weights for ', trim(this%label), '.'
627 write(message(2), '(a)') ' Spacing:'
628 do idir = 1, this%mesh%box%dim
629 write(message(2), '(a,f16.8)') trim(message(2)), this%mesh%spacing(idir)
630 end do
631 call messages_info(2)
632
633 do istencil = 1, this%stencil%size
634 select case(this%mesh%box%dim)
635 case(1)
636 write(message(1), '(a,i3,1i4,f25.10)') ' ', istencil, this%stencil%points(1:1, istencil), this%w(istencil, 1)
637 case(2)
638 write(message(1), '(a,i3,2i4,f25.10)') ' ', istencil, this%stencil%points(1:2, istencil), this%w(istencil, 1)
639 case(3)
640 write(message(1), '(a,i3,3i4,f25.10)') ' ', istencil, this%stencil%points(1:3, istencil), this%w(istencil, 1)
641 end select
642 call messages_info(1)
643 end do
644
645 end if
646
648
649 end subroutine nl_operator_output_weights
650
651 ! ---------------------------------------------------------
653 subroutine nl_operator_adjoint(op, opt, mesh, self_adjoint)
654 type(nl_operator_t), target, intent(in) :: op
655 type(nl_operator_t), target, intent(out) :: opt
656 type(mesh_t), target, intent(in) :: mesh
657 logical, intent(in) :: self_adjoint
658
659 integer :: ip, jp, kp, lp, index, ii, p1(1:mesh%box%dim)
660 real(real64), pointer, contiguous :: vol_pp(:), weights(:, :), tmp(:)
661 real(real64) :: factor
662
663 push_sub(nl_operator_adjoint)
664
665 if (self_adjoint) then
666 ! make operator self-adjoint
667 factor = m_one
668 else
669 ! make operator skew-adjoint
670 factor = -m_one
671 end if
672
673 call nl_operator_copy(opt, op)
674
675 if (mesh%parallel_in_domains) then
676 ! get ghost point values
677 safe_allocate(vol_pp(1:mesh%np_part))
678 vol_pp(1:mesh%np) = mesh%vol_pp(1:mesh%np)
679 call dpar_vec_ghost_update(mesh%pv, vol_pp)
680
681 if (.not.op%const_w) then
682 safe_allocate(weights(1:op%stencil%size, 1:mesh%np_part))
683 safe_allocate(tmp(1:mesh%np_part))
684 do ii = 1, op%stencil%size
685 tmp(1:mesh%np) = op%w(ii, 1:mesh%np)
686 call dpar_vec_ghost_update(mesh%pv, tmp)
687 weights(ii, 1:mesh%np_part) = tmp(1:mesh%np_part)
688 end do
689 safe_deallocate_p(tmp)
690 end if
691 else
692 vol_pp => mesh%vol_pp
693 weights => op%w
694 end if
695
696 if (.not.op%const_w) then
697 opt%w = m_zero
698 do ip = 1, mesh%np
699 do jp = 1, op%stencil%size
700 index = nl_operator_get_index(op, jp, ip)
701 if (index <= mesh%np) then
702 ! point is in the mesh
703 do lp = 1, op%stencil%size
704 kp = nl_operator_get_index(op, lp, index)
705 if (kp == ip) then
706 opt%w(jp, ip) = m_half*weights(jp, ip) + factor*m_half*(vol_pp(index)/vol_pp(ip))*weights(lp, index)
707 end if
708 end do
709 else if (index <= mesh%np + mesh%pv%np_ghost) then
710 ! point is in the ghost layer
711 call mesh_local_index_to_coords(mesh, index, p1)
712 do lp = 1, op%stencil%size
713 kp = mesh_local_index_from_coords(mesh, p1(1:mesh%box%dim) + &
714 op%stencil%points(1:mesh%box%dim, lp))
715 if (kp == ip) then
716 opt%w(jp, ip) = m_half*weights(jp, ip) + factor*m_half*(vol_pp(index)/vol_pp(ip))*weights(lp, index)
717 end if
718 end do
719 end if
720 end do
721 end do
722 else
723 opt%w(1:op%stencil%size, 1) = op%w(1:op%stencil%size, 1)
724 end if
725
727
728 if (mesh%parallel_in_domains) then
729 safe_deallocate_p(vol_pp)
730 safe_deallocate_p(weights)
731 end if
732
733 pop_sub(nl_operator_adjoint)
734 end subroutine nl_operator_adjoint
735
736 ! ---------------------------------------------------------
737 subroutine nl_operator_end(op)
738 type(nl_operator_t), intent(inout) :: op
739
740 push_sub(nl_operator_end)
741
742 if (accel_is_enabled() .and. op%const_w) then
743
744 call accel_release_buffer(op%buff_ri)
745 select case (function_opencl)
746 case (op_invmap)
747 call accel_release_buffer(op%buff_imin)
748 call accel_release_buffer(op%buff_imax)
749
750 case (op_map)
751 call accel_release_buffer(op%buff_map)
752 if (op%mesh%parallel_in_domains) then
753 call accel_release_buffer(op%buff_all)
754 call accel_release_buffer(op%buff_inner)
755 call accel_release_buffer(op%buff_outer)
756 end if
757
758 case (op_nomap)
759 call accel_release_buffer(op%buff_map)
760 call accel_release_buffer(op%buff_stencil)
761 call accel_release_buffer(op%buff_xyz_to_ip)
762 call accel_release_buffer(op%buff_ip_to_xyz)
763 end select
764
765 call accel_release_buffer(op%buff_weights)
766 call accel_release_buffer(op%buff_half_weights)
767 end if
768
769 safe_deallocate_a(op%inner%imin)
770 safe_deallocate_a(op%inner%imax)
771 safe_deallocate_a(op%inner%ri)
772 safe_deallocate_a(op%outer%imin)
773 safe_deallocate_a(op%outer%imax)
774 safe_deallocate_a(op%outer%ri)
775
776 safe_deallocate_a(op%w)
777
778 safe_deallocate_a(op%ri)
779 safe_deallocate_a(op%rimap)
780 safe_deallocate_a(op%rimap_inv)
781 safe_deallocate_a(op%nn)
782
783 call stencil_end(op%stencil)
784
785 pop_sub(nl_operator_end)
786 end subroutine nl_operator_end
787
788
789 ! ---------------------------------------------------------
790 integer pure function nl_operator_get_index(op, is, ip) result(res)
791 type(nl_operator_t), intent(in) :: op
792 integer, intent(in) :: is
793 integer, intent(in) :: ip
794
795 res = ip + op%ri(is, op%rimap(ip))
796 end function nl_operator_get_index
797
798 ! ---------------------------------------------------------
799
801 type(nl_operator_t), intent(inout) :: op
802
804
805 ! Update the GPU weights
806 if (accel_is_enabled() .and. op%const_w) then
807 call accel_create_buffer(op%buff_weights, accel_mem_read_only, type_float, op%stencil%size)
808 call accel_create_buffer(op%buff_half_weights, accel_mem_read_only, type_float, op%stencil%size)
809 end if
810
813
814
815 ! ---------------------------------------------------------
816
817 subroutine nl_operator_update_gpu_buffers(op)
818 type(nl_operator_t), intent(inout) :: op
819
821
822 ! Update the GPU weights
823 if (accel_is_enabled() .and. op%const_w) then
824 call accel_write_buffer(op%buff_weights, op%stencil%size, op%w(:, 1))
825 call accel_write_buffer(op%buff_half_weights, op%stencil%size, -m_half*op%w(:, 1))
826 end if
827
829 end subroutine nl_operator_update_gpu_buffers
831 ! ---------------------------------------------------------
832
833 integer pure function nl_operator_np_zero_bc(op) result(np_bc)
834 type(nl_operator_t), intent(in) :: op
835
836 integer :: jj, ii
837
838 np_bc = 0
839 do jj = 1, op%nri
840 ii = op%rimap_inv(jj + 1) + maxval(op%ri(1:op%stencil%size, jj))
841 np_bc = max(np_bc, ii)
842 end do
843
844 end function nl_operator_np_zero_bc
845
846 ! ------------------------------------------------------
847
848 logical pure function nl_operator_compact_boundaries()
849
852
853
854#include "undef.F90"
855#include "real.F90"
856#include "nl_operator_inc.F90"
857
858#include "undef.F90"
859#include "complex.F90"
860#include "nl_operator_inc.F90"
861
862end module nl_operator_oct_m
863
864!! Local Variables:
865!! mode: f90
866!! coding: utf-8
867!! End:
subroutine, public accel_release_buffer(this)
Definition: accel.F90:1246
subroutine, public accel_kernel_build(this, file_name, kernel_name, flags)
Definition: accel.F90:2067
pure logical function, public accel_is_enabled()
Definition: accel.F90:400
integer, parameter, public accel_mem_read_only
Definition: accel.F90:183
integer pure function, public accel_max_workgroup_size()
Definition: accel.F90:1472
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.
type(debug_t), save, public debug
Definition: debug.F90:154
real(real64), parameter, public m_zero
Definition: global.F90:187
real(real64), parameter, public m_half
Definition: global.F90:193
real(real64), parameter, public m_one
Definition: global.F90:188
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:522
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:935
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:947
logical pure function, public mesh_compact_boundaries(mesh)
Definition: mesh.F90:822
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1057
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
Definition: messages.F90:624
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:723
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1097
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), public type_float
Definition: types.F90:133
type(type_t), public type_integer
Definition: types.F90:135
Describes mesh distribution to nodes.
Definition: mesh.F90:186
index type for non-local operators
data type for non local operators
int true(void)