Neko  0.8.1
A portable framework for high-order spectral element flow simulations
boussinesq_source_term_cpu.f90
Go to the documentation of this file.
1 ! Copyright (c) 2023, 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 !
35  use num_types, only : rp
36  use field_list, only : field_list_t
37  use field, only : field_t
38  use math, only : add2s2, cadd
39  implicit none
40  private
41 
43 
44 contains
45 
52  subroutine boussinesq_source_term_compute_cpu(fields, s, ref_value, g, beta)
53  type(field_list_t), intent(inout) :: fields
54  type(field_t), intent(inout) :: s
55  real(kind=rp), intent(in) :: ref_value
56  real(kind=rp), intent(in) :: g(3)
57  real(kind=rp), intent(in) :: beta
58  integer :: n_fields, i, n
59 
60  n_fields = size(fields%fields)
61  n = fields%fields(1)%f%dof%size()
62 
63  do i=1, n_fields
64  call add2s2(fields%fields(i)%f%x, s%x, g(i)*beta, n)
65  call cadd(fields%fields(i)%f%x, -g(i)*beta*ref_value, n)
66  end do
68 
Implements the cpu kernel for the boussinesq_source_term_t type.
subroutine, public boussinesq_source_term_compute_cpu(fields, s, ref_value, g, beta)
Computes the Boussinesq source term on the cpu.
Defines a field.
Definition: field.f90:34
Definition: math.f90:60
subroutine, public cadd(a, s, n)
Add a scalar to vector .
Definition: math.f90:258
subroutine, public add2s2(a, b, c1, n)
Vector addition with scalar multiplication (multiplication on second argument)
Definition: math.f90:589
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
field_list_t, To be able to group fields together
Definition: field_list.f90:7