Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
pnpn_res.f90
Go to the documentation of this file.
1! Copyright (c) 2022, 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
40 use space, only : space_t
41 use mesh, only : mesh_t
42 use num_types, only : rp
43 implicit none
44 private
45
47 type, public, abstract :: pnpn_prs_res_t
48 contains
49 procedure(prs_res), nopass, deferred :: compute
50 end type pnpn_prs_res_t
51
53 type, public, abstract :: pnpn_vel_res_t
54 contains
55 procedure(vel_res), nopass, deferred :: compute
56 end type pnpn_vel_res_t
57
58 abstract interface
59 subroutine prs_res(p, p_res, u, v, w, u_e, v_e, w_e, f_x, f_y, f_z, c_xh,&
60 gs_Xh, bc_prs_surface, bc_sym_surface, Ax, bd, dt, mu, rho)
61 import field_t
62 import ax_t
63 import gs_t
64 import facet_normal_t
65 import coef_t
66 import rp
67 type(field_t), intent(inout) :: p, u, v, w
68 type(field_t), intent(in) :: u_e, v_e, w_e
69 type(field_t), intent(inout) :: p_res
71 type(field_t), intent(in) :: f_x, f_y, f_z
72 type(coef_t), intent(inout) :: c_Xh
73 type(gs_t), intent(inout) :: gs_Xh
74 type(facet_normal_t), intent(in) :: bc_prs_surface
75 type(facet_normal_t), intent(in) :: bc_sym_surface
76 class(ax_t), intent(inout) :: Ax
77 real(kind=rp), intent(in) :: bd
78 real(kind=rp), intent(in) :: dt
79 type(field_t), intent(in) :: mu
80 type(field_t), intent(in) :: rho
81 end subroutine prs_res
82 end interface
83
84 abstract interface
85 subroutine vel_res(Ax, u, v, w, u_res, v_res, w_res, &
86 p, f_x, f_y, f_z, c_Xh, msh, Xh, mu, rho, bd, dt, n)
87 import field_t
88 import ax_t
89 import gs_t
90 import facet_normal_t
91 import space_t
92 import coef_t
93 import mesh_t
94 import rp
95 class(ax_t), intent(in) :: Ax
96 type(mesh_t), intent(inout) :: msh
97 type(space_t), intent(inout) :: Xh
98 type(field_t), intent(inout) :: p, u, v, w
99 type(field_t), intent(inout) :: u_res, v_res, w_res
100 type(field_t), intent(in) :: f_x, f_y, f_z
101 type(coef_t), intent(inout) :: c_Xh
102 type(field_t), intent(in) :: mu
103 type(field_t), intent(in) :: rho
104 real(kind=rp), intent(in) :: bd
105 real(kind=rp), intent(in) :: dt
106 integer, intent(in) :: n
107 end subroutine vel_res
108
109 end interface
110
111 interface
112
117 module subroutine pnpn_prs_res_factory(object)
118 class(pnpn_prs_res_t), allocatable, intent(inout) :: object
119 end subroutine pnpn_prs_res_factory
120
125 module subroutine pnpn_vel_res_factory(object)
126 class(pnpn_vel_res_t), allocatable, intent(inout) :: object
127 end subroutine pnpn_vel_res_factory
128
133 module subroutine pnpn_prs_res_stress_factory(object)
134 class(pnpn_prs_res_t), allocatable, intent(inout) :: object
135 end subroutine pnpn_prs_res_stress_factory
136
141 module subroutine pnpn_vel_res_stress_factory(object)
142 class(pnpn_vel_res_t), allocatable, intent(inout) :: object
143 end subroutine pnpn_vel_res_stress_factory
144
145 end interface
146
147 public :: pnpn_prs_res_factory, pnpn_vel_res_factory, &
148 pnpn_prs_res_stress_factory, pnpn_vel_res_stress_factory
149
150end module pnpn_residual
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 Pressure and velocity residuals in the Pn-Pn formulation.
Definition pnpn_res.f90:34
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 pressure residual.
Definition pnpn_res.f90:47
Abstract type to compute velocity residual.
Definition pnpn_res.f90:53
The function space for the SEM solution fields.
Definition space.f90:62