Neko 1.99.1
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, &
80 use, intrinsic :: iso_c_binding, only: c_ptr
81 implicit none
82 private
83
94
95contains
96
98 subroutine field_rzero(a, n)
99 integer, intent(in), optional :: n
100 type(field_t), intent(inout) :: a
101 integer :: size
102
103 if (present(n)) then
104 size = n
105 else
106 size = a%size()
107 end if
108
109 if (neko_bcknd_device .eq. 1) then
110 call device_rzero(a%x_d, size)
111 else
112 call rzero(a%x, size)
113 end if
114 end subroutine field_rzero
115
117 subroutine field_rone(a, n)
118 integer, intent(in), optional :: n
119 type(field_t), intent(inout) :: a
120 integer :: size
121
122 if (present(n)) then
123 size = n
124 else
125 size = a%size()
126 end if
127
128 if (neko_bcknd_device .eq. 1) then
129 call device_rone(a%x_d, size)
130 else
131 call rone(a%x, size)
132 end if
133 end subroutine field_rone
134
136 subroutine field_copy(a, b, n)
137 integer, intent(in), optional :: n
138 type(field_t), intent(in) :: b
139 type(field_t), intent(inout) :: a
140 integer :: size
141
142 if (present(n)) then
143 size = n
144 else
145 size = a%size()
146 end if
147
148 if (neko_bcknd_device .eq. 1) then
149 call device_copy(a%x_d, b%x_d, size)
150 else
151 call copy(a%x, b%x, size)
152 end if
153 end subroutine field_copy
154
156 subroutine field_cmult(a, c, n)
157 integer, intent(in), optional :: n
158 type(field_t), intent(inout) :: a
159 real(kind=rp), intent(in) :: c
160 integer :: size
161
162 if (present(n)) then
163 size = n
164 else
165 size = a%size()
166 end if
167
168 if (neko_bcknd_device .eq. 1) then
169 call device_cmult(a%x_d, c, size)
170 else
171 call cmult(a%x, c, size)
172 end if
173 end subroutine field_cmult
174
176 subroutine field_cadd(a, s, n)
177 integer, intent(in), optional :: n
178 type(field_t), intent(inout) :: a
179 real(kind=rp), intent(in) :: s
180 integer :: size
181
182 if (present(n)) then
183 size = n
184 else
185 size = a%size()
186 end if
187
188 if (neko_bcknd_device .eq. 1) then
189 call device_cadd(a%x_d, s, size)
190 else
191 call cadd(a%x, s, size)
192 end if
193 end subroutine field_cadd
194
196 subroutine field_cfill(a, c, n)
197 integer, intent(in), optional :: n
198 type(field_t), intent(inout) :: a
199 real(kind=rp), intent(in) :: c
200 integer :: size
201
202 if (present(n)) then
203 size = n
204 else
205 size = a%size()
206 end if
207
208 if (neko_bcknd_device .eq. 1) then
209 call device_cfill(a%x_d, c, size)
210 else
211 call cfill(a%x, c, size)
212 end if
213 end subroutine field_cfill
214
216 subroutine field_invcol1(a, n)
217 integer, intent(in), optional :: n
218 type(field_t), intent(inout) :: a
219 integer :: size
220
221 if (present(n)) then
222 size = n
223 else
224 size = a%size()
225 end if
226
227 if (neko_bcknd_device .eq. 1) then
228 call device_invcol1(a%x_d, size)
229 else
230 call invcol1(a%x, size)
231 end if
232
233 end subroutine field_invcol1
234
236 subroutine field_invcol3(a, b, c, n)
237 integer, intent(in), optional :: n
238 type(field_t), intent(inout) :: a
239 type(field_t), intent(in) :: b
240 type(field_t), intent(in) :: c
241 integer :: size
242
243 if (present(n)) then
244 size = n
245 else
246 size = a%size()
247 end if
248
249 if (neko_bcknd_device .eq. 1) then
250 call device_invcol3(a%x_d, b%x_d, c%x_d, size)
251 else
252 call invcol3(a%x, b%x, c%x, size)
253 end if
254
255 end subroutine field_invcol3
256
259 subroutine field_vdot3(dot, u1, u2, u3, v1, v2, v3, n)
260 integer, intent(in), optional :: n
261 type(field_t), intent(in) :: u1, u2, u3
262 type(field_t), intent(in) :: v1, v2, v3
263 type(field_t), intent(out) :: dot
264 integer :: size
265
266 if (present(n)) then
267 size = n
268 else
269 size = dot%size()
270 end if
271
272 if (neko_bcknd_device .eq. 1) then
273 call device_vdot3(dot%x_d, &
274 u1%x_d, u2%x_d, u3%x_d, &
275 v1%x_d, v2%x_d, v3%x_d, &
276 size)
277 else
278 call vdot3(dot%x, &
279 u1%x, u2%x, u3%x, &
280 v1%x, v2%x, v3%x, &
281 size)
282 end if
283
284 end subroutine field_vdot3
285
287 subroutine field_add2(a, b, n)
288 integer, intent(in), optional :: n
289 type(field_t), intent(inout) :: a
290 type(field_t), intent(in) :: b
291 integer :: size
292
293 if (present(n)) then
294 size = n
295 else
296 size = a%size()
297 end if
298
299 if (neko_bcknd_device .eq. 1) then
300 call device_add2(a%x_d, b%x_d, size)
301 else
302 call add2(a%x, b%x, size)
303 end if
304
305 end subroutine field_add2
306
308 subroutine field_add3(a, b, c, n)
309 integer, intent(in), optional :: n
310 type(field_t), intent(inout) :: a
311 type(field_t), intent(in) :: b, c
312 integer :: size
313
314 if (present(n)) then
315 size = n
316 else
317 size = a%size()
318 end if
319
320 if (neko_bcknd_device .eq. 1) then
321 call device_add3(a%x_d, b%x_d, c%x_d, size)
322 else
323 call add3(a%x, b%x, c%x, size)
324 end if
325
326 end subroutine field_add3
327
329 subroutine field_add4(a, b, c, d, n)
330 integer, intent(in), optional :: n
331 type(field_t), intent(inout) :: a
332 type(field_t), intent(in) :: b, c, d
333 integer :: size
334
335 if (present(n)) then
336 size = n
337 else
338 size = a%size()
339 end if
340
341 if (neko_bcknd_device .eq. 1) then
342 call device_add4(a%x_d, b%x_d, c%x_d, d%x_d, size)
343 else
344 call add4(a%x, b%x, c%x, d%x, size)
345 end if
346
347 end subroutine field_add4
348
350 subroutine field_sub2(a, b, n)
351 integer, intent(in), optional :: n
352 type(field_t), intent(inout) :: a
353 type(field_t), intent(inout) :: b
354 integer :: size
355
356 if (present(n)) then
357 size = n
358 else
359 size = a%size()
360 end if
361
362 if (neko_bcknd_device .eq. 1) then
363 call device_sub2(a%x_d, b%x_d, size)
364 else
365 call sub2(a%x, b%x, size)
366 end if
367
368 end subroutine field_sub2
369
371 subroutine field_sub3(a, b, c, n)
372 integer, intent(in), optional :: n
373 type(field_t), intent(inout) :: a
374 type(field_t), intent(in) :: b
375 type(field_t), intent(in) :: c
376 integer :: size
377
378 if (present(n)) then
379 size = n
380 else
381 size = a%size()
382 end if
383
384 if (neko_bcknd_device .eq. 1) then
385 call device_sub3(a%x_d, b%x_d, c%x_d, size)
386 else
387 call sub3(a%x, b%x, c%x, size)
388 end if
389
390 end subroutine field_sub3
391
392
395 subroutine field_add2s1(a, b, c1, n)
396 integer, intent(in), optional :: n
397 type(field_t), intent(inout) :: a
398 type(field_t), intent(inout) :: b
399 real(kind=rp), intent(in) :: c1
400 integer :: size
401
402 if (present(n)) then
403 size = n
404 else
405 size = a%size()
406 end if
407
408 if (neko_bcknd_device .eq. 1) then
409 call device_add2s1(a%x_d, b%x_d, c1, size)
410 else
411 call add2s1(a%x, b%x, c1, size)
412 end if
413
414 end subroutine field_add2s1
415
418 subroutine field_add2s2(a, b, c1, n)
419 integer, intent(in), optional :: n
420 type(field_t), intent(inout) :: a
421 type(field_t), intent(inout) :: b
422 real(kind=rp), intent(in) :: c1
423 integer :: size
424
425 if (present(n)) then
426 size = n
427 else
428 size = a%size()
429 end if
430
431 if (neko_bcknd_device .eq. 1) then
432 call device_add2s2(a%x_d, b%x_d, c1, size)
433 else
434 call add2s2(a%x, b%x, c1, size)
435 end if
436
437 end subroutine field_add2s2
438
440 subroutine field_addsqr2s2(a, b, c1, n)
441 integer, intent(in), optional :: n
442 type(field_t), intent(inout) :: a
443 type(field_t), intent(in) :: b
444 real(kind=rp), intent(in) :: c1
445 integer :: size
446
447 if (present(n)) then
448 size = n
449 else
450 size = a%size()
451 end if
452
453 if (neko_bcknd_device .eq. 1) then
454 call device_addsqr2s2(a%x_d, b%x_d, c1, size)
455 else
456 call addsqr2s2(a%x, b%x, c1, size)
457 end if
458
459 end subroutine field_addsqr2s2
460
462 subroutine field_cmult2(a, b, c, n)
463 integer, intent(in), optional :: n
464 type(field_t), intent(inout) :: a
465 type(field_t), intent(in) :: b
466 real(kind=rp), intent(in) :: c
467 integer :: size
468
469 if (present(n)) then
470 size = n
471 else
472 size = a%size()
473 end if
474
475 if (neko_bcknd_device .eq. 1) then
476 call device_cmult2(a%x_d, b%x_d, c, size)
477 else
478 call cmult2(a%x, b%x, c, size)
479 end if
480
481 end subroutine field_cmult2
482
484 subroutine field_invcol2(a, b, n)
485 integer, intent(in), optional :: n
486 type(field_t), intent(inout) :: a
487 type(field_t), intent(in) :: b
488 integer :: size
489
490 if (present(n)) then
491 size = n
492 else
493 size = a%size()
494 end if
495
496 if (neko_bcknd_device .eq. 1) then
497 call device_invcol2(a%x_d, b%x_d, size)
498 else
499 call invcol2(a%x, b%x, size)
500 end if
501
502 end subroutine field_invcol2
503
504
506 subroutine field_col2(a, b, n)
507 integer, intent(in), optional :: n
508 type(field_t), intent(inout) :: a
509 type(field_t), intent(in) :: b
510 integer :: size
511
512 if (present(n)) then
513 size = n
514 else
515 size = a%size()
516 end if
517
518 if (neko_bcknd_device .eq. 1) then
519 call device_col2(a%x_d, b%x_d, size)
520 else
521 call col2(a%x, b%x, size)
522 end if
523
524 end subroutine field_col2
525
527 subroutine field_col3(a, b, c, n)
528 integer, intent(in), optional :: n
529 type(field_t), intent(inout) :: a
530 type(field_t), intent(in) :: b
531 type(field_t), intent(in) :: c
532 integer :: size
533
534 if (present(n)) then
535 size = n
536 else
537 size = a%size()
538 end if
539
540 if (neko_bcknd_device .eq. 1) then
541 call device_col3(a%x_d, b%x_d, c%x_d, size)
542 else
543 call col3(a%x, b%x, c%x, size)
544 end if
545
546 end subroutine field_col3
547
549 subroutine field_subcol3(a, b, c, n)
550 integer, intent(in), optional :: n
551 type(field_t), intent(inout) :: a
552 type(field_t), intent(in) :: b
553 type(field_t), intent(in) :: c
554 integer :: size
555
556 if (present(n)) then
557 size = n
558 else
559 size = a%size()
560 end if
561
562 if (neko_bcknd_device .eq. 1) then
563 call device_subcol3(a%x_d, b%x_d, c%x_d, size)
564 else
565 call subcol3(a%x, b%x, c%x, size)
566 end if
567
568 end subroutine field_subcol3
569
571 subroutine field_add3s2(a, b, c, c1, c2, n)
572 integer, intent(in), optional :: n
573 type(field_t), intent(inout) :: a
574 type(field_t), intent(in) :: b
575 type(field_t), intent(in) :: c
576 real(kind=rp), intent(in) :: c1, c2
577 integer :: size
578
579 if (present(n)) then
580 size = n
581 else
582 size = a%size()
583 end if
584
585 if (neko_bcknd_device .eq. 1) then
586 call device_add3s2(a%x_d, b%x_d, c%x_d, c1, c2, size)
587 else
588 call add3s2(a%x, b%x, c%x, c1, c2, size)
589 end if
590
591 end subroutine field_add3s2
592
594 subroutine field_addcol3(a, b, c, n)
595 integer, intent(in), optional :: n
596 type(field_t), intent(inout) :: a
597 type(field_t), intent(in) :: b
598 type(field_t), intent(in) :: c
599 integer :: size
600
601 if (present(n)) then
602 size = n
603 else
604 size = a%size()
605 end if
606
607 if (neko_bcknd_device .eq. 1) then
608 call device_addcol3(a%x_d, b%x_d, c%x_d, size)
609 else
610 call addcol3(a%x, b%x, c%x, size)
611 end if
612
613 end subroutine field_addcol3
614
616 subroutine field_addcol4(a, b, c, d, n)
617 integer, intent(in), optional :: n
618 type(field_t), intent(inout) :: a
619 type(field_t), intent(in) :: b
620 type(field_t), intent(in) :: c
621 type(field_t), intent(in) :: d
622 integer :: size
623
624 if (present(n)) then
625 size = n
626 else
627 size = a%size()
628 end if
629
630 if (neko_bcknd_device .eq. 1) then
631 call device_addcol4(a%x_d, b%x_d, c%x_d, d%x_d, size)
632 else
633 call addcol4(a%x, b%x, c%x, d%x, size)
634 end if
635
636 end subroutine field_addcol4
637
638 function field_glsum(a, n) result(sum)
639 integer, intent(in), optional :: n
640 type(field_t), intent(in) :: a
641 real(kind=rp) :: sum
642 integer :: size
643
644 if (present(n)) then
645 size = n
646 else
647 size = a%size()
648 end if
649
650 if (neko_bcknd_device .eq. 1) then
651 sum = device_glsum(a%x_d, size)
652 else
653 sum = glsum(a%x, size)
654 end if
655
656 end function field_glsum
657
658 function field_glsc2(a, b, n) result(norm)
659 integer, intent(in), optional :: n
660 type(field_t), intent(in) :: a, b
661 real(kind=rp) :: norm
662 integer :: size
663
664 if (present(n)) then
665 size = n
666 else
667 size = a%size()
668 end if
669
670 if (neko_bcknd_device .eq. 1) then
671 norm = device_glsc2(a%x_d, b%x_d, size)
672 else
673 norm = glsc2(a%x, b%x, size)
674 end if
675
676 end function field_glsc2
677
678 function field_glsc3(a, b, c, n) result(norm)
679 integer, intent(in), optional :: n
680 type(field_t), intent(in) :: a, b, c
681 real(kind=rp) :: norm
682 integer :: size
683
684 if (present(n)) then
685 size = n
686 else
687 size = a%size()
688 end if
689
690 if (neko_bcknd_device .eq. 1) then
691 norm = device_glsc3(a%x_d, b%x_d, c%x_d, size)
692 else
693 norm = glsc3(a%x, b%x, c%x, size)
694 end if
695
696 end function field_glsc3
697
698 function field_glsubnorm(a, b, n) result(norm)
699 integer, intent(in), optional :: n
700 type(field_t), intent(in) :: a, b
701 real(kind=rp) :: norm
702 integer :: size
703
704 if (present(n)) then
705 size = n
706 else
707 size = a%size()
708 end if
709
710 if (neko_bcknd_device .eq. 1) then
711 norm = device_glsubnorm(a%x_d, b%x_d, size)
712 else
713 norm = glsubnorm(a%x, b%x, size)
714 end if
715
716 end function field_glsubnorm
717
726 subroutine field_masked_gather_copy_0(a, b, mask, n, n_mask)
727 integer, intent(in) :: n, n_mask
728 real(kind=rp), dimension(n_mask), intent(inout) :: a
729 type(field_t) :: b
730 integer, dimension(0:n_mask) :: mask
731 type(c_ptr) :: mask_d, a_d
732
733 if (neko_bcknd_device .eq. 1) then
734 mask_d = device_get_ptr(mask)
735 a_d = device_get_ptr(a)
736 call device_masked_gather_copy_0(a_d, b%x_d, mask_d, n, n_mask)
737 else
738 call masked_gather_copy_0(a, b%x, mask, n, n_mask)
739 end if
740
741 end subroutine field_masked_gather_copy_0
742
751 subroutine field_masked_scatter_copy_0(a, b, mask, n, n_mask)
752 integer, intent(in) :: n, n_mask
753 real(kind=rp), dimension(n_mask), intent(in) :: b
754 type(field_t), intent(inout) :: a
755 integer, dimension(0:n_mask) :: mask
756 type(c_ptr) :: mask_d, b_d
757
758 if (neko_bcknd_device .eq. 1) then
759 mask_d = device_get_ptr(mask)
760 b_d = device_get_ptr(b)
761 call device_masked_scatter_copy_0(a%x_d, b_d, mask_d, n, n_mask)
762 else
763 call masked_scatter_copy_0(a%x, b, mask, n, n_mask)
764 end if
765
766 end subroutine field_masked_scatter_copy_0
767
770 subroutine field_pwmax2(a, b, n)
771 integer, intent(in), optional :: n
772 type(field_t), intent(inout) :: a
773 type(field_t), intent(in) :: b
774 integer :: size, i
775
776 if (present(n)) then
777 size = n
778 else
779 size = a%size()
780 end if
781
782 if (neko_bcknd_device .eq. 1) then
783 call device_pwmax2(a%x_d, b%x_d, size)
784 else
785 call pwmax2(a%x, b%x, size)
786 end if
787
788 end subroutine field_pwmax2
789
792 subroutine field_pwmax3(a, b, c, n)
793 integer, intent(in), optional :: n
794 type(field_t), intent(inout) :: a
795 type(field_t), intent(in) :: b
796 type(field_t), intent(in) :: c
797 integer :: size, i
798
799 if (present(n)) then
800 size = n
801 else
802 size = a%size()
803 end if
804
805 if (neko_bcknd_device .eq. 1) then
806 call device_pwmax3(a%x_d, b%x_d, c%x_d, size)
807 else
808 call pwmax3(a%x, b%x, c%x, size)
809 end if
810 end subroutine field_pwmax3
811
814 subroutine field_cpwmax2(a, b, n)
815 integer, intent(in), optional :: n
816 type(field_t), intent(inout) :: a
817 real(kind=rp), 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_cpwmax2(a%x_d, b, size)
828 else
829 call cpwmax2(a%x, b, size)
830 end if
831
832 end subroutine field_cpwmax2
833
836 subroutine field_cpwmax3(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 real(kind=rp), 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_cpwmax3(a%x_d, b%x_d, c, size)
851 else
852 call cpwmax3(a%x, b%x, c, size)
853 end if
854
855 end subroutine field_cpwmax3
856
859 subroutine field_pwmin2(a, b, n)
860 integer, intent(in), optional :: n
861 type(field_t), intent(inout) :: a
862 type(field_t), intent(in) :: b
863 integer :: size, i
864
865 if (present(n)) then
866 size = n
867 else
868 size = a%size()
869 end if
870
871 if (neko_bcknd_device .eq. 1) then
872 call device_pwmin2(a%x_d, b%x_d, size)
873 else
874 call pwmin2(a%x, b%x, size)
875 end if
876
877 end subroutine field_pwmin2
878
881 subroutine field_pwmin3(a, b, c, n)
882 integer, intent(in), optional :: n
883 type(field_t), intent(inout) :: a
884 type(field_t), intent(in) :: b
885 type(field_t), intent(in) :: c
886 integer :: size, i
887
888 if (present(n)) then
889 size = n
890 else
891 size = a%size()
892 end if
893
894 if (neko_bcknd_device .eq. 1) then
895 call device_pwmin3(a%x_d, b%x_d, c%x_d, size)
896 else
897 call pwmin3(a%x, b%x, c%x, size)
898 end if
899 end subroutine field_pwmin3
900
903 subroutine field_cpwmin2(a, b, n)
904 integer, intent(in), optional :: n
905 type(field_t), intent(inout) :: a
906 real(kind=rp), 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_cpwmin2(a%x_d, b, size)
917 else
918 call cpwmin2(a%x, b, size)
919 end if
920
921 end subroutine field_cpwmin2
922
925 subroutine field_cpwmin3(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 real(kind=rp), 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_cpwmin3(a%x_d, b%x_d, c, size)
940 else
941 call cpwmin3(a%x, b%x, c, size)
942 end if
943
944 end subroutine field_cpwmin3
945
946end module field_math
Return the device pointer for an associated Fortran array.
Definition device.F90:101
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)
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 .
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 .
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_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
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
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12