Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
field_math.f90
Go to the documentation of this file.
1! Copyright (c) 2008-2020, UCHICAGO ARGONNE, LLC.
2!
3! The UChicago Argonne, LLC as Operator of Argonne National
4! Laboratory holds copyright in the Software. The copyright holder
5! reserves all rights except those expressly granted to licensees,
6! and U.S. Government license rights.
7!
8! Redistribution and use in source and binary forms, with or without
9! modification, are permitted provided that the following conditions
10! are met:
11!
12! 1. Redistributions of source code must retain the above copyright
13! notice, this list of conditions and the disclaimer below.
14!
15! 2. Redistributions in binary form must reproduce the above copyright
16! notice, this list of conditions and the disclaimer (as noted below)
17! in the documentation and/or other materials provided with the
18! distribution.
19!
20! 3. Neither the name of ANL nor the names of its contributors
21! may be used to endorse or promote products derived from this software
22! without specific prior written permission.
23!
24! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
25! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
26! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
27! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL
28! UCHICAGO ARGONNE, LLC, THE U.S. DEPARTMENT OF
29! ENERGY OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
30! SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED
31! TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
32! DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
33! THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
34! (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
35! OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36!
37! Additional BSD Notice
38! ---------------------
39! 1. This notice is required to be provided under our contract with
40! the U.S. Department of Energy (DOE). This work was produced at
41! Argonne National Laboratory under Contract
42! No. DE-AC02-06CH11357 with the DOE.
43!
44! 2. Neither the United States Government nor UCHICAGO ARGONNE,
45! LLC nor any of their employees, makes any warranty,
46! express or implied, or assumes any liability or responsibility for the
47! accuracy, completeness, or usefulness of any information, apparatus,
48! product, or process disclosed, or represents that its use would not
49! infringe privately-owned rights.
50!
51! 3. Also, reference herein to any specific commercial products, process,
52! or services by trade name, trademark, manufacturer or otherwise does
53! not necessarily constitute or imply its endorsement, recommendation,
54! or favoring by the United States Government or UCHICAGO ARGONNE LLC.
55! The views and opinions of authors expressed
56! herein do not necessarily state or reflect those of the United States
57! Government or UCHICAGO ARGONNE, LLC, and shall
58! not be used for advertising or product endorsement purposes.
59!
62 use num_types, only : rp
63 use field, only : field_t
64 use device, only : device_get_ptr
65 use math, only : rzero, rone, copy, cmult, cadd, cfill, invcol1, vdot3, add2, &
81 use, intrinsic :: iso_c_binding, only : c_ptr
82 implicit none
83 private
84
96
97contains
98
100 subroutine field_rzero(a, n)
101 integer, intent(in), optional :: n
102 type(field_t), intent(inout) :: a
103 integer :: size
104
105 if (present(n)) then
106 size = n
107 else
108 size = a%size()
109 end if
110
111 if (neko_bcknd_device .eq. 1) then
112 call device_rzero(a%x_d, size)
113 else
114 call rzero(a%x, size)
115 end if
116 end subroutine field_rzero
117
119 subroutine field_rone(a, n)
120 integer, intent(in), optional :: n
121 type(field_t), intent(inout) :: a
122 integer :: size
123
124 if (present(n)) then
125 size = n
126 else
127 size = a%size()
128 end if
129
130 if (neko_bcknd_device .eq. 1) then
131 call device_rone(a%x_d, size)
132 else
133 call rone(a%x, size)
134 end if
135 end subroutine field_rone
136
138 subroutine field_copy(a, b, n)
139 integer, intent(in), optional :: n
140 type(field_t), intent(in) :: b
141 type(field_t), intent(inout) :: a
142 integer :: size
143
144 if (present(n)) then
145 size = n
146 else
147 size = a%size()
148 end if
149
150 if (neko_bcknd_device .eq. 1) then
151 call device_copy(a%x_d, b%x_d, size)
152 else
153 call copy(a%x, b%x, size)
154 end if
155 end subroutine field_copy
156
158 subroutine field_cmult(a, c, n)
159 integer, intent(in), optional :: n
160 type(field_t), intent(inout) :: a
161 real(kind=rp), intent(in) :: c
162 integer :: size
163
164 if (present(n)) then
165 size = n
166 else
167 size = a%size()
168 end if
169
170 if (neko_bcknd_device .eq. 1) then
171 call device_cmult(a%x_d, c, size)
172 else
173 call cmult(a%x, c, size)
174 end if
175 end subroutine field_cmult
176
178 subroutine field_cadd(a, s, n)
179 integer, intent(in), optional :: n
180 type(field_t), intent(inout) :: a
181 real(kind=rp), intent(in) :: s
182 integer :: size
183
184 if (present(n)) then
185 size = n
186 else
187 size = a%size()
188 end if
189
190 if (neko_bcknd_device .eq. 1) then
191 call device_cadd(a%x_d, s, size)
192 else
193 call cadd(a%x, s, size)
194 end if
195 end subroutine field_cadd
196
198 subroutine field_cfill(a, c, n)
199 integer, intent(in), optional :: n
200 type(field_t), intent(inout) :: a
201 real(kind=rp), intent(in) :: c
202 integer :: size
203
204 if (present(n)) then
205 size = n
206 else
207 size = a%size()
208 end if
209
210 if (neko_bcknd_device .eq. 1) then
211 call device_cfill(a%x_d, c, size)
212 else
213 call cfill(a%x, c, size)
214 end if
215 end subroutine field_cfill
216
218 subroutine field_invcol1(a, n)
219 integer, intent(in), optional :: n
220 type(field_t), intent(inout) :: a
221 integer :: size
222
223 if (present(n)) then
224 size = n
225 else
226 size = a%size()
227 end if
228
229 if (neko_bcknd_device .eq. 1) then
230 call device_invcol1(a%x_d, size)
231 else
232 call invcol1(a%x, size)
233 end if
234
235 end subroutine field_invcol1
236
238 subroutine field_invcol3(a, b, c, n)
239 integer, intent(in), optional :: n
240 type(field_t), intent(inout) :: a
241 type(field_t), intent(in) :: b
242 type(field_t), intent(in) :: c
243 integer :: size
244
245 if (present(n)) then
246 size = n
247 else
248 size = a%size()
249 end if
250
251 if (neko_bcknd_device .eq. 1) then
252 call device_invcol3(a%x_d, b%x_d, c%x_d, size)
253 else
254 call invcol3(a%x, b%x, c%x, size)
255 end if
256
257 end subroutine field_invcol3
258
261 subroutine field_vdot3(dot, u1, u2, u3, v1, v2, v3, n)
262 integer, intent(in), optional :: n
263 type(field_t), intent(in) :: u1, u2, u3
264 type(field_t), intent(in) :: v1, v2, v3
265 type(field_t), intent(out) :: dot
266 integer :: size
267
268 if (present(n)) then
269 size = n
270 else
271 size = dot%size()
272 end if
273
274 if (neko_bcknd_device .eq. 1) then
275 call device_vdot3(dot%x_d, &
276 u1%x_d, u2%x_d, u3%x_d, &
277 v1%x_d, v2%x_d, v3%x_d, &
278 size)
279 else
280 call vdot3(dot%x, &
281 u1%x, u2%x, u3%x, &
282 v1%x, v2%x, v3%x, &
283 size)
284 end if
285
286 end subroutine field_vdot3
287
289 subroutine field_add2(a, b, n)
290 integer, intent(in), optional :: n
291 type(field_t), intent(inout) :: a
292 type(field_t), intent(in) :: b
293 integer :: size
294
295 if (present(n)) then
296 size = n
297 else
298 size = a%size()
299 end if
300
301 if (neko_bcknd_device .eq. 1) then
302 call device_add2(a%x_d, b%x_d, size)
303 else
304 call add2(a%x, b%x, size)
305 end if
306
307 end subroutine field_add2
308
310 subroutine field_add3(a, b, c, n)
311 integer, intent(in), optional :: n
312 type(field_t), intent(inout) :: a
313 type(field_t), intent(in) :: b, c
314 integer :: size
315
316 if (present(n)) then
317 size = n
318 else
319 size = a%size()
320 end if
321
322 if (neko_bcknd_device .eq. 1) then
323 call device_add3(a%x_d, b%x_d, c%x_d, size)
324 else
325 call add3(a%x, b%x, c%x, size)
326 end if
327
328 end subroutine field_add3
329
331 subroutine field_add4(a, b, c, d, n)
332 integer, intent(in), optional :: n
333 type(field_t), intent(inout) :: a
334 type(field_t), intent(in) :: b, c, d
335 integer :: size
336
337 if (present(n)) then
338 size = n
339 else
340 size = a%size()
341 end if
342
343 if (neko_bcknd_device .eq. 1) then
344 call device_add4(a%x_d, b%x_d, c%x_d, d%x_d, size)
345 else
346 call add4(a%x, b%x, c%x, d%x, size)
347 end if
348
349 end subroutine field_add4
350
352 subroutine field_sub2(a, b, n)
353 integer, intent(in), optional :: n
354 type(field_t), intent(inout) :: a
355 type(field_t), intent(inout) :: b
356 integer :: size
357
358 if (present(n)) then
359 size = n
360 else
361 size = a%size()
362 end if
363
364 if (neko_bcknd_device .eq. 1) then
365 call device_sub2(a%x_d, b%x_d, size)
366 else
367 call sub2(a%x, b%x, size)
368 end if
369
370 end subroutine field_sub2
371
373 subroutine field_sub3(a, b, c, n)
374 integer, intent(in), optional :: n
375 type(field_t), intent(inout) :: a
376 type(field_t), intent(in) :: b
377 type(field_t), intent(in) :: c
378 integer :: size
379
380 if (present(n)) then
381 size = n
382 else
383 size = a%size()
384 end if
385
386 if (neko_bcknd_device .eq. 1) then
387 call device_sub3(a%x_d, b%x_d, c%x_d, size)
388 else
389 call sub3(a%x, b%x, c%x, size)
390 end if
391
392 end subroutine field_sub3
393
394
397 subroutine field_add2s1(a, b, c1, n)
398 integer, intent(in), optional :: n
399 type(field_t), intent(inout) :: a
400 type(field_t), intent(inout) :: b
401 real(kind=rp), intent(in) :: c1
402 integer :: size
403
404 if (present(n)) then
405 size = n
406 else
407 size = a%size()
408 end if
409
410 if (neko_bcknd_device .eq. 1) then
411 call device_add2s1(a%x_d, b%x_d, c1, size)
412 else
413 call add2s1(a%x, b%x, c1, size)
414 end if
415
416 end subroutine field_add2s1
417
420 subroutine field_add2s2(a, b, c1, n)
421 integer, intent(in), optional :: n
422 type(field_t), intent(inout) :: a
423 type(field_t), intent(inout) :: b
424 real(kind=rp), intent(in) :: c1
425 integer :: size
426
427 if (present(n)) then
428 size = n
429 else
430 size = a%size()
431 end if
432
433 if (neko_bcknd_device .eq. 1) then
434 call device_add2s2(a%x_d, b%x_d, c1, size)
435 else
436 call add2s2(a%x, b%x, c1, size)
437 end if
438
439 end subroutine field_add2s2
440
442 subroutine field_addsqr2s2(a, b, c1, n)
443 integer, intent(in), optional :: n
444 type(field_t), intent(inout) :: a
445 type(field_t), intent(in) :: b
446 real(kind=rp), intent(in) :: c1
447 integer :: size
448
449 if (present(n)) then
450 size = n
451 else
452 size = a%size()
453 end if
454
455 if (neko_bcknd_device .eq. 1) then
456 call device_addsqr2s2(a%x_d, b%x_d, c1, size)
457 else
458 call addsqr2s2(a%x, b%x, c1, size)
459 end if
460
461 end subroutine field_addsqr2s2
462
464 subroutine field_cmult2(a, b, c, n)
465 integer, intent(in), optional :: n
466 type(field_t), intent(inout) :: a
467 type(field_t), intent(in) :: b
468 real(kind=rp), intent(in) :: c
469 integer :: size
470
471 if (present(n)) then
472 size = n
473 else
474 size = a%size()
475 end if
476
477 if (neko_bcknd_device .eq. 1) then
478 call device_cmult2(a%x_d, b%x_d, c, size)
479 else
480 call cmult2(a%x, b%x, c, size)
481 end if
482
483 end subroutine field_cmult2
484
486 subroutine field_invcol2(a, b, n)
487 integer, intent(in), optional :: n
488 type(field_t), intent(inout) :: a
489 type(field_t), intent(in) :: b
490 integer :: size
491
492 if (present(n)) then
493 size = n
494 else
495 size = a%size()
496 end if
497
498 if (neko_bcknd_device .eq. 1) then
499 call device_invcol2(a%x_d, b%x_d, size)
500 else
501 call invcol2(a%x, b%x, size)
502 end if
503
504 end subroutine field_invcol2
505
506
508 subroutine field_col2(a, b, n)
509 integer, intent(in), optional :: n
510 type(field_t), intent(inout) :: a
511 type(field_t), intent(in) :: b
512 integer :: size
513
514 if (present(n)) then
515 size = n
516 else
517 size = a%size()
518 end if
519
520 if (neko_bcknd_device .eq. 1) then
521 call device_col2(a%x_d, b%x_d, size)
522 else
523 call col2(a%x, b%x, size)
524 end if
525
526 end subroutine field_col2
527
529 subroutine field_col3(a, b, c, n)
530 integer, intent(in), optional :: n
531 type(field_t), intent(inout) :: a
532 type(field_t), intent(in) :: b
533 type(field_t), intent(in) :: c
534 integer :: size
535
536 if (present(n)) then
537 size = n
538 else
539 size = a%size()
540 end if
541
542 if (neko_bcknd_device .eq. 1) then
543 call device_col3(a%x_d, b%x_d, c%x_d, size)
544 else
545 call col3(a%x, b%x, c%x, size)
546 end if
547
548 end subroutine field_col3
549
551 subroutine field_subcol3(a, b, c, n)
552 integer, intent(in), optional :: n
553 type(field_t), intent(inout) :: a
554 type(field_t), intent(in) :: b
555 type(field_t), intent(in) :: c
556 integer :: size
557
558 if (present(n)) then
559 size = n
560 else
561 size = a%size()
562 end if
563
564 if (neko_bcknd_device .eq. 1) then
565 call device_subcol3(a%x_d, b%x_d, c%x_d, size)
566 else
567 call subcol3(a%x, b%x, c%x, size)
568 end if
569
570 end subroutine field_subcol3
571
573 subroutine field_add3s2(a, b, c, c1, c2, n)
574 integer, intent(in), optional :: n
575 type(field_t), intent(inout) :: a
576 type(field_t), intent(in) :: b
577 type(field_t), intent(in) :: c
578 real(kind=rp), intent(in) :: c1, c2
579 integer :: size
580
581 if (present(n)) then
582 size = n
583 else
584 size = a%size()
585 end if
586
587 if (neko_bcknd_device .eq. 1) then
588 call device_add3s2(a%x_d, b%x_d, c%x_d, c1, c2, size)
589 else
590 call add3s2(a%x, b%x, c%x, c1, c2, size)
591 end if
592
593 end subroutine field_add3s2
594
596 subroutine field_addcol3(a, b, c, n)
597 integer, intent(in), optional :: n
598 type(field_t), intent(inout) :: a
599 type(field_t), intent(in) :: b
600 type(field_t), intent(in) :: c
601 integer :: size
602
603 if (present(n)) then
604 size = n
605 else
606 size = a%size()
607 end if
608
609 if (neko_bcknd_device .eq. 1) then
610 call device_addcol3(a%x_d, b%x_d, c%x_d, size)
611 else
612 call addcol3(a%x, b%x, c%x, size)
613 end if
614
615 end subroutine field_addcol3
616
618 subroutine field_addcol4(a, b, c, d, n)
619 integer, intent(in), optional :: n
620 type(field_t), intent(inout) :: a
621 type(field_t), intent(in) :: b
622 type(field_t), intent(in) :: c
623 type(field_t), intent(in) :: d
624 integer :: size
625
626 if (present(n)) then
627 size = n
628 else
629 size = a%size()
630 end if
631
632 if (neko_bcknd_device .eq. 1) then
633 call device_addcol4(a%x_d, b%x_d, c%x_d, d%x_d, size)
634 else
635 call addcol4(a%x, b%x, c%x, d%x, size)
636 end if
637
638 end subroutine field_addcol4
639
640 function field_glsum(a, n) result(sum)
641 integer, intent(in), optional :: n
642 type(field_t), intent(in) :: a
643 real(kind=rp) :: sum
644 integer :: size
645
646 if (present(n)) then
647 size = n
648 else
649 size = a%size()
650 end if
651
652 if (neko_bcknd_device .eq. 1) then
653 sum = device_glsum(a%x_d, size)
654 else
655 sum = glsum(a%x, size)
656 end if
657
658 end function field_glsum
659
661 function field_glmax(a, n) result(val)
662 integer, intent(in), optional :: n
663 type(field_t), intent(in) :: a
664 real(kind=rp) :: val
665 integer :: size
666
667 if (present(n)) then
668 size = n
669 else
670 size = a%size()
671 end if
672
673 if (neko_bcknd_device .eq. 1) then
674 val = device_glmax(a%x_d, size)
675 else
676 val = glmax(a%x, size)
677 end if
678
679 end function field_glmax
680
682 function field_glmin(a, n) result(val)
683 integer, intent(in), optional :: n
684 type(field_t), intent(in) :: a
685 real(kind=rp) :: val
686 integer :: size
687
688 if (present(n)) then
689 size = n
690 else
691 size = a%size()
692 end if
693
694 if (neko_bcknd_device .eq. 1) then
695 val = device_glmin(a%x_d, size)
696 else
697 val = glmin(a%x, size)
698 end if
699
700 end function field_glmin
701
702 function field_glsc2(a, b, n) result(norm)
703 integer, intent(in), optional :: n
704 type(field_t), intent(in) :: a, b
705 real(kind=rp) :: norm
706 integer :: size
707
708 if (present(n)) then
709 size = n
710 else
711 size = a%size()
712 end if
713
714 if (neko_bcknd_device .eq. 1) then
715 norm = device_glsc2(a%x_d, b%x_d, size)
716 else
717 norm = glsc2(a%x, b%x, size)
718 end if
719
720 end function field_glsc2
721
722 function field_glsc3(a, b, c, n) result(norm)
723 integer, intent(in), optional :: n
724 type(field_t), intent(in) :: a, b, c
725 real(kind=rp) :: norm
726 integer :: size
727
728 if (present(n)) then
729 size = n
730 else
731 size = a%size()
732 end if
733
734 if (neko_bcknd_device .eq. 1) then
735 norm = device_glsc3(a%x_d, b%x_d, c%x_d, size)
736 else
737 norm = glsc3(a%x, b%x, c%x, size)
738 end if
739
740 end function field_glsc3
741
742 function field_glsubnorm(a, b, n) result(norm)
743 integer, intent(in), optional :: n
744 type(field_t), intent(in) :: a, b
745 real(kind=rp) :: norm
746 integer :: size
747
748 if (present(n)) then
749 size = n
750 else
751 size = a%size()
752 end if
753
754 if (neko_bcknd_device .eq. 1) then
755 norm = device_glsubnorm(a%x_d, b%x_d, size)
756 else
757 norm = glsubnorm(a%x, b%x, size)
758 end if
759
760 end function field_glsubnorm
761
770 subroutine field_masked_gather_copy_0(a, b, mask, n, n_mask)
771 integer, intent(in) :: n, n_mask
772 real(kind=rp), dimension(n_mask), intent(inout) :: a
773 type(field_t) :: b
774 integer, dimension(0:n_mask) :: mask
775 type(c_ptr) :: mask_d, a_d
776
777 if (neko_bcknd_device .eq. 1) then
778 mask_d = device_get_ptr(mask)
779 a_d = device_get_ptr(a)
780 call device_masked_gather_copy_0(a_d, b%x_d, mask_d, n, n_mask)
781 else
782 call masked_gather_copy_0(a, b%x, mask, n, n_mask)
783 end if
784
785 end subroutine field_masked_gather_copy_0
786
795 subroutine field_masked_scatter_copy_0(a, b, mask, n, n_mask)
796 integer, intent(in) :: n, n_mask
797 real(kind=rp), dimension(n_mask), intent(in) :: b
798 type(field_t), intent(inout) :: a
799 integer, dimension(0:n_mask) :: mask
800 type(c_ptr) :: mask_d, b_d
801
802 if (neko_bcknd_device .eq. 1) then
803 mask_d = device_get_ptr(mask)
804 b_d = device_get_ptr(b)
805 call device_masked_scatter_copy_0(a%x_d, b_d, mask_d, n, n_mask)
806 else
807 call masked_scatter_copy_0(a%x, b, mask, n, n_mask)
808 end if
809
810 end subroutine field_masked_scatter_copy_0
811
814 subroutine field_pwmax2(a, b, n)
815 integer, intent(in), optional :: n
816 type(field_t), intent(inout) :: a
817 type(field_t), intent(in) :: b
818 integer :: size, i
819
820 if (present(n)) then
821 size = n
822 else
823 size = a%size()
824 end if
825
826 if (neko_bcknd_device .eq. 1) then
827 call device_pwmax2(a%x_d, b%x_d, size)
828 else
829 call pwmax2(a%x, b%x, size)
830 end if
831
832 end subroutine field_pwmax2
833
836 subroutine field_pwmax3(a, b, c, n)
837 integer, intent(in), optional :: n
838 type(field_t), intent(inout) :: a
839 type(field_t), intent(in) :: b
840 type(field_t), intent(in) :: c
841 integer :: size, i
842
843 if (present(n)) then
844 size = n
845 else
846 size = a%size()
847 end if
848
849 if (neko_bcknd_device .eq. 1) then
850 call device_pwmax3(a%x_d, b%x_d, c%x_d, size)
851 else
852 call pwmax3(a%x, b%x, c%x, size)
853 end if
854 end subroutine field_pwmax3
855
858 subroutine field_cpwmax2(a, b, n)
859 integer, intent(in), optional :: n
860 type(field_t), intent(inout) :: a
861 real(kind=rp), intent(in) :: b
862 integer :: size, i
863
864 if (present(n)) then
865 size = n
866 else
867 size = a%size()
868 end if
869
870 if (neko_bcknd_device .eq. 1) then
871 call device_cpwmax2(a%x_d, b, size)
872 else
873 call cpwmax2(a%x, b, size)
874 end if
875
876 end subroutine field_cpwmax2
877
880 subroutine field_cpwmax3(a, b, c, n)
881 integer, intent(in), optional :: n
882 type(field_t), intent(inout) :: a
883 type(field_t), intent(in) :: b
884 real(kind=rp), intent(in) :: c
885 integer :: size, i
886
887 if (present(n)) then
888 size = n
889 else
890 size = a%size()
891 end if
892
893 if (neko_bcknd_device .eq. 1) then
894 call device_cpwmax3(a%x_d, b%x_d, c, size)
895 else
896 call cpwmax3(a%x, b%x, c, size)
897 end if
898
899 end subroutine field_cpwmax3
900
903 subroutine field_pwmin2(a, b, n)
904 integer, intent(in), optional :: n
905 type(field_t), intent(inout) :: a
906 type(field_t), intent(in) :: b
907 integer :: size, i
908
909 if (present(n)) then
910 size = n
911 else
912 size = a%size()
913 end if
914
915 if (neko_bcknd_device .eq. 1) then
916 call device_pwmin2(a%x_d, b%x_d, size)
917 else
918 call pwmin2(a%x, b%x, size)
919 end if
920
921 end subroutine field_pwmin2
922
925 subroutine field_pwmin3(a, b, c, n)
926 integer, intent(in), optional :: n
927 type(field_t), intent(inout) :: a
928 type(field_t), intent(in) :: b
929 type(field_t), intent(in) :: c
930 integer :: size, i
931
932 if (present(n)) then
933 size = n
934 else
935 size = a%size()
936 end if
937
938 if (neko_bcknd_device .eq. 1) then
939 call device_pwmin3(a%x_d, b%x_d, c%x_d, size)
940 else
941 call pwmin3(a%x, b%x, c%x, size)
942 end if
943 end subroutine field_pwmin3
944
947 subroutine field_cpwmin2(a, b, n)
948 integer, intent(in), optional :: n
949 type(field_t), intent(inout) :: a
950 real(kind=rp), intent(in) :: b
951 integer :: size, i
952
953 if (present(n)) then
954 size = n
955 else
956 size = a%size()
957 end if
958
959 if (neko_bcknd_device .eq. 1) then
960 call device_cpwmin2(a%x_d, b, size)
961 else
962 call cpwmin2(a%x, b, size)
963 end if
964
965 end subroutine field_cpwmin2
966
969 subroutine field_cpwmin3(a, b, c, n)
970 integer, intent(in), optional :: n
971 type(field_t), intent(inout) :: a
972 type(field_t), intent(in) :: b
973 real(kind=rp), intent(in) :: c
974 integer :: size, i
975
976 if (present(n)) then
977 size = n
978 else
979 size = a%size()
980 end if
981
982 if (neko_bcknd_device .eq. 1) then
983 call device_cpwmin3(a%x_d, b%x_d, c, size)
984 else
985 call cpwmin3(a%x, b%x, c, size)
986 end if
987
988 end subroutine field_cpwmin3
989
990end module field_math
Return the device pointer for an associated Fortran array.
Definition device.F90:107
subroutine, public device_pwmin2(a_d, b_d, n, strm)
Compute the point-wise minimum of two vectors .
subroutine, public device_add2s1(a_d, b_d, c1, n, strm)
subroutine, public device_sub3(a_d, b_d, c_d, n, strm)
Vector subtraction .
subroutine, public device_masked_scatter_copy_0(a_d, b_d, mask_d, n, n_mask, strm)
Scatter a masked vector .
subroutine, public device_add2s2(a_d, b_d, c1, n, strm)
Vector addition with scalar multiplication (multiplication on first argument)
real(kind=rp) function, public device_glmax(a_d, n, strm)
Max of a vector of length n.
subroutine, public device_add2(a_d, b_d, n, strm)
Vector addition .
subroutine, public device_addcol3(a_d, b_d, c_d, n, strm)
Returns .
real(kind=rp) function, public device_glsum(a_d, n, strm)
Sum a vector of length n.
subroutine, public device_pwmax3(a_d, b_d, c_d, n, strm)
Compute the point-wise maximum of two vectors .
subroutine, public device_invcol1(a_d, n, strm)
Invert a vector .
subroutine, public device_add3s2(a_d, b_d, c_d, c1, c2, n, strm)
Returns .
subroutine, public device_cpwmax3(a_d, b_d, c, n, strm)
Compute the point-wise maximum of a vector and a scalar .
subroutine, public device_rzero(a_d, n, strm)
Zero a real vector.
subroutine, public device_rone(a_d, n, strm)
Set all elements to one.
subroutine, public device_cmult(a_d, c, n, strm)
Multiplication by constant c .
subroutine, public device_vdot3(dot_d, u1_d, u2_d, u3_d, v1_d, v2_d, v3_d, n, strm)
Compute a dot product (3-d version) assuming vector components etc.
real(kind=rp) function, public device_glsubnorm(a_d, b_d, n, strm)
Returns the norm of the difference of two vectors .
subroutine, public device_sub2(a_d, b_d, n, strm)
Vector substraction .
subroutine, public device_copy(a_d, b_d, n, strm)
Copy a vector .
subroutine, public device_invcol3(a_d, b_d, c_d, n, strm)
Vector division .
subroutine, public device_pwmin3(a_d, b_d, c_d, n, strm)
Compute the point-wise minimum of two vectors .
subroutine, public device_col2(a_d, b_d, n, strm)
Vector multiplication .
subroutine, public device_add4(a_d, b_d, c_d, d_d, n, strm)
subroutine, public device_cpwmax2(a_d, c, n, strm)
Compute the point-wise maximum of a vector and a scalar .
subroutine, public device_subcol3(a_d, b_d, c_d, n, strm)
Returns .
subroutine, public device_masked_gather_copy_0(a_d, b_d, mask_d, n, n_mask, strm)
Gather a masked vector .
subroutine, public device_invcol2(a_d, b_d, n, strm)
Vector division .
subroutine, public device_addsqr2s2(a_d, b_d, c1, n, strm)
Returns .
real(kind=rp) function, public device_glsc3(a_d, b_d, c_d, n, strm)
Weighted inner product .
subroutine, public device_cpwmin3(a_d, b_d, c, n, strm)
Compute the point-wise minimum of a vector and a scalar .
real(kind=rp) function, public device_glsc2(a_d, b_d, n, strm)
Weighted inner product .
subroutine, public device_cmult2(a_d, b_d, c, n, strm)
Multiplication by constant c .
subroutine, public device_col3(a_d, b_d, c_d, n, strm)
Vector multiplication with 3 vectors .
subroutine, public device_addcol4(a_d, b_d, c_d, d_d, n, strm)
Returns .
subroutine, public device_cfill(a_d, c, n, strm)
Set all elements to a constant c .
subroutine, public device_cpwmin2(a_d, c, n, strm)
Compute the point-wise minimum of a vector and a scalar .
subroutine, public device_add3(a_d, b_d, c_d, n, strm)
Vector addition .
real(kind=rp) function, public device_glmin(a_d, n, strm)
Min of a vector of length n.
subroutine, public device_pwmax2(a_d, b_d, n, strm)
Compute the point-wise maximum of two vectors .
Device abstraction, common interface for various accelerators.
Definition device.F90:34
subroutine, public field_cpwmax2(a, b, n)
Point-wise max operation for field and constant .
subroutine, public field_add2s1(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on first argument)
subroutine, public field_cpwmin3(a, b, c, n)
Point-wise min operation for field and constant .
subroutine, public field_invcol2(a, b, n)
Vector division .
subroutine field_add4(a, b, c, d, n)
Vector addition .
subroutine, public field_cadd(a, s, n)
Add a scalar to vector .
subroutine, public field_invcol1(a, n)
Invert a vector .
subroutine, public field_col2(a, b, n)
Vector multiplication .
subroutine, public field_sub2(a, b, n)
Vector substraction .
subroutine, public field_sub3(a, b, c, n)
Vector subtraction .
real(kind=rp) function, public field_glmin(a, n)
Global minimum of all elements in a field .
subroutine, public field_cmult2(a, b, c, n)
Multiplication by constant c .
real(kind=rp) function, public field_glsc2(a, b, n)
subroutine, public field_cfill(a, c, n)
Set all elements to a constant c .
subroutine, public field_add2s2(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on second argument)
real(kind=rp) function, public field_glmax(a, n)
Global maximum of all elements in a field .
real(kind=rp) function, public field_glsum(a, n)
subroutine, public field_rzero(a, n)
Zero a real vector.
subroutine, public field_addcol3(a, b, c, n)
Returns .
real(kind=rp) function, public field_glsubnorm(a, b, n)
subroutine, public field_rone(a, n)
Set all elements to one.
subroutine, public field_invcol3(a, b, c, n)
Invert a vector .
subroutine, public field_add3s2(a, b, c, c1, c2, n)
Returns .
subroutine, public field_add2(a, b, n)
Vector addition .
subroutine, public field_cpwmin2(a, b, n)
Point-wise min operation for field and constant .
subroutine, public field_pwmax2(a, b, n)
Point-wise max operation .
subroutine, public field_pwmin3(a, b, c, n)
Point-wise min operation .
subroutine, public field_col3(a, b, c, n)
Vector multiplication with 3 vectors .
subroutine, public field_copy(a, b, n)
Copy a vector .
subroutine, public field_subcol3(a, b, c, n)
Returns .
subroutine, public field_add3(a, b, c, n)
Vector addition .
subroutine, public field_pwmax3(a, b, c, n)
Point-wise max operation .
subroutine, public field_cpwmax3(a, b, c, n)
Point-wise max operation for field and constant .
subroutine, public field_cmult(a, c, n)
Multiplication by constant c .
real(kind=rp) function, public field_glsc3(a, b, c, n)
subroutine, public field_addcol4(a, b, c, d, n)
Returns .
subroutine, public field_addsqr2s2(a, b, c1, n)
Returns .
subroutine, public field_vdot3(dot, u1, u2, u3, v1, v2, v3, n)
Compute a dot product (3-d version) assuming vector components etc.
subroutine, public field_pwmin2(a, b, n)
Point-wise min operation .
subroutine, public field_masked_scatter_copy_0(a, b, mask, n, n_mask)
Gather a contigous array into a field .
subroutine, public field_masked_gather_copy_0(a, b, mask, n, n_mask)
Gather a field to reduced contigous array .
Defines a field.
Definition field.f90:34
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:411
subroutine, public cmult2(a, b, c, n)
Multiplication by constant c .
Definition math.f90:423
subroutine, public invcol2(a, b, n)
Vector division .
Definition math.f90:840
real(kind=rp) function, public glsc3(a, b, c, n)
Weighted inner product .
Definition math.f90:1067
subroutine, public cadd(a, s, n)
Add a scalar to vector .
Definition math.f90:462
subroutine, public addsqr2s2(a, b, c1, n)
Returns .
Definition math.f90:826
subroutine, public add2s1(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on first argument)
Definition math.f90:797
subroutine, public cpwmin2(a, b, n)
Point-wise minimum of scalar and vector .
Definition math.f90:1458
real(kind=rp) function, public glsc2(a, b, n)
Weighted inner product .
Definition math.f90:1048
subroutine, public masked_scatter_copy_0(a, b, mask, n, n_mask)
Scatter a contigous vector to masked positions in a target array .
Definition math.f90:358
subroutine, public subcol3(a, b, c, n)
Returns .
Definition math.f90:881
subroutine, public rone(a, n)
Set all elements to one.
Definition math.f90:238
subroutine, public cpwmin3(a, b, c, n)
Point-wise minimum of scalar and vector .
Definition math.f90:1470
subroutine, public pwmax3(a, b, c, n)
Point-wise maximum of two vectors .
Definition math.f90:1397
subroutine, public add3(a, b, c, n)
Vector addition .
Definition math.f90:739
real(kind=rp) function, public glsum(a, n)
Sum a vector of length n.
Definition math.f90:499
subroutine, public sub3(a, b, c, n)
Vector subtraction .
Definition math.f90:781
subroutine, public addcol4(a, b, c, d, n)
Returns .
Definition math.f90:972
subroutine, public add2(a, b, n)
Vector addition .
Definition math.f90:726
subroutine, public cfill(a, c, n)
Set all elements to a constant c .
Definition math.f90:487
subroutine, public invcol3(a, b, c, n)
Invert a vector .
Definition math.f90:626
subroutine, public add3s2(a, b, c, c1, c2, n)
Returns .
Definition math.f90:895
subroutine, public pwmax2(a, b, n)
Point-wise maximum of two vectors .
Definition math.f90:1385
subroutine, public pwmin2(a, b, n)
Point-wise minimum of two vectors .
Definition math.f90:1434
subroutine, public masked_gather_copy_0(a, b, mask, n, n_mask)
Gather a masked vector to reduced contigous vector .
Definition math.f90:312
subroutine, public addcol3(a, b, c, n)
Returns .
Definition math.f90:958
subroutine, public invcol1(a, n)
Invert a vector .
Definition math.f90:614
subroutine, public cpwmax3(a, b, c, n)
Point-wise maximum of scalar and vector .
Definition math.f90:1421
subroutine, public col2(a, b, n)
Vector multiplication .
Definition math.f90:854
real(kind=rp) function, public glmax(a, n)
Max of a vector of length n.
Definition math.f90:516
subroutine, public copy(a, b, n)
Copy a vector .
Definition math.f90:249
subroutine, public add4(a, b, c, d, n)
Vector addition .
Definition math.f90:753
subroutine, public col3(a, b, c, n)
Vector multiplication with 3 vectors .
Definition math.f90:867
subroutine, public vdot3(dot, u1, u2, u3, v1, v2, v3, n)
Compute a dot product (3-d version) assuming vector components etc.
Definition math.f90:684
subroutine, public rzero(a, n)
Zero a real vector.
Definition math.f90:205
real(kind=rp) function, public glsubnorm(a, b, n)
Returns the norm of the difference of two vectors .
Definition math.f90:1109
subroutine, public cpwmax2(a, b, n)
Point-wise maximum of scalar and vector .
Definition math.f90:1409
subroutine, public sub2(a, b, n)
Vector substraction .
Definition math.f90:768
subroutine, public pwmin3(a, b, c, n)
Point-wise minimum of two vectors .
Definition math.f90:1446
subroutine, public add2s2(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on second argument)
Definition math.f90:812
real(kind=rp) function, public glmin(a, n)
Min of a vector of length n.
Definition math.f90:546
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12