Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
ale_manager.f90
Go to the documentation of this file.
1! Copyright (c) 2025-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 num_types, only : rp, dp
36 use json_module, only : json_file
38 use field, only : field_t
39 use coefs, only : coef_t
40 use ax_product, only : ax_t, ax_helm_factory
41 use krylov, only : ksp_t, ksp_monitor_t, krylov_solver_factory
42 use precon, only : pc_t, precon_factory, precon_destroy
43 use bc_list, only : bc_list_t
45 use gather_scatter, only : gs_t, gs_op_add
46 use dofmap, only : dofmap_t
47 use jacobi, only : jacobi_t
48 use hsmg, only : hsmg_t
49 use phmg, only : phmg_t
51 use sx_jacobi, only : sx_jacobi_t
53 use file, only : file_t
54 use logger, only : neko_log, log_size
55 use advection, only : advection_t
64 use utils, only : neko_error
66 use mpi_f08, only : mpi_wtime, mpi_barrier
67 use comm, only : neko_comm
68 use registry, only : neko_registry
70 use time_state, only : time_state_t
71 use fld_file, only : fld_file_t
74 use math, only : glmin, pi
75 use field_math, only : field_rzero, field_add2, &
77
78 implicit none
79 private
80
81 public :: compute_stiffness_ale
83 public :: update_ale_mesh
84 public :: log_rot_angles
85 public :: log_pivot
86
87 type, public :: ale_manager_t
88 ! Default
89 logical :: active = .false.
90 logical :: has_moving_boundary = .false.
91
93 type(zero_dirichlet_t) :: bc_moving
94 type(zero_dirichlet_t) :: bc_fixed
95
96 type(ale_config_t) :: config
97
99 type(field_t), pointer :: wm_x => null()
100 type(field_t), pointer :: wm_y => null()
101 type(field_t), pointer :: wm_z => null()
102
104 type(field_series_t) :: wm_x_lag
105 type(field_series_t) :: wm_y_lag
106 type(field_series_t) :: wm_z_lag
107
109 type(field_t) :: x_ref, y_ref, z_ref
110
112 type(pivot_state_t), allocatable :: ale_pivot(:)
113 type(body_kinematics_t), allocatable :: body_kin(:)
114
117 type(field_t), allocatable :: base_shapes(:)
118
120 type(field_t) :: phi_total
121
122 real(kind=rp), pointer :: global_pivot_pos(:) => null()
123 real(kind=rp), pointer :: global_pivot_vel_lag(:, :) => null()
124
125 ! Basis Vectors for orientation
126 real(kind=rp), pointer :: global_basis_pos(:) => null()
127 ! Store history for the ghost trackers
128 real(kind=rp), pointer :: global_basis_vel_lag(:, :) => null()
129 ! Private handles to the ghost trackers (2 per body)
130 integer, allocatable :: ghost_handles(:,:)
131 ! Rotation matrices
132 real(kind=rp), allocatable :: body_rot_matrices(:,:,:)
133
134 type(point_tracker_t), allocatable :: trackers(:)
135 integer :: n_trackers = 0
136
137 procedure(user_ale_mesh_velocity_intf), nopass, pointer :: &
138 user_ale_mesh_vel => null()
139 procedure(user_ale_base_shapes_intf), nopass, pointer :: &
140 user_ale_base_shapes => null()
141 procedure(user_ale_rigid_kinematics_intf), nopass, pointer :: &
142 user_ale_rigid_kinematics => null()
143
144 contains
145 procedure, pass(this) :: init => ale_manager_init
146 procedure, pass(this) :: free => ale_manager_free
147 procedure, pass(this) :: mesh_preview
148 procedure, pass(this) :: solve_base_mesh_displacement
149 procedure, pass(this) :: advance_mesh
150 procedure, pass(this) :: update_mesh_velocity
151 procedure, pass(this) :: set_pivot_restart
152 procedure, pass(this) :: set_coef_restart
153 procedure, pass(this) :: request_tracker
154 procedure, pass(this) :: get_tracker_pos
155 procedure, pass(this) :: compute_rotation_matrix
156 procedure, pass(this) :: prep_checkpoint => set_pivot_basis_for_checkpoint
157 procedure, pass(this) :: ghost_tracker_coord_step
158 procedure, pass(this) :: log_rot_angles
159 procedure, pass(this) :: log_pivot
160 end type ale_manager_t
161
162 type(ale_manager_t), public, pointer :: neko_ale => null()
163
164contains
165
168 subroutine ale_manager_init(this, coef, json, user)
169 class(ale_manager_t), intent(inout), target :: this
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
174 type(time_state_t) :: t_init
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
191 logical :: moving_
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
196
197 if (json%valid_path('case.fluid.ale')) then
198 call json_get(json, 'case.fluid.ale.enabled', this%active)
199 end if
200 call json_get_or_default(json, 'case.numerics.oifs', oifs, .false.)
201
202 if (.not. this%active) then
203 neko_ale => null()
204 return
205 else if (this%active) then
206 if (neko_bcknd_device .eq. 1) then
207 call neko_error("ALE not currently supported with device backend.")
208 end if
209 if (oifs) then
210 call neko_error("ALE not currently supported with OIFS.")
211 end if
212 neko_ale => this
213 end if
214
215 call neko_log%section("ALE Initialization")
216
217 tmp_logical = .false.
218 n = coef%dof%size()
219
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")
223
224 this%x_ref%x = coef%dof%x
225 this%y_ref%x = coef%dof%y
226 this%z_ref%x = coef%dof%z
227
228 ! Set user function pointers.
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
232
233 ! Check user association states
234 has_user_kin = associated(this%user_ale_rigid_kinematics)
235 has_user_mesh = associated(this%user_ale_mesh_vel)
236
237 ! Enable B history (Blag, Blaglag)
238 call coef%enable_B_history()
239 call json_get(json, 'case.numerics.time_order', time_order)
240
241 ! Stuff for zone_id checks
242 n_moving_zones = 0
243 if (allocated(moving_zone_ids)) deallocate(moving_zone_ids)
244 allocate(moving_zone_ids(0))
245
246 ! Register mesh velocity fields
247 call neko_registry%add_field(coef%dof, 'wm_x')
248 call neko_registry%add_field(coef%dof, 'wm_y')
249 call neko_registry%add_field(coef%dof, 'wm_z')
250 this%wm_x => neko_registry%get_field('wm_x')
251 this%wm_y => neko_registry%get_field('wm_y')
252 this%wm_z => neko_registry%get_field('wm_z')
253
254 call get_ale_solver_params_json(this, json, ksp_solver, precon_type, &
255 precon_params, abstol, ksp_max_iter, res_monitor)
256
257 ! Mark BCs
258 call this%bc_moving%init_from_components(coef)
259 call this%bc_fixed%init_from_components(coef)
260
261 if (json%valid_path('case.fluid.boundary_conditions')) then
262 call json%info('case.fluid.boundary_conditions', n_children = n_bcs)
263
264 do i = 1, n_bcs
265 call json_extract_item(json, 'case.fluid.boundary_conditions', &
266 i, bc_subdict)
267
268 if (allocated(bc_type)) deallocate(bc_type)
269 call json_get(bc_subdict, 'type', bc_type)
270
271 if (allocated(zone_indices)) deallocate(zone_indices)
272 call json_get(bc_subdict, 'zone_indices', zone_indices)
273
274 moving_ = .false.
275 if (trim(bc_type) == 'no_slip') then
276 call json_get_or_default(bc_subdict, 'moving', moving_, .false.)
277 end if
278
279 if (moving_) then
280 do j = 1, size(zone_indices)
281 ! we append unique moving zone ids for future checks
282 call append_unique_int(moving_zone_ids, n_moving_zones, &
283 zone_indices(j))
284 call this%bc_moving%mark_zone(coef%msh%labeled_zones(&
285 zone_indices(j)))
286 end do
287 this%has_moving_boundary = .true.
288 else
289 do j = 1, size(zone_indices)
290 call this%bc_fixed%mark_zone(coef%msh%labeled_zones(&
291 zone_indices(j)))
292 end do
293 end if
294 end do
295 end if
296
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)
302
303 ! Mesh Stiffness
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'")
309 end if
310 end if
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
315 call neko_log%message(log_buf)
316 call neko_log%message('Mesh Stiffness : ' // &
317 trim(this%config%stiffness_type))
318 end if
319 call neko_log%message(' ')
320
321 ! Bodies
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))
335
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
341
342 do i = 1, n_bodies
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
346 end do
347
348 do i = 1, n_bodies
349 call json_extract_item(json, 'case.fluid.ale.bodies', i, body_sub)
350 this%config%bodies(i)%id = i
351
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
355 else
356 write(this%config%bodies(i)%name, '(A,I0)') 'body_', i
357 endif
358
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
362 else
363 call neko_error("ALE: body " // &
364 trim(this%config%bodies(i)%name) // &
365 " must have 'zone_indices'")
366 endif
367
368 ! Oscillation
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, &
373 expected_size = 3)
374 this%config%bodies(i)%osc_amp = tmp_vec
375 call json_get(body_sub, 'oscillation.frequency', tmp_vec, &
376 expected_size = 3)
377 this%config%bodies(i)%osc_freq = tmp_vec
378 end if
379
380 ! Rotation
381 if (body_sub%valid_path('rotation')) then
382 ! Check if pivot exists.
383 if (.not. body_sub%valid_path('pivot')) then
384 call neko_error("ale.bodies.pivot is missing " // &
385 "from the case file.")
386 end if
387
388 call json_get(body_sub, 'rotation.type', tmp_str)
389 this%config%bodies(i)%rotation_type = tmp_str
390
391 select case (trim(tmp_str))
392 case ('harmonic')
393 call json_get(body_sub, 'rotation.amplitude_deg', tmp_vec, &
394 expected_size = 3)
395 this%config%bodies(i)%rot_amp_degree = tmp_vec
396
397 call json_get(body_sub, 'rotation.frequency', tmp_vec, &
398 expected_size = 3)
399 this%config%bodies(i)%rot_freq = tmp_vec
400
401
402 case ('ramp')
403 call json_get(body_sub, 'rotation.ramp_t0', tmp_vec, &
404 expected_size = 3)
405 this%config%bodies(i)%ramp_t0 = tmp_vec
406
407 call json_get(body_sub, 'rotation.ramp_omega0', tmp_vec, &
408 expected_size = 3)
409 this%config%bodies(i)%ramp_omega0 = tmp_vec
410
411
412 case ('smooth_step')
413 call json_get_or_default(body_sub, 'rotation.axis', &
414 tmp_int, 3)
415 if (tmp_int >= 1 .and. tmp_int <= 3) then
416 this%config%bodies(i)%rotation_axis = tmp_int
417 else
418 call neko_error("ALE: rotation.axis must be (integer) " // &
419 "1 -> x, 2 -> y, or 3 -> z")
420 end if
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
424
425 call json_get(body_sub, 'rotation.target_angle_deg', tmp_val)
426 this%config%bodies(i)%target_rot_angle_deg = tmp_val
427
428 ! This mode can be neglected. It exists just in case!
429 ! user motion can be combined with any other type from above
430 ! using user interfaces.
431 case ('user')
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.")
437 end if
438
439 case default
440 call neko_error("ALE: rotation.type must be 'harmonic', " // &
441 "'ramp', 'smooth_step', or 'user'.")
442 end select
443 end if
444
445 ! Rotation Center
446 if (body_sub%valid_path('pivot')) then
447 call json_get_or_default(body_sub, 'pivot.type', tmp_str, &
448 'relative')
449 this%config%bodies(i)%rotation_center_type = tmp_str
450
451 call json_get(body_sub, 'pivot.value', tmp_vec, expected_size = 3)
452 this%config%bodies(i)%rot_center = tmp_vec
453
454
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'.")
460 end if
461 end if
462
463 ! Stiff Geom
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
471
472 select case (trim(this%config%bodies(i)%stiff_geom%decay_profile))
473 case ('gaussian')
474 call json_get_or_default(body_sub, &
475 'stiff_geom.cutoff_coef', &
476 this%config%bodies(i)%stiff_geom%cutoff_coef, 9.0_rp)
477 case ('tanh')
478 call json_get_or_default(body_sub, &
479 'stiff_geom.cutoff_coef', &
480 this%config%bodies(i)%stiff_geom%cutoff_coef, 3.5_rp)
481 case default
482 call neko_error("ALE: Invalid stiff_geom.decay_profile: " // &
483 trim(this%config%bodies(i)%stiff_geom%decay_profile))
484 end select
485
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, &
489 expected_size = 3)
490 this%config%bodies(i)%stiff_geom%center = tmp_vec
491
492 call json_get(body_sub, 'stiff_geom.radius', &
493 this%config%bodies(i)%stiff_geom%radius)
494 case ('cheap_dist')
495 call json_get(body_sub, 'stiff_geom.stiff_dist', &
496 this%config%bodies(i)%stiff_geom%stiff_dist)
497 case ('box')
498 call neko_error("ALE: stiff_geom.type 'box' not yet" // &
499 " implemented.")
500 case default
501 call neko_error("ALE: Invalid stiff_geom.type: " // &
502 trim(this%config%bodies(i)%stiff_geom%type))
503 end select
504 else
505 call neko_error("ALE: Body '" // &
506 trim(this%config%bodies(i)%name) // &
507 "' must have 'stiff_geom' definition.")
508 end if
509
510 ! Initialize the pivots.
511 call init_pivot_state(this%ale_pivot(i), this%config%bodies(i))
512
513 call this%base_shapes(i)%init(coef%dof, &
514 "phi_" // trim(this%config%bodies(i)%name))
515 call field_rzero(this%base_shapes(i))
516
517 ! Create Ghost Trackers for numerically forming the rotation matrix
518 ! of each body.
519 ! Basis X (Pivot + 1.0 in X)
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)
523 ! Basis Y (Pivot + 1.0 in Y)
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)
527
528 call neko_log%message('Registered Body : ' // &
529 trim(this%config%bodies(i)%name))
530
531 ! Logging Stiff Body
532 call neko_log%message(' ')
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)
536 call neko_log%message(log_buf)
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
541 call neko_log%message(log_buf)
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
546 call neko_log%message(log_buf)
547 write(log_buf, '(A,ES23.15)') ' Radius :', &
548 this%config%bodies(i)%stiff_geom%radius
549 call neko_log%message(log_buf)
550 case ('cheap_dist')
551 write(log_buf, '(A,ES23.15)') ' Stiff Dist:', &
552 this%config%bodies(i)%stiff_geom%stiff_dist
553 call neko_log%message(log_buf)
554 end select
555 end if
556 call neko_log%message(' ')
557
558 ! Logging Oscillation
559 has_builtin_osc = any(abs(this%config%bodies(i)%osc_amp) > 0.0_rp)
560
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
567 call neko_log%message(log_buf)
568 write(log_buf, '(A,3(ES18.11,1X))') ' Freq :', &
569 this%config%bodies(i)%osc_freq
570 call neko_log%message(log_buf)
571 else
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
576 call neko_log%message(log_buf)
577 write(log_buf, '(A,3(ES18.11,1X))') ' Freq :', &
578 this%config%bodies(i)%osc_freq
579 call neko_log%message(log_buf)
580 end if
581 else
582 if (has_user_kin .or. has_user_mesh) then
583 call neko_log%message(' Oscillation : User-defined')
584 else
585 call neko_log%message(' Oscillation : None')
586 end if
587 end if
588 call neko_log%message(' ')
589
590 ! Logging Rotation
591 has_builtin_rot = (trim(this%config%bodies(i)%rotation_type) &
592 /= 'user')
593
594 if (trim(this%config%bodies(i)%rotation_type) == 'user') then
595
596 call neko_log%message(' Rotation Type: User-defined')
597
598 elseif (has_builtin_rot) then
599
600 ! Check parameters active
601 is_rot_active = .false.
602 select case (trim(this%config%bodies(i)%rotation_type))
603 case ('harmonic')
604 is_rot_active = any(abs(this%config%bodies(i)%rot_amp_degree) &
605 > 0.0_rp)
606 case ('ramp')
607 is_rot_active = any(abs(this%config%bodies(i)%ramp_omega0) &
608 > 0.0_rp)
609 case ('smooth_step')
610 is_rot_active = &
611 (abs(this%config%bodies(i)%target_rot_angle_deg) > 0.0_rp)
612 end select
613
614 if (is_rot_active) then
615 ! Harmonic
616 if (trim(this%config%bodies(i)%rotation_type) &
617 == 'harmonic') then
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')
621 else
622 call neko_log%message(' Rotation : ' // &
623 'Theta(t) = Amp*sin(2*pi*Freq*t)')
624 end if
625 write(log_buf, '(A,3(ES18.11,1X))') ' Amp (deg) :', &
626 this%config%bodies(i)%rot_amp_degree
627 call neko_log%message(log_buf)
628 write(log_buf, '(A,3(ES18.11,1X))') ' Freq :', &
629 this%config%bodies(i)%rot_freq
630 call neko_log%message(log_buf)
631
632 ! Ramp
633 elseif (trim(this%config%bodies(i)%rotation_type) &
634 == 'ramp') then
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')
638 else
639 call neko_log%message(' Rotation : ' // &
640 'Omega(t) = Omega0*(1 - exp(-4.6*t/t0))')
641 end if
642 write(log_buf, '(A,3(ES18.11,1X))') ' Omega0 :', &
643 this%config%bodies(i)%ramp_omega0
644 call neko_log%message(log_buf)
645 write(log_buf, '(A,3(ES18.11,1X))') ' t0 :', &
646 this%config%bodies(i)%ramp_t0
647 call neko_log%message(log_buf)
648
649 ! Smooth Step
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')
655 else
656 call neko_log%message(' Rotation : ' // &
657 'Smooth Step Control')
658 end if
659 write(log_buf, '(A,I10)') ' Rotation Axis :', &
660 this%config%bodies(i)%rotation_axis
661 call neko_log%message(log_buf)
662 write(log_buf, '(A,ES18.11)') ' Target Rot ' // &
663 'Angle (deg) :', &
664 this%config%bodies(i)%target_rot_angle_deg
665 call neko_log%message(log_buf)
666 write(log_buf, '(A,4(ES18.11,1X))') &
667 ' Control Times [t0, t1, t2, t3] :', &
668 this%config%bodies(i)%step_control_times
669 call neko_log%message(log_buf)
670 end if
671 else
672 if (has_user_kin .or. has_user_mesh) then
673 call neko_log%message(' Rotation Type: User-defined')
674 else
675 call neko_log%message(' Rotation Type: None')
676 end if
677 end if
678
679 end if
680
681 ! Logging Pivot
682 call neko_log%message(' ')
683 call neko_log%message(' Pivot Type : ' // &
684 trim(this%config%bodies(i)%rotation_center_type))
685
686 write(log_buf, '(A,3(ES18.11,1X))') ' Init Pivot:', &
687 this%config%bodies(i)%rot_center
688 call neko_log%message(log_buf)
689 call neko_log%message(' ')
690
691 end do
692 else
693 call neko_error("ALE: No 'ale bodies' found in case file!")
694 end if
695
696 if (this%config%nbodies > 1) then
697 call this%phi_total%init(coef%dof, "phi_total")
698 call field_rzero(this%phi_total)
699 end if
700
701 ! Check to be sure moving no_slip ids belong to an ALE body
702 do i = 1, n_moving_zones
703 z = moving_zone_ids(i)
704 found_zone = .false.
705 j = 1
706 do while ((.not. found_zone) .and. (j <= this%config%nbodies))
707 if (any(this%config%bodies(j)%zone_indices == z)) then
708 found_zone = .true.
709 end if
710 j = j + 1
711 end do
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."
717 call neko_error(trim(log_buf_l))
718 end if
719 end do
720
721 ! Any id registered in ALE bodies must have no_slip with moving: true in BCs.
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)
726 found_zone = .false.
727 if (n_moving_zones > 0) then
728 if (any(moving_zone_ids(1:n_moving_zones) == z)) then
729 found_zone = .true.
730 end if
731 end if
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."
737 call neko_error(trim(log_buf_l))
738 end if
739 end do
740 end if
741 end do
742
743 ! Check no zone ID is assigned to more than one ALE body.
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)
748
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), "')."
757 call neko_error(trim(log_buf_l))
758 end if
759 end if
760 end do
761
762 end do
763 end if
764 end do
765
766 ! Find the smooth blending function for mesh displacement.
767 call this%solve_base_mesh_displacement(coef, abstol, ksp_solver, &
768 ksp_max_iter, precon_type, precon_params, res_monitor)
769
770 ! If we are restarting, we skip this. It will be handled
771 ! properly by chkp file.
772 if (.not. json%valid_path('case.restart_file')) then
773 t_init%t = 0.0_rp
774 t_init%tstep = 0
775 t_init%dt = 0.0_rp
776 call this%update_mesh_velocity(coef, t_init)
777 end if
778
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)
782 do i = 1, time_order
783 call field_rzero(this%wm_x_lag%lf(i))
784 call field_rzero(this%wm_y_lag%lf(i))
785 call field_rzero(this%wm_z_lag%lf(i))
786 end do
787
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)
793
794 ! Performing mesh_preview.
795 call this%mesh_preview(coef, json)
796
797 call neko_log%end_section()
798 end subroutine ale_manager_init
799
804 subroutine solve_base_mesh_displacement(this, coef, abstol, ksp_solver, &
805 ksp_max_iter, precon_type, precon_params, res_monitor)
806 class(ale_manager_t), intent(inout) :: this
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
816 type(file_t) :: phi_file
817 type(field_t) :: rhs_field
818 type(field_t) :: corr_field
819 type(ksp_monitor_t) :: monitor(1)
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
824 integer :: j
825 real(kind=rp), allocatable :: h1_restore(:, :, :, :)
826 real(kind=rp), allocatable :: h2_restore(:, :, :, :)
827 type(zero_dirichlet_t) :: bc_active_body
828 type(zero_dirichlet_t) :: bc_inactive_body
829 type(bc_list_t) :: bcloc
830 type(bc_list_t) :: bcloc_zeros_only
831
832
833 if (.not. this%active) return
834 if (.not. this%has_moving_boundary) return
835 if (this%config%nbodies == 0) return
836
837 call neko_log%message(" ")
838 call neko_log%message("Starting base mesh motion solve ...")
839 n = coef%dof%size()
840
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)
844 call ale_precon_factory(pc, ksp, coef, coef%dof, &
845 coef%gs_h, this%bc_list, precon_type, precon_params)
846
847 ! Save original h1/h2
848 h1_restore = coef%h1
849 h2_restore = coef%h2
850
851 call rhs_field%init(coef%dof)
852 call corr_field%init(coef%dof)
853
854
855 ! User Defined Base Shapes (Skip Solver).
856 if (associated(this%user_ale_base_shapes)) then
857 call neko_log%message(" Using user-defined base shapes " // &
858 "(skipping Laplace solve)")
859
860 ! Call User Hook (Populates this%base_shapes)
861 call this%user_ale_base_shapes(this%base_shapes)
862
863 ! Compute phi_total (Sum of all user shapes)
864 if (this%config%nbodies > 1) then
865 call field_rzero(this%phi_total)
866 do body_idx = 1, this%config%nbodies
867 call field_add2(this%phi_total, this%base_shapes(body_idx), n)
868 end do
869 end if
870
871 ! Output Shapes
872 if (this%config%if_output_phi) then
873 ! Individual Bodies
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))
878 call phi_file%free()
879 call neko_log%message(' phi_' // &
880 trim(this%config%bodies(body_idx)%name) // '.fld saved.')
881 end do
882
883 ! Total
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)
888 call phi_file%free()
889 end if
890 end if
891 else
892 ! Standard Laplace Solve (Requires Stiffness)
893
894 ! Compute Stiffness
895 call compute_stiffness_ale(coef, this%config)
896
897 ! Output Stiffness if requested (for diagnostic)
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)
902 call phi_file%free()
903 call field_rzero(rhs_field)
904 end if
905
906 ! Loop over bodies and Solve Laplace
907 do body_idx = 1, this%config%nbodies
908 call mpi_barrier(neko_comm, ierr)
909 sample_start_time = mpi_wtime()
910 call neko_log%message(" Solving laplace for body: " // &
911 trim(this%config%bodies(body_idx)%name))
912
913 call bc_active_body%init_from_components(coef)
914 call bc_inactive_body%init_from_components(coef)
915
916 ! Mark zones
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))
920 end do
921
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))
928 end do
929 end if
930 end do
931
932 call bc_active_body%finalize()
933 call bc_inactive_body%finalize()
934
935 ! The Full list for the solver (Freeze everything to 0 correction)
936 call bcloc%init()
937 call bcloc%append(this%bc_fixed)
938 call bcloc%append(bc_active_body)
939 call bcloc%append(bc_inactive_body)
940
941 ! The "Zeros Only" list for the field (Reset other boundaries)
942 call bcloc_zeros_only%init()
943 call bcloc_zeros_only%append(this%bc_fixed)
944 call bcloc_zeros_only%append(bc_inactive_body)
945
946 call field_rzero(this%base_shapes(body_idx))
947
948 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
949 ! phi = phi_corr + phi_lifted!
950 ! A*phi_corr = -A*phi_lifted !
951 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
952
953 ! Lift BC (Dirichlet = 1.0 on moving body)
954 m = bc_active_body%msk(0)
955 do i = 1, m
956 k = bc_active_body%msk(i)
957 this%base_shapes(body_idx)%x(k, 1, 1, 1) = 1.0_rp
958 end do
959
960 ! Apply Zeros to others.
961 ! This ensures fixed walls and other bodies are 0.0,
962 ! even if they share grid with a moving wall.
963 call bcloc_zeros_only%apply_scalar(this%base_shapes(body_idx)%x, n)
964
965 ! Compute RHS: RHS = -A * Phi_lifted.
966 ! The following is motivated by implementation in Nek5000.
967 call ax%compute(rhs_field%x, this%base_shapes(body_idx)%x, &
968 coef, coef%msh, coef%Xh)
969 call field_cmult(rhs_field, -1.0_rp)
970
971 ! Here we use the FULL list to apply zero Dirichlet BC
972 ! on all boundaries.
973 call bcloc%apply_scalar(rhs_field%x, n)
974 call coef%gs_h%op(rhs_field, gs_op_add)
975
976 ! Solve
977 call field_rzero(corr_field)
978 call pc%update()
979 monitor(1) = ksp%solve(ax, corr_field, &
980 rhs_field%x, n, coef, bcloc, coef%gs_h)
981
982 ! phi = phi_lifted + phi_corr
983 call field_add2(this%base_shapes(body_idx), corr_field, n)
984
985 ! Update Total Phi
986 ! phi_total should be between 0 and 1.
987 if (this%config%nbodies > 1) then
988 call field_add2(this%phi_total, this%base_shapes(body_idx), n)
989 end if
990
991 call mpi_barrier(neko_comm, ierr)
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 ", &
996 sample_time, " (s)"
997
998 call neko_log%message(log_buf)
999
1000 call bc_active_body%free()
1001 call bc_inactive_body%free()
1002 call bcloc%free()
1003 call bcloc_zeros_only%free()
1004
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()
1010 call neko_log%message(' phi_' // &
1011 trim(this%config%bodies(body_idx)%name) // '.fld saved.')
1012 end if
1013 end do
1014
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()
1020 end if
1021 end if
1022
1023 call rhs_field%free()
1024 call corr_field%free()
1025
1026 ! Restore h1/h2 to what they were before
1027 coef%h1 = h1_restore
1028 coef%h2 = h2_restore
1029
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
1034 call ksp%free()
1035 deallocate(ksp)
1036 end if
1037 if (allocated(pc)) then
1038 call precon_destroy(pc)
1039 deallocate(pc)
1040 end if
1041
1042 end subroutine solve_base_mesh_displacement
1043
1046 subroutine update_mesh_velocity(this, coef, time_s)
1047 class(ale_manager_t), intent(inout) :: this
1048 type(coef_t), intent(in) :: coef
1049 type(time_state_t), intent(in) :: time_s
1050 integer :: i
1051 type(body_kinematics_t) :: kin
1052 real(kind=rp) :: rot_mat(3,3)
1053 real(kind=rp) :: rot_center(3)
1054
1055 if (.not. this%active) return
1056 if (.not. this%has_moving_boundary) return
1057 call profiler_start_region('ALE add mesh velocity')
1058
1059 call field_rzero(this%wm_x)
1060 call field_rzero(this%wm_y)
1061 call field_rzero(this%wm_z)
1062
1063 do i = 1, this%config%nbodies
1064 ! Compute kinematics for built-in motions
1065 ! "kin" will be like solid body kinematics at current time
1067 this%config%bodies(i), time_s)
1068
1069 ! User modifier (Superposition or Override)
1070 if (associated(this%user_ale_rigid_kinematics)) then
1071 call this%user_ale_rigid_kinematics(this%config%bodies(i)%id, &
1072 time_s, &
1073 kin%vel_trans, &
1074 kin%vel_ang)
1075 end if
1076
1077 kin%center = this%ale_pivot(i)%pos
1078 this%ale_pivot(i)%vel = kin%vel_trans
1079
1080 this%body_kin(i)%center = this%ale_pivot(i)%pos
1081 this%body_kin(i)%vel_trans = kin%vel_trans
1082 this%body_kin(i)%vel_ang = kin%vel_ang
1083
1084 ! Compute rotation matrix at current time
1085 call this%compute_rotation_matrix(i, time_s)
1086 rot_mat = this%body_rot_matrices(:,:,i)
1087 rot_center = this%config%bodies(i)%rot_center
1088
1089 ! Accumulate contribution from each body and add to mesh velocity
1090 call add_kinematics_to_mesh_velocity(this%wm_x, this%wm_y, &
1091 this%wm_z, this%x_ref, this%y_ref, this%z_ref , &
1092 this%base_shapes(i), coef, kin, rot_mat, rot_center)
1093
1094 ! For checkpointing
1095 call this%prep_checkpoint(i)
1096 end do
1097
1098 ! If user has provided a custom function for mesh velocity.
1099 ! User mesh velocity will be added to the ale computed mesh velocity.
1100 ! This routine should not be used for rigid body motions!
1101 if (associated(this%user_ale_mesh_vel)) then
1102 call this%user_ale_mesh_vel(this%wm_x, this%wm_y, this%wm_z, &
1103 coef, this%x_ref, this%y_ref, this%z_ref, this%base_shapes, time_s)
1104 end if
1105 call profiler_end_region('ALE add mesh velocity')
1106
1107 end subroutine update_mesh_velocity
1108
1110 subroutine advance_mesh(this, coef, time, nadv)
1111 class(ale_manager_t), intent(inout) :: this
1112 type(coef_t), intent(inout) :: coef
1113 type(time_state_t), intent(in) :: time
1114 integer, intent(in) :: nadv
1115 integer :: i
1116
1117 if (.not. this%active) return
1118 if (.not. this%has_moving_boundary) return
1119 call profiler_start_region('ALE update mesh')
1120 do i = 1, this%config%nbodies
1121 ! Advance Point Trackers attached to this body.
1122 ! Can be used for torque calculation (simcomp) at a point distanced
1123 ! from the body.
1124 ! or other purposes like tracking movement (user_check).
1125 call this%ghost_tracker_coord_step(this%body_kin(i), time, nadv, i)
1126 ! Update Pivot Location if requested
1127 call update_pivot_location(this%ale_pivot(i), &
1128 this%ale_pivot(i)%pos, &
1129 this%ale_pivot(i)%vel, &
1130 time, &
1131 nadv, &
1132 this%config%bodies(i))
1133 end do
1134
1135 ! Update lagged B terms (geometry history)
1136 call coef%update_B_history()
1137
1138 ! Update mesh coordinates
1139 call update_ale_mesh(coef, this%wm_x, this%wm_y, this%wm_z, &
1140 this%wm_x_lag, this%wm_y_lag, this%wm_z_lag, &
1141 time, nadv, "ab")
1142
1143 ! Update internal history of mesh velocity.
1144 call this%wm_x_lag%update()
1145 call this%wm_y_lag%update()
1146 call this%wm_z_lag%update()
1147 call profiler_end_region('ALE update mesh')
1148 end subroutine advance_mesh
1149
1150 ! Compute mesh stiffness with per-body gain/decay from stiff_geom.
1151 subroutine compute_stiffness_ale(coef, params)
1152 type(coef_t), intent(inout) :: coef
1153 type(ale_config_t), intent(in) :: params
1154 if (neko_bcknd_device .eq. 1) then
1155 call compute_stiffness_ale_device(coef, params)
1156 else
1157 call compute_stiffness_ale_cpu(coef, params)
1158 end if
1159 end subroutine compute_stiffness_ale
1160
1161 ! Adds kinematics to mesh velocity.
1162 subroutine add_kinematics_to_mesh_velocity(wx, wy, wz, &
1163 x_ref, y_ref, z_ref, phi, coef, kinematics, rot_mat, initial_pivot_loc)
1164 type(field_t), intent(inout) :: wx, wy, wz
1165 type(field_t), intent(in) :: x_ref, y_ref, z_ref
1166 type(field_t), intent(in) :: phi
1167 type(coef_t), intent(in) :: coef
1168 type(body_kinematics_t), intent(in) :: kinematics
1169 real(kind=rp), intent(in) :: initial_pivot_loc(3)
1170 real(kind=rp), intent(in) :: rot_mat(3,3)
1171 if (neko_bcknd_device .eq. 1) then
1173 x_ref, y_ref, z_ref, &
1174 phi, coef, kinematics, rot_mat, initial_pivot_loc)
1175 else
1176 call add_kinematics_to_mesh_velocity_cpu(wx, wy, wz, &
1177 x_ref, y_ref, z_ref, &
1178 phi, coef, kinematics, rot_mat, initial_pivot_loc)
1179 end if
1180 end subroutine add_kinematics_to_mesh_velocity
1181
1182 ! Updates mesh position by integrating mesh velocity in time using AB scheme.
1183 subroutine update_ale_mesh(c_Xh, wm_x, wm_y, wm_z, wm_x_lag, wm_y_lag, &
1184 wm_z_lag, time, nadv, scheme_)
1185 type(coef_t), intent(inout) :: c_xh
1186 type(field_t), intent(in) :: wm_x, wm_y, wm_z
1187 type(field_series_t), intent(in) :: wm_x_lag, wm_y_lag, wm_z_lag
1188 type(time_state_t), intent(in) :: time
1189 integer, intent(in) :: nadv
1190 character(len=*), intent(in) :: scheme_
1191 if (neko_bcknd_device .eq. 1) then
1192 call update_ale_mesh_device(c_xh, wm_x, wm_y, wm_z, wm_x_lag, wm_y_lag, &
1193 wm_z_lag, time, nadv, scheme_)
1194 else
1195 call update_ale_mesh_cpu(c_xh, wm_x, wm_y, wm_z, wm_x_lag, wm_y_lag, &
1196 wm_z_lag, time, nadv, scheme_)
1197 end if
1198 end subroutine update_ale_mesh
1199
1200
1201 subroutine ale_manager_free(this)
1202 class(ale_manager_t), intent(inout), target :: this
1203 integer :: i
1204
1205 if (.not. this%active) return
1206
1207 call this%bc_moving%free()
1208 call this%bc_fixed%free()
1209 call this%bc_list%free()
1210
1211 if (allocated(this%base_shapes)) then
1212 do i = 1, size(this%base_shapes)
1213 call this%base_shapes(i)%free()
1214 end do
1215 deallocate(this%base_shapes)
1216 end if
1217 if (this%config%nbodies > 1) then
1218 call this%phi_total%free()
1219 end if
1220 call this%wm_x_lag%free()
1221 call this%wm_y_lag%free()
1222 call this%wm_z_lag%free()
1223 call this%x_ref%free()
1224 call this%y_ref%free()
1225 call this%z_ref%free()
1226
1227 if (allocated(this%ale_pivot)) deallocate(this%ale_pivot)
1228 if (allocated(this%config%bodies)) deallocate(this%config%bodies)
1229 if (allocated(this%body_kin)) deallocate(this%body_kin)
1230 if (associated(this%global_pivot_pos)) deallocate(this%global_pivot_pos)
1231 if (associated(this%global_pivot_vel_lag)) &
1232 deallocate(this%global_pivot_vel_lag)
1233 if (associated(this%global_basis_pos)) deallocate(this%global_basis_pos)
1234 if (associated(this%global_basis_vel_lag)) &
1235 deallocate(this%global_basis_vel_lag)
1236 if (allocated(this%ghost_handles)) deallocate(this%ghost_handles)
1237 if (allocated(this%body_rot_matrices)) deallocate(this%body_rot_matrices)
1238 if (allocated(this%trackers)) deallocate(this%trackers)
1239 if (associated(neko_ale, this)) nullify(neko_ale)
1240
1241 end subroutine ale_manager_free
1242
1244 subroutine ale_precon_factory(pc, ksp, coef, dof, gs, bclst, pctype, params)
1245 class(pc_t), allocatable, target, intent(inout) :: pc
1246 class(ksp_t), target, intent(inout) :: ksp
1247 type(coef_t), target, intent(in) :: coef
1248 type(dofmap_t), target, intent(in) :: dof
1249 type(gs_t), target, intent(inout) :: gs
1250 type(bc_list_t), target, intent(inout) :: bclst
1251 character(len=*), intent(in) :: pctype
1252 type(json_file), intent(inout) :: params
1253 call precon_factory(pc, pctype)
1254 select type (pcp => pc)
1255 type is (jacobi_t)
1256 call pcp%init(coef, dof, gs)
1257 type is (sx_jacobi_t)
1258 call pcp%init(coef, dof, gs)
1259 type is (device_jacobi_t)
1260 call pcp%init(coef, dof, gs)
1261 type is (hsmg_t)
1262 call pcp%init(coef, bclst, params)
1263 type is (phmg_t)
1264 call pcp%init(coef, bclst, params)
1265 end select
1266 call ksp%set_pc(pc)
1267 end subroutine ale_precon_factory
1268
1269 ! Sets the pivot state at restart.
1270 subroutine set_pivot_restart(this, time_restart)
1271 class(ale_manager_t), intent(inout) :: this
1272 real(kind=dp), intent(in) :: time_restart
1273 type(body_kinematics_t) :: kin_restart
1274 integer :: i, idx, handle_1, handle_2, offset_base
1275 type(time_state_t) :: time_state_dummy
1276 time_state_dummy%t = time_restart
1277
1278 !if (.not. allocated(this%global_pivot_pos)) return
1279
1280 do i = 1, this%config%nbodies
1281
1282 call compute_body_kinematics_built_in(kin_restart, &
1283 this%config%bodies(i), time_state_dummy)
1284
1285 ! User Modifier (Superposition or Override)
1286 if (associated(this%user_ale_rigid_kinematics)) then
1287 call this%user_ale_rigid_kinematics(this%config%bodies(i)%id, &
1288 time_state_dummy, &
1289 kin_restart%vel_trans, &
1290 kin_restart%vel_ang)
1291 end if
1292
1293 this%ale_pivot(i)%vel = kin_restart%vel_trans
1294
1295 idx = (i - 1) * 3
1296 ! Restore Position
1297 this%ale_pivot(i)%pos(1:3) = this%global_pivot_pos(idx + 1:idx + 3)
1298 this%body_kin(i)%center = this%ale_pivot(i)%pos
1299 this%body_kin(i)%vel_trans = kin_restart%vel_trans
1300 this%body_kin(i)%vel_ang = kin_restart%vel_ang
1301
1302 ! Restore Velocity History
1303 this%ale_pivot(i)%vel_lag(1:3, 1:3) = &
1304 this%global_pivot_vel_lag(idx + 1:idx + 3, :)
1305
1306
1307 offset_base = (i-1)*6
1308 handle_1 = this%ghost_handles(1, i)
1309 handle_2 = this%ghost_handles(2, i)
1310
1311 if (handle_1 > 0 .and. handle_1 <= this%n_trackers) then
1312 this%trackers(handle_1)%pos = &
1313 this%global_basis_pos(offset_base + 1 : offset_base + 3)
1314
1315 ! Restore velocity history for ghost-x
1316 this%trackers(handle_1)%vel_lag = &
1317 this%global_basis_vel_lag(offset_base + 1 : offset_base + 3, :)
1318 end if
1319
1320 if (handle_2 > 0 .and. handle_2 <= this%n_trackers) then
1321 this%trackers(handle_2)%pos = &
1322 this%global_basis_pos(offset_base + 4 : offset_base + 6)
1323
1324 ! Restore velocity history for ghost-y
1325 this%trackers(handle_2)%vel_lag = &
1326 this%global_basis_vel_lag(offset_base + 4 : offset_base + 6, :)
1327 end if
1328 end do
1329 end subroutine set_pivot_restart
1330
1331 ! Restores the current coef and related metrics
1332 ! and the pivot states at restart.
1333 subroutine set_coef_restart(this, coef, adv, time_restart)
1334 class(ale_manager_t), intent(inout) :: this
1335 class(advection_t), intent(inout) :: adv
1336 type(coef_t), intent(inout) :: coef
1337 real(kind=dp), intent(in) :: time_restart
1338
1339 if (.not. this%active) return
1340
1341 call this%set_pivot_restart(time_restart)
1342 call coef%recompute_metrics()
1343 call adv%recompute_metrics(coef, .true.)
1344
1345 end subroutine set_coef_restart
1346
1347 subroutine set_pivot_basis_for_checkpoint(this, body_idx)
1348 class(ale_manager_t), intent(inout) :: this
1349 integer, intent(in) :: body_idx
1350 integer :: idx, offset_base, h1, h2
1351
1352 if (.not. this%active) return
1353 if (.not. this%has_moving_boundary) return
1354
1355 idx = (body_idx - 1) * 3
1356 this%global_pivot_pos(idx + 1:idx + 3) = this%ale_pivot(body_idx)%pos(1:3)
1357 this%global_pivot_vel_lag(idx + 1:idx + 3, :) = &
1358 this%ale_pivot(body_idx)%vel_lag(1:3, 1:3)
1359
1360 h1 = this%ghost_handles(1, body_idx)
1361 h2 = this%ghost_handles(2, body_idx)
1362
1363 offset_base = (body_idx-1)*6
1364
1365 ! Save Positions
1366 this%global_basis_pos(offset_base + 1 : offset_base + 3) = &
1367 this%get_tracker_pos(h1)
1368 this%global_basis_pos(offset_base + 4 : offset_base + 6) = &
1369 this%get_tracker_pos(h2)
1370
1371 ! Ghost-x history
1372 this%global_basis_vel_lag(offset_base + 1 : offset_base + 3, :) = &
1373 this%trackers(h1)%vel_lag
1374
1375 ! Ghost-y history
1376 this%global_basis_vel_lag(offset_base + 4 : offset_base + 6, :) = &
1377 this%trackers(h2)%vel_lag
1378 end subroutine set_pivot_basis_for_checkpoint
1379
1380 ! Append val to arr if not already present.
1381 subroutine append_unique_int(arr, n, val)
1382 integer, allocatable, intent(inout) :: arr(:)
1383 integer, intent(inout) :: n
1384 integer, intent(in) :: val
1385 integer, allocatable :: tmp(:)
1386 integer :: k
1387
1388 do k = 1, n
1389 if (arr(k) == val) return
1390 end do
1391
1392 allocate(tmp(n + 1))
1393 if (n > 0) tmp(1:n) = arr(1:n)
1394 tmp(n + 1) = val
1395
1396 if (allocated(arr)) deallocate(arr)
1397 call move_alloc(tmp, arr)
1398
1399 n = n + 1
1400 end subroutine append_unique_int
1401
1402
1404 subroutine mesh_preview(this, coef, json)
1405 class(ale_manager_t), intent(inout) :: this
1406 type(coef_t), intent(inout) :: coef
1407 type(json_file), intent(inout) :: json
1408
1409 logical :: mesh_preview_active
1410 real(kind=rp) :: t_start
1411 real(kind=rp) :: t_end
1412 real(kind=rp) :: dt
1413 integer :: output_freq
1414 integer :: step, n_steps, file_index
1415 type(time_state_t) :: t_state
1416 type(file_t) :: out_file
1417 character(len=128) :: log_buf
1418 type(field_t) :: dummy_field
1419 integer :: nadv, nadv_sim
1420 real(kind=rp) :: min_jac
1421 integer :: n
1422
1423 mesh_preview_active = .false.
1424
1425 if (json%valid_path('case.fluid.ale.mesh_preview.enabled')) then
1426 call json%get('case.fluid.ale.mesh_preview.enabled', &
1427 mesh_preview_active)
1428 end if
1429
1430 if (.not. mesh_preview_active) return
1431
1432 call json_get_or_default(json, 'case.fluid.ale.mesh_preview.start_time', &
1433 t_start, 0.0_rp)
1434 call json_get(json, 'case.fluid.ale.mesh_preview.end_time', &
1435 t_end)
1436 call json_get(json, 'case.fluid.ale.mesh_preview.dt', &
1437 dt)
1438 call json_get(json, &
1439 'case.fluid.ale.mesh_preview.output_freq', &
1440 output_freq)
1441
1442 call neko_log%section("ALE Mesh Preview")
1443 call neko_log%message("Executing mesh motion preview...")
1444
1445
1446
1447 n_steps = int((t_end - t_start) / dt)
1448 call json_get(json, 'case.numerics.time_order', nadv_sim)
1449
1450 write(log_buf, '(A, ES23.15)') ' Start Time : ', t_start
1451 call neko_log%message(log_buf)
1452 write(log_buf, '(A, ES23.15)') ' End Time : ', t_end
1453 call neko_log%message(log_buf)
1454 write(log_buf, '(A, ES23.15)') ' dt : ', dt
1455 call neko_log%message(log_buf)
1456 write(log_buf, '(A, I0)') ' Num Steps : ', n_steps
1457 call neko_log%message(log_buf)
1458 write(log_buf, '(A, I0)') ' Output Freq: ', output_freq
1459 call neko_log%message(log_buf)
1460 call neko_log%message('')
1461
1462 ! Setup Dummy Field for Output
1463 call dummy_field%init(coef%dof, "mesh_preview")
1464 call field_rzero(dummy_field)
1465
1466 ! Time Loop Setup
1467 step = 0
1468 nadv = 1
1469 t_state%t = t_start
1470 t_state%dt = dt
1471 t_state%tstep = 0
1472 t_state%dtlag = dt
1473 n = coef%dof%size()
1474 file_index = -2
1475
1476 min_jac = glmin(coef%jac, n)
1477 call save_mesh_preview_step(coef, dummy_field, out_file, t_state, step, &
1478 file_index)
1479 write(log_buf, '(A,I0, A,ES23.15, A,ES18.11)') &
1480 "Initial Mesh and Mass matrix saved! Step: ", step, " | Time:", &
1481 t_state%t, " | Min Jac: ", min_jac
1482
1483 call neko_log%message(trim(log_buf))
1484 call this%update_mesh_velocity(coef, t_state)
1485
1486 do step = 1, n_steps
1487 t_state%tstep = step
1488 t_state%t = t_start + (step * dt)
1489 nadv = min(step, nadv_sim)
1490
1491 call this%advance_mesh(coef, t_state, nadv)
1492 call coef%recompute_metrics()
1493
1494 min_jac = glmin(coef%jac, n)
1495
1496 if (min_jac <= 0.0_rp) then
1497 write(log_buf, '(A, ES18.11, A, ES23.15)') &
1498 "Negative Jacobian detected (", min_jac, ") at t = ", &
1499 t_state%t
1500 call neko_log%message(log_buf)
1501
1502 call save_mesh_preview_step(coef, dummy_field, out_file, &
1503 t_state, step, file_index)
1504
1505 write(log_buf, '(A,I0, A,ES23.15, A,ES18.11)') &
1506 "Mesh and Mass matrix saved! Step: ", step, " | Time:", &
1507 t_state%t, " | Min Jac:", min_jac
1508 call neko_log%message(trim(log_buf))
1509
1510 call neko_error("ALE Mesh Preview Aborted: Negative Jacobian found.")
1511 end if
1512
1513 if (mod(step, output_freq) == 0) then
1514
1515 call save_mesh_preview_step(coef, dummy_field, out_file, t_state, &
1516 step, file_index)
1517 write(log_buf, '(A,I0, A,ES23.15, A,ES18.11)') &
1518 "Mesh and Mass matrix saved! Step: ", step, " | Time:", &
1519 t_state%t, " | Min Jac:", min_jac
1520 call neko_log%message(trim(log_buf))
1521
1522 end if
1523
1524 call this%update_mesh_velocity(coef, t_state)
1525
1526 end do
1527
1528 call dummy_field%free()
1529 call neko_log%end_section()
1530 call neko_log%message("Mesh preview complete.")
1531 call neko_error("ALE Mesh Preview Finished Successfully.")
1532
1533 end subroutine mesh_preview
1534
1535 subroutine save_mesh_preview_step(coef, dummy_field, out_file, t_state, &
1536 step, file_index)
1537 type(coef_t), intent(inout) :: coef
1538 type(field_t), intent(inout) :: dummy_field
1539 type(file_t), intent(inout) :: out_file
1540 type(time_state_t), intent(in) :: t_state
1541 integer, intent(in) :: step
1542 integer, intent(inout) :: file_index
1543
1544 file_index = file_index + 1
1545 dummy_field%x = coef%B
1546
1547 call out_file%init("mesh_preview.fld")
1548 select type (ft => out_file%file_type)
1549 type is (fld_file_t)
1550 ft%write_mesh = .true.
1551 end select
1552
1553 call out_file%set_counter(file_index)
1554 call out_file%write(dummy_field, t = t_state%t)
1555 call out_file%free()
1556
1557 end subroutine save_mesh_preview_step
1558
1559 ! Asign a tracker point to a body. The tracker moves with body's
1560 ! rigid motion.
1561 function request_tracker(this, initial_pos, body_id) result(handle)
1562 class(ale_manager_t), intent(inout) :: this
1563 real(kind=rp), intent(in) :: initial_pos(3)
1564 integer, intent(in) :: body_id
1565 integer :: handle
1566 type(point_tracker_t), allocatable :: tmp(:)
1567
1568 handle = -100
1569 if (.not. this%active) return
1570 if (.not. this%has_moving_boundary) return
1571
1572 if (.not. allocated(this%trackers)) then
1573 allocate(this%trackers(30))
1574 this%n_trackers = 0
1575 elseif (this%n_trackers >= size(this%trackers)) then
1576 allocate(tmp(size(this%trackers) + 30))
1577 tmp(1:size(this%trackers)) = this%trackers
1578 deallocate(this%trackers)
1579 call move_alloc(tmp, this%trackers)
1580 end if
1581 this%n_trackers = this%n_trackers + 1
1582 handle = this%n_trackers
1583
1584 this%trackers(handle)%pos = initial_pos
1585 this%trackers(handle)%body_id = body_id
1586 this%trackers(handle)%vel_lag = this%ale_pivot(body_id)%vel_lag
1587 end function request_tracker
1588
1589 function get_tracker_pos(this, handle) result(pos)
1590 class(ale_manager_t), intent(in) :: this
1591 integer, intent(in) :: handle
1592 real(kind=rp) :: pos(3)
1593
1594 if (handle > 0 .and. handle <= this%n_trackers) then
1595 pos = this%trackers(handle)%pos
1596 else
1597 pos = 0.0_rp
1598 end if
1599 end function get_tracker_pos
1600
1601
1603 subroutine compute_rotation_matrix(this, body_idx, time)
1604 class(ale_manager_t), intent(inout) :: this
1605 integer, intent(in) :: body_idx
1606 type(time_state_t), intent(in) :: time
1607 integer :: h_x, h_y
1608 real(kind=rp) :: p(3), gx(3), gy(3)
1609 real(kind=rp) :: u(3), v(3), w(3), v_temp(3)
1610
1611 if (.not. this%active) return
1612 if (.not. this%has_moving_boundary) return
1613
1614 ! Get Points
1615 h_x = this%ghost_handles(1, body_idx)
1616 h_y = this%ghost_handles(2, body_idx)
1617
1618 p = this%ale_pivot(body_idx)%pos
1619 gx = this%get_tracker_pos(h_x)
1620 gy = this%get_tracker_pos(h_y)
1621
1622 ! Construct u (New X-axis)
1623 u = gx - p
1624 u = u / sqrt(sum(u**2))
1625
1626 ! Construct w via cross product (Z-axis)
1627 v_temp = gy - p
1628 w(1) = u(2)*v_temp(3) - u(3)*v_temp(2)
1629 w(2) = u(3)*v_temp(1) - u(1)*v_temp(3)
1630 w(3) = u(1)*v_temp(2) - u(2)*v_temp(1)
1631 w = w / sqrt(sum(w**2))
1632
1633 ! Construct v via orthogonalization (Y-axis)
1634 v(1) = w(2)*u(3) - w(3)*u(2)
1635 v(2) = w(3)*u(1) - w(1)*u(3)
1636 v(3) = w(1)*u(2) - w(2)*u(1)
1637
1638 this%body_rot_matrices(:, 1, body_idx) = u
1639 this%body_rot_matrices(:, 2, body_idx) = v
1640 this%body_rot_matrices(:, 3, body_idx) = w
1641
1642 end subroutine compute_rotation_matrix
1643
1644
1648 subroutine log_rot_angles(this, time, body_idxs)
1649 class(ale_manager_t), intent(in) :: this
1650 type(time_state_t), intent(in) :: time
1651 integer, optional, intent(in) :: body_idxs(:)
1652
1653 integer :: i, idx, n_log
1654 real(kind=rp) :: roll_deg, pitch_deg, yaw_deg
1655 real(kind=rp) :: r(3,3)
1656 character(len=256) :: log_buf
1657 real(kind=rp), parameter :: rad_to_deg = 180.0_rp / pi
1658
1659 if (.not. this%active) return
1660 if (.not. this%has_moving_boundary) return
1661
1662 if (present(body_idxs)) then
1663 n_log = size(body_idxs)
1664 else
1665 n_log = this%config%nbodies
1666 end if
1667
1668 call neko_log%message(" ")
1669 call neko_log%message("---------Rotation log---------")
1670 call neko_log%message("variable, time step, time, body, " // &
1671 "x_val, y_val, z_val")
1672
1673 ! If body_idxs is provided, only log those. Otherwise, log all.
1674 do i = 1, n_log
1675
1676 if (present(body_idxs)) then
1677 idx = body_idxs(i)
1678 else
1679 idx = i
1680 end if
1681
1682 r = this%body_rot_matrices(:,:,idx)
1683
1684 ! Angles
1685 yaw_deg = atan2(r(2,1), r(1,1)) * rad_to_deg
1686 pitch_deg = atan2(-r(3,1), sqrt(r(3,2)**2 + r(3,3)**2)) * rad_to_deg
1687 roll_deg = atan2(r(3,2), r(3,3)) * rad_to_deg
1688
1689 ! Log Rotation Angles (Roll, Pitch, Yaw) -> (X, Y, Z)
1690 write(log_buf, '(A, I0, A, ES13.6, A, A, A, 3(ES17.10, :, 2X))') &
1691 "Total_Rot_deg ", time%tstep, " ", time%t, " ", &
1692 trim(this%config%bodies(idx)%name), " ", &
1693 roll_deg, pitch_deg, yaw_deg
1694 call neko_log%message(trim(log_buf))
1695
1696 end do
1697
1698 end subroutine log_rot_angles
1699
1703 subroutine log_pivot(this, time, body_idxs)
1704 class(ale_manager_t), intent(in) :: this
1705 type(time_state_t), intent(in) :: time
1706 integer, optional, intent(in) :: body_idxs(:)
1707 integer :: i, idx, n_log
1708 real(kind=rp) :: pivot_pos(3), pivot_vel(3)
1709 character(len=256) :: log_buf
1710
1711 if (.not. this%active) return
1712 if (.not. this%has_moving_boundary) return
1713
1714 if (present(body_idxs)) then
1715 n_log = size(body_idxs)
1716 else
1717 n_log = this%config%nbodies
1718 end if
1719
1720 call neko_log%message(" ")
1721 call neko_log%message("----------Pivot Log-----------")
1722 call neko_log%message("variable, time step, time, body, " // &
1723 "x_val, y_val, z_val")
1724
1725 ! If body_idxs is provided, only log those. Otherwise, log all.
1726 do i = 1, n_log
1727
1728 if (present(body_idxs)) then
1729 idx = body_idxs(i)
1730 else
1731 idx = i
1732 end if
1733
1734 pivot_pos = this%ale_pivot(idx)%pos
1735 pivot_vel = this%ale_pivot(idx)%vel
1736
1737 ! Pivot Position
1738 write(log_buf, '(A, I0, A, ES13.6, A, A, A, 3(ES17.10, :, 2X))') &
1739 "Total_Pivot_pos ", time%tstep, " ", time%t, " ", &
1740 trim(this%config%bodies(idx)%name), " ", &
1741 this%ale_pivot(idx)%pos
1742 call neko_log%message(trim(log_buf))
1743
1744 ! Pivot Velocity
1745 write(log_buf, '(A, I0, A, ES13.6, A, A, A, 3(ES17.10, :, 2X))') &
1746 "Total_Pivot_vel ", time%tstep, " ", time%t, " ", &
1747 trim(this%config%bodies(idx)%name), " ", &
1748 this%ale_pivot(idx)%vel
1749 call neko_log%message(trim(log_buf))
1750 end do
1751
1752 end subroutine log_pivot
1753
1754 subroutine ghost_tracker_coord_step(this, kin_object, time_s, nadv, body_idx)
1755 class(ale_manager_t), intent(inout) :: this
1756 type(body_kinematics_t), intent(in) :: kin_object
1757 type(time_state_t), intent(in) :: time_s
1758 integer, intent(in) :: nadv
1759 integer, intent(in) :: body_idx
1760 integer :: t
1761 real(kind=rp) :: p_vel(3), rel_pos(3), v_tan(3)
1762
1763 if (.not. this%active) return
1764 if (.not. this%has_moving_boundary) return
1765
1766 if (allocated(this%trackers)) then
1767 do t = 1, this%n_trackers
1768 if (this%trackers(t)%body_id == this%config%bodies(body_idx)%id) then
1769 if (t == this%ghost_handles(1, body_idx) .or. &
1770 t == this%ghost_handles(2, body_idx)) then
1771
1772 ! Calculate the Arm vector (r) at current step
1773 rel_pos = this%trackers(t)%pos - kin_object%center
1774
1775 ! Calculate tangential velocity (Omega \cross r)
1776 v_tan(1) = kin_object%vel_ang(2) * rel_pos(3) - &
1777 kin_object%vel_ang(3) * rel_pos(2)
1778 v_tan(2) = kin_object%vel_ang(3) * rel_pos(1) - &
1779 kin_object%vel_ang(1) * rel_pos(3)
1780 v_tan(3) = kin_object%vel_ang(1) * rel_pos(2) - &
1781 kin_object%vel_ang(2) * rel_pos(1)
1782
1783 ! Total velocity
1784 p_vel = kin_object%vel_trans + v_tan
1785
1786 if (time_s%tstep > 0) then
1787 call ab_integrate_point_pos(this%trackers(t)%pos, &
1788 this%trackers(t)%vel_lag, p_vel, time_s, nadv)
1789 end if
1790
1791 end if
1792
1793 end if
1794 end do
1795
1796 end if
1797 end subroutine ghost_tracker_coord_step
1798
1799 subroutine get_ale_solver_params_json(this, json, ksp_solver, precon_type, &
1800 precon_params, abstol, ksp_max_iter, res_monitor)
1801 class(ale_manager_t), intent(inout) :: this
1802 type(json_file), intent(inout) :: json
1803 character(len=:), allocatable, intent(inout) :: ksp_solver
1804 character(len=:), allocatable, intent(inout) :: precon_type
1805 type(json_file), intent(inout) :: precon_params
1806 real(kind=rp), intent(out) :: abstol
1807 integer, intent(out) :: ksp_max_iter
1808 logical, intent(out) :: res_monitor
1809 logical :: tmp_logical
1810 character(len=:), allocatable :: tmp_str
1811
1812 if (allocated(ksp_solver)) deallocate(ksp_solver)
1813 if (allocated(precon_type)) deallocate(precon_type)
1814
1815 call json_get_or_default(json, 'case.fluid.ale.solver.type', &
1816 ksp_solver, 'cg')
1817
1818 call json_get_or_default(json, &
1819 'case.fluid.ale.solver.preconditioner.type', precon_type, 'jacobi')
1820
1821 if (json%valid_path('case.fluid.ale.solver.preconditioner')) then
1822 call json_get(json, 'case.fluid.ale.solver.preconditioner', &
1823 precon_params)
1824 end if
1825
1826 call json_get_or_default(json, 'case.fluid.ale.solver.absolute_tolerance', &
1827 abstol, 1.0e-10_rp)
1828 call json_get_or_default(json, 'case.fluid.ale.solver.monitor', &
1829 res_monitor, .false.)
1830 call json_get_or_default(json, 'case.fluid.ale.solver.max_iterations', &
1831 ksp_max_iter, 10000)
1832
1833 if (json%valid_path('case.fluid.ale.solver.output_base_shape')) then
1834 call json%get('case.fluid.ale.solver.output_base_shape', tmp_logical)
1835 this%config%if_output_phi = tmp_logical
1836 end if
1837 if (json%valid_path('case.fluid.ale.solver.output_stiffness')) then
1838 call json%get('case.fluid.ale.solver.output_stiffness', tmp_logical)
1839 this%config%if_output_stiffness = tmp_logical
1840 end if
1841
1842 ! Mesh Stiffness
1843 if (json%valid_path('case.fluid.ale.solver.mesh_stiffness.type')) then
1844 call json%get('case.fluid.ale.solver.mesh_stiffness.type', tmp_str)
1845 this%config%stiffness_type = tmp_str
1846 if (.not. (trim(tmp_str) == 'built-in')) then
1847 call neko_error("ALE: stiffness_type must be 'built-in'")
1848 end if
1849 end if
1850 end subroutine get_ale_solver_params_json
1851end module ale_manager
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Retrieves a parameter by name or throws an error.
Abstract interface for user defined ALE base shapes.
Abstract interface for user defined ALE mesh velocity.
Abstract interface for user defined ALE rigid body kinematics.
Subroutines to add advection terms to the RHS of a transport equation.
Definition advection.f90:34
ALE Manager: Handles Mesh Motion.
type(ale_manager_t), pointer, public neko_ale
subroutine ale_manager_free(this)
subroutine, public log_pivot(this, time, body_idxs)
Logs pivot positions for all or selected bodies. can be called in usercompute. eg: call neko_alelog_p...
subroutine ale_manager_init(this, coef, json, user)
Initialize ALE Manager Sets up solver, registers fields, solves for base shape, etc.
subroutine ale_precon_factory(pc, ksp, coef, dof, gs, bclst, pctype, params)
Factory for ALE Preconditioner.
real(kind=rp) function, dimension(3) get_tracker_pos(this, handle)
subroutine set_pivot_basis_for_checkpoint(this, body_idx)
subroutine get_ale_solver_params_json(this, json, ksp_solver, precon_type, precon_params, abstol, ksp_max_iter, res_monitor)
subroutine solve_base_mesh_displacement(this, coef, abstol, ksp_solver, ksp_max_iter, precon_type, precon_params, res_monitor)
Solves the Laplace equation to determine the base shape (phi) for each body. It finds a smooth blendi...
subroutine compute_rotation_matrix(this, body_idx, time)
Computes Rotation Matrix.
subroutine set_coef_restart(this, coef, adv, time_restart)
subroutine, public update_ale_mesh(c_xh, wm_x, wm_y, wm_z, wm_x_lag, wm_y_lag, wm_z_lag, time, nadv, scheme_)
subroutine, public add_kinematics_to_mesh_velocity(wx, wy, wz, x_ref, y_ref, z_ref, phi, coef, kinematics, rot_mat, initial_pivot_loc)
subroutine update_mesh_velocity(this, coef, time_s)
Updates the mesh velocity field based on current time and kinematics Sums contributions from all bodi...
subroutine set_pivot_restart(this, time_restart)
subroutine, public compute_stiffness_ale(coef, params)
subroutine mesh_preview(this, coef, json)
Performs a preview of the mesh motion to verify quality/topology.
subroutine append_unique_int(arr, n, val)
integer function request_tracker(this, initial_pos, body_id)
subroutine, public log_rot_angles(this, time, body_idxs)
Logs rotation angles for all or selected bodies. can be called in usercompute. eg: call neko_alelog_r...
subroutine ghost_tracker_coord_step(this, kin_object, time_s, nadv, body_idx)
subroutine save_mesh_preview_step(coef, dummy_field, out_file, t_state, step, file_index)
subroutine advance_mesh(this, coef, time, nadv)
Main routine to advance the mesh in time.
Defines data structures and algorithms for configuring, calculating, and time-integrating the rigid-b...
subroutine, public compute_body_kinematics_built_in(kinematics, body_conf, time)
Compute built-in kinematics for a body. Uses inputs from JSON. CPU-only.
subroutine, public ab_integrate_point_pos(pos, vel_lag, current_vel, time, nadv)
Advance a single point position (x,y,z) from the point's velocity using AB time-integration.
subroutine, public init_pivot_state(pivot, body_conf)
Initialize pivot state.
subroutine, public update_pivot_location(pivot, pivot_loc, pivot_vel, time, nadv, body_conf)
Updates pivot location.
subroutine, public add_kinematics_to_mesh_velocity_cpu(wx, wy, wz, x_ref, y_ref, z_ref, phi, coef, kinematics, rot_mat, inital_pivot_loc)
Adds kinematics to mesh velocity (CPU)
subroutine, public compute_stiffness_ale_cpu(coef, params)
Compute mesh stiffness with per-body gain/decay from stiff_geom.
subroutine, public update_ale_mesh_cpu(c_xh, wm_x, wm_y, wm_z, wm_x_lag, wm_y_lag, wm_z_lag, time, nadv, scheme_type)
Updates mesh position by integrating mesh velocity in time using AB (CPU)
subroutine, public add_kinematics_to_mesh_velocity_device(wx, wy, wz, x_ref, y_ref, z_ref, phi, coef, kinematics, rot_mat, inital_pivot_loc)
subroutine, public update_ale_mesh_device(c_xh, wm_x, wm_y, wm_z, wm_x_lag, wm_y_lag, wm_z_lag, time, nadv, scheme_type)
subroutine, public compute_stiffness_ale_device(coef, params)
Defines a Matrix-vector product.
Definition ax.f90:34
Defines a list of bc_t.
Definition bc_list.f90:34
Coefficients.
Definition coef.f90:34
Definition comm.F90:1
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:43
Jacobi preconditioner accelerator backend.
Defines a mapping of the degrees of freedom.
Definition dofmap.f90:35
subroutine, public field_rzero(a, n)
Zero a real vector.
subroutine, public field_add2(a, b, n)
Vector addition .
subroutine, public field_cmult(a, c, n)
Multiplication by constant c .
Contains the field_serties_t type.
Defines a field.
Definition field.f90:34
Module for file I/O operations.
Definition file.f90:34
NEKTON fld file format.
Definition fld_file.f90:35
Gather-scatter.
Krylov preconditioner.
Definition pc_hsmg.f90:61
Jacobi preconditioner.
Definition pc_jacobi.f90:34
Utilities for retrieving parameters from the case files.
Implements the base abstract type for Krylov solvers plus helper types.
Definition krylov.f90:34
Logging routines.
Definition log.f90:34
type(log_t), public neko_log
Global log stream.
Definition log.f90:79
integer, parameter, public log_size
Definition log.f90:45
Definition math.f90:60
real(kind=rp), parameter, public pi
Definition math.f90:75
real(kind=rp) function, public glmin(a, n)
Min of a vector of length n.
Definition math.f90:547
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public dp
Definition num_types.f90:9
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Hybrid ph-multigrid preconditioner.
Definition phmg.f90:34
Krylov preconditioner.
Definition precon.f90:34
Profiling interface.
Definition profiler.F90:34
subroutine, public profiler_start_region(name, region_id)
Started a named (name) profiler region.
Definition profiler.F90:78
subroutine, public profiler_end_region(name, region_id)
End the most recently started profiler region.
Definition profiler.F90:107
Defines a registry for storing solution fields.
Definition registry.f90:34
type(registry_t), target, public neko_registry
Global field registry.
Definition registry.f90:144
Jacobi preconditioner SX-Aurora backend.
Module with things related to the simulation time.
Interfaces for user interaction with NEKO.
Definition user_intf.f90:34
Utilities.
Definition utils.f90:35
Defines a zero-valued Dirichlet boundary condition.
Base abstract type for computing the advection operator.
Definition advection.f90:46
Calculated Kinematics for a body at current time.
State history for time-integration of pivots.
Type for a tracked point linked to a body.
Base type for a matrix-vector product providing .
Definition ax.f90:43
A list of allocatable `bc_t`. Follows the standard interface of lists.
Definition bc_list.f90:48
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:62
Defines a jacobi preconditioner.
Stores a series (sequence) of fields, logically connected to a base field, and arranged according to ...
A wrapper around a polymorphic generic_file_t that handles its init. This is essentially a factory fo...
Definition file.f90:56
Interface for NEKTON fld files.
Definition fld_file.f90:65
Defines a jacobi preconditioner.
Definition pc_jacobi.f90:45
Type for storing initial and final residuals in a Krylov solver.
Definition krylov.f90:56
Base abstract type for a canonical Krylov method, solving .
Definition krylov.f90:73
Defines a canonical Krylov preconditioner.
Definition precon.f90:40
Defines a jacobi preconditioner for SX-Aurora.
A struct that contains all info about the time, expand as needed.
A type collecting all the overridable user routines and flag to suppress type injection from custom m...
Zero-valued Dirichlet boundary condition. Used for no-slip walls, but also for various auxillary cond...