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 :: &
70
76 private
77 integer :: nri = 0
78 integer, allocatable :: imin(:)
79 integer, allocatable :: imax(:)
80 integer, allocatable :: ri(:, :)
81 end type nl_operator_index_t
82
84 type nl_operator_t
85 private
86 type(stencil_t), public :: stencil
87 type(mesh_t), pointer :: mesh => null()
88 integer, allocatable :: nn(:)
89 integer, public :: np = 0
90 ! When running in parallel mode, the next three arrays are unique on each node.
91 real(real64), allocatable, public :: w(:,:)
92
93 logical, public :: const_w = .true.
94
95 type(accel_mem_t), public :: buff_weights
96 type(accel_mem_t), public :: buff_half_weights
97
98 character(len=40) :: label
99
101 integer, public :: nri = 0
102 integer, allocatable, public :: ri(:,:)
103 integer, allocatable, public :: rimap(:)
104 integer, allocatable, public :: rimap_inv(:)
105
106 integer :: ninner = 0
107 integer :: nouter = 0
108
109 type(nl_operator_index_t) :: inner
110 type(nl_operator_index_t) :: outer
111
112 type(accel_kernel_t) :: kernel
113 type(accel_mem_t) :: buff_imin
114 type(accel_mem_t) :: buff_imax
115 type(accel_mem_t) :: buff_ri
116 type(accel_mem_t) :: buff_map
117 type(accel_mem_t) :: buff_all
118 type(accel_mem_t) :: buff_inner
119 type(accel_mem_t) :: buff_outer
120 type(accel_mem_t) :: buff_stencil
121 type(accel_mem_t) :: buff_ip_to_xyz
122 type(accel_mem_t) :: buff_xyz_to_ip
123
124 ! For multigrid solvers
125 type(nl_operator_t), public, pointer :: coarser => null()
126
127 end type nl_operator_t
128
129 integer, parameter :: &
130 OP_FORTRAN = 0, &
131 op_vec = 1, &
132 op_min = op_fortran, &
133 op_max = op_vec
134
135 integer, parameter :: &
136 OP_INVMAP = 1, &
137 op_map = 2, &
138 op_nomap = 3
139
140 integer, public, parameter :: OP_ALL = 3, op_inner = 1, op_outer = 2
141
142 interface
143 integer function op_is_available(opid, type)
144 implicit none
145 integer, intent(in) :: opid, type
146 end function op_is_available
147 end interface
148
149 integer :: dfunction_global = -1
150 integer :: zfunction_global = -1
151 integer :: function_opencl
152
153contains
154
155 ! ---------------------------------------------------------
158 subroutine nl_operator_global_init(namespace)
159 type(namespace_t), intent(in) :: namespace
160
161 integer :: default
162
164
165 !%Variable OperateDouble
166 !%Type integer
167 !%Section Execution::Optimization
168 !%Default optimized
169 !%Description
170 !% This variable selects the subroutine used to apply non-local
171 !% operators over the grid for real functions.
172 !%Option fortran 0
173 !% The standard Fortran function.
174 !%Option optimized 1
175 !% This version is optimized using vector primitives (if available).
176 !%End
177
178 !%Variable OperateComplex
179 !%Type integer
180 !%Section Execution::Optimization
181 !%Default optimized
182 !%Description
183 !% This variable selects the subroutine used to apply non-local
184 !% operators over the grid for complex functions.
185 !%Option fortran 0
186 !% The standard Fortran function.
187 !%Option optimized 1
188 !% This version is optimized using vector primitives (if available).
189 !%End
191 default = op_vec
192
193 call parse_variable(namespace, 'OperateDouble', default, dfunction_global)
194 if (.not. varinfo_valid_option('OperateDouble', dfunction_global)) call messages_input_error(namespace, 'OperateDouble')
195
196 call parse_variable(namespace, 'OperateComplex', default, zfunction_global)
197 if (.not. varinfo_valid_option('OperateComplex', zfunction_global)) call messages_input_error(namespace, 'OperateComplex')
200 if (accel_is_enabled()) then
202 !%Variable OperateAccel
203 !%Type integer
204 !%Default map
205 !%Section Execution::Optimization
206 !%Description
207 !% This variable selects the subroutine used to apply non-local
208 !% operators over the grid when an accelerator device is used.
209 !%Option invmap 1
210 !% The standard implementation ported to OpenCL.
211 !%Option map 2
212 !% A different version, more suitable for GPUs.
213 !%End
214 call parse_variable(namespace, 'OperateAccel', op_map, function_opencl)
216 call messages_obsolete_variable(namespace, 'OperateOpenCL', 'OperateAccel')
218 end if
219
221 end subroutine nl_operator_global_init
222
223 ! ---------------------------------------------------------
225 subroutine nl_operator_global_end()
226 push_sub(nl_operator_global_end)
227
229 end subroutine nl_operator_global_end
231 ! ---------------------------------------------------------
235 subroutine nl_operator_init(op, label)
236 type(nl_operator_t), intent(inout) :: op
237 character(len=*), intent(in) :: label
239 push_sub(nl_operator_init)
240
241 op%label = label
242
243 pop_sub(nl_operator_init)
244 end subroutine nl_operator_init
247 ! ---------------------------------------------------------
248 subroutine nl_operator_copy(opo, opi)
249 type(nl_operator_t), intent(inout) :: opo
250 type(nl_operator_t), target, intent(in) :: opi
251
252 push_sub(nl_operator_copy)
254 ! We cannot currently copy the GPU kernel for the nl_operator
255 assert(.not. accel_is_enabled())
256
257 call nl_operator_end(opo)
258 call nl_operator_init(opo, opi%label)
259
260 call stencil_copy(opi%stencil, opo%stencil)
261
262 opo%np = opi%np
263 opo%mesh => opi%mesh
264
265 safe_allocate_source_a(opo%nn, opi%nn)
266 safe_allocate_source_a(opo%w, opi%w)
267
268 opo%const_w = opi%const_w
269
270 opo%nri = opi%nri
271 assert(allocated(opi%ri))
272
273 safe_allocate_source_a(opo%ri, opi%ri)
274 safe_allocate_source_a(opo%rimap, opi%rimap)
275 safe_allocate_source_a(opo%rimap_inv, opi%rimap_inv)
276
277 if (opi%mesh%parallel_in_domains) then
278 opo%inner%nri = opi%inner%nri
279 safe_allocate_source_a(opo%inner%imin, opi%inner%imin)
280 safe_allocate_source_a(opo%inner%imax, opi%inner%imax)
281 safe_allocate_source_a(opo%inner%ri, opi%inner%ri)
282
283 opo%outer%nri = opi%outer%nri
284 safe_allocate_source_a(opo%outer%imin, opi%outer%imin)
285 safe_allocate_source_a(opo%outer%imax, opi%outer%imax)
286 safe_allocate_source_a(opo%outer%ri, opi%outer%ri)
287 end if
288
289 ! We create the corresponding GPU buffers
290 if (accel_is_enabled() .and. opo%const_w) then
291 call accel_create_buffer(opo%buff_weights, accel_mem_read_only, type_float, opo%stencil%size)
292 call accel_write_buffer(opo%buff_weights, opo%stencil%size, opo%w(:, 1))
293 call accel_create_buffer(opo%buff_half_weights, accel_mem_read_only, type_float, opo%stencil%size)
294 call accel_write_buffer(opo%buff_half_weights, opo%stencil%size, -m_half*opo%w(:, 1))
295 end if
296
297
298 pop_sub(nl_operator_copy)
299 end subroutine nl_operator_copy
300
301
302 ! ---------------------------------------------------------
312 subroutine nl_operator_build(space, mesh, op, np, const_w, regenerate)
313 class(space_t), intent(in) :: space
314 type(mesh_t), target, intent(in) :: mesh
315 type(nl_operator_t), intent(inout) :: op
316 integer, intent(in) :: np
317 logical, optional, intent(in) :: const_w
318 logical, optional, intent(in) :: regenerate
319
320 integer :: ii, jj, p1(space%dim), time, current
321 integer, allocatable :: st1(:), st2(:), st1r(:)
322 integer :: ir, maxp, iinner, iouter
323 logical :: change, force_change
324 character(len=200) :: flags
325 integer, allocatable :: inner_points(:), outer_points(:), all_points(:)
326
327 push_sub(nl_operator_build)
328
329 if (mesh%parallel_in_domains .and. .not. const_w) then
330 call messages_experimental('Domain parallelization with curvilinear coordinates')
331 end if
332
333 assert(np > 0)
334
335 ! store values in structure
336 op%np = np
337 op%mesh => mesh
338 if (.not. optional_default(regenerate, .false.)) then
339 op%const_w = optional_default(const_w, .false.)
340 end if
341
342 ! allocate weights op%w
343 if (op%const_w) then
344 safe_allocate(op%w(1:op%stencil%size, 1))
345 message(1) = 'Debug: nl_operator_build: working with constant weights.'
346 call messages_info(1, debug_only=.true.)
347 else
348 safe_allocate(op%w(1:op%stencil%size, 1:op%np))
349 message(1) = 'Debug: nl_operator_build: working with non-constant weights.'
350 call messages_info(1, debug_only=.true.)
351 end if
352
353 ! set initially to zero
354 op%w = m_zero
355
356 ! Build lookup table
357 safe_allocate(st1(1:op%stencil%size))
358 safe_allocate(st1r(1:op%stencil%size))
359 safe_allocate(st2(1:op%stencil%size))
360
361 op%nri = 0
362 do time = 1, 2
363 st2 = 0
364 do ii = 1, np
365 p1 = 0
366 call mesh_local_index_to_coords(mesh, ii, p1)
367
368 do jj = 1, op%stencil%size
369 ! Get local index of p1 plus current stencil point.
370 st1(jj) = mesh_local_index_from_coords(mesh, p1 + op%stencil%points(:, jj))
371
372 assert(st1(jj) > 0)
373 end do
374
375 st1(1:op%stencil%size) = st1(1:op%stencil%size) - ii
376
377 change = any(st1 /= st2)
378
379 !the next is to detect when we move from a point that does not
380 !have boundary points as neighbours to one that has
381 force_change = any(st1 + ii > mesh%np) .and. all(st2 + ii - 1 <= mesh%np)
382
383 ! if the stencil changes
384 if (change .or. force_change) then
385 !store it
386 st2(:) = st1(:)
387
388 !first time, just count
389 if (time == 1) op%nri = op%nri + 1
390
391 !second time, store
392 if (time == 2) then
393 current = current + 1
394 op%ri(1:op%stencil%size, current) = st1(1:op%stencil%size)
395 end if
396 end if
397
398 if (time == 2) op%rimap(ii) = current
399
400 end do
401
402 !after counting, allocate
403 if (time == 1) then
404 safe_deallocate_a(op%ri)
405 safe_deallocate_a(op%rimap)
406 safe_deallocate_a(op%rimap_inv)
408 safe_allocate(op%ri(1:op%stencil%size, 1:op%nri))
409 safe_allocate(op%rimap(1:op%np))
410 safe_allocate(op%rimap_inv(1:op%nri + 1))
411 op%ri = 0
412 op%rimap = 0
413 op%rimap_inv = 0
414 current = 0
415
416 ! the sizes
417 if (mesh%use_curvilinear) then
418 safe_allocate(op%nn(1:op%nri))
419 ! for the moment all the sizes are the same
420 op%nn = op%stencil%size
421 end if
422 end if
423
424 end do
425
426 !the inverse mapping
427 op%rimap_inv(1) = 0
428 do jj = 1, op%np
429 op%rimap_inv(op%rimap(jj) + 1) = jj
430 end do
431 op%rimap_inv(op%nri + 1) = op%np
432
433 safe_deallocate_a(st1)
434 safe_deallocate_a(st1r)
435 safe_deallocate_a(st2)
436
437 if (op%mesh%parallel_in_domains) then
438 !now build the arrays required to apply the nl_operator by parts
439
440 !count points
441 op%inner%nri = 0
442 op%outer%nri = 0
443 do ir = 1, op%nri
444 maxp = op%rimap_inv(ir + 1) + maxval(op%ri(1:op%stencil%size, ir))
445 if (maxp <= np) then
446 !inner point
447 op%inner%nri = op%inner%nri + 1
448 assert(op%inner%nri <= op%nri)
449 else
450 !outer point
451 op%outer%nri = op%outer%nri + 1
452 assert(op%outer%nri <= op%nri)
453 end if
454 end do
455
456 assert(op%inner%nri + op%outer%nri == op%nri)
457
458 if (optional_default(regenerate, .false.)) then
459 safe_deallocate_a(op%inner%imin)
460 safe_deallocate_a(op%inner%imax)
461 safe_deallocate_a(op%inner%ri)
462 safe_deallocate_a(op%outer%imin)
463 safe_deallocate_a(op%outer%imax)
464 safe_deallocate_a(op%outer%ri)
465 end if
466 safe_allocate(op%inner%imin(1:op%inner%nri + 1))
467 safe_allocate(op%inner%imax(1:op%inner%nri))
468 safe_allocate(op%inner%ri(1:op%stencil%size, 1:op%inner%nri))
469
470 safe_allocate(op%outer%imin(1:op%outer%nri + 1))
471 safe_allocate(op%outer%imax(1:op%outer%nri))
472 safe_allocate(op%outer%ri(1:op%stencil%size, 1:op%outer%nri))
473
474 !now populate the arrays
475 iinner = 0
476 iouter = 0
477 do ir = 1, op%nri
478 maxp = op%rimap_inv(ir + 1) + maxval(op%ri(1:op%stencil%size, ir))
479 if (maxp <= np) then
480 !inner point
481 iinner = iinner + 1
482 op%inner%imin(iinner) = op%rimap_inv(ir)
483 op%inner%imax(iinner) = op%rimap_inv(ir + 1)
484 op%inner%ri(1:op%stencil%size, iinner) = op%ri(1:op%stencil%size, ir)
485 else
486 !outer point
487 iouter = iouter + 1
488 op%outer%imin(iouter) = op%rimap_inv(ir)
489 op%outer%imax(iouter) = op%rimap_inv(ir + 1)
490 op%outer%ri(1:op%stencil%size, iouter) = op%ri(1:op%stencil%size, ir)
491 end if
492 end do
493
494 !verify that all points in the inner operator are actually inner
495 do ir = 1, op%inner%nri
496 do ii = op%inner%imin(ir) + 1, op%inner%imax(ir)
497 assert(all(ii + op%inner%ri(1:op%stencil%size, ir) <= mesh%np))
498 end do
499 end do
500
501 end if
502
503 if (accel_is_enabled() .and. op%const_w) then
504
505 write(flags, '(i5)') op%stencil%size
506 flags='-DNDIM=3 -DSTENCIL_SIZE='//trim(adjustl(flags))
507
508 if (op%mesh%parallel_in_domains) flags = '-DINDIRECT '//trim(flags)
509
510 select case (function_opencl)
511 case (op_invmap)
512 call accel_kernel_build(op%kernel, 'operate.cl', 'operate', flags)
513 case (op_map)
514 call accel_kernel_build(op%kernel, 'operate.cl', 'operate_map', flags)
515 end select
516
517 ! conversion to i8 needed to avoid integer overflow
518 call accel_create_buffer(op%buff_ri, accel_mem_read_only, type_integer, int(op%nri, int64)*op%stencil%size)
519 call accel_write_buffer(op%buff_ri, op%stencil%size, op%nri, op%ri)
520
521 select case (function_opencl)
522 case (op_invmap)
523 call accel_create_buffer(op%buff_imin, accel_mem_read_only, type_integer, op%nri)
524 call accel_write_buffer(op%buff_imin, op%nri, op%rimap_inv(1:))
525 call accel_create_buffer(op%buff_imax, accel_mem_read_only, type_integer, op%nri)
526 call accel_write_buffer(op%buff_imax, op%nri, op%rimap_inv(2:))
527
528 case (op_map)
529
531 call accel_write_buffer(op%buff_map, op%mesh%np, (op%rimap - 1)*op%stencil%size)
532
533 if (op%mesh%parallel_in_domains) then
534
535 safe_allocate(inner_points(1:op%mesh%np))
536 safe_allocate(outer_points(1:op%mesh%np))
537 safe_allocate(all_points(1:op%mesh%np))
538
539 op%ninner = 0
540 op%nouter = 0
541
542 do ii = 1, op%mesh%np
543 all_points(ii) = ii - 1
544 maxp = ii + maxval(op%ri(1:op%stencil%size, op%rimap(ii)))
545 if (maxp <= op%mesh%np) then
546 op%ninner = op%ninner + 1
547 inner_points(op%ninner) = ii - 1
548 else
549 op%nouter = op%nouter + 1
550 outer_points(op%nouter) = ii - 1
551 end if
552 end do
553
555 call accel_write_buffer(op%buff_all, op%mesh%np, all_points)
556
558 call accel_write_buffer(op%buff_inner, op%ninner, inner_points)
559
561 call accel_write_buffer(op%buff_outer, op%nouter, outer_points)
562
563 safe_deallocate_a(inner_points)
564 safe_deallocate_a(outer_points)
565 safe_deallocate_a(all_points)
566
567 end if
568 end select
569 end if
570
571 pop_sub(nl_operator_build)
572
573 end subroutine nl_operator_build
574
575 ! ---------------------------------------------------------
576 subroutine nl_operator_output_weights(this)
577 type(nl_operator_t), intent(inout) :: this
578
579 integer :: istencil, idir
580
582
583 write(message(1), '(3a)') 'Debug info: Finite difference weights for ', trim(this%label), '.'
584 write(message(2), '(a)') ' Spacing:'
585 do idir = 1, this%mesh%box%dim
586 write(message(2), '(a,f16.8)') trim(message(2)), this%mesh%spacing(idir)
587 end do
588 call messages_info(2, debug_only=.true.)
589
590 do istencil = 1, this%stencil%size
591 select case(this%mesh%box%dim)
592 case(1)
593 write(message(1), '(a,i3,1i4,f25.10)') ' ', istencil, this%stencil%points(1:1, istencil), this%w(istencil, 1)
594 case(2)
595 write(message(1), '(a,i3,2i4,f25.10)') ' ', istencil, this%stencil%points(1:2, istencil), this%w(istencil, 1)
596 case(3)
597 write(message(1), '(a,i3,3i4,f25.10)') ' ', istencil, this%stencil%points(1:3, istencil), this%w(istencil, 1)
598 end select
599 call messages_info(1, debug_only=.true.)
600 end do
601
603
604 end subroutine nl_operator_output_weights
605
606 ! ---------------------------------------------------------
607 subroutine nl_operator_end(op)
608 type(nl_operator_t), intent(inout) :: op
609
610 push_sub(nl_operator_end)
611
612 if (accel_is_enabled() .and. op%const_w) then
614 end if
615
616 safe_deallocate_a(op%inner%imin)
617 safe_deallocate_a(op%inner%imax)
618 safe_deallocate_a(op%inner%ri)
619 safe_deallocate_a(op%outer%imin)
620 safe_deallocate_a(op%outer%imax)
621 safe_deallocate_a(op%outer%ri)
622
623 safe_deallocate_a(op%w)
624
625 safe_deallocate_a(op%ri)
626 safe_deallocate_a(op%rimap)
627 safe_deallocate_a(op%rimap_inv)
628 safe_deallocate_a(op%nn)
629
630 call stencil_end(op%stencil)
631
632 pop_sub(nl_operator_end)
633 end subroutine nl_operator_end
634
635
636 subroutine nl_operator_clear_gpu_buffers(op)
637 type(nl_operator_t), intent(inout) :: op
638
640
641 call accel_release_buffer(op%buff_ri)
642 select case (function_opencl)
643 case (op_invmap)
644 call accel_release_buffer(op%buff_imin)
645 call accel_release_buffer(op%buff_imax)
646
647 case (op_map)
648 call accel_release_buffer(op%buff_map)
649 if (op%mesh%parallel_in_domains) then
650 call accel_release_buffer(op%buff_all)
651 call accel_release_buffer(op%buff_inner)
652 call accel_release_buffer(op%buff_outer)
653 end if
654
655 case (op_nomap)
656 call accel_release_buffer(op%buff_map)
657 call accel_release_buffer(op%buff_stencil)
658 call accel_release_buffer(op%buff_xyz_to_ip)
659 call accel_release_buffer(op%buff_ip_to_xyz)
660 end select
661
662 call accel_release_buffer(op%buff_weights)
663 call accel_release_buffer(op%buff_half_weights)
664
666 end subroutine nl_operator_clear_gpu_buffers
667
668 ! ---------------------------------------------------------
669 integer pure function nl_operator_get_index(op, is, ip) result(res)
670 type(nl_operator_t), intent(in) :: op
671 integer, intent(in) :: is
672 integer, intent(in) :: ip
673
674 res = ip + op%ri(is, op%rimap(ip))
675 end function nl_operator_get_index
676
677 ! ---------------------------------------------------------
678
680 type(nl_operator_t), intent(inout) :: op
681
683
684 ! Update the GPU weights
685 if (accel_is_enabled() .and. op%const_w) then
686 call accel_create_buffer(op%buff_weights, accel_mem_read_only, type_float, op%stencil%size)
687 call accel_create_buffer(op%buff_half_weights, accel_mem_read_only, type_float, op%stencil%size)
688 end if
689
692
693
694 ! ---------------------------------------------------------
695
696 subroutine nl_operator_update_gpu_buffers(op)
697 type(nl_operator_t), intent(inout) :: op
698
700
701 ! Update the GPU weights
702 if (accel_is_enabled() .and. op%const_w) then
703 call accel_write_buffer(op%buff_weights, op%stencil%size, op%w(:, 1))
704 call accel_write_buffer(op%buff_half_weights, op%stencil%size, -m_half*op%w(:, 1))
705 end if
706
708 end subroutine nl_operator_update_gpu_buffers
709
710 ! ---------------------------------------------------------
711
712 integer pure function nl_operator_np_zero_bc(op) result(np_bc)
713 type(nl_operator_t), intent(in) :: op
714
715 integer :: jj, ii
716
717 np_bc = 0
718 do jj = 1, op%nri
719 ii = op%rimap_inv(jj + 1) + maxval(op%ri(1:op%stencil%size, jj))
720 np_bc = max(np_bc, ii)
721 end do
722
723 end function nl_operator_np_zero_bc
724
725
726 ! ---------------------------------------------------------
728 subroutine nl_operator_remove_zero_weight_points(op, space, mesh)
729 type(nl_operator_t), intent(inout) :: op
730 type(space_t), intent(in) :: space
731 class(mesh_t), intent(in) :: mesh
732
733 integer :: ip, size
734 real(real64), parameter :: tol = 1.0e-14_real64
735 real(real64) :: max_weight, new_w(op%stencil%size)
736 integer :: new_points(space%dim, op%stencil%size)
737
738 if (.not. op%const_w) return
739
741
742 max_weight = maxval(abs(op%w(:, 1)))
743 size = 0
744 do ip = 1, op%stencil%size
745 if (abs(op%w(ip, 1)) > tol * max_weight) then
746 size = size +1
747 new_w(size) = op%w(ip, 1)
748 new_points(:,size) = op%stencil%points(:, ip)
749 end if
750 end do
751
752 ! We regenerate the stencil without the zero-weight points
753 op%stencil%size = size
754 safe_deallocate_a(op%stencil%points)
755 safe_allocate(op%stencil%points(space%dim, op%stencil%size))
756 op%stencil%points(:, :) = new_points(:, 1:size)
757 safe_deallocate_a(op%w)
759 call nl_operator_build(space, mesh, op, mesh%np, const_w=op%const_w, regenerate=.true.)
760 op%w(1:size, 1) = new_w(1:size)
761
765#include "undef.F90"
766#include "real.F90"
767#include "nl_operator_inc.F90"
768
769#include "undef.F90"
770#include "complex.F90"
771#include "nl_operator_inc.F90"
772
773end module nl_operator_oct_m
775!! Local Variables:
776!! mode: f90
777!! coding: utf-8
778!! End:
subroutine, public accel_kernel_build(this, file_name, kernel_name, flags)
Definition: accel.F90:1340
subroutine, public accel_release_buffer(this, async)
Definition: accel.F90:918
pure logical function, public accel_is_enabled()
Definition: accel.F90:418
integer, parameter, public accel_mem_read_only
Definition: accel.F90:195
integer pure function, public accel_max_workgroup_size()
Definition: accel.F90:1098
This module implements batches of mesh functions.
Definition: batch.F90:135
Module implementing boundary conditions in Octopus.
Definition: boundaries.F90:124
real(real64), parameter, public m_zero
Definition: global.F90:190
real(real64), parameter, public m_half
Definition: global.F90:196
This module implements the index, used for the mesh points.
Definition: index.F90:124
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:117
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
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
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1029
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:162
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:697
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1069
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:600
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:147
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 nl_operator_clear_gpu_buffers(op)
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)
subroutine, public nl_operator_remove_zero_weight_points(op, space, mesh)
Removes the zero-weight points for constant weight stencils.
subroutine, public nl_operator_global_end()
integer, parameter, public op_outer
subroutine, public nl_operator_build(space, mesh, op, np, const_w, regenerate)
Creates the nonlocal operators for the stencils used for finite differences.
subroutine, public znl_operator_operate_batch(op, fi, fo, ghost_update, profile, points, factor, async)
integer, parameter op_nomap
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:119
Some general things and nomenclature:
Definition: par_vec.F90:173
This module defines stencils used in Octopus.
Definition: stencil.F90:137
subroutine, public stencil_end(this)
Definition: stencil.F90:216
subroutine, public stencil_copy(input, output)
Definition: stencil.F90:198
type(type_t), public type_float
Definition: types.F90:135
type(type_t), public type_integer
Definition: types.F90:137
Describes mesh distribution to nodes.
Definition: mesh.F90:187
index type for non-local operators
data type for non local operators
int true(void)