Neko  0.8.99
A portable framework for high-order spectral element flow simulations
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 !
60 module math
61  use num_types, only : rp, dp, sp, qp, i4
62  use comm
63  implicit none
64  private
65 
67  real(kind=rp), public, parameter :: neko_eps = epsilon(1.0_rp)
68 
70  real(kind=rp), public, parameter :: neko_m_ln2 = log(2.0_rp)
71 
73  real(kind=rp), public, parameter :: pi = 4._rp*atan(1._rp)
74 
75  interface abscmp
76  module procedure sabscmp, dabscmp, qabscmp
77  end interface abscmp
78 
79  interface sort
80  module procedure sortrp, sorti4
81  end interface sort
82 
83  interface swap
84  module procedure swapdp, swapi4
85  end interface swap
86 
87  interface reord
88  module procedure reorddp, reordi4
89  end interface reord
90 
91  interface flipv
92  module procedure flipvdp, flipvi4
93  end interface flipv
94 
95  interface relcmp
96  module procedure srelcmp, drelcmp, qrelcmp
97  end interface relcmp
98 
99  public :: abscmp, rzero, izero, row_zero, rone, copy, cmult, cadd, cfill, &
105  swap, reord, flipv, cadd2
106 
107 contains
108 
110  pure function sabscmp(x, y)
111  real(kind=sp), intent(in) :: x
112  real(kind=sp), intent(in) :: y
113  logical :: sabscmp
114 
115  sabscmp = abs(x - y) .lt. neko_eps
116 
117  end function sabscmp
118 
120  pure function dabscmp(x, y)
121  real(kind=dp), intent(in) :: x
122  real(kind=dp), intent(in) :: y
123  logical :: dabscmp
124 
125  dabscmp = abs(x - y) .lt. neko_eps
126 
127  end function dabscmp
128 
130  pure function qabscmp(x, y)
131  real(kind=qp), intent(in) :: x
132  real(kind=qp), intent(in) :: y
133  logical :: qabscmp
134 
135  qabscmp = abs(x - y) .lt. neko_eps
136 
137  end function qabscmp
138 
140  pure function srelcmp(x, y, eps)
141  real(kind=sp), intent(in) :: x
142  real(kind=sp), intent(in) :: y
143  real(kind=sp), intent(in), optional :: eps
144  logical :: srelcmp
145  if (present(eps)) then
146  srelcmp = abs(x - y) .le. eps*abs(y)
147  else
148  srelcmp = abs(x - y) .le. neko_eps*abs(y)
149  end if
150 
151  end function srelcmp
152 
154  pure function drelcmp(x, y, eps)
155  real(kind=dp), intent(in) :: x
156  real(kind=dp), intent(in) :: y
157  real(kind=dp), intent(in), optional :: eps
158  logical :: drelcmp
159  if (present(eps)) then
160  drelcmp = abs(x - y) .le. eps*abs(y)
161  else
162  drelcmp = abs(x - y) .le. neko_eps*abs(y)
163  end if
164 
165  end function drelcmp
166 
167 
169  pure function qrelcmp(x, y, eps)
170  real(kind=qp), intent(in) :: x
171  real(kind=qp), intent(in) :: y
172  real(kind=qp), intent(in), optional :: eps
173  logical :: qrelcmp
174  if (present(eps)) then
175  qrelcmp = abs(x - y)/abs(y) .lt. eps
176  else
177  qrelcmp = abs(x - y)/abs(y) .lt. neko_eps
178  end if
179 
180  end function qrelcmp
181 
183  subroutine rzero(a, n)
184  integer, intent(in) :: n
185  real(kind=rp), dimension(n), intent(inout) :: a
186  integer :: i
187 
188  do i = 1, n
189  a(i) = 0.0_rp
190  end do
191  end subroutine rzero
192 
194  subroutine izero(a, n)
195  integer, intent(in) :: n
196  integer, dimension(n), intent(inout) :: a
197  integer :: i
198 
199  do i = 1, n
200  a(i) = 0
201  end do
202  end subroutine izero
203 
205  subroutine row_zero(a, m, n, e)
206  integer, intent(in) :: m, n, e
207  real(kind=rp), intent(inout) :: a(m,n)
208  integer :: j
209 
210  do j = 1,n
211  a(e,j) = 0.0_rp
212  end do
213  end subroutine row_zero
214 
216  subroutine rone(a, n)
217  integer, intent(in) :: n
218  real(kind=rp), dimension(n), intent(inout) :: a
219  integer :: i
220 
221  do i = 1, n
222  a(i) = 1.0_rp
223  end do
224  end subroutine rone
225 
227  subroutine copy(a, b, n)
228  integer, intent(in) :: n
229  real(kind=rp), dimension(n), intent(in) :: b
230  real(kind=rp), dimension(n), intent(inout) :: a
231  integer :: i
232 
233  do i = 1, n
234  a(i) = b(i)
235  end do
236 
237  end subroutine copy
238 
246  subroutine masked_copy(a, b, mask, n, m)
247  integer, intent(in) :: n, m
248  real(kind=rp), dimension(n), intent(in) :: b
249  real(kind=rp), dimension(n), intent(inout) :: a
250  integer, dimension(0:m) :: mask
251  integer :: i, j
252 
253  do i = 1, m
254  j = mask(i)
255  a(j) = b(j)
256  end do
257 
258  end subroutine masked_copy
259 
262  subroutine cfill_mask(a, c, size, mask, mask_size)
263  integer, intent(in) :: size, mask_size
264  real(kind=rp), dimension(size), intent(inout) :: a
265  real(kind=rp), intent(in) :: c
266  integer, dimension(mask_size), intent(in) :: mask
267  integer :: i
268 
269  do i = 1, mask_size
270  a(mask(i)) = c
271  end do
272 
273  end subroutine cfill_mask
274 
276  subroutine cmult(a, c, n)
277  integer, intent(in) :: n
278  real(kind=rp), dimension(n), intent(inout) :: a
279  real(kind=rp), intent(in) :: c
280  integer :: i
281 
282  do i = 1, n
283  a(i) = c * a(i)
284  end do
285  end subroutine cmult
286 
288  subroutine cadd(a, s, n)
289  integer, intent(in) :: n
290  real(kind=rp), dimension(n), intent(inout) :: a
291  real(kind=rp), intent(in) :: s
292  integer :: i
293 
294  do i = 1, n
295  a(i) = a(i) + s
296  end do
297  end subroutine cadd
298 
300  subroutine cadd2(a, b, s, n)
301  integer, intent(in) :: n
302  real(kind=rp), dimension(n), intent(inout) :: a
303  real(kind=rp), dimension(n), intent(in) :: b
304  real(kind=rp), intent(in) :: s
305  integer :: i
306 
307  do i = 1, n
308  a(i) = b(i) + s
309  end do
310  end subroutine cadd2
311 
313  subroutine cfill(a, c, n)
314  integer, intent(in) :: n
315  real(kind=rp), dimension(n), intent(inout) :: a
316  real(kind=rp), intent(in) :: c
317  integer :: i
318 
319  do i = 1, n
320  a(i) = c
321  end do
322  end subroutine cfill
323 
325  function glsum(a, n)
326  integer, intent(in) :: n
327  real(kind=rp), dimension(n) :: a
328  real(kind=rp) :: tmp, glsum
329  integer :: i, ierr
330  tmp = 0.0_rp
331  do i = 1, n
332  tmp = tmp + a(i)
333  end do
334  call mpi_allreduce(tmp, glsum, 1, &
335  mpi_real_precision, mpi_sum, neko_comm, ierr)
336 
337  end function glsum
338 
340  function glmax(a, n)
341  integer, intent(in) :: n
342  real(kind=rp), dimension(n) :: a
343  real(kind=rp) :: tmp, glmax
344  integer :: i, ierr
345  tmp = a(1)
346  do i = 2, n
347  tmp = max(tmp,a(i))
348  end do
349  call mpi_allreduce(tmp, glmax, 1, &
350  mpi_real_precision, mpi_max, neko_comm, ierr)
351  end function glmax
352 
354  function glimax(a, n)
355  integer, intent(in) :: n
356  integer, dimension(n) :: a
357  integer :: tmp, glimax
358  integer :: i, ierr
359  tmp = a(1)
360  do i = 2, n
361  tmp = max(tmp,a(i))
362  end do
363  call mpi_allreduce(tmp, glimax, 1, &
364  mpi_integer, mpi_max, neko_comm, ierr)
365  end function glimax
366 
368  function glmin(a, n)
369  integer, intent(in) :: n
370  real(kind=rp), dimension(n) :: a
371  real(kind=rp) :: tmp, glmin
372  integer :: i, ierr
373  tmp = a(1)
374  do i = 2, n
375  tmp = min(tmp,a(i))
376  end do
377  call mpi_allreduce(tmp, glmin, 1, &
378  mpi_real_precision, mpi_min, neko_comm, ierr)
379  end function glmin
380 
382  function glimin(a, n)
383  integer, intent(in) :: n
384  integer, dimension(n) :: a
385  integer :: tmp, glimin
386  integer :: i, ierr
387  tmp = a(1)
388  do i = 2, n
389  tmp = min(tmp,a(i))
390  end do
391  call mpi_allreduce(tmp, glimin, 1, &
392  mpi_integer, mpi_min, neko_comm, ierr)
393  end function glimin
394 
395 
396 
397 
399  subroutine chsign(a, n)
400  integer, intent(in) :: n
401  real(kind=rp), dimension(n), intent(inout) :: a
402  integer :: i
403 
404  do i = 1, n
405  a(i) = -a(i)
406  end do
407 
408  end subroutine chsign
409 
411  function vlmax(vec,n) result(tmax)
412  integer :: n, i
413  real(kind=rp), intent(in) :: vec(n)
414  real(kind=rp) :: tmax
415  tmax = real(-99d20, rp)
416  do i = 1, n
417  tmax = max(tmax, vec(i))
418  end do
419  end function vlmax
420 
422  function vlmin(vec,n) result(tmin)
423  integer, intent(in) :: n
424  real(kind=rp), intent(in) :: vec(n)
425  real(kind=rp) :: tmin
426  integer :: i
427  tmin = real(99.0e20, rp)
428  do i = 1, n
429  tmin = min(tmin, vec(i))
430  end do
431  end function vlmin
432 
434  subroutine invcol1(a, n)
435  integer, intent(in) :: n
436  real(kind=rp), dimension(n), intent(inout) :: a
437  integer :: i
438 
439  do i = 1, n
440  a(i) = 1.0_rp / a(i)
441  end do
442 
443  end subroutine invcol1
444 
446  subroutine invcol3(a, b, c, n)
447  integer, intent(in) :: n
448  real(kind=rp), dimension(n), intent(inout) :: a
449  real(kind=rp), dimension(n), intent(in) :: b,c
450  integer :: i
451 
452  do i = 1, n
453  a(i) = b(i) / c(i)
454  end do
455 
456  end subroutine invcol3
457 
459  subroutine invers2(a, b, n)
460  integer, intent(in) :: n
461  real(kind=rp), dimension(n), intent(inout) :: a
462  real(kind=rp), dimension(n), intent(in) :: b
463  integer :: i
464 
465  do i = 1, n
466  a(i) = 1.0_rp / b(i)
467  end do
468 
469  end subroutine invers2
470 
473  subroutine vcross(u1, u2, u3, v1, v2, v3, w1, w2, w3, n)
474  integer, intent(in) :: n
475  real(kind=rp), dimension(n), intent(in) :: v1, v2, v3
476  real(kind=rp), dimension(n), intent(in) :: w1, w2, w3
477  real(kind=rp), dimension(n), intent(out) :: u1, u2, u3
478  integer :: i
479 
480  do i = 1, n
481  u1(i) = v2(i)*w3(i) - v3(i)*w2(i)
482  u2(i) = v3(i)*w1(i) - v1(i)*w3(i)
483  u3(i) = v1(i)*w2(i) - v2(i)*w1(i)
484  end do
485 
486  end subroutine vcross
487 
490  subroutine vdot2(dot, u1, u2, v1, v2, n)
491  integer, intent(in) :: n
492  real(kind=rp), dimension(n), intent(in) :: u1, u2
493  real(kind=rp), dimension(n), intent(in) :: v1, v2
494  real(kind=rp), dimension(n), intent(out) :: dot
495  integer :: i
496  do i = 1, n
497  dot(i) = u1(i)*v1(i) + u2(i)*v2(i)
498  end do
499 
500  end subroutine vdot2
501 
504  subroutine vdot3(dot, u1, u2, u3, v1, v2, v3, n)
505  integer, intent(in) :: n
506  real(kind=rp), dimension(n), intent(in) :: u1, u2, u3
507  real(kind=rp), dimension(n), intent(in) :: v1, v2, v3
508  real(kind=rp), dimension(n), intent(out) :: dot
509  integer :: i
510 
511  do i = 1, n
512  dot(i) = u1(i)*v1(i) + u2(i)*v2(i) + u3(i)*v3(i)
513  end do
514 
515  end subroutine vdot3
516 
518  function vlsc3(u, v, w, n) result(s)
519  integer, intent(in) :: n
520  real(kind=rp), dimension(n), intent(in) :: u, v, w
521  real(kind=rp) :: s
522  integer :: i
523 
524  s = 0.0_rp
525  do i = 1, n
526  s = s + u(i)*v(i)*w(i)
527  end do
528 
529  end function vlsc3
530 
532  function vlsc2(u, v, n) result(s)
533  integer, intent(in) :: n
534  real(kind=rp), dimension(n), intent(in) :: u, v
535  real(kind=rp) :: s
536  integer :: i
537 
538  s = 0.0_rp
539  do i = 1, n
540  s = s + u(i)*v(i)
541  end do
542 
543  end function vlsc2
544 
546  subroutine add2(a, b, n)
547  integer, intent(in) :: n
548  real(kind=rp), dimension(n), intent(inout) :: a
549  real(kind=rp), dimension(n), intent(in) :: b
550  integer :: i
551 
552  do i = 1, n
553  a(i) = a(i) + b(i)
554  end do
555 
556  end subroutine add2
557 
559  subroutine add3(a, b, c, n)
560  integer, intent(in) :: n
561  real(kind=rp), dimension(n), intent(inout) :: a
562  real(kind=rp), dimension(n), intent(in) :: b
563  real(kind=rp), dimension(n), intent(in) :: c
564  integer :: i
565 
566  do i = 1, n
567  a(i) = b(i) + c(i)
568  end do
569 
570  end subroutine add3
571 
573  subroutine add4(a, b, c, d, n)
574  integer, intent(in) :: n
575  real(kind=rp), dimension(n), intent(inout) :: d
576  real(kind=rp), dimension(n), intent(inout) :: c
577  real(kind=rp), dimension(n), intent(inout) :: b
578  real(kind=rp), dimension(n), intent(out) :: a
579  integer :: i
580 
581  do i = 1, n
582  a(i) = b(i) + c(i) + d(i)
583  end do
584 
585  end subroutine add4
586 
588  subroutine sub2(a, b, n)
589  integer, intent(in) :: n
590  real(kind=rp), dimension(n), intent(inout) :: a
591  real(kind=rp), dimension(n), intent(inout) :: b
592  integer :: i
593 
594  do i = 1, n
595  a(i) = a(i) - b(i)
596  end do
597 
598  end subroutine sub2
599 
601  subroutine sub3(a, b, c, n)
602  integer, intent(in) :: n
603  real(kind=rp), dimension(n), intent(inout) :: a
604  real(kind=rp), dimension(n), intent(in) :: b
605  real(kind=rp), dimension(n), intent(in) :: c
606  integer :: i
607 
608  do i = 1, n
609  a(i) = b(i) - c(i)
610  end do
611 
612  end subroutine sub3
613 
614 
617  subroutine add2s1(a, b, c1, n)
618  integer, intent(in) :: n
619  real(kind=rp), dimension(n), intent(inout) :: a
620  real(kind=rp), dimension(n), intent(inout) :: b
621  real(kind=rp), intent(in) :: c1
622  integer :: i
623 
624  do i = 1, n
625  a(i) = c1 * a(i) + b(i)
626  end do
627 
628  end subroutine add2s1
629 
632  subroutine add2s2(a, b, c1, n)
633  integer, intent(in) :: n
634  real(kind=rp), dimension(n), intent(inout) :: a
635  real(kind=rp), dimension(n), intent(inout) :: b
636  real(kind=rp), intent(in) :: c1
637  integer :: i
638 
639  do i = 1, n
640  a(i) = a(i) + c1 * b(i)
641  end do
642 
643  end subroutine add2s2
644 
646  subroutine addsqr2s2(a, b, c1, n)
647  integer, intent(in) :: n
648  real(kind=rp), dimension(n), intent(inout) :: a
649  real(kind=rp), dimension(n), intent(in) :: b
650  real(kind=rp), intent(in) :: c1
651  integer :: i
652 
653  do i = 1,n
654  a(i) = a(i) + c1 * ( b(i) * b(i) )
655  end do
656 
657  end subroutine addsqr2s2
658 
660  subroutine cmult2(a, b, c, n)
661  integer, intent(in) :: n
662  real(kind=rp), dimension(n), intent(inout) :: a
663  real(kind=rp), dimension(n), intent(in) :: b
664  real(kind=rp), intent(in) :: c
665  integer :: i
666 
667  do i = 1, n
668  a(i) = c * b(i)
669  end do
670 
671  end subroutine cmult2
672 
674  subroutine invcol2(a, b, n)
675  integer, intent(in) :: n
676  real(kind=rp), dimension(n), intent(inout) :: a
677  real(kind=rp), dimension(n), intent(in) :: b
678  integer :: i
679 
680  do i = 1, n
681  a(i) = a(i) /b(i)
682  end do
683 
684  end subroutine invcol2
685 
686 
688  subroutine col2(a, b, n)
689  integer, intent(in) :: n
690  real(kind=rp), dimension(n), intent(inout) :: a
691  real(kind=rp), dimension(n), intent(in) :: b
692  integer :: i
693 
694  do i = 1, n
695  a(i) = a(i) * b(i)
696  end do
697 
698  end subroutine col2
699 
701  subroutine col3(a, b, c, n)
702  integer, intent(in) :: n
703  real(kind=rp), dimension(n), intent(inout) :: a
704  real(kind=rp), dimension(n), intent(in) :: b
705  real(kind=rp), dimension(n), intent(in) :: c
706  integer :: i
707 
708  do i = 1, n
709  a(i) = b(i) * c(i)
710  end do
711 
712  end subroutine col3
713 
715  subroutine subcol3(a, b, c, n)
716  integer, intent(in) :: n
717  real(kind=rp), dimension(n), intent(inout) :: a
718  real(kind=rp), dimension(n), intent(in) :: b
719  real(kind=rp), dimension(n), intent(in) :: c
720  integer :: i
721 
722  do i = 1,n
723  a(i) = a(i) - b(i) * c(i)
724  end do
725 
726  end subroutine subcol3
727 
729  subroutine add3s2(a, b, c, c1, c2 ,n)
730  integer, intent(in) :: n
731  real(kind=rp), dimension(n), intent(inout) :: a
732  real(kind=rp), dimension(n), intent(in) :: b
733  real(kind=rp), dimension(n), intent(in) :: c
734  real(kind=rp), intent(in) :: c1, c2
735  integer :: i
736 
737  do i = 1,n
738  a(i) = c1 * b(i) + c2 * c(i)
739  end do
740 
741  end subroutine add3s2
742 
743 
745  subroutine subcol4(a, b, c, d, n)
746  integer, intent(in) :: n
747  real(kind=rp), dimension(n), intent(inout) :: a
748  real(kind=rp), dimension(n), intent(in) :: b
749  real(kind=rp), dimension(n), intent(in) :: c
750  real(kind=rp), dimension(n), intent(in) :: d
751  integer :: i
752 
753  do i = 1,n
754  a(i) = a(i) - b(i) * c(i) * d(i)
755  end do
756 
757  end subroutine subcol4
758 
760  subroutine addcol3(a, b, c, n)
761  integer, intent(in) :: n
762  real(kind=rp), dimension(n), intent(inout) :: a
763  real(kind=rp), dimension(n), intent(in) :: b
764  real(kind=rp), dimension(n), intent(in) :: c
765  integer :: i
766 
767  do i = 1,n
768  a(i) = a(i) + b(i) * c(i)
769  end do
770 
771  end subroutine addcol3
772 
774  subroutine addcol4(a, b, c, d, n)
775  integer, intent(in) :: n
776  real(kind=rp), dimension(n), intent(inout) :: a
777  real(kind=rp), dimension(n), intent(in) :: b
778  real(kind=rp), dimension(n), intent(in) :: c
779  real(kind=rp), dimension(n), intent(in) :: d
780  integer :: i
781 
782  do i = 1,n
783  a(i) = a(i) + b(i) * c(i) * d(i)
784  end do
785 
786  end subroutine addcol4
787 
789  subroutine ascol5(a, b, c, d, e, n)
790  integer, intent(in) :: n
791  real(kind=rp), dimension(n), intent(inout) :: a
792  real(kind=rp), dimension(n), intent(in) :: b
793  real(kind=rp), dimension(n), intent(in) :: c
794  real(kind=rp), dimension(n), intent(in) :: d
795  real(kind=rp), dimension(n), intent(in) :: e
796  integer :: i
797 
798  do i = 1,n
799  a(i) = b(i)*c(i)-d(i)*e(i)
800  end do
801 
802  end subroutine ascol5
803 
805  subroutine p_update(a, b, c, c1, c2, n)
806  integer, intent(in) :: n
807  real(kind=rp), dimension(n), intent(inout) :: a
808  real(kind=rp), dimension(n), intent(in) :: b
809  real(kind=rp), dimension(n), intent(in) :: c
810  real(kind=rp), intent(in) :: c1, c2
811  integer :: i
812 
813  do i = 1,n
814  a(i) = b(i) + c1*(a(i)-c2*c(i))
815  end do
816 
817  end subroutine p_update
818 
820  subroutine x_update(a, b, c, c1, c2, n)
821  integer, intent(in) :: n
822  real(kind=rp), dimension(n), intent(inout) :: a
823  real(kind=rp), dimension(n), intent(in) :: b
824  real(kind=rp), dimension(n), intent(in) :: c
825  real(kind=rp), intent(in) :: c1, c2
826  integer :: i
827 
828  do i = 1,n
829  a(i) = a(i) + c1*b(i)+c2*c(i)
830  end do
831 
832  end subroutine x_update
833 
835  function glsc2(a, b, n)
836  integer, intent(in) :: n
837  real(kind=rp), dimension(n), intent(in) :: a
838  real(kind=rp), dimension(n), intent(in) :: b
839  real(kind=rp) :: glsc2, tmp
840  integer :: i, ierr
841 
842  tmp = 0.0_rp
843  do i = 1, n
844  tmp = tmp + a(i) * b(i)
845  end do
846 
847  call mpi_allreduce(tmp, glsc2, 1, &
848  mpi_real_precision, mpi_sum, neko_comm, ierr)
849 
850  end function glsc2
851 
853  function glsc3(a, b, c, n)
854  integer, intent(in) :: n
855  real(kind=rp), dimension(n), intent(in) :: a
856  real(kind=rp), dimension(n), intent(in) :: b
857  real(kind=rp), dimension(n), intent(in) :: c
858  real(kind=rp) :: glsc3, tmp
859  integer :: i, ierr
860 
861  tmp = 0.0_rp
862  do i = 1, n
863  tmp = tmp + a(i) * b(i) * c(i)
864  end do
865 
866  call mpi_allreduce(tmp, glsc3, 1, &
867  mpi_real_precision, mpi_sum, neko_comm, ierr)
868 
869  end function glsc3
870  function glsc4(a, b, c, d, n)
871  integer, intent(in) :: n
872  real(kind=rp), dimension(n), intent(in) :: a
873  real(kind=rp), dimension(n), intent(in) :: b
874  real(kind=rp), dimension(n), intent(in) :: c
875  real(kind=rp), dimension(n), intent(in) :: d
876  real(kind=rp) :: glsc4, tmp
877  integer :: i, ierr
878 
879  tmp = 0.0_rp
880  do i = 1, n
881  tmp = tmp + a(i) * b(i) * c(i) * d(i)
882  end do
883 
884  call mpi_allreduce(tmp, glsc4, 1, &
885  mpi_real_precision, mpi_sum, neko_comm, ierr)
886 
887  end function glsc4
888 
894  subroutine sortrp(a, ind, n)
895  integer, intent(in) :: n
896  real(kind=rp), intent(inout) :: a(n)
897  integer, intent(out) :: ind(n)
898  real(kind=rp) :: aa
899  integer :: j, ir, i, ii, l
900 
901  do j = 1, n
902  ind(j) = j
903  end do
904 
905  if (n .le. 1) return
906 
907 
908  l = n/2+1
909  ir = n
910  do while (.true.)
911  if (l .gt. 1) then
912  l = l-1
913  aa = a(l)
914  ii = ind(l)
915  else
916  aa = a(ir)
917  ii = ind(ir)
918  a(ir) = a(1)
919  ind(ir) = ind(1)
920  ir = ir - 1
921  if (ir .eq. 1) then
922  a(1) = aa
923  ind(1) = ii
924  return
925  end if
926  end if
927  i = l
928  j = l+l
929  do while (j .le. ir)
930  if (j .lt. ir) then
931  if ( a(j) .lt. a(j+1) ) j = j + 1
932  end if
933  if (aa .lt. a(j)) then
934  a(i) = a(j)
935  ind(i) = ind(j)
936  i = j
937  j = j+j
938  else
939  j = ir+1
940  end if
941  end do
942  a(i) = aa
943  ind(i) = ii
944  end do
945  end subroutine sortrp
946 
952  subroutine sorti4(a, ind, n)
953  integer, intent(in) :: n
954  integer(i4), intent(inout) :: a(n)
955  integer, intent(out) :: ind(n)
956  integer(i4) :: aa
957  integer :: j, ir, i, ii, l
958 
959  do j = 1, n
960  ind(j) = j
961  end do
962 
963  if (n .le. 1) return
964 
965  l = n/2+1
966  ir = n
967  do while (.true.)
968  if (l .gt. 1) then
969  l = l - 1
970  aa = a(l)
971  ii = ind(l)
972  else
973  aa = a(ir)
974  ii = ind(ir)
975  a(ir) = a( 1)
976  ind(ir) = ind( 1)
977  ir = ir - 1
978  if (ir .eq. 1) then
979  a(1) = aa
980  ind(1) = ii
981  return
982  end if
983  end if
984  i = l
985  j = l + l
986  do while (j .le. ir)
987  if (j .lt. ir) then
988  if ( a(j) .lt. a(j + 1) ) j = j + 1
989  end if
990  if (aa .lt. a(j)) then
991  a(i) = a(j)
992  ind(i) = ind(j)
993  i = j
994  j = j + j
995  else
996  j = ir + 1
997  end if
998  end do
999  a(i) = aa
1000  ind(i) = ii
1001  end do
1002  end subroutine sorti4
1003 
1008  subroutine swapdp(b, ind, n)
1009  integer, intent(in) :: n
1010  real(kind=rp), intent(inout) :: b(n)
1011  integer, intent(in) :: ind(n)
1012  real(kind=rp) :: temp(n)
1013  integer :: i, jj
1014 
1015  do i = 1, n
1016  temp(i) = b(i)
1017  end do
1018  do i = 1, n
1019  jj = ind(i)
1020  b(i) = temp(jj)
1021  end do
1022  end subroutine swapdp
1023 
1028  subroutine swapi4(b, ind, n)
1029  integer, intent(in) :: n
1030  integer(i4), intent(inout) :: b(n)
1031  integer, intent(in) :: ind(n)
1032  integer(i4) :: temp(n)
1033  integer :: i, jj
1034 
1035  do i = 1, n
1036  temp(i) = b(i)
1037  end do
1038  do i = 1, n
1039  jj = ind(i)
1040  b(i) = temp(jj)
1041  end do
1042  end subroutine swapi4
1043 
1048  subroutine reorddp(b, ind, n)
1049  integer, intent(in) :: n
1050  real(kind=rp), intent(inout) :: b(n)
1051  integer, intent(in) :: ind(n)
1052  real(kind=rp) :: temp(n)
1053  integer :: i, jj
1054 
1055  do i = 1, n
1056  temp(i) = b(i)
1057  end do
1058  do i = 1, n
1059  jj = ind(i)
1060  b(jj) = temp(i)
1061  end do
1062  end subroutine reorddp
1063 
1068  subroutine reordi4(b, ind, n)
1069  integer, intent(in) :: n
1070  integer(i4), intent(inout) :: b(n)
1071  integer, intent(in) :: ind(n)
1072  integer(i4) :: temp(n)
1073  integer :: i, jj
1074 
1075  do i = 1, n
1076  temp(i) = b(i)
1077  end do
1078  do i = 1, n
1079  jj = ind(i)
1080  b(jj) = temp(i)
1081  end do
1082  end subroutine reordi4
1083 
1088  subroutine flipvdp(b, ind, n)
1089  integer, intent(in) :: n
1090  real(kind=rp), intent(inout) :: b(n)
1091  integer, intent(inout) :: ind(n)
1092  real(kind=rp) :: temp(n)
1093  integer :: tempind(n)
1094  integer :: i, jj
1095 
1096  do i = 1, n
1097  jj = n+1-i
1098  temp(jj) = b(i)
1099  tempind(jj) = ind(i)
1100  end do
1101  do i = 1,n
1102  b(i) = temp(i)
1103  ind(i) = tempind(i)
1104  end do
1105  end subroutine flipvdp
1106 
1111  subroutine flipvi4(b, ind, n)
1112  integer, intent(in) :: n
1113  integer(i4), intent(inout) :: b(n)
1114  integer, intent(inout) :: ind(n)
1115  integer(i4) :: temp(n)
1116  integer :: tempind(n)
1117  integer :: i, jj
1118 
1119  do i = 1, n
1120  jj = n+1-i
1121  temp(jj) = b(i)
1122  tempind(jj) = ind(i)
1123  end do
1124  do i = 1,n
1125  b(i) = temp(i)
1126  ind(i) = tempind(i)
1127  end do
1128  end subroutine flipvi4
1129 
1130 end module math
double real
Definition: device_config.h:12
Definition: comm.F90:1
type(mpi_comm) neko_comm
MPI communicator.
Definition: comm.F90:16
type(mpi_datatype) mpi_real_precision
MPI type for working precision of REAL types.
Definition: comm.F90:22
Definition: math.f90:60
subroutine, public cmult(a, c, n)
Multiplication by constant c .
Definition: math.f90:277
subroutine, public cmult2(a, b, c, n)
Multiplication by constant c .
Definition: math.f90:661
subroutine, public row_zero(a, m, n, e)
Sets row e to 0 in matrix a.
Definition: math.f90:206
subroutine, public invcol2(a, b, n)
Vector division .
Definition: math.f90:675
real(kind=rp) function, public vlsc2(u, v, n)
Compute multiplication sum .
Definition: math.f90:533
real(kind=rp), parameter, public pi
Definition: math.f90:73
real(kind=rp) function, public glsc3(a, b, c, n)
Weighted inner product .
Definition: math.f90:854
pure logical function qabscmp(x, y)
Return double precision absolute comparison .
Definition: math.f90:131
subroutine, public ascol5(a, b, c, d, e, n)
Returns .
Definition: math.f90:790
subroutine, public invers2(a, b, n)
Compute inverted vector .
Definition: math.f90:460
subroutine, public cadd2(a, b, s, n)
Add a scalar to vector .
Definition: math.f90:301
subroutine, public cadd(a, s, n)
Add a scalar to vector .
Definition: math.f90:289
subroutine reorddp(b, ind, n)
reorder double precision array - inverse of swap
Definition: math.f90:1049
subroutine, public addsqr2s2(a, b, c1, n)
Returns .
Definition: math.f90:647
real(kind=rp) function, public glsc4(a, b, c, d, n)
Definition: math.f90:871
subroutine swapdp(b, ind, n)
sort double precision array acording to ind vector
Definition: math.f90:1009
subroutine flipvi4(b, ind, n)
Flip single integer vector b and ind.
Definition: math.f90:1112
subroutine, public add2s1(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on first argument)
Definition: math.f90:618
real(kind=rp) function, public glsc2(a, b, n)
Weighted inner product .
Definition: math.f90:836
subroutine, public subcol3(a, b, c, n)
Returns .
Definition: math.f90:716
subroutine, public rone(a, n)
Set all elements to one.
Definition: math.f90:217
subroutine flipvdp(b, ind, n)
Flip double precision vector b and ind.
Definition: math.f90:1089
subroutine, public x_update(a, b, c, c1, c2, n)
Returns .
Definition: math.f90:821
subroutine, public add3(a, b, c, n)
Vector addition .
Definition: math.f90:560
subroutine swapi4(b, ind, n)
sort single integer array acording to ind vector
Definition: math.f90:1029
integer function, public glimin(a, n)
Min of an integer vector of length n.
Definition: math.f90:383
real(kind=rp) function, public glsum(a, n)
Sum a vector of length n.
Definition: math.f90:326
subroutine, public sub3(a, b, c, n)
Vector subtraction .
Definition: math.f90:602
subroutine, public addcol4(a, b, c, d, n)
Returns .
Definition: math.f90:775
subroutine, public add2(a, b, n)
Vector addition .
Definition: math.f90:547
subroutine, public cfill(a, c, n)
Set all elements to a constant c .
Definition: math.f90:314
subroutine, public invcol3(a, b, c, n)
Invert a vector .
Definition: math.f90:447
pure logical function sabscmp(x, y)
Return single precision absolute comparison .
Definition: math.f90:111
subroutine, public add3s2(a, b, c, c1, c2, n)
Returns .
Definition: math.f90:730
subroutine, public subcol4(a, b, c, d, n)
Returns .
Definition: math.f90:746
subroutine sorti4(a, ind, n)
Heap Sort for single integer arrays.
Definition: math.f90:953
subroutine, public addcol3(a, b, c, n)
Returns .
Definition: math.f90:761
pure logical function dabscmp(x, y)
Return double precision absolute comparison .
Definition: math.f90:121
subroutine, public invcol1(a, n)
Invert a vector .
Definition: math.f90:435
subroutine, public masked_copy(a, b, mask, n, m)
Copy a masked vector .
Definition: math.f90:247
real(kind=rp), parameter, public neko_m_ln2
Definition: math.f90:70
subroutine, public chsign(a, n)
Change sign of vector .
Definition: math.f90:400
subroutine, public col2(a, b, n)
Vector multiplication .
Definition: math.f90:689
subroutine, public izero(a, n)
Zero an integer vector.
Definition: math.f90:195
real(kind=rp) function, public glmax(a, n)
Max of a vector of length n.
Definition: math.f90:341
subroutine, public copy(a, b, n)
Copy a vector .
Definition: math.f90:228
subroutine, public add4(a, b, c, d, n)
Vector addition .
Definition: math.f90:574
subroutine, public col3(a, b, c, n)
Vector multiplication with 3 vectors .
Definition: math.f90:702
pure logical function drelcmp(x, y, eps)
Return double precision relative comparison .
Definition: math.f90:155
real(kind=rp), parameter, public neko_eps
Machine epsilon .
Definition: math.f90:67
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:505
pure logical function qrelcmp(x, y, eps)
Return quad precision relative comparison .
Definition: math.f90:170
subroutine, public rzero(a, n)
Zero a real vector.
Definition: math.f90:184
subroutine, public vdot2(dot, u1, u2, v1, v2, n)
Compute a dot product (2-d version) assuming vector components etc.
Definition: math.f90:491
real(kind=rp) function, public vlmin(vec, n)
minimun value of a vector of length n
Definition: math.f90:423
real(kind=rp) function, public vlmax(vec, n)
maximum value of a vector of length n
Definition: math.f90:412
integer function, public glimax(a, n)
Max of an integer vector of length n.
Definition: math.f90:355
subroutine sortrp(a, ind, n)
Heap Sort for double precision arrays.
Definition: math.f90:895
subroutine, public sub2(a, b, n)
Vector substraction .
Definition: math.f90:589
subroutine, public cfill_mask(a, c, size, mask, mask_size)
Fill a constant to a masked vector. .
Definition: math.f90:263
subroutine, public add2s2(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on second argument)
Definition: math.f90:633
real(kind=rp) function, public glmin(a, n)
Min of a vector of length n.
Definition: math.f90:369
subroutine, public vcross(u1, u2, u3, v1, v2, v3, w1, w2, w3, n)
Compute a cross product assuming vector components etc.
Definition: math.f90:474
pure logical function srelcmp(x, y, eps)
Return single precision relative comparison .
Definition: math.f90:141
real(kind=rp) function, public vlsc3(u, v, w, n)
Compute multiplication sum .
Definition: math.f90:519
subroutine reordi4(b, ind, n)
reorder single integer array - inverse of swap
Definition: math.f90:1069
subroutine, public p_update(a, b, c, c1, c2, n)
Returns .
Definition: math.f90:806
integer, parameter, public qp
Definition: num_types.f90:10
integer, parameter, public i4
Definition: num_types.f90:6
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