Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
force_torque.f90
Go to the documentation of this file.
1! Copyright (c) 2023, The Neko Authors
2! All rights reserved.
3!
4! Redistribution and use in source and binary forms, with or without
5! modification, are permitted provided that the following conditions
6! are met:
7!
8! * Redistributions of source code must retain the above copyright
9! notice, this list of conditions and the following disclaimer.
10!
11! * Redistributions in binary form must reproduce the above
12! copyright notice, this list of conditions and the following
13! disclaimer in the documentation and/or other materials provided
14! with the distribution.
15!
16! * Neither the name of the authors nor the names of its
17! contributors may be used to endorse or promote products derived
18! from this software without specific prior written permission.
19!
20! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31! POSSIBILITY OF SUCH DAMAGE.
32!
33!
35
37 use num_types, only : rp, dp, sp
39 use json_module, only : json_file
41 use registry, only : neko_registry
43 use time_state, only : time_state_t
44 use field, only : field_t
45 use operators, only : curl
46 use case, only : case_t
49 use coefs, only : coef_t
50 use operators, only : strain_rate
51 use vector, only : vector_t
52 use dirichlet, only : dirichlet_t
55 use logger, only : log_size, neko_log
60 use mpi_f08, only : mpi_integer, mpi_sum, mpi_allreduce
61 use comm, only : neko_comm
63 use ale_manager, only : neko_ale
65 use utils, only : neko_error
66
67 implicit none
68 private
69
72 type, public, extends(simulation_component_t) :: force_torque_t
74 type(field_t), pointer :: u => null()
76 type(field_t), pointer :: v => null()
78 type(field_t), pointer :: w => null()
80 type(field_t), pointer :: p => null()
82 type(field_t), pointer :: mu => null()
83
84 ! Masked working arrays
85 type(vector_t) :: n1, n2, n3
86 type(vector_t) :: r1, r2, r3
87 type(vector_t) :: force1, force2, force3
88 type(vector_t) :: force4, force5, force6
89 type(vector_t) :: pmsk
90 type(vector_t) :: mu_msk
91 type(vector_t) :: s11msk, s22msk, s33msk, s12msk, s13msk, s23msk
92 real(kind=rp) :: center(3) = 0.0_rp
93 real(kind=rp) :: scale
94 integer :: zone_id
95 character(len=20) :: zone_name
96 type(coef_t), pointer :: coef => null()
97 type(dirichlet_t) :: bc
98 character(len=80) :: print_format
99 ! Pointer to the live pivot state inside ale_manager
100 type(pivot_state_t), pointer :: pivot_link => null()
101 logical :: moving_center = .false.
102 logical :: update_normals = .false.
103 character(len=64) :: linked_body_name = 'NOT_LINKED'
104 ! Stores the Time=0 offset from the initial pivot
105 real(kind=rp) :: local_offset(3) = 0.0_rp
106 ! Current Pivot Position
107 real(kind=rp), pointer :: body_p(:) => null()
108 ! Current Rotation Matrix
109 real(kind=rp), pointer :: body_r(:,:) => null()
110
111 contains
113 procedure, pass(this) :: init => force_torque_init_from_json
115 procedure, private, pass(this) :: init_common => force_torque_init_common
117 generic :: init_from_components => &
118 init_from_controllers, init_from_controllers_properties
120 procedure, pass(this) :: init_from_controllers => &
124 procedure, pass(this) :: init_from_controllers_properties => &
127 procedure, pass(this) :: free => force_torque_free
129 procedure, pass(this) :: compute_ => force_torque_compute
131 procedure, private, pass(this) :: ale_link => setup_ale_link
132 end type force_torque_t
133
134contains
135
137 subroutine force_torque_init_from_json(this, json, case)
138 class(force_torque_t), intent(inout), target :: this
139 type(json_file), intent(inout) :: json
140 class(case_t), intent(inout), target :: case
141 integer :: zone_id
142 real(kind=rp), allocatable :: center(:)
143 character(len=:), allocatable :: zone_name, fluid_name, center_type
144 character(len=:), allocatable :: name
145 real(kind=rp) :: scale
146 logical :: long_print
147
148 call json_get_or_default(json, "name", name, "force_torque")
149 call this%init_base(json, case)
150
151 call json_get_or_default(json, 'fluid_name', fluid_name, 'fluid')
152 call json_get_or_lookup(json, 'zone_id', zone_id)
153 call json_get_or_default(json, 'zone_name', zone_name, ' ')
154 call json_get_or_lookup_or_default(json, 'scale', scale, 1.0_rp)
155 call json_get_or_default(json, 'long_print', long_print, .false.)
156 call json_get_or_lookup(json, 'center', center)
157 call json_get_or_default(json, 'center_type', center_type, 'fixed')
158 if (trim(center_type) /= 'fixed' .and. &
159 trim(center_type) /= 'pivot' .and. &
160 trim(center_type) /= 'body_attached') then
161 call neko_error("force_torque: center_type must be 'fixed'"// &
162 ", 'pivot', or 'body_attached'.")
163 end if
164
165 if (allocated(center)) this%center = center
166
167 call this%init_common(name, fluid_name, zone_id, zone_name, this%center, &
168 scale, case%fluid%c_xh, long_print, center_type = center_type)
169 end subroutine force_torque_init_from_json
170
185 subroutine force_torque_init_from_controllers(this, name, case, order, &
186 preprocess_controller, compute_controller, output_controller, &
187 fluid_name, zone_id, zone_name, center, scale, coef, long_print)
188 class(force_torque_t), intent(inout) :: this
189 character(len=*), intent(in) :: name
190 class(case_t), intent(inout), target :: case
191 integer :: order
192 type(time_based_controller_t), intent(in) :: preprocess_controller
193 type(time_based_controller_t), intent(in) :: compute_controller
194 type(time_based_controller_t), intent(in) :: output_controller
195 character(len=*), intent(in) :: fluid_name
196 character(len=*), intent(in) :: zone_name
197 integer, intent(in) :: zone_id
198 real(kind=rp), intent(in) :: center(3)
199 real(kind=rp), intent(in) :: scale
200 type(coef_t), target, intent(in) :: coef
201 logical, intent(in) :: long_print
202
203 call this%init_base_from_components(case, order, preprocess_controller, &
204 compute_controller, output_controller)
205 call this%init_common(name, fluid_name, zone_id, zone_name, center, scale, &
206 coef, long_print)
207
209
229 case, order, preprocess_control, preprocess_value, compute_control, &
230 compute_value, output_control, output_value, fluid_name, zone_name, &
231 zone_id, center, scale, coef, long_print)
232 class(force_torque_t), intent(inout) :: this
233 character(len=*), intent(in) :: name
234 class(case_t), intent(inout), target :: case
235 integer :: order
236 character(len=*), intent(in) :: preprocess_control
237 real(kind=rp), intent(in) :: preprocess_value
238 character(len=*), intent(in) :: compute_control
239 real(kind=rp), intent(in) :: compute_value
240 character(len=*), intent(in) :: output_control
241 real(kind=rp), intent(in) :: output_value
242 character(len=*), intent(in) :: fluid_name
243 character(len=*), intent(in) :: zone_name
244 integer, intent(in) :: zone_id
245 real(kind=rp), intent(in) :: center(3)
246 real(kind=rp), intent(in) :: scale
247 type(coef_t), target, intent(in) :: coef
248 logical, intent(in) :: long_print
249
250 call this%init_base_from_components(case, order, preprocess_control, &
251 preprocess_value, compute_control, compute_value, output_control, &
252 output_value)
253 call this%init_common(name, fluid_name, zone_id, zone_name, center, scale, &
254 coef, long_print)
255
257
268 subroutine force_torque_init_common(this, name, fluid_name, zone_id, &
269 zone_name, center, scale, coef, long_print, center_type)
270 class(force_torque_t), intent(inout) :: this
271 character(len=*), intent(in) :: name
272 real(kind=rp), intent(in) :: center(3)
273 real(kind=rp), intent(in) :: scale
274 character(len=*), intent(in) :: fluid_name
275 character(len=*), intent(in) :: zone_name
276 integer, intent(in) :: zone_id
277 type(coef_t), target, intent(in) :: coef
278 logical, intent(in) :: long_print
279 character(len=*), intent(in), optional :: center_type
280 character(len=:), allocatable :: ctype_str
281 integer :: n_pts, glb_n_pts, ierr
282 real(kind=rp) :: avg_r(3)
283 character(len=1000) :: log_buf
284 this%name = name
285 this%coef => coef
286 this%zone_id = zone_id
287 this%scale = scale
288 this%zone_name = zone_name
289
290 if (present(center_type)) then
291 ctype_str = center_type
292 else
293 ctype_str = 'fixed' ! Default behavior
294 end if
295
296 call this%ale_link(zone_id, ctype_str, center)
297
298 if (ctype_str /= 'fixed' .and. .not. this%moving_center) then
299 ctype_str = 'fixed (reverted from ' // ctype_str // ')'
300 end if
301
302 ! Set fixed center if not linked to an ALE body
303 if (.not. this%moving_center) then
304 this%center = center
305 end if
306
307
308 if (long_print) then
309 this%print_format = '(I7,E20.10,E20.10,E20.10,E20.10,A)'
310 else
311 this%print_format = '(I7,E13.5,E13.5,E13.5,E13.5,A)'
312 end if
313
314 this%u => neko_registry%get_field_by_name("u")
315 this%v => neko_registry%get_field_by_name("v")
316 this%w => neko_registry%get_field_by_name("w")
317 this%p => neko_registry%get_field_by_name("p")
318 this%mu => neko_registry%get_field_by_name(fluid_name // '_mu_tot')
319
320
321 call this%bc%init_base(this%coef)
322 call this%bc%mark_zone(this%case%msh%labeled_zones(this%zone_id))
323 call this%bc%finalize()
324 n_pts = this%bc%msk(0)
325 if (n_pts .gt. 0) then
326 call this%n1%init(n_pts)
327 call this%n2%init(n_pts)
328 call this%n3%init(n_pts)
329 call this%r1%init(n_pts)
330 call this%r2%init(n_pts)
331 call this%r3%init(n_pts)
332 call this%force1%init(n_pts)
333 call this%force2%init(n_pts)
334 call this%force3%init(n_pts)
335 call this%force4%init(n_pts)
336 call this%force5%init(n_pts)
337 call this%force6%init(n_pts)
338 call this%s11msk%init(n_pts)
339 call this%s22msk%init(n_pts)
340 call this%s33msk%init(n_pts)
341 call this%s12msk%init(n_pts)
342 call this%s13msk%init(n_pts)
343 call this%s23msk%init(n_pts)
344 call this%pmsk%init(n_pts)
345 call this%mu_msk%init(n_pts)
346 end if
347
348 call setup_normals(this%coef, this%bc%msk, this%bc%facet, &
349 this%n1%x, this%n2%x, this%n3%x, n_pts)
350 call masked_gather_copy_0(this%r1%x, this%coef%dof%x, this%bc%msk, &
351 this%u%size(), n_pts)
352 call masked_gather_copy_0(this%r2%x, this%coef%dof%y, this%bc%msk, &
353 this%u%size(), n_pts)
354 call masked_gather_copy_0(this%r3%x, this%coef%dof%z, this%bc%msk, &
355 this%u%size(), n_pts)
356
357 call mpi_allreduce(n_pts, glb_n_pts, 1, &
358 mpi_integer, mpi_sum, neko_comm, ierr)
359 ! Calculating avg pos here
360 avg_r(1) = glsum(this%r1%x, n_pts)/glb_n_pts
361 avg_r(2) = glsum(this%r2%x, n_pts)/glb_n_pts
362 avg_r(3) = glsum(this%r3%x, n_pts)/glb_n_pts
363 ! Print some information
364 call neko_log%section('Force/torque calculation')
365 write(log_buf, '(A,I4,A,A)') 'Zone ', zone_id, ' ', trim(zone_name)
366 call neko_log%message(log_buf)
367
368 write(log_buf, '(A,I6, I6)') 'Global number of GLL points in zone: ', &
369 glb_n_pts
370 call neko_log%message(log_buf)
371
372 write(log_buf, '(A,A)') 'Center Type: ', trim(ctype_str)
373 call neko_log%message(log_buf)
374
375 if (trim(ctype_str) == 'pivot' .or. &
376 trim(ctype_str) == 'body_attached') then
377 write(log_buf, '(A,A)') 'Linked to ALE Body movement: ', &
378 trim(this%linked_body_name)
379 call neko_log%message(log_buf)
380 write(log_buf, '(A,E15.7,E15.7,E15.7)') &
381 'Initial center for torque calculation: ', this%center
382 call neko_log%message(log_buf)
383 else
384 write(log_buf, '(A,E15.7,E15.7,E15.7)') &
385 'Fixed center for torque calculation: ', this%center
386 call neko_log%message(log_buf)
387 end if
388
389 write(log_buf, '(A,E15.7,E15.7,E15.7)') &
390 'Average of zone''s coordinates: ', avg_r
391 call neko_log%message(log_buf)
392
393 write(log_buf, '(A,E15.7)') 'Scale: ', scale
394 call neko_log%message(log_buf)
395 call neko_log%end_section()
396
397
398 call cadd(this%r1%x, -this%center(1), n_pts)
399 call cadd(this%r2%x, -this%center(2), n_pts)
400 call cadd(this%r3%x, -this%center(3), n_pts)
401 if (neko_bcknd_device .eq. 1 .and. n_pts .gt. 0) then
402 call device_memcpy(this%n1%x, this%n1%x_d, n_pts, host_to_device, &
403 .false.)
404 call device_memcpy(this%n2%x, this%n2%x_d, n_pts, host_to_device, &
405 .false.)
406 call device_memcpy(this%n3%x, this%n3%x_d, n_pts, host_to_device, &
407 .true.)
408 call device_memcpy(this%r1%x, this%r1%x_d, n_pts, host_to_device, &
409 .false.)
410 call device_memcpy(this%r2%x, this%r2%x_d, n_pts, host_to_device, &
411 .false.)
412 call device_memcpy(this%r3%x, this%r3%x_d, n_pts, host_to_device, &
413 .true.)
414 end if
415
416 end subroutine force_torque_init_common
417
419 subroutine force_torque_free(this)
420 class(force_torque_t), intent(inout) :: this
421 call this%free_base()
422
423 call this%n1%free()
424 call this%n2%free()
425 call this%n3%free()
426
427 call this%r1%free()
428 call this%r2%free()
429 call this%r3%free()
430
431 call this%force1%free()
432 call this%force2%free()
433 call this%force3%free()
434
435 call this%force4%free()
436 call this%force5%free()
437 call this%force6%free()
438
439 call this%pmsk%free()
440 call this%mu_msk%free()
441 call this%s11msk%free()
442 call this%s22msk%free()
443 call this%s33msk%free()
444 call this%s12msk%free()
445 call this%s13msk%free()
446 call this%s23msk%free()
447
448 nullify(this%u)
449 nullify(this%v)
450 nullify(this%w)
451 nullify(this%p)
452 nullify(this%coef)
453 nullify(this%mu)
454 nullify(this%pivot_link)
455 end subroutine force_torque_free
456
459 subroutine force_torque_compute(this, time)
460 class(force_torque_t), intent(inout) :: this
461 type(time_state_t), intent(in) :: time
462
463 real(kind=rp) :: dgtq(12) = 0.0_rp
464 integer :: n_pts, temp_indices(6)
465 type(field_t), pointer :: s11, s22, s33, s12, s13, s23
466 character(len=1000) :: log_buf
467 real(kind=rp) :: rot_offset(3)
468 n_pts = this%bc%msk(0)
469
470
471 ! body_attached
472 if (this%moving_center .and. associated(this%body_P)) then
473 if (associated(this%body_R)) then
474 ! R * Offset
475 rot_offset(1) = this%body_R(1,1)*this%local_offset(1) + &
476 this%body_R(1,2)*this%local_offset(2) + &
477 this%body_R(1,3)*this%local_offset(3)
478 rot_offset(2) = this%body_R(2,1)*this%local_offset(1) + &
479 this%body_R(2,2)*this%local_offset(2) + &
480 this%body_R(2,3)*this%local_offset(3)
481 rot_offset(3) = this%body_R(3,1)*this%local_offset(1) + &
482 this%body_R(3,2)*this%local_offset(2) + &
483 this%body_R(3,3)*this%local_offset(3)
484
485 this%center = this%body_P + rot_offset
486 end if
487 end if
488
489 if (this%update_normals) then
490
491 call setup_normals(this%coef, this%bc%msk, this%bc%facet, &
492 this%n1%x, this%n2%x, this%n3%x, n_pts)
493
494 call masked_gather_copy_0(this%r1%x, this%coef%dof%x, this%bc%msk, &
495 this%u%size(), n_pts)
496 call masked_gather_copy_0(this%r2%x, this%coef%dof%y, this%bc%msk, &
497 this%u%size(), n_pts)
498 call masked_gather_copy_0(this%r3%x, this%coef%dof%z, this%bc%msk, &
499 this%u%size(), n_pts)
500
501 call cadd(this%r1%x, -this%center(1), n_pts)
502 call cadd(this%r2%x, -this%center(2), n_pts)
503 call cadd(this%r3%x, -this%center(3), n_pts)
504
505 if (neko_bcknd_device .eq. 1 .and. n_pts .gt. 0) then
506 call device_memcpy(this%n1%x, this%n1%x_d, n_pts, &
507 host_to_device, .false.)
508 call device_memcpy(this%n2%x, this%n2%x_d, n_pts, &
509 host_to_device, .false.)
510 call device_memcpy(this%n3%x, this%n3%x_d, n_pts, &
511 host_to_device, .true.)
512 call device_memcpy(this%r1%x, this%r1%x_d, n_pts, &
513 host_to_device, .false.)
514 call device_memcpy(this%r2%x, this%r2%x_d, n_pts, &
515 host_to_device, .false.)
516 call device_memcpy(this%r3%x, this%r3%x_d, n_pts, &
517 host_to_device, .true.)
518 end if
519 end if
520
521 call neko_scratch_registry%request_field(s11, temp_indices(1), .false.)
522 call neko_scratch_registry%request_field(s12, temp_indices(2), .false.)
523 call neko_scratch_registry%request_field(s13, temp_indices(3), .false.)
524 call neko_scratch_registry%request_field(s22, temp_indices(4), .false.)
525 call neko_scratch_registry%request_field(s23, temp_indices(5), .false.)
526 call neko_scratch_registry%request_field(s33, temp_indices(6), .false.)
527
528 call strain_rate(s11%x, s22%x, s33%x, s12%x, &
529 s13%x, s23%x, this%u, this%v, this%w, this%coef)
530
531 ! On the CPU we can actually just use the original subroutines...
532 if (neko_bcknd_device .eq. 0) then
533 call masked_gather_copy_0(this%s11msk%x, s11%x, this%bc%msk, &
534 this%u%size(), n_pts)
535 call masked_gather_copy_0(this%s22msk%x, s22%x, this%bc%msk, &
536 this%u%size(), n_pts)
537 call masked_gather_copy_0(this%s33msk%x, s33%x, this%bc%msk, &
538 this%u%size(), n_pts)
539 call masked_gather_copy_0(this%s12msk%x, s12%x, this%bc%msk, &
540 this%u%size(), n_pts)
541 call masked_gather_copy_0(this%s13msk%x, s13%x, this%bc%msk, &
542 this%u%size(), n_pts)
543 call masked_gather_copy_0(this%s23msk%x, s23%x, this%bc%msk, &
544 this%u%size(), n_pts)
545 call masked_gather_copy_0(this%pmsk%x, this%p%x, this%bc%msk, &
546 this%u%size(), n_pts)
547 call masked_gather_copy_0(this%mu_msk%x, this%mu%x, this%bc%msk, &
548 this%u%size(), n_pts)
549 call calc_force_array(this%force1%x, this%force2%x, this%force3%x, &
550 this%force4%x, this%force5%x, this%force6%x, &
551 this%s11msk%x, &
552 this%s22msk%x, &
553 this%s33msk%x, &
554 this%s12msk%x, &
555 this%s13msk%x, &
556 this%s23msk%x, &
557 this%pmsk%x, &
558 this%n1%x, &
559 this%n2%x, &
560 this%n3%x, &
561 this%mu_msk%x, &
562 n_pts)
563 dgtq(1) = glsum(this%force1%x, n_pts)
564 dgtq(2) = glsum(this%force2%x, n_pts)
565 dgtq(3) = glsum(this%force3%x, n_pts)
566 dgtq(4) = glsum(this%force4%x, n_pts)
567 dgtq(5) = glsum(this%force5%x, n_pts)
568 dgtq(6) = glsum(this%force6%x, n_pts)
569 call vcross(this%s11msk%x, this%s22msk%x, this%s33msk%x, &
570 this%r1%x, this%r2%x, this%r3%x, &
571 this%force1%x, this%force2%x, this%force3%x, n_pts)
572
573 dgtq(7) = glsum(this%s11msk%x, n_pts)
574 dgtq(8) = glsum(this%s22msk%x, n_pts)
575 dgtq(9) = glsum(this%s33msk%x, n_pts)
576 call vcross(this%s11msk%x, this%s22msk%x, this%s33msk%x, &
577 this%r1%x, this%r2%x, this%r3%x, &
578 this%force4%x, this%force5%x, this%force6%x, n_pts)
579 dgtq(10) = glsum(this%s11msk%x, n_pts)
580 dgtq(11) = glsum(this%s22msk%x, n_pts)
581 dgtq(12) = glsum(this%s33msk%x, n_pts)
582 else
583 if (n_pts .gt. 0) then
584 call device_masked_gather_copy_0(this%s11msk%x_d, s11%x_d, &
585 this%bc%msk_d, this%u%size(), n_pts)
586 call device_masked_gather_copy_0(this%s22msk%x_d, s22%x_d, &
587 this%bc%msk_d, this%u%size(), n_pts)
588 call device_masked_gather_copy_0(this%s33msk%x_d, s33%x_d, &
589 this%bc%msk_d, this%u%size(), n_pts)
590 call device_masked_gather_copy_0(this%s12msk%x_d, s12%x_d, &
591 this%bc%msk_d, this%u%size(), n_pts)
592 call device_masked_gather_copy_0(this%s13msk%x_d, s13%x_d, &
593 this%bc%msk_d, this%u%size(), n_pts)
594 call device_masked_gather_copy_0(this%s23msk%x_d, s23%x_d, &
595 this%bc%msk_d, this%u%size(), n_pts)
596 call device_masked_gather_copy_0(this%pmsk%x_d, this%p%x_d, &
597 this%bc%msk_d, this%u%size(), n_pts)
598 call device_masked_gather_copy_0(this%mu_msk%x_d, this%mu%x_d, &
599 this%bc%msk_d, this%u%size(), n_pts)
600
601 call device_calc_force_array(this%force1%x_d, this%force2%x_d, &
602 this%force3%x_d, &
603 this%force4%x_d, &
604 this%force5%x_d, &
605 this%force6%x_d, &
606 this%s11msk%x_d, &
607 this%s22msk%x_d, &
608 this%s33msk%x_d, &
609 this%s12msk%x_d, &
610 this%s13msk%x_d, &
611 this%s23msk%x_d, &
612 this%pmsk%x_d, &
613 this%n1%x_d, &
614 this%n2%x_d, &
615 this%n3%x_d, &
616 this%mu_msk%x_d, &
617 n_pts)
618 ! Overwriting masked s11, s22, s33 as they are no longer needed
619 call device_vcross(this%s11msk%x_d, this%s22msk%x_d, &
620 this%s33msk%x_d, &
621 this%r1%x_d, this%r2%x_d, this%r3%x_d, &
622 this%force1%x_d, this%force2%x_d, &
623 this%force3%x_d, n_pts)
624 call device_vcross(this%s12msk%x_d, this%s13msk%x_d, this%s23msk%x_d,&
625 this%r1%x_d, this%r2%x_d, this%r3%x_d, &
626 this%force4%x_d, this%force5%x_d, this%force6%x_d, n_pts)
627 end if
628 dgtq(1) = device_glsum(this%force1%x_d, n_pts)
629 dgtq(2) = device_glsum(this%force2%x_d, n_pts)
630 dgtq(3) = device_glsum(this%force3%x_d, n_pts)
631 dgtq(4) = device_glsum(this%force4%x_d, n_pts)
632 dgtq(5) = device_glsum(this%force5%x_d, n_pts)
633 dgtq(6) = device_glsum(this%force6%x_d, n_pts)
634 dgtq(7) = device_glsum(this%s11msk%x_d, n_pts)
635 dgtq(8) = device_glsum(this%s22msk%x_d, n_pts)
636 dgtq(9) = device_glsum(this%s33msk%x_d, n_pts)
637 dgtq(10) = device_glsum(this%s12msk%x_d, n_pts)
638 dgtq(11) = device_glsum(this%s13msk%x_d, n_pts)
639 dgtq(12) = device_glsum(this%s23msk%x_d, n_pts)
640 end if
641 dgtq = this%scale*dgtq
642 write(log_buf, '(A, I4, A, A)') 'Force and torque on zone ', &
643 this%zone_id, ' ', this%zone_name
644 call neko_log%message(log_buf)
645 write(log_buf, '(A)') &
646 'Time step, time, total force/torque, pressure, viscous, direction'
647 call neko_log%message(log_buf)
648 write(log_buf, this%print_format) &
649 time%tstep, time%t, dgtq(1) + dgtq(4), dgtq(1), dgtq(4), ', forcex'
650 call neko_log%message(log_buf)
651 write(log_buf, this%print_format) &
652 time%tstep, time%t, dgtq(2) + dgtq(5), dgtq(2), dgtq(5), ', forcey'
653 call neko_log%message(log_buf)
654 write(log_buf, this%print_format) &
655 time%tstep, time%t, dgtq(3) + dgtq(6), dgtq(3), dgtq(6), ', forcez'
656 call neko_log%message(log_buf)
657 write(log_buf, this%print_format) &
658 time%tstep, time%t, dgtq(7) + dgtq(10), dgtq(7), dgtq(10), ', torquex'
659 call neko_log%message(log_buf)
660 write(log_buf, this%print_format) &
661 time%tstep, time%t, dgtq(8) + dgtq(11), dgtq(8), dgtq(11), ', torquey'
662 call neko_log%message(log_buf)
663 write(log_buf, this%print_format) &
664 time%tstep, time%t, dgtq(9) + dgtq(12), dgtq(9), dgtq(12), ', torquez'
665 call neko_log%message(log_buf)
666 call neko_scratch_registry%relinquish_field(temp_indices)
667
668 end subroutine force_torque_compute
669
671 subroutine setup_ale_link(this, zone_id, center_type, center_in)
672 class(force_torque_t), intent(inout) :: this
673 integer, intent(in) :: zone_id
674 character(len=*), intent(in) :: center_type
675 real(kind=rp), intent(in) :: center_in(3)
676
677 character(len=128) :: log_buf
678 integer :: i, j, nbodies, nindices
679 logical :: body_found
680 logical :: ale_active
681
682 ale_active = .false.
683 this%moving_center = .false.
684 this%update_normals = .false.
685 this%local_offset = 0.0_rp
686 nullify(this%body_P)
687 nullify(this%body_R)
688
689 ! Check Global Pointer for ALE
690 if (associated(neko_ale)) then
691 if (neko_ale%active) then
692 ale_active = .true.
693 this%update_normals = .true.
694
695 if (trim(center_type) == 'pivot' .or. &
696 trim(center_type) == 'body_attached') then
697
698 body_found = .false.
699 nbodies = neko_ale%config%nbodies
700 i = 1
701 do while (i <= nbodies .and. .not. body_found)
702 if (allocated(neko_ale%config%bodies(i)%zone_indices)) then
703 nindices = size(neko_ale%config%bodies(i)%zone_indices)
704 j = 1
705 do while (j <= nindices .and. .not. body_found)
706 if (neko_ale%config%bodies(i)%zone_indices(j) &
707 == zone_id) then
708
709 ! Body found
710 this%moving_center = .true.
711 this%pivot_link => neko_ale%ale_pivot(i)
712 this%linked_body_name = neko_ale%config%bodies(i)%name
713 ! body_id = neko_ale%config%bodies(i)%id
714 body_found = .true.
715
716 ! Point to Live Data
717 this%body_P => neko_ale%ale_pivot(i)%pos
718 this%body_R => neko_ale%body_rot_matrices(:, :, i)
719
720 ! Calculate local offset (Using Time=0 data)
721 if (trim(center_type) == 'pivot') then
722 ! Attached to center -> Offset is 0
723 this%local_offset = 0.0_rp
724 this%center = this%body_P ! Initial Position
725 else
726 ! Detached from pivot -->
727 ! offset = JSON_Input - Initial_Pivot
728 ! Note: we assume center_in is the Time=0 global coord
729 this%local_offset = center_in - &
730 neko_ale%config%bodies(i)%rot_center
731 ! Set initial position for init_common
732 this%center = center_in
733 end if
734
735 end if
736 j = j + 1
737 end do
738 end if
739 i = i + 1
740 end do
741
742 if (.not. body_found) then
743 call neko_log%message(' ')
744 write(log_buf, '(A,I0,A)') 'Warning: Zone ', zone_id, &
745 ' requested "' // trim(center_type) // &
746 '" center, but is not registered as an ALE body.'
747
748 call neko_log%message(log_buf)
749
750 call neko_log%message('Reverting to FIXED center ' // &
751 'using JSON coordinates.')
752
753 this%moving_center = .false.
754 end if
755
756 end if
757 end if
758 end if
759
760 if ((.not. ale_active) .and. (trim(center_type) == 'pivot' .or. &
761 trim(center_type) == 'body_attached')) then
762
763 call neko_log%message(' ')
764
765 write(log_buf, '(A,I0,A)') "Warning: Zone ", zone_id, &
766 " requested '" // trim(center_type) // "' center, " // &
767 "but ALE is not active."
768 call neko_log%message(log_buf)
769
770 call neko_log%message("pivot and body_attached work " // &
771 "only for ALE simulations.")
772
773 call neko_log%message("Reverting to 'fixed' center " // &
774 "using JSON coordinates.")
775
776 end if
777
778 end subroutine setup_ale_link
779
780end module force_torque
Copy data between host and device (or device and device)
Definition device.F90:71
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Retrieves a parameter by name or throws an error.
ALE Manager: Handles Mesh Motion.
type(ale_manager_t), pointer, public neko_ale
Defines data structures and algorithms for configuring, calculating, and time-integrating the rigid-b...
Defines a boundary condition.
Definition bc.f90:34
Defines a simulation case.
Definition case.f90:34
Coefficients.
Definition coef.f90:34
Definition comm.F90:1
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:43
real(kind=rp) function, public device_glsum(a_d, n, strm)
Sum a vector of length n.
subroutine, public device_vcross(u1_d, u2_d, u3_d, v1_d, v2_d, v3_d, w1_d, w2_d, w3_d, n, strm)
Compute a cross product (3-d version) assuming vector components etc.
subroutine, public device_masked_gather_copy_0(a_d, b_d, mask_d, n, n_mask, strm)
Gather a masked vector .
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public host_to_device
Definition device.F90:47
Defines a dirichlet boundary condition.
Definition dirichlet.f90:34
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 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.
Defines a field.
Definition field.f90:34
Implements the force_torque_t type.
subroutine force_torque_compute(this, time)
Compute the force_torque field.
subroutine force_torque_init_from_controllers_properties(this, name, case, order, preprocess_control, preprocess_value, compute_control, compute_value, output_control, output_value, fluid_name, zone_name, zone_id, center, scale, coef, long_print)
Constructor from components, passing properties to the time_based_controller` components in the base ...
subroutine force_torque_free(this)
Destructor.
subroutine setup_ale_link(this, zone_id, center_type, center_in)
Routine to configure ALE connectivity for force/torque module.
subroutine force_torque_init_from_controllers(this, name, case, order, preprocess_controller, compute_controller, output_controller, fluid_name, zone_id, zone_name, center, scale, coef, long_print)
Constructor from components, passing controllers.
subroutine force_torque_init_common(this, name, fluid_name, zone_id, zone_name, center, scale, coef, long_print, center_type)
Common part of constructors.
subroutine force_torque_init_from_json(this, json, case)
Constructor from json.
Utilities for retrieving parameters from the case files.
Logging routines.
Definition log.f90:34
type(log_t), public neko_log
Global log stream.
Definition log.f90:79
integer, parameter, public log_size
Definition log.f90:45
Definition math.f90:60
subroutine, public cadd(a, s, n)
Add a scalar to vector .
Definition math.f90:463
real(kind=rp) function, public glsum(a, n)
Sum a vector of length n.
Definition math.f90:500
subroutine, public masked_gather_copy_0(a, b, mask, n, n_mask)
Gather a masked vector to reduced contigous vector .
Definition math.f90:313
subroutine, public vcross(u1, u2, u3, v1, v2, v3, w1, w2, w3, n)
Compute a cross product assuming vector components etc.
Definition math.f90:654
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public dp
Definition num_types.f90:9
integer, parameter, public sp
Definition num_types.f90:8
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Operators.
Definition operators.f90:34
subroutine, public curl(w1, w2, w3, u1, u2, u3, work1, work2, coef, event)
subroutine, public strain_rate(s11, s22, s33, s12, s13, s23, u, v, w, coef)
Compute the strain rate tensor, i.e 0.5 * du_i/dx_j + du_j/dx_i.
Implements output_controller_t
Defines a registry for storing solution fields.
Definition registry.f90:34
type(registry_t), target, public neko_registry
Global field registry.
Definition registry.f90:144
Defines a registry for storing and requesting temporary objects This can be used when you have a func...
type(scratch_registry_t), target, public neko_scratch_registry
Global scratch registry.
Simulation components are objects that encapsulate functionality that can be fit to a particular comp...
subroutine compute_(this, time)
Dummy compute function.
Contains the time_based_controller_t type.
Module with things related to the simulation time.
Utilities.
Definition utils.f90:35
Defines a vector.
Definition vector.f90:34
State history for time-integration of pivots.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:62
Generic Dirichlet boundary condition on .
Definition dirichlet.f90:49
A simulation component that computes the force and torque on a given boundary zone.
Base abstract class for simulation components.
A utility type for determining whether an action should be executed based on the current time value....
A struct that contains all info about the time, expand as needed.