Neko  0.8.99
A portable framework for high-order spectral element flow simulations
scalar_residual.f90
Go to the documentation of this file.
1 ! Copyright (c) 2022-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 gather_scatter, only : gs_t
36  use ax_product, only : ax_t
37  use field, only : field_t
38  use coefs, only : coef_t
39  use source_scalar, only : source_scalar_t
40  use facet_normal, only : facet_normal_t
41  use space, only : space_t
42  use mesh, only : mesh_t
43  use num_types, only : rp
44  implicit none
45  private
46 
48  type, public, abstract :: scalar_residual_t
49  contains
50  procedure(scalar_residual_interface), nopass, deferred :: compute
51  end type scalar_residual_t
52 
53  abstract interface
54 
67  subroutine scalar_residual_interface(Ax, s, s_res, f_Xh, c_Xh, msh, Xh, &
68  lambda, rhocp, bd, dt, n)
69  import field_t
70  import ax_t
71  import gs_t
72  import facet_normal_t
73  import source_scalar_t
74  import space_t
75  import coef_t
76  import mesh_t
77  import rp
78  class(ax_t), intent(in) :: Ax
79  type(mesh_t), intent(inout) :: msh
80  type(space_t), intent(inout) :: Xh
81  type(field_t), intent(inout) :: s
82  type(field_t), intent(inout) :: s_res
83  type(field_t), intent(inout) :: f_Xh
84  type(coef_t), intent(inout) :: c_Xh
85  type(field_t), intent(in) :: lambda
86  real(kind=rp), intent(in) :: rhocp
87  real(kind=rp), intent(in) :: bd
88  real(kind=rp), intent(in) :: dt
89  integer, intent(in) :: n
90  end subroutine scalar_residual_interface
91  end interface
92 
93  interface
94 
97  module subroutine scalar_residual_factory(object)
98  class(scalar_residual_t), allocatable, intent(inout) :: object
99  end subroutine scalar_residual_factory
100  end interface
101 
102  public :: scalar_residual_factory
103 
104 end module scalar_residual
Interface for computing the residual of a scalar transport equation.
Defines a Matrix-vector product.
Definition: ax.f90:34
Coefficients.
Definition: coef.f90:34
Dirichlet condition applied in the facet normal direction.
Defines a field.
Definition: field.f90:34
Gather-scatter.
Defines a mesh.
Definition: mesh.f90:34
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Defines the residual for the scalar transport equation.
Source terms for scalars.
Defines a function space.
Definition: space.f90:34
Base type for a matrix-vector product providing .
Definition: ax.f90:43
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition: coef.f90:55
Dirichlet condition in facet normal direction.
Abstract type to compute scalar residual.
Defines a source term for the scalar transport equation term .
The function space for the SEM solution fields.
Definition: space.f90:62