Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
adjoint_pnpn_res.f90
Go to the documentation of this file.
1
34module adjoint_pnpn_residual
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 facet_normal, only : facet_normal_t
40 use space, only : space_t
41 use mesh, only : mesh_t
42 use num_types, only : rp
43 use, intrinsic :: iso_c_binding, only : c_ptr
44 implicit none
45 private
46
48 type, public, abstract :: adjoint_pnpn_prs_res_t
49 contains
50 procedure(adjoint_prs_res), nopass, deferred :: compute
52
54 type, public, abstract :: adjoint_pnpn_vel_res_t
55 contains
56 procedure(adjoint_vel_res), nopass, deferred :: compute
58
59 abstract interface
60
79 subroutine adjoint_prs_res(p, p_res, u, v, w, f_x, f_y, f_z, c_xh,&
80 gs_Xh, bc_prs_surface, bc_sym_surface, Ax, bd, dt, mu, rho, event)
81 import field_t
82 import ax_t
83 import gs_t
84 import facet_normal_t
85 import coef_t
86 import rp
87 import c_ptr
88 type(field_t), intent(inout) :: p
89 type(field_t), intent(inout) :: u
90 type(field_t), intent(inout) :: v
91 type(field_t), intent(inout) :: w
92 type(field_t), intent(inout) :: p_res
93 type(field_t), intent(in) :: f_x
94 type(field_t), intent(in) :: f_y
95 type(field_t), intent(in) :: f_z
96 type(coef_t), intent(inout) :: c_Xh
97
98 type(gs_t), intent(inout) :: gs_Xh
99
100 type(facet_normal_t), intent(in) :: bc_prs_surface
101
102 type(facet_normal_t), intent(in) :: bc_sym_surface
103
104 class(ax_t), intent(inout) :: Ax
105 real(kind=rp), intent(in) :: bd
106 real(kind=rp), intent(in) :: dt
107 type(field_t), intent(in) :: mu
108 type(field_t), intent(in) :: rho
109 type(c_ptr), intent(inout) :: event
110 end subroutine adjoint_prs_res
111 end interface
112
113 abstract interface
114
134 subroutine adjoint_vel_res(Ax, u, v, w, u_res, v_res, w_res, &
135 p, f_x, f_y, f_z, c_Xh, msh, Xh, mu, rho, bd, dt, n)
136 import field_t
137 import ax_t
138 import space_t
139 import coef_t
140 import mesh_t
141 import rp
142 class(ax_t), intent(in) :: Ax
143 type(mesh_t), intent(inout) :: msh
144 type(space_t), intent(inout) :: Xh
145 type(field_t), intent(inout) :: p
146 type(field_t), intent(inout) :: u
147 type(field_t), intent(inout) :: v
148 type(field_t), intent(inout) :: w
149 type(field_t), intent(inout) :: u_res
150
151 type(field_t), intent(inout) :: v_res
152
153 type(field_t), intent(inout) :: w_res
154
155 type(field_t), intent(in) :: f_x
156 type(field_t), intent(in) :: f_y
157 type(field_t), intent(in) :: f_z
158 type(coef_t), intent(inout) :: c_Xh
159
160 type(field_t), intent(in) :: mu
161 type(field_t), intent(in) :: rho
162 real(kind=rp), intent(in) :: bd
163 real(kind=rp), intent(in) :: dt
164 integer, intent(in) :: n
165 end subroutine adjoint_vel_res
166
167 end interface
168
169 interface
170
175 module subroutine adjoint_pnpn_prs_res_factory(object)
176 class(adjoint_pnpn_prs_res_t), allocatable, intent(inout) :: object
177 end subroutine adjoint_pnpn_prs_res_factory
178
183 module subroutine adjoint_pnpn_vel_res_factory(object)
184 class(adjoint_pnpn_vel_res_t), allocatable, intent(inout) :: object
185 end subroutine adjoint_pnpn_vel_res_factory
186
187 end interface
188
189 public :: adjoint_pnpn_prs_res_factory, adjoint_pnpn_vel_res_factory
190
191end module adjoint_pnpn_residual
Compute adjoint pressure residual.
Abstract type to compute pressure residual.
Abstract type to compute velocity residual.