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
139 cp = transfer(case_iptr, c_null_ptr)
140 if (c_associated(cp))
then
141 call c_f_pointer(cp, c)
145 case_iptr = transfer(cp, 0_c_intptr_t)
150 if (c_associated(case_json))
then
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)
165 call case_init(c, json_case)
170 call neko_simcomps%init(c)
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
277 type(json_file) :: dt_params
279 real(kind=dp),
save :: step_loop_start = 0.0_dp
281 cptr = transfer(case_iptr, c_null_ptr)
282 if (c_associated(cptr))
then
283 call c_f_pointer(cptr, c)
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)
291 if (c%time%tstep .eq. 0)
then
292 call simulation_init(c, dt_controller)
294 step_loop_start = mpi_wtime()
297 if (.not. c%time%is_done())
then
298 call simulation_step(c, dt_controller, step_loop_start)
301 if (c%time%is_done())
then
302 call simulation_finalize(c)
303 if (
allocated(dt_controller))
then
304 deallocate(dt_controller)
309 call neko_error(
'Invalid Neko case')
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
447 if (field_name(len+1) .eq. c_null_char)
exit
449 name(len:len) = field_name(len)
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)
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
469 type(c_ptr),
intent(inout) :: dof_ptr, x_ptr, y_ptr, z_ptr
470 integer,
intent(inout) :: size
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()
478 call neko_error(
'Invalid Neko case')
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
521 if (field_name(len+1) .eq. c_null_char)
exit
523 name(len:len) = field_name(len)
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)
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
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
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)
562 call neko_error(
'Invalid Neko case')
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
577 dr_inv = c_loc(xh%dr_inv)
578 ds_inv = c_loc(xh%ds_inv)
579 dt_inv = c_loc(xh%dt_inv)
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
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
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)
644 call neko_error(
'Invalid Neko case')
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
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)
672 call neko_error(
'Invalid Neko case')
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
725 logical(c_bool) :: same_name
729 if (field_name(len+1) .eq. c_null_char)
exit
731 name(len:len) = field_name(len)
734 f1 => neko_api_user_cb_get_field(field_idx)
735 f2 => neko_api_user_cb_get_field(trim(name(1:len)))
737 same_name = trim(f1%name) .eq. trim(f2%name)