170 type(
coef_t),
intent(inout) :: coef
171 type(json_file),
intent(inout) :: json
172 type(json_file) :: body_sub, bc_subdict
173 type(json_file) :: precon_params
175 type(
user_t),
intent(in) :: user
176 integer,
allocatable :: zone_indices(:)
177 integer :: time_order
178 integer :: n_moving_zones
179 integer :: z, tmp_int, ksp_max_iter
180 integer,
allocatable :: moving_zone_ids(:)
181 integer :: i, j, k, n_bcs, n, n_bodies
182 real(kind=
rp),
allocatable :: tmp_vec(:)
183 real(kind=
rp) :: tmp_val, abstol
184 character(len=128) :: log_buf
185 character(len=256) :: log_buf_l
186 character(len=:),
allocatable :: bc_type
187 character(len=:),
allocatable :: tmp_str
188 character(len=:),
allocatable :: ksp_solver
189 character(len=:),
allocatable :: precon_type
190 logical :: tmp_logical, oifs
192 logical :: found_zone
193 logical :: has_user_kin, has_user_mesh
194 logical :: has_builtin_osc, has_builtin_rot, is_rot_active
195 logical :: res_monitor
197 if (json%valid_path(
'case.fluid.ale'))
then
198 call json_get(json,
'case.fluid.ale.enabled', this%active)
202 if (.not. this%active)
then
205 else if (this%active)
then
207 call neko_error(
"ALE not currently supported with device backend.")
210 call neko_error(
"ALE not currently supported with OIFS.")
215 call neko_log%section(
"ALE Initialization")
217 tmp_logical = .false.
220 call this%x_ref%init(coef%dof,
"x_ref")
221 call this%y_ref%init(coef%dof,
"y_ref")
222 call this%z_ref%init(coef%dof,
"z_ref")
224 this%x_ref%x = coef%dof%x
225 this%y_ref%x = coef%dof%y
226 this%z_ref%x = coef%dof%z
229 this%user_ale_mesh_vel =>
user%ale_mesh_velocity
230 this%user_ale_base_shapes =>
user%ale_base_shapes
231 this%user_ale_rigid_kinematics =>
user%ale_rigid_kinematics
234 has_user_kin =
associated(this%user_ale_rigid_kinematics)
235 has_user_mesh =
associated(this%user_ale_mesh_vel)
238 call coef%enable_B_history()
239 call json_get(json,
'case.numerics.time_order', time_order)
243 if (
allocated(moving_zone_ids))
deallocate(moving_zone_ids)
244 allocate(moving_zone_ids(0))
255 precon_params, abstol, ksp_max_iter, res_monitor)
258 call this%bc_moving%init_from_components(coef)
259 call this%bc_fixed%init_from_components(coef)
261 if (json%valid_path(
'case.fluid.boundary_conditions'))
then
262 call json%info(
'case.fluid.boundary_conditions', n_children = n_bcs)
268 if (
allocated(bc_type))
deallocate(bc_type)
269 call json_get(bc_subdict,
'type', bc_type)
271 if (
allocated(zone_indices))
deallocate(zone_indices)
272 call json_get(bc_subdict,
'zone_indices', zone_indices)
275 if (trim(bc_type) ==
'no_slip')
then
280 do j = 1,
size(zone_indices)
284 call this%bc_moving%mark_zone(coef%msh%labeled_zones(&
287 this%has_moving_boundary = .true.
289 do j = 1,
size(zone_indices)
290 call this%bc_fixed%mark_zone(coef%msh%labeled_zones(&
297 call this%bc_moving%finalize()
298 call this%bc_fixed%finalize()
299 call this%bc_list%init()
300 call this%bc_list%append(this%bc_moving)
301 call this%bc_list%append(this%bc_fixed)
304 if (json%valid_path(
'case.fluid.ale.solver.mesh_stiffness.type'))
then
305 call json%get(
'case.fluid.ale.solver.mesh_stiffness.type', tmp_str)
306 this%config%stiffness_type = tmp_str
307 if (.not. (trim(tmp_str) ==
'built-in'))
then
308 call neko_error(
"ALE: stiffness_type must be 'built-in'")
311 if (.not.
associated(this%user_ale_base_shapes))
then
312 call neko_log%message(
'Solver Type : (' // &
313 trim(ksp_solver) //
', ' // trim(precon_type) //
')')
314 write(log_buf,
'(A,ES13.6)')
'Abs tol :', abstol
316 call neko_log%message(
'Mesh Stiffness : ' // &
317 trim(this%config%stiffness_type))
322 if (json%valid_path(
'case.fluid.ale.bodies'))
then
323 call json%info(
'case.fluid.ale.bodies', n_children = n_bodies)
324 this%config%nbodies = n_bodies
325 allocate(this%config%bodies(n_bodies))
326 allocate(this%ale_pivot(n_bodies))
327 allocate(this%body_kin(n_bodies))
328 allocate(this%base_shapes(n_bodies))
329 allocate(this%global_pivot_pos(3 * this%config%nbodies))
330 allocate(this%global_pivot_vel_lag(3 * this%config%nbodies, 3))
331 allocate(this%global_basis_pos(6 * this%config%nbodies))
332 allocate(this%ghost_handles(2, this%config%nbodies))
333 allocate(this%global_basis_vel_lag(6 * this%config%nbodies, 3))
334 allocate(this%body_rot_matrices(3, 3, this%config%nbodies))
336 this%global_pivot_pos = 0.0_rp
337 this%global_pivot_vel_lag = 0.0_rp
338 this%global_basis_pos = 0.0_rp
339 this%global_basis_vel_lag = 0.0_rp
340 this%body_rot_matrices = 0.0_rp
343 this%body_rot_matrices(1, 1, i) = 1.0_rp
344 this%body_rot_matrices(2, 2, i) = 1.0_rp
345 this%body_rot_matrices(3, 3, i) = 1.0_rp
350 this%config%bodies(i)%id = i
352 if (body_sub%valid_path(
'name'))
then
353 call json_get(body_sub,
'name', tmp_str)
354 this%config%bodies(i)%name = tmp_str
356 write(this%config%bodies(i)%name,
'(A,I0)')
'body_', i
359 if (body_sub%valid_path(
'zone_indices'))
then
360 call json_get(body_sub,
'zone_indices', zone_indices)
361 this%config%bodies(i)%zone_indices = zone_indices
364 trim(this%config%bodies(i)%name) // &
365 " must have 'zone_indices'")
369 this%config%bodies(i)%osc_amp = 0.0_rp
370 this%config%bodies(i)%osc_freq = 0.0_rp
371 if (body_sub%valid_path(
'oscillation'))
then
372 call json_get(body_sub,
'oscillation.amplitude', tmp_vec, &
374 this%config%bodies(i)%osc_amp = tmp_vec
375 call json_get(body_sub,
'oscillation.frequency', tmp_vec, &
377 this%config%bodies(i)%osc_freq = tmp_vec
381 if (body_sub%valid_path(
'rotation'))
then
383 if (.not. body_sub%valid_path(
'pivot'))
then
384 call neko_error(
"ale.bodies.pivot is missing " // &
385 "from the case file.")
388 call json_get(body_sub,
'rotation.type', tmp_str)
389 this%config%bodies(i)%rotation_type = tmp_str
391 select case (trim(tmp_str))
393 call json_get(body_sub,
'rotation.amplitude_deg', tmp_vec, &
395 this%config%bodies(i)%rot_amp_degree = tmp_vec
397 call json_get(body_sub,
'rotation.frequency', tmp_vec, &
399 this%config%bodies(i)%rot_freq = tmp_vec
403 call json_get(body_sub,
'rotation.ramp_t0', tmp_vec, &
405 this%config%bodies(i)%ramp_t0 = tmp_vec
407 call json_get(body_sub,
'rotation.ramp_omega0', tmp_vec, &
409 this%config%bodies(i)%ramp_omega0 = tmp_vec
415 if (tmp_int >= 1 .and. tmp_int <= 3)
then
416 this%config%bodies(i)%rotation_axis = tmp_int
418 call neko_error(
"ALE: rotation.axis must be (integer) " // &
419 "1 -> x, 2 -> y, or 3 -> z")
421 call json_get(body_sub,
'rotation.step_control_times', &
422 tmp_vec, expected_size = 4)
423 this%config%bodies(i)%step_control_times = tmp_vec
425 call json_get(body_sub,
'rotation.target_angle_deg', tmp_val)
426 this%config%bodies(i)%target_rot_angle_deg = tmp_val
432 if (.not.
associated(this%user_ale_mesh_vel) .and. &
433 .not.
associated(this%user_ale_rigid_kinematics))
then
434 call neko_error(
"'user' rotation is chosen, but " // &
435 "neither 'user_ale_rigid_kinematics' nor " // &
436 "'user_ale_mesh_velocity' is provided.")
440 call neko_error(
"ALE: rotation.type must be 'harmonic', " // &
441 "'ramp', 'smooth_step', or 'user'.")
446 if (body_sub%valid_path(
'pivot'))
then
449 this%config%bodies(i)%rotation_center_type = tmp_str
451 call json_get(body_sub,
'pivot.value', tmp_vec, expected_size = 3)
452 this%config%bodies(i)%rot_center = tmp_vec
455 tmp_str = this%config%bodies(i)%rotation_center_type
456 if (trim(tmp_str) /=
'relative' .and. &
457 trim(tmp_str) /=
'relative_sin')
then
458 call neko_error(
"ALE: pivot.type must be " // &
459 "'relative', or 'relative_sin'.")
464 if (body_sub%valid_path(
'stiff_geom'))
then
465 call json_get(body_sub,
'stiff_geom.type', tmp_str)
466 this%config%bodies(i)%stiff_geom%type = tmp_str
467 call json_get(body_sub,
'stiff_geom.gain', &
468 this%config%bodies(i)%stiff_geom%gain)
469 call json_get(body_sub,
'stiff_geom.decay_profile', tmp_str)
470 this%config%bodies(i)%stiff_geom%decay_profile = tmp_str
472 select case (trim(this%config%bodies(i)%stiff_geom%decay_profile))
475 'stiff_geom.cutoff_coef', &
476 this%config%bodies(i)%stiff_geom%cutoff_coef, 9.0_rp)
479 'stiff_geom.cutoff_coef', &
480 this%config%bodies(i)%stiff_geom%cutoff_coef, 3.5_rp)
482 call neko_error(
"ALE: Invalid stiff_geom.decay_profile: " // &
483 trim(this%config%bodies(i)%stiff_geom%decay_profile))
486 select case (trim(this%config%bodies(i)%stiff_geom%type))
487 case (
'cylinder',
'sphere')
488 call json_get(body_sub,
'stiff_geom.center', tmp_vec, &
490 this%config%bodies(i)%stiff_geom%center = tmp_vec
492 call json_get(body_sub,
'stiff_geom.radius', &
493 this%config%bodies(i)%stiff_geom%radius)
495 call json_get(body_sub,
'stiff_geom.stiff_dist', &
496 this%config%bodies(i)%stiff_geom%stiff_dist)
498 call neko_error(
"ALE: stiff_geom.type 'box' not yet" // &
501 call neko_error(
"ALE: Invalid stiff_geom.type: " // &
502 trim(this%config%bodies(i)%stiff_geom%type))
506 trim(this%config%bodies(i)%name) // &
507 "' must have 'stiff_geom' definition.")
513 call this%base_shapes(i)%init(coef%dof, &
514 "phi_" // trim(this%config%bodies(i)%name))
520 this%ghost_handles(1, i) = this%request_tracker( &
521 this%config%bodies(i)%rot_center + [1.0_rp, 0.0_rp, 0.0_rp], &
522 this%config%bodies(i)%id)
524 this%ghost_handles(2, i) = this%request_tracker( &
525 this%config%bodies(i)%rot_center + [0.0_rp, 1.0_rp, 0.0_rp], &
526 this%config%bodies(i)%id)
528 call neko_log%message(
'Registered Body : ' // &
529 trim(this%config%bodies(i)%name))
533 if (.not.
associated(this%user_ale_base_shapes))
then
534 write(log_buf,
'(A,A)')
' Stiff Type : ', &
535 trim(this%config%bodies(i)%stiff_geom%type)
537 write(log_buf,
'(A,ES18.11,A,A,A,ES10.4)')
' Gain : ', &
538 this%config%bodies(i)%stiff_geom%gain,
' | Profile: ', &
539 trim(this%config%bodies(i)%stiff_geom%decay_profile), &
540 ' | Cutoff: ', this%config%bodies(i)%stiff_geom%cutoff_coef
542 select case (trim(this%config%bodies(i)%stiff_geom%type))
543 case (
'cylinder',
'sphere')
544 write(log_buf,
'(A,3(ES23.15,1X))')
' Center :', &
545 this%config%bodies(i)%stiff_geom%center
547 write(log_buf,
'(A,ES23.15)')
' Radius :', &
548 this%config%bodies(i)%stiff_geom%radius
551 write(log_buf,
'(A,ES23.15)')
' Stiff Dist:', &
552 this%config%bodies(i)%stiff_geom%stiff_dist
559 has_builtin_osc = any(abs(this%config%bodies(i)%osc_amp) > 0.0_rp)
561 if (has_builtin_osc)
then
562 if (has_user_kin .or. has_user_mesh)
then
563 call neko_log%message(
' Oscillation : ' // &
564 'X(t) = Amp*sin(2*pi*Freq*t) + User')
565 write(log_buf,
'(A,3(ES18.11,1X))')
' Amp :', &
566 this%config%bodies(i)%osc_amp
568 write(log_buf,
'(A,3(ES18.11,1X))')
' Freq :', &
569 this%config%bodies(i)%osc_freq
572 call neko_log%message(
' Oscillation : ' // &
573 'X(t) = Amp*sin(2*pi*Freq*t)')
574 write(log_buf,
'(A,3(ES18.11,1X))')
' Amp :', &
575 this%config%bodies(i)%osc_amp
577 write(log_buf,
'(A,3(ES18.11,1X))')
' Freq :', &
578 this%config%bodies(i)%osc_freq
582 if (has_user_kin .or. has_user_mesh)
then
583 call neko_log%message(
' Oscillation : User-defined')
585 call neko_log%message(
' Oscillation : None')
591 has_builtin_rot = (trim(this%config%bodies(i)%rotation_type) &
594 if (trim(this%config%bodies(i)%rotation_type) ==
'user')
then
596 call neko_log%message(
' Rotation Type: User-defined')
598 elseif (has_builtin_rot)
then
601 is_rot_active = .false.
602 select case (trim(this%config%bodies(i)%rotation_type))
604 is_rot_active = any(abs(this%config%bodies(i)%rot_amp_degree) &
607 is_rot_active = any(abs(this%config%bodies(i)%ramp_omega0) &
611 (abs(this%config%bodies(i)%target_rot_angle_deg) > 0.0_rp)
614 if (is_rot_active)
then
616 if (trim(this%config%bodies(i)%rotation_type) &
618 if (has_user_kin .or. has_user_mesh)
then
619 call neko_log%message(
' Rotation : ' // &
620 'Theta(t) = Amp*sin(2*pi*Freq*t) + User')
622 call neko_log%message(
' Rotation : ' // &
623 'Theta(t) = Amp*sin(2*pi*Freq*t)')
625 write(log_buf,
'(A,3(ES18.11,1X))')
' Amp (deg) :', &
626 this%config%bodies(i)%rot_amp_degree
628 write(log_buf,
'(A,3(ES18.11,1X))')
' Freq :', &
629 this%config%bodies(i)%rot_freq
633 elseif (trim(this%config%bodies(i)%rotation_type) &
635 if (has_user_kin .or. has_user_mesh)
then
636 call neko_log%message(
' Rotation : ' // &
637 'Omega(t) = Omega0*(1 - exp(-4.6*t/t0)) + User')
639 call neko_log%message(
' Rotation : ' // &
640 'Omega(t) = Omega0*(1 - exp(-4.6*t/t0))')
642 write(log_buf,
'(A,3(ES18.11,1X))')
' Omega0 :', &
643 this%config%bodies(i)%ramp_omega0
645 write(log_buf,
'(A,3(ES18.11,1X))')
' t0 :', &
646 this%config%bodies(i)%ramp_t0
650 elseif (trim(this%config%bodies(i)%rotation_type) &
651 ==
'smooth_step')
then
652 if (has_user_kin .or. has_user_mesh)
then
653 call neko_log%message(
' Rotation : ' // &
654 'Smooth Step Control + User')
656 call neko_log%message(
' Rotation : ' // &
657 'Smooth Step Control')
659 write(log_buf,
'(A,I10)')
' Rotation Axis :', &
660 this%config%bodies(i)%rotation_axis
662 write(log_buf,
'(A,ES18.11)')
' Target Rot ' // &
664 this%config%bodies(i)%target_rot_angle_deg
666 write(log_buf,
'(A,4(ES18.11,1X))') &
667 ' Control Times [t0, t1, t2, t3] :', &
668 this%config%bodies(i)%step_control_times
672 if (has_user_kin .or. has_user_mesh)
then
673 call neko_log%message(
' Rotation Type: User-defined')
675 call neko_log%message(
' Rotation Type: None')
683 call neko_log%message(
' Pivot Type : ' // &
684 trim(this%config%bodies(i)%rotation_center_type))
686 write(log_buf,
'(A,3(ES18.11,1X))')
' Init Pivot:', &
687 this%config%bodies(i)%rot_center
693 call neko_error(
"ALE: No 'ale bodies' found in case file!")
696 if (this%config%nbodies > 1)
then
697 call this%phi_total%init(coef%dof,
"phi_total")
702 do i = 1, n_moving_zones
703 z = moving_zone_ids(i)
706 do while ((.not. found_zone) .and. (j <= this%config%nbodies))
707 if (any(this%config%bodies(j)%zone_indices == z))
then
712 if (.not. found_zone)
then
713 write(log_buf_l,
'(A,I0,A)') &
714 "ALE: zone index ", z, &
715 " has BC no_slip with moving: true, " // &
716 "but it is not registered in ALE bodies."
722 do j = 1, this%config%nbodies
723 if (
allocated(this%config%bodies(j)%zone_indices))
then
724 do i = 1,
size(this%config%bodies(j)%zone_indices)
725 z = this%config%bodies(j)%zone_indices(i)
727 if (n_moving_zones > 0)
then
728 if (any(moving_zone_ids(1:n_moving_zones) == z))
then
732 if (.not. found_zone)
then
733 write(log_buf_l,
'(A,I0,A,A)') &
734 "ALE: zone index ", z, &
735 " is registered in ALE bodies, ", &
736 "but the BC is not no_slip with moving: true."
744 do j = 1, this%config%nbodies
745 if (
allocated(this%config%bodies(j)%zone_indices))
then
746 do i = 1,
size(this%config%bodies(j)%zone_indices)
747 z = this%config%bodies(j)%zone_indices(i)
749 do k = j + 1, this%config%nbodies
750 if (
allocated(this%config%bodies(k)%zone_indices))
then
751 if (any(this%config%bodies(k)%zone_indices == z))
then
752 write(log_buf_l,
'(A,I0,A,A,A,A,A)') &
753 "ALE: zone index ", z, &
754 " is assigned to multiple bodies ('", &
755 trim(this%config%bodies(j)%name),
"' and '", &
756 trim(this%config%bodies(k)%name),
"')."
767 call this%solve_base_mesh_displacement(coef, abstol, ksp_solver, &
768 ksp_max_iter, precon_type, precon_params, res_monitor)
772 if (.not. json%valid_path(
'case.restart_file'))
then
776 call this%update_mesh_velocity(coef, t_init)
779 call this%wm_x_lag%init(this%wm_x, time_order)
780 call this%wm_y_lag%init(this%wm_y, time_order)
781 call this%wm_z_lag%init(this%wm_z, time_order)
788 if (
allocated(moving_zone_ids))
deallocate(moving_zone_ids)
789 if (
allocated(bc_type))
deallocate(bc_type)
790 if (
allocated(zone_indices))
deallocate(zone_indices)
791 if (
allocated(ksp_solver))
deallocate(ksp_solver)
792 if (
allocated(precon_type))
deallocate(precon_type)
795 call this%mesh_preview(coef, json)
805 ksp_max_iter, precon_type, precon_params, res_monitor)
807 class(
ax_t),
allocatable :: Ax
808 class(
ksp_t),
allocatable :: ksp
809 class(
pc_t),
allocatable :: pc
810 type(
coef_t),
intent(inout) :: coef
811 real(kind=
rp),
intent(in) :: abstol
812 logical,
intent(in) :: res_monitor
813 character(len=*),
intent(in) :: ksp_solver, precon_type
814 integer,
intent(in) :: ksp_max_iter
815 type(json_file),
intent(inout) :: precon_params
820 real(kind=
rp) :: sample_start_time, sample_end_time
821 real(kind=
rp) :: sample_time
822 character(len=LOG_SIZE) :: log_buf
823 integer :: n, i, m, k, ierr, body_idx, z_idx
825 real(kind=
rp),
allocatable :: h1_restore(:, :, :, :)
826 real(kind=
rp),
allocatable :: h2_restore(:, :, :, :)
833 if (.not. this%active)
return
834 if (.not. this%has_moving_boundary)
return
835 if (this%config%nbodies == 0)
return
838 call neko_log%message(
"Starting base mesh motion solve ...")
841 call ax_helm_factory(ax, full_formulation = .false.)
842 call krylov_solver_factory(ksp, n, ksp_solver, &
843 ksp_max_iter, abstol, monitor = res_monitor)
845 coef%gs_h, this%bc_list, precon_type, precon_params)
851 call rhs_field%init(coef%dof)
852 call corr_field%init(coef%dof)
856 if (
associated(this%user_ale_base_shapes))
then
857 call neko_log%message(
" Using user-defined base shapes " // &
858 "(skipping Laplace solve)")
861 call this%user_ale_base_shapes(this%base_shapes)
864 if (this%config%nbodies > 1)
then
866 do body_idx = 1, this%config%nbodies
867 call field_add2(this%phi_total, this%base_shapes(body_idx), n)
872 if (this%config%if_output_phi)
then
874 do body_idx = 1, this%config%nbodies
875 call phi_file%init(
'phi_' // &
876 trim(this%config%bodies(body_idx)%name) //
'.fld')
877 call phi_file%write(this%base_shapes(body_idx))
880 trim(this%config%bodies(body_idx)%name) //
'.fld saved.')
884 if (this%config%nbodies > 1)
then
885 call neko_log%message(
" phi_total.fld saved.")
886 call phi_file%init(
'phi_total.fld')
887 call phi_file%write(this%phi_total)
898 if (this%config%if_output_stiffness)
then
899 rhs_field%x = coef%h1
900 call phi_file%init(
'stiffness.fld')
901 call phi_file%write(rhs_field)
907 do body_idx = 1, this%config%nbodies
909 sample_start_time = mpi_wtime()
910 call neko_log%message(
" Solving laplace for body: " // &
911 trim(this%config%bodies(body_idx)%name))
913 call bc_active_body%init_from_components(coef)
914 call bc_inactive_body%init_from_components(coef)
917 do j = 1,
size(this%config%bodies(body_idx)%zone_indices)
918 z_idx = this%config%bodies(body_idx)%zone_indices(j)
919 call bc_active_body%mark_zone(coef%msh%labeled_zones(z_idx))
922 do i = 1, this%config%nbodies
923 if (i /= body_idx)
then
924 do j = 1,
size(this%config%bodies(i)%zone_indices)
925 z_idx = this%config%bodies(i)%zone_indices(j)
926 call bc_inactive_body%mark_zone(&
927 coef%msh%labeled_zones(z_idx))
932 call bc_active_body%finalize()
933 call bc_inactive_body%finalize()
937 call bcloc%append(this%bc_fixed)
938 call bcloc%append(bc_active_body)
939 call bcloc%append(bc_inactive_body)
942 call bcloc_zeros_only%init()
943 call bcloc_zeros_only%append(this%bc_fixed)
944 call bcloc_zeros_only%append(bc_inactive_body)
954 m = bc_active_body%msk(0)
956 k = bc_active_body%msk(i)
957 this%base_shapes(body_idx)%x(k, 1, 1, 1) = 1.0_rp
963 call bcloc_zeros_only%apply_scalar(this%base_shapes(body_idx)%x, n)
967 call ax%compute(rhs_field%x, this%base_shapes(body_idx)%x, &
968 coef, coef%msh, coef%Xh)
973 call bcloc%apply_scalar(rhs_field%x, n)
974 call coef%gs_h%op(rhs_field, gs_op_add)
979 monitor(1) = ksp%solve(ax, corr_field, &
980 rhs_field%x, n, coef, bcloc, coef%gs_h)
983 call field_add2(this%base_shapes(body_idx), corr_field, n)
987 if (this%config%nbodies > 1)
then
988 call field_add2(this%phi_total, this%base_shapes(body_idx), n)
992 sample_end_time = mpi_wtime()
993 sample_time = sample_end_time - sample_start_time
994 write(log_buf,
'(A, A, A, ES11.4, A)')
" Laplace solve for '", &
995 trim(this%config%bodies(body_idx)%name),
"' took ", &
1000 call bc_active_body%free()
1001 call bc_inactive_body%free()
1003 call bcloc_zeros_only%free()
1005 if (this%config%if_output_phi)
then
1006 call phi_file%init(
'phi_' // &
1007 trim(this%config%bodies(body_idx)%name) //
'.fld')
1008 call phi_file%write(this%base_shapes(body_idx))
1009 call phi_file%free()
1011 trim(this%config%bodies(body_idx)%name) //
'.fld saved.')
1015 if (this%config%if_output_phi .and. (this%config%nbodies > 1))
then
1016 call neko_log%message(
" phi_total.fld saved.")
1017 call phi_file%init(
'phi_total.fld')
1018 call phi_file%write(this%phi_total)
1019 call phi_file%free()
1023 call rhs_field%free()
1024 call corr_field%free()
1027 coef%h1 = h1_restore
1028 coef%h2 = h2_restore
1030 if (
allocated(h1_restore))
deallocate(h1_restore)
1031 if (
allocated(h2_restore))
deallocate(h2_restore)
1032 if (
allocated(ax))
deallocate(ax)
1033 if (
allocated(ksp))
then
1037 if (
allocated(pc))
then
1038 call precon_destroy(pc)