Neko 1.99.2
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 device
65 use math, only : rzero, rone, copy, cmult, cadd, cfill, invcol1, vdot3, &
77 use, intrinsic :: iso_c_binding, only : c_ptr
78 implicit none
79 private
80
90
91contains
92
94 subroutine vector_rzero(a, n)
95 integer, intent(in), optional :: n
96 type(vector_t), intent(inout) :: a
97 integer :: size
98
99 if (present(n)) then
100 size = n
101 else
102 size = a%size()
103 end if
104
105 if (neko_bcknd_device .eq. 1) then
106 call device_rzero(a%x_d, size)
107 else
108 call rzero(a%x, size)
109 end if
110 end subroutine vector_rzero
111
113 subroutine vector_rone(a, n)
114 integer, intent(in), optional :: n
115 type(vector_t), intent(inout) :: a
116 integer :: size
117
118 if (present(n)) then
119 size = n
120 else
121 size = a%size()
122 end if
123
124 if (neko_bcknd_device .eq. 1) then
125 call device_rone(a%x_d, size)
126 else
127 call rone(a%x, size)
128 end if
129 end subroutine vector_rone
130
132 subroutine vector_copy(a, b, n)
133 integer, intent(in), optional :: n
134 type(vector_t), intent(in) :: b
135 type(vector_t), intent(inout) :: a
136 integer :: size
137
138 if (present(n)) then
139 size = n
140 else
141 size = a%size()
142 end if
143
144 if (neko_bcknd_device .eq. 1) then
145 call device_copy(a%x_d, b%x_d, size)
146 else
147 call copy(a%x, b%x, size)
148 end if
149 end subroutine vector_copy
150
152 subroutine vector_cmult(a, c, n)
153 integer, intent(in), optional :: n
154 type(vector_t), intent(inout) :: a
155 real(kind=rp), intent(in) :: c
156 integer :: size
157
158 if (present(n)) then
159 size = n
160 else
161 size = a%size()
162 end if
163
164 if (neko_bcknd_device .eq. 1) then
165 call device_cmult(a%x_d, c, size)
166 else
167 call cmult(a%x, c, size)
168 end if
169 end subroutine vector_cmult
170
172 subroutine vector_cadd(a, s, n)
173 integer, intent(in), optional :: n
174 type(vector_t), intent(inout) :: a
175 real(kind=rp), intent(in) :: s
176 integer :: size
177
178 if (present(n)) then
179 size = n
180 else
181 size = a%size()
182 end if
183
184 if (neko_bcknd_device .eq. 1) then
185 call device_cadd(a%x_d, s, size)
186 else
187 call cadd(a%x, s, size)
188 end if
189 end subroutine vector_cadd
190
192 subroutine vector_cfill(a, c, n)
193 integer, intent(in), optional :: n
194 type(vector_t), intent(inout) :: a
195 real(kind=rp), intent(in) :: c
196 integer :: size
197
198 if (present(n)) then
199 size = n
200 else
201 size = a%size()
202 end if
203
204 if (neko_bcknd_device .eq. 1) then
205 call device_cfill(a%x_d, c, size)
206 else
207 call cfill(a%x, c, size)
208 end if
209 end subroutine vector_cfill
210
212 subroutine vector_invcol1(a, n)
213 integer, intent(in), optional :: n
214 type(vector_t), intent(inout) :: a
215 integer :: size
216
217 if (present(n)) then
218 size = n
219 else
220 size = a%size()
221 end if
222
223 if (neko_bcknd_device .eq. 1) then
224 call device_invcol1(a%x_d, size)
225 else
226 call invcol1(a%x, size)
227 end if
228
229 end subroutine vector_invcol1
230
232 subroutine vector_invcol3(a, b, c, n)
233 integer, intent(in), optional :: n
234 type(vector_t), intent(inout) :: a
235 type(vector_t), intent(in) :: b
236 type(vector_t), intent(in) :: c
237 integer :: size
238
239 if (present(n)) then
240 size = n
241 else
242 size = a%size()
243 end if
244
245 if (neko_bcknd_device .eq. 1) then
246 call device_invcol3(a%x_d, b%x_d, c%x_d, size)
247 else
248 call invcol3(a%x, b%x, c%x, size)
249 end if
250
251 end subroutine vector_invcol3
252
255 subroutine vector_vdot3(dot, u1, u2, u3, v1, v2, v3, n)
256 integer, intent(in), optional :: n
257 type(vector_t), intent(in) :: u1, u2, u3
258 type(vector_t), intent(in) :: v1, v2, v3
259 type(vector_t), intent(out) :: dot
260 integer :: size
261
262 if (present(n)) then
263 size = n
264 else
265 size = dot%size()
266 end if
267
268 if (neko_bcknd_device .eq. 1) then
269 call device_vdot3(dot%x_d, &
270 u1%x_d, u2%x_d, u3%x_d, &
271 v1%x_d, v2%x_d, v3%x_d, &
272 size)
273 else
274 call vdot3(dot%x, &
275 u1%x, u2%x, u3%x, &
276 v1%x, v2%x, v3%x, &
277 size)
278 end if
279
280 end subroutine vector_vdot3
281
283 subroutine vector_add2(a, b, n)
284 integer, intent(in), optional :: n
285 type(vector_t), intent(inout) :: a
286 type(vector_t), intent(in) :: b
287 integer :: size
288
289 if (present(n)) then
290 size = n
291 else
292 size = a%size()
293 end if
294
295 if (neko_bcknd_device .eq. 1) then
296 call device_add2(a%x_d, b%x_d, size)
297 else
298 call add2(a%x, b%x, size)
299 end if
300
301 end subroutine vector_add2
302
304 subroutine vector_add3(a, b, c, n)
305 integer, intent(in), optional :: n
306 type(vector_t), intent(inout) :: a
307 type(vector_t), intent(in) :: b, c
308 integer :: size
309
310 if (present(n)) then
311 size = n
312 else
313 size = a%size()
314 end if
315
316 if (neko_bcknd_device .eq. 1) then
317 call device_add3(a%x_d, b%x_d, c%x_d, size)
318 else
319 call add3(a%x, b%x, c%x, size)
320 end if
321
322 end subroutine vector_add3
323
325 subroutine vector_add4(a, b, c, d, n)
326 integer, intent(in), optional :: n
327 type(vector_t), intent(inout) :: a
328 type(vector_t), intent(in) :: b, c, d
329 integer :: size
330
331 if (present(n)) then
332 size = n
333 else
334 size = a%size()
335 end if
336
337 if (neko_bcknd_device .eq. 1) then
338 call device_add4(a%x_d, b%x_d, c%x_d, d%x_d, size)
339 else
340 call add4(a%x, b%x, c%x, d%x, size)
341 end if
342
343 end subroutine vector_add4
344
346 subroutine vector_sub2(a, b, n)
347 integer, intent(in), optional :: n
348 type(vector_t), intent(inout) :: a
349 type(vector_t), intent(inout) :: b
350 integer :: size
351
352 if (present(n)) then
353 size = n
354 else
355 size = a%size()
356 end if
357
358 if (neko_bcknd_device .eq. 1) then
359 call device_sub2(a%x_d, b%x_d, size)
360 else
361 call sub2(a%x, b%x, size)
362 end if
363
364 end subroutine vector_sub2
365
367 subroutine vector_sub3(a, b, c, n)
368 integer, intent(in), optional :: n
369 type(vector_t), intent(inout) :: a
370 type(vector_t), intent(in) :: b
371 type(vector_t), intent(in) :: c
372 integer :: size
373
374 if (present(n)) then
375 size = n
376 else
377 size = a%size()
378 end if
379
380 if (neko_bcknd_device .eq. 1) then
381 call device_sub3(a%x_d, b%x_d, c%x_d, size)
382 else
383 call sub3(a%x, b%x, c%x, size)
384 end if
385
386 end subroutine vector_sub3
387
388
391 subroutine vector_add2s1(a, b, c1, n)
392 integer, intent(in), optional :: n
393 type(vector_t), intent(inout) :: a
394 type(vector_t), intent(inout) :: b
395 real(kind=rp), intent(in) :: c1
396 integer :: size
397
398 if (present(n)) then
399 size = n
400 else
401 size = a%size()
402 end if
403
404 if (neko_bcknd_device .eq. 1) then
405 call device_add2s1(a%x_d, b%x_d, c1, size)
406 else
407 call add2s1(a%x, b%x, c1, size)
408 end if
409
410 end subroutine vector_add2s1
411
414 subroutine vector_add2s2(a, b, c1, n)
415 integer, intent(in), optional :: n
416 type(vector_t), intent(inout) :: a
417 type(vector_t), intent(inout) :: b
418 real(kind=rp), intent(in) :: c1
419 integer :: size
420
421 if (present(n)) then
422 size = n
423 else
424 size = a%size()
425 end if
426
427 if (neko_bcknd_device .eq. 1) then
428 call device_add2s2(a%x_d, b%x_d, c1, size)
429 else
430 call add2s2(a%x, b%x, c1, size)
431 end if
432
433 end subroutine vector_add2s2
434
436 subroutine vector_addsqr2s2(a, b, c1, n)
437 integer, intent(in), optional :: n
438 type(vector_t), intent(inout) :: a
439 type(vector_t), intent(in) :: b
440 real(kind=rp), intent(in) :: c1
441 integer :: size
442
443 if (present(n)) then
444 size = n
445 else
446 size = a%size()
447 end if
448
449 if (neko_bcknd_device .eq. 1) then
450 call device_addsqr2s2(a%x_d, b%x_d, c1, size)
451 else
452 call addsqr2s2(a%x, b%x, c1, size)
453 end if
454
455 end subroutine vector_addsqr2s2
456
458 subroutine vector_cmult2(a, b, c, n)
459 integer, intent(in), optional :: n
460 type(vector_t), intent(inout) :: a
461 type(vector_t), intent(in) :: b
462 real(kind=rp), intent(in) :: c
463 integer :: size
464
465 if (present(n)) then
466 size = n
467 else
468 size = a%size()
469 end if
470
471 if (neko_bcknd_device .eq. 1) then
472 call device_cmult2(a%x_d, b%x_d, c, size)
473 else
474 call cmult2(a%x, b%x, c, size)
475 end if
476
477 end subroutine vector_cmult2
478
480 subroutine vector_invcol2(a, b, n)
481 integer, intent(in), optional :: n
482 type(vector_t), intent(inout) :: a
483 type(vector_t), intent(in) :: b
484 integer :: size
485
486 if (present(n)) then
487 size = n
488 else
489 size = a%size()
490 end if
491
492 if (neko_bcknd_device .eq. 1) then
493 call device_invcol2(a%x_d, b%x_d, size)
494 else
495 call invcol2(a%x, b%x, size)
496 end if
497
498 end subroutine vector_invcol2
499
500
502 subroutine vector_col2(a, b, n)
503 integer, intent(in), optional :: n
504 type(vector_t), intent(inout) :: a
505 type(vector_t), intent(in) :: b
506 integer :: size
507
508 if (present(n)) then
509 size = n
510 else
511 size = a%size()
512 end if
513
514 if (neko_bcknd_device .eq. 1) then
515 call device_col2(a%x_d, b%x_d, size)
516 else
517 call col2(a%x, b%x, size)
518 end if
519
520 end subroutine vector_col2
521
523 subroutine vector_col3(a, b, c, n)
524 integer, intent(in), optional :: n
525 type(vector_t), intent(inout) :: a
526 type(vector_t), intent(in) :: b
527 type(vector_t), intent(in) :: c
528 integer :: size
529
530 if (present(n)) then
531 size = n
532 else
533 size = a%size()
534 end if
535
536 if (neko_bcknd_device .eq. 1) then
537 call device_col3(a%x_d, b%x_d, c%x_d, size)
538 else
539 call col3(a%x, b%x, c%x, size)
540 end if
541
542 end subroutine vector_col3
543
545 subroutine vector_subcol3(a, b, c, n)
546 integer, intent(in), optional :: n
547 type(vector_t), intent(inout) :: a
548 type(vector_t), intent(in) :: b
549 type(vector_t), intent(in) :: c
550 integer :: size
551
552 if (present(n)) then
553 size = n
554 else
555 size = a%size()
556 end if
557
558 if (neko_bcknd_device .eq. 1) then
559 call device_subcol3(a%x_d, b%x_d, c%x_d, size)
560 else
561 call subcol3(a%x, b%x, c%x, size)
562 end if
563
564 end subroutine vector_subcol3
565
567 subroutine vector_add3s2(a, b, c, c1, c2, n)
568 integer, intent(in), optional :: n
569 type(vector_t), intent(inout) :: a
570 type(vector_t), intent(in) :: b
571 type(vector_t), intent(in) :: c
572 real(kind=rp), intent(in) :: c1, c2
573 integer :: size
574
575 if (present(n)) then
576 size = n
577 else
578 size = a%size()
579 end if
580
581 if (neko_bcknd_device .eq. 1) then
582 call device_add3s2(a%x_d, b%x_d, c%x_d, c1, c2, size)
583 else
584 call add3s2(a%x, b%x, c%x, c1, c2, size)
585 end if
586
587 end subroutine vector_add3s2
588
590 subroutine vector_addcol3(a, b, c, n)
591 integer, intent(in), optional :: n
592 type(vector_t), intent(inout) :: a
593 type(vector_t), intent(in) :: b
594 type(vector_t), intent(in) :: c
595 integer :: size
596
597 if (present(n)) then
598 size = n
599 else
600 size = a%size()
601 end if
602
603 if (neko_bcknd_device .eq. 1) then
604 call device_addcol3(a%x_d, b%x_d, c%x_d, size)
605 else
606 call addcol3(a%x, b%x, c%x, size)
607 end if
608
609 end subroutine vector_addcol3
610
612 subroutine vector_addcol4(a, b, c, d, n)
613 integer, intent(in), optional :: n
614 type(vector_t), intent(inout) :: a
615 type(vector_t), intent(in) :: b
616 type(vector_t), intent(in) :: c
617 type(vector_t), intent(in) :: d
618 integer :: size
619
620 if (present(n)) then
621 size = n
622 else
623 size = a%size()
624 end if
625
626 if (neko_bcknd_device .eq. 1) then
627 call device_addcol4(a%x_d, b%x_d, c%x_d, d%x_d, size)
628 else
629 call addcol4(a%x, b%x, c%x, d%x, size)
630 end if
631
632 end subroutine vector_addcol4
633
634 function vector_glsum(a, n) result(sum)
635 integer, intent(in), optional :: n
636 type(vector_t), intent(in) :: a
637 real(kind=rp) :: sum
638 integer :: size
639
640 if (present(n)) then
641 size = n
642 else
643 size = a%size()
644 end if
645
646 if (neko_bcknd_device .eq. 1) then
647 sum = device_glsum(a%x_d, size)
648 else
649 sum = glsum(a%x, size)
650 end if
651
652 end function vector_glsum
653
654 function vector_glsc2(a, b, n) result(norm)
655 integer, intent(in), optional :: n
656 type(vector_t), intent(in) :: a, b
657 real(kind=rp) :: norm
658 integer :: size
659
660 if (present(n)) then
661 size = n
662 else
663 size = a%size()
664 end if
665
666 if (neko_bcknd_device .eq. 1) then
667 norm = device_glsc2(a%x_d, b%x_d, size)
668 else
669 norm = glsc2(a%x, b%x, size)
670 end if
671
672 end function vector_glsc2
673
674 function vector_glsc3(a, b, c, n) result(norm)
675 integer, intent(in), optional :: n
676 type(vector_t), intent(in) :: a, b, c
677 real(kind=rp) :: norm
678 integer :: size
679
680 if (present(n)) then
681 size = n
682 else
683 size = a%size()
684 end if
685
686 if (neko_bcknd_device .eq. 1) then
687 norm = device_glsc3(a%x_d, b%x_d, c%x_d, size)
688 else
689 norm = glsc3(a%x, b%x, c%x, size)
690 end if
691
692 end function vector_glsc3
693
694 function vector_glsubnorm(a, b, n) result(norm)
695 integer, intent(in), optional :: n
696 type(vector_t), intent(in) :: a, b
697 real(kind=rp) :: norm
698 integer :: size
699
700 if (present(n)) then
701 size = n
702 else
703 size = a%size()
704 end if
705
706 if (neko_bcknd_device .eq. 1) then
707 norm = device_glsubnorm(a%x_d, b%x_d, size)
708 else
709 norm = glsubnorm(a%x, b%x, size)
710 end if
711
712 end function vector_glsubnorm
713
722 subroutine vector_masked_gather_copy_0(a, b, mask, n, n_mask)
723 integer, intent(in) :: n, n_mask
724 real(kind=rp), dimension(n_mask), intent(inout) :: a
725 type(vector_t) :: b
726 integer, dimension(0:n_mask) :: mask
727 type(c_ptr) :: mask_d, a_d
728
729 if (neko_bcknd_device .eq. 1) then
730 mask_d = device_get_ptr(mask)
731 a_d = device_get_ptr(a)
732 call device_masked_gather_copy_0(a_d, b%x_d, mask_d, n, n_mask)
733 else
734 call masked_gather_copy_0(a, b%x, mask, n, n_mask)
735 end if
736
737 end subroutine vector_masked_gather_copy_0
738
747 subroutine vector_masked_scatter_copy_0(a, b, mask, n, n_mask)
748 integer, intent(in) :: n, n_mask
749 real(kind=rp), dimension(n_mask), intent(in) :: b
750 type(vector_t), intent(inout) :: a
751 integer, dimension(0:n_mask) :: mask
752 type(c_ptr) :: mask_d, b_d
753
754 if (neko_bcknd_device .eq. 1) then
755 mask_d = device_get_ptr(mask)
756 b_d = device_get_ptr(b)
757 call device_masked_scatter_copy_0(a%x_d, b_d, mask_d, n, n_mask)
758 else
759 call masked_scatter_copy_0(a%x, b, mask, n, n_mask)
760 end if
761
762 end subroutine vector_masked_scatter_copy_0
763
764
765
766end module vector_math
Return the device pointer for an associated Fortran array.
Definition device.F90:101
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)
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_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_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 .
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 .
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:411
subroutine, public cmult2(a, b, c, n)
Multiplication by constant c .
Definition math.f90:423
subroutine, public invcol2(a, b, n)
Vector division .
Definition math.f90:840
real(kind=rp) function, public glsc3(a, b, c, n)
Weighted inner product .
Definition math.f90:1067
subroutine, public cadd(a, s, n)
Add a scalar to vector .
Definition math.f90:462
subroutine, public addsqr2s2(a, b, c1, n)
Returns .
Definition math.f90:826
subroutine, public add2s1(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on first argument)
Definition math.f90:797
real(kind=rp) function, public glsc2(a, b, n)
Weighted inner product .
Definition math.f90:1048
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:358
subroutine, public subcol3(a, b, c, n)
Returns .
Definition math.f90:881
subroutine, public rone(a, n)
Set all elements to one.
Definition math.f90:238
subroutine, public add3(a, b, c, n)
Vector addition .
Definition math.f90:739
real(kind=rp) function, public glsum(a, n)
Sum a vector of length n.
Definition math.f90:499
subroutine, public sub3(a, b, c, n)
Vector subtraction .
Definition math.f90:781
subroutine, public addcol4(a, b, c, d, n)
Returns .
Definition math.f90:972
subroutine, public add2(a, b, n)
Vector addition .
Definition math.f90:726
subroutine, public cfill(a, c, n)
Set all elements to a constant c .
Definition math.f90:487
subroutine, public invcol3(a, b, c, n)
Invert a vector .
Definition math.f90:626
subroutine, public add3s2(a, b, c, c1, c2, n)
Returns .
Definition math.f90:895
subroutine, public masked_gather_copy_0(a, b, mask, n, n_mask)
Gather a masked vector to reduced contigous vector .
Definition math.f90:312
subroutine, public addcol3(a, b, c, n)
Returns .
Definition math.f90:958
subroutine, public invcol1(a, n)
Invert a vector .
Definition math.f90:614
subroutine, public col2(a, b, n)
Vector multiplication .
Definition math.f90:854
subroutine, public copy(a, b, n)
Copy a vector .
Definition math.f90:249
subroutine, public add4(a, b, c, d, n)
Vector addition .
Definition math.f90:753
subroutine, public col3(a, b, c, n)
Vector multiplication with 3 vectors .
Definition math.f90:867
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:684
subroutine, public rzero(a, n)
Zero a real vector.
Definition math.f90:205
real(kind=rp) function, public glsubnorm(a, b, n)
Returns the norm of the difference of two vectors .
Definition math.f90:1109
subroutine, public sub2(a, b, n)
Vector substraction .
Definition math.f90:768
subroutine, public add2s2(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on second argument)
Definition math.f90:812
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
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_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_sub2(a, b, n)
Vector substraction .
subroutine, public vector_addsqr2s2(a, b, c1, n)
Returns .
subroutine, public vector_col2(a, b, n)
Vector multiplication .
subroutine, public vector_cfill(a, c, n)
Set all elements to a constant c .
real(kind=rp) function, public vector_glsum(a, n)
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 .
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_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)
Gather a contigous array into a vector .
subroutine, public vector_invcol3(a, b, c, n)
Invert a vector .
Defines a vector.
Definition vector.f90:34