Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
drag_torque.f90
Go to the documentation of this file.
1! Copyright (c) 2008-2020, UCHICAGO ARGONNE, LLC.
2! Copyright (c) 2025, The Neko Authors
3!
4! The UChicago Argonne, LLC as Operator of Argonne National
5! Laboratory holds copyright in the Software. The copyright holder
6! reserves all rights except those expressly granted to licensees,
7! and U.S. Government license rights.
8!
9! Redistribution and use in source and binary forms, with or without
10! modification, are permitted provided that the following conditions
11! are met:
12!
13! 1. Redistributions of source code must retain the above copyright
14! notice, this list of conditions and the disclaimer below.
15!
16! 2. Redistributions in binary form must reproduce the above copyright
17! notice, this list of conditions and the disclaimer (as noted below)
18! in the documentation and/or other materials provided with the
19! distribution.
20!
21! 3. Neither the name of ANL nor the names of its contributors
22! may be used to endorse or promote products derived from this software
23! without specific prior written permission.
24!
25! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
26! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
27! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
28! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL
29! UCHICAGO ARGONNE, LLC, THE U.S. DEPARTMENT OF
30! ENERGY OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
31! SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED
32! TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
33! DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
34! THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
35! (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
36! OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37!
38! Additional BSD Notice
39! ---------------------
40! 1. This notice is required to be provided under our contract with
41! the U.S. Department of Energy (DOE). This work was produced at
42! Argonne National Laboratory under Contract
43! No. DE-AC02-06CH11357 with the DOE.
44!
45! 2. Neither the United States Government nor UCHICAGO ARGONNE,
46! LLC nor any of their employees, makes any warranty,
47! express or implied, or assumes any liability or responsibility for the
48! accuracy, completeness, or usefulness of any information, apparatus,
49! product, or process disclosed, or represents that its use would not
50! infringe privately-owned rights.
51!
52! 3. Also, reference herein to any specific commercial products, process,
53! or services by trade name, trademark, manufacturer or otherwise does
54! not necessarily constitute or imply its endorsement, recommendation,
55! or favoring by the United States Government or UCHICAGO ARGONNE LLC.
56! The views and opinions of authors expressed
57! herein do not necessarily state or reflect those of the United States
58! Government or UCHICAGO ARGONNE, LLC, and shall
59! not be used for advertising or product endorsement purposes.
60!
61
63 use field, only : field_t
64 use coefs, only : coef_t
65 use facet_zone, only : facet_zone_t
66 use math, only : rzero, col3, vdot3, col2
67 use space, only : space_t
68 use num_types, only : rp
69 use utils, only : nonlinear_index
70 use iso_c_binding, only : c_ptr
75 use mpi_f08, only : mpi_allreduce, mpi_in_place, mpi_sum
76 use utils, only : neko_error
77 implicit none
78 private
84
85contains
95 subroutine drag_torque_zone(dgtq, tstep, zone, center, s11, s22, s33, s12, &
96 s13, s23, p, coef, visc)
97 integer, intent(in) :: tstep
98 type(facet_zone_t) :: zone
99 type(coef_t), intent(inout) :: coef
100 real(kind=rp), intent(inout) :: s11(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
101 coef%msh%nelv)
102 real(kind=rp), intent(inout) :: s22(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
103 coef%msh%nelv)
104 real(kind=rp), intent(inout) :: s33(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
105 coef%msh%nelv)
106 real(kind=rp), intent(inout) :: s12(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
107 coef%msh%nelv)
108 real(kind=rp), intent(inout) :: s13(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
109 coef%msh%nelv)
110 real(kind=rp), intent(inout) :: s23(coef%Xh%lx, coef%Xh%ly, coef%Xh%lz, &
111 coef%msh%nelv)
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 ! pressure
116 real(kind=rp) :: dragpy = 0.0_rp
117 real(kind=rp) :: dragpz = 0.0_rp
118 real(kind=rp) :: dragvx = 0.0_rp ! viscous
119 real(kind=rp) :: dragvy = 0.0_rp
120 real(kind=rp) :: dragvz = 0.0_rp
121 real(kind=rp) :: torqpx = 0.0_rp ! pressure
122 real(kind=rp) :: torqpy = 0.0_rp
123 real(kind=rp) :: torqpz = 0.0_rp
124 real(kind=rp) :: torqvx = 0.0_rp ! viscous
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
129 dragx = 0.0
130 dragy = 0.0
131 dragz = 0.0
132
133!
134! Fill up viscous array w/ default
135!
136 dragpx = 0.0
137 dragpy = 0.0
138 dragpz = 0.0
139 dragvx = 0.0
140 dragvy = 0.0
141 dragvz = 0.0
142 do mem = 1, zone%size
143 ie = zone%facet_el(mem)%x(2)
144 ifc = zone%facet_el(mem)%x(1)
145 call drag_torque_facet(dgtq, coef%dof%x, coef%dof%y, coef%dof%z, &
146 center, s11, s22, s33, s12, s13, s23, p%x, visc, ifc, ie, coef, &
147 coef%Xh)
148
149 dragpx = dragpx + dgtq(1,1) ! pressure
150 dragpy = dragpy + dgtq(2,1)
151 dragpz = dragpz + dgtq(3,1)
152
153 dragvx = dragvx + dgtq(1,2) ! viscous
154 dragvy = dragvy + dgtq(2,2)
155 dragvz = dragvz + dgtq(3,2)
156
157 torqpx = torqpx + dgtq(1,3) ! pressure
158 torqpy = torqpy + dgtq(2,3)
159 torqpz = torqpz + dgtq(3,3)
160
161 torqvx = torqvx + dgtq(1,4) ! viscous
162 torqvy = torqvy + dgtq(2,4)
163 torqvz = torqvz + dgtq(3,4)
164 end do
165!
166! Sum contributions from all processors
167!
168 call mpi_allreduce(mpi_in_place, dragpx, 1, &
169 mpi_real_precision, mpi_sum, neko_comm, ierr)
170 call mpi_allreduce(mpi_in_place, dragpy, 1, &
171 mpi_real_precision, mpi_sum, neko_comm, ierr)
172 call mpi_allreduce(mpi_in_place, dragpz, 1, &
173 mpi_real_precision, mpi_sum, neko_comm, ierr)
174 call mpi_allreduce(mpi_in_place, dragvx, 1, &
175 mpi_real_precision, mpi_sum, neko_comm, ierr)
176 call mpi_allreduce(mpi_in_place, dragvy, 1, &
177 mpi_real_precision, mpi_sum, neko_comm, ierr)
178 call mpi_allreduce(mpi_in_place, dragvz, 1, &
179 mpi_real_precision, mpi_sum, neko_comm, ierr)
180 !Torque
181 call mpi_allreduce(mpi_in_place, torqpx, 1, &
182 mpi_real_precision, mpi_sum, neko_comm, ierr)
183 call mpi_allreduce(mpi_in_place, torqpy, 1, &
184 mpi_real_precision, mpi_sum, neko_comm, ierr)
185 call mpi_allreduce(mpi_in_place, torqpz, 1, &
186 mpi_real_precision, mpi_sum, neko_comm, ierr)
187 call mpi_allreduce(mpi_in_place, torqvx, 1, &
188 mpi_real_precision, mpi_sum, neko_comm, ierr)
189 call mpi_allreduce(mpi_in_place, torqvy, 1, &
190 mpi_real_precision, mpi_sum, neko_comm, ierr)
191 call mpi_allreduce(mpi_in_place, torqvz, 1, &
192 mpi_real_precision, mpi_sum, neko_comm, ierr)
193
194 dgtq(1,1) = dragpx ! pressure
195 dgtq(2,1) = dragpy
196 dgtq(3,1) = dragpz
197
198 dgtq(1,2) = dragvx ! viscous
199 dgtq(2,2) = dragvy
200 dgtq(3,2) = dragvz
201
202 dgtq(1,3) = torqpx ! pressure
203 dgtq(2,3) = torqpy
204 dgtq(3,3) = torqpz
205
206 dgtq(1,4) = torqvx ! viscous
207 dgtq(2,4) = torqvy
208 dgtq(3,4) = torqvz
209
210 end subroutine drag_torque_zone
211
223 subroutine drag_torque_facet(dgtq, xm0, ym0, zm0, center, s11, s22, s33, &
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
244 integer :: js1
245 integer :: jf1
246 integer :: jskip1
247 integer :: js2
248 integer :: jf2
249 integer :: jskip2
250 real(kind=rp) :: s11_, s12_, s22_, s13_, s23_, s33_
251
252
253 nx = xh%lx
254 ny = xh%ly
255 nz = xh%lz
256 skpdat(1,1) = 1
257 skpdat(2,1) = nx*(ny-1) + 1
258 skpdat(3,1) = nx
259 skpdat(4,1) = 1
260 skpdat(5,1) = ny*(nz-1) + 1
261 skpdat(6,1) = ny
262
263 skpdat(1,2) = 1 + (nx-1)
264 skpdat(2,2) = nx*(ny-1)+1 + (nx-1)
265 skpdat(3,2) = nx
266 skpdat(4,2) = 1
267 skpdat(5,2) = ny*(nz-1)+1
268 skpdat(6,2) = ny
269
270 skpdat(1,3) = 1
271 skpdat(2,3) = nx
272 skpdat(3,3) = 1
273 skpdat(4,3) = 1
274 skpdat(5,3) = ny*(nz-1)+1
275 skpdat(6,3) = ny
276
277 skpdat(1,4) = 1 + nx*(ny-1)
278 skpdat(2,4) = nx + nx*(ny-1)
279 skpdat(3,4) = 1
280 skpdat(4,4) = 1
281 skpdat(5,4) = ny*(nz-1)+1
282 skpdat(6,4) = ny
283
284 skpdat(1,5) = 1
285 skpdat(2,5) = nx
286 skpdat(3,5) = 1
287 skpdat(4,5) = 1
288 skpdat(5,5) = ny
289 skpdat(6,5) = 1
290
291 skpdat(1,6) = 1 + nx*ny*(nz-1)
292 skpdat(2,6) = nx + nx*ny*(nz-1)
293 skpdat(3,6) = 1
294 skpdat(4,6) = 1
295 skpdat(5,6) = ny
296 skpdat(6,6) = 1
297 pf = f
298 js1 = skpdat(1, pf)
299 jf1 = skpdat(2, pf)
300 jskip1 = skpdat(3, pf)
301 js2 = skpdat(4, pf)
302 jf2 = skpdat(5, pf)
303 jskip2 = skpdat(6, pf)
304 call rzero(dgtq,12)
305 i = 0
306 a = 0
307 do j2 = js2, jf2, jskip2
308 do j1 = js1, jf1, jskip1
309 i = i + 1
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)
314 v = visc
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)
324 dgtq = dgtq + dgtq_i
325 end do
326 end do
327 end subroutine drag_torque_facet
328
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
352
353 call rzero(dgtq, 12)
354 s21 = s12
355 s32 = s23
356 s31 = s13
357 !pressure drag
358 dgtq(1,1) = p*n1
359 dgtq(2,1) = p*n2
360 dgtq(3,1) = p*n3
361 ! viscous drag
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)
365 r1 = x - center(1)
366 r2 = y - center(2)
367 r3 = z - center(3)
368 !pressure torque
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)
372 !viscous torque
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)
376 end subroutine drag_torque_pt
377
387 subroutine calc_force_array(force1, force2, force3, force4, force5, force6,&
388 s11, s22, s33, s12, s13, s23,&
389 p, n1, n2, n3, mu, n_pts)
390 integer :: 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)
401 !pressure force
402 call col3(force1, p, n1, n_pts)
403 call col3(force2, p, n2, n_pts)
404 call col3(force3, p, n3, n_pts)
405 ! viscous force
406 v2 = -2.0_rp*mu
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)
413 end subroutine calc_force_array
414
424 subroutine device_calc_force_array(force1, force2, force3,&
425 force4, force5, force6,&
426 s11, s22, s33, s12, s13, s23,&
427 p, n1, n2, n3, mu, n_pts)
428 integer :: 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
433 type(c_ptr) :: mu
434 if (neko_bcknd_device .eq. 1) then
435 call device_rzero(force4, n_pts)
436 call device_rzero(force5, n_pts)
437 call device_rzero(force6, n_pts)
438 !pressure force
439 call device_col3(force1, p, n1, n_pts)
440 call device_col3(force2, p, n2, n_pts)
441 call device_col3(force3, p, n3, n_pts)
442
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)
446 call device_cmult(force4, -2.0_rp, n_pts)
447 call device_cmult(force5, -2.0_rp, n_pts)
448 call device_cmult(force6, -2.0_rp, n_pts)
449 call device_col2(force4, mu, n_pts)
450 call device_col2(force5, mu, n_pts)
451 call device_col2(force6, mu, n_pts)
452 else
453 call neko_error('error in drag_torque, no device bcklnd configured')
454 end if
455 end subroutine device_calc_force_array
456
457
467 subroutine setup_normals(coef, mask, facets, n1, n2, n3, n_pts)
468 type(coef_t) :: coef
469 integer :: n_pts
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)
473 integer :: i
474
475 do i = 1, n_pts
476 fid = facets(i)
477 idx = nonlinear_index(mask(i), coef%Xh%lx, coef%Xh%lx, coef%Xh%lx)
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)
483 end do
484 end subroutine setup_normals
485
486end module drag_torque
__inline__ __device__ void nonlinear_index(const int idx, const int lx, int *index)
Definition bc_utils.h:44
Coefficients.
Definition coef.f90:34
Definition comm.F90:1
type(mpi_datatype), public mpi_real_precision
MPI type for working precision of REAL types.
Definition comm.F90:51
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:43
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.
Defines a field.
Definition field.f90:34
Object for handling masks in Neko.
Definition mask.f90:34
Definition math.f90:60
subroutine, public col2(a, b, n)
Vector multiplication .
Definition math.f90:854
subroutine, public col3(a, b, c, n)
Vector multiplication with 3 vectors .
Definition math.f90:867
subroutine, public vdot3(dot, u1, u2, u3, v1, v2, v3, n)
Compute a dot product (3-d version) assuming vector components etc.
Definition math.f90:684
subroutine, public rzero(a, n)
Zero a real vector.
Definition math.f90:205
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Defines a function space.
Definition space.f90:34
Utilities.
Definition utils.f90:35
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:56
The function space for the SEM solution fields.
Definition space.f90:63