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_job_info() bind(c, name="neko_job_info")
98 logical :: initialized
99
100 call mpi_initialized(initialized)
101
102 if (.not.initialized) then
103 call neko_warning('Neko has not been initialised')
104 else
105 call neko_job_info()
106 call neko_log%newline()
107 end if
108
109 end subroutine neko_api_job_info
110
113 subroutine neko_api_case_allocate(case_iptr) &
114 bind(c, name="neko_case_allocate")
115 integer(c_intptr_t), intent(inout) :: case_iptr
116 type(case_t), pointer :: C
117 type(c_ptr) :: cp
118
119 allocate(c)
120 cp = c_loc(c)
121 case_iptr = transfer(cp, 0_c_intptr_t)
122
123 end subroutine neko_api_case_allocate
124
128 subroutine neko_api_case_init(case_json, case_len, case_iptr) &
129 bind(c, name="neko_case_init")
130 type(c_ptr), intent(in) :: case_json
131 integer(c_int), value :: case_len
132 integer(c_intptr_t), intent(inout) :: case_iptr
133 type(json_file) :: json_case
134 type(case_t), pointer :: C
135 type(c_ptr) :: cp
136
137 ! Check if the case has already been allocated
138 ! e.g. if a user callback has been injected
139 cp = transfer(case_iptr, c_null_ptr)
140 if (c_associated(cp)) then
141 call c_f_pointer(cp, c)
142 else
143 allocate(c)
144 cp = c_loc(c)
145 case_iptr = transfer(cp, 0_c_intptr_t)
146 end if
147
148 ! Convert passed in serialised JSON object into a Fortran
149 ! character string and create a json_file object
150 if (c_associated(case_json)) then
151 block
152 character(kind=c_char,len=case_len+1),pointer :: s
153 character(len=:), allocatable :: fcase_json
154 call c_f_pointer(case_json, s)
155 fcase_json = s(1:case_len)
156 call json_case%load_from_string(fcase_json)
157 deallocate(fcase_json)
158 nullify(s)
159 end block
160 end if
161
162 !
163 ! Create case
164 !
165 call case_init(c, json_case)
166
167 !
168 ! Create simulation components
169 !
170 call neko_simcomps%init(c)
171
172 end subroutine neko_api_case_init
173
174
177 subroutine neko_api_case_free(case_iptr) bind(c, name="neko_case_free")
178 integer(c_intptr_t), intent(inout) :: case_iptr
179 type(case_t), pointer :: C
180 type(c_ptr) :: cp
181
182 cp = transfer(case_iptr, c_null_ptr)
183 if (c_associated(cp)) then
184 call c_f_pointer(cp, c)
185 call case_free(c)
186 else
187 call neko_error('Invalid Neko case')
188 end if
189
190 end subroutine neko_api_case_free
191
195 function neko_api_case_time(case_iptr) result(time) &
196 bind(c, name="neko_case_time")
197 integer(c_intptr_t), intent(inout) :: case_iptr
198 real(kind=c_rp) :: time
199 type(case_t), pointer :: c
200 type(c_ptr) :: cptr
201
202 cptr = transfer(case_iptr, c_null_ptr)
203 if (c_associated(cptr)) then
204 call c_f_pointer(cptr, c)
205 time = c%time%t
206 else
207 call neko_error('Invalid Neko case')
208 end if
209
210 end function neko_api_case_time
211
215 function neko_api_case_end_time(case_iptr) result(end_time) &
216 bind(c, name="neko_case_end_time")
217 integer(c_intptr_t), intent(inout) :: case_iptr
218 real(kind=c_rp) :: end_time
219 type(case_t), pointer :: c
220 type(c_ptr) :: cptr
221
222 cptr = transfer(case_iptr, c_null_ptr)
223 if (c_associated(cptr)) then
224 call c_f_pointer(cptr, c)
225 end_time = c%time%end_time
226 else
227 call neko_error('Invalid Neko case')
228 end if
229
230 end function neko_api_case_end_time
231
235 function neko_api_case_tstep(case_iptr) result(tstep) &
236 bind(c, name="neko_case_tstep")
237 integer(c_intptr_t), intent(inout) :: case_iptr
238 type(case_t), pointer :: c
239 type(c_ptr) :: cptr
240 integer(c_int) :: tstep
241
242 cptr = transfer(case_iptr, c_null_ptr)
243 if (c_associated(cptr)) then
244 call c_f_pointer(cptr, c)
245 tstep = c%time%tstep
246 else
247 call neko_error('Invalid Neko case')
248 end if
249
250 end function neko_api_case_tstep
251
254 subroutine neko_api_solve(case_iptr) bind(c, name="neko_solve")
255 integer(c_intptr_t), intent(inout) :: case_iptr
256 type(case_t), pointer :: C
257 type(c_ptr) :: cp
258
259 cp = transfer(case_iptr, c_null_ptr)
260 if (c_associated(cp)) then
261 call c_f_pointer(cp, c)
262 call neko_solve(c)
263 else
264 call neko_error('Invalid Neko case')
265 end if
266
267 end subroutine neko_api_solve
268
271 subroutine neko_api_step(case_iptr) &
272 bind(c, name="neko_step")
274 integer(c_intptr_t), intent(inout) :: case_iptr
275 type(case_t), pointer :: C
276 type(c_ptr) :: cptr
277 type(json_file) :: dt_params
278 type(time_step_controller_t), save, allocatable :: dt_controller
279 real(kind=dp), save :: step_loop_start = 0.0_dp
280
281 cptr = transfer(case_iptr, c_null_ptr)
282 if (c_associated(cptr)) then
283 call c_f_pointer(cptr, c)
284
285 if (.not. allocated(dt_controller)) then
286 allocate(dt_controller)
287 call json_extract_object(c%params, 'case.time', dt_params)
288 call dt_controller%init(dt_params)
289 end if
290
291 if (c%time%tstep .eq. 0) then
292 call simulation_init(c, dt_controller)
293
294 step_loop_start = mpi_wtime()
295 end if
296
297 if (.not. c%time%is_done()) then
298 call simulation_step(c, dt_controller, step_loop_start)
299 end if
300
301 if (c%time%is_done()) then
302 call simulation_finalize(c)
303 if (allocated(dt_controller)) then
304 deallocate(dt_controller)
305 end if
306 end if
307
308 else
309 call neko_error('Invalid Neko case')
310 end if
311
312 end subroutine neko_api_step
313
318 subroutine neko_api_output_ctrl_execute(case_iptr, force_output) &
319 bind(c, name="neko_output_ctrl_execute")
320 integer(c_intptr_t), intent(inout) :: case_iptr
321 logical(kind=c_bool), value :: force_output
322 logical :: f_force_output
323 type(case_t), pointer :: C
324 type(c_ptr) :: cp
325 type(time_state_t) :: f_time
326
327 cp = transfer(case_iptr, c_null_ptr)
328 if (c_associated(cp)) then
329 call c_f_pointer(cp, c)
330 else
331 call neko_error('Invalid Neko case')
332 end if
333
334 f_force_output = transfer(force_output, f_force_output)
335 call c%output_controller%execute(c%time, f_force_output)
336
337 end subroutine neko_api_output_ctrl_execute
338
341 function neko_api_field(field_name) result(field_ptr) &
342 bind(c, name='neko_field')
343 character(kind=c_char), dimension(*), intent(in) :: field_name
344 character(len=8192) :: name
345 type(field_t), pointer :: field
346 type(c_ptr) :: field_ptr
347 integer :: len
348
349 len = 0
350 do
351 if (field_name(len+1) .eq. c_null_char) exit
352 len = len + 1
353 name(len:len) = field_name(len)
354 end do
355
356 field => neko_field_registry%get_field(trim(name(1:len)))
357
358 field_ptr = c_loc(field%x)
359
360 end function neko_api_field
361
364 function neko_api_field_order(field_name) result(field_lx) &
365 bind(c, name='neko_field_order')
366 character(kind=c_char), dimension(*), intent(in) :: field_name
367 character(len=8192) :: name
368 type(field_t), pointer :: field
369 integer(c_int) :: field_lx
370 integer :: len
371
372 len = 0
373 do
374 if (field_name(len+1) .eq. c_null_char) exit
375 len = len + 1
376 name(len:len) = field_name(len)
377 end do
378
379 field => neko_field_registry%get_field(trim(name(1:len)))
380
381 field_lx = field%Xh%lx
382
383 end function neko_api_field_order
384
387 function neko_api_field_nelements(field_name) result(field_nelv) &
388 bind(c, name='neko_field_nelements')
389 character(kind=c_char), dimension(*), intent(in) :: field_name
390 character(len=8192) :: name
391 type(field_t), pointer :: field
392 integer(c_int) :: field_nelv
393 integer :: len
394
395 len = 0
396 do
397 if (field_name(len+1) .eq. c_null_char) exit
398 len = len + 1
399 name(len:len) = field_name(len)
400 end do
401
402 field => neko_field_registry%get_field(trim(name(1:len)))
403
404 field_nelv = field%msh%nelv
405
406 end function neko_api_field_nelements
407
410 function neko_api_field_size(field_name) result(field_size) &
411 bind(c, name='neko_field_size')
412 character(kind=c_char), dimension(*), intent(in) :: field_name
413 character(len=8192) :: name
414 type(field_t), pointer :: field
415 integer(c_int) :: field_size
416 integer :: len
417
418 len = 0
419 do
420 if (field_name(len+1) .eq. c_null_char) exit
421 len = len + 1
422 name(len:len) = field_name(len)
423 end do
424
425 field => neko_field_registry%get_field(trim(name(1:len)))
426
427 field_size = field%dof%size()
428
429 end function neko_api_field_size
430
437 subroutine neko_api_field_dofmap(field_name, dof_ptr, x_ptr, y_ptr, z_ptr) &
438 bind(c, name='neko_field_dofmap')
439 character(kind=c_char), dimension(*), intent(in) :: field_name
440 type(c_ptr), intent(inout) :: dof_ptr, x_ptr, y_ptr, z_ptr
441 character(len=8192) :: name
442 type(field_t), pointer :: field
443 integer :: len
444
445 len = 0
446 do
447 if (field_name(len+1) .eq. c_null_char) exit
448 len = len + 1
449 name(len:len) = field_name(len)
450 end do
451
452 field => neko_field_registry%get_field(trim(name(1:len)))
453 call neko_api_wrap_dofmap(field%dof, dof_ptr, x_ptr, y_ptr, z_ptr)
454
455 end subroutine neko_api_field_dofmap
456
464 subroutine neko_api_case_fluid_dofmap(case_iptr, dof_ptr, &
465 x_ptr, y_ptr, z_ptr, size) bind(c, name='neko_case_fluid_dofmap')
466 integer(c_intptr_t), intent(inout) :: case_iptr
467 type(case_t), pointer :: C
468 type(c_ptr) :: cptr
469 type(c_ptr), intent(inout) :: dof_ptr, x_ptr, y_ptr, z_ptr
470 integer, intent(inout) :: size
471
472 cptr = transfer(case_iptr, c_null_ptr)
473 if (c_associated(cptr)) then
474 call c_f_pointer(cptr, c)
475 call neko_api_wrap_dofmap(c%fluid%dm_Xh, dof_ptr, x_ptr, y_ptr, z_ptr)
476 size = c%fluid%dm_Xh%size()
477 else
478 call neko_error('Invalid Neko case')
479 end if
480
481 end subroutine neko_api_case_fluid_dofmap
482
484 subroutine neko_api_wrap_dofmap(dm, dof_ptr, x_ptr, y_ptr, z_ptr)
485 type(dofmap_t), target, intent(in) :: dm
486 type(c_ptr), intent(inout) :: dof_ptr, x_ptr, y_ptr, z_ptr
487
488 dof_ptr = c_loc(dm%dof)
489 x_ptr = c_loc(dm%x)
490 y_ptr = c_loc(dm%y)
491 z_ptr = c_loc(dm%z)
492
493 end subroutine neko_api_wrap_dofmap
494
508 subroutine neko_api_field_space(field_name, lx, zg, &
509 dr_inv, ds_inv, dt_inv, wx, wy, wz, dx, dy, dz) &
510 bind(c, name='neko_field_space')
511 character(kind=c_char), dimension(*), intent(in) :: field_name
512 integer, intent(inout) :: lx
513 type(c_ptr), intent(inout) :: zg, dr_inv, ds_inv, dt_inv
514 type(c_ptr), intent(inout) :: wx, wy, wz, dx, dy, dz
515 character(len=8192) :: name
516 type(field_t), pointer :: field
517 integer :: len
518
519 len = 0
520 do
521 if (field_name(len+1) .eq. c_null_char) exit
522 len = len + 1
523 name(len:len) = field_name(len)
524 end do
525
526 field => neko_field_registry%get_field(trim(name(1:len)))
527 call neko_api_wrap_space(field%Xh, lx, zg, dr_inv, ds_inv, dt_inv, &
528 wx, wy, wz, dx, dy, dz)
529
530 end subroutine neko_api_field_space
531
545 subroutine neko_api_case_fluid_space(case_iptr, lx, zg, &
546 dr_inv, ds_inv, dt_inv, wx, wy, wz, dx, dy, dz) &
547 bind(c, name='neko_case_fluid_space')
548 integer(c_intptr_t), intent(inout) :: case_iptr
549 type(case_t), pointer :: C
550 type(c_ptr) :: cptr
551 integer, intent(inout) :: lx
552 type(c_ptr), intent(inout) :: zg, dr_inv, ds_inv, dt_inv
553 type(c_ptr), intent(inout) :: wx, wy, wz, dx, dy, dz
554
555
556 cptr = transfer(case_iptr, c_null_ptr)
557 if (c_associated(cptr)) then
558 call c_f_pointer(cptr, c)
559 call neko_api_wrap_space(c%fluid%Xh, lx, zg, dr_inv, ds_inv, dt_inv, &
560 wx, wy, wz, dx, dy, dz)
561 else
562 call neko_error('Invalid Neko case')
563 end if
564
565 end subroutine neko_api_case_fluid_space
566
568 subroutine neko_api_wrap_space(Xh, lx, zg, dr_inv, ds_inv, dt_inv, &
569 wx, wy, wz, dx, dy, dz)
570 type(space_t), target, intent(in) :: Xh
571 integer, intent(inout) :: lx
572 type(c_ptr), intent(inout) :: zg, dr_inv, ds_inv, dt_inv
573 type(c_ptr), intent(inout) :: wx, wy, wz, dx, dy, dz
574
575 lx = xh%lx
576 zg = c_loc(xh%zg)
577 dr_inv = c_loc(xh%dr_inv)
578 ds_inv = c_loc(xh%ds_inv)
579 dt_inv = c_loc(xh%dt_inv)
580 wx = c_loc(xh%wx)
581 wy = c_loc(xh%wy)
582 wz = c_loc(xh%wz)
583 dx = c_loc(xh%dx)
584 dy = c_loc(xh%dy)
585 dz = c_loc(xh%dz)
586
587 end subroutine neko_api_wrap_space
588
591 subroutine neko_api_case_fluid_coef(case_iptr, G11, G22, G33, G12, G13, G23, &
592 mult, dxdr, dydr, dzdr, dxds, dyds, dzds, dxdt, dydt, dzdt, &
593 drdx, drdy, drdz, dsdx, dsdy, dsdz, dtdx, dtdy, dtdz, &
594 jac, B, area, nx, ny, nz) bind(c, name='neko_case_fluid_coef')
595 integer(c_intptr_t), intent(inout) :: case_iptr
596 type(case_t), pointer :: C
597 type(c_ptr) :: cptr
598 type(c_ptr), intent(inout) :: G11, G22, G33, G12, G13, G23
599 type(c_ptr), intent(inout) :: mult
600 type(c_ptr), intent(inout) :: dxdr, dydr, dzdr
601 type(c_ptr), intent(inout) :: dxds, dyds, dzds
602 type(c_ptr), intent(inout) :: dxdt, dydt, dzdt
603 type(c_ptr), intent(inout) :: drdx, drdy, drdz
604 type(c_ptr), intent(inout) :: dsdx, dsdy, dsdz
605 type(c_ptr), intent(inout) :: dtdx, dtdy, dtdz
606 type(c_ptr), intent(inout) :: jac, B, area, nx, ny, nz
607
608
609 cptr = transfer(case_iptr, c_null_ptr)
610 if (c_associated(cptr)) then
611 call c_f_pointer(cptr, c)
612 g11 = c_loc(c%fluid%c_Xh%G11)
613 g22 = c_loc(c%fluid%c_Xh%G22)
614 g33 = c_loc(c%fluid%c_Xh%G33)
615 g12 = c_loc(c%fluid%c_Xh%G12)
616 g13 = c_loc(c%fluid%c_Xh%G13)
617 g23 = c_loc(c%fluid%c_Xh%G23)
618 mult = c_loc(c%fluid%c_Xh%mult)
619 dxdr = c_loc(c%fluid%c_Xh%dxdr)
620 dydr = c_loc(c%fluid%c_Xh%dydr)
621 dzdr = c_loc(c%fluid%c_Xh%dzdr)
622 dxds = c_loc(c%fluid%c_Xh%dxds)
623 dyds = c_loc(c%fluid%c_Xh%dyds)
624 dzds = c_loc(c%fluid%c_Xh%dzds)
625 dxdt = c_loc(c%fluid%c_Xh%dxdt)
626 dydt = c_loc(c%fluid%c_Xh%dydt)
627 dzdt = c_loc(c%fluid%c_Xh%dzdt)
628 drdx = c_loc(c%fluid%c_Xh%drdx)
629 drdy = c_loc(c%fluid%c_Xh%drdy)
630 drdz = c_loc(c%fluid%c_Xh%drdz)
631 dsdx = c_loc(c%fluid%c_Xh%dsdx)
632 dsdy = c_loc(c%fluid%c_Xh%dsdy)
633 dsdz = c_loc(c%fluid%c_Xh%dsdz)
634 dtdx = c_loc(c%fluid%c_Xh%dtdx)
635 dtdy = c_loc(c%fluid%c_Xh%dtdy)
636 dtdz = c_loc(c%fluid%c_Xh%dtdz)
637 jac = c_loc(c%fluid%c_Xh%jac)
638 b = c_loc(c%fluid%c_Xh%B)
639 area = c_loc(c%fluid%c_Xh%area)
640 nx = c_loc(c%fluid%c_Xh%nx)
641 ny = c_loc(c%fluid%c_Xh%ny)
642 nz = c_loc(c%fluid%c_Xh%nz)
643 else
644 call neko_error('Invalid Neko case')
645 end if
646
647 end subroutine neko_api_case_fluid_coef
648
657 subroutine neko_api_user_setup(case_iptr, initial_cb, preprocess_cb, &
658 compute_cb, dirichlet_cb, material_cb, source_cb) &
659 bind(c, name='neko_user_setup')
660 integer(c_intptr_t), intent(inout) :: case_iptr
661 type(c_funptr), value :: initial_cb, preprocess_cb, compute_cb
662 type(c_funptr), value :: dirichlet_cb, material_cb, source_cb
663 type(case_t), pointer :: C
664 type(c_ptr) :: cptr
665
666 cptr = transfer(case_iptr, c_null_ptr)
667 if (c_associated(cptr)) then
668 call c_f_pointer(cptr, c)
669 call neko_api_user_cb_register(c%user, initial_cb, preprocess_cb, &
670 compute_cb, dirichlet_cb, material_cb, source_cb)
671 else
672 call neko_error('Invalid Neko case')
673 end if
674
675 end subroutine neko_api_user_setup
676
679 function neko_api_user_cb_field_by_name(field_name) result(field_ptr) &
680 bind(c, name='neko_cb_field_by_name')
681 character(kind=c_char), dimension(*), intent(in) :: field_name
682 character(len=8192) :: name
683 type(field_t), pointer :: field
684 type(c_ptr) :: field_ptr
685 integer :: len
686
687 len = 0
688 do
689 if (field_name(len+1) .eq. c_null_char) exit
690 len = len + 1
691 name(len:len) = field_name(len)
692 end do
693
694 field => neko_api_user_cb_get_field(trim(name(1:len)))
695
696 field_ptr = c_loc(field%x)
697
698 end function neko_api_user_cb_field_by_name
699
702 function neko_api_user_cb_field_by_index(field_idx) &
703 result(field_ptr) bind(c, name='neko_cb_field_by_index')
704 integer, intent(in) :: field_idx
705 type(field_t), pointer :: field
706 type(c_ptr) :: field_ptr
707
708 field => neko_api_user_cb_get_field(field_idx)
709
710 field_ptr = c_loc(field%x)
711
712 end function neko_api_user_cb_field_by_index
713
717 function neko_api_user_cb_field_name_at_index(field_idx, field_name) &
718 result(same_name) bind(c, name='neko_cb_field_name_at_index')
719 integer, intent(in) :: field_idx
720 character(kind=c_char), dimension(*), intent(in) :: field_name
721 character(len=8192) :: name
722 type(field_t), pointer :: f1, f2
723 type(c_ptr) :: field_ptr
724 integer :: len
725 logical(c_bool) :: same_name
726
727 len = 0
728 do
729 if (field_name(len+1) .eq. c_null_char) exit
730 len = len + 1
731 name(len:len) = field_name(len)
732 end do
733
734 f1 => neko_api_user_cb_get_field(field_idx)
735 f2 => neko_api_user_cb_get_field(trim(name(1:len)))
736
737 same_name = trim(f1%name) .eq. trim(f2%name)
738
739 end function neko_api_user_cb_field_name_at_index
740
741
742end 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)