Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
matrix_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 matrix, only : matrix_t
64 use device
65 use math, only : rzero, rone, copy, cmult, cadd, cfill, invcol1, vdot3, &
79 use, intrinsic :: iso_c_binding, only : c_ptr
80 implicit none
81 private
82
91
92contains
93
95 subroutine matrix_rzero(a, n)
96 type(matrix_t), intent(inout) :: a
97 integer, intent(in), optional :: n
98 integer :: size
99
100 if (present(n)) then
101 size = n
102 else
103 size = a%size()
104 end if
105
106 if (neko_bcknd_device .eq. 1) then
107 call device_rzero(a%x_d, size)
108 else
109 call rzero(a%x, size)
110 end if
111 end subroutine matrix_rzero
112
114 subroutine matrix_rone(a, n)
115 type(matrix_t), intent(inout) :: a
116 integer, intent(in), optional :: n
117 integer :: size
118
119 if (present(n)) then
120 size = n
121 else
122 size = a%size()
123 end if
124
125 if (neko_bcknd_device .eq. 1) then
126 call device_rone(a%x_d, size)
127 else
128 call rone(a%x, size)
129 end if
130 end subroutine matrix_rone
131
133 subroutine matrix_copy(a, b, n)
134 type(matrix_t), intent(inout) :: a
135 type(matrix_t), intent(in) :: b
136 integer, intent(in), optional :: n
137 integer :: size
138
139 if (present(n)) then
140 size = n
141 else
142 size = a%size()
143 end if
144
145 if (neko_bcknd_device .eq. 1) then
146 call device_copy(a%x_d, b%x_d, size)
147 else
148 call copy(a%x, b%x, size)
149 end if
150 end subroutine matrix_copy
151
153 subroutine matrix_cmult(a, c, n)
154 type(matrix_t), intent(inout) :: a
155 real(kind=rp), intent(in) :: c
156 integer, intent(in), optional :: n
157 integer :: size
158
159 if (present(n)) then
160 size = n
161 else
162 size = a%size()
163 end if
164
165 if (neko_bcknd_device .eq. 1) then
166 call device_cmult(a%x_d, c, size)
167 else
168 call cmult(a%x, c, size)
169 end if
170 end subroutine matrix_cmult
171
173 subroutine matrix_cadd(a, s, n)
174 type(matrix_t), intent(inout) :: a
175 real(kind=rp), intent(in) :: s
176 integer, intent(in), optional :: n
177 integer :: size
178
179 if (present(n)) then
180 size = n
181 else
182 size = a%size()
183 end if
184
185 if (neko_bcknd_device .eq. 1) then
186 call device_cadd(a%x_d, s, size)
187 else
188 call cadd(a%x, s, size)
189 end if
190 end subroutine matrix_cadd
191
193 subroutine matrix_cfill(a, c, n)
194 type(matrix_t), intent(inout) :: a
195 real(kind=rp), intent(in) :: c
196 integer, intent(in), optional :: n
197 integer :: size
198
199 if (present(n)) then
200 size = n
201 else
202 size = a%size()
203 end if
204
205 if (neko_bcknd_device .eq. 1) then
206 call device_cfill(a%x_d, c, size)
207 else
208 call cfill(a%x, c, size)
209 end if
210 end subroutine matrix_cfill
211
213 subroutine matrix_invcol1(a, n)
214 type(matrix_t), intent(inout) :: a
215 integer, intent(in), optional :: n
216 integer :: size
217
218 if (present(n)) then
219 size = n
220 else
221 size = a%size()
222 end if
223
224 if (neko_bcknd_device .eq. 1) then
225 call device_invcol1(a%x_d, size)
226 else
227 call invcol1(a%x, size)
228 end if
229
230 end subroutine matrix_invcol1
231
233 subroutine matrix_invcol3(a, b, c, n)
234 type(matrix_t), intent(inout) :: a
235 type(matrix_t), intent(in) :: b
236 type(matrix_t), intent(in) :: c
237 integer, intent(in), optional :: n
238 integer :: size
239
240 if (present(n)) then
241 size = n
242 else
243 size = a%size()
244 end if
245
246 if (neko_bcknd_device .eq. 1) then
247 call device_invcol3(a%x_d, b%x_d, c%x_d, size)
248 else
249 call invcol3(a%x, b%x, c%x, size)
250 end if
251
252 end subroutine matrix_invcol3
253
255 subroutine matrix_add2(a, b, n)
256 type(matrix_t), intent(inout) :: a
257 type(matrix_t), intent(in) :: b
258 integer, intent(in), optional :: n
259 integer :: size
260
261 if (present(n)) then
262 size = n
263 else
264 size = a%size()
265 end if
266
267 if (neko_bcknd_device .eq. 1) then
268 call device_add2(a%x_d, b%x_d, size)
269 else
270 call add2(a%x, b%x, size)
271 end if
272
273 end subroutine matrix_add2
274
276 subroutine matrix_add3(a, b, c, n)
277 type(matrix_t), intent(inout) :: a
278 type(matrix_t), intent(in) :: b, c
279 integer, intent(in), optional :: n
280 integer :: size
281
282 if (present(n)) then
283 size = n
284 else
285 size = a%size()
286 end if
287
288 if (neko_bcknd_device .eq. 1) then
289 call device_add3(a%x_d, b%x_d, c%x_d, size)
290 else
291 call add3(a%x, b%x, c%x, size)
292 end if
293
294 end subroutine matrix_add3
295
297 subroutine matrix_add4(a, b, c, d, n)
298 type(matrix_t), intent(inout) :: a
299 type(matrix_t), intent(in) :: b, c, d
300 integer, intent(in), optional :: n
301 integer :: size
302
303 if (present(n)) then
304 size = n
305 else
306 size = a%size()
307 end if
308
309 if (neko_bcknd_device .eq. 1) then
310 call device_add4(a%x_d, b%x_d, c%x_d, d%x_d, size)
311 else
312 call add4(a%x, b%x, c%x, d%x, size)
313 end if
314
315 end subroutine matrix_add4
316
318 subroutine matrix_sub2(a, b, n)
319 type(matrix_t), intent(inout) :: a
320 type(matrix_t), intent(in) :: b
321 integer, intent(in), optional :: n
322 integer :: size
323
324 if (present(n)) then
325 size = n
326 else
327 size = a%size()
328 end if
329
330 if (neko_bcknd_device .eq. 1) then
331 call device_sub2(a%x_d, b%x_d, size)
332 else
333 call sub2(a%x, b%x, size)
334 end if
335
336 end subroutine matrix_sub2
337
339 subroutine matrix_sub3(a, b, c, n)
340 type(matrix_t), intent(inout) :: a
341 type(matrix_t), intent(in) :: b
342 type(matrix_t), intent(in) :: c
343 integer, intent(in), optional :: n
344 integer :: size
345
346 if (present(n)) then
347 size = n
348 else
349 size = a%size()
350 end if
351
352 if (neko_bcknd_device .eq. 1) then
353 call device_sub3(a%x_d, b%x_d, c%x_d, size)
354 else
355 call sub3(a%x, b%x, c%x, size)
356 end if
357
358 end subroutine matrix_sub3
359
360
363 subroutine matrix_add2s1(a, b, c1, n)
364 type(matrix_t), intent(inout) :: a
365 type(matrix_t), intent(in) :: b
366 real(kind=rp), intent(in) :: c1
367 integer, intent(in), optional :: n
368 integer :: size
369
370 if (present(n)) then
371 size = n
372 else
373 size = a%size()
374 end if
375
376 if (neko_bcknd_device .eq. 1) then
377 call device_add2s1(a%x_d, b%x_d, c1, size)
378 else
379 call add2s1(a%x, b%x, c1, size)
380 end if
381
382 end subroutine matrix_add2s1
383
386 subroutine matrix_add2s2(a, b, c1, n)
387 type(matrix_t), intent(inout) :: a
388 type(matrix_t), intent(in) :: b
389 real(kind=rp), intent(in) :: c1
390 integer, intent(in), optional :: n
391 integer :: size
392
393 if (present(n)) then
394 size = n
395 else
396 size = a%size()
397 end if
398
399 if (neko_bcknd_device .eq. 1) then
400 call device_add2s2(a%x_d, b%x_d, c1, size)
401 else
402 call add2s2(a%x, b%x, c1, size)
403 end if
404
405 end subroutine matrix_add2s2
406
408 subroutine matrix_addsqr2s2(a, b, c1, n)
409 type(matrix_t), intent(inout) :: a
410 type(matrix_t), intent(in) :: b
411 real(kind=rp), intent(in) :: c1
412 integer, intent(in), optional :: n
413 integer :: size
414
415 if (present(n)) then
416 size = n
417 else
418 size = a%size()
419 end if
420
421 if (neko_bcknd_device .eq. 1) then
422 call device_addsqr2s2(a%x_d, b%x_d, c1, size)
423 else
424 call addsqr2s2(a%x, b%x, c1, size)
425 end if
426
427 end subroutine matrix_addsqr2s2
428
430 subroutine matrix_cmult2(a, b, c, n)
431 type(matrix_t), intent(inout) :: a
432 type(matrix_t), intent(in) :: b
433 real(kind=rp), intent(in) :: c
434 integer, intent(in), optional :: n
435 integer :: size
436
437 if (present(n)) then
438 size = n
439 else
440 size = a%size()
441 end if
442
443 if (neko_bcknd_device .eq. 1) then
444 call device_cmult2(a%x_d, b%x_d, c, size)
445 else
446 call cmult2(a%x, b%x, c, size)
447 end if
448
449 end subroutine matrix_cmult2
450
452 subroutine matrix_invcol2(a, b, n)
453 type(matrix_t), intent(inout) :: a
454 type(matrix_t), intent(in) :: b
455 integer, intent(in), optional :: n
456 integer :: size
457
458 if (present(n)) then
459 size = n
460 else
461 size = a%size()
462 end if
463
464 if (neko_bcknd_device .eq. 1) then
465 call device_invcol2(a%x_d, b%x_d, size)
466 else
467 call invcol2(a%x, b%x, size)
468 end if
469
470 end subroutine matrix_invcol2
471
472
474 subroutine matrix_col2(a, b, n)
475 type(matrix_t), intent(inout) :: a
476 type(matrix_t), intent(in) :: b
477 integer, intent(in), optional :: n
478 integer :: size
479
480 if (present(n)) then
481 size = n
482 else
483 size = a%size()
484 end if
485
486 if (neko_bcknd_device .eq. 1) then
487 call device_col2(a%x_d, b%x_d, size)
488 else
489 call col2(a%x, b%x, size)
490 end if
491
492 end subroutine matrix_col2
493
495 subroutine matrix_col3(a, b, c, n)
496 type(matrix_t), intent(inout) :: a
497 type(matrix_t), intent(in) :: b
498 type(matrix_t), intent(in) :: c
499 integer, intent(in), optional :: n
500 integer :: size
501
502 if (present(n)) then
503 size = n
504 else
505 size = a%size()
506 end if
507
508 if (neko_bcknd_device .eq. 1) then
509 call device_col3(a%x_d, b%x_d, c%x_d, size)
510 else
511 call col3(a%x, b%x, c%x, size)
512 end if
513
514 end subroutine matrix_col3
515
517 subroutine matrix_subcol3(a, b, c, n)
518 type(matrix_t), intent(inout) :: a
519 type(matrix_t), intent(in) :: b
520 type(matrix_t), intent(in) :: c
521 integer, intent(in), optional :: n
522 integer :: size
523
524 if (present(n)) then
525 size = n
526 else
527 size = a%size()
528 end if
529
530 if (neko_bcknd_device .eq. 1) then
531 call device_subcol3(a%x_d, b%x_d, c%x_d, size)
532 else
533 call subcol3(a%x, b%x, c%x, size)
534 end if
535
536 end subroutine matrix_subcol3
537
539 subroutine matrix_add3s2(a, b, c, c1, c2, n)
540 type(matrix_t), intent(inout) :: a
541 type(matrix_t), intent(in) :: b
542 type(matrix_t), intent(in) :: c
543 real(kind=rp), intent(in) :: c1, c2
544 integer, intent(in), optional :: n
545 integer :: size
546
547 if (present(n)) then
548 size = n
549 else
550 size = a%size()
551 end if
552
553 if (neko_bcknd_device .eq. 1) then
554 call device_add3s2(a%x_d, b%x_d, c%x_d, c1, c2, size)
555 else
556 call add3s2(a%x, b%x, c%x, c1, c2, size)
557 end if
558
559 end subroutine matrix_add3s2
560
562 subroutine matrix_addcol3(a, b, c, n)
563 type(matrix_t), intent(inout) :: a
564 type(matrix_t), intent(in) :: b
565 type(matrix_t), intent(in) :: c
566 integer, intent(in), optional :: n
567 integer :: size
568
569 if (present(n)) then
570 size = n
571 else
572 size = a%size()
573 end if
574
575 if (neko_bcknd_device .eq. 1) then
576 call device_addcol3(a%x_d, b%x_d, c%x_d, size)
577 else
578 call addcol3(a%x, b%x, c%x, size)
579 end if
580
581 end subroutine matrix_addcol3
582
584 subroutine matrix_addcol4(a, b, c, d, n)
585 type(matrix_t), intent(inout) :: a
586 type(matrix_t), intent(in) :: b
587 type(matrix_t), intent(in) :: c
588 type(matrix_t), intent(in) :: d
589 integer, intent(in), optional :: n
590 integer :: size
591
592 if (present(n)) then
593 size = n
594 else
595 size = a%size()
596 end if
597
598 if (neko_bcknd_device .eq. 1) then
599 call device_addcol4(a%x_d, b%x_d, c%x_d, d%x_d, size)
600 else
601 call addcol4(a%x, b%x, c%x, d%x, size)
602 end if
603
604 end subroutine matrix_addcol4
605
607 function matrix_glsum(a, n) result(sum)
608 integer, intent(in), optional :: n
609 type(matrix_t), intent(in) :: a
610 real(kind=rp) :: sum
611 integer :: size
612
613 if (present(n)) then
614 size = n
615 else
616 size = a%size()
617 end if
618
619 if (neko_bcknd_device .eq. 1) then
620 sum = device_glsum(a%x_d, size)
621 else
622 sum = glsum(a%x, size)
623 end if
624
625 end function matrix_glsum
626
628 function matrix_glmax(a, n) result(val)
629 integer, intent(in), optional :: n
630 type(matrix_t), intent(in) :: a
631 real(kind=rp) :: val
632 integer :: size
633
634 if (present(n)) then
635 size = n
636 else
637 size = a%size()
638 end if
639
640 if (neko_bcknd_device .eq. 1) then
641 val = device_glmax(a%x_d, size)
642 else
643 val = glmax(a%x, size)
644 end if
645
646 end function matrix_glmax
647
649 function matrix_glmin(a, n) result(val)
650 integer, intent(in), optional :: n
651 type(matrix_t), intent(in) :: a
652 real(kind=rp) :: val
653 integer :: size
654
655 if (present(n)) then
656 size = n
657 else
658 size = a%size()
659 end if
660
661 if (neko_bcknd_device .eq. 1) then
662 val = device_glmin(a%x_d, size)
663 else
664 val = glmin(a%x, size)
665 end if
666
667 end function matrix_glmin
668
670 function matrix_glsc2(a, b, n) result(ip)
671 integer, intent(in), optional :: n
672 type(matrix_t), intent(in) :: a, b
673 real(kind=rp) :: ip
674 integer :: size
675
676 if (present(n)) then
677 size = n
678 else
679 size = a%size()
680 end if
681
682 if (neko_bcknd_device .eq. 1) then
683 ip = device_glsc2(a%x_d, b%x_d, size)
684 else
685 ip = glsc2(a%x, b%x, size)
686 end if
687
688 end function matrix_glsc2
689
691 function matrix_glsc3(a, b, c, n) result(ip)
692 integer, intent(in), optional :: n
693 type(matrix_t), intent(in) :: a, b, c
694 real(kind=rp) :: ip
695 integer :: size
696
697 if (present(n)) then
698 size = n
699 else
700 size = a%size()
701 end if
702
703 if (neko_bcknd_device .eq. 1) then
704 ip = device_glsc3(a%x_d, b%x_d, c%x_d, size)
705 else
706 ip = glsc3(a%x, b%x, c%x, size)
707 end if
708
709 end function matrix_glsc3
710
713 function matrix_glsubnorm(a, b, n) result(norm)
714 integer, intent(in), optional :: n
715 type(matrix_t), intent(in) :: a, b
716 real(kind=rp) :: norm
717 integer :: size
718
719 if (present(n)) then
720 size = n
721 else
722 size = a%size()
723 end if
724
725 if (neko_bcknd_device .eq. 1) then
726 norm = device_glsubnorm(a%x_d, b%x_d, size)
727 else
728 norm = glsubnorm(a%x, b%x, size)
729 end if
730
731 end function matrix_glsubnorm
732
733end module matrix_math
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_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 .
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
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
real(kind=rp) function, public glmax(a, n)
Max of a vector of length n.
Definition math.f90:516
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
real(kind=rp) function, public glmin(a, n)
Min of a vector of length n.
Definition math.f90:546
real(kind=rp) function, public matrix_glmax(a, n)
Global maximum of all elements in a matrix .
real(kind=rp) function, public matrix_glsum(a, n)
Global sum of all elements in a matrix .
subroutine, public matrix_copy(a, b, n)
Copy a matrix .
subroutine, public matrix_invcol2(a, b, n)
Vector division .
subroutine, public matrix_add2(a, b, n)
Vector addition .
subroutine, public matrix_addcol4(a, b, c, d, n)
Returns .
subroutine, public matrix_cadd(a, s, n)
Add a scalar to matrix .
subroutine, public matrix_cmult(a, c, n)
Multiplication by constant c .
subroutine, public matrix_add3s2(a, b, c, c1, c2, n)
Returns .
subroutine, public matrix_cmult2(a, b, c, n)
Multiplication by constant c .
subroutine matrix_add4(a, b, c, d, n)
Vector addition .
subroutine, public matrix_col3(a, b, c, n)
Vector multiplication with 3 vectors .
subroutine, public matrix_cfill(a, c, n)
Set all elements to a constant c .
subroutine, public matrix_invcol1(a, n)
Invert elements of a matrix .
real(kind=rp) function, public matrix_glsc3(a, b, c, n)
Global inner product of three matrices .
subroutine, public matrix_sub2(a, b, n)
Vector substraction .
subroutine, public matrix_subcol3(a, b, c, n)
Returns .
subroutine, public matrix_addsqr2s2(a, b, c1, n)
Returns .
real(kind=rp) function, public matrix_glmin(a, n)
Global minimum of all elements in a matrix .
subroutine, public matrix_add3(a, b, c, n)
Vector addition .
subroutine, public matrix_rzero(a, n)
Zero a real matrix .
real(kind=rp) function, public matrix_glsc2(a, b, n)
Global inner product of two matrices .
subroutine, public matrix_add2s2(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on second argument)
subroutine, public matrix_addcol3(a, b, c, n)
Returns .
subroutine, public matrix_add2s1(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on first argument)
subroutine, public matrix_invcol3(a, b, c, n)
Element division of two matrices .
real(kind=rp) function, public matrix_glsubnorm(a, b, n)
Global subtracted norm of two matrices .
subroutine, public matrix_col2(a, b, n)
Vector multiplication .
subroutine, public matrix_rone(a, n)
Set all elements to one.
subroutine, public matrix_sub3(a, b, c, n)
Vector subtraction .
Defines a matrix.
Definition matrix.f90:34
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12