Neko  0.8.1
A portable framework for high-order spectral element flow simulations
les_model.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 !
33 !
35 module les_model
36  use num_types, only : rp
37  use field, only : field_t, field_ptr_t
38  use json_module, only : json_file
40  use dofmap, only : dofmap_t
41  use coefs, only : coef_t
42  use gs_ops, only : gs_op_add
43  use neko_config, only : neko_bcknd_device
45  implicit none
46  private
47 
49  type, abstract, public :: les_model_t
51  type(field_t), pointer :: nut => null()
53  type(field_t), pointer :: delta => null()
55  type(coef_t), pointer :: coef => null()
56  contains
58  procedure, pass(this) :: init_base => les_model_init_base
60  procedure, pass(this) :: free_base => les_model_free_base
62  procedure, pass(this) :: compute_delta => les_model_compute_delta
64  procedure(les_model_init), pass(this), deferred :: init
66  procedure(les_model_free), pass(this), deferred :: free
68  procedure(les_model_compute), pass(this), deferred :: compute
69  end type les_model_t
70 
71  abstract interface
72 
75  subroutine les_model_compute(this, t, tstep)
76  import les_model_t, rp
77  class(les_model_t), intent(inout) :: this
78  real(kind=rp), intent(in) :: t
79  integer, intent(in) :: tstep
80  end subroutine les_model_compute
81  end interface
82 
83  abstract interface
84 
88  subroutine les_model_init(this, dofmap, coef, json)
89  import les_model_t, json_file, dofmap_t, coef_t
90  class(les_model_t), intent(inout) :: this
91  type(coef_t), intent(in) :: coef
92  type(dofmap_t), intent(in) :: dofmap
93  type(json_file), intent(inout) :: json
94  end subroutine les_model_init
95  end interface
96 
97  abstract interface
98 
99  subroutine les_model_free(this)
100  import les_model_t
101  class(les_model_t), intent(inout) :: this
102  end subroutine les_model_free
103  end interface
104 
105 contains
110  subroutine les_model_init_base(this, dofmap, coef, nut_name)
111  class(les_model_t), intent(inout) :: this
112  type(dofmap_t), intent(in) :: dofmap
113  type(coef_t), target, intent(in) :: coef
114  character(len=*), intent(in) :: nut_name
115 
116  if (.not. neko_field_registry%field_exists(trim(nut_name))) then
117  call neko_field_registry%add_field(dofmap, trim(nut_name))
118  end if
119  if (.not. neko_field_registry%field_exists("les_delta")) then
120  call neko_field_registry%add_field(dofmap, "les_delta")
121  end if
122  this%nut => neko_field_registry%get_field(trim(nut_name))
123  this%delta => neko_field_registry%get_field("les_delta")
124  this%coef => coef
125 
126  call this%compute_delta()
127  end subroutine les_model_init_base
128 
130  subroutine les_model_free_base(this)
131  class(les_model_t), intent(inout) :: this
132 
133  nullify(this%nut)
134  nullify(this%coef)
135  end subroutine les_model_free_base
136 
143  subroutine les_model_compute_delta(this)
144  class(les_model_t), intent(inout) :: this
145  integer :: e, i, j, k
146  integer :: im, ip, jm, jp, km, kp
147  real(kind=rp) :: di, dj, dk, ndim_inv
148 
149 
150  do e=1, this%coef%msh%nelv
151  do k=1, this%coef%Xh%lz
152  km = max(1, k-1)
153  kp = min(this%coef%Xh%lz, k+1)
154 
155  do j=1, this%coef%Xh%ly
156  jm = max(1, j-1)
157  jp = min(this%coef%Xh%ly, j+1)
158 
159  do i=1, this%coef%Xh%lx
160  im = max(1, i-1)
161  ip = min(this%coef%Xh%lx, i+1)
162 
163  di = (this%coef%dof%x(ip,j,k,e) - this%coef%dof%x(im,j,k,e))**2 &
164  + (this%coef%dof%y(ip,j,k,e) - this%coef%dof%y(im,j,k,e))**2 &
165  + (this%coef%dof%z(ip,j,k,e) - this%coef%dof%z(im,j,k,e))**2
166 
167  dj = (this%coef%dof%x(i,jp,k,e) - this%coef%dof%x(i,jm,k,e))**2 &
168  + (this%coef%dof%y(i,jp,k,e) - this%coef%dof%y(i,jm,k,e))**2 &
169  + (this%coef%dof%z(i,jp,k,e) - this%coef%dof%z(i,jm,k,e))**2
170 
171  dk = (this%coef%dof%x(i,j,kp,e) - this%coef%dof%x(i,j,km,e))**2 &
172  + (this%coef%dof%y(i,j,kp,e) - this%coef%dof%y(i,j,km,e))**2 &
173  + (this%coef%dof%z(i,j,kp,e) - this%coef%dof%z(i,j,km,e))**2
174 
175  di = sqrt(di) / (ip - im)
176  dj = sqrt(dj) / (jp - jm)
177  dk = sqrt(dk) / (kp - km)
178  this%delta%x(i,j,k,e) = (di * dj * dk)**(1.0_rp / 3.0_rp)
179 
180  enddo
181  enddo
182  enddo
183  enddo
184 
185  call this%coef%gs_h%op(this%delta, gs_op_add)
186 
187  if (neko_bcknd_device .eq. 1) then
188  call device_memcpy(this%delta%x, this%delta%x_d, this%delta%dof%size(),&
189  host_to_device, sync=.false.)
190  end if
191 
192  end subroutine les_model_compute_delta
193 
194 end module les_model
Copy data between host and device (or device and device)
Definition: device.F90:51
Compute eddy viscosity.
Definition: les_model.f90:75
Common constructor.
Definition: les_model.f90:88
Coefficients.
Definition: coef.f90:34
Device abstraction, common interface for various accelerators.
Definition: device.F90:34
integer, parameter, public host_to_device
Definition: device.F90:47
Defines a mapping of the degrees of freedom.
Definition: dofmap.f90:35
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
Defines Gather-scatter operations.
Definition: gs_ops.f90:34
integer, parameter, public gs_op_add
Definition: gs_ops.f90:36
Implements les_model_t.
Definition: les_model.f90:35
subroutine les_model_init_base(this, dofmap, coef, nut_name)
Constructor for the les_model_t (base) class.
Definition: les_model.f90:111
subroutine les_model_compute_delta(this)
Compute the LES lengthscale. For each GLL point, we take the distance between its neighbours in all 3...
Definition: les_model.f90:144
subroutine les_model_free_base(this)
Destructor for the les_model_t (base) class.
Definition: les_model.f90:131
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
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition: coef.f90:54
field_ptr_t, To easily obtain a pointer to a field
Definition: field.f90:80
Base abstract type for LES models based on the Boussinesq approximation.
Definition: les_model.f90:49
#define max(a, b)
Definition: tensor.cu:40