Neko 1.0.1
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
92 end subroutine neko_api_registry_init
93
95 subroutine neko_api_registry_free() bind(c, name="neko_registry_free")
96
97 call neko_registry%free()
98
99 end subroutine neko_api_registry_free
100
103 subroutine neko_api_case_allocate(case_iptr) &
104 bind(c, name="neko_case_allocate")
105 integer(c_intptr_t), intent(inout) :: case_iptr
106 type(case_t), pointer :: C
107 type(c_ptr) :: cp
108
109 allocate(c)
110 cp = c_loc(c)
111 case_iptr = transfer(cp, 0_c_intptr_t)
112
113 end subroutine neko_api_case_allocate
114
118 subroutine neko_api_case_init(case_json, case_len, case_iptr) &
119 bind(c, name="neko_case_init")
120 type(c_ptr), intent(in) :: case_json
121 integer(c_int), value :: case_len
122 integer(c_intptr_t), intent(inout) :: case_iptr
123 type(json_file) :: json_case
124 type(case_t), pointer :: C
125 type(c_ptr) :: cp
126
127 ! Check if the case has already been allocated
128 ! e.g. if a user callback has been injected
129 cp = transfer(case_iptr, c_null_ptr)
130 if (c_associated(cp)) then
131 call c_f_pointer(cp, c)
132 else
133 allocate(c)
134 cp = c_loc(c)
135 case_iptr = transfer(cp, 0_c_intptr_t)
136 end if
137
138 ! Convert passed in serialised JSON object into a Fortran
139 ! character string and create a json_file object
140 if (c_associated(case_json)) then
141 block
142 character(kind=c_char,len=case_len+1),pointer :: s
143 character(len=:), allocatable :: fcase_json
144 call c_f_pointer(case_json, s)
145 fcase_json = s(1:case_len)
146 call json_case%load_from_string(fcase_json)
147 deallocate(fcase_json)
148 nullify(s)
149 end block
150 end if
151
152 !
153 ! Create case
154 !
155 call c%init(json_case)
156
157 !
158 ! Create simulation components
159 !
160 call neko_simcomps%init(c)
161
162 end subroutine neko_api_case_init
163
164
167 subroutine neko_api_case_free(case_iptr) bind(c, name="neko_case_free")
168 integer(c_intptr_t), intent(inout) :: case_iptr
169 type(case_t), pointer :: C
170 type(c_ptr) :: cp
171
172 cp = transfer(case_iptr, c_null_ptr)
173 if (c_associated(cp)) then
174 call c_f_pointer(cp, c)
175 call c%free()
176 else
177 call neko_error('Invalid Neko case')
178 end if
179
180 end subroutine neko_api_case_free
181
185 function neko_api_case_time(case_iptr) result(time) &
186 bind(c, name="neko_case_time")
187 integer(c_intptr_t), intent(inout) :: case_iptr
188 real(kind=c_rp) :: time
189 type(case_t), pointer :: c
190 type(c_ptr) :: cptr
191
192 cptr = transfer(case_iptr, c_null_ptr)
193 if (c_associated(cptr)) then
194 call c_f_pointer(cptr, c)
195 time = c%time%t
196 else
197 call neko_error('Invalid Neko case')
198 end if
199
200 end function neko_api_case_time
201
205 function neko_api_case_end_time(case_iptr) result(end_time) &
206 bind(c, name="neko_case_end_time")
207 integer(c_intptr_t), intent(inout) :: case_iptr
208 real(kind=c_rp) :: end_time
209 type(case_t), pointer :: c
210 type(c_ptr) :: cptr
211
212 cptr = transfer(case_iptr, c_null_ptr)
213 if (c_associated(cptr)) then
214 call c_f_pointer(cptr, c)
215 end_time = c%time%end_time
216 else
217 call neko_error('Invalid Neko case')
218 end if
219
220 end function neko_api_case_end_time
221
225 function neko_api_case_tstep(case_iptr) result(tstep) &
226 bind(c, name="neko_case_tstep")
227 integer(c_intptr_t), intent(inout) :: case_iptr
228 type(case_t), pointer :: c
229 type(c_ptr) :: cptr
230 integer(c_int) :: tstep
231
232 cptr = transfer(case_iptr, c_null_ptr)
233 if (c_associated(cptr)) then
234 call c_f_pointer(cptr, c)
235 tstep = c%time%tstep
236 else
237 call neko_error('Invalid Neko case')
238 end if
239
240 end function neko_api_case_tstep
241
244 subroutine neko_api_solve(case_iptr) bind(c, name="neko_solve")
245 integer(c_intptr_t), intent(inout) :: case_iptr
246 type(case_t), pointer :: C
247 type(c_ptr) :: cp
248
249 cp = transfer(case_iptr, c_null_ptr)
250 if (c_associated(cp)) then
251 call c_f_pointer(cp, c)
252 call neko_solve(c)
253 else
254 call neko_error('Invalid Neko case')
255 end if
256
257 end subroutine neko_api_solve
258
261 subroutine neko_api_step(case_iptr) &
262 bind(c, name="neko_step")
264 integer(c_intptr_t), intent(inout) :: case_iptr
265 type(case_t), pointer :: C
266 type(c_ptr) :: cptr
267 type(json_file) :: dt_params
268 type(time_step_controller_t), save, allocatable :: dt_controller
269 real(kind=dp), save :: step_loop_start = 0.0_dp
270
271 cptr = transfer(case_iptr, c_null_ptr)
272 if (c_associated(cptr)) then
273 call c_f_pointer(cptr, c)
274
275 if (.not. allocated(dt_controller)) then
276 allocate(dt_controller)
277 call json_get(c%params, 'case.time', dt_params)
278 call dt_controller%init(dt_params)
279 end if
280
281 if (c%time%tstep .eq. 0) then
282 call simulation_init(c, dt_controller)
283
284 step_loop_start = mpi_wtime()
285 end if
286
287 if (.not. c%time%is_done()) then
288 call simulation_step(c, dt_controller, step_loop_start)
289 end if
290
291 if (c%time%is_done()) then
292 call simulation_finalize(c)
293 if (allocated(dt_controller)) then
294 deallocate(dt_controller)
295 end if
296 end if
297
298 else
299 call neko_error('Invalid Neko case')
300 end if
301
302 end subroutine neko_api_step
303
308 subroutine neko_api_output_ctrl_execute(case_iptr, force_output) &
309 bind(c, name="neko_output_ctrl_execute")
310 integer(c_intptr_t), intent(inout) :: case_iptr
311 logical(kind=c_bool), value :: force_output
312 logical :: f_force_output
313 type(case_t), pointer :: C
314 type(c_ptr) :: cp
315 type(time_state_t) :: f_time
316
317 cp = transfer(case_iptr, c_null_ptr)
318 if (c_associated(cp)) then
319 call c_f_pointer(cp, c)
320 else
321 call neko_error('Invalid Neko case')
322 end if
323
324 f_force_output = transfer(force_output, f_force_output)
325 call c%output_controller%execute(c%time, f_force_output)
326
327 end subroutine neko_api_output_ctrl_execute
328
331 function neko_api_field(field_name) result(field_ptr) &
332 bind(c, name='neko_field')
333 character(kind=c_char), dimension(*), intent(in) :: field_name
334 character(len=8192) :: name
335 type(field_t), pointer :: field
336 type(c_ptr) :: field_ptr
337 integer :: len
338
339 len = 0
340 do
341 if (field_name(len+1) .eq. c_null_char) exit
342 len = len + 1
343 name(len:len) = field_name(len)
344 end do
345
346 field => neko_registry%get_field(trim(name(1:len)))
347
348 field_ptr = c_loc(field%x)
349
350 end function neko_api_field
351
354 function neko_api_field_order(field_name) result(field_lx) &
355 bind(c, name='neko_field_order')
356 character(kind=c_char), dimension(*), intent(in) :: field_name
357 character(len=8192) :: name
358 type(field_t), pointer :: field
359 integer(c_int) :: field_lx
360 integer :: len
361
362 len = 0
363 do
364 if (field_name(len+1) .eq. c_null_char) exit
365 len = len + 1
366 name(len:len) = field_name(len)
367 end do
368
369 field => neko_registry%get_field(trim(name(1:len)))
370
371 field_lx = field%Xh%lx
372
373 end function neko_api_field_order
374
377 function neko_api_field_nelements(field_name) result(field_nelv) &
378 bind(c, name='neko_field_nelements')
379 character(kind=c_char), dimension(*), intent(in) :: field_name
380 character(len=8192) :: name
381 type(field_t), pointer :: field
382 integer(c_int) :: field_nelv
383 integer :: len
384
385 len = 0
386 do
387 if (field_name(len+1) .eq. c_null_char) exit
388 len = len + 1
389 name(len:len) = field_name(len)
390 end do
391
392 field => neko_registry%get_field(trim(name(1:len)))
393
394 field_nelv = field%msh%nelv
395
396 end function neko_api_field_nelements
397
400 function neko_api_field_size(field_name) result(field_size) &
401 bind(c, name='neko_field_size')
402 character(kind=c_char), dimension(*), intent(in) :: field_name
403 character(len=8192) :: name
404 type(field_t), pointer :: field
405 integer(c_int) :: field_size
406 integer :: len
407
408 len = 0
409 do
410 if (field_name(len+1) .eq. c_null_char) exit
411 len = len + 1
412 name(len:len) = field_name(len)
413 end do
414
415 field => neko_registry%get_field(trim(name(1:len)))
416
417 field_size = field%dof%size()
418
419 end function neko_api_field_size
420
427 subroutine neko_api_field_dofmap(field_name, dof_ptr, x_ptr, y_ptr, z_ptr) &
428 bind(c, name='neko_field_dofmap')
429 character(kind=c_char), dimension(*), intent(in) :: field_name
430 type(c_ptr), intent(inout) :: dof_ptr, x_ptr, y_ptr, z_ptr
431 character(len=8192) :: name
432 type(field_t), pointer :: field
433 integer :: len
434
435 len = 0
436 do
437 if (field_name(len+1) .eq. c_null_char) exit
438 len = len + 1
439 name(len:len) = field_name(len)
440 end do
441
442 field => neko_registry%get_field(trim(name(1:len)))
443 call neko_api_wrap_dofmap(field%dof, dof_ptr, x_ptr, y_ptr, z_ptr)
444
445 end subroutine neko_api_field_dofmap
446
454 subroutine neko_api_case_fluid_dofmap(case_iptr, dof_ptr, &
455 x_ptr, y_ptr, z_ptr, size) bind(c, name='neko_case_fluid_dofmap')
456 integer(c_intptr_t), intent(inout) :: case_iptr
457 type(case_t), pointer :: C
458 type(c_ptr) :: cptr
459 type(c_ptr), intent(inout) :: dof_ptr, x_ptr, y_ptr, z_ptr
460 integer, intent(inout) :: size
461
462 cptr = transfer(case_iptr, c_null_ptr)
463 if (c_associated(cptr)) then
464 call c_f_pointer(cptr, c)
465 call neko_api_wrap_dofmap(c%fluid%dm_Xh, dof_ptr, x_ptr, y_ptr, z_ptr)
466 size = c%fluid%dm_Xh%size()
467 else
468 call neko_error('Invalid Neko case')
469 end if
470
471 end subroutine neko_api_case_fluid_dofmap
472
474 subroutine neko_api_wrap_dofmap(dm, dof_ptr, x_ptr, y_ptr, z_ptr)
475 type(dofmap_t), target, intent(in) :: dm
476 type(c_ptr), intent(inout) :: dof_ptr, x_ptr, y_ptr, z_ptr
477
478 dof_ptr = c_loc(dm%dof)
479 x_ptr = c_loc(dm%x)
480 y_ptr = c_loc(dm%y)
481 z_ptr = c_loc(dm%z)
482
483 end subroutine neko_api_wrap_dofmap
484
498 subroutine neko_api_field_space(field_name, lx, zg, &
499 dr_inv, ds_inv, dt_inv, wx, wy, wz, dx, dy, dz) &
500 bind(c, name='neko_field_space')
501 character(kind=c_char), dimension(*), intent(in) :: field_name
502 integer, intent(inout) :: lx
503 type(c_ptr), intent(inout) :: zg, dr_inv, ds_inv, dt_inv
504 type(c_ptr), intent(inout) :: wx, wy, wz, dx, dy, dz
505 character(len=8192) :: name
506 type(field_t), pointer :: field
507 integer :: len
508
509 len = 0
510 do
511 if (field_name(len+1) .eq. c_null_char) exit
512 len = len + 1
513 name(len:len) = field_name(len)
514 end do
515
516 field => neko_registry%get_field(trim(name(1:len)))
517 call neko_api_wrap_space(field%Xh, lx, zg, dr_inv, ds_inv, dt_inv, &
518 wx, wy, wz, dx, dy, dz)
519
520 end subroutine neko_api_field_space
521
535 subroutine neko_api_case_fluid_space(case_iptr, lx, zg, &
536 dr_inv, ds_inv, dt_inv, wx, wy, wz, dx, dy, dz) &
537 bind(c, name='neko_case_fluid_space')
538 integer(c_intptr_t), intent(inout) :: case_iptr
539 type(case_t), pointer :: C
540 type(c_ptr) :: cptr
541 integer, intent(inout) :: lx
542 type(c_ptr), intent(inout) :: zg, dr_inv, ds_inv, dt_inv
543 type(c_ptr), intent(inout) :: wx, wy, wz, dx, dy, dz
544
545
546 cptr = transfer(case_iptr, c_null_ptr)
547 if (c_associated(cptr)) then
548 call c_f_pointer(cptr, c)
549 call neko_api_wrap_space(c%fluid%Xh, lx, zg, dr_inv, ds_inv, dt_inv, &
550 wx, wy, wz, dx, dy, dz)
551 else
552 call neko_error('Invalid Neko case')
553 end if
554
555 end subroutine neko_api_case_fluid_space
556
558 subroutine neko_api_wrap_space(Xh, lx, zg, dr_inv, ds_inv, dt_inv, &
559 wx, wy, wz, dx, dy, dz)
560 type(space_t), target, intent(in) :: Xh
561 integer, intent(inout) :: lx
562 type(c_ptr), intent(inout) :: zg, dr_inv, ds_inv, dt_inv
563 type(c_ptr), intent(inout) :: wx, wy, wz, dx, dy, dz
564
565 lx = xh%lx
566 zg = c_loc(xh%zg)
567 dr_inv = c_loc(xh%dr_inv)
568 ds_inv = c_loc(xh%ds_inv)
569 dt_inv = c_loc(xh%dt_inv)
570 wx = c_loc(xh%wx)
571 wy = c_loc(xh%wy)
572 wz = c_loc(xh%wz)
573 dx = c_loc(xh%dx)
574 dy = c_loc(xh%dy)
575 dz = c_loc(xh%dz)
576
577 end subroutine neko_api_wrap_space
578
581 subroutine neko_api_case_fluid_coef(case_iptr, G11, G22, G33, G12, G13, G23, &
582 mult, dxdr, dydr, dzdr, dxds, dyds, dzds, dxdt, dydt, dzdt, &
583 drdx, drdy, drdz, dsdx, dsdy, dsdz, dtdx, dtdy, dtdz, &
584 jac, B, area, nx, ny, nz) bind(c, name='neko_case_fluid_coef')
585 integer(c_intptr_t), intent(inout) :: case_iptr
586 type(case_t), pointer :: C
587 type(c_ptr) :: cptr
588 type(c_ptr), intent(inout) :: G11, G22, G33, G12, G13, G23
589 type(c_ptr), intent(inout) :: mult
590 type(c_ptr), intent(inout) :: dxdr, dydr, dzdr
591 type(c_ptr), intent(inout) :: dxds, dyds, dzds
592 type(c_ptr), intent(inout) :: dxdt, dydt, dzdt
593 type(c_ptr), intent(inout) :: drdx, drdy, drdz
594 type(c_ptr), intent(inout) :: dsdx, dsdy, dsdz
595 type(c_ptr), intent(inout) :: dtdx, dtdy, dtdz
596 type(c_ptr), intent(inout) :: jac, B, area, nx, ny, nz
597
598
599 cptr = transfer(case_iptr, c_null_ptr)
600 if (c_associated(cptr)) then
601 call c_f_pointer(cptr, c)
602 g11 = c_loc(c%fluid%c_Xh%G11)
603 g22 = c_loc(c%fluid%c_Xh%G22)
604 g33 = c_loc(c%fluid%c_Xh%G33)
605 g12 = c_loc(c%fluid%c_Xh%G12)
606 g13 = c_loc(c%fluid%c_Xh%G13)
607 g23 = c_loc(c%fluid%c_Xh%G23)
608 mult = c_loc(c%fluid%c_Xh%mult)
609 dxdr = c_loc(c%fluid%c_Xh%dxdr)
610 dydr = c_loc(c%fluid%c_Xh%dydr)
611 dzdr = c_loc(c%fluid%c_Xh%dzdr)
612 dxds = c_loc(c%fluid%c_Xh%dxds)
613 dyds = c_loc(c%fluid%c_Xh%dyds)
614 dzds = c_loc(c%fluid%c_Xh%dzds)
615 dxdt = c_loc(c%fluid%c_Xh%dxdt)
616 dydt = c_loc(c%fluid%c_Xh%dydt)
617 dzdt = c_loc(c%fluid%c_Xh%dzdt)
618 drdx = c_loc(c%fluid%c_Xh%drdx)
619 drdy = c_loc(c%fluid%c_Xh%drdy)
620 drdz = c_loc(c%fluid%c_Xh%drdz)
621 dsdx = c_loc(c%fluid%c_Xh%dsdx)
622 dsdy = c_loc(c%fluid%c_Xh%dsdy)
623 dsdz = c_loc(c%fluid%c_Xh%dsdz)
624 dtdx = c_loc(c%fluid%c_Xh%dtdx)
625 dtdy = c_loc(c%fluid%c_Xh%dtdy)
626 dtdz = c_loc(c%fluid%c_Xh%dtdz)
627 jac = c_loc(c%fluid%c_Xh%jac)
628 b = c_loc(c%fluid%c_Xh%B)
629 area = c_loc(c%fluid%c_Xh%area)
630 nx = c_loc(c%fluid%c_Xh%nx)
631 ny = c_loc(c%fluid%c_Xh%ny)
632 nz = c_loc(c%fluid%c_Xh%nz)
633 else
634 call neko_error('Invalid Neko case')
635 end if
636
637 end subroutine neko_api_case_fluid_coef
638
647 subroutine neko_api_user_setup(case_iptr, initial_cb, preprocess_cb, &
648 compute_cb, dirichlet_cb, material_cb, source_cb) &
649 bind(c, name='neko_user_setup')
650 integer(c_intptr_t), intent(inout) :: case_iptr
651 type(c_funptr), value :: initial_cb, preprocess_cb, compute_cb
652 type(c_funptr), value :: dirichlet_cb, material_cb, source_cb
653 type(case_t), pointer :: C
654 type(c_ptr) :: cptr
655
656 cptr = transfer(case_iptr, c_null_ptr)
657 if (c_associated(cptr)) then
658 call c_f_pointer(cptr, c)
659 call neko_api_user_cb_register(c%user, initial_cb, preprocess_cb, &
660 compute_cb, dirichlet_cb, material_cb, source_cb)
661 else
662 call neko_error('Invalid Neko case')
663 end if
664
665 end subroutine neko_api_user_setup
666
669 function neko_api_user_cb_field_by_name(field_name) result(field_ptr) &
670 bind(c, name='neko_cb_field_by_name')
671 character(kind=c_char), dimension(*), intent(in) :: field_name
672 character(len=8192) :: name
673 type(field_t), pointer :: field
674 type(c_ptr) :: field_ptr
675 integer :: len
676
677 len = 0
678 do
679 if (field_name(len+1) .eq. c_null_char) exit
680 len = len + 1
681 name(len:len) = field_name(len)
682 end do
683
684 field => neko_api_user_cb_get_field(trim(name(1:len)))
685
686 field_ptr = c_loc(field%x)
687
689
693 result(field_ptr) bind(c, name='neko_cb_field_by_index')
694 integer, intent(in) :: field_idx
695 type(field_t), pointer :: field
696 type(c_ptr) :: field_ptr
697
698 field => neko_api_user_cb_get_field(field_idx)
699
700 field_ptr = c_loc(field%x)
701
703
707 function neko_api_user_cb_field_name_at_index(field_idx, field_name) &
708 result(same_name) bind(c, name='neko_cb_field_name_at_index')
709 integer, intent(in) :: field_idx
710 character(kind=c_char), dimension(*), intent(in) :: field_name
711 character(len=8192) :: name
712 type(field_t), pointer :: f1, f2
713 type(c_ptr) :: field_ptr
714 integer :: len
715 logical(c_bool) :: same_name
716
717 len = 0
718 do
719 if (field_name(len+1) .eq. c_null_char) exit
720 len = len + 1
721 name(len:len) = field_name(len)
722 end do
723
724 f1 => neko_api_user_cb_get_field(field_idx)
725 f2 => neko_api_user_cb_get_field(trim(name(1:len)))
726
727 same_name = trim(f1%name) .eq. trim(f2%name)
728
730
731
732end 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:475
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:560
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:650
subroutine neko_api_step(case_iptr)
Compute a time-step for a neko case.
Definition neko_api.f90:263
integer(c_int) function neko_api_field_nelements(field_name)
Retrive the number of elements in a field.
Definition neko_api.f90:379
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:501
subroutine neko_api_case_free(case_iptr)
Destroy a Neko case.
Definition neko_api.f90:168
subroutine neko_api_case_allocate(case_iptr)
Allocate memory for a Neko case.
Definition neko_api.f90:105
subroutine neko_api_output_ctrl_execute(case_iptr, force_output)
Execute the Case's output controller.
Definition neko_api.f90:310
subroutine neko_api_case_init(case_json, case_len, case_iptr)
Initalise a Neko case.
Definition neko_api.f90:120
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:429
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:538
real(kind=c_rp) function neko_api_case_end_time(case_iptr)
Retrive the end time of a case.
Definition neko_api.f90:207
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:694
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:709
type(c_ptr) function neko_api_field(field_name)
Retrive a pointer to a flow field.
Definition neko_api.f90:333
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:402
subroutine neko_api_solve(case_iptr)
Solve a neko case.
Definition neko_api.f90:245
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:585
subroutine neko_api_registry_free()
Destroy a Neko field registry.
Definition neko_api.f90:96
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:671
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:456
integer(c_int) function neko_api_field_order(field_name)
Retrive the order of a field.
Definition neko_api.f90:356
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:187
integer(c_int) function neko_api_case_tstep(case_iptr)
Retrive the time-step of a case.
Definition neko_api.f90:227
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)