70 use iso_c_binding,
only : c_ptr
75 use mpi_f08,
only : mpi_allreduce, mpi_in_place, mpi_sum
96 s13, s23, p, coef, visc)
97 integer,
intent(in) :: tstep
99 type(
coef_t),
intent(inout) :: coef
100 real(kind=
rp),
intent(inout) :: s11(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
102 real(kind=
rp),
intent(inout) :: s22(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
104 real(kind=
rp),
intent(inout) :: s33(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
106 real(kind=
rp),
intent(inout) :: s12(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
108 real(kind=
rp),
intent(inout) :: s13(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
110 real(kind=
rp),
intent(inout) :: s23(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
112 type(
field_t),
intent(inout) :: p
113 real(kind=
rp),
intent(in) :: visc, center(3)
114 real(kind=
rp) :: dgtq(3,4)
115 real(kind=
rp) :: dragpx = 0.0_rp
116 real(kind=
rp) :: dragpy = 0.0_rp
117 real(kind=
rp) :: dragpz = 0.0_rp
118 real(kind=
rp) :: dragvx = 0.0_rp
119 real(kind=
rp) :: dragvy = 0.0_rp
120 real(kind=
rp) :: dragvz = 0.0_rp
121 real(kind=
rp) :: torqpx = 0.0_rp
122 real(kind=
rp) :: torqpy = 0.0_rp
123 real(kind=
rp) :: torqpz = 0.0_rp
124 real(kind=
rp) :: torqvx = 0.0_rp
125 real(kind=
rp) :: torqvy = 0.0_rp
126 real(kind=
rp) :: torqvz = 0.0_rp
127 real(kind=
rp) :: dragx, dragy, dragz
128 integer :: ie, ifc, mem, ierr
142 do mem = 1, zone%size
143 ie = zone%facet_el(mem)%x(2)
144 ifc = zone%facet_el(mem)%x(1)
146 center, s11, s22, s33, s12, s13, s23, p%x, visc, ifc, ie, coef, &
149 dragpx = dragpx + dgtq(1,1)
150 dragpy = dragpy + dgtq(2,1)
151 dragpz = dragpz + dgtq(3,1)
153 dragvx = dragvx + dgtq(1,2)
154 dragvy = dragvy + dgtq(2,2)
155 dragvz = dragvz + dgtq(3,2)
157 torqpx = torqpx + dgtq(1,3)
158 torqpy = torqpy + dgtq(2,3)
159 torqpz = torqpz + dgtq(3,3)
161 torqvx = torqvx + dgtq(1,4)
162 torqvy = torqvy + dgtq(2,4)
163 torqvz = torqvz + dgtq(3,4)
168 call mpi_allreduce(mpi_in_place, dragpx, 1, &
170 call mpi_allreduce(mpi_in_place, dragpy, 1, &
172 call mpi_allreduce(mpi_in_place, dragpz, 1, &
174 call mpi_allreduce(mpi_in_place, dragvx, 1, &
176 call mpi_allreduce(mpi_in_place, dragvy, 1, &
178 call mpi_allreduce(mpi_in_place, dragvz, 1, &
181 call mpi_allreduce(mpi_in_place, torqpx, 1, &
183 call mpi_allreduce(mpi_in_place, torqpy, 1, &
185 call mpi_allreduce(mpi_in_place, torqpz, 1, &
187 call mpi_allreduce(mpi_in_place, torqvx, 1, &
189 call mpi_allreduce(mpi_in_place, torqvy, 1, &
191 call mpi_allreduce(mpi_in_place, torqvz, 1, &
224 s12, s13, s23, pm1, visc, f, e, coef, Xh)
225 type(
coef_t),
intent(in) :: coef
226 type(
space_t),
intent(in) :: xh
227 real(kind=
rp),
intent(out) :: dgtq(3,4)
228 real(kind=
rp),
intent(in) :: center(3)
229 real(kind=
rp),
intent(in) :: xm0(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
230 real(kind=
rp),
intent(in) :: ym0(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
231 real(kind=
rp),
intent(in) :: zm0(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
232 real(kind=
rp),
intent(in) :: s11(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
233 real(kind=
rp),
intent(in) :: s22(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
234 real(kind=
rp),
intent(in) :: s33(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
235 real(kind=
rp),
intent(in) :: s12(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
236 real(kind=
rp),
intent(in) :: s13(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
237 real(kind=
rp),
intent(in) :: s23(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
238 real(kind=
rp),
intent(in) :: pm1(xh%lx, xh%ly, xh%lz, coef%msh%nelv)
239 real(kind=
rp),
intent(in) :: visc
240 integer,
intent(in) :: f, e
241 integer :: pf, i, j1, j2
242 real(kind=
rp) :: n1, n2, n3, a, v, dgtq_i(3,4)
243 integer :: skpdat(6,6), nx, ny, nz
250 real(kind=
rp) :: s11_, s12_, s22_, s13_, s23_, s33_
257 skpdat(2,1) = nx*(ny-1) + 1
260 skpdat(5,1) = ny*(nz-1) + 1
263 skpdat(1,2) = 1 + (nx-1)
264 skpdat(2,2) = nx*(ny-1)+1 + (nx-1)
267 skpdat(5,2) = ny*(nz-1)+1
274 skpdat(5,3) = ny*(nz-1)+1
277 skpdat(1,4) = 1 + nx*(ny-1)
278 skpdat(2,4) = nx + nx*(ny-1)
281 skpdat(5,4) = ny*(nz-1)+1
291 skpdat(1,6) = 1 + nx*ny*(nz-1)
292 skpdat(2,6) = nx + nx*ny*(nz-1)
300 jskip1 = skpdat(3, pf)
303 jskip2 = skpdat(6, pf)
307 do j2 = js2, jf2, jskip2
308 do j1 = js1, jf1, jskip1
310 n1 = coef%nx(i,1,f,e) * coef%area(i,1,f,e)
311 n2 = coef%ny(i,1,f,e) * coef%area(i,1,f,e)
312 n3 = coef%nz(i,1,f,e) * coef%area(i,1,f,e)
313 a = a + coef%area(i,1,f,e)
315 s11_ = s11(j1, j2, 1, e)
316 s12_ = s12(j1, j2, 1, e)
317 s22_ = s22(j1, j2, 1, e)
318 s13_ = s13(j1, j2, 1, e)
319 s23_ = s23(j1, j2, 1, e)
320 s33_ = s33(j1, j2, 1, e)
321 call drag_torque_pt(dgtq_i, xm0(j1, j2, 1, e), ym0(j1, j2, 1, e), &
322 zm0(j1, j2, 1, e), center, s11_, s22_, s33_, s12_, s13_, s23_, &
323 pm1(j1, j2, 1, e), n1, n2, n3, v)
341 subroutine drag_torque_pt(dgtq, x, y, z, center, s11, s22, s33, s12, s13, &
342 s23, p, n1, n2, n3, v)
343 real(kind=
rp),
intent(inout) :: dgtq(3, 4)
344 real(kind=
rp),
intent(in) :: x
345 real(kind=
rp),
intent(in) :: y
346 real(kind=
rp),
intent(in) :: z
347 real(kind=
rp),
intent(in) :: p
348 real(kind=
rp),
intent(in) :: v
349 real(kind=
rp),
intent(in) :: n1, n2, n3, center(3)
350 real(kind=
rp),
intent(in) :: s11, s12, s22, s13, s23, s33
351 real(kind=
rp) :: s21, s31, s32, r1, r2, r3
362 dgtq(1,2) = -2*v*(s11*n1 + s12*n2 + s13*n3)
363 dgtq(2,2) = -2*v*(s21*n1 + s22*n2 + s23*n3)
364 dgtq(3,2) = -2*v*(s31*n1 + s32*n2 + s33*n3)
369 dgtq(1,3) = r2*dgtq(3,1) - r3*dgtq(2,1)
370 dgtq(2,3) = r3*dgtq(1,1) - r1*dgtq(3,1)
371 dgtq(3,3) = r1*dgtq(2,1) - r2*dgtq(1,1)
373 dgtq(1,4) = r2*dgtq(3,2) - r3*dgtq(2,2)
374 dgtq(2,4) = r3*dgtq(1,2) - r1*dgtq(3,2)
375 dgtq(3,4) = r1*dgtq(2,2) - r2*dgtq(1,2)
388 s11, s22, s33, s12, s13, s23,&
389 p, n1, n2, n3, mu, n_pts)
391 real(kind=
rp),
intent(inout),
dimension(n_pts) :: force1, force2, force3
392 real(kind=
rp),
intent(inout),
dimension(n_pts) :: force4, force5, force6
393 real(kind=
rp),
intent(in) :: p(n_pts)
394 real(kind=
rp),
intent(in) :: mu(n_pts)
395 real(kind=
rp),
intent(in) :: n1(n_pts), n2(n_pts), n3(n_pts)
396 real(kind=
rp),
intent(in),
dimension(n_pts) :: s11, s12, s22, s13, s23, s33
397 real(kind=
rp) :: v2(n_pts)
398 call rzero(force4, n_pts)
399 call rzero(force5, n_pts)
400 call rzero(force6, n_pts)
402 call col3(force1, p, n1, n_pts)
403 call col3(force2, p, n2, n_pts)
404 call col3(force3, p, n3, n_pts)
407 call vdot3(force4, s11, s12, s13, n1, n2, n3, n_pts)
408 call vdot3(force5, s12, s22, s23, n1, n2, n3, n_pts)
409 call vdot3(force6, s13, s23, s33, n1, n2, n3, n_pts)
410 call col2(force4, v2, n_pts)
411 call col2(force5, v2, n_pts)
412 call col2(force6, v2, n_pts)
425 force4, force5, force6,&
426 s11, s22, s33, s12, s13, s23,&
427 p, n1, n2, n3, mu, n_pts)
429 type(c_ptr) :: force1, force2, force3
430 type(c_ptr) :: force4, force5, force6
431 type(c_ptr) :: p, n1, n2, n3
432 type(c_ptr) :: s11, s12, s22, s13, s23, s33
443 call device_vdot3(force4, s11, s12, s13, n1, n2, n3, n_pts)
444 call device_vdot3(force5, s12, s22, s23, n1, n2, n3, n_pts)
445 call device_vdot3(force6, s13, s23, s33, n1, n2, n3, n_pts)
453 call neko_error(
'error in drag_torque, no device bcklnd configured')
470 real(kind=
rp),
dimension(n_pts) :: n1, n2, n3
471 integer ::
mask(0:n_pts), facets(0:n_pts), fid, idx(4)
472 real(kind=
rp) :: normal(3), area(3)
478 normal = coef%get_normal(idx(1), idx(2), idx(3), idx(4), fid)
479 area = coef%get_area(idx(1), idx(2), idx(3), idx(4), fid)
480 n1(i) = normal(1)*area(1)
481 n2(i) = normal(2)*area(2)
482 n3(i) = normal(3)*area(3)
__inline__ __device__ void nonlinear_index(const int idx, const int lx, int *index)
type(mpi_datatype), public mpi_real_precision
MPI type for working precision of REAL types.
type(mpi_comm), public neko_comm
MPI communicator.
subroutine, public device_rzero(a_d, n, strm)
Zero a real vector.
subroutine, public device_cmult(a_d, c, n, strm)
Multiplication by constant c .
subroutine, public device_vdot3(dot_d, u1_d, u2_d, u3_d, v1_d, v2_d, v3_d, n, strm)
Compute a dot product (3-d version) assuming vector components etc.
subroutine, public device_col2(a_d, b_d, n, strm)
Vector multiplication .
subroutine, public device_col3(a_d, b_d, c_d, n, strm)
Vector multiplication with 3 vectors .
subroutine, public drag_torque_zone(dgtq, tstep, zone, center, s11, s22, s33, s12, s13, s23, p, coef, visc)
Some functions to calculate the lift/drag and torque Calculation can be done on a zone,...
subroutine, public setup_normals(coef, mask, facets, n1, n2, n3, n_pts)
Computes the normals for a given set of boundary points accessed by the mask.
subroutine, public drag_torque_facet(dgtq, xm0, ym0, zm0, center, s11, s22, s33, s12, s13, s23, pm1, visc, f, e, coef, xh)
Calculate drag and torque over a facet.
subroutine, public calc_force_array(force1, force2, force3, force4, force5, force6, s11, s22, s33, s12, s13, s23, p, n1, n2, n3, mu, n_pts)
Calculate drag and torque from array of points.
subroutine, public device_calc_force_array(force1, force2, force3, force4, force5, force6, s11, s22, s33, s12, s13, s23, p, n1, n2, n3, mu, n_pts)
Calculate drag and torque from array of points.
subroutine, public drag_torque_pt(dgtq, x, y, z, center, s11, s22, s33, s12, s13, s23, p, n1, n2, n3, v)
Calculate drag and torque from one point.
Defines a zone as a subset of facets in a mesh.
Object for handling masks in Neko.
subroutine, public col2(a, b, n)
Vector multiplication .
subroutine, public col3(a, b, c, n)
Vector multiplication with 3 vectors .
subroutine, public vdot3(dot, u1, u2, u3, v1, v2, v3, n)
Compute a dot product (3-d version) assuming vector components etc.
subroutine, public rzero(a, n)
Zero a real vector.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Defines a function space.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
The function space for the SEM solution fields.