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