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