Neko 1.99.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-2025, 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
36 use, intrinsic :: iso_c_binding
37 implicit none
38 private
39
40 interface
41
49 module subroutine neko_api_user_cb_register(user, initial_cb, &
50 preprocess_cb, compute_cb, dirichlet_cb, material_cb, source_cb)
51 type(user_t), intent(inout) :: user
52 type(c_funptr), value :: initial_cb, preprocess_cb, compute_cb
53 type(c_funptr), value :: dirichlet_cb, material_cb, source_cb
54 end subroutine neko_api_user_cb_register
55 end interface
56
57 interface
58
60 module function neko_api_user_cb_get_field_by_name(field_name) result(f)
61 character(len=*), intent(in) :: field_name
62 type(field_t), pointer :: f
63 end function neko_api_user_cb_get_field_by_name
64 end interface
65
66 interface
67
69 module function neko_api_user_cb_get_field_by_index(field_idx) result(f)
70 integer, intent(in) :: field_idx
71 type(field_t), pointer :: f
72 end function neko_api_user_cb_get_field_by_index
73 end interface
74
75 interface neko_api_user_cb_get_field
76 module procedure neko_api_user_cb_get_field_by_name, &
77 neko_api_user_cb_get_field_by_index
78 end interface neko_api_user_cb_get_field
79
80contains
81
83 subroutine neko_api_init() bind(c, name="neko_init")
84
85 call neko_init()
86
87 end subroutine neko_api_init
88
90 subroutine neko_api_finalize() bind(c, name="neko_finalize")
91
92 call neko_finalize()
93
94 end subroutine neko_api_finalize
95
97 subroutine neko_api_device_init() bind(c, name="neko_device_init")
98
99 call device_init
100
101 end subroutine neko_api_device_init
102
104 subroutine neko_api_device_finalize() bind(c, name="neko_device_finalize")
105
106 call device_finalize
107
108 end subroutine neko_api_device_finalize
109
111 subroutine neko_api_job_info() bind(c, name="neko_job_info")
112 logical :: initialized
113
114 call mpi_initialized(initialized)
115
116 if (.not.initialized) then
117 call neko_warning('Neko has not been initialised')
118 else
119 call neko_job_info()
120 call neko_log%newline()
121 end if
122
123 end subroutine neko_api_job_info
124
126 subroutine neko_api_field_registry_init() bind(c, name="neko_field_registry_init")
127
128 call neko_field_registry%init()
129
130 end subroutine neko_api_field_registry_init
131
133 subroutine neko_api_field_registry_free() bind(c, name="neko_field_registry_free")
134
135 call neko_field_registry%free()
136
137 end subroutine neko_api_field_registry_free
138
141 subroutine neko_api_case_allocate(case_iptr) &
142 bind(c, name="neko_case_allocate")
143 integer(c_intptr_t), intent(inout) :: case_iptr
144 type(case_t), pointer :: C
145 type(c_ptr) :: cp
146
147 allocate(c)
148 cp = c_loc(c)
149 case_iptr = transfer(cp, 0_c_intptr_t)
150
151 end subroutine neko_api_case_allocate
152
156 subroutine neko_api_case_init(case_json, case_len, case_iptr) &
157 bind(c, name="neko_case_init")
158 type(c_ptr), intent(in) :: case_json
159 integer(c_int), value :: case_len
160 integer(c_intptr_t), intent(inout) :: case_iptr
161 type(json_file) :: json_case
162 type(case_t), pointer :: C
163 type(c_ptr) :: cp
164
165 ! Check if the case has already been allocated
166 ! e.g. if a user callback has been injected
167 cp = transfer(case_iptr, c_null_ptr)
168 if (c_associated(cp)) then
169 call c_f_pointer(cp, c)
170 else
171 allocate(c)
172 cp = c_loc(c)
173 case_iptr = transfer(cp, 0_c_intptr_t)
174 end if
175
176 ! Convert passed in serialised JSON object into a Fortran
177 ! character string and create a json_file object
178 if (c_associated(case_json)) then
179 block
180 character(kind=c_char,len=case_len+1),pointer :: s
181 character(len=:), allocatable :: fcase_json
182 call c_f_pointer(case_json, s)
183 fcase_json = s(1:case_len)
184 call json_case%load_from_string(fcase_json)
185 deallocate(fcase_json)
186 nullify(s)
187 end block
188 end if
189
190 !
191 ! Create case
192 !
193 call case_init(c, json_case)
194
195 !
196 ! Create simulation components
197 !
198 call neko_simcomps%init(c)
199
200 end subroutine neko_api_case_init
201
202
205 subroutine neko_api_case_free(case_iptr) bind(c, name="neko_case_free")
206 integer(c_intptr_t), intent(inout) :: case_iptr
207 type(case_t), pointer :: C
208 type(c_ptr) :: cp
209
210 cp = transfer(case_iptr, c_null_ptr)
211 if (c_associated(cp)) then
212 call c_f_pointer(cp, c)
213 call case_free(c)
214 else
215 call neko_error('Invalid Neko case')
216 end if
217
218 end subroutine neko_api_case_free
219
223 function neko_api_case_time(case_iptr) result(time) &
224 bind(c, name="neko_case_time")
225 integer(c_intptr_t), intent(inout) :: case_iptr
226 real(kind=c_rp) :: time
227 type(case_t), pointer :: c
228 type(c_ptr) :: cptr
229
230 cptr = transfer(case_iptr, c_null_ptr)
231 if (c_associated(cptr)) then
232 call c_f_pointer(cptr, c)
233 time = c%time%t
234 else
235 call neko_error('Invalid Neko case')
236 end if
237
238 end function neko_api_case_time
239
243 function neko_api_case_end_time(case_iptr) result(end_time) &
244 bind(c, name="neko_case_end_time")
245 integer(c_intptr_t), intent(inout) :: case_iptr
246 real(kind=c_rp) :: end_time
247 type(case_t), pointer :: c
248 type(c_ptr) :: cptr
249
250 cptr = transfer(case_iptr, c_null_ptr)
251 if (c_associated(cptr)) then
252 call c_f_pointer(cptr, c)
253 end_time = c%time%end_time
254 else
255 call neko_error('Invalid Neko case')
256 end if
257
258 end function neko_api_case_end_time
259
263 function neko_api_case_tstep(case_iptr) result(tstep) &
264 bind(c, name="neko_case_tstep")
265 integer(c_intptr_t), intent(inout) :: case_iptr
266 type(case_t), pointer :: c
267 type(c_ptr) :: cptr
268 integer(c_int) :: tstep
269
270 cptr = transfer(case_iptr, c_null_ptr)
271 if (c_associated(cptr)) then
272 call c_f_pointer(cptr, c)
273 tstep = c%time%tstep
274 else
275 call neko_error('Invalid Neko case')
276 end if
277
278 end function neko_api_case_tstep
279
282 subroutine neko_api_solve(case_iptr) bind(c, name="neko_solve")
283 integer(c_intptr_t), intent(inout) :: case_iptr
284 type(case_t), pointer :: C
285 type(c_ptr) :: cp
286
287 cp = transfer(case_iptr, c_null_ptr)
288 if (c_associated(cp)) then
289 call c_f_pointer(cp, c)
290 call neko_solve(c)
291 else
292 call neko_error('Invalid Neko case')
293 end if
294
295 end subroutine neko_api_solve
296
299 subroutine neko_api_step(case_iptr) &
300 bind(c, name="neko_step")
302 integer(c_intptr_t), intent(inout) :: case_iptr
303 type(case_t), pointer :: C
304 type(c_ptr) :: cptr
305 type(json_file) :: dt_params
306 type(time_step_controller_t), save, allocatable :: dt_controller
307 real(kind=dp), save :: step_loop_start = 0.0_dp
308
309 cptr = transfer(case_iptr, c_null_ptr)
310 if (c_associated(cptr)) then
311 call c_f_pointer(cptr, c)
312
313 if (.not. allocated(dt_controller)) then
314 allocate(dt_controller)
315 call json_get(c%params, 'case.time', dt_params)
316 call dt_controller%init(dt_params)
317 end if
318
319 if (c%time%tstep .eq. 0) then
320 call simulation_init(c, dt_controller)
321
322 step_loop_start = mpi_wtime()
323 end if
324
325 if (.not. c%time%is_done()) then
326 call simulation_step(c, dt_controller, step_loop_start)
327 end if
328
329 if (c%time%is_done()) then
330 call simulation_finalize(c)
331 if (allocated(dt_controller)) then
332 deallocate(dt_controller)
333 end if
334 end if
335
336 else
337 call neko_error('Invalid Neko case')
338 end if
339
340 end subroutine neko_api_step
341
346 subroutine neko_api_output_ctrl_execute(case_iptr, force_output) &
347 bind(c, name="neko_output_ctrl_execute")
348 integer(c_intptr_t), intent(inout) :: case_iptr
349 logical(kind=c_bool), value :: force_output
350 logical :: f_force_output
351 type(case_t), pointer :: C
352 type(c_ptr) :: cp
353 type(time_state_t) :: f_time
354
355 cp = transfer(case_iptr, c_null_ptr)
356 if (c_associated(cp)) then
357 call c_f_pointer(cp, c)
358 else
359 call neko_error('Invalid Neko case')
360 end if
361
362 f_force_output = transfer(force_output, f_force_output)
363 call c%output_controller%execute(c%time, f_force_output)
364
365 end subroutine neko_api_output_ctrl_execute
366
369 function neko_api_field(field_name) result(field_ptr) &
370 bind(c, name='neko_field')
371 character(kind=c_char), dimension(*), intent(in) :: field_name
372 character(len=8192) :: name
373 type(field_t), pointer :: field
374 type(c_ptr) :: field_ptr
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_field_registry%get_field(trim(name(1:len)))
385
386 field_ptr = c_loc(field%x)
387
388 end function neko_api_field
389
392 function neko_api_field_order(field_name) result(field_lx) &
393 bind(c, name='neko_field_order')
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_lx
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_field_registry%get_field(trim(name(1:len)))
408
409 field_lx = field%Xh%lx
410
411 end function neko_api_field_order
412
415 function neko_api_field_nelements(field_name) result(field_nelv) &
416 bind(c, name='neko_field_nelements')
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_nelv
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_field_registry%get_field(trim(name(1:len)))
431
432 field_nelv = field%msh%nelv
433
434 end function neko_api_field_nelements
435
438 function neko_api_field_size(field_name) result(field_size) &
439 bind(c, name='neko_field_size')
440 character(kind=c_char), dimension(*), intent(in) :: field_name
441 character(len=8192) :: name
442 type(field_t), pointer :: field
443 integer(c_int) :: field_size
444 integer :: len
445
446 len = 0
447 do
448 if (field_name(len+1) .eq. c_null_char) exit
449 len = len + 1
450 name(len:len) = field_name(len)
451 end do
452
453 field => neko_field_registry%get_field(trim(name(1:len)))
454
455 field_size = field%dof%size()
456
457 end function neko_api_field_size
458
465 subroutine neko_api_field_dofmap(field_name, dof_ptr, x_ptr, y_ptr, z_ptr) &
466 bind(c, name='neko_field_dofmap')
467 character(kind=c_char), dimension(*), intent(in) :: field_name
468 type(c_ptr), intent(inout) :: dof_ptr, x_ptr, y_ptr, z_ptr
469 character(len=8192) :: name
470 type(field_t), pointer :: field
471 integer :: len
472
473 len = 0
474 do
475 if (field_name(len+1) .eq. c_null_char) exit
476 len = len + 1
477 name(len:len) = field_name(len)
478 end do
479
480 field => neko_field_registry%get_field(trim(name(1:len)))
481 call neko_api_wrap_dofmap(field%dof, dof_ptr, x_ptr, y_ptr, z_ptr)
482
483 end subroutine neko_api_field_dofmap
484
492 subroutine neko_api_case_fluid_dofmap(case_iptr, dof_ptr, &
493 x_ptr, y_ptr, z_ptr, size) bind(c, name='neko_case_fluid_dofmap')
494 integer(c_intptr_t), intent(inout) :: case_iptr
495 type(case_t), pointer :: C
496 type(c_ptr) :: cptr
497 type(c_ptr), intent(inout) :: dof_ptr, x_ptr, y_ptr, z_ptr
498 integer, intent(inout) :: size
499
500 cptr = transfer(case_iptr, c_null_ptr)
501 if (c_associated(cptr)) then
502 call c_f_pointer(cptr, c)
503 call neko_api_wrap_dofmap(c%fluid%dm_Xh, dof_ptr, x_ptr, y_ptr, z_ptr)
504 size = c%fluid%dm_Xh%size()
505 else
506 call neko_error('Invalid Neko case')
507 end if
508
509 end subroutine neko_api_case_fluid_dofmap
510
512 subroutine neko_api_wrap_dofmap(dm, dof_ptr, x_ptr, y_ptr, z_ptr)
513 type(dofmap_t), target, intent(in) :: dm
514 type(c_ptr), intent(inout) :: dof_ptr, x_ptr, y_ptr, z_ptr
515
516 dof_ptr = c_loc(dm%dof)
517 x_ptr = c_loc(dm%x)
518 y_ptr = c_loc(dm%y)
519 z_ptr = c_loc(dm%z)
520
521 end subroutine neko_api_wrap_dofmap
522
536 subroutine neko_api_field_space(field_name, lx, zg, &
537 dr_inv, ds_inv, dt_inv, wx, wy, wz, dx, dy, dz) &
538 bind(c, name='neko_field_space')
539 character(kind=c_char), dimension(*), intent(in) :: field_name
540 integer, intent(inout) :: lx
541 type(c_ptr), intent(inout) :: zg, dr_inv, ds_inv, dt_inv
542 type(c_ptr), intent(inout) :: wx, wy, wz, dx, dy, dz
543 character(len=8192) :: name
544 type(field_t), pointer :: field
545 integer :: len
546
547 len = 0
548 do
549 if (field_name(len+1) .eq. c_null_char) exit
550 len = len + 1
551 name(len:len) = field_name(len)
552 end do
553
554 field => neko_field_registry%get_field(trim(name(1:len)))
555 call neko_api_wrap_space(field%Xh, lx, zg, dr_inv, ds_inv, dt_inv, &
556 wx, wy, wz, dx, dy, dz)
557
558 end subroutine neko_api_field_space
559
573 subroutine neko_api_case_fluid_space(case_iptr, lx, zg, &
574 dr_inv, ds_inv, dt_inv, wx, wy, wz, dx, dy, dz) &
575 bind(c, name='neko_case_fluid_space')
576 integer(c_intptr_t), intent(inout) :: case_iptr
577 type(case_t), pointer :: C
578 type(c_ptr) :: cptr
579 integer, intent(inout) :: lx
580 type(c_ptr), intent(inout) :: zg, dr_inv, ds_inv, dt_inv
581 type(c_ptr), intent(inout) :: wx, wy, wz, dx, dy, dz
582
583
584 cptr = transfer(case_iptr, c_null_ptr)
585 if (c_associated(cptr)) then
586 call c_f_pointer(cptr, c)
587 call neko_api_wrap_space(c%fluid%Xh, lx, zg, dr_inv, ds_inv, dt_inv, &
588 wx, wy, wz, dx, dy, dz)
589 else
590 call neko_error('Invalid Neko case')
591 end if
592
593 end subroutine neko_api_case_fluid_space
594
596 subroutine neko_api_wrap_space(Xh, lx, zg, dr_inv, ds_inv, dt_inv, &
597 wx, wy, wz, dx, dy, dz)
598 type(space_t), target, intent(in) :: Xh
599 integer, intent(inout) :: lx
600 type(c_ptr), intent(inout) :: zg, dr_inv, ds_inv, dt_inv
601 type(c_ptr), intent(inout) :: wx, wy, wz, dx, dy, dz
602
603 lx = xh%lx
604 zg = c_loc(xh%zg)
605 dr_inv = c_loc(xh%dr_inv)
606 ds_inv = c_loc(xh%ds_inv)
607 dt_inv = c_loc(xh%dt_inv)
608 wx = c_loc(xh%wx)
609 wy = c_loc(xh%wy)
610 wz = c_loc(xh%wz)
611 dx = c_loc(xh%dx)
612 dy = c_loc(xh%dy)
613 dz = c_loc(xh%dz)
614
615 end subroutine neko_api_wrap_space
616
619 subroutine neko_api_case_fluid_coef(case_iptr, G11, G22, G33, G12, G13, G23, &
620 mult, dxdr, dydr, dzdr, dxds, dyds, dzds, dxdt, dydt, dzdt, &
621 drdx, drdy, drdz, dsdx, dsdy, dsdz, dtdx, dtdy, dtdz, &
622 jac, B, area, nx, ny, nz) bind(c, name='neko_case_fluid_coef')
623 integer(c_intptr_t), intent(inout) :: case_iptr
624 type(case_t), pointer :: C
625 type(c_ptr) :: cptr
626 type(c_ptr), intent(inout) :: G11, G22, G33, G12, G13, G23
627 type(c_ptr), intent(inout) :: mult
628 type(c_ptr), intent(inout) :: dxdr, dydr, dzdr
629 type(c_ptr), intent(inout) :: dxds, dyds, dzds
630 type(c_ptr), intent(inout) :: dxdt, dydt, dzdt
631 type(c_ptr), intent(inout) :: drdx, drdy, drdz
632 type(c_ptr), intent(inout) :: dsdx, dsdy, dsdz
633 type(c_ptr), intent(inout) :: dtdx, dtdy, dtdz
634 type(c_ptr), intent(inout) :: jac, B, area, nx, ny, nz
635
636
637 cptr = transfer(case_iptr, c_null_ptr)
638 if (c_associated(cptr)) then
639 call c_f_pointer(cptr, c)
640 g11 = c_loc(c%fluid%c_Xh%G11)
641 g22 = c_loc(c%fluid%c_Xh%G22)
642 g33 = c_loc(c%fluid%c_Xh%G33)
643 g12 = c_loc(c%fluid%c_Xh%G12)
644 g13 = c_loc(c%fluid%c_Xh%G13)
645 g23 = c_loc(c%fluid%c_Xh%G23)
646 mult = c_loc(c%fluid%c_Xh%mult)
647 dxdr = c_loc(c%fluid%c_Xh%dxdr)
648 dydr = c_loc(c%fluid%c_Xh%dydr)
649 dzdr = c_loc(c%fluid%c_Xh%dzdr)
650 dxds = c_loc(c%fluid%c_Xh%dxds)
651 dyds = c_loc(c%fluid%c_Xh%dyds)
652 dzds = c_loc(c%fluid%c_Xh%dzds)
653 dxdt = c_loc(c%fluid%c_Xh%dxdt)
654 dydt = c_loc(c%fluid%c_Xh%dydt)
655 dzdt = c_loc(c%fluid%c_Xh%dzdt)
656 drdx = c_loc(c%fluid%c_Xh%drdx)
657 drdy = c_loc(c%fluid%c_Xh%drdy)
658 drdz = c_loc(c%fluid%c_Xh%drdz)
659 dsdx = c_loc(c%fluid%c_Xh%dsdx)
660 dsdy = c_loc(c%fluid%c_Xh%dsdy)
661 dsdz = c_loc(c%fluid%c_Xh%dsdz)
662 dtdx = c_loc(c%fluid%c_Xh%dtdx)
663 dtdy = c_loc(c%fluid%c_Xh%dtdy)
664 dtdz = c_loc(c%fluid%c_Xh%dtdz)
665 jac = c_loc(c%fluid%c_Xh%jac)
666 b = c_loc(c%fluid%c_Xh%B)
667 area = c_loc(c%fluid%c_Xh%area)
668 nx = c_loc(c%fluid%c_Xh%nx)
669 ny = c_loc(c%fluid%c_Xh%ny)
670 nz = c_loc(c%fluid%c_Xh%nz)
671 else
672 call neko_error('Invalid Neko case')
673 end if
674
675 end subroutine neko_api_case_fluid_coef
676
685 subroutine neko_api_user_setup(case_iptr, initial_cb, preprocess_cb, &
686 compute_cb, dirichlet_cb, material_cb, source_cb) &
687 bind(c, name='neko_user_setup')
688 integer(c_intptr_t), intent(inout) :: case_iptr
689 type(c_funptr), value :: initial_cb, preprocess_cb, compute_cb
690 type(c_funptr), value :: dirichlet_cb, material_cb, source_cb
691 type(case_t), pointer :: C
692 type(c_ptr) :: cptr
693
694 cptr = transfer(case_iptr, c_null_ptr)
695 if (c_associated(cptr)) then
696 call c_f_pointer(cptr, c)
697 call neko_api_user_cb_register(c%user, initial_cb, preprocess_cb, &
698 compute_cb, dirichlet_cb, material_cb, source_cb)
699 else
700 call neko_error('Invalid Neko case')
701 end if
702
703 end subroutine neko_api_user_setup
704
707 function neko_api_user_cb_field_by_name(field_name) result(field_ptr) &
708 bind(c, name='neko_cb_field_by_name')
709 character(kind=c_char), dimension(*), intent(in) :: field_name
710 character(len=8192) :: name
711 type(field_t), pointer :: field
712 type(c_ptr) :: field_ptr
713 integer :: len
714
715 len = 0
716 do
717 if (field_name(len+1) .eq. c_null_char) exit
718 len = len + 1
719 name(len:len) = field_name(len)
720 end do
721
722 field => neko_api_user_cb_get_field(trim(name(1:len)))
723
724 field_ptr = c_loc(field%x)
725
726 end function neko_api_user_cb_field_by_name
727
730 function neko_api_user_cb_field_by_index(field_idx) &
731 result(field_ptr) bind(c, name='neko_cb_field_by_index')
732 integer, intent(in) :: field_idx
733 type(field_t), pointer :: field
734 type(c_ptr) :: field_ptr
735
736 field => neko_api_user_cb_get_field(field_idx)
737
738 field_ptr = c_loc(field%x)
739
740 end function neko_api_user_cb_field_by_index
741
745 function neko_api_user_cb_field_name_at_index(field_idx, field_name) &
746 result(same_name) bind(c, name='neko_cb_field_name_at_index')
747 integer, intent(in) :: field_idx
748 character(kind=c_char), dimension(*), intent(in) :: field_name
749 character(len=8192) :: name
750 type(field_t), pointer :: f1, f2
751 type(c_ptr) :: field_ptr
752 integer :: len
753 logical(c_bool) :: same_name
754
755 len = 0
756 do
757 if (field_name(len+1) .eq. c_null_char) exit
758 len = len + 1
759 name(len:len) = field_name(len)
760 end do
761
762 f1 => neko_api_user_cb_get_field(field_idx)
763 f2 => neko_api_user_cb_get_field(trim(name(1:len)))
764
765 same_name = trim(f1%name) .eq. trim(f2%name)
766
767 end function neko_api_user_cb_field_name_at_index
768
769
770end module neko_api
Defines a field.
Definition field.f90:34
Neko C API.
Definition neko_api.f90:34
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)