Neko  0.8.99
A portable framework for high-order spectral element flow simulations
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 mesh, only : mesh_t
65  use facet_zone, only : facet_zone_t
66  use comm
67  use math, only : rzero
68  use space, only: space_t
69  use num_types, only: rp
70  implicit none
71  private
76 
77 contains
87  subroutine drag_torque_zone(dgtq, tstep, zone, center, s11, s22, s33, s12, s13, s23,&
88  p, coef, visc)
89  integer, intent(in) :: tstep
90  type(facet_zone_t) :: zone
91  type(coef_t), intent(inout) :: coef
92  real(kind=rp), intent(inout) :: s11(coef%Xh%lx,coef%Xh%lx,coef%Xh%lz,coef%msh%nelv)
93  real(kind=rp), intent(inout) :: s22(coef%Xh%lx,coef%Xh%lx,coef%Xh%lz,coef%msh%nelv)
94  real(kind=rp), intent(inout) :: s33(coef%Xh%lx,coef%Xh%lx,coef%Xh%lz,coef%msh%nelv)
95  real(kind=rp), intent(inout) :: s12(coef%Xh%lx,coef%Xh%lx,coef%Xh%lz,coef%msh%nelv)
96  real(kind=rp), intent(inout) :: s13(coef%Xh%lx,coef%Xh%lx,coef%Xh%lz,coef%msh%nelv)
97  real(kind=rp), intent(inout) :: s23(coef%Xh%lx,coef%Xh%lx,coef%Xh%lz,coef%msh%nelv)
98  type(field_t), intent(inout) :: p
99  real(kind=rp), intent(in) :: visc, center(3)
100  real(kind=rp) :: dgtq(3,4)
101  real(kind=rp) :: dragpx = 0.0_rp ! pressure
102  real(kind=rp) :: dragpy = 0.0_rp
103  real(kind=rp) :: dragpz = 0.0_rp
104  real(kind=rp) :: dragvx = 0.0_rp ! viscous
105  real(kind=rp) :: dragvy = 0.0_rp
106  real(kind=rp) :: dragvz = 0.0_rp
107  real(kind=rp) :: torqpx = 0.0_rp ! pressure
108  real(kind=rp) :: torqpy = 0.0_rp
109  real(kind=rp) :: torqpz = 0.0_rp
110  real(kind=rp) :: torqvx = 0.0_rp ! viscous
111  real(kind=rp) :: torqvy = 0.0_rp
112  real(kind=rp) :: torqvz = 0.0_rp
113  real(kind=rp) :: dragx, dragy, dragz
114  integer :: ie, ifc, mem, ierr
115  dragx = 0.0
116  dragy = 0.0
117  dragz = 0.0
118 
119 !
120 ! Fill up viscous array w/ default
121 !
122  dragpx = 0.0
123  dragpy = 0.0
124  dragpz = 0.0
125  dragvx = 0.0
126  dragvy = 0.0
127  dragvz = 0.0
128  do mem = 1,zone%size
129  ie = zone%facet_el(mem)%x(2)
130  ifc = zone%facet_el(mem)%x(1)
131  call drag_torque_facet(dgtq,coef%dof%x,coef%dof%y,coef%dof%z,&
132  center,&
133  s11, s22, s33, s12, s13, s23,&
134  p%x,visc,ifc,ie, coef, coef%Xh)
135 
136  dragpx = dragpx + dgtq(1,1) ! pressure
137  dragpy = dragpy + dgtq(2,1)
138  dragpz = dragpz + dgtq(3,1)
139 
140  dragvx = dragvx + dgtq(1,2) ! viscous
141  dragvy = dragvy + dgtq(2,2)
142  dragvz = dragvz + dgtq(3,2)
143 
144  torqpx = torqpx + dgtq(1,3) ! pressure
145  torqpy = torqpy + dgtq(2,3)
146  torqpz = torqpz + dgtq(3,3)
147 
148  torqvx = torqvx + dgtq(1,4) ! viscous
149  torqvy = torqvy + dgtq(2,4)
150  torqvz = torqvz + dgtq(3,4)
151  end do
152 !
153 ! Sum contributions from all processors
154 !
155  call mpi_allreduce(mpi_in_place,dragpx, 1, &
156  mpi_real_precision, mpi_sum, neko_comm, ierr)
157  call mpi_allreduce(mpi_in_place,dragpy, 1, &
158  mpi_real_precision, mpi_sum, neko_comm, ierr)
159  call mpi_allreduce(mpi_in_place,dragpz, 1, &
160  mpi_real_precision, mpi_sum, neko_comm, ierr)
161  call mpi_allreduce(mpi_in_place,dragvx, 1, &
162  mpi_real_precision, mpi_sum, neko_comm, ierr)
163  call mpi_allreduce(mpi_in_place,dragvy, 1, &
164  mpi_real_precision, mpi_sum, neko_comm, ierr)
165  call mpi_allreduce(mpi_in_place,dragvz, 1, &
166  mpi_real_precision, mpi_sum, neko_comm, ierr)
167  !Torque
168  call mpi_allreduce(mpi_in_place,torqpx, 1, &
169  mpi_real_precision, mpi_sum, neko_comm, ierr)
170  call mpi_allreduce(mpi_in_place,torqpy, 1, &
171  mpi_real_precision, mpi_sum, neko_comm, ierr)
172  call mpi_allreduce(mpi_in_place,torqpz, 1, &
173  mpi_real_precision, mpi_sum, neko_comm, ierr)
174  call mpi_allreduce(mpi_in_place,torqvx, 1, &
175  mpi_real_precision, mpi_sum, neko_comm, ierr)
176  call mpi_allreduce(mpi_in_place,torqvy, 1, &
177  mpi_real_precision, mpi_sum, neko_comm, ierr)
178  call mpi_allreduce(mpi_in_place,torqvz, 1, &
179  mpi_real_precision, mpi_sum, neko_comm, ierr)
180 
181  dgtq(1,1) = dragpx ! pressure
182  dgtq(2,1) = dragpy
183  dgtq(3,1) = dragpz
184 
185  dgtq(1,2) = dragvx ! viscous
186  dgtq(2,2) = dragvy
187  dgtq(3,2) = dragvz
188 
189  dgtq(1,3) = torqpx ! pressure
190  dgtq(2,3) = torqpy
191  dgtq(3,3) = torqpz
192 
193  dgtq(1,4) = torqvx ! viscous
194  dgtq(2,4) = torqvy
195  dgtq(3,4) = torqvz
196 
197  end subroutine drag_torque_zone
198 
210  subroutine drag_torque_facet(dgtq,xm0,ym0,zm0, center,&
211  s11, s22, s33, s12, s13, s23,&
212  pm1,visc,f,e, coef, Xh)
213  type(coef_t), intent(in) :: coef
214  type(space_t), intent(in) :: xh
215  real(kind=rp), intent(out) :: dgtq(3,4)
216  real(kind=rp), intent(in) :: center(3)
217  real(kind=rp), intent(in) :: xm0(xh%lx,xh%ly,xh%lz,coef%msh%nelv)
218  real(kind=rp), intent(in) :: ym0(xh%lx,xh%ly,xh%lz,coef%msh%nelv)
219  real(kind=rp), intent(in) :: zm0(xh%lx,xh%ly,xh%lz,coef%msh%nelv)
220  real(kind=rp), intent(in) :: s11(xh%lx,xh%ly,xh%lz,coef%msh%nelv)
221  real(kind=rp), intent(in) :: s22(xh%lx,xh%ly,xh%lz,coef%msh%nelv)
222  real(kind=rp), intent(in) :: s33(xh%lx,xh%ly,xh%lz,coef%msh%nelv)
223  real(kind=rp), intent(in) :: s12(xh%lx,xh%ly,xh%lz,coef%msh%nelv)
224  real(kind=rp), intent(in) :: s13(xh%lx,xh%ly,xh%lz,coef%msh%nelv)
225  real(kind=rp), intent(in) :: s23(xh%lx,xh%ly,xh%lz,coef%msh%nelv)
226  real(kind=rp), intent(in) :: pm1(xh%lx,xh%ly,xh%lz,coef%msh%nelv)
227  real(kind=rp), intent(in) :: visc
228  integer, intent(in) :: f,e
229  integer :: pf, i, j1, j2
230  real(kind=rp) :: n1,n2,n3, a, v, dgtq_i(3,4)
231  integer :: skpdat(6,6), nx, ny, nz
232  integer :: js1
233  integer :: jf1
234  integer :: jskip1
235  integer :: js2
236  integer :: jf2
237  integer :: jskip2
238  real(kind=rp) :: s11_, s12_, s22_, s13_, s23_, s33_
239 
240 
241  nx = xh%lx
242  ny = xh%ly
243  nz = xh%lz
244  skpdat(1,1)=1
245  skpdat(2,1)=nx*(ny-1)+1
246  skpdat(3,1)=nx
247  skpdat(4,1)=1
248  skpdat(5,1)=ny*(nz-1)+1
249  skpdat(6,1)=ny
250 
251  skpdat(1,2)=1 + (nx-1)
252  skpdat(2,2)=nx*(ny-1)+1 + (nx-1)
253  skpdat(3,2)=nx
254  skpdat(4,2)=1
255  skpdat(5,2)=ny*(nz-1)+1
256  skpdat(6,2)=ny
257 
258  skpdat(1,3)=1
259  skpdat(2,3)=nx
260  skpdat(3,3)=1
261  skpdat(4,3)=1
262  skpdat(5,3)=ny*(nz-1)+1
263  skpdat(6,3)=ny
264 
265  skpdat(1,4)=1 + nx*(ny-1)
266  skpdat(2,4)=nx + nx*(ny-1)
267  skpdat(3,4)=1
268  skpdat(4,4)=1
269  skpdat(5,4)=ny*(nz-1)+1
270  skpdat(6,4)=ny
271 
272  skpdat(1,5)=1
273  skpdat(2,5)=nx
274  skpdat(3,5)=1
275  skpdat(4,5)=1
276  skpdat(5,5)=ny
277  skpdat(6,5)=1
278 
279  skpdat(1,6)=1 + nx*ny*(nz-1)
280  skpdat(2,6)=nx + nx*ny*(nz-1)
281  skpdat(3,6)=1
282  skpdat(4,6)=1
283  skpdat(5,6)=ny
284  skpdat(6,6)=1
285  pf = f
286  js1 = skpdat(1,pf)
287  jf1 = skpdat(2,pf)
288  jskip1 = skpdat(3,pf)
289  js2 = skpdat(4,pf)
290  jf2 = skpdat(5,pf)
291  jskip2 = skpdat(6,pf)
292  call rzero(dgtq,12)
293  i = 0
294  a = 0
295  do j2=js2,jf2,jskip2
296  do j1=js1,jf1,jskip1
297  i = i+1
298  n1 = coef%nx(i,1,f,e)*coef%area(i,1,f,e)
299  n2 = coef%ny(i,1,f,e)*coef%area(i,1,f,e)
300  n3 = coef%nz(i,1,f,e)*coef%area(i,1,f,e)
301  a = a + coef%area(i,1,f,e)
302  v = visc
303  s11_ = s11(j1,j2,1,e)
304  s12_ = s12(j1,j2,1,e)
305  s22_ = s22(j1,j2,1,e)
306  s13_ = s13(j1,j2,1,e)
307  s23_ = s23(j1,j2,1,e)
308  s33_ = s33(j1,j2,1,e)
309  call drag_torque_pt(dgtq_i,xm0(j1,j2,1,e), ym0(j1,j2,1,e),zm0(j1,j2,1,e), center,&
310  s11_, s22_, s33_, s12_, s13_, s23_,&
311  pm1(j1,j2,1,e), n1, n2, n3, v)
312  dgtq = dgtq + dgtq_i
313  end do
314  end do
315  end subroutine drag_torque_facet
316 
329  subroutine drag_torque_pt(dgtq, x, y, z, center, s11, s22, s33, s12, s13, s23,&
330  p, n1, n2, n3, v)
331  real(kind=rp), intent(inout) :: dgtq(3,4)
332  real(kind=rp), intent(in) :: x
333  real(kind=rp), intent(in) :: y
334  real(kind=rp), intent(in) :: z
335  real(kind=rp), intent(in) :: p
336  real(kind=rp), intent(in) :: v
337  real(kind=rp), intent(in) :: n1, n2, n3, center(3)
338  real(kind=rp), intent(in) :: s11, s12, s22, s13, s23, s33
339  real(kind=rp) :: s21, s31, s32, r1, r2, r3
340  call rzero(dgtq,12)
341  s21 = s12
342  s32 = s23
343  s31 = s13
344  !pressure drag
345  dgtq(1,1) = p*n1
346  dgtq(2,1) = p*n2
347  dgtq(3,1) = p*n3
348  ! viscous drag
349  dgtq(1,2) = -2*v*(s11*n1 + s12*n2 + s13*n3)
350  dgtq(2,2) = -2*v*(s21*n1 + s22*n2 + s23*n3)
351  dgtq(3,2) = -2*v*(s31*n1 + s32*n2 + s33*n3)
352  r1 = x-center(1)
353  r2 = y-center(2)
354  r3 = z-center(3)
355  !pressure torque
356  dgtq(1,3) = (r2*dgtq(3,1)-r3*dgtq(2,1))
357  dgtq(2,3) = (r3*dgtq(1,1)-r1*dgtq(3,1))
358  dgtq(3,3) = (r1*dgtq(2,1)-r2*dgtq(1,1))
359  !viscous torque
360  dgtq(1,4) = (r2*dgtq(3,2)-r3*dgtq(2,2))
361  dgtq(2,4) = (r3*dgtq(1,2)-r1*dgtq(3,2))
362  dgtq(3,4) = (r1*dgtq(2,2)-r2*dgtq(1,2))
363  end subroutine drag_torque_pt
364 
365 end module drag_torque
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:22
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,...
Definition: drag_torque.f90:89
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 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.
Definition: facet_zone.f90:34
Defines a field.
Definition: field.f90:34
Definition: math.f90:60
subroutine, public rzero(a, n)
Zero a real vector.
Definition: math.f90:184
Defines a mesh.
Definition: mesh.f90:34
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Defines a function space.
Definition: space.f90:34
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