Neko  0.8.1
A portable framework for high-order spectral element flow simulations
boussinesq_source_term.f90
Go to the documentation of this file.
1 ! Copyright (c) 2024, 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 json_module, only : json_file
40  use source_term, only : source_term_t
41  use coefs, only : coef_t
42  use neko_config, only : neko_bcknd_device
43  use utils, only : neko_error
47  implicit none
48  private
49 
58  type, public, extends(source_term_t) :: boussinesq_source_term_t
60  type(field_t), pointer :: s => null()
62  real(kind=rp) :: ref_value = 0
64  real(kind=rp) :: g(3)
66  real(kind=rp) :: beta
67  contains
69  procedure, pass(this) :: init => boussinesq_source_term_init_from_json
71  procedure, pass(this) :: init_from_compenents => &
74  procedure, pass(this) :: free => boussinesq_source_term_free
76  procedure, pass(this) :: compute_ => boussinesq_source_term_compute
78 
79 contains
84  subroutine boussinesq_source_term_init_from_json(this, json, fields, coef)
85  class(boussinesq_source_term_t), intent(inout) :: this
86  type(json_file), intent(inout) :: json
87  type(field_list_t), intent(inout), target :: fields
88  type(coef_t), intent(inout) :: coef
89  real(kind=rp), allocatable :: values(:)
90  real(kind=rp) :: start_time, end_time, ref_value
91  character(len=:), allocatable :: scalar_name
92  real(kind=rp), allocatable :: g(:)
93  real(kind=rp) :: beta
94 
95  if (.not. size(fields%fields) == 3) then
96  call neko_error("Boussinesq term expects 3 fields to work on.")
97  end if
98 
99  call json_get_or_default(json, "scalar_field", start_time, 0.0_rp)
100  call json_get_or_default(json, "end_time", end_time, huge(0.0_rp))
101 
102  call json_get_or_default(json, "scalar_field", scalar_name, "s")
103  call json_get(json, "g", g)
104 
105  if (.not. size(g) == 3) then
106  call neko_error("The gravity vector should have 3 components")
107  end if
108 
109  call json_get(json, "reference_value", ref_value)
110  call json_get_or_default(json, "beta", beta, 1.0_rp/ref_value)
111 
112  call boussinesq_source_term_init_from_components(this, fields, scalar_name,&
113  ref_value, g, beta, coef, start_time, end_time)
114 
116 
127  scalar_name, ref_value, g, beta, coef, start_time, end_time)
128  class(boussinesq_source_term_t), intent(inout) :: this
129  class(field_list_t), intent(inout), target :: fields
130  character(len=*), intent(in) :: scalar_name
131  real(kind=rp), intent(in) :: ref_value
132  real(kind=rp), intent(in) :: g(3)
133  real(kind=rp), intent(in) :: beta
134  type(coef_t) :: coef
135  real(kind=rp), intent(in) :: start_time
136  real(kind=rp), intent(in) :: end_time
137 
138  call this%free()
139  call this%init_base(fields, coef, start_time, end_time)
140 
141  if (.not. neko_field_registry%field_exists(scalar_name)) then
142  call neko_field_registry%add_field(this%fields%fields(1)%f%dof, "s")
143  end if
144  this%s => neko_field_registry%get_field("s")
145 
146  this%ref_value = ref_value
147  this%g = g
148  this%beta = beta
150 
153  class(boussinesq_source_term_t), intent(inout) :: this
154 
155  call this%free_base()
156  nullify(this%s)
157  end subroutine boussinesq_source_term_free
158 
162  subroutine boussinesq_source_term_compute(this, t, tstep)
163  class(boussinesq_source_term_t), intent(inout) :: this
164  real(kind=rp), intent(in) :: t
165  integer, intent(in) :: tstep
166  integer :: n_fields, i, n
167 
168  n_fields = size(this%fields%fields)
169  n = this%fields%fields(1)%f%dof%size()
170 
171  if (neko_bcknd_device .eq. 1) then
172  call boussinesq_source_term_compute_device(this%fields, this%s,&
173  this%ref_value, this%g, this%beta)
174  else
175  call boussinesq_source_term_compute_cpu(this%fields, this%s,&
176  this%ref_value, this%g, this%beta)
177  end if
178  end subroutine boussinesq_source_term_compute
179 
180 end module boussinesq_source_term
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Definition: json_utils.f90:53
Retrieves a parameter by name or throws an error.
Definition: json_utils.f90:44
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.
Implements the device kernel for the boussinesq_source_term_t type.
subroutine, public boussinesq_source_term_compute_device(fields, s, ref_value, g, beta)
Computes the Boussinesq source term on the device.
Implements the boussinesq_source_term_t type.
subroutine boussinesq_source_term_compute(this, t, tstep)
Computes the source term and adds the result to fields.
subroutine boussinesq_source_term_init_from_json(this, json, fields, coef)
The common constructor using a JSON object.
subroutine boussinesq_source_term_init_from_components(this, fields, scalar_name, ref_value, g, beta, coef, start_time, end_time)
The constructor from type components.
subroutine boussinesq_source_term_free(this)
Destructor.
Coefficients.
Definition: coef.f90:34
Defines a registry for storing solution fields.
type(field_registry_t), target, public neko_field_registry
Global field registry.
Defines a field.
Definition: field.f90:34
Utilities for retrieving parameters from the case files.
Definition: json_utils.f90:34
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
Implements the source_term_t type and a wrapper source_term_wrapper_t.
Definition: source_term.f90:34
Utilities.
Definition: utils.f90:35
Bouyancy source term accroding to the Boussinesq approximation.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition: coef.f90:54
field_list_t, To be able to group fields together
Definition: field_list.f90:7
Base abstract type for source terms.
Definition: source_term.f90:44