Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
vector_math.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!
62 use num_types, only : rp
63 use vector, only : vector_t
64 use mask, only : mask_t
65 use device
66 use utils, only : neko_error
67 use math, only : rzero, rone, copy, cmult, cadd, cfill, invcol1, vdot3, &
85 use, intrinsic :: iso_c_binding, only : c_ptr
86 implicit none
87 private
88
100
101contains
102
104 subroutine vector_rzero(a, n)
105 integer, intent(in), optional :: n
106 type(vector_t), intent(inout) :: a
107 integer :: size
108
109 if (present(n)) then
110 size = n
111 else
112 size = a%size()
113 end if
114
115 if (neko_bcknd_device .eq. 1) then
116 call device_rzero(a%x_d, size)
117 else
118 call rzero(a%x, size)
119 end if
120 end subroutine vector_rzero
121
123 subroutine vector_rone(a, n)
124 integer, intent(in), optional :: n
125 type(vector_t), intent(inout) :: a
126 integer :: size
127
128 if (present(n)) then
129 size = n
130 else
131 size = a%size()
132 end if
133
134 if (neko_bcknd_device .eq. 1) then
135 call device_rone(a%x_d, size)
136 else
137 call rone(a%x, size)
138 end if
139 end subroutine vector_rone
140
142 subroutine vector_copy(a, b, n)
143 integer, intent(in), optional :: n
144 type(vector_t), intent(in) :: b
145 type(vector_t), intent(inout) :: a
146 integer :: size
147
148 if (present(n)) then
149 size = n
150 else
151 size = a%size()
152 end if
153
154 if (neko_bcknd_device .eq. 1) then
155 call device_copy(a%x_d, b%x_d, size)
156 else
157 call copy(a%x, b%x, size)
158 end if
159 end subroutine vector_copy
160
162 subroutine vector_cmult(a, c, n)
163 integer, intent(in), optional :: n
164 type(vector_t), intent(inout) :: a
165 real(kind=rp), intent(in) :: c
166 integer :: size
167
168 if (present(n)) then
169 size = n
170 else
171 size = a%size()
172 end if
173
174 if (neko_bcknd_device .eq. 1) then
175 call device_cmult(a%x_d, c, size)
176 else
177 call cmult(a%x, c, size)
178 end if
179 end subroutine vector_cmult
180
182 subroutine vector_cadd(a, s, n)
183 integer, intent(in), optional :: n
184 type(vector_t), intent(inout) :: a
185 real(kind=rp), intent(in) :: s
186 integer :: size
187
188 if (present(n)) then
189 size = n
190 else
191 size = a%size()
192 end if
193
194 if (neko_bcknd_device .eq. 1) then
195 call device_cadd(a%x_d, s, size)
196 else
197 call cadd(a%x, s, size)
198 end if
199 end subroutine vector_cadd
200
202 subroutine vector_cadd2(a, b, s, n)
203 integer, intent(in), optional :: n
204 type(vector_t), intent(inout) :: a
205 type(vector_t), intent(in) :: b
206 real(kind=rp), intent(in) :: s
207 integer :: size
208
209 if (present(n)) then
210 size = n
211 else
212 size = a%size()
213 end if
214
215 if (neko_bcknd_device .eq. 1) then
216 call device_cadd2(a%x_d, b%x_d, s, size)
217 else
218 call cadd2(a%x, b%x, s, size)
219 end if
220 end subroutine vector_cadd2
221
222
224 subroutine vector_cfill(a, c, n)
225 integer, intent(in), optional :: n
226 type(vector_t), intent(inout) :: a
227 real(kind=rp), intent(in) :: c
228 integer :: size
229
230 if (present(n)) then
231 size = n
232 else
233 size = a%size()
234 end if
235
236 if (neko_bcknd_device .eq. 1) then
237 call device_cfill(a%x_d, c, size)
238 else
239 call cfill(a%x, c, size)
240 end if
241 end subroutine vector_cfill
242
244 subroutine vector_invcol1(a, n)
245 integer, intent(in), optional :: n
246 type(vector_t), intent(inout) :: a
247 integer :: size
248
249 if (present(n)) then
250 size = n
251 else
252 size = a%size()
253 end if
254
255 if (neko_bcknd_device .eq. 1) then
256 call device_invcol1(a%x_d, size)
257 else
258 call invcol1(a%x, size)
259 end if
260
261 end subroutine vector_invcol1
262
264 subroutine vector_invcol3(a, b, c, n)
265 integer, intent(in), optional :: n
266 type(vector_t), intent(inout) :: a
267 type(vector_t), intent(in) :: b
268 type(vector_t), intent(in) :: c
269 integer :: size
270
271 if (present(n)) then
272 size = n
273 else
274 size = a%size()
275 end if
276
277 if (neko_bcknd_device .eq. 1) then
278 call device_invcol3(a%x_d, b%x_d, c%x_d, size)
279 else
280 call invcol3(a%x, b%x, c%x, size)
281 end if
282
283 end subroutine vector_invcol3
284
287 subroutine vector_vdot3(dot, u1, u2, u3, v1, v2, v3, n)
288 integer, intent(in), optional :: n
289 type(vector_t), intent(in) :: u1, u2, u3
290 type(vector_t), intent(in) :: v1, v2, v3
291 type(vector_t), intent(out) :: dot
292 integer :: size
293
294 if (present(n)) then
295 size = n
296 else
297 size = dot%size()
298 end if
299
300 if (neko_bcknd_device .eq. 1) then
301 call device_vdot3(dot%x_d, &
302 u1%x_d, u2%x_d, u3%x_d, &
303 v1%x_d, v2%x_d, v3%x_d, &
304 size)
305 else
306 call vdot3(dot%x, &
307 u1%x, u2%x, u3%x, &
308 v1%x, v2%x, v3%x, &
309 size)
310 end if
311
312 end subroutine vector_vdot3
313
315 subroutine vector_add2(a, b, n)
316 integer, intent(in), optional :: n
317 type(vector_t), intent(inout) :: a
318 type(vector_t), intent(in) :: b
319 integer :: size
320
321 if (present(n)) then
322 size = n
323 else
324 size = a%size()
325 end if
326
327 if (neko_bcknd_device .eq. 1) then
328 call device_add2(a%x_d, b%x_d, size)
329 else
330 call add2(a%x, b%x, size)
331 end if
332
333 end subroutine vector_add2
334
336 subroutine vector_add3(a, b, c, n)
337 integer, intent(in), optional :: n
338 type(vector_t), intent(inout) :: a
339 type(vector_t), intent(in) :: b, c
340 integer :: size
341
342 if (present(n)) then
343 size = n
344 else
345 size = a%size()
346 end if
347
348 if (neko_bcknd_device .eq. 1) then
349 call device_add3(a%x_d, b%x_d, c%x_d, size)
350 else
351 call add3(a%x, b%x, c%x, size)
352 end if
353
354 end subroutine vector_add3
355
357 subroutine vector_add4(a, b, c, d, n)
358 integer, intent(in), optional :: n
359 type(vector_t), intent(inout) :: a
360 type(vector_t), intent(in) :: b, c, d
361 integer :: size
362
363 if (present(n)) then
364 size = n
365 else
366 size = a%size()
367 end if
368
369 if (neko_bcknd_device .eq. 1) then
370 call device_add4(a%x_d, b%x_d, c%x_d, d%x_d, size)
371 else
372 call add4(a%x, b%x, c%x, d%x, size)
373 end if
374
375 end subroutine vector_add4
376
378 subroutine vector_sub2(a, b, n)
379 integer, intent(in), optional :: n
380 type(vector_t), intent(inout) :: a
381 type(vector_t), intent(inout) :: b
382 integer :: size
383
384 if (present(n)) then
385 size = n
386 else
387 size = a%size()
388 end if
389
390 if (neko_bcknd_device .eq. 1) then
391 call device_sub2(a%x_d, b%x_d, size)
392 else
393 call sub2(a%x, b%x, size)
394 end if
395
396 end subroutine vector_sub2
397
399 subroutine vector_sub3(a, b, c, n)
400 integer, intent(in), optional :: n
401 type(vector_t), intent(inout) :: a
402 type(vector_t), intent(in) :: b
403 type(vector_t), intent(in) :: c
404 integer :: size
405
406 if (present(n)) then
407 size = n
408 else
409 size = a%size()
410 end if
411
412 if (neko_bcknd_device .eq. 1) then
413 call device_sub3(a%x_d, b%x_d, c%x_d, size)
414 else
415 call sub3(a%x, b%x, c%x, size)
416 end if
417
418 end subroutine vector_sub3
419
420
423 subroutine vector_add2s1(a, b, c1, n)
424 integer, intent(in), optional :: n
425 type(vector_t), intent(inout) :: a
426 type(vector_t), intent(inout) :: b
427 real(kind=rp), intent(in) :: c1
428 integer :: size
429
430 if (present(n)) then
431 size = n
432 else
433 size = a%size()
434 end if
435
436 if (neko_bcknd_device .eq. 1) then
437 call device_add2s1(a%x_d, b%x_d, c1, size)
438 else
439 call add2s1(a%x, b%x, c1, size)
440 end if
441
442 end subroutine vector_add2s1
443
446 subroutine vector_add2s2(a, b, c1, n)
447 integer, intent(in), optional :: n
448 type(vector_t), intent(inout) :: a
449 type(vector_t), intent(inout) :: b
450 real(kind=rp), intent(in) :: c1
451 integer :: size
452
453 if (present(n)) then
454 size = n
455 else
456 size = a%size()
457 end if
458
459 if (neko_bcknd_device .eq. 1) then
460 call device_add2s2(a%x_d, b%x_d, c1, size)
461 else
462 call add2s2(a%x, b%x, c1, size)
463 end if
464
465 end subroutine vector_add2s2
466
468 subroutine vector_addsqr2s2(a, b, c1, n)
469 integer, intent(in), optional :: n
470 type(vector_t), intent(inout) :: a
471 type(vector_t), intent(in) :: b
472 real(kind=rp), intent(in) :: c1
473 integer :: size
474
475 if (present(n)) then
476 size = n
477 else
478 size = a%size()
479 end if
480
481 if (neko_bcknd_device .eq. 1) then
482 call device_addsqr2s2(a%x_d, b%x_d, c1, size)
483 else
484 call addsqr2s2(a%x, b%x, c1, size)
485 end if
486
487 end subroutine vector_addsqr2s2
488
490 subroutine vector_cmult2(a, b, c, n)
491 integer, intent(in), optional :: n
492 type(vector_t), intent(inout) :: a
493 type(vector_t), intent(in) :: b
494 real(kind=rp), intent(in) :: c
495 integer :: size
496
497 if (present(n)) then
498 size = n
499 else
500 size = a%size()
501 end if
502
503 if (neko_bcknd_device .eq. 1) then
504 call device_cmult2(a%x_d, b%x_d, c, size)
505 else
506 call cmult2(a%x, b%x, c, size)
507 end if
508
509 end subroutine vector_cmult2
510
512 subroutine vector_invcol2(a, b, n)
513 integer, intent(in), optional :: n
514 type(vector_t), intent(inout) :: a
515 type(vector_t), intent(in) :: b
516 integer :: size
517
518 if (present(n)) then
519 size = n
520 else
521 size = a%size()
522 end if
523
524 if (neko_bcknd_device .eq. 1) then
525 call device_invcol2(a%x_d, b%x_d, size)
526 else
527 call invcol2(a%x, b%x, size)
528 end if
529
530 end subroutine vector_invcol2
531
532
534 subroutine vector_col2(a, b, n)
535 integer, intent(in), optional :: n
536 type(vector_t), intent(inout) :: a
537 type(vector_t), intent(in) :: b
538 integer :: size
539
540 if (present(n)) then
541 size = n
542 else
543 size = a%size()
544 end if
545
546 if (neko_bcknd_device .eq. 1) then
547 call device_col2(a%x_d, b%x_d, size)
548 else
549 call col2(a%x, b%x, size)
550 end if
551
552 end subroutine vector_col2
553
555 subroutine vector_col3(a, b, c, n)
556 integer, intent(in), optional :: n
557 type(vector_t), intent(inout) :: a
558 type(vector_t), intent(in) :: b
559 type(vector_t), intent(in) :: c
560 integer :: size
561
562 if (present(n)) then
563 size = n
564 else
565 size = a%size()
566 end if
567
568 if (neko_bcknd_device .eq. 1) then
569 call device_col3(a%x_d, b%x_d, c%x_d, size)
570 else
571 call col3(a%x, b%x, c%x, size)
572 end if
573
574 end subroutine vector_col3
575
577 subroutine vector_subcol3(a, b, c, n)
578 integer, intent(in), optional :: n
579 type(vector_t), intent(inout) :: a
580 type(vector_t), intent(in) :: b
581 type(vector_t), intent(in) :: c
582 integer :: size
583
584 if (present(n)) then
585 size = n
586 else
587 size = a%size()
588 end if
589
590 if (neko_bcknd_device .eq. 1) then
591 call device_subcol3(a%x_d, b%x_d, c%x_d, size)
592 else
593 call subcol3(a%x, b%x, c%x, size)
594 end if
595
596 end subroutine vector_subcol3
597
599 subroutine vector_add3s2(a, b, c, c1, c2, n)
600 integer, intent(in), optional :: n
601 type(vector_t), intent(inout) :: a
602 type(vector_t), intent(in) :: b
603 type(vector_t), intent(in) :: c
604 real(kind=rp), intent(in) :: c1, c2
605 integer :: size
606
607 if (present(n)) then
608 size = n
609 else
610 size = a%size()
611 end if
612
613 if (neko_bcknd_device .eq. 1) then
614 call device_add3s2(a%x_d, b%x_d, c%x_d, c1, c2, size)
615 else
616 call add3s2(a%x, b%x, c%x, c1, c2, size)
617 end if
618
619 end subroutine vector_add3s2
620
622 subroutine vector_addcol3(a, b, c, n)
623 integer, intent(in), optional :: n
624 type(vector_t), intent(inout) :: a
625 type(vector_t), intent(in) :: b
626 type(vector_t), intent(in) :: c
627 integer :: size
628
629 if (present(n)) then
630 size = n
631 else
632 size = a%size()
633 end if
634
635 if (neko_bcknd_device .eq. 1) then
636 call device_addcol3(a%x_d, b%x_d, c%x_d, size)
637 else
638 call addcol3(a%x, b%x, c%x, size)
639 end if
640
641 end subroutine vector_addcol3
642
644 subroutine vector_addcol4(a, b, c, d, n)
645 integer, intent(in), optional :: n
646 type(vector_t), intent(inout) :: a
647 type(vector_t), intent(in) :: b
648 type(vector_t), intent(in) :: c
649 type(vector_t), intent(in) :: d
650 integer :: size
651
652 if (present(n)) then
653 size = n
654 else
655 size = a%size()
656 end if
657
658 if (neko_bcknd_device .eq. 1) then
659 call device_addcol4(a%x_d, b%x_d, c%x_d, d%x_d, size)
660 else
661 call addcol4(a%x, b%x, c%x, d%x, size)
662 end if
663
664 end subroutine vector_addcol4
665
666 function vector_glsum(a, n) result(sum)
667 integer, intent(in), optional :: n
668 type(vector_t), intent(in) :: a
669 real(kind=rp) :: sum
670 integer :: size
671
672 if (present(n)) then
673 size = n
674 else
675 size = a%size()
676 end if
677
678 if (neko_bcknd_device .eq. 1) then
679 sum = device_glsum(a%x_d, size)
680 else
681 sum = glsum(a%x, size)
682 end if
683
684 end function vector_glsum
685
687 function vector_glmax(a, n) result(val)
688 integer, intent(in), optional :: n
689 type(vector_t), intent(in) :: a
690 real(kind=rp) :: val
691 integer :: size
692
693 if (present(n)) then
694 size = n
695 else
696 size = a%size()
697 end if
698
699 if (neko_bcknd_device .eq. 1) then
700 val = device_glmax(a%x_d, size)
701 else
702 val = glmax(a%x, size)
703 end if
704
705 end function vector_glmax
706
708 function vector_glmin(a, n) result(val)
709 integer, intent(in), optional :: n
710 type(vector_t), intent(in) :: a
711 real(kind=rp) :: val
712 integer :: size
713
714 if (present(n)) then
715 size = n
716 else
717 size = a%size()
718 end if
719
720 if (neko_bcknd_device .eq. 1) then
721 val = device_glmin(a%x_d, size)
722 else
723 val = glmin(a%x, size)
724 end if
725
726 end function vector_glmin
727
728 function vector_glsc2(a, b, n) result(norm)
729 integer, intent(in), optional :: n
730 type(vector_t), intent(in) :: a, b
731 real(kind=rp) :: norm
732 integer :: size
733
734 if (present(n)) then
735 size = n
736 else
737 size = a%size()
738 end if
739
740 if (neko_bcknd_device .eq. 1) then
741 norm = device_glsc2(a%x_d, b%x_d, size)
742 else
743 norm = glsc2(a%x, b%x, size)
744 end if
745
746 end function vector_glsc2
747
748 function vector_glsc3(a, b, c, n) result(norm)
749 integer, intent(in), optional :: n
750 type(vector_t), intent(in) :: a, b, c
751 real(kind=rp) :: norm
752 integer :: size
753
754 if (present(n)) then
755 size = n
756 else
757 size = a%size()
758 end if
759
760 if (neko_bcknd_device .eq. 1) then
761 norm = device_glsc3(a%x_d, b%x_d, c%x_d, size)
762 else
763 norm = glsc3(a%x, b%x, c%x, size)
764 end if
765
766 end function vector_glsc3
767
768 function vector_glsubnorm(a, b, n) result(norm)
769 integer, intent(in), optional :: n
770 type(vector_t), intent(in) :: a, b
771 real(kind=rp) :: norm
772 integer :: size
773
774 if (present(n)) then
775 size = n
776 else
777 size = a%size()
778 end if
779
780 if (neko_bcknd_device .eq. 1) then
781 norm = device_glsubnorm(a%x_d, b%x_d, size)
782 else
783 norm = glsubnorm(a%x, b%x, size)
784 end if
785
786 end function vector_glsubnorm
787
796 subroutine vector_masked_gather_copy_0(a, b, mask, n, n_mask)
797 integer, intent(in) :: n, n_mask
798 type(vector_t), intent(inout) :: a
799 real(kind=rp), dimension(n), intent(in) :: b
800 integer, dimension(0:n_mask) :: mask
801 type(c_ptr) :: mask_d, b_d
802
803 if (n .lt. 1 .or. n_mask .lt. 1) return !Avoid getting null pointers
804
805 if (neko_bcknd_device .eq. 1) then
806 mask_d = device_get_ptr(mask)
807 b_d = device_get_ptr(b)
808 call device_masked_gather_copy_0(a%x_d, b_d, mask_d, n, n_mask)
809 else
810 call masked_gather_copy_0(a%x, b, mask, n, n_mask)
811 end if
812
813 end subroutine vector_masked_gather_copy_0
814
824 subroutine vector_face_masked_gather_copy_0(a, b, mask, facet, lx, ly, lz, &
825 n_mask)
826 integer, intent(in) :: lx, ly, lz, n_mask
827 type(vector_t), intent(inout) :: a
828 real(kind=rp), dimension(:, :, :, :), intent(in) :: b
829 integer, dimension(0:n_mask), intent(in) :: mask
830 integer, dimension(0:n_mask), intent(in) :: facet
831 type(c_ptr) :: mask_d, facet_d, b_d
832
833 if (neko_bcknd_device .eq. 1) then
834 mask_d = device_get_ptr(mask)
835 facet_d = device_get_ptr(facet)
836 b_d = device_get_ptr(b)
837 call device_face_masked_gather_copy_0(a%x_d, b_d, mask_d, facet_d, &
838 size(b, 1), size(b, 2), lx, ly, lz, n_mask)
839 else
840 call face_masked_gather_copy_0(a%x, b, mask, facet, lx, ly, lz, n_mask)
841 end if
842
844
851 subroutine vector_masked_gather_copy(a, b, mask, n)
852 type(vector_t), intent(inout) :: a
853 real(kind=rp), dimension(:), intent(in) :: b
854 type(mask_t), intent(in) :: mask
855 integer, intent(in) :: n
856 type(c_ptr) :: mask_d, b_d
857
858 if (n .lt. 1 .or. mask%size() .lt. 1) return !Avoid getting null pointers
859
860 if (neko_bcknd_device .eq. 1) then
861 mask_d = mask%get_d()
862 b_d = device_get_ptr(b)
863 call device_masked_gather_copy_aligned(a%x_d, b_d, mask_d, n, &
864 mask%size())
865 else
866 call masked_gather_copy(a%x, b, mask%get(), n, mask%size())
867 end if
868
869 end subroutine vector_masked_gather_copy
870
879 subroutine vector_masked_scatter_copy_0(a, b, mask, n, n_mask)
880 integer, intent(in) :: n, n_mask
881 real(kind=rp), dimension(n), intent(inout) :: a
882 type(vector_t), intent(in) :: b
883 integer, dimension(0:n_mask) :: mask
884 type(c_ptr) :: mask_d, a_d
885
886 if (n .lt. 1 .or. n_mask .lt. 1) return !Avoid getting null pointers
887
888 if (neko_bcknd_device .eq. 1) then
889 a_d = device_get_ptr(a)
890 mask_d = device_get_ptr(mask)
891 call device_masked_scatter_copy_0(a_d, b%x_d, mask_d, n, n_mask)
892 else
893 call masked_scatter_copy_0(a, b%x, mask, n, n_mask)
894 end if
895
896 end subroutine vector_masked_scatter_copy_0
897
904 subroutine vector_masked_scatter_copy(a, b, mask, n)
905 real(kind=rp), dimension(:), intent(inout) :: a
906 type(vector_t), intent(in) :: b
907 type(mask_t), intent(in) :: mask
908 integer, intent(in) :: n
909 type(c_ptr) :: mask_d, a_d
910
911 if (n .lt. 1 .or. mask%size() .lt. 1) return !Avoid getting null pointers
912
913 if (neko_bcknd_device .eq. 1) then
914 a_d = device_get_ptr(a)
915 mask_d = mask%get_d()
916 call device_masked_scatter_copy_aligned(a_d, b%x_d, mask_d, n, &
917 mask%size())
918 else
919 call masked_scatter_copy(a, b%x, mask%get(), n, mask%size())
920 end if
921
922 end subroutine vector_masked_scatter_copy
923
925 subroutine vector_cwrap(a, min_value, max_value, n)
926 integer, intent(in), optional :: n
927 type(vector_t), intent(inout) :: a
928 real(kind=rp), intent(in) :: min_value, max_value
929 integer :: size
930
931 if (present(n)) then
932 size = n
933 else
934 size = a%size()
935 end if
936
937 if (a%size() .lt. 1) return ! Avoid getting null pointers
938
939 if (neko_bcknd_device .eq. 1) then
940 call device_cwrap(a%x_d, min_value, max_value, size)
941 else
942 call cwrap(a%x, min_value, max_value, size)
943 end if
944
945 end subroutine vector_cwrap
946
947
948
949end module vector_math
Return the device pointer for an associated Fortran array.
Definition device.F90:107
subroutine, public device_masked_scatter_copy_aligned(a_d, b_d, mask_d, n, n_mask, strm)
Scatter a masked vector .
subroutine, public device_add2s1(a_d, b_d, c1, n, strm)
subroutine, public device_sub3(a_d, b_d, c_d, n, strm)
Vector subtraction .
subroutine, public device_masked_scatter_copy_0(a_d, b_d, mask_d, n, n_mask, strm)
Scatter a masked vector .
subroutine, public device_add2s2(a_d, b_d, c1, n, strm)
Vector addition with scalar multiplication (multiplication on first argument)
real(kind=rp) function, public device_glmax(a_d, n, strm)
Max of a vector of length n.
subroutine, public device_add2(a_d, b_d, n, strm)
Vector addition .
subroutine, public device_addcol3(a_d, b_d, c_d, n, strm)
Returns .
real(kind=rp) function, public device_glsum(a_d, n, strm)
Sum a vector of length n.
subroutine, public device_invcol1(a_d, n, strm)
Invert a vector .
subroutine, public device_add3s2(a_d, b_d, c_d, c1, c2, n, strm)
Returns .
subroutine, public device_rzero(a_d, n, strm)
Zero a real vector.
subroutine, public device_rone(a_d, n, strm)
Set all elements to one.
subroutine, public device_cmult(a_d, c, n, strm)
Multiplication by constant c .
subroutine, public device_vdot3(dot_d, u1_d, u2_d, u3_d, v1_d, v2_d, v3_d, n, strm)
Compute a dot product (3-d version) assuming vector components etc.
real(kind=rp) function, public device_glsubnorm(a_d, b_d, n, strm)
Returns the norm of the difference of two vectors .
subroutine, public device_cadd2(a_d, b_d, c, n, strm)
Add a scalar to vector .
subroutine, public device_sub2(a_d, b_d, n, strm)
Vector substraction .
subroutine, public device_copy(a_d, b_d, n, strm)
Copy a vector .
subroutine, public device_invcol3(a_d, b_d, c_d, n, strm)
Vector division .
subroutine, public device_col2(a_d, b_d, n, strm)
Vector multiplication .
subroutine, public device_masked_gather_copy_aligned(a_d, b_d, mask_d, n, n_mask, strm)
Gather a masked vector .
subroutine, public device_add4(a_d, b_d, c_d, d_d, n, strm)
subroutine, public device_subcol3(a_d, b_d, c_d, n, strm)
Returns .
subroutine, public device_masked_gather_copy_0(a_d, b_d, mask_d, n, n_mask, strm)
Gather a masked vector .
subroutine, public device_invcol2(a_d, b_d, n, strm)
Vector division .
subroutine, public device_addsqr2s2(a_d, b_d, c1, n, strm)
Returns .
real(kind=rp) function, public device_glsc3(a_d, b_d, c_d, n, strm)
Weighted inner product .
subroutine, public device_cwrap(a_d, min_val, max_val, n, strm)
Wrap value around a range (min, max)
subroutine, public device_face_masked_gather_copy_0(a_d, b_d, mask_d, facet_d, n1, n2, lx, ly, lz, n_mask, strm)
Gather a face-local SEM field .
real(kind=rp) function, public device_glsc2(a_d, b_d, n, strm)
Weighted inner product .
subroutine, public device_cmult2(a_d, b_d, c, n, strm)
Multiplication by constant c .
subroutine, public device_col3(a_d, b_d, c_d, n, strm)
Vector multiplication with 3 vectors .
subroutine, public device_addcol4(a_d, b_d, c_d, d_d, n, strm)
Returns .
subroutine, public device_cfill(a_d, c, n, strm)
Set all elements to a constant c .
subroutine, public device_add3(a_d, b_d, c_d, n, strm)
Vector addition .
real(kind=rp) function, public device_glmin(a_d, n, strm)
Min of a vector of length n.
Device abstraction, common interface for various accelerators.
Definition device.F90:34
Object for handling masks in Neko.
Definition mask.f90:34
Definition math.f90:60
subroutine, public cmult(a, c, n)
Multiplication by constant c .
Definition math.f90:447
subroutine, public cmult2(a, b, c, n)
Multiplication by constant c .
Definition math.f90:459
subroutine, public invcol2(a, b, n)
Vector division .
Definition math.f90:890
real(kind=rp) function, public glsc3(a, b, c, n)
Weighted inner product .
Definition math.f90:1117
subroutine, public masked_scatter_copy(a, b, mask, n, n_mask)
Scatter a contigous vector to masked positions in a target array .
Definition math.f90:417
subroutine, public cadd2(a, b, s, n)
Add a scalar to vector .
Definition math.f90:510
subroutine, public face_masked_gather_copy_0(a, b, mask, facet, lx, ly, lz, n_mask)
Gather values from a face-local SEM field to a reduced contiguous vector.
Definition math.f90:338
subroutine, public cadd(a, s, n)
Add a scalar to vector .
Definition math.f90:498
subroutine, public addsqr2s2(a, b, c1, n)
Returns .
Definition math.f90:876
subroutine, public cwrap(a, min_val, max_val, n)
Wrap value around a range [min, max)
Definition math.f90:535
subroutine, public add2s1(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on first argument)
Definition math.f90:847
real(kind=rp) function, public glsc2(a, b, n)
Weighted inner product .
Definition math.f90:1098
subroutine, public masked_scatter_copy_0(a, b, mask, n, n_mask)
Scatter a contigous vector to masked positions in a target array .
Definition math.f90:394
subroutine, public subcol3(a, b, c, n)
Returns .
Definition math.f90:931
subroutine, public rone(a, n)
Set all elements to one.
Definition math.f90:241
subroutine, public masked_gather_copy(a, b, mask, n, n_mask)
Gather a masked vector to reduced contigous vector .
Definition math.f90:371
subroutine, public add3(a, b, c, n)
Vector addition .
Definition math.f90:789
real(kind=rp) function, public glsum(a, n)
Sum a vector of length n.
Definition math.f90:549
subroutine, public sub3(a, b, c, n)
Vector subtraction .
Definition math.f90:831
subroutine, public addcol4(a, b, c, d, n)
Returns .
Definition math.f90:1022
subroutine, public add2(a, b, n)
Vector addition .
Definition math.f90:776
subroutine, public cfill(a, c, n)
Set all elements to a constant c .
Definition math.f90:523
subroutine, public invcol3(a, b, c, n)
Invert a vector .
Definition math.f90:676
subroutine, public add3s2(a, b, c, c1, c2, n)
Returns .
Definition math.f90:945
subroutine, public masked_gather_copy_0(a, b, mask, n, n_mask)
Gather a masked vector to reduced contigous vector .
Definition math.f90:315
subroutine, public addcol3(a, b, c, n)
Returns .
Definition math.f90:1008
subroutine, public invcol1(a, n)
Invert a vector .
Definition math.f90:664
subroutine, public col2(a, b, n)
Vector multiplication .
Definition math.f90:904
real(kind=rp) function, public glmax(a, n)
Max of a vector of length n.
Definition math.f90:566
subroutine, public copy(a, b, n)
Copy a vector .
Definition math.f90:252
subroutine, public add4(a, b, c, d, n)
Vector addition .
Definition math.f90:803
subroutine, public col3(a, b, c, n)
Vector multiplication with 3 vectors .
Definition math.f90:917
subroutine, public vdot3(dot, u1, u2, u3, v1, v2, v3, n)
Compute a dot product (3-d version) assuming vector components etc.
Definition math.f90:734
subroutine, public rzero(a, n)
Zero a real vector.
Definition math.f90:208
real(kind=rp) function, public glsubnorm(a, b, n)
Returns the norm of the difference of two vectors .
Definition math.f90:1159
subroutine, public sub2(a, b, n)
Vector substraction .
Definition math.f90:818
subroutine, public add2s2(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on second argument)
Definition math.f90:862
real(kind=rp) function, public glmin(a, n)
Min of a vector of length n.
Definition math.f90:596
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Utilities.
Definition utils.f90:35
subroutine, public vector_masked_gather_copy_0(a, b, mask, n, n_mask)
Gather a vector to reduced contigous array .
real(kind=rp) function, public vector_glmin(a, n)
Global minimum of all elements in a vector .
real(kind=rp) function, public vector_glsc3(a, b, c, n)
subroutine, public vector_add2(a, b, n)
Vector addition .
subroutine, public vector_vdot3(dot, u1, u2, u3, v1, v2, v3, n)
Compute a dot product (3-d version) assuming vector components etc.
subroutine, public vector_copy(a, b, n)
Copy a vector .
subroutine, public vector_sub3(a, b, c, n)
Vector subtraction .
subroutine, public vector_cmult(a, c, n)
Multiplication by constant c .
subroutine, public vector_invcol1(a, n)
Invert a vector .
subroutine, public vector_add3s2(a, b, c, c1, c2, n)
Returns .
subroutine, public vector_masked_gather_copy(a, b, mask, n)
Gather a vector to reduced contigous array .
subroutine, public vector_sub2(a, b, n)
Vector substraction .
subroutine, public vector_face_masked_gather_copy_0(a, b, mask, facet, lx, ly, lz, n_mask)
Gather a face-local SEM field to a reduced contiguous vector.
subroutine, public vector_addsqr2s2(a, b, c1, n)
Returns .
subroutine, public vector_col2(a, b, n)
Vector multiplication .
subroutine, public vector_cadd2(a, b, s, n)
Add a scalar to vector .
subroutine, public vector_cfill(a, c, n)
Set all elements to a constant c .
real(kind=rp) function, public vector_glsum(a, n)
subroutine, public vector_cwrap(a, min_value, max_value, n)
Wrap vector elements into the range [min_value, max_value)
real(kind=rp) function, public vector_glsc2(a, b, n)
subroutine, public vector_addcol3(a, b, c, n)
Returns .
real(kind=rp) function, public vector_glsubnorm(a, b, n)
subroutine, public vector_invcol2(a, b, n)
Vector division .
subroutine, public vector_subcol3(a, b, c, n)
Returns .
subroutine, public vector_cadd(a, s, n)
Add a scalar to vector .
subroutine vector_add4(a, b, c, d, n)
Vector addition .
subroutine, public vector_add3(a, b, c, n)
Vector addition .
real(kind=rp) function, public vector_glmax(a, n)
Global maximum of all elements in a vector .
subroutine, public vector_add2s1(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on first argument)
subroutine, public vector_add2s2(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on second argument)
subroutine, public vector_addcol4(a, b, c, d, n)
Returns .
subroutine, public vector_masked_scatter_copy(a, b, mask, n)
Scatter a contiguous vector into an array .
subroutine, public vector_cmult2(a, b, c, n)
Multiplication by constant c .
subroutine, public vector_rone(a, n)
Set all elements to one.
subroutine, public vector_rzero(a, n)
Zero a real vector.
subroutine, public vector_col3(a, b, c, n)
Vector multiplication with 3 vectors .
subroutine, public vector_masked_scatter_copy_0(a, b, mask, n, n_mask)
Scatter a contiguous vector into an array .
subroutine, public vector_invcol3(a, b, c, n)
Invert a vector .
Defines a vector.
Definition vector.f90:34
Type for consistently handling masks in Neko. This type encapsulates the mask array and its associate...
Definition mask.f90:51