Neko 1.99.1
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, add2, &
77 use, intrinsic :: iso_c_binding, only: c_ptr
78 implicit none
79 private
80
89
90contains
91
93 subroutine matrix_rzero(a, n)
94 type(matrix_t), intent(inout) :: a
95 integer, intent(in), optional :: n
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 matrix_rzero
110
112 subroutine matrix_rone(a, n)
113 type(matrix_t), intent(inout) :: a
114 integer, intent(in), optional :: n
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 matrix_rone
129
131 subroutine matrix_copy(a, b, n)
132 type(matrix_t), intent(inout) :: a
133 type(matrix_t), intent(in) :: b
134 integer, intent(in), optional :: n
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 matrix_copy
149
151 subroutine matrix_cmult(a, c, n)
152 type(matrix_t), intent(inout) :: a
153 real(kind=rp), intent(in) :: c
154 integer, intent(in), optional :: n
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 matrix_cmult
169
171 subroutine matrix_cadd(a, s, n)
172 type(matrix_t), intent(inout) :: a
173 real(kind=rp), intent(in) :: s
174 integer, intent(in), optional :: n
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 matrix_cadd
189
191 subroutine matrix_cfill(a, c, n)
192 type(matrix_t), intent(inout) :: a
193 real(kind=rp), intent(in) :: c
194 integer, intent(in), optional :: n
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 matrix_cfill
209
211 subroutine matrix_invcol1(a, n)
212 type(matrix_t), intent(inout) :: a
213 integer, intent(in), optional :: n
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 matrix_invcol1
229
231 subroutine matrix_invcol3(a, b, c, n)
232 type(matrix_t), intent(inout) :: a
233 type(matrix_t), intent(in) :: b
234 type(matrix_t), intent(in) :: c
235 integer, intent(in), optional :: n
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 matrix_invcol3
251
253 subroutine matrix_add2(a, b, n)
254 type(matrix_t), intent(inout) :: a
255 type(matrix_t), intent(in) :: b
256 integer, intent(in), optional :: n
257 integer :: size
258
259 if (present(n)) then
260 size = n
261 else
262 size = a%size()
263 end if
264
265 if (neko_bcknd_device .eq. 1) then
266 call device_add2(a%x_d, b%x_d, size)
267 else
268 call add2(a%x, b%x, size)
269 end if
270
271 end subroutine matrix_add2
272
274 subroutine matrix_add3(a, b, c, n)
275 type(matrix_t), intent(inout) :: a
276 type(matrix_t), intent(in) :: b, c
277 integer, intent(in), optional :: n
278 integer :: size
279
280 if (present(n)) then
281 size = n
282 else
283 size = a%size()
284 end if
285
286 if (neko_bcknd_device .eq. 1) then
287 call device_add3(a%x_d, b%x_d, c%x_d, size)
288 else
289 call add3(a%x, b%x, c%x, size)
290 end if
291
292 end subroutine matrix_add3
293
295 subroutine matrix_add4(a, b, c, d, n)
296 type(matrix_t), intent(inout) :: a
297 type(matrix_t), intent(in) :: b, c, d
298 integer, intent(in), optional :: n
299 integer :: size
300
301 if (present(n)) then
302 size = n
303 else
304 size = a%size()
305 end if
306
307 if (neko_bcknd_device .eq. 1) then
308 call device_add4(a%x_d, b%x_d, c%x_d, d%x_d, size)
309 else
310 call add4(a%x, b%x, c%x, d%x, size)
311 end if
312
313 end subroutine matrix_add4
314
316 subroutine matrix_sub2(a, b, n)
317 type(matrix_t), intent(inout) :: a
318 type(matrix_t), intent(in) :: b
319 integer, intent(in), optional :: n
320 integer :: size
321
322 if (present(n)) then
323 size = n
324 else
325 size = a%size()
326 end if
327
328 if (neko_bcknd_device .eq. 1) then
329 call device_sub2(a%x_d, b%x_d, size)
330 else
331 call sub2(a%x, b%x, size)
332 end if
333
334 end subroutine matrix_sub2
335
337 subroutine matrix_sub3(a, b, c, n)
338 type(matrix_t), intent(inout) :: a
339 type(matrix_t), intent(in) :: b
340 type(matrix_t), intent(in) :: c
341 integer, intent(in), optional :: n
342 integer :: size
343
344 if (present(n)) then
345 size = n
346 else
347 size = a%size()
348 end if
349
350 if (neko_bcknd_device .eq. 1) then
351 call device_sub3(a%x_d, b%x_d, c%x_d, size)
352 else
353 call sub3(a%x, b%x, c%x, size)
354 end if
355
356 end subroutine matrix_sub3
357
358
361 subroutine matrix_add2s1(a, b, c1, n)
362 type(matrix_t), intent(inout) :: a
363 type(matrix_t), intent(in) :: b
364 real(kind=rp), intent(in) :: c1
365 integer, intent(in), optional :: n
366 integer :: size
367
368 if (present(n)) then
369 size = n
370 else
371 size = a%size()
372 end if
373
374 if (neko_bcknd_device .eq. 1) then
375 call device_add2s1(a%x_d, b%x_d, c1, size)
376 else
377 call add2s1(a%x, b%x, c1, size)
378 end if
379
380 end subroutine matrix_add2s1
381
384 subroutine matrix_add2s2(a, b, c1, n)
385 type(matrix_t), intent(inout) :: a
386 type(matrix_t), intent(in) :: b
387 real(kind=rp), intent(in) :: c1
388 integer, intent(in), optional :: n
389 integer :: size
390
391 if (present(n)) then
392 size = n
393 else
394 size = a%size()
395 end if
396
397 if (neko_bcknd_device .eq. 1) then
398 call device_add2s2(a%x_d, b%x_d, c1, size)
399 else
400 call add2s2(a%x, b%x, c1, size)
401 end if
402
403 end subroutine matrix_add2s2
404
406 subroutine matrix_addsqr2s2(a, b, c1, n)
407 type(matrix_t), intent(inout) :: a
408 type(matrix_t), intent(in) :: b
409 real(kind=rp), intent(in) :: c1
410 integer, intent(in), optional :: n
411 integer :: size
412
413 if (present(n)) then
414 size = n
415 else
416 size = a%size()
417 end if
418
419 if (neko_bcknd_device .eq. 1) then
420 call device_addsqr2s2(a%x_d, b%x_d, c1, size)
421 else
422 call addsqr2s2(a%x, b%x, c1, size)
423 end if
424
425 end subroutine matrix_addsqr2s2
426
428 subroutine matrix_cmult2(a, b, c, n)
429 type(matrix_t), intent(inout) :: a
430 type(matrix_t), intent(in) :: b
431 real(kind=rp), intent(in) :: c
432 integer, intent(in), optional :: n
433 integer :: size
434
435 if (present(n)) then
436 size = n
437 else
438 size = a%size()
439 end if
440
441 if (neko_bcknd_device .eq. 1) then
442 call device_cmult2(a%x_d, b%x_d, c, size)
443 else
444 call cmult2(a%x, b%x, c, size)
445 end if
446
447 end subroutine matrix_cmult2
448
450 subroutine matrix_invcol2(a, b, n)
451 type(matrix_t), intent(inout) :: a
452 type(matrix_t), intent(in) :: b
453 integer, intent(in), optional :: n
454 integer :: size
455
456 if (present(n)) then
457 size = n
458 else
459 size = a%size()
460 end if
461
462 if (neko_bcknd_device .eq. 1) then
463 call device_invcol2(a%x_d, b%x_d, size)
464 else
465 call invcol2(a%x, b%x, size)
466 end if
467
468 end subroutine matrix_invcol2
469
470
472 subroutine matrix_col2(a, b, n)
473 type(matrix_t), intent(inout) :: a
474 type(matrix_t), intent(in) :: b
475 integer, intent(in), optional :: n
476 integer :: size
477
478 if (present(n)) then
479 size = n
480 else
481 size = a%size()
482 end if
483
484 if (neko_bcknd_device .eq. 1) then
485 call device_col2(a%x_d, b%x_d, size)
486 else
487 call col2(a%x, b%x, size)
488 end if
489
490 end subroutine matrix_col2
491
493 subroutine matrix_col3(a, b, c, n)
494 type(matrix_t), intent(inout) :: a
495 type(matrix_t), intent(in) :: b
496 type(matrix_t), intent(in) :: c
497 integer, intent(in), optional :: n
498 integer :: size
499
500 if (present(n)) then
501 size = n
502 else
503 size = a%size()
504 end if
505
506 if (neko_bcknd_device .eq. 1) then
507 call device_col3(a%x_d, b%x_d, c%x_d, size)
508 else
509 call col3(a%x, b%x, c%x, size)
510 end if
511
512 end subroutine matrix_col3
513
515 subroutine matrix_subcol3(a, b, c, n)
516 type(matrix_t), intent(inout) :: a
517 type(matrix_t), intent(in) :: b
518 type(matrix_t), intent(in) :: c
519 integer, intent(in), optional :: n
520 integer :: size
521
522 if (present(n)) then
523 size = n
524 else
525 size = a%size()
526 end if
527
528 if (neko_bcknd_device .eq. 1) then
529 call device_subcol3(a%x_d, b%x_d, c%x_d, size)
530 else
531 call subcol3(a%x, b%x, c%x, size)
532 end if
533
534 end subroutine matrix_subcol3
535
537 subroutine matrix_add3s2(a, b, c, c1, c2, n)
538 type(matrix_t), intent(inout) :: a
539 type(matrix_t), intent(in) :: b
540 type(matrix_t), intent(in) :: c
541 real(kind=rp), intent(in) :: c1, c2
542 integer, intent(in), optional :: n
543 integer :: size
544
545 if (present(n)) then
546 size = n
547 else
548 size = a%size()
549 end if
550
551 if (neko_bcknd_device .eq. 1) then
552 call device_add3s2(a%x_d, b%x_d, c%x_d, c1, c2, size)
553 else
554 call add3s2(a%x, b%x, c%x, c1, c2, size)
555 end if
556
557 end subroutine matrix_add3s2
558
560 subroutine matrix_addcol3(a, b, c, n)
561 type(matrix_t), intent(inout) :: a
562 type(matrix_t), intent(in) :: b
563 type(matrix_t), intent(in) :: c
564 integer, intent(in), optional :: n
565 integer :: size
566
567 if (present(n)) then
568 size = n
569 else
570 size = a%size()
571 end if
572
573 if (neko_bcknd_device .eq. 1) then
574 call device_addcol3(a%x_d, b%x_d, c%x_d, size)
575 else
576 call addcol3(a%x, b%x, c%x, size)
577 end if
578
579 end subroutine matrix_addcol3
580
582 subroutine matrix_addcol4(a, b, c, d, n)
583 type(matrix_t), intent(inout) :: a
584 type(matrix_t), intent(in) :: b
585 type(matrix_t), intent(in) :: c
586 type(matrix_t), intent(in) :: d
587 integer, intent(in), optional :: n
588 integer :: size
589
590 if (present(n)) then
591 size = n
592 else
593 size = a%size()
594 end if
595
596 if (neko_bcknd_device .eq. 1) then
597 call device_addcol4(a%x_d, b%x_d, c%x_d, d%x_d, size)
598 else
599 call addcol4(a%x, b%x, c%x, d%x, size)
600 end if
601
602 end subroutine matrix_addcol4
603
605 function matrix_glsum(a, n) result(sum)
606 integer, intent(in), optional :: n
607 type(matrix_t), intent(in) :: a
608 real(kind=rp) :: sum
609 integer :: size
610
611 if (present(n)) then
612 size = n
613 else
614 size = a%size()
615 end if
616
617 if (neko_bcknd_device .eq. 1) then
618 sum = device_glsum(a%x_d, size)
619 else
620 sum = glsum(a%x, size)
621 end if
622
623 end function matrix_glsum
624
626 function matrix_glsc2(a, b, n) result(ip)
627 integer, intent(in), optional :: n
628 type(matrix_t), intent(in) :: a, b
629 real(kind=rp) :: ip
630 integer :: size
631
632 if (present(n)) then
633 size = n
634 else
635 size = a%size()
636 end if
637
638 if (neko_bcknd_device .eq. 1) then
639 ip = device_glsc2(a%x_d, b%x_d, size)
640 else
641 ip = glsc2(a%x, b%x, size)
642 end if
643
644 end function matrix_glsc2
645
647 function matrix_glsc3(a, b, c, n) result(ip)
648 integer, intent(in), optional :: n
649 type(matrix_t), intent(in) :: a, b, c
650 real(kind=rp) :: ip
651 integer :: size
652
653 if (present(n)) then
654 size = n
655 else
656 size = a%size()
657 end if
658
659 if (neko_bcknd_device .eq. 1) then
660 ip = device_glsc3(a%x_d, b%x_d, c%x_d, size)
661 else
662 ip = glsc3(a%x, b%x, c%x, size)
663 end if
664
665 end function matrix_glsc3
666
669 function matrix_glsubnorm(a, b, n) result(norm)
670 integer, intent(in), optional :: n
671 type(matrix_t), intent(in) :: a, b
672 real(kind=rp) :: norm
673 integer :: size
674
675 if (present(n)) then
676 size = n
677 else
678 size = a%size()
679 end if
680
681 if (neko_bcknd_device .eq. 1) then
682 norm = device_glsubnorm(a%x_d, b%x_d, size)
683 else
684 norm = glsubnorm(a%x, b%x, size)
685 end if
686
687 end function matrix_glsubnorm
688
689end 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)
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
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
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 .
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