Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
cuda_mma_math.f90
Go to the documentation of this file.
1! Copyright (c) 2025, The Neko-TOP 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!
34 use num_types, only: rp, c_rp
35 use, intrinsic :: iso_c_binding, only: c_int, c_ptr
36
37 implicit none
38 public
39
40 interface
41 subroutine mma_gensub1_cuda(low_d, upp_d, x_d, xmin_d, xmax_d, asyinit, n)&
42 bind(c, name = 'mma_gensub1_cuda')
43 import c_rp, c_int, c_ptr
44 type(c_ptr), value :: low_d, upp_d, x_d, xmin_d, xmax_d
45 real(c_rp) :: asyinit
46 integer(c_int) :: n
47 end subroutine mma_gensub1_cuda
48
49 subroutine mma_gensub2_cuda(low_d, upp_d, x_d, xold1_d, xold2_d, xmin_d, &
50 xmax_d, asydecr, asyincr, n) bind(c, name = 'mma_gensub2_cuda')
51 import c_rp, c_int, c_ptr
52 type(c_ptr), value :: low_d, upp_d, x_d, xold1_d, xold2_d, xmin_d, xmax_d
53 real(c_rp) :: asydecr, asyincr
54 integer(c_int) :: n
55 end subroutine mma_gensub2_cuda
56
57 subroutine mma_gensub3_cuda(x_d, df0dx_d, dfdx_d, low_d, upp_d, min_d, &
58 max_d, alpha_d, beta_d, p0j_d, q0j_d, pij_d, qij_d, n, m) &
59 bind(c, name = 'mma_gensub3_cuda')
60 import c_int, c_ptr
61 type(c_ptr), value :: x_d, df0dx_d, dfdx_d, low_d, upp_d, min_d, max_d, &
62 alpha_d, beta_d, p0j_d, q0j_d, pij_d, qij_d
63 integer(c_int) :: n, m
64 end subroutine mma_gensub3_cuda
65
66 subroutine mma_gensub4_cuda(x_d, low_d, upp_d, pij_d, qij_d, n, m, bi_d) &
67 bind(c, name = 'mma_gensub4_cuda')
68 import c_int, c_ptr
69 type(c_ptr), value :: x_d, low_d, upp_d, pij_d, qij_d, bi_d
70 integer(c_int) :: n, m
71 end subroutine mma_gensub4_cuda
72
73 subroutine cuda_mma_max(xsi_d, x_d, alpha_d, n) &
74 bind(c, name = 'cuda_mma_max')
75 import c_int, c_ptr
76 type(c_ptr), value :: xsi_d, x_d, alpha_d
77 integer(c_int) :: n
78 end subroutine cuda_mma_max
79
80 subroutine cuda_rex(rex_d, x_d, low_d, upp_d, pij_d, p0j_d, qij_d, q0j_d, &
81 lambda_d, xsi_d, eta_d, n, m) bind(c, name = 'cuda_rex')
82 import c_int, c_ptr
83 type(c_ptr), value :: rex_d, x_d, low_d, upp_d, pij_d, p0j_d, qij_d, &
84 q0j_d, lambda_d, xsi_d, eta_d
85 integer(c_int) :: n, m
86 end subroutine cuda_rex
87
88 subroutine cuda_relambda(relambda_d, x_d, upp_d, low_d, pij_d, qij_d, n, &
89 m) bind(c, name = 'cuda_relambda')
90 import c_int, c_ptr
91 type(c_ptr), value :: relambda_d, x_d, upp_d, low_d, pij_d, qij_d
92 integer(c_int) :: n, m
93 end subroutine cuda_relambda
94
95 subroutine cuda_sub2cons2(rexsi_d, xsi_d, x_d, alpha_d, epsi, n) &
96 bind(c, name = 'cuda_sub2cons2')
97 import c_rp, c_int, c_ptr
98 type(c_ptr), value :: rexsi_d, xsi_d, x_d, alpha_d
99 real(c_rp) :: epsi
100 integer(c_int) :: n
101 end subroutine cuda_sub2cons2
102
103 real(c_rp) function cuda_maxval(rex_d, n) bind(c, name = 'cuda_maxval')
104 import c_rp, c_int, c_ptr
105 type(c_ptr), value :: rex_d
106 integer(c_int) :: n
107 end function cuda_maxval
108
109 real(c_rp) function cuda_norm(rex_d, n) bind(c, name = 'cuda_norm')
110 import c_rp, c_int, c_ptr
111 type(c_ptr), value :: rex_d
112 integer(c_int) :: n
113 end function cuda_norm
114
115 subroutine cuda_delx(delx_d, x_d, low_d, upp_d, pij_d, qij_d, p0j_d, &
116 q0j_d, alpha_d, beta_d, lambda_d, epsi, n, m) &
117 bind(c, name = 'cuda_delx')
118 import c_rp, c_int, c_ptr
119 type(c_ptr), value :: delx_d, x_d, low_d, upp_d, pij_d, qij_d, p0j_d, &
120 q0j_d, alpha_d, beta_d, lambda_d
121 real(c_rp) :: epsi
122 integer(c_int) :: n, m
123 end subroutine cuda_delx
124
125
126
127 subroutine cuda_gg(GG_d, x_d, low_d, upp_d, pij_d, qij_d, n, m) &
128 bind(c, name = 'cuda_GG')
129 import c_int, c_ptr
130 type(c_ptr), value :: GG_d, x_d, low_d, upp_d, pij_d, qij_d
131 integer(c_int) :: n, m
132 end subroutine cuda_gg
133
134 subroutine cuda_diagx(diagx_d, x_d, xsi_d, low_d, upp_d, p0j_d, q0j_d, &
135 pij_d, qij_d, alpha_d, beta_d, eta_d, lambda_d, n, m) &
136 bind(c, name = 'cuda_diagx')
137 import c_int, c_ptr
138 type(c_ptr), value :: diagx_d, x_d, xsi_d, low_d, upp_d, p0j_d, q0j_d, &
139 pij_d, qij_d, alpha_d, beta_d, eta_d, lambda_d
140 integer(c_int) :: n, m
141 end subroutine cuda_diagx
142
143 subroutine cuda_bb(bb_d, GG_d, delx_d, diagx_d, n, m) &
144 bind(c, name = 'cuda_bb')
145 import c_int, c_ptr
146 type(c_ptr), value :: bb_d, GG_d, delx_d, diagx_d
147 integer(c_int) :: n, m
148 end subroutine cuda_bb
149
150 subroutine cuda_aa(AA_d, GG_d, diagx_d, n, m) bind(c, name = 'cuda_AA')
151 import c_int, c_ptr
152 type(c_ptr), value :: AA_d, GG_d, diagx_d
153 integer(c_int) :: n, m
154 end subroutine cuda_aa
155
156 subroutine cuda_dx(dx_d, delx_d, diagx_d, GG_d, dlambda_d, n, m) &
157 bind(c, name = 'cuda_dx')
158 import c_int, c_ptr
159 type(c_ptr), value :: dx_d, delx_d, diagx_d, GG_d, dlambda_d
160 integer(c_int) :: n, m
161 end subroutine cuda_dx
162
163 subroutine cuda_dxsi(dxsi_d, xsi_d, dx_d, x_d, alpha_d, epsi, n) &
164 bind(c, name = 'cuda_dxsi')
165 import c_rp, c_int, c_ptr
166 type(c_ptr), value :: dxsi_d, xsi_d, dx_d, x_d, alpha_d
167 real(c_rp) :: epsi
168 integer(c_int) :: n
169 end subroutine cuda_dxsi
170
171 subroutine cuda_deta(deta_d, eta_d, dx_d, x_d, beta_d, epsi, n) &
172 bind(c, name = 'cuda_deta')
173 import c_rp, c_int, c_ptr
174 type(c_ptr), value :: deta_d, eta_d, dx_d, x_d, beta_d
175 real(c_rp) :: epsi
176 integer(c_int) :: n
177 end subroutine cuda_deta
178
179 real(c_rp) function cuda_maxval2(dxx_d, xx_d, cons, n) &
180 bind(c, name = 'cuda_maxval2')
181 import c_rp, c_int, c_ptr
182 type(c_ptr), value :: dxx_d, xx_d
183 real(c_rp) :: cons
184 integer(c_int) :: n
185 end function cuda_maxval2
186
187 real(c_rp) function cuda_maxval3(dx_d, x_d, alpha_d, cons, n) &
188 bind(c, name = 'cuda_maxval3')
189 import c_rp, c_int, c_ptr
190 type(c_ptr), value :: dx_d, x_d, alpha_d
191 real(c_rp) :: cons
192 integer(c_int) :: n
193 end function cuda_maxval3
194
195 subroutine cuda_kkt_rex(rex_d, df0dx_d, dfdx_d, xsi_d, eta_d, lambda_d, &
196 n, m) bind(c, name = 'cuda_kkt_rex')
197 import c_int, c_ptr
198 type(c_ptr), value :: rex_d, df0dx_d, dfdx_d, xsi_d, eta_d, lambda_d
199 integer(c_int) :: n, m
200 end subroutine cuda_kkt_rex
201
202
203 subroutine cuda_maxcons(a_d, b, c, d_d, n) bind(c, name = 'cuda_maxcons')
204 import c_rp, c_int, c_ptr
205 type(c_ptr), value :: a_d, d_d
206 real(c_rp) :: b, c
207 integer(c_int) :: n
208 end subroutine cuda_maxcons
209
210
211 real(c_rp) function cuda_lcsc2(a_d, b_d, n) bind(c, name = 'cuda_lcsc2')
212 import c_rp, c_int, c_ptr
213 type(c_ptr), value :: a_d, b_d
214 integer(c_int) :: n
215 end function cuda_lcsc2
216
217 subroutine cuda_mpisum(a_d, n) bind(c, name = 'cuda_mpisum')
218 import c_int, c_ptr
219 type(c_ptr), value :: a_d
220 integer(c_int) :: n
221 end subroutine cuda_mpisum
222
223 subroutine cuda_add2inv2(a_d, b_d, c, n) bind(c, name = 'cuda_add2inv2')
224 import c_rp, c_int, c_ptr
225 type(c_ptr), value :: a_d, b_d
226 integer(c_int) :: n
227 real(c_rp) :: c
228 end subroutine cuda_add2inv2
229
230 subroutine cuda_max2(a_d, b, c_d, d, n) bind(c, name = 'cuda_max2')
231 import c_rp, c_int, c_ptr
232 type(c_ptr), value :: a_d, c_d
233 integer(c_int) :: n
234 real(c_rp) :: b, d
235 end subroutine cuda_max2
236
237 subroutine cuda_updatebb(bb_d, dellambda_d, dely_d, d_d, mu_d, y_d, delz, &
238 m) bind(c, name = 'cuda_updatebb')
239 import c_rp, c_int, c_ptr
240 type(c_ptr), value :: bb_d, dellambda_d, dely_d, d_d, mu_d, y_d
241 integer(c_int) :: m
242 real(c_rp) :: delz
243 end subroutine cuda_updatebb
244
245 subroutine cuda_updateaa(AA_d, globaltmp_mm_d, s_d, lambda_d, d_d, mu_d, &
246 y_d, a_d, zeta, z, m) bind(c, name = 'cuda_updateAA')
247 import c_rp, c_int, c_ptr
248 type(c_ptr), value :: AA_d, globaltmp_mm_d, s_d, lambda_d, d_d, mu_d, &
249 y_d, a_d
250 integer(c_int) :: m
251 real(c_rp) :: zeta, z
252 end subroutine cuda_updateaa
253
254 subroutine cuda_dy(dy_d, dely_d, dlambda_d, d_d, mu_d, y_d, n) &
255 bind(c, name = 'cuda_dy')
256 import c_int, c_ptr
257 type(c_ptr), value :: dy_d, dely_d, dlambda_d, d_d, mu_d, y_d
258 integer(c_int) :: n
259 end subroutine cuda_dy
260
261 end interface
262
263end module cuda_mma_math