Neko  0.8.99
A portable framework for high-order spectral element flow simulations
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 !
60 module field_math
62  use num_types, only: rp
63  use field, only: field_t
64  use math, only: rzero, rone, copy, cmult, cadd, cfill, invcol1, vdot3, add2, &
73  implicit none
74  private
75 
83 
84 contains
85 
87  subroutine field_rzero(a, n)
88  integer, intent(in), optional :: n
89  type(field_t), intent(inout) :: a
90  integer :: size
91 
92  if (present(n)) then
93  size = n
94  else
95  size = a%size()
96  end if
97 
98  if (neko_bcknd_device .eq. 1) then
99  call device_rzero(a%x_d, size)
100  else
101  call rzero(a%x, size)
102  end if
103  end subroutine field_rzero
104 
106  subroutine field_rone(a, n)
107  integer, intent(in), optional :: n
108  type(field_t), intent(inout) :: a
109  integer :: size
110 
111  if (present(n)) then
112  size = n
113  else
114  size = a%size()
115  end if
116 
117  if (neko_bcknd_device .eq. 1) then
118  call device_rone(a%x_d, size)
119  else
120  call rone(a%x, size)
121  end if
122  end subroutine field_rone
123 
125  subroutine field_copy(a, b, n)
126  integer, intent(in), optional :: n
127  type(field_t), intent(in) :: b
128  type(field_t), intent(inout) :: a
129  integer :: size
130 
131  if (present(n)) then
132  size = n
133  else
134  size = a%size()
135  end if
136 
137  if (neko_bcknd_device .eq. 1) then
138  call device_copy(a%x_d, b%x_d, size)
139  else
140  call copy(a%x, b%x, size)
141  end if
142  end subroutine field_copy
143 
145  subroutine field_cmult(a, c, n)
146  integer, intent(in), optional :: n
147  type(field_t), intent(inout) :: a
148  real(kind=rp), intent(in) :: c
149  integer :: size
150 
151  if (present(n)) then
152  size = n
153  else
154  size = a%size()
155  end if
156 
157  if (neko_bcknd_device .eq. 1) then
158  call device_cmult(a%x_d, c, size)
159  else
160  call cmult(a%x, c, size)
161  end if
162  end subroutine field_cmult
163 
165  subroutine field_cadd(a, s, n)
166  integer, intent(in), optional :: n
167  type(field_t), intent(inout) :: a
168  real(kind=rp), intent(in) :: s
169  integer :: size
170 
171  if (present(n)) then
172  size = n
173  else
174  size = a%size()
175  end if
176 
177  if (neko_bcknd_device .eq. 1) then
178  call device_cadd(a%x_d, s, size)
179  else
180  call cadd(a%x, s, size)
181  end if
182  end subroutine field_cadd
183 
185  subroutine field_cfill(a, c, n)
186  integer, intent(in), optional :: n
187  type(field_t), intent(inout) :: a
188  real(kind=rp), intent(in) :: c
189  integer :: size
190 
191  if (present(n)) then
192  size = n
193  else
194  size = a%size()
195  end if
196 
197  if (neko_bcknd_device .eq. 1) then
198  call device_cfill(a%x_d, c, size)
199  else
200  call cfill(a%x, c, size)
201  end if
202  end subroutine field_cfill
203 
205  subroutine field_invcol1(a, n)
206  integer, intent(in), optional :: n
207  type(field_t), intent(inout) :: a
208  integer :: size
209 
210  if (present(n)) then
211  size = n
212  else
213  size = a%size()
214  end if
215 
216  if (neko_bcknd_device .eq. 1) then
217  call device_invcol1(a%x_d, size)
218  else
219  call invcol1(a%x, size)
220  end if
221 
222  end subroutine field_invcol1
223 
226  subroutine field_vdot3(dot, u1, u2, u3, v1, v2, v3, n)
227  integer, intent(in), optional :: n
228  type(field_t), intent(in) :: u1, u2, u3
229  type(field_t), intent(in) :: v1, v2, v3
230  type(field_t), intent(out) :: dot
231  integer :: size
232 
233  if (present(n)) then
234  size = n
235  else
236  size = dot%size()
237  end if
238 
239  if (neko_bcknd_device .eq. 1) then
240  call device_vdot3(dot%x_d, &
241  u1%x_d, u2%x_d, u3%x_d, &
242  v1%x_d, v2%x_d, v3%x_d, &
243  size)
244  else
245  call vdot3(dot%x, &
246  u1%x, u2%x, u3%x, &
247  v1%x, v2%x, v3%x, &
248  size)
249  end if
250 
251  end subroutine field_vdot3
252 
254  subroutine field_add2(a, b, n)
255  integer, intent(in), optional :: n
256  type(field_t), intent(inout) :: a
257  type(field_t), intent(in) :: b
258  integer :: size
259 
260  if (present(n)) then
261  size = n
262  else
263  size = a%size()
264  end if
265 
266  if (neko_bcknd_device .eq. 1) then
267  call device_add2(a%x_d, b%x_d, size)
268  else
269  call add2(a%x, b%x, size)
270  end if
271 
272  end subroutine field_add2
273 
275  subroutine field_sub2(a, b, n)
276  integer, intent(in), optional :: n
277  type(field_t), intent(inout) :: a
278  type(field_t), intent(inout) :: b
279  integer :: size
280 
281  if (present(n)) then
282  size = n
283  else
284  size = a%size()
285  end if
286 
287  if (neko_bcknd_device .eq. 1) then
288  call device_sub2(a%x_d, b%x_d, size)
289  else
290  call sub2(a%x, b%x, size)
291  end if
292 
293  end subroutine field_sub2
294 
296  subroutine field_sub3(a, b, c, n)
297  integer, intent(in), optional :: n
298  type(field_t), intent(inout) :: c
299  type(field_t), intent(inout) :: b
300  type(field_t), intent(out) :: a
301  integer :: size
302 
303  if (present(n)) then
304  size = n
305  else
306  size = a%size()
307  end if
308 
309  if (neko_bcknd_device .eq. 1) then
310  call device_sub3(a%x_d, b%x_d, c%x_d, size)
311  else
312  call sub3(a%x, b%x, c%x, size)
313  end if
314 
315  end subroutine field_sub3
316 
317 
320  subroutine field_add2s1(a, b, c1, n)
321  integer, intent(in), optional :: n
322  type(field_t), intent(inout) :: a
323  type(field_t), intent(inout) :: b
324  real(kind=rp), intent(in) :: c1
325  integer :: size
326 
327  if (present(n)) then
328  size = n
329  else
330  size = a%size()
331  end if
332 
333  if (neko_bcknd_device .eq. 1) then
334  call device_add2s1(a%x_d, b%x_d, c1, size)
335  else
336  call add2s1(a%x, b%x, c1, size)
337  end if
338 
339  end subroutine field_add2s1
340 
343  subroutine field_add2s2(a, b, c1, n)
344  integer, intent(in), optional :: n
345  type(field_t), intent(inout) :: a
346  type(field_t), intent(inout) :: b
347  real(kind=rp), intent(in) :: c1
348  integer :: size
349 
350  if (present(n)) then
351  size = n
352  else
353  size = a%size()
354  end if
355 
356  if (neko_bcknd_device .eq. 1) then
357  call device_add2s2(a%x_d, b%x_d, c1, size)
358  else
359  call add2s2(a%x, b%x, c1, size)
360  end if
361 
362  end subroutine field_add2s2
363 
365  subroutine field_addsqr2s2(a, b, c1, n)
366  integer, intent(in), optional :: n
367  type(field_t), intent(inout) :: a
368  type(field_t), intent(in) :: b
369  real(kind=rp), intent(in) :: c1
370  integer :: size
371 
372  if (present(n)) then
373  size = n
374  else
375  size = a%size()
376  end if
377 
378  if (neko_bcknd_device .eq. 1) then
379  call device_addsqr2s2(a%x_d, b%x_d, c1, size)
380  else
381  call addsqr2s2(a%x, b%x, c1, size)
382  end if
383 
384  end subroutine field_addsqr2s2
385 
387  subroutine field_cmult2(a, b, c, n)
388  integer, intent(in), optional :: n
389  type(field_t), intent(inout) :: a
390  type(field_t), intent(in) :: b
391  real(kind=rp), intent(in) :: c
392  integer :: size
393 
394  if (present(n)) then
395  size = n
396  else
397  size = a%size()
398  end if
399 
400  if (neko_bcknd_device .eq. 1) then
401  call device_cmult2(a%x_d, b%x_d, c, size)
402  else
403  call cmult2(a%x, b%x, c, size)
404  end if
405 
406  end subroutine field_cmult2
407 
409  subroutine field_invcol2(a, b, n)
410  integer, intent(in), optional :: n
411  type(field_t), intent(inout) :: a
412  type(field_t), intent(in) :: b
413  integer :: size
414 
415  if (present(n)) then
416  size = n
417  else
418  size = a%size()
419  end if
420 
421  if (neko_bcknd_device .eq. 1) then
422  call device_invcol2(a%x_d, b%x_d, size)
423  else
424  call invcol2(a%x, b%x, size)
425  end if
426 
427  end subroutine field_invcol2
428 
429 
431  subroutine field_col2(a, b, n)
432  integer, intent(in), optional :: n
433  type(field_t), intent(inout) :: a
434  type(field_t), intent(in) :: b
435  integer :: size
436 
437  if (present(n)) then
438  size = n
439  else
440  size = a%size()
441  end if
442 
443  if (neko_bcknd_device .eq. 1) then
444  call device_col2(a%x_d, b%x_d, size)
445  else
446  call col2(a%x, b%x, size)
447  end if
448 
449  end subroutine field_col2
450 
452  subroutine field_col3(a, b, c, n)
453  integer, intent(in), optional :: n
454  type(field_t), intent(inout) :: a
455  type(field_t), intent(in) :: b
456  type(field_t), intent(in) :: c
457  integer :: size
458 
459  if (present(n)) then
460  size = n
461  else
462  size = a%size()
463  end if
464 
465  if (neko_bcknd_device .eq. 1) then
466  call device_col3(a%x_d, b%x_d, c%x_d, size)
467  else
468  call col3(a%x, b%x, c%x, size)
469  end if
470 
471  end subroutine field_col3
472 
474  subroutine field_subcol3(a, b, c, n)
475  integer, intent(in), optional :: n
476  type(field_t), intent(inout) :: a
477  type(field_t), intent(in) :: b
478  type(field_t), intent(in) :: c
479  integer :: size
480 
481  if (present(n)) then
482  size = n
483  else
484  size = a%size()
485  end if
486 
487  if (neko_bcknd_device .eq. 1) then
488  call device_subcol3(a%x_d, b%x_d, c%x_d, size)
489  else
490  call subcol3(a%x, b%x, c%x, size)
491  end if
492 
493  end subroutine field_subcol3
494 
496  subroutine field_add3s2(a, b, c, c1, c2, n)
497  integer, intent(in), optional :: n
498  type(field_t), intent(inout) :: a
499  type(field_t), intent(in) :: b
500  type(field_t), intent(in) :: c
501  real(kind=rp), intent(in) :: c1, c2
502  integer :: size
503 
504  if (present(n)) then
505  size = n
506  else
507  size = a%size()
508  end if
509 
510  if (neko_bcknd_device .eq. 1) then
511  call device_add3s2(a%x_d, b%x_d, c%x_d, c1, c2, size)
512  else
513  call add3s2(a%x, b%x, c%x, c1, c2, size)
514  end if
515 
516  end subroutine field_add3s2
517 
519  subroutine field_addcol3(a, b, c, n)
520  integer, intent(in), optional :: n
521  type(field_t), intent(inout) :: a
522  type(field_t), intent(in) :: b
523  type(field_t), intent(in) :: c
524  integer :: size
525 
526  if (present(n)) then
527  size = n
528  else
529  size = a%size()
530  end if
531 
532  if (neko_bcknd_device .eq. 1) then
533  call device_addcol3(a%x_d, b%x_d, c%x_d, size)
534  else
535  call addcol3(a%x, b%x, c%x, size)
536  end if
537 
538  end subroutine field_addcol3
539 
541  subroutine field_addcol4(a, b, c, d, n)
542  integer, intent(in), optional :: n
543  type(field_t), intent(inout) :: a
544  type(field_t), intent(in) :: b
545  type(field_t), intent(in) :: c
546  type(field_t), intent(in) :: d
547  integer :: size
548 
549  if (present(n)) then
550  size = n
551  else
552  size = a%size()
553  end if
554 
555  if (neko_bcknd_device .eq. 1) then
556  call device_addcol4(a%x_d, b%x_d, c%x_d, d%x_d, size)
557  else
558  call addcol4(a%x, b%x, c%x, d%x, size)
559  end if
560 
561  end subroutine field_addcol4
562 
563  function field_glsum(a, n) result(sum)
564  integer, intent(in), optional :: n
565  type(field_t), intent(in) :: a
566  real(kind=rp) :: sum
567  integer :: size
568 
569  if (present(n)) then
570  size = n
571  else
572  size = a%size()
573  end if
574 
575  if (neko_bcknd_device .eq. 1) then
576  sum = device_glsum(a%x_d, size)
577  else
578  sum = glsum(a%x, size)
579  end if
580 
581  end function field_glsum
582 
583  function field_glsc2(a, b, n) result(norm)
584  integer, intent(in), optional :: n
585  type(field_t), intent(in) :: a, b
586  real(kind=rp) :: norm
587  integer :: size
588 
589  if (present(n)) then
590  size = n
591  else
592  size = a%size()
593  end if
594 
595  if (neko_bcknd_device .eq. 1) then
596  norm = device_glsc2(a%x_d, b%x_d, size)
597  else
598  norm = glsc2(a%x, b%x, size)
599  end if
600 
601  end function field_glsc2
602 
603  function field_glsc3(a, b, c, n) result(norm)
604  integer, intent(in), optional :: n
605  type(field_t), intent(in) :: a, b, c
606  real(kind=rp) :: norm
607  integer :: size
608 
609  if (present(n)) then
610  size = n
611  else
612  size = a%size()
613  end if
614 
615  if (neko_bcknd_device .eq. 1) then
616  norm = device_glsc3(a%x_d, b%x_d, c%x_d, size)
617  else
618  norm = glsc3(a%x, b%x, c%x, size)
619  end if
620 
621  end function field_glsc3
622 
623 end module field_math
subroutine, public device_add2(a_d, b_d, n)
Vector addition .
subroutine, public device_addcol3(a_d, b_d, c_d, n)
Returns .
subroutine, public device_col2(a_d, b_d, n)
Vector multiplication .
subroutine, public device_add2s1(a_d, b_d, c1, n)
subroutine, public device_rzero(a_d, n)
Zero a real vector.
subroutine, public device_rone(a_d, n)
Set all elements to one.
subroutine, public device_add2s2(a_d, b_d, c1, n)
Vector addition with scalar multiplication (multiplication on first argument)
subroutine, public device_invcol1(a_d, n)
Invert a vector .
subroutine, public device_col3(a_d, b_d, c_d, n)
Vector multiplication with 3 vectors .
subroutine, public device_cadd(a_d, c, n)
Add a scalar to vector .
subroutine, public device_vdot3(dot_d, u1_d, u2_d, u3_d, v1_d, v2_d, v3_d, n)
Compute a dot product (3-d version) assuming vector components etc.
subroutine, public device_cmult2(a_d, b_d, c, n)
Multiplication by constant c .
subroutine, public device_cmult(a_d, c, n)
Multiplication by constant c .
real(kind=rp) function, public device_glsc2(a_d, b_d, n)
Weighted inner product .
subroutine, public device_sub3(a_d, b_d, c_d, n)
Vector subtraction .
real(kind=rp) function, public device_glsc3(a_d, b_d, c_d, n)
Weighted inner product .
real(kind=rp) function, public device_glsum(a_d, n)
Sum a vector of length n.
subroutine, public device_copy(a_d, b_d, n)
Copy a vector .
subroutine, public device_add3s2(a_d, b_d, c_d, c1, c2, n)
Returns .
subroutine, public device_subcol3(a_d, b_d, c_d, n)
Returns .
subroutine, public device_sub2(a_d, b_d, n)
Vector substraction .
subroutine, public device_cfill(a_d, c, n)
Set all elements to a constant c .
subroutine, public device_addcol4(a_d, b_d, c_d, d_d, n)
Returns .
subroutine, public device_invcol2(a_d, b_d, n)
Vector division .
subroutine, public device_addsqr2s2(a_d, b_d, c1, n)
Returns .
subroutine, public field_add2s1(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on first argument)
Definition: field_math.f90:321
subroutine, public field_invcol2(a, b, n)
Vector division .
Definition: field_math.f90:410
subroutine, public field_cadd(a, s, n)
Add a scalar to vector .
Definition: field_math.f90:166
subroutine, public field_invcol1(a, n)
Invert a vector .
Definition: field_math.f90:206
subroutine, public field_col2(a, b, n)
Vector multiplication .
Definition: field_math.f90:432
subroutine, public field_sub2(a, b, n)
Vector substraction .
Definition: field_math.f90:276
subroutine, public field_sub3(a, b, c, n)
Vector subtraction .
Definition: field_math.f90:297
subroutine, public field_cmult2(a, b, c, n)
Multiplication by constant c .
Definition: field_math.f90:388
real(kind=rp) function, public field_glsc2(a, b, n)
Definition: field_math.f90:584
subroutine, public field_cfill(a, c, n)
Set all elements to a constant c .
Definition: field_math.f90:186
subroutine, public field_add2s2(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on second argument)
Definition: field_math.f90:344
real(kind=rp) function, public field_glsum(a, n)
Definition: field_math.f90:564
subroutine, public field_rzero(a, n)
Zero a real vector.
Definition: field_math.f90:88
subroutine, public field_addcol3(a, b, c, n)
Returns .
Definition: field_math.f90:520
subroutine, public field_rone(a, n)
Set all elements to one.
Definition: field_math.f90:107
subroutine, public field_add3s2(a, b, c, c1, c2, n)
Returns .
Definition: field_math.f90:497
subroutine, public field_add2(a, b, n)
Vector addition .
Definition: field_math.f90:255
subroutine, public field_col3(a, b, c, n)
Vector multiplication with 3 vectors .
Definition: field_math.f90:453
subroutine, public field_copy(a, b, n)
Copy a vector .
Definition: field_math.f90:126
subroutine, public field_subcol3(a, b, c, n)
Returns .
Definition: field_math.f90:475
subroutine, public field_cmult(a, c, n)
Multiplication by constant c .
Definition: field_math.f90:146
real(kind=rp) function, public field_glsc3(a, b, c, n)
Definition: field_math.f90:604
subroutine, public field_addcol4(a, b, c, d, n)
Returns .
Definition: field_math.f90:542
subroutine, public field_addsqr2s2(a, b, c1, n)
Returns .
Definition: field_math.f90:366
subroutine, public field_vdot3(dot, u1, u2, u3, v1, v2, v3, n)
Compute a dot product (3-d version) assuming vector components etc.
Definition: field_math.f90:227
Defines a field.
Definition: field.f90:34
Definition: math.f90:60
subroutine, public cmult(a, c, n)
Multiplication by constant c .
Definition: math.f90:277
subroutine, public cmult2(a, b, c, n)
Multiplication by constant c .
Definition: math.f90:661
subroutine, public invcol2(a, b, n)
Vector division .
Definition: math.f90:675
real(kind=rp) function, public glsc3(a, b, c, n)
Weighted inner product .
Definition: math.f90:854
subroutine, public cadd(a, s, n)
Add a scalar to vector .
Definition: math.f90:289
subroutine, public addsqr2s2(a, b, c1, n)
Returns .
Definition: math.f90:647
subroutine, public add2s1(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on first argument)
Definition: math.f90:618
real(kind=rp) function, public glsc2(a, b, n)
Weighted inner product .
Definition: math.f90:836
subroutine, public subcol3(a, b, c, n)
Returns .
Definition: math.f90:716
subroutine, public rone(a, n)
Set all elements to one.
Definition: math.f90:217
real(kind=rp) function, public glsum(a, n)
Sum a vector of length n.
Definition: math.f90:326
subroutine, public sub3(a, b, c, n)
Vector subtraction .
Definition: math.f90:602
subroutine, public addcol4(a, b, c, d, n)
Returns .
Definition: math.f90:775
subroutine, public add2(a, b, n)
Vector addition .
Definition: math.f90:547
subroutine, public cfill(a, c, n)
Set all elements to a constant c .
Definition: math.f90:314
subroutine, public add3s2(a, b, c, c1, c2, n)
Returns .
Definition: math.f90:730
subroutine, public addcol3(a, b, c, n)
Returns .
Definition: math.f90:761
subroutine, public invcol1(a, n)
Invert a vector .
Definition: math.f90:435
subroutine, public col2(a, b, n)
Vector multiplication .
Definition: math.f90:689
subroutine, public copy(a, b, n)
Copy a vector .
Definition: math.f90:228
subroutine, public col3(a, b, c, n)
Vector multiplication with 3 vectors .
Definition: math.f90:702
subroutine, public vdot3(dot, u1, u2, u3, v1, v2, v3, n)
Compute a dot product (3-d version) assuming vector components etc.
Definition: math.f90:505
subroutine, public rzero(a, n)
Zero a real vector.
Definition: math.f90:184
subroutine, public sub2(a, b, n)
Vector substraction .
Definition: math.f90:589
subroutine, public add2s2(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on second argument)
Definition: math.f90:633
Build configurations.
Definition: neko_config.f90:34
integer, parameter neko_bcknd_device
Definition: neko_config.f90:44
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12