Neko  0.9.99
A portable framework for high-order spectral element flow simulations
source_term_fctry.f90
Go to the documentation of this file.
1 ! Copyright (c) 2023-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 !
33 !
35 submodule(source_term) source_term_fctry
40  use json_utils, only : json_get
42  implicit none
43 
44  ! List of all possible types created by the factory routine
45  character(len=20) :: SOURCE_KNOWN_TYPES(4) = [character(len=20) :: &
46  "constant", &
47  "boussinesq", &
48  "coriolis", &
49  "brinkman"]
50 
51 contains
52 
57  module subroutine source_term_factory(object, json, fields, coef)
58  class(source_term_t), allocatable, intent(inout) :: object
59  type(json_file), intent(inout) :: json
60  type(field_list_t), intent(inout) :: fields
61  type(coef_t), intent(inout) :: coef
62  character(len=:), allocatable :: type_name
63  character(len=:), allocatable :: type_string
64 
65  call json_get(json, "type", type_name)
66 
67  if (trim(type_name) .eq. "constant") then
68  allocate(const_source_term_t::object)
69  else if (trim(type_name) .eq. "boussinesq") then
70  allocate(boussinesq_source_term_t::object)
71  else if (trim(type_name) .eq. "coriolis") then
72  allocate(coriolis_source_term_t::object)
73  else if (trim(type_name) .eq. "brinkman") then
74  allocate(brinkman_source_term_t::object)
75  else
76  type_string = concat_string_array(source_known_types, &
77  new_line('A') // "- ", .true.)
78  call neko_error("Unknown source term type: " &
79  // trim(type_name) // ". Known types are: " &
80  // type_string)
81  end if
82 
83  ! Initialize
84  call object%init(json, fields, coef)
85 
86  end subroutine source_term_factory
87 
88 end submodule source_term_fctry
Retrieves a parameter by name or throws an error.
Definition: json_utils.f90:45
Implements the boussinesq_source_term_t type.
Implements the brinkman_source_term_t type.
Implements the const_source_term_t type.
Implements the coriolis_source_term_t type. Maintainer: Timofey Mukha.
Utilities for retrieving parameters from the case files.
Definition: json_utils.f90:34
Implements the source_term_t type and a wrapper source_term_wrapper_t.
Definition: source_term.f90:34
Utilities.
Definition: utils.f90:35
character(:) function, allocatable, public concat_string_array(array, sep, prepend)
Concatenate an array of strings into one string with array items separated by spaces.
Definition: utils.f90:255
Bouyancy source term accroding to the Boussinesq approximation.
A Brinkman source term. The region and strength are controlled by assigning regions types and brinkma...
A constant source term. The strength is specified with the values keyword, which should be an array,...
This source term adds the Coriolis force.