Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
ale_routines_device.F90
Go to the documentation of this file.
2 use num_types, only : rp
3 use field, only : field_t
4 use coefs, only : coef_t
6 use time_state, only : time_state_t
8 use mesh, only : mesh_t
9 use utils, only : neko_error
10 use device_math, only : device_add2s2
11 use math, only : rzero
16 implicit none
17 private
18
23
24contains
25
26 subroutine compute_stiffness_ale_device(coef, params)
27 type(coef_t), intent(inout) :: coef
28 type(ale_config_t), intent(in) :: params
29 call neko_error("ALE: compute_stiffness_ale_device not implemented yet")
30 end subroutine compute_stiffness_ale_device
31
32 subroutine compute_cheap_dist_device(d, coef, msh, zone_indices)
33 real(kind=rp), intent(inout), target :: d(:)
34 type(coef_t), intent(in) :: coef
35 type(mesh_t), intent(in) :: msh
36 integer, intent(in) :: zone_indices(:)
37 call neko_error("ALE: compute_cheap_dist_device not implemented yet")
38 end subroutine compute_cheap_dist_device
39
41 x_ref, y_ref, z_ref, phi, coef, kinematics, rot_mat, inital_pivot_loc)
42 type(field_t), intent(inout) :: wx, wy, wz
43 type(field_t), intent(in) :: x_ref, y_ref, z_ref
44 type(field_t), intent(in) :: phi
45 type(coef_t), intent(in) :: coef
46 type(body_kinematics_t), intent(in) :: kinematics
47 real(kind=rp), intent(in) :: inital_pivot_loc(3)
48 real(kind=rp), intent(in) :: rot_mat(3,3)
49 call neko_error("ALE: " // &
50 "add_kinematics_to_mesh_velocity_device not implemented yet")
52
53
54 subroutine update_ale_mesh_device(c_Xh, wm_x, wm_y, wm_z, wm_x_lag, &
55 wm_y_lag, wm_z_lag, time, nadv, scheme_type)
56
57 type(coef_t), intent(inout) :: c_xh
58 type(field_t), intent(in) :: wm_x, wm_y, wm_z
59 type(field_series_t), intent(in) :: wm_x_lag, wm_y_lag, wm_z_lag
60 type(time_state_t), intent(in) :: time
61 type(ab_time_scheme_t) :: ab_scheme_obj
62 integer, intent(in) :: nadv
63 integer :: j, n
64 character(len=*), intent(in) :: scheme_type
65 real(kind=rp) :: ab_coeffs(4), dt_history(10), factor
66
67 call rzero(ab_coeffs, 4)
68 if (trim(scheme_type) .eq. 'ab') then
69 dt_history(1) = time%dt
70 dt_history(2) = time%dtlag(1)
71 dt_history(3) = time%dtlag(2)
72 call ab_scheme_obj%compute_coeffs(ab_coeffs, dt_history, nadv)
73 else
74 call neko_error("ALE: Unknown mesh time-integration scheme")
75 end if
76
77 n = c_xh%dof%size()
78
79 factor = time%dt * ab_coeffs(1)
80 ! for now we leave it like this. Can be replaced by a single routine
81 ! that does all three components at once
82 call device_add2s2(c_xh%dof%x_d, wm_x%x_d, factor, n)
83 call device_add2s2(c_xh%dof%y_d, wm_y%x_d, factor, n)
84 call device_add2s2(c_xh%dof%z_d, wm_z%x_d, factor, n)
85
86 ! History Terms
87 do j = 2, nadv
88 factor = time%dt * ab_coeffs(j)
89 call device_add2s2(c_xh%dof%x_d, wm_x_lag%lf(j - 1)%x_d, factor, n)
90 call device_add2s2(c_xh%dof%y_d, wm_y_lag%lf(j - 1)%x_d, factor, n)
91 call device_add2s2(c_xh%dof%z_d, wm_z_lag%lf(j - 1)%x_d, factor, n)
92 end do
93
94 end subroutine update_ale_mesh_device
95
96end module ale_routines_device
Adam-Bashforth scheme for time integration.
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 compute_cheap_dist_device(d, coef, msh, zone_indices)
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)
Coefficients.
Definition coef.f90:34
subroutine, public device_add2s2(a_d, b_d, c1, n, strm)
Vector addition with scalar multiplication (multiplication on first argument)
Contains the field_serties_t type.
Defines a field.
Definition field.f90:34
Definition math.f90:60
subroutine, public rzero(a, n)
Zero a real vector.
Definition math.f90:206
Defines a mesh.
Definition mesh.f90:34
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Module with things related to the simulation time.
Utilities.
Definition utils.f90:35
Explicit Adam-Bashforth scheme for time integration.
Calculated Kinematics for a body at current time.
State history for time-integration of pivots.
Type for a tracked point linked to a body.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:62
Stores a series (sequence) of fields, logically connected to a base field, and arranged according to ...
A struct that contains all info about the time, expand as needed.