Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
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!
60module math
61 use num_types, only : rp, dp, sp, qp, i4, xp
63 use mpi_f08, only : mpi_min, mpi_max, mpi_sum, mpi_in_place, mpi_integer, &
64 mpi_allreduce
65 implicit none
66 private
67
69 real(kind=rp), public, parameter :: neko_eps = epsilon(1.0_rp)
70
72 real(kind=rp), public, parameter :: neko_m_ln2 = log(2.0_rp)
73
75 real(kind=rp), public, parameter :: pi = 4._rp*atan(1._rp)
76
77 interface abscmp
78 module procedure sabscmp, dabscmp, qabscmp
79 end interface abscmp
80
81 interface sort
82 module procedure sortrp, sorti4
83 end interface sort
84
85 interface swap
86 module procedure swapdp, swapi4
87 end interface swap
88
89 interface reord
90 module procedure reorddp, reordi4
91 end interface reord
92
93 interface flipv
94 module procedure flipvdp, flipvi4
95 end interface flipv
96
97 interface relcmp
98 module procedure srelcmp, drelcmp, qrelcmp
99 end interface relcmp
100
101 public :: abscmp, rzero, izero, row_zero, rone, copy, cmult, cadd, cfill, &
113
114contains
115
117 pure function sabscmp(x, y, tol)
118 real(kind=sp), intent(in) :: x
119 real(kind=sp), intent(in) :: y
120 real(kind=sp), intent(in), optional :: tol
121 logical :: sabscmp
122
123 if (present(tol)) then
124 sabscmp = abs(x - y) .lt. tol
125 else
126 sabscmp = abs(x - y) .lt. neko_eps
127 end if
128
129 end function sabscmp
130
132 pure function dabscmp(x, y, tol)
133 real(kind=dp), intent(in) :: x
134 real(kind=dp), intent(in) :: y
135 real(kind=dp), intent(in), optional :: tol
136 logical :: dabscmp
137
138 if (present(tol)) then
139 dabscmp = abs(x - y) .lt. tol
140 else
141 dabscmp = abs(x - y) .lt. neko_eps
142 end if
143
144 end function dabscmp
145
147 pure function qabscmp(x, y, tol)
148 real(kind=qp), intent(in) :: x
149 real(kind=qp), intent(in) :: y
150 real(kind=qp), intent(in), optional :: tol
151 logical :: qabscmp
152
153 if (present(tol)) then
154 qabscmp = abs(x - y) .lt. tol
155 else
156 qabscmp = abs(x - y) .lt. neko_eps
157 end if
158
159 end function qabscmp
160
162 pure function srelcmp(x, y, eps)
163 real(kind=sp), intent(in) :: x
164 real(kind=sp), intent(in) :: y
165 real(kind=sp), intent(in), optional :: eps
166 logical :: srelcmp
167 if (present(eps)) then
168 srelcmp = abs(x - y) .le. eps*abs(y)
169 else
170 srelcmp = abs(x - y) .le. neko_eps*abs(y)
171 end if
172
173 end function srelcmp
174
176 pure function drelcmp(x, y, eps)
177 real(kind=dp), intent(in) :: x
178 real(kind=dp), intent(in) :: y
179 real(kind=dp), intent(in), optional :: eps
180 logical :: drelcmp
181 if (present(eps)) then
182 drelcmp = abs(x - y) .le. eps*abs(y)
183 else
184 drelcmp = abs(x - y) .le. neko_eps*abs(y)
185 end if
186
187 end function drelcmp
188
189
191 pure function qrelcmp(x, y, eps)
192 real(kind=qp), intent(in) :: x
193 real(kind=qp), intent(in) :: y
194 real(kind=qp), intent(in), optional :: eps
195 logical :: qrelcmp
196 if (present(eps)) then
197 qrelcmp = abs(x - y)/abs(y) .lt. eps
198 else
199 qrelcmp = abs(x - y)/abs(y) .lt. neko_eps
200 end if
201
202 end function qrelcmp
203
205 subroutine rzero(a, n)
206 integer, intent(in) :: n
207 real(kind=rp), dimension(n), intent(inout) :: a
208 integer :: i
209
210 do concurrent(i = 1:n)
211 a(i) = 0.0_rp
212 end do
213 end subroutine rzero
214
216 subroutine izero(a, n)
217 integer, intent(in) :: n
218 integer, dimension(n), intent(inout) :: a
219 integer :: i
220
221 do concurrent(i = 1:n)
222 a(i) = 0
223 end do
224 end subroutine izero
225
227 subroutine row_zero(a, m, n, e)
228 integer, intent(in) :: m, n, e
229 real(kind=rp), intent(inout) :: a(m,n)
230 integer :: j
231
232 do concurrent(j = 1:n)
233 a(e,j) = 0.0_rp
234 end do
235 end subroutine row_zero
236
238 subroutine rone(a, n)
239 integer, intent(in) :: n
240 real(kind=rp), dimension(n), intent(inout) :: a
241 integer :: i
242
243 do concurrent(i = 1:n)
244 a(i) = 1.0_rp
245 end do
246 end subroutine rone
247
249 subroutine copy(a, b, n)
250 integer, intent(in) :: n
251 real(kind=rp), dimension(n), intent(in) :: b
252 real(kind=rp), dimension(n), intent(inout) :: a
253 integer :: i
254
255 do concurrent(i = 1:n)
256 a(i) = b(i)
257 end do
258
259 end subroutine copy
260
268 subroutine masked_copy_0(a, b, mask, n, n_mask)
269 integer, intent(in) :: n, n_mask
270 real(kind=rp), dimension(n), intent(in) :: b
271 real(kind=rp), dimension(n), intent(inout) :: a
272 integer, dimension(0:n_mask) :: mask
273 integer :: i, j
274
275 do concurrent(i = 1:n_mask)
276 j = mask(i)
277 a(j) = b(j)
278 end do
279
280 end subroutine masked_copy_0
281
289 subroutine masked_copy(a, b, mask, n, n_mask)
290 integer, intent(in) :: n, n_mask
291 real(kind=rp), dimension(n), intent(in) :: b
292 real(kind=rp), dimension(n), intent(inout) :: a
293 integer, dimension(n_mask) :: mask
294 integer :: i, j
295
296 do concurrent(i = 1:n_mask)
297 j = mask(i)
298 a(j) = b(j)
299 end do
300
301 end subroutine masked_copy
302
312 subroutine masked_gather_copy_0(a, b, mask, n, n_mask)
313 integer, intent(in) :: n, n_mask
314 real(kind=rp), dimension(n), intent(in) :: b
315 real(kind=rp), dimension(n_mask), intent(inout) :: a
316 integer, dimension(0:n_mask) :: mask
317 integer :: i, j
318
319 do concurrent(i = 1:n_mask)
320 j = mask(i)
321 a(i) = b(j)
322 end do
323
324 end subroutine masked_gather_copy_0
325
335 subroutine masked_gather_copy(a, b, mask, n, n_mask)
336 integer, intent(in) :: n, n_mask
337 real(kind=rp), dimension(n), intent(in) :: b
338 real(kind=rp), dimension(n_mask), intent(inout) :: a
339 integer, dimension(n_mask) :: mask
340 integer :: i, j
341
342 do concurrent(i = 1:n_mask)
343 j = mask(i)
344 a(i) = b(j)
345 end do
346
347 end subroutine masked_gather_copy
348
358 subroutine masked_scatter_copy_0(a, b, mask, n, n_mask)
359 integer, intent(in) :: n, n_mask
360 real(kind=rp), dimension(n_mask), intent(in) :: b
361 real(kind=rp), dimension(n), intent(inout) :: a
362 integer, dimension(0:n_mask) :: mask
363 integer :: i, j
364
365 do concurrent(i = 1:n_mask)
366 j = mask(i)
367 a(j) = b(i)
368 end do
369
370 end subroutine masked_scatter_copy_0
371
381 subroutine masked_scatter_copy(a, b, mask, n, n_mask)
382 integer, intent(in) :: n, n_mask
383 real(kind=rp), dimension(n_mask), intent(in) :: b
384 real(kind=rp), dimension(n), intent(inout) :: a
385 integer, dimension(n_mask) :: mask
386 integer :: i, j
387
388 do concurrent(i = 1:n_mask)
389 j = mask(i)
390 a(j) = b(i)
391 end do
392
393 end subroutine masked_scatter_copy
394
397 subroutine cfill_mask(a, c, n, mask, n_mask)
398 integer, intent(in) :: n, n_mask
399 real(kind=rp), dimension(n), intent(inout) :: a
400 real(kind=rp), intent(in) :: c
401 integer, dimension(n_mask), intent(in) :: mask
402 integer :: i
403
404 do concurrent(i = 1:n_mask)
405 a(mask(i)) = c
406 end do
407
408 end subroutine cfill_mask
409
411 subroutine cmult(a, c, n)
412 integer, intent(in) :: n
413 real(kind=rp), dimension(n), intent(inout) :: a
414 real(kind=rp), intent(in) :: c
415 integer :: i
416
417 do concurrent(i = 1:n)
418 a(i) = c * a(i)
419 end do
420 end subroutine cmult
421
423 subroutine cmult2(a, b, c, n)
424 integer, intent(in) :: n
425 real(kind=rp), dimension(n), intent(inout) :: a
426 real(kind=rp), dimension(n), intent(in) :: b
427 real(kind=rp), intent(in) :: c
428 integer :: i
429
430 do concurrent(i = 1:n)
431 a(i) = c * b(i)
432 end do
433
434 end subroutine cmult2
435
437 subroutine cdiv(a, c, n)
438 integer, intent(in) :: n
439 real(kind=rp), dimension(n), intent(inout) :: a
440 real(kind=rp), intent(in) :: c
441 integer :: i
442
443 do concurrent(i = 1:n)
444 a(i) = c / a(i)
445 end do
446 end subroutine cdiv
447
449 subroutine cdiv2(a, b, c, n)
450 integer, intent(in) :: n
451 real(kind=rp), dimension(n), intent(inout) :: a
452 real(kind=rp), dimension(n), intent(in) :: b
453 real(kind=rp), intent(in) :: c
454 integer :: i
455
456 do concurrent(i = 1:n)
457 a(i) = c / b(i)
458 end do
459 end subroutine cdiv2
460
462 subroutine cadd(a, s, n)
463 integer, intent(in) :: n
464 real(kind=rp), dimension(n), intent(inout) :: a
465 real(kind=rp), intent(in) :: s
466 integer :: i
467
468 do concurrent(i = 1:n)
469 a(i) = a(i) + s
470 end do
471 end subroutine cadd
472
474 subroutine cadd2(a, b, s, n)
475 integer, intent(in) :: n
476 real(kind=rp), dimension(n), intent(inout) :: a
477 real(kind=rp), dimension(n), intent(in) :: b
478 real(kind=rp), intent(in) :: s
479 integer :: i
480
481 do concurrent(i = 1:n)
482 a(i) = b(i) + s
483 end do
484 end subroutine cadd2
485
487 subroutine cfill(a, c, n)
488 integer, intent(in) :: n
489 real(kind=rp), dimension(n), intent(inout) :: a
490 real(kind=rp), intent(in) :: c
491 integer :: i
492
493 do concurrent(i = 1:n)
494 a(i) = c
495 end do
496 end subroutine cfill
497
499 function glsum(a, n)
500 integer, intent(in) :: n
501 real(kind=rp), dimension(n) :: a
502 real(kind=rp) :: glsum
503 real(kind=xp) :: tmp
504 integer :: i, ierr
505 tmp = 0.0_rp
506 do i = 1, n
507 tmp = tmp + a(i)
508 end do
509 call mpi_allreduce(mpi_in_place, tmp, 1, &
510 mpi_extra_precision, mpi_sum, neko_comm, ierr)
511 glsum = tmp
512
513 end function glsum
514
516 function glmax(a, n)
517 integer, intent(in) :: n
518 real(kind=rp), dimension(n) :: a
519 real(kind=rp) :: tmp, glmax
520 integer :: i, ierr
521
522 tmp = -huge(0.0_rp)
523 do i = 1, n
524 tmp = max(tmp,a(i))
525 end do
526 call mpi_allreduce(tmp, glmax, 1, &
527 mpi_real_precision, mpi_max, neko_comm, ierr)
528 end function glmax
529
531 function glimax(a, n)
532 integer, intent(in) :: n
533 integer, dimension(n) :: a
534 integer :: tmp, glimax
535 integer :: i, ierr
536
537 tmp = -huge(0)
538 do i = 1, n
539 tmp = max(tmp,a(i))
540 end do
541 call mpi_allreduce(tmp, glimax, 1, &
542 mpi_integer, mpi_max, neko_comm, ierr)
543 end function glimax
544
546 function glmin(a, n)
547 integer, intent(in) :: n
548 real(kind=rp), dimension(n) :: a
549 real(kind=rp) :: tmp, glmin
550 integer :: i, ierr
551
552 tmp = huge(0.0_rp)
553 do i = 1, n
554 tmp = min(tmp,a(i))
555 end do
556 call mpi_allreduce(tmp, glmin, 1, &
557 mpi_real_precision, mpi_min, neko_comm, ierr)
558 end function glmin
559
561 function glimin(a, n)
562 integer, intent(in) :: n
563 integer, dimension(n) :: a
564 integer :: tmp, glimin
565 integer :: i, ierr
566
567 tmp = huge(0)
568 do i = 1, n
569 tmp = min(tmp,a(i))
570 end do
571 call mpi_allreduce(tmp, glimin, 1, &
572 mpi_integer, mpi_min, neko_comm, ierr)
573 end function glimin
574
575
576
577
579 subroutine chsign(a, n)
580 integer, intent(in) :: n
581 real(kind=rp), dimension(n), intent(inout) :: a
582 integer :: i
583
584 do concurrent(i = 1:n)
585 a(i) = -a(i)
586 end do
587
588 end subroutine chsign
589
591 function vlmax(vec,n) result(tmax)
592 integer :: n, i
593 real(kind=rp), intent(in) :: vec(n)
594 real(kind=rp) :: tmax
595 tmax = real(-99d20, rp)
596 do i = 1, n
597 tmax = max(tmax, vec(i))
598 end do
599 end function vlmax
600
602 function vlmin(vec,n) result(tmin)
603 integer, intent(in) :: n
604 real(kind=rp), intent(in) :: vec(n)
605 real(kind=rp) :: tmin
606 integer :: i
607 tmin = real(99.0e20, rp)
608 do i = 1, n
609 tmin = min(tmin, vec(i))
610 end do
611 end function vlmin
612
614 subroutine invcol1(a, n)
615 integer, intent(in) :: n
616 real(kind=rp), dimension(n), intent(inout) :: a
617 integer :: i
618
619 do concurrent(i = 1:n)
620 a(i) = 1.0_xp / real(a(i), xp)
621 end do
622
623 end subroutine invcol1
624
626 subroutine invcol3(a, b, c, n)
627 integer, intent(in) :: n
628 real(kind=rp), dimension(n), intent(inout) :: a
629 real(kind=rp), dimension(n), intent(in) :: b, c
630 integer :: i
631
632 do concurrent(i = 1:n)
633 a(i) = real(b(i), xp) / c(i)
634 end do
635
636 end subroutine invcol3
637
639 subroutine invers2(a, b, n)
640 integer, intent(in) :: n
641 real(kind=rp), dimension(n), intent(inout) :: a
642 real(kind=rp), dimension(n), intent(in) :: b
643 integer :: i
644
645 do concurrent(i = 1:n)
646 a(i) = 1.0_xp / real(b(i), xp)
647 end do
648
649 end subroutine invers2
650
653 subroutine vcross(u1, u2, u3, v1, v2, v3, w1, w2, w3, n)
654 integer, intent(in) :: n
655 real(kind=rp), dimension(n), intent(in) :: v1, v2, v3
656 real(kind=rp), dimension(n), intent(in) :: w1, w2, w3
657 real(kind=rp), dimension(n), intent(out) :: u1, u2, u3
658 integer :: i
659
660 do concurrent(i = 1:n)
661 u1(i) = v2(i)*w3(i) - v3(i)*w2(i)
662 u2(i) = v3(i)*w1(i) - v1(i)*w3(i)
663 u3(i) = v1(i)*w2(i) - v2(i)*w1(i)
664 end do
665
666 end subroutine vcross
667
670 subroutine vdot2(dot, u1, u2, v1, v2, n)
671 integer, intent(in) :: n
672 real(kind=rp), dimension(n), intent(in) :: u1, u2
673 real(kind=rp), dimension(n), intent(in) :: v1, v2
674 real(kind=rp), dimension(n), intent(out) :: dot
675 integer :: i
676 do concurrent(i = 1:n)
677 dot(i) = u1(i)*v1(i) + u2(i)*v2(i)
678 end do
679
680 end subroutine vdot2
681
684 subroutine vdot3(dot, u1, u2, u3, v1, v2, v3, n)
685 integer, intent(in) :: n
686 real(kind=rp), dimension(n), intent(in) :: u1, u2, u3
687 real(kind=rp), dimension(n), intent(in) :: v1, v2, v3
688 real(kind=rp), dimension(n), intent(out) :: dot
689 integer :: i
690
691 do concurrent(i = 1:n)
692 dot(i) = u1(i)*v1(i) + u2(i)*v2(i) + u3(i)*v3(i)
693 end do
694
695 end subroutine vdot3
696
698 function vlsc3(u, v, w, n) result(s)
699 integer, intent(in) :: n
700 real(kind=rp), dimension(n), intent(in) :: u, v, w
701 real(kind=rp) :: s
702 integer :: i
703
704 s = 0.0_rp
705 do i = 1, n
706 s = s + u(i)*v(i)*w(i)
707 end do
708
709 end function vlsc3
710
712 function vlsc2(u, v, n) result(s)
713 integer, intent(in) :: n
714 real(kind=rp), dimension(n), intent(in) :: u, v
715 real(kind=rp) :: s
716 integer :: i
717
718 s = 0.0_rp
719 do i = 1, n
720 s = s + u(i)*v(i)
721 end do
722
723 end function vlsc2
724
726 subroutine add2(a, b, n)
727 integer, intent(in) :: n
728 real(kind=rp), dimension(n), intent(inout) :: a
729 real(kind=rp), dimension(n), intent(in) :: b
730 integer :: i
731
732 do concurrent(i = 1:n)
733 a(i) = a(i) + b(i)
734 end do
735
736 end subroutine add2
737
739 subroutine add3(a, b, c, n)
740 integer, intent(in) :: n
741 real(kind=rp), dimension(n), intent(inout) :: a
742 real(kind=rp), dimension(n), intent(in) :: b
743 real(kind=rp), dimension(n), intent(in) :: c
744 integer :: i
745
746 do concurrent(i = 1:n)
747 a(i) = b(i) + c(i)
748 end do
749
750 end subroutine add3
751
753 subroutine add4(a, b, c, d, n)
754 integer, intent(in) :: n
755 real(kind=rp), dimension(n), intent(out) :: a
756 real(kind=rp), dimension(n), intent(in) :: d
757 real(kind=rp), dimension(n), intent(in) :: c
758 real(kind=rp), dimension(n), intent(in) :: b
759 integer :: i
760
761 do concurrent(i = 1:n)
762 a(i) = b(i) + c(i) + d(i)
763 end do
764
765 end subroutine add4
766
768 subroutine sub2(a, b, n)
769 integer, intent(in) :: n
770 real(kind=rp), dimension(n), intent(inout) :: a
771 real(kind=rp), dimension(n), intent(in) :: b
772 integer :: i
773
774 do concurrent(i = 1:n)
775 a(i) = a(i) - b(i)
776 end do
777
778 end subroutine sub2
779
781 subroutine sub3(a, b, c, n)
782 integer, intent(in) :: n
783 real(kind=rp), dimension(n), intent(inout) :: a
784 real(kind=rp), dimension(n), intent(in) :: b
785 real(kind=rp), dimension(n), intent(in) :: c
786 integer :: i
787
788 do concurrent(i = 1:n)
789 a(i) = b(i) - c(i)
790 end do
791
792 end subroutine sub3
793
794
797 subroutine add2s1(a, b, c1, n)
798 integer, intent(in) :: n
799 real(kind=rp), dimension(n), intent(inout) :: a
800 real(kind=rp), dimension(n), intent(in) :: b
801 real(kind=rp), intent(in) :: c1
802 integer :: i
803
804 do concurrent(i = 1:n)
805 a(i) = c1 * a(i) + b(i)
806 end do
807
808 end subroutine add2s1
809
812 subroutine add2s2(a, b, c1, n)
813 integer, intent(in) :: n
814 real(kind=rp), dimension(n), intent(inout) :: a
815 real(kind=rp), dimension(n), intent(in) :: b
816 real(kind=rp), intent(in) :: c1
817 integer :: i
818
819 do concurrent(i = 1:n)
820 a(i) = a(i) + c1 * b(i)
821 end do
822
823 end subroutine add2s2
824
826 subroutine addsqr2s2(a, b, c1, n)
827 integer, intent(in) :: n
828 real(kind=rp), dimension(n), intent(inout) :: a
829 real(kind=rp), dimension(n), intent(in) :: b
830 real(kind=rp), intent(in) :: c1
831 integer :: i
832
833 do concurrent(i = 1:n)
834 a(i) = a(i) + c1 * ( b(i) * b(i) )
835 end do
836
837 end subroutine addsqr2s2
838
840 subroutine invcol2(a, b, n)
841 integer, intent(in) :: n
842 real(kind=rp), dimension(n), intent(inout) :: a
843 real(kind=rp), dimension(n), intent(in) :: b
844 integer :: i
845
846 do concurrent(i = 1:n)
847 a(i) = real(a(i), xp) / b(i)
848 end do
849
850 end subroutine invcol2
851
852
854 subroutine col2(a, b, n)
855 integer, intent(in) :: n
856 real(kind=rp), dimension(n), intent(inout) :: a
857 real(kind=rp), dimension(n), intent(in) :: b
858 integer :: i
859
860 do concurrent(i = 1:n)
861 a(i) = a(i) * b(i)
862 end do
863
864 end subroutine col2
865
867 subroutine col3(a, b, c, n)
868 integer, intent(in) :: n
869 real(kind=rp), dimension(n), intent(inout) :: a
870 real(kind=rp), dimension(n), intent(in) :: b
871 real(kind=rp), dimension(n), intent(in) :: c
872 integer :: i
873
874 do concurrent(i = 1:n)
875 a(i) = b(i) * c(i)
876 end do
877
878 end subroutine col3
879
881 subroutine subcol3(a, b, c, n)
882 integer, intent(in) :: n
883 real(kind=rp), dimension(n), intent(inout) :: a
884 real(kind=rp), dimension(n), intent(in) :: b
885 real(kind=rp), dimension(n), intent(in) :: c
886 integer :: i
887
888 do concurrent(i = 1:n)
889 a(i) = a(i) - b(i) * c(i)
890 end do
891
892 end subroutine subcol3
893
895 subroutine add3s2(a, b, c, c1, c2 ,n)
896 integer, intent(in) :: n
897 real(kind=rp), dimension(n), intent(inout) :: a
898 real(kind=rp), dimension(n), intent(in) :: b
899 real(kind=rp), dimension(n), intent(in) :: c
900 real(kind=rp), intent(in) :: c1, c2
901 integer :: i
902
903 do concurrent(i = 1:n)
904 a(i) = c1 * b(i) + c2 * c(i)
905 end do
906
907 end subroutine add3s2
908
910 subroutine add4s3(a, b, c, d, c1, c2, c3, n)
911 integer, intent(in) :: n
912 real(kind=rp), dimension(n), intent(inout) :: a
913 real(kind=rp), dimension(n), intent(in) :: b
914 real(kind=rp), dimension(n), intent(in) :: c
915 real(kind=rp), dimension(n), intent(in) :: d
916 real(kind=rp), intent(in) :: c1, c2, c3
917 integer :: i
918
919 do concurrent(i = 1:n)
920 a(i) = c1 * b(i) + c2 * c(i) + c3 * d(i)
921 end do
922
923 end subroutine add4s3
924
926 subroutine add5s4(a, b, c, d, e, c1, c2, c3, c4, n)
927 integer, intent(in) :: n
928 real(kind=rp), dimension(n), intent(inout) :: a
929 real(kind=rp), dimension(n), intent(in) :: b
930 real(kind=rp), dimension(n), intent(in) :: c
931 real(kind=rp), dimension(n), intent(in) :: d
932 real(kind=rp), dimension(n), intent(in) :: e
933 real(kind=rp), intent(in) :: c1, c2, c3, c4
934 integer :: i
935
936 do concurrent(i = 1:n)
937 a(i) = a(i) + c1 * b(i) + c2 * c(i) + c3 * d(i) + c4 * e(i)
938 end do
939
940 end subroutine add5s4
941
943 subroutine subcol4(a, b, c, d, n)
944 integer, intent(in) :: n
945 real(kind=rp), dimension(n), intent(inout) :: a
946 real(kind=rp), dimension(n), intent(in) :: b
947 real(kind=rp), dimension(n), intent(in) :: c
948 real(kind=rp), dimension(n), intent(in) :: d
949 integer :: i
950
951 do concurrent(i = 1:n)
952 a(i) = a(i) - b(i) * c(i) * d(i)
953 end do
954
955 end subroutine subcol4
956
958 subroutine addcol3(a, b, c, n)
959 integer, intent(in) :: n
960 real(kind=rp), dimension(n), intent(inout) :: a
961 real(kind=rp), dimension(n), intent(in) :: b
962 real(kind=rp), dimension(n), intent(in) :: c
963 integer :: i
964
965 do concurrent(i = 1:n)
966 a(i) = a(i) + b(i) * c(i)
967 end do
968
969 end subroutine addcol3
970
972 subroutine addcol4(a, b, c, d, n)
973 integer, intent(in) :: n
974 real(kind=rp), dimension(n), intent(inout) :: a
975 real(kind=rp), dimension(n), intent(in) :: b
976 real(kind=rp), dimension(n), intent(in) :: c
977 real(kind=rp), dimension(n), intent(in) :: d
978 integer :: i
979
980 do concurrent(i = 1:n)
981 a(i) = a(i) + b(i) * c(i) * d(i)
982 end do
983
984 end subroutine addcol4
985
987 subroutine addcol3s2(a, b, c, s, n)
988 integer, intent(in) :: n
989 real(kind=rp), dimension(n), intent(inout) :: a
990 real(kind=rp), dimension(n), intent(in) :: b
991 real(kind=rp), dimension(n), intent(in) :: c
992 real(kind=rp), intent(in) :: s
993 integer :: i
994
995 do concurrent(i = 1:n)
996 a(i) = a(i) + s * b(i) * c(i)
997 end do
998
999 end subroutine addcol3s2
1000
1002 subroutine ascol5(a, b, c, d, e, n)
1003 integer, intent(in) :: n
1004 real(kind=rp), dimension(n), intent(inout) :: a
1005 real(kind=rp), dimension(n), intent(in) :: b
1006 real(kind=rp), dimension(n), intent(in) :: c
1007 real(kind=rp), dimension(n), intent(in) :: d
1008 real(kind=rp), dimension(n), intent(in) :: e
1009 integer :: i
1010
1011 do concurrent(i = 1:n)
1012 a(i) = b(i)*c(i) - d(i)*e(i)
1013 end do
1014
1015 end subroutine ascol5
1016
1018 subroutine p_update(a, b, c, c1, c2, n)
1019 integer, intent(in) :: n
1020 real(kind=rp), dimension(n), intent(inout) :: a
1021 real(kind=rp), dimension(n), intent(in) :: b
1022 real(kind=rp), dimension(n), intent(in) :: c
1023 real(kind=rp), intent(in) :: c1, c2
1024 integer :: i
1025
1026 do concurrent(i = 1:n)
1027 a(i) = b(i) + c1*(a(i)-c2*c(i))
1028 end do
1029
1030 end subroutine p_update
1031
1033 subroutine x_update(a, b, c, c1, c2, n)
1034 integer, intent(in) :: n
1035 real(kind=rp), dimension(n), intent(inout) :: a
1036 real(kind=rp), dimension(n), intent(in) :: b
1037 real(kind=rp), dimension(n), intent(in) :: c
1038 real(kind=rp), intent(in) :: c1, c2
1039 integer :: i
1040
1041 do concurrent(i = 1:n)
1042 a(i) = a(i) + c1*b(i)+c2*c(i)
1043 end do
1044
1045 end subroutine x_update
1046
1048 function glsc2(a, b, n)
1049 integer, intent(in) :: n
1050 real(kind=rp), dimension(n), intent(in) :: a
1051 real(kind=rp), dimension(n), intent(in) :: b
1052 real(kind=rp) :: glsc2
1053 real(kind=xp) :: tmp
1054 integer :: i, ierr
1055
1056 tmp = 0.0_xp
1057 do i = 1, n
1058 tmp = tmp + a(i) * b(i)
1059 end do
1060
1061 call mpi_allreduce(mpi_in_place, tmp, 1, &
1062 mpi_extra_precision, mpi_sum, neko_comm, ierr)
1063 glsc2 = tmp
1064 end function glsc2
1065
1067 function glsc3(a, b, c, n)
1068 integer, intent(in) :: n
1069 real(kind=rp), dimension(n), intent(in) :: a
1070 real(kind=rp), dimension(n), intent(in) :: b
1071 real(kind=rp), dimension(n), intent(in) :: c
1072 real(kind=rp) :: glsc3
1073 real(kind=xp) :: tmp
1074 integer :: i, ierr
1075
1076 tmp = 0.0_xp
1077 do i = 1, n
1078 tmp = tmp + a(i) * b(i) * c(i)
1079 end do
1080
1081 call mpi_allreduce(mpi_in_place, tmp, 1, &
1082 mpi_extra_precision, mpi_sum, neko_comm, ierr)
1083 glsc3 = tmp
1084
1085 end function glsc3
1086 function glsc4(a, b, c, d, n)
1087 integer, intent(in) :: n
1088 real(kind=rp), dimension(n), intent(in) :: a
1089 real(kind=rp), dimension(n), intent(in) :: b
1090 real(kind=rp), dimension(n), intent(in) :: c
1091 real(kind=rp), dimension(n), intent(in) :: d
1092 real(kind=rp) :: glsc4
1093 real(kind=xp) :: tmp
1094 integer :: i, ierr
1095
1096 tmp = 0.0_xp
1097 do i = 1, n
1098 tmp = tmp + a(i) * b(i) * c(i) * d(i)
1099 end do
1100
1101 call mpi_allreduce(mpi_in_place, tmp, 1, &
1102 mpi_extra_precision, mpi_sum, neko_comm, ierr)
1103 glsc4 = tmp
1104
1105 end function glsc4
1106
1109 function glsubnorm(a, b, n)
1110 integer, intent(in) :: n
1111 real(kind=rp), dimension(n), intent(in) :: a
1112 real(kind=rp), dimension(n), intent(in) :: b
1113 real(kind=rp) :: glsubnorm
1114 real(kind=xp) :: tmp
1115 integer :: i, ierr
1116
1117 tmp = 0.0_xp
1118 do i = 1, n
1119 tmp = tmp + (a(i) - b(i))**2
1120 end do
1121
1122 call mpi_allreduce(mpi_in_place, tmp, 1, &
1123 mpi_extra_precision, mpi_sum, neko_comm, ierr)
1124 glsubnorm = sqrt(tmp)
1125
1126 end function glsubnorm
1127
1133 subroutine sortrp(a, ind, n)
1134 integer, intent(in) :: n
1135 real(kind=rp), intent(inout) :: a(n)
1136 integer, intent(out) :: ind(n)
1137 real(kind=rp) :: aa
1138 integer :: j, ir, i, ii, l
1139
1140 do j = 1, n
1141 ind(j) = j
1142 end do
1143
1144 if (n .le. 1) return
1145
1146
1147 l = n/2+1
1148 ir = n
1149 do while (.true.)
1150 if (l .gt. 1) then
1151 l = l-1
1152 aa = a(l)
1153 ii = ind(l)
1154 else
1155 aa = a(ir)
1156 ii = ind(ir)
1157 a(ir) = a(1)
1158 ind(ir) = ind(1)
1159 ir = ir - 1
1160 if (ir .eq. 1) then
1161 a(1) = aa
1162 ind(1) = ii
1163 return
1164 end if
1165 end if
1166 i = l
1167 j = l+l
1168 do while (j .le. ir)
1169 if (j .lt. ir) then
1170 if ( a(j) .lt. a(j+1) ) j = j + 1
1171 end if
1172 if (aa .lt. a(j)) then
1173 a(i) = a(j)
1174 ind(i) = ind(j)
1175 i = j
1176 j = j+j
1177 else
1178 j = ir+1
1179 end if
1180 end do
1181 a(i) = aa
1182 ind(i) = ii
1183 end do
1184 end subroutine sortrp
1185
1191 subroutine sorti4(a, ind, n)
1192 integer, intent(in) :: n
1193 integer(i4), intent(inout) :: a(n)
1194 integer, intent(out) :: ind(n)
1195 integer(i4) :: aa
1196 integer :: j, ir, i, ii, l
1197
1198 do j = 1, n
1199 ind(j) = j
1200 end do
1201
1202 if (n .le. 1) return
1203
1204 l = n/2+1
1205 ir = n
1206 do while (.true.)
1207 if (l .gt. 1) then
1208 l = l - 1
1209 aa = a(l)
1210 ii = ind(l)
1211 else
1212 aa = a(ir)
1213 ii = ind(ir)
1214 a(ir) = a( 1)
1215 ind(ir) = ind( 1)
1216 ir = ir - 1
1217 if (ir .eq. 1) then
1218 a(1) = aa
1219 ind(1) = ii
1220 return
1221 end if
1222 end if
1223 i = l
1224 j = l + l
1225 do while (j .le. ir)
1226 if (j .lt. ir) then
1227 if ( a(j) .lt. a(j + 1) ) j = j + 1
1228 end if
1229 if (aa .lt. a(j)) then
1230 a(i) = a(j)
1231 ind(i) = ind(j)
1232 i = j
1233 j = j + j
1234 else
1235 j = ir + 1
1236 end if
1237 end do
1238 a(i) = aa
1239 ind(i) = ii
1240 end do
1241 end subroutine sorti4
1242
1247 subroutine swapdp(b, ind, n)
1248 integer, intent(in) :: n
1249 real(kind=rp), intent(inout) :: b(n)
1250 integer, intent(in) :: ind(n)
1251 real(kind=rp) :: temp(n)
1252 integer :: i, jj
1253
1254 do i = 1, n
1255 temp(i) = b(i)
1256 end do
1257 do i = 1, n
1258 jj = ind(i)
1259 b(i) = temp(jj)
1260 end do
1261 end subroutine swapdp
1262
1267 subroutine swapi4(b, ind, n)
1268 integer, intent(in) :: n
1269 integer(i4), intent(inout) :: b(n)
1270 integer, intent(in) :: ind(n)
1271 integer(i4) :: temp(n)
1272 integer :: i, jj
1273
1274 do i = 1, n
1275 temp(i) = b(i)
1276 end do
1277 do i = 1, n
1278 jj = ind(i)
1279 b(i) = temp(jj)
1280 end do
1281 end subroutine swapi4
1282
1287 subroutine reorddp(b, ind, n)
1288 integer, intent(in) :: n
1289 real(kind=rp), intent(inout) :: b(n)
1290 integer, intent(in) :: ind(n)
1291 real(kind=rp) :: temp(n)
1292 integer :: i, jj
1293
1294 do i = 1, n
1295 temp(i) = b(i)
1296 end do
1297 do i = 1, n
1298 jj = ind(i)
1299 b(jj) = temp(i)
1300 end do
1301 end subroutine reorddp
1302
1307 subroutine reordi4(b, ind, n)
1308 integer, intent(in) :: n
1309 integer(i4), intent(inout) :: b(n)
1310 integer, intent(in) :: ind(n)
1311 integer(i4) :: temp(n)
1312 integer :: i, jj
1313
1314 do i = 1, n
1315 temp(i) = b(i)
1316 end do
1317 do i = 1, n
1318 jj = ind(i)
1319 b(jj) = temp(i)
1320 end do
1321 end subroutine reordi4
1322
1327 subroutine flipvdp(b, ind, n)
1328 integer, intent(in) :: n
1329 real(kind=rp), intent(inout) :: b(n)
1330 integer, intent(inout) :: ind(n)
1331 real(kind=rp) :: temp(n)
1332 integer :: tempind(n)
1333 integer :: i, jj
1334
1335 do i = 1, n
1336 jj = n+1-i
1337 temp(jj) = b(i)
1338 tempind(jj) = ind(i)
1339 end do
1340 do i = 1,n
1341 b(i) = temp(i)
1342 ind(i) = tempind(i)
1343 end do
1344 end subroutine flipvdp
1345
1350 subroutine flipvi4(b, ind, n)
1351 integer, intent(in) :: n
1352 integer(i4), intent(inout) :: b(n)
1353 integer, intent(inout) :: ind(n)
1354 integer(i4) :: temp(n)
1355 integer :: tempind(n)
1356 integer :: i, jj
1357
1358 do i = 1, n
1359 jj = n+1-i
1360 temp(jj) = b(i)
1361 tempind(jj) = ind(i)
1362 end do
1363 do i = 1,n
1364 b(i) = temp(i)
1365 ind(i) = tempind(i)
1366 end do
1367 end subroutine flipvi4
1368
1372 subroutine absval(a, n)
1373 integer, intent(in) :: n
1374 real(kind=rp), dimension(n), intent(inout) :: a
1375 integer :: i
1376 do i = 1, n
1377 a(i) = abs(a(i))
1378 end do
1379 end subroutine absval
1380
1381 ! ========================================================================== !
1382 ! Point-wise operations
1383
1385 subroutine pwmax2(a, b, n)
1386 integer, intent(in) :: n
1387 real(kind=rp), dimension(n), intent(inout) :: a
1388 real(kind=rp), dimension(n), intent(in) :: b
1389 integer :: i
1390
1391 do i = 1, n
1392 a(i) = max(a(i), b(i))
1393 end do
1394 end subroutine pwmax2
1395
1397 subroutine pwmax3(a, b, c, n)
1398 integer, intent(in) :: n
1399 real(kind=rp), dimension(n), intent(inout) :: a
1400 real(kind=rp), dimension(n), intent(in) :: b, c
1401 integer :: i
1402
1403 do i = 1, n
1404 a(i) = max(b(i), c(i))
1405 end do
1406 end subroutine pwmax3
1407
1409 subroutine cpwmax2(a, b, n)
1410 integer, intent(in) :: n
1411 real(kind=rp), dimension(n), intent(inout) :: a
1412 real(kind=rp), intent(in) :: b
1413 integer :: i
1414
1415 do i = 1, n
1416 a(i) = max(a(i), b)
1417 end do
1418 end subroutine cpwmax2
1419
1421 subroutine cpwmax3(a, b, c, n)
1422 integer, intent(in) :: n
1423 real(kind=rp), dimension(n), intent(inout) :: a
1424 real(kind=rp), dimension(n), intent(in) :: b
1425 real(kind=rp), intent(in) :: c
1426 integer :: i
1427
1428 do i = 1, n
1429 a(i) = max(b(i), c)
1430 end do
1431 end subroutine cpwmax3
1432
1434 subroutine pwmin2(a, b, n)
1435 integer, intent(in) :: n
1436 real(kind=rp), dimension(n), intent(inout) :: a
1437 real(kind=rp), dimension(n), intent(in) :: b
1438 integer :: i
1439
1440 do i = 1, n
1441 a(i) = min(a(i), b(i))
1442 end do
1443 end subroutine pwmin2
1444
1446 subroutine pwmin3(a, b, c, n)
1447 integer, intent(in) :: n
1448 real(kind=rp), dimension(n), intent(inout) :: a
1449 real(kind=rp), dimension(n), intent(in) :: b, c
1450 integer :: i
1451
1452 do i = 1, n
1453 a(i) = min(b(i), c(i))
1454 end do
1455 end subroutine pwmin3
1456
1458 subroutine cpwmin2(a, b, n)
1459 integer, intent(in) :: n
1460 real(kind=rp), dimension(n), intent(inout) :: a
1461 real(kind=rp), intent(in) :: b
1462 integer :: i
1463
1464 do i = 1, n
1465 a(i) = min(a(i), b)
1466 end do
1467 end subroutine cpwmin2
1468
1470 subroutine cpwmin3(a, b, c, n)
1471 integer, intent(in) :: n
1472 real(kind=rp), dimension(n), intent(inout) :: a
1473 real(kind=rp), dimension(n), intent(in) :: b
1474 real(kind=rp), intent(in) :: c
1475 integer :: i
1476
1477 do i = 1, n
1478 a(i) = min(b(i), c)
1479 end do
1480 end subroutine cpwmin3
1481
1482 ! M33INV and M44INV by David G. Simpson pure function version from
1483 ! https://fortranwiki.org/fortran/show/Matrix+inversion
1484 ! Invert 3x3 matrix
1485 function matinv39(a11, a12, a13, a21, a22, a23, a31, a32, a33) &
1486 result(b)
1487 real(kind=rp), intent(in) :: a11, a12, a13, a21, a22, a23, a31, a32, a33
1488 real(xp) :: a(3,3) !! Matrix
1489 real(rp) :: b(3,3) !! Inverse matrix
1490 a(1,1) = a11
1491 a(1,2) = a12
1492 a(1,3) = a13
1493 a(2,1) = a21
1494 a(2,2) = a22
1495 a(2,3) = a23
1496 a(3,1) = a31
1497 a(3,2) = a32
1498 a(3,3) = a33
1499 b = matinv3(a)
1500 end function matinv39
1501
1506 function matinv3(A) result(B)
1507 !! Performs a direct calculation of the inverse of a 3×3 matrix.
1508 real(kind=xp), intent(in) :: a(3,3) !! Matrix
1509 real(kind=xp) :: b(3,3) !! Inverse matrix
1510 real(kind=xp) :: detinv
1511
1512 ! Calculate the inverse determinant of the matrix
1513 ! first index x,y,z, second r, s, t
1514 detinv = 1.0_xp / real(a(1,1)*a(2,2)*a(3,3) - a(1,1)*a(2,3)*a(3,2) &
1515 - a(1,2)*a(2,1)*a(3,3) + a(1,2)*a(2,3)*a(3,1)&
1516 + a(1,3)*a(2,1)*a(3,2) - a(1,3)*a(2,2)*a(3,1), xp)
1517 ! Calculate the inverse of the matrix
1518 ! first index r, s, t, second x, y, z
1519 b(1,1) = +detinv * (a(2,2)*a(3,3) - a(2,3)*a(3,2))
1520 b(2,1) = -detinv * (a(2,1)*a(3,3) - a(2,3)*a(3,1))
1521 b(3,1) = +detinv * (a(2,1)*a(3,2) - a(2,2)*a(3,1))
1522 b(1,2) = -detinv * (a(1,2)*a(3,3) - a(1,3)*a(3,2))
1523 b(2,2) = +detinv * (a(1,1)*a(3,3) - a(1,3)*a(3,1))
1524 b(3,2) = -detinv * (a(1,1)*a(3,2) - a(1,2)*a(3,1))
1525 b(1,3) = +detinv * (a(1,2)*a(2,3) - a(1,3)*a(2,2))
1526 b(2,3) = -detinv * (a(1,1)*a(2,3) - a(1,3)*a(2,1))
1527 b(3,3) = +detinv * (a(1,1)*a(2,2) - a(1,2)*a(2,1))
1528 end function matinv3
1529
1532 function math_stepf(x) result(val)
1533 real(kind=rp), intent(in) :: x
1534 real(kind=rp) :: val
1535 real(kind=rp), parameter :: xdmin = 0.0001_rp
1536 real(kind=rp), parameter :: xdmax = 0.9999_rp
1537 real(kind=rp) :: g
1538
1539 if (x <= xdmin) then
1540 ! Below the lower bound, the function is 0
1541 val = 0.0_rp
1542 else if (x >= xdmax) then
1543 ! Above the upper bound, the function is 1
1544 val = 1.0_rp
1545 else
1546 ! g(x) = 1/(x-1) + 1/x
1547 g = (1.0_rp / (x - 1.0_rp)) + (1.0_rp / x)
1548
1549 ! The sigmoid: S(x) = 1 / (1 + exp(g))
1550 val = 1.0_rp / (1.0_rp + exp(g))
1551 end if
1552 end function math_stepf
1553
1555 function math_dstepf(x) result(val)
1556 real(kind=rp), intent(in) :: x
1557 real(kind=rp) :: val
1558 real(kind=rp), parameter :: xdmin = 0.0001_rp
1559 real(kind=rp), parameter :: xdmax = 0.9999_rp
1560 real(kind=rp) :: arg, g, dg, s_val
1561
1562 if (x <= xdmin .or. x >= xdmax) then
1563 val = 0.0_rp
1564 else
1565 ! The step function is S(x) = 1 / (1 + exp(g(x)))
1566 ! where g(x) = 1/(x-1) + 1/x
1567 ! S'(x) = -S(x) * (1 - S(x)) * g'(x)
1568
1569 g = (1.0_rp / (x - 1.0_rp)) + (1.0_rp / x)
1570
1571 ! Derivative of g(x)
1572 dg = -(1.0_rp / ((x - 1.0_rp)**2)) - (1.0_rp / (x**2))
1573
1574 ! Recompute S(x) locally
1575 s_val = 1.0_rp / (1.0_rp + exp(g))
1576
1577 val = -s_val * (1.0_rp - s_val) * dg
1578 end if
1579 end function math_dstepf
1580
1581end module math
double real
Definition comm.F90:1
type(mpi_datatype), public mpi_real_precision
MPI type for working precision of REAL types.
Definition comm.F90:51
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:43
type(mpi_datatype), public mpi_extra_precision
Definition comm.F90:52
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:412
subroutine, public cmult2(a, b, c, n)
Multiplication by constant c .
Definition math.f90:424
subroutine, public row_zero(a, m, n, e)
Sets row e to 0 in matrix a.
Definition math.f90:228
subroutine, public invcol2(a, b, n)
Vector division .
Definition math.f90:841
real(kind=rp) function, public vlsc2(u, v, n)
Compute multiplication sum .
Definition math.f90:713
pure logical function, public dabscmp(x, y, tol)
Return double precision absolute comparison .
Definition math.f90:133
real(kind=rp), parameter, public pi
Definition math.f90:75
pure logical function qabscmp(x, y, tol)
Return double precision absolute comparison .
Definition math.f90:148
real(kind=rp) function, public glsc3(a, b, c, n)
Weighted inner product .
Definition math.f90:1068
subroutine, public ascol5(a, b, c, d, e, n)
Returns .
Definition math.f90:1003
subroutine, public addcol3s2(a, b, c, s, n)
Returns .
Definition math.f90:988
subroutine, public masked_scatter_copy(a, b, mask, n, n_mask)
Scatter a contigous vector to masked positions in a target array .
Definition math.f90:382
subroutine, public invers2(a, b, n)
Compute inverted vector .
Definition math.f90:640
subroutine, public cadd2(a, b, s, n)
Add a scalar to vector .
Definition math.f90:475
real(rp) function, dimension(3, 3), public matinv39(a11, a12, a13, a21, a22, a23, a31, a32, a33)
Definition math.f90:1487
subroutine, public cadd(a, s, n)
Add a scalar to vector .
Definition math.f90:463
subroutine, public masked_copy(a, b, mask, n, n_mask)
Copy a masked vector .
Definition math.f90:290
subroutine reorddp(b, ind, n)
reorder double precision array - inverse of swap
Definition math.f90:1288
subroutine, public addsqr2s2(a, b, c1, n)
Returns .
Definition math.f90:827
real(kind=rp) function, public glsc4(a, b, c, d, n)
Definition math.f90:1087
subroutine, public cdiv2(a, b, c, n)
Division of constant c by elements of a .
Definition math.f90:450
real(kind=rp) function, public math_stepf(x)
Smooth step function S(x) Returns 0 for x <= 0, 1 for x >= 1, and smooth transition in between.
Definition math.f90:1533
subroutine swapdp(b, ind, n)
sort double precision array acording to ind vector
Definition math.f90:1248
subroutine flipvi4(b, ind, n)
Flip single integer vector b and ind.
Definition math.f90:1351
subroutine, public add2s1(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on first argument)
Definition math.f90:798
subroutine, public cpwmin2(a, b, n)
Point-wise minimum of scalar and vector .
Definition math.f90:1459
real(kind=rp) function, public glsc2(a, b, n)
Weighted inner product .
Definition math.f90:1049
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:359
subroutine, public subcol3(a, b, c, n)
Returns .
Definition math.f90:882
subroutine, public rone(a, n)
Set all elements to one.
Definition math.f90:239
subroutine flipvdp(b, ind, n)
Flip double precision vector b and ind.
Definition math.f90:1328
subroutine, public cpwmin3(a, b, c, n)
Point-wise minimum of scalar and vector .
Definition math.f90:1471
subroutine, public pwmax3(a, b, c, n)
Point-wise maximum of two vectors .
Definition math.f90:1398
subroutine, public masked_gather_copy(a, b, mask, n, n_mask)
Gather a masked vector to reduced contigous vector .
Definition math.f90:336
subroutine, public x_update(a, b, c, c1, c2, n)
Returns .
Definition math.f90:1034
subroutine, public add3(a, b, c, n)
Vector addition .
Definition math.f90:740
subroutine swapi4(b, ind, n)
sort single integer array acording to ind vector
Definition math.f90:1268
integer function, public glimin(a, n)
Min of an integer vector of length n.
Definition math.f90:562
real(kind=rp) function, public glsum(a, n)
Sum a vector of length n.
Definition math.f90:500
subroutine, public sub3(a, b, c, n)
Vector subtraction .
Definition math.f90:782
subroutine, public addcol4(a, b, c, d, n)
Returns .
Definition math.f90:973
subroutine, public add2(a, b, n)
Vector addition .
Definition math.f90:727
subroutine, public cfill(a, c, n)
Set all elements to a constant c .
Definition math.f90:488
subroutine, public absval(a, n)
Take the absolute value of an array.
Definition math.f90:1373
subroutine, public invcol3(a, b, c, n)
Invert a vector .
Definition math.f90:627
subroutine, public add3s2(a, b, c, c1, c2, n)
Returns .
Definition math.f90:896
real(kind=xp) function, dimension(3, 3), public matinv3(a)
Performs a direct calculation of the inverse of a 3×3 matrix. M33INV and M44INV by David G....
Definition math.f90:1507
subroutine, public pwmax2(a, b, n)
Point-wise maximum of two vectors .
Definition math.f90:1386
subroutine, public pwmin2(a, b, n)
Point-wise minimum of two vectors .
Definition math.f90:1435
subroutine, public masked_gather_copy_0(a, b, mask, n, n_mask)
Gather a masked vector to reduced contigous vector .
Definition math.f90:313
subroutine, public subcol4(a, b, c, d, n)
Returns .
Definition math.f90:944
subroutine sorti4(a, ind, n)
Heap Sort for single integer arrays.
Definition math.f90:1192
subroutine, public addcol3(a, b, c, n)
Returns .
Definition math.f90:959
subroutine, public invcol1(a, n)
Invert a vector .
Definition math.f90:615
subroutine, public cdiv(a, c, n)
Division of constant c by elements of a .
Definition math.f90:438
real(kind=rp), parameter, public neko_m_ln2
Definition math.f90:72
subroutine, public chsign(a, n)
Change sign of vector .
Definition math.f90:580
subroutine, public cpwmax3(a, b, c, n)
Point-wise maximum of scalar and vector .
Definition math.f90:1422
subroutine, public col2(a, b, n)
Vector multiplication .
Definition math.f90:855
subroutine, public izero(a, n)
Zero an integer vector.
Definition math.f90:217
real(kind=rp) function, public glmax(a, n)
Max of a vector of length n.
Definition math.f90:517
subroutine, public copy(a, b, n)
Copy a vector .
Definition math.f90:250
subroutine, public add4s3(a, b, c, d, c1, c2, c3, n)
Returns .
Definition math.f90:911
subroutine, public add4(a, b, c, d, n)
Vector addition .
Definition math.f90:754
subroutine, public col3(a, b, c, n)
Vector multiplication with 3 vectors .
Definition math.f90:868
subroutine, public add5s4(a, b, c, d, e, c1, c2, c3, c4, n)
Returns .
Definition math.f90:927
real(kind=rp) function, public math_dstepf(x)
Derivative of math_stepf with respect to x: d(stepf)/dx.
Definition math.f90:1556
pure logical function drelcmp(x, y, eps)
Return double precision relative comparison .
Definition math.f90:177
real(kind=rp), parameter, public neko_eps
Machine epsilon .
Definition math.f90:69
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:685
pure logical function, public sabscmp(x, y, tol)
Return single precision absolute comparison .
Definition math.f90:118
pure logical function qrelcmp(x, y, eps)
Return quad precision relative comparison .
Definition math.f90:192
subroutine, public rzero(a, n)
Zero a real vector.
Definition math.f90:206
real(kind=rp) function, public glsubnorm(a, b, n)
Returns the norm of the difference of two vectors .
Definition math.f90:1110
subroutine, public vdot2(dot, u1, u2, v1, v2, n)
Compute a dot product (2-d version) assuming vector components etc.
Definition math.f90:671
subroutine, public cpwmax2(a, b, n)
Point-wise maximum of scalar and vector .
Definition math.f90:1410
real(kind=rp) function, public vlmin(vec, n)
minimun value of a vector of length n
Definition math.f90:603
subroutine, public cfill_mask(a, c, n, mask, n_mask)
Fill a constant to a masked vector. .
Definition math.f90:398
real(kind=rp) function, public vlmax(vec, n)
maximum value of a vector of length n
Definition math.f90:592
integer function, public glimax(a, n)
Max of an integer vector of length n.
Definition math.f90:532
subroutine sortrp(a, ind, n)
Heap Sort for double precision arrays.
Definition math.f90:1134
subroutine, public sub2(a, b, n)
Vector substraction .
Definition math.f90:769
subroutine, public pwmin3(a, b, c, n)
Point-wise minimum of two vectors .
Definition math.f90:1447
subroutine, public add2s2(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on second argument)
Definition math.f90:813
real(kind=rp) function, public glmin(a, n)
Min of a vector of length n.
Definition math.f90:547
subroutine, public masked_copy_0(a, b, mask, n, n_mask)
Copy a masked vector .
Definition math.f90:269
subroutine, public vcross(u1, u2, u3, v1, v2, v3, w1, w2, w3, n)
Compute a cross product assuming vector components etc.
Definition math.f90:654
pure logical function srelcmp(x, y, eps)
Return single precision relative comparison .
Definition math.f90:163
real(kind=rp) function, public vlsc3(u, v, w, n)
Compute multiplication sum .
Definition math.f90:699
subroutine reordi4(b, ind, n)
reorder single integer array - inverse of swap
Definition math.f90:1308
subroutine, public p_update(a, b, c, c1, c2, n)
Returns .
Definition math.f90:1019
integer, parameter, public qp
Definition num_types.f90:10
integer, parameter, public i4
Definition num_types.f90:6
integer, parameter, public xp
Definition num_types.f90:14
integer, parameter, public dp
Definition num_types.f90:9
integer, parameter, public sp
Definition num_types.f90:8
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
#define max(a, b)
Definition tensor.cu:40