Neko  0.9.99
A portable framework for high-order spectral element flow simulations
wall_model_fctry.f90
Go to the documentation of this file.
1 ! Copyright (c) 2021-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 submodule(wall_model) wall_model_fctry
34  use vreman, only : vreman_t
35  use spalding, only : spalding_t
36  use rough_log_law, only : rough_log_law_t
37  use utils, only : concat_string_array
38  use json_utils, only : json_get
39  implicit none
40 
41  ! List of all possible types created by the factory routine
42  character(len=20) :: WALLM_KNOWN_TYPES(2) = [character(len=20) :: &
43  "spalding", &
44  "rough_log_law"]
45 
46 contains
47 
56  module subroutine wall_model_factory(object, coef, msk, facet, nu, &
57  json)
58  class(wall_model_t), allocatable, intent(inout) :: object
59  type(coef_t), intent(in) :: coef
60  integer, intent(in) :: msk(:)
61  integer, intent(in) :: facet(:)
62  real(kind=rp), intent(in) :: nu
63  type(json_file), intent(inout) :: json
64  character(len=:), allocatable :: type_name
65  character(len=:), allocatable :: type_string
66  integer :: h_index
67 
68  type_string = concat_string_array(wallm_known_types, &
69  new_line('A') // "- ", prepend = .true.)
70 
71  call json_get(json, "type", type_name)
72 
73  if (trim(type_name) .eq. "spalding") then
74  allocate(spalding_t::object)
75  else if (trim(type_name) .eq. "rough_log_law") then
76  allocate(rough_log_law_t::object)
77  else
78  call neko_error("Unknown wall model type: " // trim(type_name) // &
79  trim(type_name) // ". Known types are: " // type_string)
80  stop
81  end if
82 
83  call json_get(json, "h_index", h_index)
84 
85  ! Initialize
86  call object%init(coef, msk, facet, nu, h_index, json)
87 
88  end subroutine wall_model_factory
89 
90 end submodule wall_model_fctry
Retrieves a parameter by name or throws an error.
Definition: json_utils.f90:45
Utilities for retrieving parameters from the case files.
Definition: json_utils.f90:34
Implements rough_log_law_t.
Implements spalding_t.
Definition: spalding.f90:35
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
Implements vreman_t.
Definition: vreman.f90:35
Implements wall_model_t.
Definition: wall_model.f90:35
Wall model based on the log-law for a rough wall. The formula defining the law is ....
Wall model based on Spalding's law of the wall. Reference: http://dx.doi.org/10.1115/1....
Definition: spalding.f90:53
Implements the Vreman LES model.
Definition: vreman.f90:52