Neko 1.99.1
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 interface pwmax
102 module procedure pwmax_vec2, pwmax_vec3, pwmax_scal2, pwmax_scal3
103 end interface pwmax
104
105 interface pwmin
106 module procedure pwmin_vec2, pwmin_vec3, pwmin_sca2, pwmin_sca3
107 end interface pwmin
108
109 public :: abscmp, rzero, izero, row_zero, rone, copy, cmult, cadd, cfill, &
118
119contains
120
122 pure function sabscmp(x, y, tol)
123 real(kind=sp), intent(in) :: x
124 real(kind=sp), intent(in) :: y
125 real(kind=sp), intent(in), optional :: tol
126 logical :: sabscmp
127
128 if (present(tol)) then
129 sabscmp = abs(x - y) .lt. tol
130 else
131 sabscmp = abs(x - y) .lt. neko_eps
132 end if
133
134 end function sabscmp
135
137 pure function dabscmp(x, y, tol)
138 real(kind=dp), intent(in) :: x
139 real(kind=dp), intent(in) :: y
140 real(kind=dp), intent(in), optional :: tol
141 logical :: dabscmp
142
143 if (present(tol)) then
144 dabscmp = abs(x - y) .lt. tol
145 else
146 dabscmp = abs(x - y) .lt. neko_eps
147 end if
148
149 end function dabscmp
150
152 pure function qabscmp(x, y, tol)
153 real(kind=qp), intent(in) :: x
154 real(kind=qp), intent(in) :: y
155 real(kind=qp), intent(in), optional :: tol
156 logical :: qabscmp
157
158 if (present(tol)) then
159 qabscmp = abs(x - y) .lt. tol
160 else
161 qabscmp = abs(x - y) .lt. neko_eps
162 end if
163
164 end function qabscmp
165
167 pure function srelcmp(x, y, eps)
168 real(kind=sp), intent(in) :: x
169 real(kind=sp), intent(in) :: y
170 real(kind=sp), intent(in), optional :: eps
171 logical :: srelcmp
172 if (present(eps)) then
173 srelcmp = abs(x - y) .le. eps*abs(y)
174 else
175 srelcmp = abs(x - y) .le. neko_eps*abs(y)
176 end if
177
178 end function srelcmp
179
181 pure function drelcmp(x, y, eps)
182 real(kind=dp), intent(in) :: x
183 real(kind=dp), intent(in) :: y
184 real(kind=dp), intent(in), optional :: eps
185 logical :: drelcmp
186 if (present(eps)) then
187 drelcmp = abs(x - y) .le. eps*abs(y)
188 else
189 drelcmp = abs(x - y) .le. neko_eps*abs(y)
190 end if
191
192 end function drelcmp
193
194
196 pure function qrelcmp(x, y, eps)
197 real(kind=qp), intent(in) :: x
198 real(kind=qp), intent(in) :: y
199 real(kind=qp), intent(in), optional :: eps
200 logical :: qrelcmp
201 if (present(eps)) then
202 qrelcmp = abs(x - y)/abs(y) .lt. eps
203 else
204 qrelcmp = abs(x - y)/abs(y) .lt. neko_eps
205 end if
206
207 end function qrelcmp
208
210 subroutine rzero(a, n)
211 integer, intent(in) :: n
212 real(kind=rp), dimension(n), intent(inout) :: a
213 integer :: i
214
215 do i = 1, n
216 a(i) = 0.0_rp
217 end do
218 end subroutine rzero
219
221 subroutine izero(a, n)
222 integer, intent(in) :: n
223 integer, dimension(n), intent(inout) :: a
224 integer :: i
225
226 do i = 1, n
227 a(i) = 0
228 end do
229 end subroutine izero
230
232 subroutine row_zero(a, m, n, e)
233 integer, intent(in) :: m, n, e
234 real(kind=rp), intent(inout) :: a(m,n)
235 integer :: j
236
237 do j = 1,n
238 a(e,j) = 0.0_rp
239 end do
240 end subroutine row_zero
241
243 subroutine rone(a, n)
244 integer, intent(in) :: n
245 real(kind=rp), dimension(n), intent(inout) :: a
246 integer :: i
247
248 do i = 1, n
249 a(i) = 1.0_rp
250 end do
251 end subroutine rone
252
254 subroutine copy(a, b, n)
255 integer, intent(in) :: n
256 real(kind=rp), dimension(n), intent(in) :: b
257 real(kind=rp), dimension(n), intent(inout) :: a
258 integer :: i
259
260 do i = 1, n
261 a(i) = b(i)
262 end do
263
264 end subroutine copy
265
273 subroutine masked_copy_0(a, b, mask, n, n_mask)
274 integer, intent(in) :: n, n_mask
275 real(kind=rp), dimension(n), intent(in) :: b
276 real(kind=rp), dimension(n), intent(inout) :: a
277 integer, dimension(0:n_mask) :: mask
278 integer :: i, j
279
280 do i = 1, n_mask
281 j = mask(i)
282 a(j) = b(j)
283 end do
284
285 end subroutine masked_copy_0
286
294 subroutine masked_copy(a, b, mask, n, n_mask)
295 integer, intent(in) :: n, n_mask
296 real(kind=rp), dimension(n), intent(in) :: b
297 real(kind=rp), dimension(n), intent(inout) :: a
298 integer, dimension(n_mask) :: mask
299 integer :: i, j
300
301 do i = 1, n_mask
302 j = mask(i)
303 a(j) = b(j)
304 end do
305
306 end subroutine masked_copy
307
317 subroutine masked_gather_copy_0(a, b, mask, n, n_mask)
318 integer, intent(in) :: n, n_mask
319 real(kind=rp), dimension(n), intent(in) :: b
320 real(kind=rp), dimension(n_mask), intent(inout) :: a
321 integer, dimension(0:n_mask) :: mask
322 integer :: i, j
323
324 do i = 1, n_mask
325 j = mask(i)
326 a(i) = b(j)
327 end do
328
329 end subroutine masked_gather_copy_0
330
340 subroutine masked_gather_copy(a, b, mask, n, n_mask)
341 integer, intent(in) :: n, n_mask
342 real(kind=rp), dimension(n), intent(in) :: b
343 real(kind=rp), dimension(n_mask), intent(inout) :: a
344 integer, dimension(n_mask) :: mask
345 integer :: i, j
346
347 do i = 1, n_mask
348 j = mask(i)
349 a(i) = b(j)
350 end do
351
352 end subroutine masked_gather_copy
353
363 subroutine masked_scatter_copy_0(a, b, mask, n, n_mask)
364 integer, intent(in) :: n, n_mask
365 real(kind=rp), dimension(n_mask), intent(in) :: b
366 real(kind=rp), dimension(n), intent(inout) :: a
367 integer, dimension(0:n_mask) :: mask
368 integer :: i, j
369
370 do i = 1, n_mask
371 j = mask(i)
372 a(j) = b(i)
373 end do
374
375 end subroutine masked_scatter_copy_0
376
386 subroutine masked_scatter_copy(a, b, mask, n, n_mask)
387 integer, intent(in) :: n, n_mask
388 real(kind=rp), dimension(n_mask), intent(in) :: b
389 real(kind=rp), dimension(n), intent(inout) :: a
390 integer, dimension(n_mask) :: mask
391 integer :: i, j
392
393 do i = 1, n_mask
394 j = mask(i)
395 a(j) = b(i)
396 end do
397
398 end subroutine masked_scatter_copy
399
402 subroutine cfill_mask(a, c, n, mask, n_mask)
403 integer, intent(in) :: n, n_mask
404 real(kind=rp), dimension(n), intent(inout) :: a
405 real(kind=rp), intent(in) :: c
406 integer, dimension(n_mask), intent(in) :: mask
407 integer :: i
408
409 do i = 1, n_mask
410 a(mask(i)) = c
411 end do
412
413 end subroutine cfill_mask
414
416 subroutine cmult(a, c, n)
417 integer, intent(in) :: n
418 real(kind=rp), dimension(n), intent(inout) :: a
419 real(kind=rp), intent(in) :: c
420 integer :: i
421
422 do i = 1, n
423 a(i) = c * a(i)
424 end do
425 end subroutine cmult
426
428 subroutine cmult2(a, b, c, n)
429 integer, intent(in) :: n
430 real(kind=rp), dimension(n), intent(inout) :: a
431 real(kind=rp), dimension(n), intent(in) :: b
432 real(kind=rp), intent(in) :: c
433 integer :: i
434
435 do i = 1, n
436 a(i) = c * b(i)
437 end do
438
439 end subroutine cmult2
440
442 subroutine cdiv(a, c, n)
443 integer, intent(in) :: n
444 real(kind=rp), dimension(n), intent(inout) :: a
445 real(kind=rp), intent(in) :: c
446 integer :: i
447
448 do i = 1, n
449 a(i) = c / a(i)
450 end do
451 end subroutine cdiv
452
454 subroutine cdiv2(a, b, c, n)
455 integer, intent(in) :: n
456 real(kind=rp), dimension(n), intent(inout) :: a
457 real(kind=rp), dimension(n), intent(in) :: b
458 real(kind=rp), intent(in) :: c
459 integer :: i
460
461 do i = 1, n
462 a(i) = c / b(i)
463 end do
464 end subroutine cdiv2
465
467 subroutine cadd(a, s, n)
468 integer, intent(in) :: n
469 real(kind=rp), dimension(n), intent(inout) :: a
470 real(kind=rp), intent(in) :: s
471 integer :: i
472
473 do i = 1, n
474 a(i) = a(i) + s
475 end do
476 end subroutine cadd
477
479 subroutine cadd2(a, b, s, n)
480 integer, intent(in) :: n
481 real(kind=rp), dimension(n), intent(inout) :: a
482 real(kind=rp), dimension(n), intent(in) :: b
483 real(kind=rp), intent(in) :: s
484 integer :: i
485
486 do i = 1, n
487 a(i) = b(i) + s
488 end do
489 end subroutine cadd2
490
492 subroutine cfill(a, c, n)
493 integer, intent(in) :: n
494 real(kind=rp), dimension(n), intent(inout) :: a
495 real(kind=rp), intent(in) :: c
496 integer :: i
497
498 do i = 1, n
499 a(i) = c
500 end do
501 end subroutine cfill
502
504 function glsum(a, n)
505 integer, intent(in) :: n
506 real(kind=rp), dimension(n) :: a
507 real(kind=rp) :: glsum
508 real(kind=xp) :: tmp
509 integer :: i, ierr
510 tmp = 0.0_rp
511 do i = 1, n
512 tmp = tmp + a(i)
513 end do
514 call mpi_allreduce(mpi_in_place, tmp, 1, &
515 mpi_extra_precision, mpi_sum, neko_comm, ierr)
516 glsum = tmp
517
518 end function glsum
519
521 function glmax(a, n)
522 integer, intent(in) :: n
523 real(kind=rp), dimension(n) :: a
524 real(kind=rp) :: tmp, glmax
525 integer :: i, ierr
526
527 tmp = -huge(0.0_rp)
528 do i = 1, n
529 tmp = max(tmp,a(i))
530 end do
531 call mpi_allreduce(tmp, glmax, 1, &
532 mpi_real_precision, mpi_max, neko_comm, ierr)
533 end function glmax
534
536 function glimax(a, n)
537 integer, intent(in) :: n
538 integer, dimension(n) :: a
539 integer :: tmp, glimax
540 integer :: i, ierr
541
542 tmp = -huge(0)
543 do i = 1, n
544 tmp = max(tmp,a(i))
545 end do
546 call mpi_allreduce(tmp, glimax, 1, &
547 mpi_integer, mpi_max, neko_comm, ierr)
548 end function glimax
549
551 function glmin(a, n)
552 integer, intent(in) :: n
553 real(kind=rp), dimension(n) :: a
554 real(kind=rp) :: tmp, glmin
555 integer :: i, ierr
556
557 tmp = huge(0.0_rp)
558 do i = 1, n
559 tmp = min(tmp,a(i))
560 end do
561 call mpi_allreduce(tmp, glmin, 1, &
562 mpi_real_precision, mpi_min, neko_comm, ierr)
563 end function glmin
564
566 function glimin(a, n)
567 integer, intent(in) :: n
568 integer, dimension(n) :: a
569 integer :: tmp, glimin
570 integer :: i, ierr
571
572 tmp = huge(0)
573 do i = 1, n
574 tmp = min(tmp,a(i))
575 end do
576 call mpi_allreduce(tmp, glimin, 1, &
577 mpi_integer, mpi_min, neko_comm, ierr)
578 end function glimin
579
580
581
582
584 subroutine chsign(a, n)
585 integer, intent(in) :: n
586 real(kind=rp), dimension(n), intent(inout) :: a
587 integer :: i
588
589 do i = 1, n
590 a(i) = -a(i)
591 end do
592
593 end subroutine chsign
594
596 function vlmax(vec,n) result(tmax)
597 integer :: n, i
598 real(kind=rp), intent(in) :: vec(n)
599 real(kind=rp) :: tmax
600 tmax = real(-99d20, rp)
601 do i = 1, n
602 tmax = max(tmax, vec(i))
603 end do
604 end function vlmax
605
607 function vlmin(vec,n) result(tmin)
608 integer, intent(in) :: n
609 real(kind=rp), intent(in) :: vec(n)
610 real(kind=rp) :: tmin
611 integer :: i
612 tmin = real(99.0e20, rp)
613 do i = 1, n
614 tmin = min(tmin, vec(i))
615 end do
616 end function vlmin
617
619 subroutine invcol1(a, n)
620 integer, intent(in) :: n
621 real(kind=rp), dimension(n), intent(inout) :: a
622 integer :: i
623
624 do i = 1, n
625 a(i) = 1.0_xp / real(a(i),xp)
626 end do
627
628 end subroutine invcol1
629
631 subroutine invcol3(a, b, c, n)
632 integer, intent(in) :: n
633 real(kind=rp), dimension(n), intent(inout) :: a
634 real(kind=rp), dimension(n), intent(in) :: b,c
635 integer :: i
636
637 do i = 1, n
638 a(i) = real(b(i),xp) / c(i)
639 end do
640
641 end subroutine invcol3
642
644 subroutine invers2(a, b, n)
645 integer, intent(in) :: n
646 real(kind=rp), dimension(n), intent(inout) :: a
647 real(kind=rp), dimension(n), intent(in) :: b
648 integer :: i
649
650 do i = 1, n
651 a(i) = 1.0_xp / real(b(i),xp)
652 end do
653
654 end subroutine invers2
655
658 subroutine vcross(u1, u2, u3, v1, v2, v3, w1, w2, w3, n)
659 integer, intent(in) :: n
660 real(kind=rp), dimension(n), intent(in) :: v1, v2, v3
661 real(kind=rp), dimension(n), intent(in) :: w1, w2, w3
662 real(kind=rp), dimension(n), intent(out) :: u1, u2, u3
663 integer :: i
664
665 do i = 1, n
666 u1(i) = v2(i)*w3(i) - v3(i)*w2(i)
667 u2(i) = v3(i)*w1(i) - v1(i)*w3(i)
668 u3(i) = v1(i)*w2(i) - v2(i)*w1(i)
669 end do
670
671 end subroutine vcross
672
675 subroutine vdot2(dot, u1, u2, v1, v2, n)
676 integer, intent(in) :: n
677 real(kind=rp), dimension(n), intent(in) :: u1, u2
678 real(kind=rp), dimension(n), intent(in) :: v1, v2
679 real(kind=rp), dimension(n), intent(out) :: dot
680 integer :: i
681 do i = 1, n
682 dot(i) = u1(i)*v1(i) + u2(i)*v2(i)
683 end do
684
685 end subroutine vdot2
686
689 subroutine vdot3(dot, u1, u2, u3, v1, v2, v3, n)
690 integer, intent(in) :: n
691 real(kind=rp), dimension(n), intent(in) :: u1, u2, u3
692 real(kind=rp), dimension(n), intent(in) :: v1, v2, v3
693 real(kind=rp), dimension(n), intent(out) :: dot
694 integer :: i
695
696 do i = 1, n
697 dot(i) = u1(i)*v1(i) + u2(i)*v2(i) + u3(i)*v3(i)
698 end do
699
700 end subroutine vdot3
701
703 function vlsc3(u, v, w, n) result(s)
704 integer, intent(in) :: n
705 real(kind=rp), dimension(n), intent(in) :: u, v, w
706 real(kind=rp) :: s
707 integer :: i
708
709 s = 0.0_rp
710 do i = 1, n
711 s = s + u(i)*v(i)*w(i)
712 end do
713
714 end function vlsc3
715
717 function vlsc2(u, v, n) result(s)
718 integer, intent(in) :: n
719 real(kind=rp), dimension(n), intent(in) :: u, v
720 real(kind=rp) :: s
721 integer :: i
722
723 s = 0.0_rp
724 do i = 1, n
725 s = s + u(i)*v(i)
726 end do
727
728 end function vlsc2
729
731 subroutine add2(a, b, n)
732 integer, intent(in) :: n
733 real(kind=rp), dimension(n), intent(inout) :: a
734 real(kind=rp), dimension(n), intent(in) :: b
735 integer :: i
736
737 do i = 1, n
738 a(i) = a(i) + b(i)
739 end do
740
741 end subroutine add2
742
744 subroutine add3(a, b, c, n)
745 integer, intent(in) :: n
746 real(kind=rp), dimension(n), intent(inout) :: a
747 real(kind=rp), dimension(n), intent(in) :: b
748 real(kind=rp), dimension(n), intent(in) :: c
749 integer :: i
750
751 do i = 1, n
752 a(i) = b(i) + c(i)
753 end do
754
755 end subroutine add3
756
758 subroutine add4(a, b, c, d, n)
759 integer, intent(in) :: n
760 real(kind=rp), dimension(n), intent(out) :: a
761 real(kind=rp), dimension(n), intent(in) :: d
762 real(kind=rp), dimension(n), intent(in) :: c
763 real(kind=rp), dimension(n), intent(in) :: b
764 integer :: i
765
766 do i = 1, n
767 a(i) = b(i) + c(i) + d(i)
768 end do
769
770 end subroutine add4
771
773 subroutine sub2(a, b, n)
774 integer, intent(in) :: n
775 real(kind=rp), dimension(n), intent(inout) :: a
776 real(kind=rp), dimension(n), intent(inout) :: b
777 integer :: i
778
779 do i = 1, n
780 a(i) = a(i) - b(i)
781 end do
782
783 end subroutine sub2
784
786 subroutine sub3(a, b, c, n)
787 integer, intent(in) :: n
788 real(kind=rp), dimension(n), intent(inout) :: a
789 real(kind=rp), dimension(n), intent(in) :: b
790 real(kind=rp), dimension(n), intent(in) :: c
791 integer :: i
792
793 do i = 1, n
794 a(i) = b(i) - c(i)
795 end do
796
797 end subroutine sub3
798
799
802 subroutine add2s1(a, b, c1, n)
803 integer, intent(in) :: n
804 real(kind=rp), dimension(n), intent(inout) :: a
805 real(kind=rp), dimension(n), intent(inout) :: b
806 real(kind=rp), intent(in) :: c1
807 integer :: i
808
809 do i = 1, n
810 a(i) = c1 * a(i) + b(i)
811 end do
812
813 end subroutine add2s1
814
817 subroutine add2s2(a, b, c1, n)
818 integer, intent(in) :: n
819 real(kind=rp), dimension(n), intent(inout) :: a
820 real(kind=rp), dimension(n), intent(in) :: b
821 real(kind=rp), intent(in) :: c1
822 integer :: i
823
824 do i = 1, n
825 a(i) = a(i) + c1 * b(i)
826 end do
827
828 end subroutine add2s2
829
831 subroutine addsqr2s2(a, b, c1, n)
832 integer, intent(in) :: n
833 real(kind=rp), dimension(n), intent(inout) :: a
834 real(kind=rp), dimension(n), intent(in) :: b
835 real(kind=rp), intent(in) :: c1
836 integer :: i
837
838 do i = 1,n
839 a(i) = a(i) + c1 * ( b(i) * b(i) )
840 end do
841
842 end subroutine addsqr2s2
843
845 subroutine invcol2(a, b, n)
846 integer, intent(in) :: n
847 real(kind=rp), dimension(n), intent(inout) :: a
848 real(kind=rp), dimension(n), intent(in) :: b
849 integer :: i
850
851 do i = 1, n
852 a(i) = real(a(i),xp) /b(i)
853 end do
854
855 end subroutine invcol2
856
857
859 subroutine col2(a, b, n)
860 integer, intent(in) :: n
861 real(kind=rp), dimension(n), intent(inout) :: a
862 real(kind=rp), dimension(n), intent(in) :: b
863 integer :: i
864
865 do i = 1, n
866 a(i) = a(i) * b(i)
867 end do
868
869 end subroutine col2
870
872 subroutine col3(a, b, c, n)
873 integer, intent(in) :: n
874 real(kind=rp), dimension(n), intent(inout) :: a
875 real(kind=rp), dimension(n), intent(in) :: b
876 real(kind=rp), dimension(n), intent(in) :: c
877 integer :: i
878
879 do i = 1, n
880 a(i) = b(i) * c(i)
881 end do
882
883 end subroutine col3
884
886 subroutine subcol3(a, b, c, n)
887 integer, intent(in) :: n
888 real(kind=rp), dimension(n), intent(inout) :: a
889 real(kind=rp), dimension(n), intent(in) :: b
890 real(kind=rp), dimension(n), intent(in) :: c
891 integer :: i
892
893 do i = 1,n
894 a(i) = a(i) - b(i) * c(i)
895 end do
896
897 end subroutine subcol3
898
900 subroutine add3s2(a, b, c, c1, c2 ,n)
901 integer, intent(in) :: n
902 real(kind=rp), dimension(n), intent(inout) :: a
903 real(kind=rp), dimension(n), intent(in) :: b
904 real(kind=rp), dimension(n), intent(in) :: c
905 real(kind=rp), intent(in) :: c1, c2
906 integer :: i
907
908 do i = 1,n
909 a(i) = c1 * b(i) + c2 * c(i)
910 end do
911
912 end subroutine add3s2
913
914
916 subroutine subcol4(a, b, c, d, n)
917 integer, intent(in) :: n
918 real(kind=rp), dimension(n), intent(inout) :: a
919 real(kind=rp), dimension(n), intent(in) :: b
920 real(kind=rp), dimension(n), intent(in) :: c
921 real(kind=rp), dimension(n), intent(in) :: d
922 integer :: i
923
924 do i = 1,n
925 a(i) = a(i) - b(i) * c(i) * d(i)
926 end do
927
928 end subroutine subcol4
929
931 subroutine addcol3(a, b, c, n)
932 integer, intent(in) :: n
933 real(kind=rp), dimension(n), intent(inout) :: a
934 real(kind=rp), dimension(n), intent(in) :: b
935 real(kind=rp), dimension(n), intent(in) :: c
936 integer :: i
937
938 do i = 1,n
939 a(i) = a(i) + b(i) * c(i)
940 end do
941
942 end subroutine addcol3
943
945 subroutine addcol4(a, b, c, d, n)
946 integer, intent(in) :: n
947 real(kind=rp), dimension(n), intent(inout) :: a
948 real(kind=rp), dimension(n), intent(in) :: b
949 real(kind=rp), dimension(n), intent(in) :: c
950 real(kind=rp), dimension(n), intent(in) :: d
951 integer :: i
952
953 do i = 1,n
954 a(i) = a(i) + b(i) * c(i) * d(i)
955 end do
956
957 end subroutine addcol4
958
960 subroutine ascol5(a, b, c, d, e, n)
961 integer, intent(in) :: n
962 real(kind=rp), dimension(n), intent(inout) :: a
963 real(kind=rp), dimension(n), intent(in) :: b
964 real(kind=rp), dimension(n), intent(in) :: c
965 real(kind=rp), dimension(n), intent(in) :: d
966 real(kind=rp), dimension(n), intent(in) :: e
967 integer :: i
968
969 do i = 1,n
970 a(i) = b(i)*c(i)-d(i)*e(i)
971 end do
972
973 end subroutine ascol5
974
976 subroutine p_update(a, b, c, c1, c2, n)
977 integer, intent(in) :: n
978 real(kind=rp), dimension(n), intent(inout) :: a
979 real(kind=rp), dimension(n), intent(in) :: b
980 real(kind=rp), dimension(n), intent(in) :: c
981 real(kind=rp), intent(in) :: c1, c2
982 integer :: i
983
984 do i = 1,n
985 a(i) = b(i) + c1*(a(i)-c2*c(i))
986 end do
987
988 end subroutine p_update
989
991 subroutine x_update(a, b, c, c1, c2, n)
992 integer, intent(in) :: n
993 real(kind=rp), dimension(n), intent(inout) :: a
994 real(kind=rp), dimension(n), intent(in) :: b
995 real(kind=rp), dimension(n), intent(in) :: c
996 real(kind=rp), intent(in) :: c1, c2
997 integer :: i
998
999 do i = 1,n
1000 a(i) = a(i) + c1*b(i)+c2*c(i)
1001 end do
1002
1003 end subroutine x_update
1004
1006 function glsc2(a, b, n)
1007 integer, intent(in) :: n
1008 real(kind=rp), dimension(n), intent(in) :: a
1009 real(kind=rp), dimension(n), intent(in) :: b
1010 real(kind=rp) :: glsc2
1011 real(kind=xp) :: tmp
1012 integer :: i, ierr
1013
1014 tmp = 0.0_xp
1015 do i = 1, n
1016 tmp = tmp + a(i) * b(i)
1017 end do
1018
1019 call mpi_allreduce(mpi_in_place, tmp, 1, &
1020 mpi_extra_precision, mpi_sum, neko_comm, ierr)
1021 glsc2 = tmp
1022 end function glsc2
1023
1025 function glsc3(a, b, c, n)
1026 integer, intent(in) :: n
1027 real(kind=rp), dimension(n), intent(in) :: a
1028 real(kind=rp), dimension(n), intent(in) :: b
1029 real(kind=rp), dimension(n), intent(in) :: c
1030 real(kind=rp) :: glsc3
1031 real(kind=xp) :: tmp
1032 integer :: i, ierr
1033
1034 tmp = 0.0_xp
1035 do i = 1, n
1036 tmp = tmp + a(i) * b(i) * c(i)
1037 end do
1038
1039 call mpi_allreduce(mpi_in_place, tmp, 1, &
1040 mpi_extra_precision, mpi_sum, neko_comm, ierr)
1041 glsc3 = tmp
1042
1043 end function glsc3
1044 function glsc4(a, b, c, d, n)
1045 integer, intent(in) :: n
1046 real(kind=rp), dimension(n), intent(in) :: a
1047 real(kind=rp), dimension(n), intent(in) :: b
1048 real(kind=rp), dimension(n), intent(in) :: c
1049 real(kind=rp), dimension(n), intent(in) :: d
1050 real(kind=rp) :: glsc4
1051 real(kind=xp) :: tmp
1052 integer :: i, ierr
1053
1054 tmp = 0.0_xp
1055 do i = 1, n
1056 tmp = tmp + a(i) * b(i) * c(i) * d(i)
1057 end do
1058
1059 call mpi_allreduce(mpi_in_place, tmp, 1, &
1060 mpi_extra_precision, mpi_sum, neko_comm, ierr)
1061 glsc4 = tmp
1062
1063 end function glsc4
1064
1067 function glsubnorm(a, b, 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) :: glsubnorm
1072 real(kind=xp) :: tmp
1073 integer :: i, ierr
1074
1075 tmp = 0.0_xp
1076 do i = 1, n
1077 tmp = tmp + (a(i) - b(i))**2
1078 end do
1079
1080 call mpi_allreduce(mpi_in_place, tmp, 1, &
1081 mpi_extra_precision, mpi_sum, neko_comm, ierr)
1082 glsubnorm = sqrt(tmp)
1083
1084 end function glsubnorm
1085
1091 subroutine sortrp(a, ind, n)
1092 integer, intent(in) :: n
1093 real(kind=rp), intent(inout) :: a(n)
1094 integer, intent(out) :: ind(n)
1095 real(kind=rp) :: aa
1096 integer :: j, ir, i, ii, l
1097
1098 do j = 1, n
1099 ind(j) = j
1100 end do
1101
1102 if (n .le. 1) return
1103
1104
1105 l = n/2+1
1106 ir = n
1107 do while (.true.)
1108 if (l .gt. 1) then
1109 l = l-1
1110 aa = a(l)
1111 ii = ind(l)
1112 else
1113 aa = a(ir)
1114 ii = ind(ir)
1115 a(ir) = a(1)
1116 ind(ir) = ind(1)
1117 ir = ir - 1
1118 if (ir .eq. 1) then
1119 a(1) = aa
1120 ind(1) = ii
1121 return
1122 end if
1123 end if
1124 i = l
1125 j = l+l
1126 do while (j .le. ir)
1127 if (j .lt. ir) then
1128 if ( a(j) .lt. a(j+1) ) j = j + 1
1129 end if
1130 if (aa .lt. a(j)) then
1131 a(i) = a(j)
1132 ind(i) = ind(j)
1133 i = j
1134 j = j+j
1135 else
1136 j = ir+1
1137 end if
1138 end do
1139 a(i) = aa
1140 ind(i) = ii
1141 end do
1142 end subroutine sortrp
1143
1149 subroutine sorti4(a, ind, n)
1150 integer, intent(in) :: n
1151 integer(i4), intent(inout) :: a(n)
1152 integer, intent(out) :: ind(n)
1153 integer(i4) :: aa
1154 integer :: j, ir, i, ii, l
1155
1156 do j = 1, n
1157 ind(j) = j
1158 end do
1159
1160 if (n .le. 1) return
1161
1162 l = n/2+1
1163 ir = n
1164 do while (.true.)
1165 if (l .gt. 1) then
1166 l = l - 1
1167 aa = a(l)
1168 ii = ind(l)
1169 else
1170 aa = a(ir)
1171 ii = ind(ir)
1172 a(ir) = a( 1)
1173 ind(ir) = ind( 1)
1174 ir = ir - 1
1175 if (ir .eq. 1) then
1176 a(1) = aa
1177 ind(1) = ii
1178 return
1179 end if
1180 end if
1181 i = l
1182 j = l + l
1183 do while (j .le. ir)
1184 if (j .lt. ir) then
1185 if ( a(j) .lt. a(j + 1) ) j = j + 1
1186 end if
1187 if (aa .lt. a(j)) then
1188 a(i) = a(j)
1189 ind(i) = ind(j)
1190 i = j
1191 j = j + j
1192 else
1193 j = ir + 1
1194 end if
1195 end do
1196 a(i) = aa
1197 ind(i) = ii
1198 end do
1199 end subroutine sorti4
1200
1205 subroutine swapdp(b, ind, n)
1206 integer, intent(in) :: n
1207 real(kind=rp), intent(inout) :: b(n)
1208 integer, intent(in) :: ind(n)
1209 real(kind=rp) :: temp(n)
1210 integer :: i, jj
1211
1212 do i = 1, n
1213 temp(i) = b(i)
1214 end do
1215 do i = 1, n
1216 jj = ind(i)
1217 b(i) = temp(jj)
1218 end do
1219 end subroutine swapdp
1220
1225 subroutine swapi4(b, ind, n)
1226 integer, intent(in) :: n
1227 integer(i4), intent(inout) :: b(n)
1228 integer, intent(in) :: ind(n)
1229 integer(i4) :: temp(n)
1230 integer :: i, jj
1231
1232 do i = 1, n
1233 temp(i) = b(i)
1234 end do
1235 do i = 1, n
1236 jj = ind(i)
1237 b(i) = temp(jj)
1238 end do
1239 end subroutine swapi4
1240
1245 subroutine reorddp(b, ind, n)
1246 integer, intent(in) :: n
1247 real(kind=rp), intent(inout) :: b(n)
1248 integer, intent(in) :: ind(n)
1249 real(kind=rp) :: temp(n)
1250 integer :: i, jj
1251
1252 do i = 1, n
1253 temp(i) = b(i)
1254 end do
1255 do i = 1, n
1256 jj = ind(i)
1257 b(jj) = temp(i)
1258 end do
1259 end subroutine reorddp
1260
1265 subroutine reordi4(b, ind, n)
1266 integer, intent(in) :: n
1267 integer(i4), intent(inout) :: b(n)
1268 integer, intent(in) :: ind(n)
1269 integer(i4) :: temp(n)
1270 integer :: i, jj
1271
1272 do i = 1, n
1273 temp(i) = b(i)
1274 end do
1275 do i = 1, n
1276 jj = ind(i)
1277 b(jj) = temp(i)
1278 end do
1279 end subroutine reordi4
1280
1285 subroutine flipvdp(b, ind, n)
1286 integer, intent(in) :: n
1287 real(kind=rp), intent(inout) :: b(n)
1288 integer, intent(inout) :: ind(n)
1289 real(kind=rp) :: temp(n)
1290 integer :: tempind(n)
1291 integer :: i, jj
1292
1293 do i = 1, n
1294 jj = n+1-i
1295 temp(jj) = b(i)
1296 tempind(jj) = ind(i)
1297 end do
1298 do i = 1,n
1299 b(i) = temp(i)
1300 ind(i) = tempind(i)
1301 end do
1302 end subroutine flipvdp
1303
1308 subroutine flipvi4(b, ind, n)
1309 integer, intent(in) :: n
1310 integer(i4), intent(inout) :: b(n)
1311 integer, intent(inout) :: ind(n)
1312 integer(i4) :: temp(n)
1313 integer :: tempind(n)
1314 integer :: i, jj
1315
1316 do i = 1, n
1317 jj = n+1-i
1318 temp(jj) = b(i)
1319 tempind(jj) = ind(i)
1320 end do
1321 do i = 1,n
1322 b(i) = temp(i)
1323 ind(i) = tempind(i)
1324 end do
1325 end subroutine flipvi4
1326
1330 subroutine absval(a, n)
1331 integer, intent(in) :: n
1332 real(kind=rp), dimension(n), intent(inout) :: a
1333 integer :: i
1334 do i = 1, n
1335 a(i) = abs(a(i))
1336 end do
1337 end subroutine absval
1338
1339 ! ========================================================================== !
1340 ! Point-wise operations
1341
1343 subroutine pwmax_vec2(a, b, n)
1344 integer, intent(in) :: n
1345 real(kind=rp), dimension(n), intent(inout) :: a
1346 real(kind=rp), dimension(n), intent(in) :: b
1347 integer :: i
1348
1349 do i = 1, n
1350 a(i) = max(a(i), b(i))
1351 end do
1352 end subroutine pwmax_vec2
1353
1355 subroutine pwmax_vec3(a, b, c, n)
1356 integer, intent(in) :: n
1357 real(kind=rp), dimension(n), intent(inout) :: a
1358 real(kind=rp), dimension(n), intent(in) :: b, c
1359 integer :: i
1360
1361 do i = 1, n
1362 a(i) = max(b(i), c(i))
1363 end do
1364 end subroutine pwmax_vec3
1365
1367 subroutine pwmax_scal2(a, b, n)
1368 integer, intent(in) :: n
1369 real(kind=rp), dimension(n), intent(inout) :: a
1370 real(kind=rp), intent(in) :: b
1371 integer :: i
1372
1373 do i = 1, n
1374 a(i) = max(a(i), b)
1375 end do
1376 end subroutine pwmax_scal2
1377
1379 subroutine pwmax_scal3(a, b, c, n)
1380 integer, intent(in) :: n
1381 real(kind=rp), dimension(n), intent(inout) :: a
1382 real(kind=rp), dimension(n), intent(in) :: b
1383 real(kind=rp), intent(in) :: c
1384 integer :: i
1385
1386 do i = 1, n
1387 a(i) = max(b(i), c)
1388 end do
1389 end subroutine pwmax_scal3
1390
1392 subroutine pwmin_vec2(a, b, n)
1393 integer, intent(in) :: n
1394 real(kind=rp), dimension(n), intent(inout) :: a
1395 real(kind=rp), dimension(n), intent(in) :: b
1396 integer :: i
1397
1398 do i = 1, n
1399 a(i) = min(a(i), b(i))
1400 end do
1401 end subroutine pwmin_vec2
1402
1404 subroutine pwmin_vec3(a, b, c, n)
1405 integer, intent(in) :: n
1406 real(kind=rp), dimension(n), intent(inout) :: a
1407 real(kind=rp), dimension(n), intent(in) :: b, c
1408 integer :: i
1409
1410 do i = 1, n
1411 a(i) = min(b(i), c(i))
1412 end do
1413 end subroutine pwmin_vec3
1414
1416 subroutine pwmin_sca2(a, b, n)
1417 integer, intent(in) :: n
1418 real(kind=rp), dimension(n), intent(inout) :: a
1419 real(kind=rp), intent(in) :: b
1420 integer :: i
1421
1422 do i = 1, n
1423 a(i) = min(a(i), b)
1424 end do
1425 end subroutine pwmin_sca2
1426
1428 subroutine pwmin_sca3(a, b, c, n)
1429 integer, intent(in) :: n
1430 real(kind=rp), dimension(n), intent(inout) :: a
1431 real(kind=rp), dimension(n), intent(in) :: b
1432 real(kind=rp), intent(in) :: c
1433 integer :: i
1434
1435 do i = 1, n
1436 a(i) = min(b(i), c)
1437 end do
1438 end subroutine pwmin_sca3
1439
1440end 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:50
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:42
type(mpi_datatype), public mpi_extra_precision
Definition comm.F90:51
Object for handling masks in Neko.
Definition mask.f90:34
Definition math.f90:60
subroutine, public cmult(a, c, n)
Multiplication by constant c .
Definition math.f90:417
subroutine, public cmult2(a, b, c, n)
Multiplication by constant c .
Definition math.f90:429
subroutine, public row_zero(a, m, n, e)
Sets row e to 0 in matrix a.
Definition math.f90:233
subroutine, public invcol2(a, b, n)
Vector division .
Definition math.f90:846
real(kind=rp) function, public vlsc2(u, v, n)
Compute multiplication sum .
Definition math.f90:718
pure logical function sabscmp(x, y, tol)
Return single precision absolute comparison .
Definition math.f90:123
subroutine pwmin_sca3(a, b, c, n)
Point-wise minimum of scalar and vector .
Definition math.f90:1429
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:153
subroutine pwmax_scal2(a, b, n)
Point-wise maximum of scalar and vector .
Definition math.f90:1368
real(kind=rp) function, public glsc3(a, b, c, n)
Weighted inner product .
Definition math.f90:1026
subroutine, public ascol5(a, b, c, d, e, n)
Returns .
Definition math.f90:961
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:387
subroutine, public invers2(a, b, n)
Compute inverted vector .
Definition math.f90:645
subroutine, public cadd2(a, b, s, n)
Add a scalar to vector .
Definition math.f90:480
subroutine pwmin_sca2(a, b, n)
Point-wise minimum of scalar and vector .
Definition math.f90:1417
subroutine, public cadd(a, s, n)
Add a scalar to vector .
Definition math.f90:468
subroutine, public masked_copy(a, b, mask, n, n_mask)
Copy a masked vector .
Definition math.f90:295
subroutine reorddp(b, ind, n)
reorder double precision array - inverse of swap
Definition math.f90:1246
subroutine, public addsqr2s2(a, b, c1, n)
Returns .
Definition math.f90:832
real(kind=rp) function, public glsc4(a, b, c, d, n)
Definition math.f90:1045
subroutine, public cdiv2(a, b, c, n)
Division of constant c by elements of a .
Definition math.f90:455
subroutine pwmax_vec3(a, b, c, n)
Point-wise maximum of two vectors .
Definition math.f90:1356
subroutine swapdp(b, ind, n)
sort double precision array acording to ind vector
Definition math.f90:1206
subroutine flipvi4(b, ind, n)
Flip single integer vector b and ind.
Definition math.f90:1309
subroutine, public add2s1(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on first argument)
Definition math.f90:803
real(kind=rp) function, public glsc2(a, b, n)
Weighted inner product .
Definition math.f90:1007
subroutine, public masked_scatter_copy_0(a, b, mask, n, n_mask)
Scatter a contigous vector to masked positions in a target array .
Definition math.f90:364
subroutine, public subcol3(a, b, c, n)
Returns .
Definition math.f90:887
subroutine, public rone(a, n)
Set all elements to one.
Definition math.f90:244
subroutine flipvdp(b, ind, n)
Flip double precision vector b and ind.
Definition math.f90:1286
subroutine, public masked_gather_copy(a, b, mask, n, n_mask)
Gather a masked vector to reduced contigous vector .
Definition math.f90:341
subroutine, public x_update(a, b, c, c1, c2, n)
Returns .
Definition math.f90:992
subroutine, public add3(a, b, c, n)
Vector addition .
Definition math.f90:745
subroutine swapi4(b, ind, n)
sort single integer array acording to ind vector
Definition math.f90:1226
integer function, public glimin(a, n)
Min of an integer vector of length n.
Definition math.f90:567
real(kind=rp) function, public glsum(a, n)
Sum a vector of length n.
Definition math.f90:505
subroutine, public sub3(a, b, c, n)
Vector subtraction .
Definition math.f90:787
subroutine, public addcol4(a, b, c, d, n)
Returns .
Definition math.f90:946
subroutine, public add2(a, b, n)
Vector addition .
Definition math.f90:732
subroutine, public cfill(a, c, n)
Set all elements to a constant c .
Definition math.f90:493
subroutine, public absval(a, n)
Take the absolute value of an array.
Definition math.f90:1331
subroutine, public invcol3(a, b, c, n)
Invert a vector .
Definition math.f90:632
subroutine, public add3s2(a, b, c, c1, c2, n)
Returns .
Definition math.f90:901
subroutine, public masked_gather_copy_0(a, b, mask, n, n_mask)
Gather a masked vector to reduced contigous vector .
Definition math.f90:318
subroutine, public subcol4(a, b, c, d, n)
Returns .
Definition math.f90:917
subroutine sorti4(a, ind, n)
Heap Sort for single integer arrays.
Definition math.f90:1150
subroutine, public addcol3(a, b, c, n)
Returns .
Definition math.f90:932
subroutine, public invcol1(a, n)
Invert a vector .
Definition math.f90:620
subroutine, public cdiv(a, c, n)
Division of constant c by elements of a .
Definition math.f90:443
real(kind=rp), parameter, public neko_m_ln2
Definition math.f90:72
subroutine, public chsign(a, n)
Change sign of vector .
Definition math.f90:585
subroutine, public col2(a, b, n)
Vector multiplication .
Definition math.f90:860
subroutine, public izero(a, n)
Zero an integer vector.
Definition math.f90:222
real(kind=rp) function, public glmax(a, n)
Max of a vector of length n.
Definition math.f90:522
subroutine, public copy(a, b, n)
Copy a vector .
Definition math.f90:255
subroutine pwmax_scal3(a, b, c, n)
Point-wise maximum of scalar and vector .
Definition math.f90:1380
subroutine, public add4(a, b, c, d, n)
Vector addition .
Definition math.f90:759
pure logical function dabscmp(x, y, tol)
Return double precision absolute comparison .
Definition math.f90:138
subroutine, public col3(a, b, c, n)
Vector multiplication with 3 vectors .
Definition math.f90:873
subroutine pwmin_vec2(a, b, n)
Point-wise minimum of two vectors .
Definition math.f90:1393
subroutine pwmin_vec3(a, b, c, n)
Point-wise minimum of two vectors .
Definition math.f90:1405
subroutine pwmax_vec2(a, b, n)
Point-wise maximum of two vectors .
Definition math.f90:1344
pure logical function drelcmp(x, y, eps)
Return double precision relative comparison .
Definition math.f90:182
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:690
pure logical function qrelcmp(x, y, eps)
Return quad precision relative comparison .
Definition math.f90:197
subroutine, public rzero(a, n)
Zero a real vector.
Definition math.f90:211
real(kind=rp) function, public glsubnorm(a, b, n)
Returns the norm of the difference of two vectors .
Definition math.f90:1068
subroutine, public vdot2(dot, u1, u2, v1, v2, n)
Compute a dot product (2-d version) assuming vector components etc.
Definition math.f90:676
real(kind=rp) function, public vlmin(vec, n)
minimun value of a vector of length n
Definition math.f90:608
subroutine, public cfill_mask(a, c, n, mask, n_mask)
Fill a constant to a masked vector. .
Definition math.f90:403
real(kind=rp) function, public vlmax(vec, n)
maximum value of a vector of length n
Definition math.f90:597
integer function, public glimax(a, n)
Max of an integer vector of length n.
Definition math.f90:537
subroutine sortrp(a, ind, n)
Heap Sort for double precision arrays.
Definition math.f90:1092
subroutine, public sub2(a, b, n)
Vector substraction .
Definition math.f90:774
subroutine, public add2s2(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on second argument)
Definition math.f90:818
real(kind=rp) function, public glmin(a, n)
Min of a vector of length n.
Definition math.f90:552
subroutine, public masked_copy_0(a, b, mask, n, n_mask)
Copy a masked vector .
Definition math.f90:274
subroutine, public vcross(u1, u2, u3, v1, v2, v3, w1, w2, w3, n)
Compute a cross product assuming vector components etc.
Definition math.f90:659
pure logical function srelcmp(x, y, eps)
Return single precision relative comparison .
Definition math.f90:168
real(kind=rp) function, public vlsc3(u, v, w, n)
Compute multiplication sum .
Definition math.f90:704
subroutine reordi4(b, ind, n)
reorder single integer array - inverse of swap
Definition math.f90:1266
subroutine, public p_update(a, b, c, c1, c2, n)
Returns .
Definition math.f90:977
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