Neko 0.9.99
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!
3! The UChicago Argonne, LLC as Operator of Argonne National
4! Laboratory holds copyright in the Software. The copyright holder
5! reserves all rights except those expressly granted to licensees,
6! and U.S. Government license rights.
7!
8! Redistribution and use in source and binary forms, with or without
9! modification, are permitted provided that the following conditions
10! are met:
11!
12! 1. Redistributions of source code must retain the above copyright
13! notice, this list of conditions and the disclaimer below.
14!
15! 2. Redistributions in binary form must reproduce the above copyright
16! notice, this list of conditions and the disclaimer (as noted below)
17! in the documentation and/or other materials provided with the
18! distribution.
19!
20! 3. Neither the name of ANL nor the names of its contributors
21! may be used to endorse or promote products derived from this software
22! without specific prior written permission.
23!
24! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
25! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
26! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
27! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL
28! UCHICAGO ARGONNE, LLC, THE U.S. DEPARTMENT OF
29! ENERGY OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
30! SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED
31! TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
32! DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
33! THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
34! (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
35! OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36!
37! Additional BSD Notice
38! ---------------------
39! 1. This notice is required to be provided under our contract with
40! the U.S. Department of Energy (DOE). This work was produced at
41! Argonne National Laboratory under Contract
42! No. DE-AC02-06CH11357 with the DOE.
43!
44! 2. Neither the United States Government nor UCHICAGO ARGONNE,
45! LLC nor any of their employees, makes any warranty,
46! express or implied, or assumes any liability or responsibility for the
47! accuracy, completeness, or usefulness of any information, apparatus,
48! product, or process disclosed, or represents that its use would not
49! infringe privately-owned rights.
50!
51! 3. Also, reference herein to any specific commercial products, process,
52! or services by trade name, trademark, manufacturer or otherwise does
53! not necessarily constitute or imply its endorsement, recommendation,
54! or favoring by the United States Government or UCHICAGO ARGONNE LLC.
55! The views and opinions of authors expressed
56! herein do not necessarily state or reflect those of the United States
57! Government or UCHICAGO ARGONNE, LLC, and shall
58! not be used for advertising or product endorsement purposes.
59!
60
62 use field, only : field_t
63 use coefs, only : coef_t
64 use facet_zone, only : facet_zone_t
65 use comm
67 use space, only : space_t
68 use num_types, only : rp
69 use utils, only : nonlinear_index
70 use device
73 implicit none
74 private
80
81contains
91 subroutine drag_torque_zone(dgtq, tstep, zone, center, s11, s22, s33, s12, s13, s23,&
92 p, coef, visc)
93 integer, intent(in) :: tstep
94 type(facet_zone_t) :: zone
95 type(coef_t), intent(inout) :: coef
96 real(kind=rp), intent(inout) :: s11(coef%Xh%lx,coef%Xh%lx,coef%Xh%lz,coef%msh%nelv)
97 real(kind=rp), intent(inout) :: s22(coef%Xh%lx,coef%Xh%lx,coef%Xh%lz,coef%msh%nelv)
98 real(kind=rp), intent(inout) :: s33(coef%Xh%lx,coef%Xh%lx,coef%Xh%lz,coef%msh%nelv)
99 real(kind=rp), intent(inout) :: s12(coef%Xh%lx,coef%Xh%lx,coef%Xh%lz,coef%msh%nelv)
100 real(kind=rp), intent(inout) :: s13(coef%Xh%lx,coef%Xh%lx,coef%Xh%lz,coef%msh%nelv)
101 real(kind=rp), intent(inout) :: s23(coef%Xh%lx,coef%Xh%lx,coef%Xh%lz,coef%msh%nelv)
102 type(field_t), intent(inout) :: p
103 real(kind=rp), intent(in) :: visc, center(3)
104 real(kind=rp) :: dgtq(3,4)
105 real(kind=rp) :: dragpx = 0.0_rp ! pressure
106 real(kind=rp) :: dragpy = 0.0_rp
107 real(kind=rp) :: dragpz = 0.0_rp
108 real(kind=rp) :: dragvx = 0.0_rp ! viscous
109 real(kind=rp) :: dragvy = 0.0_rp
110 real(kind=rp) :: dragvz = 0.0_rp
111 real(kind=rp) :: torqpx = 0.0_rp ! pressure
112 real(kind=rp) :: torqpy = 0.0_rp
113 real(kind=rp) :: torqpz = 0.0_rp
114 real(kind=rp) :: torqvx = 0.0_rp ! viscous
115 real(kind=rp) :: torqvy = 0.0_rp
116 real(kind=rp) :: torqvz = 0.0_rp
117 real(kind=rp) :: dragx, dragy, dragz
118 integer :: ie, ifc, mem, ierr
119 dragx = 0.0
120 dragy = 0.0
121 dragz = 0.0
122
123!
124! Fill up viscous array w/ default
125!
126 dragpx = 0.0
127 dragpy = 0.0
128 dragpz = 0.0
129 dragvx = 0.0
130 dragvy = 0.0
131 dragvz = 0.0
132 do mem = 1,zone%size
133 ie = zone%facet_el(mem)%x(2)
134 ifc = zone%facet_el(mem)%x(1)
135 call drag_torque_facet(dgtq,coef%dof%x,coef%dof%y,coef%dof%z,&
136 center,&
137 s11, s22, s33, s12, s13, s23,&
138 p%x,visc,ifc,ie, coef, coef%Xh)
139
140 dragpx = dragpx + dgtq(1,1) ! pressure
141 dragpy = dragpy + dgtq(2,1)
142 dragpz = dragpz + dgtq(3,1)
143
144 dragvx = dragvx + dgtq(1,2) ! viscous
145 dragvy = dragvy + dgtq(2,2)
146 dragvz = dragvz + dgtq(3,2)
147
148 torqpx = torqpx + dgtq(1,3) ! pressure
149 torqpy = torqpy + dgtq(2,3)
150 torqpz = torqpz + dgtq(3,3)
151
152 torqvx = torqvx + dgtq(1,4) ! viscous
153 torqvy = torqvy + dgtq(2,4)
154 torqvz = torqvz + dgtq(3,4)
155 end do
156!
157! Sum contributions from all processors
158!
159 call mpi_allreduce(mpi_in_place,dragpx, 1, &
160 mpi_real_precision, mpi_sum, neko_comm, ierr)
161 call mpi_allreduce(mpi_in_place,dragpy, 1, &
162 mpi_real_precision, mpi_sum, neko_comm, ierr)
163 call mpi_allreduce(mpi_in_place,dragpz, 1, &
164 mpi_real_precision, mpi_sum, neko_comm, ierr)
165 call mpi_allreduce(mpi_in_place,dragvx, 1, &
166 mpi_real_precision, mpi_sum, neko_comm, ierr)
167 call mpi_allreduce(mpi_in_place,dragvy, 1, &
168 mpi_real_precision, mpi_sum, neko_comm, ierr)
169 call mpi_allreduce(mpi_in_place,dragvz, 1, &
170 mpi_real_precision, mpi_sum, neko_comm, ierr)
171 !Torque
172 call mpi_allreduce(mpi_in_place,torqpx, 1, &
173 mpi_real_precision, mpi_sum, neko_comm, ierr)
174 call mpi_allreduce(mpi_in_place,torqpy, 1, &
175 mpi_real_precision, mpi_sum, neko_comm, ierr)
176 call mpi_allreduce(mpi_in_place,torqpz, 1, &
177 mpi_real_precision, mpi_sum, neko_comm, ierr)
178 call mpi_allreduce(mpi_in_place,torqvx, 1, &
179 mpi_real_precision, mpi_sum, neko_comm, ierr)
180 call mpi_allreduce(mpi_in_place,torqvy, 1, &
181 mpi_real_precision, mpi_sum, neko_comm, ierr)
182 call mpi_allreduce(mpi_in_place,torqvz, 1, &
183 mpi_real_precision, mpi_sum, neko_comm, ierr)
184
185 dgtq(1,1) = dragpx ! pressure
186 dgtq(2,1) = dragpy
187 dgtq(3,1) = dragpz
188
189 dgtq(1,2) = dragvx ! viscous
190 dgtq(2,2) = dragvy
191 dgtq(3,2) = dragvz
192
193 dgtq(1,3) = torqpx ! pressure
194 dgtq(2,3) = torqpy
195 dgtq(3,3) = torqpz
196
197 dgtq(1,4) = torqvx ! viscous
198 dgtq(2,4) = torqvy
199 dgtq(3,4) = torqvz
200
201 end subroutine drag_torque_zone
202
214 subroutine drag_torque_facet(dgtq,xm0,ym0,zm0, center,&
215 s11, s22, s33, s12, s13, s23,&
216 pm1,visc,f,e, coef, Xh)
217 type(coef_t), intent(in) :: coef
218 type(space_t), intent(in) :: xh
219 real(kind=rp), intent(out) :: dgtq(3,4)
220 real(kind=rp), intent(in) :: center(3)
221 real(kind=rp), intent(in) :: xm0(xh%lx,xh%ly,xh%lz,coef%msh%nelv)
222 real(kind=rp), intent(in) :: ym0(xh%lx,xh%ly,xh%lz,coef%msh%nelv)
223 real(kind=rp), intent(in) :: zm0(xh%lx,xh%ly,xh%lz,coef%msh%nelv)
224 real(kind=rp), intent(in) :: s11(xh%lx,xh%ly,xh%lz,coef%msh%nelv)
225 real(kind=rp), intent(in) :: s22(xh%lx,xh%ly,xh%lz,coef%msh%nelv)
226 real(kind=rp), intent(in) :: s33(xh%lx,xh%ly,xh%lz,coef%msh%nelv)
227 real(kind=rp), intent(in) :: s12(xh%lx,xh%ly,xh%lz,coef%msh%nelv)
228 real(kind=rp), intent(in) :: s13(xh%lx,xh%ly,xh%lz,coef%msh%nelv)
229 real(kind=rp), intent(in) :: s23(xh%lx,xh%ly,xh%lz,coef%msh%nelv)
230 real(kind=rp), intent(in) :: pm1(xh%lx,xh%ly,xh%lz,coef%msh%nelv)
231 real(kind=rp), intent(in) :: visc
232 integer, intent(in) :: f,e
233 integer :: pf, i, j1, j2
234 real(kind=rp) :: n1,n2,n3, a, v, dgtq_i(3,4)
235 integer :: skpdat(6,6), nx, ny, nz
236 integer :: js1
237 integer :: jf1
238 integer :: jskip1
239 integer :: js2
240 integer :: jf2
241 integer :: jskip2
242 real(kind=rp) :: s11_, s12_, s22_, s13_, s23_, s33_
243
244
245 nx = xh%lx
246 ny = xh%ly
247 nz = xh%lz
248 skpdat(1,1)=1
249 skpdat(2,1)=nx*(ny-1)+1
250 skpdat(3,1)=nx
251 skpdat(4,1)=1
252 skpdat(5,1)=ny*(nz-1)+1
253 skpdat(6,1)=ny
254
255 skpdat(1,2)=1 + (nx-1)
256 skpdat(2,2)=nx*(ny-1)+1 + (nx-1)
257 skpdat(3,2)=nx
258 skpdat(4,2)=1
259 skpdat(5,2)=ny*(nz-1)+1
260 skpdat(6,2)=ny
261
262 skpdat(1,3)=1
263 skpdat(2,3)=nx
264 skpdat(3,3)=1
265 skpdat(4,3)=1
266 skpdat(5,3)=ny*(nz-1)+1
267 skpdat(6,3)=ny
268
269 skpdat(1,4)=1 + nx*(ny-1)
270 skpdat(2,4)=nx + nx*(ny-1)
271 skpdat(3,4)=1
272 skpdat(4,4)=1
273 skpdat(5,4)=ny*(nz-1)+1
274 skpdat(6,4)=ny
275
276 skpdat(1,5)=1
277 skpdat(2,5)=nx
278 skpdat(3,5)=1
279 skpdat(4,5)=1
280 skpdat(5,5)=ny
281 skpdat(6,5)=1
282
283 skpdat(1,6)=1 + nx*ny*(nz-1)
284 skpdat(2,6)=nx + nx*ny*(nz-1)
285 skpdat(3,6)=1
286 skpdat(4,6)=1
287 skpdat(5,6)=ny
288 skpdat(6,6)=1
289 pf = f
290 js1 = skpdat(1,pf)
291 jf1 = skpdat(2,pf)
292 jskip1 = skpdat(3,pf)
293 js2 = skpdat(4,pf)
294 jf2 = skpdat(5,pf)
295 jskip2 = skpdat(6,pf)
296 call rzero(dgtq,12)
297 i = 0
298 a = 0
299 do j2=js2,jf2,jskip2
300 do j1=js1,jf1,jskip1
301 i = i+1
302 n1 = coef%nx(i,1,f,e)*coef%area(i,1,f,e)
303 n2 = coef%ny(i,1,f,e)*coef%area(i,1,f,e)
304 n3 = coef%nz(i,1,f,e)*coef%area(i,1,f,e)
305 a = a + coef%area(i,1,f,e)
306 v = visc
307 s11_ = s11(j1,j2,1,e)
308 s12_ = s12(j1,j2,1,e)
309 s22_ = s22(j1,j2,1,e)
310 s13_ = s13(j1,j2,1,e)
311 s23_ = s23(j1,j2,1,e)
312 s33_ = s33(j1,j2,1,e)
313 call drag_torque_pt(dgtq_i,xm0(j1,j2,1,e), ym0(j1,j2,1,e),zm0(j1,j2,1,e), center,&
314 s11_, s22_, s33_, s12_, s13_, s23_,&
315 pm1(j1,j2,1,e), n1, n2, n3, v)
316 dgtq = dgtq + dgtq_i
317 end do
318 end do
319 end subroutine drag_torque_facet
320
333 subroutine drag_torque_pt(dgtq, x, y, z, center, s11, s22, s33, s12, s13, s23,&
334 p, n1, n2, n3, v)
335 real(kind=rp), intent(inout) :: dgtq(3,4)
336 real(kind=rp), intent(in) :: x
337 real(kind=rp), intent(in) :: y
338 real(kind=rp), intent(in) :: z
339 real(kind=rp), intent(in) :: p
340 real(kind=rp), intent(in) :: v
341 real(kind=rp), intent(in) :: n1, n2, n3, center(3)
342 real(kind=rp), intent(in) :: s11, s12, s22, s13, s23, s33
343 real(kind=rp) :: s21, s31, s32, r1, r2, r3
344 call rzero(dgtq,12)
345 s21 = s12
346 s32 = s23
347 s31 = s13
348 !pressure drag
349 dgtq(1,1) = p*n1
350 dgtq(2,1) = p*n2
351 dgtq(3,1) = p*n3
352 ! viscous drag
353 dgtq(1,2) = -2*v*(s11*n1 + s12*n2 + s13*n3)
354 dgtq(2,2) = -2*v*(s21*n1 + s22*n2 + s23*n3)
355 dgtq(3,2) = -2*v*(s31*n1 + s32*n2 + s33*n3)
356 r1 = x-center(1)
357 r2 = y-center(2)
358 r3 = z-center(3)
359 !pressure torque
360 dgtq(1,3) = (r2*dgtq(3,1)-r3*dgtq(2,1))
361 dgtq(2,3) = (r3*dgtq(1,1)-r1*dgtq(3,1))
362 dgtq(3,3) = (r1*dgtq(2,1)-r2*dgtq(1,1))
363 !viscous torque
364 dgtq(1,4) = (r2*dgtq(3,2)-r3*dgtq(2,2))
365 dgtq(2,4) = (r3*dgtq(1,2)-r1*dgtq(3,2))
366 dgtq(3,4) = (r1*dgtq(2,2)-r2*dgtq(1,2))
367 end subroutine drag_torque_pt
368
378 subroutine calc_force_array(force1, force2, force3, force4, force5, force6,&
379 s11, s22, s33, s12, s13, s23,&
380 p, n1, n2, n3, v, n_pts)
381 integer :: n_pts
382 real(kind=rp), intent(inout),dimension(n_pts) :: force1, force2, force3
383 real(kind=rp), intent(inout),dimension(n_pts) :: force4, force5, force6
384 real(kind=rp), intent(in) :: p(n_pts)
385 real(kind=rp), intent(in) :: v
386 real(kind=rp), intent(in) :: n1(n_pts), n2(n_pts), n3(n_pts)
387 real(kind=rp), intent(in), dimension(n_pts) :: s11, s12, s22, s13, s23, s33
388 real(kind=rp) :: v2
389 call rzero(force4,n_pts)
390 call rzero(force5,n_pts)
391 call rzero(force6,n_pts)
392 !pressure force
393 call col3(force1,p,n1,n_pts)
394 call col3(force2,p,n2,n_pts)
395 call col3(force3,p,n3,n_pts)
396 ! viscous force
397 v2 = -2.0_rp*v
398 call vdot3(force4,s11,s12,s13,n1,n2,n3,n_pts)
399 call vdot3(force5,s12,s22,s23,n1,n2,n3,n_pts)
400 call vdot3(force6,s13,s23,s33,n1,n2,n3,n_pts)
401 call cmult(force4,v2,n_pts)
402 call cmult(force5,v2,n_pts)
403 call cmult(force6,v2,n_pts)
404 end subroutine calc_force_array
414 subroutine device_calc_force_array(force1, force2, force3,&
415 force4, force5, force6,&
416 s11, s22, s33, s12, s13, s23,&
417 p, n1, n2, n3, v, n_pts)
418 integer :: n_pts
419 type(c_ptr) :: force1, force2, force3
420 type(c_ptr) :: force4, force5, force6
421 type(c_ptr) :: p, n1, n2, n3
422 type(c_ptr) :: s11, s12, s22, s13, s23, s33
423 real(kind=rp) :: v
424 real(kind=rp) :: v2
425 if (neko_bcknd_device .eq. 1) then
426 call device_rzero(force4,n_pts)
427 call device_rzero(force5,n_pts)
428 call device_rzero(force6,n_pts)
429 !pressure force
430 call device_col3(force1,p,n1,n_pts)
431 call device_col3(force2,p,n2,n_pts)
432 call device_col3(force3,p,n3,n_pts)
433 ! viscous force
434 v2 = -2.0_rp*v
435 call device_vdot3(force4,s11,s12,s13,n1,n2,n3,n_pts)
436 call device_vdot3(force5,s12,s22,s23,n1,n2,n3,n_pts)
437 call device_vdot3(force6,s13,s23,s33,n1,n2,n3,n_pts)
438 call device_cmult(force4,v2,n_pts)
439 call device_cmult(force5,v2,n_pts)
440 call device_cmult(force6,v2,n_pts)
441 else
442 call neko_error('error in drag_torque, no device bcklnd configured')
443 end if
444 end subroutine device_calc_force_array
445
446
447 subroutine setup_normals(coef,mask,facets,n1,n2,n3,n_pts)
448 type(coef_t) :: coef
449 integer :: n_pts
450 real(kind=rp), dimension(n_pts) :: n1, n2, n3
451 integer :: mask(0:n_pts), facets(0:n_pts), fid, idx(4)
452 real(kind=rp) :: normal(3), area(3)
453 integer :: i
454
455 do i = 1, n_pts
456 fid = facets(i)
457 idx = nonlinear_index(mask(i), coef%Xh%lx, coef%Xh%lx,&
458 coef%Xh%lx)
459 normal = coef%get_normal(idx(1), idx(2), idx(3), idx(4), fid)
460 area = coef%get_area(idx(1), idx(2), idx(3), idx(4), fid)
461 n1(i) = normal(1)*area(1)
462 n2(i) = normal(2)*area(2)
463 n3(i) = normal(3)*area(3)
464
465 end do
466
467 end subroutine setup_normals
468
469end module drag_torque
__device__ void nonlinear_index(const int idx, const int lx, int *index)
Coefficients.
Definition coef.f90:34
Definition comm.F90:1
type(mpi_comm) neko_comm
MPI communicator.
Definition comm.F90:16
type(mpi_datatype) mpi_real_precision
MPI type for working precision of REAL types.
Definition comm.F90:23
subroutine, public device_rzero(a_d, n)
Zero a real vector.
subroutine, public device_col3(a_d, b_d, c_d, n)
Vector multiplication with 3 vectors .
subroutine, public device_masked_red_copy(a_d, b_d, mask_d, n, m)
subroutine, public device_vdot3(dot_d, u1_d, u2_d, u3_d, v1_d, v2_d, v3_d, n)
Compute a dot product (3-d version) assuming vector components etc.
subroutine, public device_cmult(a_d, c, n)
Multiplication by constant c .
Device abstraction, common interface for various accelerators.
Definition device.F90:34
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)
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 device_calc_force_array(force1, force2, force3, force4, force5, force6, s11, s22, s33, s12, s13, s23, p, n1, n2, n3, v, n_pts)
Calculate drag and torque from array of points.
subroutine, public calc_force_array(force1, force2, force3, force4, force5, force6, s11, s22, s33, s12, s13, s23, p, n1, n2, n3, v, 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
Definition math.f90:60
subroutine, public cmult(a, c, n)
Multiplication by constant c .
Definition math.f90:310
subroutine, public masked_red_copy(a, b, mask, n, m)
Copy a masked vector to reduced contigous vector .
Definition math.f90:279
subroutine, public col3(a, b, c, n)
Vector multiplication with 3 vectors .
Definition math.f90:741
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:544
subroutine, public rzero(a, n)
Zero a real vector.
Definition math.f90:194
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:55
The function space for the SEM solution fields.
Definition space.f90:62