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 mask, only : mask_t
65 use device, only : device_get_ptr
66 use math, only : rzero, rone, copy, cmult, cadd, cfill, invcol1, vdot3, add2, &
84 use, intrinsic :: iso_c_binding, only : c_ptr
85 implicit none
86 private
87
100
101contains
102
104 subroutine field_rzero(a, n)
105 integer, intent(in), optional :: n
106 type(field_t), intent(inout) :: a
107 integer :: size
108
109 if (present(n)) then
110 size = n
111 else
112 size = a%size()
113 end if
114
115 if (neko_bcknd_device .eq. 1) then
116 call device_rzero(a%x_d, size)
117 else
118 call rzero(a%x, size)
119 end if
120 end subroutine field_rzero
121
123 subroutine field_rone(a, n)
124 integer, intent(in), optional :: n
125 type(field_t), intent(inout) :: a
126 integer :: size
127
128 if (present(n)) then
129 size = n
130 else
131 size = a%size()
132 end if
133
134 if (neko_bcknd_device .eq. 1) then
135 call device_rone(a%x_d, size)
136 else
137 call rone(a%x, size)
138 end if
139 end subroutine field_rone
140
142 subroutine field_copy(a, b, n)
143 integer, intent(in), optional :: n
144 type(field_t), intent(in) :: b
145 type(field_t), intent(inout) :: a
146 integer :: size
147
148 if (present(n)) then
149 size = n
150 else
151 size = a%size()
152 end if
153
154 if (neko_bcknd_device .eq. 1) then
155 call device_copy(a%x_d, b%x_d, size)
156 else
157 call copy(a%x, b%x, size)
158 end if
159 end subroutine field_copy
160
162 subroutine field_cmult(a, c, n)
163 integer, intent(in), optional :: n
164 type(field_t), intent(inout) :: a
165 real(kind=rp), intent(in) :: c
166 integer :: size
167
168 if (present(n)) then
169 size = n
170 else
171 size = a%size()
172 end if
173
174 if (neko_bcknd_device .eq. 1) then
175 call device_cmult(a%x_d, c, size)
176 else
177 call cmult(a%x, c, size)
178 end if
179 end subroutine field_cmult
180
182 subroutine field_cadd(a, s, n)
183 integer, intent(in), optional :: n
184 type(field_t), intent(inout) :: a
185 real(kind=rp), intent(in) :: s
186 integer :: size
187
188 if (present(n)) then
189 size = n
190 else
191 size = a%size()
192 end if
193
194 if (neko_bcknd_device .eq. 1) then
195 call device_cadd(a%x_d, s, size)
196 else
197 call cadd(a%x, s, size)
198 end if
199 end subroutine field_cadd
200
202 subroutine field_cfill(a, c, n)
203 integer, intent(in), optional :: n
204 type(field_t), intent(inout) :: a
205 real(kind=rp), intent(in) :: c
206 integer :: size
207
208 if (present(n)) then
209 size = n
210 else
211 size = a%size()
212 end if
213
214 if (neko_bcknd_device .eq. 1) then
215 call device_cfill(a%x_d, c, size)
216 else
217 call cfill(a%x, c, size)
218 end if
219 end subroutine field_cfill
220
222 subroutine field_invcol1(a, n)
223 integer, intent(in), optional :: n
224 type(field_t), intent(inout) :: a
225 integer :: size
226
227 if (present(n)) then
228 size = n
229 else
230 size = a%size()
231 end if
232
233 if (neko_bcknd_device .eq. 1) then
234 call device_invcol1(a%x_d, size)
235 else
236 call invcol1(a%x, size)
237 end if
238
239 end subroutine field_invcol1
240
242 subroutine field_invcol3(a, b, c, n)
243 integer, intent(in), optional :: n
244 type(field_t), intent(inout) :: a
245 type(field_t), intent(in) :: b
246 type(field_t), intent(in) :: c
247 integer :: size
248
249 if (present(n)) then
250 size = n
251 else
252 size = a%size()
253 end if
254
255 if (neko_bcknd_device .eq. 1) then
256 call device_invcol3(a%x_d, b%x_d, c%x_d, size)
257 else
258 call invcol3(a%x, b%x, c%x, size)
259 end if
260
261 end subroutine field_invcol3
262
265 subroutine field_vdot3(dot, u1, u2, u3, v1, v2, v3, n)
266 integer, intent(in), optional :: n
267 type(field_t), intent(in) :: u1, u2, u3
268 type(field_t), intent(in) :: v1, v2, v3
269 type(field_t), intent(out) :: dot
270 integer :: size
271
272 if (present(n)) then
273 size = n
274 else
275 size = dot%size()
276 end if
277
278 if (neko_bcknd_device .eq. 1) then
279 call device_vdot3(dot%x_d, &
280 u1%x_d, u2%x_d, u3%x_d, &
281 v1%x_d, v2%x_d, v3%x_d, &
282 size)
283 else
284 call vdot3(dot%x, &
285 u1%x, u2%x, u3%x, &
286 v1%x, v2%x, v3%x, &
287 size)
288 end if
289
290 end subroutine field_vdot3
291
293 subroutine field_add2(a, b, n)
294 integer, intent(in), optional :: n
295 type(field_t), intent(inout) :: a
296 type(field_t), intent(in) :: b
297 integer :: size
298
299 if (present(n)) then
300 size = n
301 else
302 size = a%size()
303 end if
304
305 if (neko_bcknd_device .eq. 1) then
306 call device_add2(a%x_d, b%x_d, size)
307 else
308 call add2(a%x, b%x, size)
309 end if
310
311 end subroutine field_add2
312
314 subroutine field_add3(a, b, c, n)
315 integer, intent(in), optional :: n
316 type(field_t), intent(inout) :: a
317 type(field_t), intent(in) :: b, c
318 integer :: size
319
320 if (present(n)) then
321 size = n
322 else
323 size = a%size()
324 end if
325
326 if (neko_bcknd_device .eq. 1) then
327 call device_add3(a%x_d, b%x_d, c%x_d, size)
328 else
329 call add3(a%x, b%x, c%x, size)
330 end if
331
332 end subroutine field_add3
333
335 subroutine field_add4(a, b, c, d, n)
336 integer, intent(in), optional :: n
337 type(field_t), intent(inout) :: a
338 type(field_t), intent(in) :: b, c, d
339 integer :: size
340
341 if (present(n)) then
342 size = n
343 else
344 size = a%size()
345 end if
346
347 if (neko_bcknd_device .eq. 1) then
348 call device_add4(a%x_d, b%x_d, c%x_d, d%x_d, size)
349 else
350 call add4(a%x, b%x, c%x, d%x, size)
351 end if
352
353 end subroutine field_add4
354
356 subroutine field_sub2(a, b, n)
357 integer, intent(in), optional :: n
358 type(field_t), intent(inout) :: a
359 type(field_t), intent(inout) :: b
360 integer :: size
361
362 if (present(n)) then
363 size = n
364 else
365 size = a%size()
366 end if
367
368 if (neko_bcknd_device .eq. 1) then
369 call device_sub2(a%x_d, b%x_d, size)
370 else
371 call sub2(a%x, b%x, size)
372 end if
373
374 end subroutine field_sub2
375
377 subroutine field_sub3(a, b, c, n)
378 integer, intent(in), optional :: n
379 type(field_t), intent(inout) :: a
380 type(field_t), intent(in) :: b
381 type(field_t), intent(in) :: c
382 integer :: size
383
384 if (present(n)) then
385 size = n
386 else
387 size = a%size()
388 end if
389
390 if (neko_bcknd_device .eq. 1) then
391 call device_sub3(a%x_d, b%x_d, c%x_d, size)
392 else
393 call sub3(a%x, b%x, c%x, size)
394 end if
395
396 end subroutine field_sub3
397
398
401 subroutine field_add2s1(a, b, c1, n)
402 integer, intent(in), optional :: n
403 type(field_t), intent(inout) :: a
404 type(field_t), intent(inout) :: b
405 real(kind=rp), intent(in) :: c1
406 integer :: size
407
408 if (present(n)) then
409 size = n
410 else
411 size = a%size()
412 end if
413
414 if (neko_bcknd_device .eq. 1) then
415 call device_add2s1(a%x_d, b%x_d, c1, size)
416 else
417 call add2s1(a%x, b%x, c1, size)
418 end if
419
420 end subroutine field_add2s1
421
424 subroutine field_add2s2(a, b, c1, n)
425 integer, intent(in), optional :: n
426 type(field_t), intent(inout) :: a
427 type(field_t), intent(inout) :: b
428 real(kind=rp), intent(in) :: c1
429 integer :: size
430
431 if (present(n)) then
432 size = n
433 else
434 size = a%size()
435 end if
436
437 if (neko_bcknd_device .eq. 1) then
438 call device_add2s2(a%x_d, b%x_d, c1, size)
439 else
440 call add2s2(a%x, b%x, c1, size)
441 end if
442
443 end subroutine field_add2s2
444
446 subroutine field_addsqr2s2(a, b, c1, n)
447 integer, intent(in), optional :: n
448 type(field_t), intent(inout) :: a
449 type(field_t), intent(in) :: b
450 real(kind=rp), intent(in) :: c1
451 integer :: size
452
453 if (present(n)) then
454 size = n
455 else
456 size = a%size()
457 end if
458
459 if (neko_bcknd_device .eq. 1) then
460 call device_addsqr2s2(a%x_d, b%x_d, c1, size)
461 else
462 call addsqr2s2(a%x, b%x, c1, size)
463 end if
464
465 end subroutine field_addsqr2s2
466
468 subroutine field_cmult2(a, b, c, n)
469 integer, intent(in), optional :: n
470 type(field_t), intent(inout) :: a
471 type(field_t), intent(in) :: b
472 real(kind=rp), intent(in) :: c
473 integer :: size
474
475 if (present(n)) then
476 size = n
477 else
478 size = a%size()
479 end if
480
481 if (neko_bcknd_device .eq. 1) then
482 call device_cmult2(a%x_d, b%x_d, c, size)
483 else
484 call cmult2(a%x, b%x, c, size)
485 end if
486
487 end subroutine field_cmult2
488
490 subroutine field_invcol2(a, b, n)
491 integer, intent(in), optional :: n
492 type(field_t), intent(inout) :: a
493 type(field_t), intent(in) :: b
494 integer :: size
495
496 if (present(n)) then
497 size = n
498 else
499 size = a%size()
500 end if
501
502 if (neko_bcknd_device .eq. 1) then
503 call device_invcol2(a%x_d, b%x_d, size)
504 else
505 call invcol2(a%x, b%x, size)
506 end if
507
508 end subroutine field_invcol2
509
510
512 subroutine field_col2(a, b, n)
513 integer, intent(in), optional :: n
514 type(field_t), intent(inout) :: a
515 type(field_t), intent(in) :: b
516 integer :: size
517
518 if (present(n)) then
519 size = n
520 else
521 size = a%size()
522 end if
523
524 if (neko_bcknd_device .eq. 1) then
525 call device_col2(a%x_d, b%x_d, size)
526 else
527 call col2(a%x, b%x, size)
528 end if
529
530 end subroutine field_col2
531
533 subroutine field_col3(a, b, c, n)
534 integer, intent(in), optional :: n
535 type(field_t), intent(inout) :: a
536 type(field_t), intent(in) :: b
537 type(field_t), intent(in) :: c
538 integer :: size
539
540 if (present(n)) then
541 size = n
542 else
543 size = a%size()
544 end if
545
546 if (neko_bcknd_device .eq. 1) then
547 call device_col3(a%x_d, b%x_d, c%x_d, size)
548 else
549 call col3(a%x, b%x, c%x, size)
550 end if
551
552 end subroutine field_col3
553
555 subroutine field_subcol3(a, b, c, n)
556 integer, intent(in), optional :: n
557 type(field_t), intent(inout) :: a
558 type(field_t), intent(in) :: b
559 type(field_t), intent(in) :: c
560 integer :: size
561
562 if (present(n)) then
563 size = n
564 else
565 size = a%size()
566 end if
567
568 if (neko_bcknd_device .eq. 1) then
569 call device_subcol3(a%x_d, b%x_d, c%x_d, size)
570 else
571 call subcol3(a%x, b%x, c%x, size)
572 end if
573
574 end subroutine field_subcol3
575
577 subroutine field_add3s2(a, b, c, c1, c2, n)
578 integer, intent(in), optional :: n
579 type(field_t), intent(inout) :: a
580 type(field_t), intent(in) :: b
581 type(field_t), intent(in) :: c
582 real(kind=rp), intent(in) :: c1, c2
583 integer :: size
584
585 if (present(n)) then
586 size = n
587 else
588 size = a%size()
589 end if
590
591 if (neko_bcknd_device .eq. 1) then
592 call device_add3s2(a%x_d, b%x_d, c%x_d, c1, c2, size)
593 else
594 call add3s2(a%x, b%x, c%x, c1, c2, size)
595 end if
596
597 end subroutine field_add3s2
598
600 subroutine field_addcol3(a, b, c, n)
601 integer, intent(in), optional :: n
602 type(field_t), intent(inout) :: a
603 type(field_t), intent(in) :: b
604 type(field_t), intent(in) :: c
605 integer :: size
606
607 if (present(n)) then
608 size = n
609 else
610 size = a%size()
611 end if
612
613 if (neko_bcknd_device .eq. 1) then
614 call device_addcol3(a%x_d, b%x_d, c%x_d, size)
615 else
616 call addcol3(a%x, b%x, c%x, size)
617 end if
618
619 end subroutine field_addcol3
620
622 subroutine field_addcol4(a, b, c, d, n)
623 integer, intent(in), optional :: n
624 type(field_t), intent(inout) :: a
625 type(field_t), intent(in) :: b
626 type(field_t), intent(in) :: c
627 type(field_t), intent(in) :: d
628 integer :: size
629
630 if (present(n)) then
631 size = n
632 else
633 size = a%size()
634 end if
635
636 if (neko_bcknd_device .eq. 1) then
637 call device_addcol4(a%x_d, b%x_d, c%x_d, d%x_d, size)
638 else
639 call addcol4(a%x, b%x, c%x, d%x, size)
640 end if
641
642 end subroutine field_addcol4
643
644 function field_glsum(a, n) result(sum)
645 integer, intent(in), optional :: n
646 type(field_t), intent(in) :: a
647 real(kind=rp) :: sum
648 integer :: size
649
650 if (present(n)) then
651 size = n
652 else
653 size = a%size()
654 end if
655
656 if (neko_bcknd_device .eq. 1) then
657 sum = device_glsum(a%x_d, size)
658 else
659 sum = glsum(a%x, size)
660 end if
661
662 end function field_glsum
663
665 function field_glmax(a, n) result(val)
666 integer, intent(in), optional :: n
667 type(field_t), intent(in) :: a
668 real(kind=rp) :: val
669 integer :: size
670
671 if (present(n)) then
672 size = n
673 else
674 size = a%size()
675 end if
676
677 if (neko_bcknd_device .eq. 1) then
678 val = device_glmax(a%x_d, size)
679 else
680 val = glmax(a%x, size)
681 end if
682
683 end function field_glmax
684
686 function field_glmin(a, n) result(val)
687 integer, intent(in), optional :: n
688 type(field_t), intent(in) :: a
689 real(kind=rp) :: val
690 integer :: size
691
692 if (present(n)) then
693 size = n
694 else
695 size = a%size()
696 end if
697
698 if (neko_bcknd_device .eq. 1) then
699 val = device_glmin(a%x_d, size)
700 else
701 val = glmin(a%x, size)
702 end if
703
704 end function field_glmin
705
706 function field_glsc2(a, b, n) result(norm)
707 integer, intent(in), optional :: n
708 type(field_t), intent(in) :: a, b
709 real(kind=rp) :: norm
710 integer :: size
711
712 if (present(n)) then
713 size = n
714 else
715 size = a%size()
716 end if
717
718 if (neko_bcknd_device .eq. 1) then
719 norm = device_glsc2(a%x_d, b%x_d, size)
720 else
721 norm = glsc2(a%x, b%x, size)
722 end if
723
724 end function field_glsc2
725
726 function field_glsc3(a, b, c, n) result(norm)
727 integer, intent(in), optional :: n
728 type(field_t), intent(in) :: a, b, c
729 real(kind=rp) :: norm
730 integer :: size
731
732 if (present(n)) then
733 size = n
734 else
735 size = a%size()
736 end if
737
738 if (neko_bcknd_device .eq. 1) then
739 norm = device_glsc3(a%x_d, b%x_d, c%x_d, size)
740 else
741 norm = glsc3(a%x, b%x, c%x, size)
742 end if
743
744 end function field_glsc3
745
746 function field_glsubnorm(a, b, n) result(norm)
747 integer, intent(in), optional :: n
748 type(field_t), intent(in) :: a, b
749 real(kind=rp) :: norm
750 integer :: size
751
752 if (present(n)) then
753 size = n
754 else
755 size = a%size()
756 end if
757
758 if (neko_bcknd_device .eq. 1) then
759 norm = device_glsubnorm(a%x_d, b%x_d, size)
760 else
761 norm = glsubnorm(a%x, b%x, size)
762 end if
763
764 end function field_glsubnorm
765
774 subroutine field_masked_gather_copy_0(a, b, mask, n, n_mask)
775 integer, intent(in) :: n, n_mask
776 real(kind=rp), dimension(n_mask), intent(inout) :: a
777 type(field_t) :: b
778 integer, dimension(0:n_mask) :: mask
779 type(c_ptr) :: mask_d, a_d
780
781 if (neko_bcknd_device .eq. 1) then
782 mask_d = device_get_ptr(mask)
783 a_d = device_get_ptr(a)
784 call device_masked_gather_copy_0(a_d, b%x_d, mask_d, n, n_mask)
785 else
786 call masked_gather_copy_0(a, b%x, mask, n, n_mask)
787 end if
788
789 end subroutine field_masked_gather_copy_0
790
797 subroutine field_masked_gather_copy(a, b, mask, n)
798 type(field_t), intent(inout) :: a
799 real(kind=rp), dimension(:), intent(inout) :: b
800 type(mask_t), intent(in) :: mask
801 integer, intent(in) :: n
802 type(c_ptr) :: mask_d, b_d
803
804 if (neko_bcknd_device .eq. 1) then
805 mask_d = mask%get_d()
806 b_d = device_get_ptr(b)
807 call device_masked_gather_copy_aligned(a%x_d, b_d, mask_d, n, &
808 mask%size())
809 else
810 call masked_gather_copy(a%x, b, mask%get(), n, mask%size())
811 end if
812
813 end subroutine field_masked_gather_copy
814
823 subroutine field_masked_scatter_copy_0(a, b, mask, n, n_mask)
824 integer, intent(in) :: n, n_mask
825 real(kind=rp), dimension(n_mask), intent(in) :: b
826 type(field_t), intent(inout) :: a
827 integer, dimension(0:n_mask) :: mask
828 type(c_ptr) :: mask_d, b_d
829
830 if (neko_bcknd_device .eq. 1) then
831 mask_d = device_get_ptr(mask)
832 b_d = device_get_ptr(b)
833 call device_masked_scatter_copy_0(a%x_d, b_d, mask_d, n, n_mask)
834 else
835 call masked_scatter_copy_0(a%x, b, mask, n, n_mask)
836 end if
837
838 end subroutine field_masked_scatter_copy_0
839
842 subroutine field_pwmax2(a, b, n)
843 integer, intent(in), optional :: n
844 type(field_t), intent(inout) :: a
845 type(field_t), intent(in) :: b
846 integer :: size, i
847
848 if (present(n)) then
849 size = n
850 else
851 size = a%size()
852 end if
853
854 if (neko_bcknd_device .eq. 1) then
855 call device_pwmax2(a%x_d, b%x_d, size)
856 else
857 call pwmax2(a%x, b%x, size)
858 end if
859
860 end subroutine field_pwmax2
861
864 subroutine field_pwmax3(a, b, c, n)
865 integer, intent(in), optional :: n
866 type(field_t), intent(inout) :: a
867 type(field_t), intent(in) :: b
868 type(field_t), intent(in) :: c
869 integer :: size, i
870
871 if (present(n)) then
872 size = n
873 else
874 size = a%size()
875 end if
876
877 if (neko_bcknd_device .eq. 1) then
878 call device_pwmax3(a%x_d, b%x_d, c%x_d, size)
879 else
880 call pwmax3(a%x, b%x, c%x, size)
881 end if
882 end subroutine field_pwmax3
883
886 subroutine field_cpwmax2(a, b, n)
887 integer, intent(in), optional :: n
888 type(field_t), intent(inout) :: a
889 real(kind=rp), intent(in) :: b
890 integer :: size, i
891
892 if (present(n)) then
893 size = n
894 else
895 size = a%size()
896 end if
897
898 if (neko_bcknd_device .eq. 1) then
899 call device_cpwmax2(a%x_d, b, size)
900 else
901 call cpwmax2(a%x, b, size)
902 end if
903
904 end subroutine field_cpwmax2
905
908 subroutine field_cpwmax3(a, b, c, n)
909 integer, intent(in), optional :: n
910 type(field_t), intent(inout) :: a
911 type(field_t), intent(in) :: b
912 real(kind=rp), intent(in) :: c
913 integer :: size, i
914
915 if (present(n)) then
916 size = n
917 else
918 size = a%size()
919 end if
920
921 if (neko_bcknd_device .eq. 1) then
922 call device_cpwmax3(a%x_d, b%x_d, c, size)
923 else
924 call cpwmax3(a%x, b%x, c, size)
925 end if
926
927 end subroutine field_cpwmax3
928
931 subroutine field_pwmin2(a, b, n)
932 integer, intent(in), optional :: n
933 type(field_t), intent(inout) :: a
934 type(field_t), intent(in) :: b
935 integer :: size, i
936
937 if (present(n)) then
938 size = n
939 else
940 size = a%size()
941 end if
942
943 if (neko_bcknd_device .eq. 1) then
944 call device_pwmin2(a%x_d, b%x_d, size)
945 else
946 call pwmin2(a%x, b%x, size)
947 end if
948
949 end subroutine field_pwmin2
950
953 subroutine field_pwmin3(a, b, c, n)
954 integer, intent(in), optional :: n
955 type(field_t), intent(inout) :: a
956 type(field_t), intent(in) :: b
957 type(field_t), intent(in) :: c
958 integer :: size, i
959
960 if (present(n)) then
961 size = n
962 else
963 size = a%size()
964 end if
965
966 if (neko_bcknd_device .eq. 1) then
967 call device_pwmin3(a%x_d, b%x_d, c%x_d, size)
968 else
969 call pwmin3(a%x, b%x, c%x, size)
970 end if
971 end subroutine field_pwmin3
972
975 subroutine field_cpwmin2(a, b, n)
976 integer, intent(in), optional :: n
977 type(field_t), intent(inout) :: a
978 real(kind=rp), intent(in) :: b
979 integer :: size, i
980
981 if (present(n)) then
982 size = n
983 else
984 size = a%size()
985 end if
986
987 if (neko_bcknd_device .eq. 1) then
988 call device_cpwmin2(a%x_d, b, size)
989 else
990 call cpwmin2(a%x, b, size)
991 end if
992
993 end subroutine field_cpwmin2
994
997 subroutine field_cpwmin3(a, b, c, n)
998 integer, intent(in), optional :: n
999 type(field_t), intent(inout) :: a
1000 type(field_t), intent(in) :: b
1001 real(kind=rp), intent(in) :: c
1002 integer :: size, i
1003
1004 if (present(n)) then
1005 size = n
1006 else
1007 size = a%size()
1008 end if
1009
1010 if (neko_bcknd_device .eq. 1) then
1011 call device_cpwmin3(a%x_d, b%x_d, c, size)
1012 else
1013 call cpwmin3(a%x, b%x, c, size)
1014 end if
1015
1016 end subroutine field_cpwmin3
1017
1018end 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_masked_gather_copy_aligned(a_d, b_d, mask_d, n, n_mask, strm)
Gather a masked vector .
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 .
subroutine, public field_masked_gather_copy(a, b, mask, n)
Gather a field to reduced contiguous array .
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 contiguous array into a field .
subroutine, public field_masked_gather_copy_0(a, b, mask, n, n_mask)
Gather a field to reduced contiguous 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:412
subroutine, public cmult2(a, b, c, n)
Multiplication by constant c .
Definition math.f90:424
subroutine, public invcol2(a, b, n)
Vector division .
Definition math.f90:841
real(kind=rp) function, public glsc3(a, b, c, n)
Weighted inner product .
Definition math.f90:1068
subroutine, public cadd(a, s, n)
Add a scalar to vector .
Definition math.f90:463
subroutine, public addsqr2s2(a, b, c1, n)
Returns .
Definition math.f90:827
subroutine, public add2s1(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on first argument)
Definition math.f90:798
subroutine, public cpwmin2(a, b, n)
Point-wise minimum of scalar and vector .
Definition math.f90:1459
real(kind=rp) function, public glsc2(a, b, n)
Weighted inner product .
Definition math.f90:1049
subroutine, public masked_scatter_copy_0(a, b, mask, n, n_mask)
Scatter a contigous vector to masked positions in a target array .
Definition math.f90:359
subroutine, public subcol3(a, b, c, n)
Returns .
Definition math.f90:882
subroutine, public rone(a, n)
Set all elements to one.
Definition math.f90:239
subroutine, public cpwmin3(a, b, c, n)
Point-wise minimum of scalar and vector .
Definition math.f90:1471
subroutine, public pwmax3(a, b, c, n)
Point-wise maximum of two vectors .
Definition math.f90:1398
subroutine, public masked_gather_copy(a, b, mask, n, n_mask)
Gather a masked vector to reduced contigous vector .
Definition math.f90:336
subroutine, public add3(a, b, c, n)
Vector addition .
Definition math.f90:740
real(kind=rp) function, public glsum(a, n)
Sum a vector of length n.
Definition math.f90:500
subroutine, public sub3(a, b, c, n)
Vector subtraction .
Definition math.f90:782
subroutine, public addcol4(a, b, c, d, n)
Returns .
Definition math.f90:973
subroutine, public add2(a, b, n)
Vector addition .
Definition math.f90:727
subroutine, public cfill(a, c, n)
Set all elements to a constant c .
Definition math.f90:488
subroutine, public invcol3(a, b, c, n)
Invert a vector .
Definition math.f90:627
subroutine, public add3s2(a, b, c, c1, c2, n)
Returns .
Definition math.f90:896
subroutine, public pwmax2(a, b, n)
Point-wise maximum of two vectors .
Definition math.f90:1386
subroutine, public pwmin2(a, b, n)
Point-wise minimum of two vectors .
Definition math.f90:1435
subroutine, public masked_gather_copy_0(a, b, mask, n, n_mask)
Gather a masked vector to reduced contigous vector .
Definition math.f90:313
subroutine, public addcol3(a, b, c, n)
Returns .
Definition math.f90:959
subroutine, public invcol1(a, n)
Invert a vector .
Definition math.f90:615
subroutine, public cpwmax3(a, b, c, n)
Point-wise maximum of scalar and vector .
Definition math.f90:1422
subroutine, public col2(a, b, n)
Vector multiplication .
Definition math.f90:855
real(kind=rp) function, public glmax(a, n)
Max of a vector of length n.
Definition math.f90:517
subroutine, public copy(a, b, n)
Copy a vector .
Definition math.f90:250
subroutine, public add4(a, b, c, d, n)
Vector addition .
Definition math.f90:754
subroutine, public col3(a, b, c, n)
Vector multiplication with 3 vectors .
Definition math.f90:868
subroutine, public vdot3(dot, u1, u2, u3, v1, v2, v3, n)
Compute a dot product (3-d version) assuming vector components etc.
Definition math.f90:685
subroutine, public rzero(a, n)
Zero a real vector.
Definition math.f90:206
real(kind=rp) function, public glsubnorm(a, b, n)
Returns the norm of the difference of two vectors .
Definition math.f90:1110
subroutine, public cpwmax2(a, b, n)
Point-wise maximum of scalar and vector .
Definition math.f90:1410
subroutine, public sub2(a, b, n)
Vector substraction .
Definition math.f90:769
subroutine, public pwmin3(a, b, c, n)
Point-wise minimum of two vectors .
Definition math.f90:1447
subroutine, public add2s2(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on second argument)
Definition math.f90:813
real(kind=rp) function, public glmin(a, n)
Min of a vector of length n.
Definition math.f90:547
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Type for consistently handling masks in Neko. This type encapsulates the mask array and its associate...
Definition mask.f90:51