Neko 1.99.1
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
intf.py
Go to the documentation of this file.
1from ctypes import (
2 CDLL,
3 util,
4 c_char_p,
5 c_int,
6 c_double,
7 c_float,
8 c_bool,
9 c_longlong,
10 POINTER,
11 byref,
12 cast,
13 addressof,
14 create_string_buffer,
15 CFUNCTYPE
16)
17from dataclasses import dataclass
18import json
19import numpy as np
20
21initial_condition = CFUNCTYPE(None, POINTER(c_char_p), c_int)
22
23preprocess = CFUNCTYPE(None, c_double, c_int)
24
25compute = CFUNCTYPE(None, c_double, c_int)
26
27dirichlet_condition = CFUNCTYPE(None, POINTER(c_int), c_int, c_double, c_int)
28
29material_properties = CFUNCTYPE(None, POINTER(c_char_p), c_int, c_double, c_int)
30
31source_term = CFUNCTYPE(None, POINTER(c_char_p), c_int, c_double, c_int)
32
33libneko = CDLL(util.find_library("neko"))
34
38
41
44
45
48
49libneko.neko_case_init.argtypes = [POINTER(c_char_p), c_int, POINTER(c_int)]
51
52libneko.neko_case_free.argtypes = [POINTER(c_int)]
54
55libneko.neko_case_time.argtypes = [POINTER(c_int)]
57
60
61libneko.neko_case_tstep.argtypes = [POINTER(c_int)]
63
64libneko.neko_solve.argtypes = [POINTER(c_int)]
66
67libneko.neko_step.argtypes = [POINTER(c_int)]
69
71libneko.neko_field.restype = POINTER(c_double)
72
75
78
81
83 POINTER(POINTER(c_longlong)),
84 POINTER(POINTER(c_double)),
85 POINTER(POINTER(c_double)),
86 POINTER(POINTER(c_double))]
88
90 POINTER(POINTER(c_longlong)),
91 POINTER(POINTER(c_double)),
92 POINTER(POINTER(c_double)),
93 POINTER(POINTER(c_double)),
94 POINTER(c_int)]
96
98 POINTER(c_int),
99 POINTER(POINTER(c_double)),
100 POINTER(POINTER(c_double)),
101 POINTER(POINTER(c_double)),
102 POINTER(POINTER(c_double)),
103 POINTER(POINTER(c_double)),
104 POINTER(POINTER(c_double)),
105 POINTER(POINTER(c_double)),
106 POINTER(POINTER(c_double)),
107 POINTER(POINTER(c_double)),
108 POINTER(POINTER(c_double))]
110
112 POINTER(c_int),
113 POINTER(POINTER(c_double)),
114 POINTER(POINTER(c_double)),
115 POINTER(POINTER(c_double)),
116 POINTER(POINTER(c_double)),
117 POINTER(POINTER(c_double)),
118 POINTER(POINTER(c_double)),
119 POINTER(POINTER(c_double)),
120 POINTER(POINTER(c_double)),
121 POINTER(POINTER(c_double)),
122 POINTER(POINTER(c_double))]
124
126 POINTER(POINTER(c_double)),
127 POINTER(POINTER(c_double)),
128 POINTER(POINTER(c_double)),
129 POINTER(POINTER(c_double)),
130 POINTER(POINTER(c_double)),
131 POINTER(POINTER(c_double)),
132 POINTER(POINTER(c_double)),
133 POINTER(POINTER(c_double)),
134 POINTER(POINTER(c_double)),
135 POINTER(POINTER(c_double)),
136 POINTER(POINTER(c_double)),
137 POINTER(POINTER(c_double)),
138 POINTER(POINTER(c_double)),
139 POINTER(POINTER(c_double)),
140 POINTER(POINTER(c_double)),
141 POINTER(POINTER(c_double)),
142 POINTER(POINTER(c_double)),
143 POINTER(POINTER(c_double)),
144 POINTER(POINTER(c_double)),
145 POINTER(POINTER(c_double)),
146 POINTER(POINTER(c_double)),
147 POINTER(POINTER(c_double)),
148 POINTER(POINTER(c_double)),
149 POINTER(POINTER(c_double)),
150 POINTER(POINTER(c_double)),
151 POINTER(POINTER(c_double)),
152 POINTER(POINTER(c_double)),
153 POINTER(POINTER(c_double)),
154 POINTER(POINTER(c_double)),
155 POINTER(POINTER(c_double)),
156 POINTER(POINTER(c_double))]
158
159libneko.neko_output_ctrl_execute.argtypes = [POINTER(c_int), c_bool]
160
161libneko.neko_user_setup.argtypes = [POINTER(c_int),
162 initial_condition,
163 preprocess,
164 compute,
165 dirichlet_condition,
166 material_properties,
167 source_term]
169
172
173libneko.neko_cb_field_name_at_index.argtypes = [POINTER(c_int), c_char_p]
175
176
177
178@dataclass
180 """ Defines a Neko mapping of the degrees of freedom """
185 ntot: int
186
187@dataclass
189 """ Defines a Neko space. """
190 lx: int
201
202@dataclass
204 """ Defines a Neko field. """
206 dm: neko_dofmap_t
207 Xh: neko_space_t
208
209@dataclass
244def init():
245 """ Initialise Neko. """
247
249 """ Finalize Neko. """
251
253 """ Initialise Neko device layer. """
255
257 """ Finalize Neko device lyaer. """
259
261 """ Display job information. """
263
265 """ Initialise a field registry. """
267
269 """ Destroy a field registry. """
271
272def case_init(case_json,
273 cb_initial_condition = initial_condition(0),
274 cb_preprocess = preprocess(0),
275 cb_compute = compute(0),
276 cb_dirichlet_condition = dirichlet_condition(0),
277 cb_material_properties = material_properties(0),
278 cb_source_term = source_term(0)):
279 """Initialise a Neko case
280
281 A case is created based on a JSON file and a couple of (optional)
282 callback functions.
283
284 """
285 cp = python_dict_to_fortran(case_json)
286 case_descr = c_int()
287
288 if (cb_initial_condition != initial_condition(0) or
289 cb_preprocess != preprocess(0) or
290 cb_compute != compute(0) or
291 cb_dirichlet_condition != dirichlet_condition(0) or
292 cb_material_properties != material_properties(0) or
293 cb_source_term != source_term(0)):
294
295 libneko.neko_case_allocate(byref(case_descr))
296 user_setup(case_descr, cb_initial_condition, cb_preprocess,
297 cb_compute, cb_dirichlet_condition, cb_material_properties,
298 cb_source_term)
299
300 libneko.neko_case_init(byref(cp), len(json.dumps(case_json)),
301 byref(case_descr))
302 return case_descr
303
304def case_free(case_descr):
305 """ Destroy a Neko case. """
306 libneko.neko_case_free(byref(case_descr))
307
308def time(case_descr):
309 """ Retrive the current time of a case. """
310 time = libneko.neko_case_time(byref(case_descr))
311 return time
312
313def end_time(case_descr):
314 """ Retrive the end time of a case. """
315 end_time = libneko.neko_case_end_time(byref(case_descr))
316 return end_time
317
318def tstep(case_descr):
319 """ Retrive the current time-step of a case. """
320 tstep = libneko.neko_case_tstep(byref(case_descr))
321 return tstep
322
323def solve(case_descr):
324 """ Solve a Neko case. """
325 libneko.neko_solve(byref(case_descr))
326
327def step(case_descr):
328 """ Compute a time-step for a Neko case. """
329 libneko.neko_step(byref(case_descr))
330
331def field(name):
332 """ Retrive a Neko field as a `neko_field_t`. """
333 a = neko_field_t(
334 np.frombuffer((c_double *
336 addressof(libneko.neko_field(name).contents)),
337 dtype=np.dtype(c_double)),
338 field_dofmap(name),
339 field_space(name)
340 )
341 return a
342
343def field_order(name):
344 """ Retrive the order of a Neko field. """
345 order = libneko.neko_field_order(name)
346 return order
347
348def field_nelements(name):
349 """ Retrive the number of elements in a Neko field. """
350 nelements = libneko.neko_field_nelements(name)
351 return nelements
352
353def field_size(name):
354 """ Retrive the total number of degrees of freedom of a Neko field. """
355 size = libneko.neko_field_size(name)
356 return size
357
358def field_dofmap(name):
359 """ Retrive a dofmap of a Neko field as a `neko_dofmap_t`. """
360 size = libneko.neko_field_size(name)
361 dof_ptr = POINTER(c_longlong)()
362 x_ptr = POINTER(c_double)()
363 y_ptr = POINTER(c_double)()
364 z_ptr = POINTER(c_double)()
365 libneko.neko_field_dofmap(name, byref(dof_ptr),
366 byref(x_ptr), byref(y_ptr), byref(z_ptr))
367
368 dofmap = wrap_dofmap(dof_ptr, x_ptr, y_ptr, z_ptr, size)
369
370 return dofmap
371
372def case_fluid_dofmap(case_descr):
373 """ Retrive a dofmap of a Neko case's fluid solver as a `neko_dofmap_t`. """
374 size = c_int()
375 dof_ptr = POINTER(c_longlong)()
376 x_ptr = POINTER(c_double)()
377 y_ptr = POINTER(c_double)()
378 z_ptr = POINTER(c_double)()
379 libneko.neko_case_fluid_dofmap(byref(case_descr), byref(dof_ptr),
380 byref(x_ptr), byref(y_ptr), byref(z_ptr),
381 byref(size))
382
383 dofmap = wrap_dofmap(dof_ptr, x_ptr, y_ptr, z_ptr, size.value)
384
385 return dofmap
386
387def wrap_dofmap(dof_ptr, x_ptr, y_ptr, z_ptr, size):
388 dofmap = neko_dofmap_t(
389 np.frombuffer((c_longlong * size).from_address(
390 addressof(dof_ptr.contents)), dtype=np.dtype(int)),
391 np.frombuffer((c_double * size).from_address(
392 addressof(x_ptr.contents)), dtype=np.dtype(c_double)),
393 np.frombuffer((c_double * size).from_address(
394 addressof(y_ptr.contents)), dtype=np.dtype(c_double)),
395 np.frombuffer((c_double * size).from_address(
396 addressof(z_ptr.contents)), dtype=np.dtype(c_double)),
397 size
398 )
399 return dofmap
400
401def field_space(name):
402 """ Retrive a space of a Neko field as a `neko_space_t`. """
403 lx = c_int()
404 zg = POINTER(c_double)()
405 dr_inv = POINTER(c_double)()
406 ds_inv = POINTER(c_double)()
407 dt_inv = POINTER(c_double)()
408 wx = POINTER(c_double)()
409 wy = POINTER(c_double)()
410 wz = POINTER(c_double)()
411 dx = POINTER(c_double)()
412 dy = POINTER(c_double)()
413 dz = POINTER(c_double)()
414 libneko.neko_field_space(name, byref(lx), byref(zg),
415 byref(dr_inv), byref(ds_inv), byref(dt_inv),
416 byref(wx), byref(wy), byref(wz),
417 byref(dx), byref(dy), byref(dz))
418
419 space = wrap_space(lx.value, zg,
420 dr_inv, ds_inv, dt_inv,
421 wx, wy, wz, dx, dy, dz)
422
423 return space
424
425def case_fluid_space(case_descr):
426 """ Retrive a space of a Neko case's fluid solver as a `neko_space_t`. """
427 lx = c_int()
428 zg = POINTER(c_double)()
429 dr_inv = POINTER(c_double)()
430 ds_inv = POINTER(c_double)()
431 dt_inv = POINTER(c_double)()
432 wx = POINTER(c_double)()
433 wy = POINTER(c_double)()
434 wz = POINTER(c_double)()
435 dx = POINTER(c_double)()
436 dy = POINTER(c_double)()
437 dz = POINTER(c_double)()
438 libneko.neko_case_fluid_space(byref(case_descr), byref(lx), byref(zg),
439 byref(dr_inv), byref(ds_inv), byref(dt_inv),
440 byref(wx), byref(wy), byref(wz),
441 byref(dx), byref(dy), byref(dz))
442
443 space = wrap_space(lx.value, zg,
444 dr_inv, ds_inv, dt_inv,
445 wx, wy, wz, dx, dy, dz)
446
447 return space
448
449def wrap_space(lx, zg, dr_inv, ds_inv, dt_inv, wx, wy, wz, dx, dy, dz):
450 space = neko_space_t(
451 lx,
452 np.frombuffer((c_double * lx * 3).from_address(
453 addressof(zg.contents))),
454 np.frombuffer((c_double * lx).from_address(
455 addressof(dr_inv.contents))),
456 np.frombuffer((c_double * lx).from_address(
457 addressof(ds_inv.contents))),
458 np.frombuffer((c_double * lx).from_address(
459 addressof(dt_inv.contents))),
460 np.frombuffer((c_double * lx).from_address(
461 addressof(wx.contents))),
462 np.frombuffer((c_double * lx).from_address(
463 addressof(wy.contents))),
464 np.frombuffer((c_double * lx).from_address(
465 addressof(wz.contents))),
466 np.frombuffer((c_double * lx * lx).from_address(
467 addressof(dx.contents))),
468 np.frombuffer((c_double * lx * lx).from_address(
469 addressof(dy.contents))),
470 np.frombuffer((c_double * lx * lx).from_address(
471 addressof(dz.contents)))
472 )
473
474 return space
475
476def case_fluid_coef(case_descr):
477 """ Retrive coefficients of a Neko case's fluid solver as a `neko_coef_t`. """
478 G11_ptr = POINTER(c_double)()
479 G22_ptr = POINTER(c_double)()
480 G33_ptr = POINTER(c_double)()
481 G12_ptr = POINTER(c_double)()
482 G13_ptr = POINTER(c_double)()
483 G23_ptr = POINTER(c_double)()
484 mult_ptr = POINTER(c_double)()
485 dxdr_ptr = POINTER(c_double)()
486 dydr_ptr = POINTER(c_double)()
487 dzdr_ptr = POINTER(c_double)()
488 dxds_ptr = POINTER(c_double)()
489 dyds_ptr = POINTER(c_double)()
490 dzds_ptr = POINTER(c_double)()
491 dxdt_ptr = POINTER(c_double)()
492 dydt_ptr = POINTER(c_double)()
493 dzdt_ptr = POINTER(c_double)()
494 drdx_ptr = POINTER(c_double)()
495 drdy_ptr = POINTER(c_double)()
496 drdz_ptr = POINTER(c_double)()
497 dsdx_ptr = POINTER(c_double)()
498 dsdy_ptr = POINTER(c_double)()
499 dsdz_ptr = POINTER(c_double)()
500 dtdx_ptr = POINTER(c_double)()
501 dtdy_ptr = POINTER(c_double)()
502 dtdz_ptr = POINTER(c_double)()
503 jac_ptr = POINTER(c_double)()
504 B_ptr = POINTER(c_double)()
505 area_ptr = POINTER(c_double)()
506 nx_ptr = POINTER(c_double)()
507 ny_ptr = POINTER(c_double)()
508 nz_ptr = POINTER(c_double)()
509
510 libneko.neko_case_fluid_coef(byref(case_descr),
511 byref(G11_ptr),
512 byref(G22_ptr),
513 byref(G33_ptr),
514 byref(G12_ptr),
515 byref(G13_ptr),
516 byref(G23_ptr),
517 byref(mult_ptr),
518 byref(dxdr_ptr),
519 byref(dydr_ptr),
520 byref(dzdr_ptr),
521 byref(dxds_ptr),
522 byref(dyds_ptr),
523 byref(dzds_ptr),
524 byref(dxdt_ptr),
525 byref(dydt_ptr),
526 byref(dzdt_ptr),
527 byref(drdx_ptr),
528 byref(drdy_ptr),
529 byref(drdz_ptr),
530 byref(dsdx_ptr),
531 byref(dsdy_ptr),
532 byref(dsdz_ptr),
533 byref(dtdx_ptr),
534 byref(dtdy_ptr),
535 byref(dtdz_ptr),
536 byref(jac_ptr),
537 byref(B_ptr),
538 byref(area_ptr),
539 byref(nx_ptr),
540 byref(ny_ptr),
541 byref(nz_ptr))
542
543 dm = case_fluid_dofmap(case_descr)
544
545 coef = neko_coef_t(
547 addressof(G11_ptr.contents))),
549 addressof(G22_ptr.contents))),
551 addressof(G33_ptr.contents))),
553 addressof(G12_ptr.contents))),
555 addressof(G13_ptr.contents))),
557 addressof(G23_ptr.contents))),
559 addressof(mult_ptr.contents))),
561 addressof(dxdr_ptr.contents))),
563 addressof(dydr_ptr.contents))),
565 addressof(dzdr_ptr.contents))),
567 addressof(dxds_ptr.contents))),
569 addressof(dyds_ptr.contents))),
571 addressof(dzds_ptr.contents))),
573 addressof(dxdt_ptr.contents))),
575 addressof(dydt_ptr.contents))),
577 addressof(dzdt_ptr.contents))),
579 addressof(drdx_ptr.contents))),
581 addressof(drdy_ptr.contents))),
583 addressof(drdz_ptr.contents))),
585 addressof(dsdx_ptr.contents))),
587 addressof(dsdy_ptr.contents))),
589 addressof(dsdz_ptr.contents))),
591 addressof(dtdx_ptr.contents))),
593 addressof(dtdy_ptr.contents))),
595 addressof(dtdz_ptr.contents))),
597 addressof(jac_ptr.contents))),
599 addressof(B_ptr.contents))),
601 addressof(area_ptr.contents))),
603 addressof(nx_ptr.contents))),
605 addressof(ny_ptr.contents))),
607 addressof(nz_ptr.contents)))
608 )
609
610 return coef
611
612def user_setup(case_descr, cb_initial_condition, cb_preprocess, cb_compute,
613 cb_dirichlet_condition, cb_material_properties, cb_source_term):
614 """ Setup user-provided callbacks. """
615 libneko.neko_user_setup(byref(case_descr),
616 cb_initial_condition,
617 cb_preprocess,
618 cb_compute,
619 cb_dirichlet_condition,
620 cb_material_properties,
621 cb_source_term)
622
623def output(case_descr, force = False):
624 libneko.neko_output_ctrl_execute(byref(case_descr), force)
625
626def callback_field(name):
627 """ Retrive a callback's Neko field as a `neko_field_t`. """
628 a = neko_field_t(
629 np.frombuffer((c_double *
631 addressof(libneko.neko_cb_field_by_name(name).contents)),
632 dtype=np.dtype(c_double)),
633 field_dofmap(name),
634 field_space(name)
635 )
636 return a
637
638def callback_field_name(index, name):
639 """ Check the name of a field at a certain index in the callback's field list. """
640 field_index = c_int(index)
641 same_field = libneko.neko_cb_field_name_at_index(byref(field_index), name)
642 return same_field > 0
643
644def bc_mask(msk, msk_size):
645 mask = np.frombuffer((c_int * msk_size).from_address(
646 addressof(msk.contents)), dtype=np.dtype(np.int32))
647 return mask
648
649#
650# https://degenerateconic.com/fortran-json-python.html
651#
652
654 return cast(create_string_buffer(s.encode()),c_char_p)
655
__global__ void dirichlet_apply_scalar_kernel(const int *__restrict__ msk, T *__restrict__ x, const T g, const int m)
Defines a field.
Definition field.f90:34
Defines an output.
Definition output.f90:34
user_setup(case_descr, cb_initial_condition, cb_preprocess, cb_compute, cb_dirichlet_condition, cb_material_properties, cb_source_term)
Definition intf.py:613
python_str_to_fortran(s)
Definition intf.py:653
field_registry_free()
Definition intf.py:268
wrap_dofmap(dof_ptr, x_ptr, y_ptr, z_ptr, size)
Definition intf.py:387
job_info()
Definition intf.py:260
wrap_space(lx, zg, dr_inv, ds_inv, dt_inv, wx, wy, wz, dx, dy, dz)
Definition intf.py:449
python_dict_to_fortran(d)
Definition intf.py:656
device_init()
Definition intf.py:252
material_properties
Definition intf.py:29
initial_condition
Definition intf.py:21
field_registry_init()
Definition intf.py:264
device_finalize()
Definition intf.py:256
dirichlet_condition
Definition intf.py:27
finalize()
Definition intf.py:248
Implements the source_term_t type and a wrapper source_term_wrapper_t.