Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
neko_api.f90
Go to the documentation of this file.
1! Copyright (c) 2022-2026, 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 neko
38 use, intrinsic :: iso_c_binding
39 implicit none
40 private
41
42contains
43
45 subroutine neko_api_init() bind(c, name="neko_init")
46
47 call neko_init()
48
49 end subroutine neko_api_init
50
52 subroutine neko_api_finalize() bind(c, name="neko_finalize")
53
54 call neko_finalize()
55
56 end subroutine neko_api_finalize
57
59 subroutine neko_api_device_init() bind(c, name="neko_device_init")
60
61 call device_init
62
63 end subroutine neko_api_device_init
64
66 subroutine neko_api_device_finalize() bind(c, name="neko_device_finalize")
67
68 call device_finalize
69
70 end subroutine neko_api_device_finalize
71
73 subroutine neko_api_job_info() bind(c, name="neko_job_info")
74 logical :: initialized
75
76 call mpi_initialized(initialized)
77
78 if (.not.initialized) then
79 call neko_warning('Neko has not been initialised')
80 else
81 call neko_job_info()
82 call neko_log%newline()
83 end if
84
85 end subroutine neko_api_job_info
86
88 subroutine neko_api_registry_init() bind(c, name="neko_registry_init")
89
90 call neko_registry%init()
91 call neko_const_registry%init()
92
93 end subroutine neko_api_registry_init
94
96 subroutine neko_api_registry_free() bind(c, name="neko_registry_free")
97
98 call neko_registry%free()
99
100 end subroutine neko_api_registry_free
101
103 subroutine neko_api_scratch_registry_init() bind(c, name="neko_scratch_registry_init")
104
105 call neko_scratch_registry%init()
106
107 end subroutine neko_api_scratch_registry_init
108
110 subroutine neko_api_scratch_registry_free() bind(c, name="neko_scratch_registry_free")
111
112 call neko_scratch_registry%free()
113
114 end subroutine neko_api_scratch_registry_free
115
118 subroutine neko_api_case_allocate(case_iptr) &
119 bind(c, name="neko_case_allocate")
120 integer(c_intptr_t), intent(inout) :: case_iptr
121 type(case_t), pointer :: C
122 type(c_ptr) :: cp
123
124 allocate(c)
125 cp = c_loc(c)
126 case_iptr = transfer(cp, 0_c_intptr_t)
127
128 end subroutine neko_api_case_allocate
129
133 subroutine neko_api_case_init(case_json, case_len, case_iptr) &
134 bind(c, name="neko_case_init")
135 type(c_ptr), intent(in) :: case_json
136 integer(c_int), value :: case_len
137 integer(c_intptr_t), intent(inout) :: case_iptr
138 type(json_file) :: json_case
139 type(case_t), pointer :: C
140 type(c_ptr) :: cp
141
142 ! Check if the case has already been allocated
143 ! e.g. if a user callback has been injected
144 cp = transfer(case_iptr, c_null_ptr)
145 if (c_associated(cp)) then
146 call c_f_pointer(cp, c)
147 else
148 allocate(c)
149 cp = c_loc(c)
150 case_iptr = transfer(cp, 0_c_intptr_t)
151 end if
152
153 ! Convert passed in serialised JSON object into a Fortran
154 ! character string and create a json_file object
155 if (c_associated(case_json)) then
156 block
157 character(kind=c_char,len=case_len+1),pointer :: s
158 character(len=:), allocatable :: fcase_json
159 call c_f_pointer(case_json, s)
160 fcase_json = s(1:case_len)
161 call json_case%load_from_string(fcase_json)
162 deallocate(fcase_json)
163 nullify(s)
164 end block
165 end if
166
167 !
168 ! Create case
169 !
170 call c%init(json_case)
171
172 !
173 ! Create simulation components
174 !
175 call neko_simcomps%init(c)
176
177 end subroutine neko_api_case_init
178
179
182 subroutine neko_api_case_free(case_iptr) bind(c, name="neko_case_free")
183 integer(c_intptr_t), intent(inout) :: case_iptr
184 type(case_t), pointer :: C
185 type(c_ptr) :: cp
186
187 cp = transfer(case_iptr, c_null_ptr)
188 if (c_associated(cp)) then
189 call c_f_pointer(cp, c)
190 call c%free()
191 else
192 call neko_error('Invalid Neko case')
193 end if
194
195 end subroutine neko_api_case_free
196
200 function neko_api_case_time(case_iptr) result(time) &
201 bind(c, name="neko_case_time")
202 integer(c_intptr_t), intent(inout) :: case_iptr
203 real(kind=c_rp) :: time
204 type(case_t), pointer :: c
205 type(c_ptr) :: cptr
206
207 cptr = transfer(case_iptr, c_null_ptr)
208 if (c_associated(cptr)) then
209 call c_f_pointer(cptr, c)
210 time = c%time%t
211 else
212 call neko_error('Invalid Neko case')
213 end if
214
215 end function neko_api_case_time
216
220 function neko_api_case_end_time(case_iptr) result(end_time) &
221 bind(c, name="neko_case_end_time")
222 integer(c_intptr_t), intent(inout) :: case_iptr
223 real(kind=c_rp) :: end_time
224 type(case_t), pointer :: c
225 type(c_ptr) :: cptr
226
227 cptr = transfer(case_iptr, c_null_ptr)
228 if (c_associated(cptr)) then
229 call c_f_pointer(cptr, c)
230 end_time = c%time%end_time
231 else
232 call neko_error('Invalid Neko case')
233 end if
234
235 end function neko_api_case_end_time
236
240 function neko_api_case_tstep(case_iptr) result(tstep) &
241 bind(c, name="neko_case_tstep")
242 integer(c_intptr_t), intent(inout) :: case_iptr
243 type(case_t), pointer :: c
244 type(c_ptr) :: cptr
245 integer(c_int) :: tstep
246
247 cptr = transfer(case_iptr, c_null_ptr)
248 if (c_associated(cptr)) then
249 call c_f_pointer(cptr, c)
250 tstep = c%time%tstep
251 else
252 call neko_error('Invalid Neko case')
253 end if
254
255 end function neko_api_case_tstep
256
259 subroutine neko_api_solve(case_iptr) bind(c, name="neko_solve")
260 integer(c_intptr_t), intent(inout) :: case_iptr
261 type(case_t), pointer :: C
262 type(c_ptr) :: cp
263
264 cp = transfer(case_iptr, c_null_ptr)
265 if (c_associated(cp)) then
266 call c_f_pointer(cp, c)
267 call neko_solve(c)
268 else
269 call neko_error('Invalid Neko case')
270 end if
271
272 end subroutine neko_api_solve
273
276 subroutine neko_api_step(case_iptr) &
277 bind(c, name="neko_step")
279 integer(c_intptr_t), intent(inout) :: case_iptr
280 type(case_t), pointer :: C
281 type(c_ptr) :: cptr
282 type(json_file) :: dt_params
283 type(time_step_controller_t), save, allocatable :: dt_controller
284 real(kind=dp), save :: step_loop_start = 0.0_dp
285
286 cptr = transfer(case_iptr, c_null_ptr)
287 if (c_associated(cptr)) then
288 call c_f_pointer(cptr, c)
289
290 if (.not. allocated(dt_controller)) then
291 allocate(dt_controller)
292 call json_get(c%params, 'case.time', dt_params)
293 call dt_controller%init(dt_params)
294 end if
295
296 if (c%time%tstep .eq. 0) then
297 call simulation_init(c, dt_controller)
298
299 step_loop_start = mpi_wtime()
300 end if
301
302 if (.not. c%time%is_done()) then
303 call simulation_step(c, dt_controller, step_loop_start)
304 end if
305
306 if (c%time%is_done()) then
307 call simulation_finalize(c)
308 if (allocated(dt_controller)) then
309 deallocate(dt_controller)
310 end if
311 end if
312
313 else
314 call neko_error('Invalid Neko case')
315 end if
316
317 end subroutine neko_api_step
318
323 subroutine neko_api_output_ctrl_execute(case_iptr, force_output) &
324 bind(c, name="neko_output_ctrl_execute")
325 integer(c_intptr_t), intent(inout) :: case_iptr
326 logical(kind=c_bool), value :: force_output
327 logical :: f_force_output
328 type(case_t), pointer :: C
329 type(c_ptr) :: cp
330 type(time_state_t) :: f_time
331
332 cp = transfer(case_iptr, c_null_ptr)
333 if (c_associated(cp)) then
334 call c_f_pointer(cp, c)
335 else
336 call neko_error('Invalid Neko case')
337 end if
338
339 f_force_output = transfer(force_output, f_force_output)
340 call c%output_controller%execute(c%time, f_force_output)
341
342 end subroutine neko_api_output_ctrl_execute
343
346 function neko_api_field(field_name) result(field_ptr) &
347 bind(c, name='neko_field')
348 character(kind=c_char), dimension(*), intent(in) :: field_name
349 character(len=8192) :: name
350 type(field_t), pointer :: field
351 type(c_ptr) :: field_ptr
352 integer :: len
353
354 len = 0
355 do
356 if (field_name(len+1) .eq. c_null_char) exit
357 len = len + 1
358 name(len:len) = field_name(len)
359 end do
360
361 field => neko_registry%get_field(trim(name(1:len)))
362
363 field_ptr = c_loc(field%x)
364
365 end function neko_api_field
366
369 function neko_api_field_order(field_name) result(field_lx) &
370 bind(c, name='neko_field_order')
371 character(kind=c_char), dimension(*), intent(in) :: field_name
372 character(len=8192) :: name
373 type(field_t), pointer :: field
374 integer(c_int) :: field_lx
375 integer :: len
376
377 len = 0
378 do
379 if (field_name(len+1) .eq. c_null_char) exit
380 len = len + 1
381 name(len:len) = field_name(len)
382 end do
383
384 field => neko_registry%get_field(trim(name(1:len)))
385
386 field_lx = field%Xh%lx
387
388 end function neko_api_field_order
389
392 function neko_api_field_nelements(field_name) result(field_nelv) &
393 bind(c, name='neko_field_nelements')
394 character(kind=c_char), dimension(*), intent(in) :: field_name
395 character(len=8192) :: name
396 type(field_t), pointer :: field
397 integer(c_int) :: field_nelv
398 integer :: len
399
400 len = 0
401 do
402 if (field_name(len+1) .eq. c_null_char) exit
403 len = len + 1
404 name(len:len) = field_name(len)
405 end do
406
407 field => neko_registry%get_field(trim(name(1:len)))
408
409 field_nelv = field%msh%nelv
410
411 end function neko_api_field_nelements
412
415 function neko_api_field_size(field_name) result(field_size) &
416 bind(c, name='neko_field_size')
417 character(kind=c_char), dimension(*), intent(in) :: field_name
418 character(len=8192) :: name
419 type(field_t), pointer :: field
420 integer(c_int) :: field_size
421 integer :: len
422
423 len = 0
424 do
425 if (field_name(len+1) .eq. c_null_char) exit
426 len = len + 1
427 name(len:len) = field_name(len)
428 end do
429
430 field => neko_registry%get_field(trim(name(1:len)))
431
432 field_size = field%dof%size()
433
434 end function neko_api_field_size
435
442 subroutine neko_api_field_dofmap(field_name, dof_ptr, x_ptr, y_ptr, z_ptr) &
443 bind(c, name='neko_field_dofmap')
444 character(kind=c_char), dimension(*), intent(in) :: field_name
445 type(c_ptr), intent(inout) :: dof_ptr, x_ptr, y_ptr, z_ptr
446 character(len=8192) :: name
447 type(field_t), pointer :: field
448 integer :: len
449
450 len = 0
451 do
452 if (field_name(len+1) .eq. c_null_char) exit
453 len = len + 1
454 name(len:len) = field_name(len)
455 end do
456
457 field => neko_registry%get_field(trim(name(1:len)))
458 call neko_api_wrap_dofmap(field%dof, dof_ptr, x_ptr, y_ptr, z_ptr)
459
460 end subroutine neko_api_field_dofmap
461
469 subroutine neko_api_case_fluid_dofmap(case_iptr, dof_ptr, &
470 x_ptr, y_ptr, z_ptr, size) bind(c, name='neko_case_fluid_dofmap')
471 integer(c_intptr_t), intent(inout) :: case_iptr
472 type(case_t), pointer :: C
473 type(c_ptr) :: cptr
474 type(c_ptr), intent(inout) :: dof_ptr, x_ptr, y_ptr, z_ptr
475 integer, intent(inout) :: size
476
477 cptr = transfer(case_iptr, c_null_ptr)
478 if (c_associated(cptr)) then
479 call c_f_pointer(cptr, c)
480 call neko_api_wrap_dofmap(c%fluid%dm_Xh, dof_ptr, x_ptr, y_ptr, z_ptr)
481 size = c%fluid%dm_Xh%size()
482 else
483 call neko_error('Invalid Neko case')
484 end if
485
486 end subroutine neko_api_case_fluid_dofmap
487
489 subroutine neko_api_wrap_dofmap(dm, dof_ptr, x_ptr, y_ptr, z_ptr)
490 type(dofmap_t), target, intent(in) :: dm
491 type(c_ptr), intent(inout) :: dof_ptr, x_ptr, y_ptr, z_ptr
492
493 dof_ptr = c_loc(dm%dof)
494 x_ptr = c_loc(dm%x)
495 y_ptr = c_loc(dm%y)
496 z_ptr = c_loc(dm%z)
497
498 end subroutine neko_api_wrap_dofmap
499
513 subroutine neko_api_field_space(field_name, lx, zg, &
514 dr_inv, ds_inv, dt_inv, wx, wy, wz, dx, dy, dz) &
515 bind(c, name='neko_field_space')
516 character(kind=c_char), dimension(*), intent(in) :: field_name
517 integer, intent(inout) :: lx
518 type(c_ptr), intent(inout) :: zg, dr_inv, ds_inv, dt_inv
519 type(c_ptr), intent(inout) :: wx, wy, wz, dx, dy, dz
520 character(len=8192) :: name
521 type(field_t), pointer :: field
522 integer :: len
523
524 len = 0
525 do
526 if (field_name(len+1) .eq. c_null_char) exit
527 len = len + 1
528 name(len:len) = field_name(len)
529 end do
530
531 field => neko_registry%get_field(trim(name(1:len)))
532 call neko_api_wrap_space(field%Xh, lx, zg, dr_inv, ds_inv, dt_inv, &
533 wx, wy, wz, dx, dy, dz)
534
535 end subroutine neko_api_field_space
536
550 subroutine neko_api_case_fluid_space(case_iptr, lx, zg, &
551 dr_inv, ds_inv, dt_inv, wx, wy, wz, dx, dy, dz) &
552 bind(c, name='neko_case_fluid_space')
553 integer(c_intptr_t), intent(inout) :: case_iptr
554 type(case_t), pointer :: C
555 type(c_ptr) :: cptr
556 integer, intent(inout) :: lx
557 type(c_ptr), intent(inout) :: zg, dr_inv, ds_inv, dt_inv
558 type(c_ptr), intent(inout) :: wx, wy, wz, dx, dy, dz
559
560
561 cptr = transfer(case_iptr, c_null_ptr)
562 if (c_associated(cptr)) then
563 call c_f_pointer(cptr, c)
564 call neko_api_wrap_space(c%fluid%Xh, lx, zg, dr_inv, ds_inv, dt_inv, &
565 wx, wy, wz, dx, dy, dz)
566 else
567 call neko_error('Invalid Neko case')
568 end if
569
570 end subroutine neko_api_case_fluid_space
571
573 subroutine neko_api_wrap_space(Xh, lx, zg, dr_inv, ds_inv, dt_inv, &
574 wx, wy, wz, dx, dy, dz)
575 type(space_t), target, intent(in) :: Xh
576 integer, intent(inout) :: lx
577 type(c_ptr), intent(inout) :: zg, dr_inv, ds_inv, dt_inv
578 type(c_ptr), intent(inout) :: wx, wy, wz, dx, dy, dz
579
580 lx = xh%lx
581 zg = c_loc(xh%zg)
582 dr_inv = c_loc(xh%dr_inv)
583 ds_inv = c_loc(xh%ds_inv)
584 dt_inv = c_loc(xh%dt_inv)
585 wx = c_loc(xh%wx)
586 wy = c_loc(xh%wy)
587 wz = c_loc(xh%wz)
588 dx = c_loc(xh%dx)
589 dy = c_loc(xh%dy)
590 dz = c_loc(xh%dz)
591
592 end subroutine neko_api_wrap_space
593
596 subroutine neko_api_case_fluid_coef(case_iptr, G11, G22, G33, G12, G13, G23, &
597 mult, dxdr, dydr, dzdr, dxds, dyds, dzds, dxdt, dydt, dzdt, &
598 drdx, drdy, drdz, dsdx, dsdy, dsdz, dtdx, dtdy, dtdz, &
599 jac, B, area, nx, ny, nz) bind(c, name='neko_case_fluid_coef')
600 integer(c_intptr_t), intent(inout) :: case_iptr
601 type(case_t), pointer :: C
602 type(c_ptr) :: cptr
603 type(c_ptr), intent(inout) :: G11, G22, G33, G12, G13, G23
604 type(c_ptr), intent(inout) :: mult
605 type(c_ptr), intent(inout) :: dxdr, dydr, dzdr
606 type(c_ptr), intent(inout) :: dxds, dyds, dzds
607 type(c_ptr), intent(inout) :: dxdt, dydt, dzdt
608 type(c_ptr), intent(inout) :: drdx, drdy, drdz
609 type(c_ptr), intent(inout) :: dsdx, dsdy, dsdz
610 type(c_ptr), intent(inout) :: dtdx, dtdy, dtdz
611 type(c_ptr), intent(inout) :: jac, B, area, nx, ny, nz
612
613
614 cptr = transfer(case_iptr, c_null_ptr)
615 if (c_associated(cptr)) then
616 call c_f_pointer(cptr, c)
617 g11 = c_loc(c%fluid%c_Xh%G11)
618 g22 = c_loc(c%fluid%c_Xh%G22)
619 g33 = c_loc(c%fluid%c_Xh%G33)
620 g12 = c_loc(c%fluid%c_Xh%G12)
621 g13 = c_loc(c%fluid%c_Xh%G13)
622 g23 = c_loc(c%fluid%c_Xh%G23)
623 mult = c_loc(c%fluid%c_Xh%mult)
624 dxdr = c_loc(c%fluid%c_Xh%dxdr)
625 dydr = c_loc(c%fluid%c_Xh%dydr)
626 dzdr = c_loc(c%fluid%c_Xh%dzdr)
627 dxds = c_loc(c%fluid%c_Xh%dxds)
628 dyds = c_loc(c%fluid%c_Xh%dyds)
629 dzds = c_loc(c%fluid%c_Xh%dzds)
630 dxdt = c_loc(c%fluid%c_Xh%dxdt)
631 dydt = c_loc(c%fluid%c_Xh%dydt)
632 dzdt = c_loc(c%fluid%c_Xh%dzdt)
633 drdx = c_loc(c%fluid%c_Xh%drdx)
634 drdy = c_loc(c%fluid%c_Xh%drdy)
635 drdz = c_loc(c%fluid%c_Xh%drdz)
636 dsdx = c_loc(c%fluid%c_Xh%dsdx)
637 dsdy = c_loc(c%fluid%c_Xh%dsdy)
638 dsdz = c_loc(c%fluid%c_Xh%dsdz)
639 dtdx = c_loc(c%fluid%c_Xh%dtdx)
640 dtdy = c_loc(c%fluid%c_Xh%dtdy)
641 dtdz = c_loc(c%fluid%c_Xh%dtdz)
642 jac = c_loc(c%fluid%c_Xh%jac)
643 b = c_loc(c%fluid%c_Xh%B)
644 area = c_loc(c%fluid%c_Xh%area)
645 nx = c_loc(c%fluid%c_Xh%nx)
646 ny = c_loc(c%fluid%c_Xh%ny)
647 nz = c_loc(c%fluid%c_Xh%nz)
648 else
649 call neko_error('Invalid Neko case')
650 end if
651
652 end subroutine neko_api_case_fluid_coef
653
662 subroutine neko_api_user_setup(case_iptr, initial_cb, preprocess_cb, &
663 compute_cb, dirichlet_cb, material_cb, source_cb) &
664 bind(c, name='neko_user_setup')
665 integer(c_intptr_t), intent(inout) :: case_iptr
666 type(c_funptr), value :: initial_cb, preprocess_cb, compute_cb
667 type(c_funptr), value :: dirichlet_cb, material_cb, source_cb
668 type(case_t), pointer :: C
669 type(c_ptr) :: cptr
670
671 cptr = transfer(case_iptr, c_null_ptr)
672 if (c_associated(cptr)) then
673 call c_f_pointer(cptr, c)
674 call neko_api_user_cb_register(c%user, initial_cb, preprocess_cb, &
675 compute_cb, dirichlet_cb, material_cb, source_cb)
676 else
677 call neko_error('Invalid Neko case')
678 end if
679
680 end subroutine neko_api_user_setup
681
684 function neko_api_user_cb_field_by_name(field_name) result(field_ptr) &
685 bind(c, name='neko_cb_field_by_name')
686 character(kind=c_char), dimension(*), intent(in) :: field_name
687 character(len=8192) :: name
688 type(field_t), pointer :: field
689 type(c_ptr) :: field_ptr
690 integer :: len
691
692 len = 0
693 do
694 if (field_name(len+1) .eq. c_null_char) exit
695 len = len + 1
696 name(len:len) = field_name(len)
697 end do
698
699 field => neko_api_user_cb_get_field(trim(name(1:len)))
700
701 field_ptr = c_loc(field%x)
702
704
708 result(field_ptr) bind(c, name='neko_cb_field_by_index')
709 integer, intent(in) :: field_idx
710 type(field_t), pointer :: field
711 type(c_ptr) :: field_ptr
712
713 field => neko_api_user_cb_get_field(field_idx)
714
715 field_ptr = c_loc(field%x)
716
718
722 function neko_api_user_cb_field_name_at_index(field_idx, field_name) &
723 result(same_name) bind(c, name='neko_cb_field_name_at_index')
724 integer, intent(in) :: field_idx
725 character(kind=c_char), dimension(*), intent(in) :: field_name
726 character(len=8192) :: name
727 type(field_t), pointer :: f1, f2
728 type(c_ptr) :: field_ptr
729 integer :: len
730 logical(c_bool) :: same_name
731
732 len = 0
733 do
734 if (field_name(len+1) .eq. c_null_char) exit
735 len = len + 1
736 name(len:len) = field_name(len)
737 end do
738
739 f1 => neko_api_user_cb_get_field(field_idx)
740 f2 => neko_api_user_cb_get_field(trim(name(1:len)))
741
742 same_name = trim(f1%name) .eq. trim(f2%name)
743
745
746
747end module neko_api
Defines a field.
Definition field.f90:34
Neko API user callbacks.
subroutine, public neko_api_user_cb_register(user, initial_cb, preprocess_cb, compute_cb, dirichlet_cb, material_cb, source_cb)
Register callbacks.
Neko C API.
Definition neko_api.f90:34
subroutine neko_api_wrap_dofmap(dm, dof_ptr, x_ptr, y_ptr, z_ptr)
Helper function to assign pointers to a dofmap's data.
Definition neko_api.f90:490
subroutine neko_api_finalize()
Finalize Neko.
Definition neko_api.f90:53
subroutine neko_api_wrap_space(xh, lx, zg, dr_inv, ds_inv, dt_inv, wx, wy, wz, dx, dy, dz)
Helper function to assign pointers to a space's data.
Definition neko_api.f90:575
subroutine neko_api_user_setup(case_iptr, initial_cb, preprocess_cb, compute_cb, dirichlet_cb, material_cb, source_cb)
Setup user-provided callbacks.
Definition neko_api.f90:665
subroutine neko_api_step(case_iptr)
Compute a time-step for a neko case.
Definition neko_api.f90:278
integer(c_int) function neko_api_field_nelements(field_name)
Retrive the number of elements in a field.
Definition neko_api.f90:394
subroutine neko_api_device_finalize()
Finalize Neko device layer.
Definition neko_api.f90:67
subroutine neko_api_job_info()
Display job information.
Definition neko_api.f90:74
subroutine neko_api_registry_init()
Initialise a Neko field registry.
Definition neko_api.f90:89
subroutine neko_api_field_space(field_name, lx, zg, dr_inv, ds_inv, dt_inv, wx, wy, wz, dx, dy, dz)
Retrive the space associated with a field.
Definition neko_api.f90:516
subroutine neko_api_case_free(case_iptr)
Destroy a Neko case.
Definition neko_api.f90:183
subroutine neko_api_case_allocate(case_iptr)
Allocate memory for a Neko case.
Definition neko_api.f90:120
subroutine neko_api_output_ctrl_execute(case_iptr, force_output)
Execute the Case's output controller.
Definition neko_api.f90:325
subroutine neko_api_case_init(case_json, case_len, case_iptr)
Initalise a Neko case.
Definition neko_api.f90:135
subroutine neko_api_field_dofmap(field_name, dof_ptr, x_ptr, y_ptr, z_ptr)
Retrive the dofmap associated with a field.
Definition neko_api.f90:444
subroutine neko_api_case_fluid_space(case_iptr, lx, zg, dr_inv, ds_inv, dt_inv, wx, wy, wz, dx, dy, dz)
Retrive the space associated with a case's fluid solver.
Definition neko_api.f90:553
real(kind=c_rp) function neko_api_case_end_time(case_iptr)
Retrive the end time of a case.
Definition neko_api.f90:222
type(c_ptr) function neko_api_user_cb_field_by_index(field_idx)
Retrive a pointer to a user callback field.
Definition neko_api.f90:709
logical(c_bool) function neko_api_user_cb_field_name_at_index(field_idx, field_name)
Check if the user callback field at a given index has a given name.
Definition neko_api.f90:724
subroutine neko_api_scratch_registry_init()
Initialise a Neko scratch registry.
Definition neko_api.f90:104
type(c_ptr) function neko_api_field(field_name)
Retrive a pointer to a flow field.
Definition neko_api.f90:348
integer(c_int) function neko_api_field_size(field_name)
Retrive the total number of degrees of freedom of a field.
Definition neko_api.f90:417
subroutine neko_api_scratch_registry_free()
Destroy a Neko scratch registry.
Definition neko_api.f90:111
subroutine neko_api_solve(case_iptr)
Solve a neko case.
Definition neko_api.f90:260
subroutine neko_api_case_fluid_coef(case_iptr, g11, g22, g33, g12, g13, g23, mult, dxdr, dydr, dzdr, dxds, dyds, dzds, dxdt, dydt, dzdt, drdx, drdy, drdz, dsdx, dsdy, dsdz, dtdx, dtdy, dtdz, jac, b, area, nx, ny, nz)
Retrive the coefficient associated with a case's fluid solver.
Definition neko_api.f90:600
subroutine neko_api_registry_free()
Destroy a Neko field registry.
Definition neko_api.f90:97
type(c_ptr) function neko_api_user_cb_field_by_name(field_name)
Retrive a pointer to a user callback field.
Definition neko_api.f90:686
subroutine neko_api_case_fluid_dofmap(case_iptr, dof_ptr, x_ptr, y_ptr, z_ptr, size)
Retrive the dofmap associated with a case's fluid solver.
Definition neko_api.f90:471
integer(c_int) function neko_api_field_order(field_name)
Retrive the order of a field.
Definition neko_api.f90:371
subroutine neko_api_device_init()
Initialise Neko device layer.
Definition neko_api.f90:60
real(kind=c_rp) function neko_api_case_time(case_iptr)
Retrive the current time of a case.
Definition neko_api.f90:202
integer(c_int) function neko_api_case_tstep(case_iptr)
Retrive the time-step of a case.
Definition neko_api.f90:242
subroutine neko_api_init()
Initialise Neko.
Definition neko_api.f90:46
Master module.
Definition neko.f90:34
Implements type time_step_controller.
void neko_finalize()
void neko_init()
void neko_job_info()
void neko_solve(int **case_iptr)