Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
mean_sqr_field.f90
Go to the documentation of this file.
1! Copyright (c) 2021, The Neko Authors
2! All rights reserved.
3!
4! Redistribution and use in source and binary forms, with or without
5! modification, are permitted provided that the following conditions
6! are met:
7!
8! * Redistributions of source code must retain the above copyright
9! notice, this list of conditions and the following disclaimer.
10!
11! * Redistributions in binary form must reproduce the above
12! copyright notice, this list of conditions and the following
13! disclaimer in the documentation and/or other materials provided
14! with the distribution.
15!
16! * Neither the name of the authors nor the names of its
17! contributors may be used to endorse or promote products derived
18! from this software without specific prior written permission.
19!
20! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31! POSSIBILITY OF SUCH DAMAGE.
32!
34!
36 use num_types, only : rp
38 use mean_field, only : mean_field_t
40 use math, only : addsqr2s2
41 implicit none
42 private
43
44 type, public, extends(mean_field_t) :: mean_sqr_field_t
45 contains
46 procedure, pass(this) :: update => mean_sqr_field_update
47 end type mean_sqr_field_t
48
49contains
50
52 subroutine mean_sqr_field_update(this, k)
53 class(mean_sqr_field_t), intent(inout) :: this
54 real(kind=rp), intent(in) :: k
55
56 if (neko_bcknd_device .eq. 1) then
57 call device_cmult(this%mf%x_d, this%time, size(this%mf%x))
58 call device_addsqr2s2(this%mf%x_d, this%f%x_d, k, size(this%mf%x))
59 this%time = this%time + k
60 call device_cmult(this%mf%x_d, 1.0_rp / this%time, size(this%mf%x))
61 else
62 this%mf%x = this%mf%x * this%time
63 call addsqr2s2(this%mf%x, this%f%x, k, this%mf%dof%size())
64 this%time = this%time + k
65 this%mf%x = this%mf%x / this%time
66 end if
67
68 end subroutine mean_sqr_field_update
69
70end module mean_sqr_field
subroutine, public device_cmult(a_d, c, n)
Multiplication by constant c .
subroutine, public device_addsqr2s2(a_d, b_d, c1, n)
Returns .
Definition math.f90:60
subroutine, public addsqr2s2(a, b, c1, n)
Returns .
Definition math.f90:686
Implements mean_field_t.
Defines a mean square field.
subroutine mean_sqr_field_update(this, k)
Update a mean sqr field.
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Computes the temporal mean of a field.