Neko 1.99.2
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
64 implicit none
65 private
66
69 type, public, extends(simulation_component_t) :: force_torque_t
71 type(field_t), pointer :: u => null()
73 type(field_t), pointer :: v => null()
75 type(field_t), pointer :: w => null()
77 type(field_t), pointer :: p => null()
79 type(field_t), pointer :: mu => null()
80
81 !Masked working arrays
82 type(vector_t) :: n1, n2, n3
83 type(vector_t) :: r1, r2, r3
84 type(vector_t) :: force1, force2, force3
85 type(vector_t) :: force4, force5, force6
86 type(vector_t) :: pmsk
87 type(vector_t) :: mu_msk
88 type(vector_t) :: s11msk, s22msk, s33msk, s12msk, s13msk, s23msk
89 real(kind=rp) :: center(3) = 0.0_rp
90 real(kind=rp) :: scale
91 integer :: zone_id
92 character(len=20) :: zone_name
93 type(coef_t), pointer :: coef => null()
94 type(dirichlet_t) :: bc
95 character(len=80) :: print_format
96
97 contains
99 procedure, pass(this) :: init => force_torque_init_from_json
101 procedure, private, pass(this) :: init_common => force_torque_init_common
103 generic :: init_from_components => &
104 init_from_controllers, init_from_controllers_properties
106 procedure, pass(this) :: init_from_controllers => &
110 procedure, pass(this) :: init_from_controllers_properties => &
113 procedure, pass(this) :: free => force_torque_free
115 procedure, pass(this) :: compute_ => force_torque_compute
116 end type force_torque_t
117
118contains
119
121 subroutine force_torque_init_from_json(this, json, case)
122 class(force_torque_t), intent(inout), target :: this
123 type(json_file), intent(inout) :: json
124 class(case_t), intent(inout), target :: case
125 integer :: zone_id
126 real(kind=rp), allocatable :: center(:)
127 character(len=:), allocatable :: zone_name, fluid_name
128 character(len=:), allocatable :: name
129 real(kind=rp) :: scale
130 logical :: long_print
131
132 call json_get_or_default(json, "name", name, "force_torque")
133 call this%init_base(json, case)
134
135 call json_get_or_default(json, 'fluid_name', fluid_name, 'fluid')
136 call json_get_or_lookup(json, 'zone_id', zone_id)
137 call json_get_or_default(json, 'zone_name', zone_name, ' ')
138 call json_get_or_lookup_or_default(json, 'scale', scale, 1.0_rp)
139 call json_get_or_default(json, 'long_print', long_print, .false.)
140 call json_get_or_lookup(json, 'center', center)
141 call this%init_common(name, fluid_name, zone_id, zone_name, center, scale, &
142 case%fluid%c_xh, long_print)
143 end subroutine force_torque_init_from_json
144
159 subroutine force_torque_init_from_controllers(this, name, case, order, &
160 preprocess_controller, compute_controller, output_controller, &
161 fluid_name, zone_id, zone_name, center, scale, coef, long_print)
162 class(force_torque_t), intent(inout) :: this
163 character(len=*), intent(in) :: name
164 class(case_t), intent(inout), target :: case
165 integer :: order
166 type(time_based_controller_t), intent(in) :: preprocess_controller
167 type(time_based_controller_t), intent(in) :: compute_controller
168 type(time_based_controller_t), intent(in) :: output_controller
169 character(len=*), intent(in) :: fluid_name
170 character(len=*), intent(in) :: zone_name
171 integer, intent(in) :: zone_id
172 real(kind=rp), intent(in) :: center(3)
173 real(kind=rp), intent(in) :: scale
174 type(coef_t), target, intent(in) :: coef
175 logical, intent(in) :: long_print
176
177 call this%init_base_from_components(case, order, preprocess_controller, &
178 compute_controller, output_controller)
179 call this%init_common(name, fluid_name, zone_id, zone_name, center, scale, &
180 coef, long_print)
181
183
203 case, order, preprocess_control, preprocess_value, compute_control, &
204 compute_value, output_control, output_value, fluid_name, zone_name, &
205 zone_id, center, scale, coef, long_print)
206 class(force_torque_t), intent(inout) :: this
207 character(len=*), intent(in) :: name
208 class(case_t), intent(inout), target :: case
209 integer :: order
210 character(len=*), intent(in) :: preprocess_control
211 real(kind=rp), intent(in) :: preprocess_value
212 character(len=*), intent(in) :: compute_control
213 real(kind=rp), intent(in) :: compute_value
214 character(len=*), intent(in) :: output_control
215 real(kind=rp), intent(in) :: output_value
216 character(len=*), intent(in) :: fluid_name
217 character(len=*), intent(in) :: zone_name
218 integer, intent(in) :: zone_id
219 real(kind=rp), intent(in) :: center(3)
220 real(kind=rp), intent(in) :: scale
221 type(coef_t), target, intent(in) :: coef
222 logical, intent(in) :: long_print
223
224 call this%init_base_from_components(case, order, preprocess_control, &
225 preprocess_value, compute_control, compute_value, output_control, &
226 output_value)
227 call this%init_common(name, fluid_name, zone_id, zone_name, center, scale, &
228 coef, long_print)
229
231
241 subroutine force_torque_init_common(this, name, fluid_name, zone_id, &
242 zone_name, center, scale, coef, long_print)
243 class(force_torque_t), intent(inout) :: this
244 character(len=*), intent(in) :: name
245 real(kind=rp), intent(in) :: center(3)
246 real(kind=rp), intent(in) :: scale
247 character(len=*), intent(in) :: fluid_name
248 character(len=*), intent(in) :: zone_name
249 integer, intent(in) :: zone_id
250 type(coef_t), target, intent(in) :: coef
251 logical, intent(in) :: long_print
252 integer :: n_pts, glb_n_pts, ierr
253 character(len=1000) :: log_buf
254
255 this%name = name
256 this%coef => coef
257 this%zone_id = zone_id
258 this%center = center
259 this%scale = scale
260 this%zone_name = zone_name
261
262 if (long_print ) then
263 this%print_format = '(I7,E20.10,E20.10,E20.10,E20.10,A)'
264 else
265 this%print_format = '(I7,E13.5,E13.5,E13.5,E13.5,A)'
266 end if
267
268 this%u => neko_registry%get_field_by_name("u")
269 this%v => neko_registry%get_field_by_name("v")
270 this%w => neko_registry%get_field_by_name("w")
271 this%p => neko_registry%get_field_by_name("p")
272 this%mu => neko_registry%get_field_by_name(fluid_name // '_mu_tot')
273
274
275 call this%bc%init_base(this%coef)
276 call this%bc%mark_zone(this%case%msh%labeled_zones(this%zone_id))
277 call this%bc%finalize()
278 n_pts = this%bc%msk(0)
279 if (n_pts .gt. 0) then
280 call this%n1%init(n_pts)
281 call this%n2%init(n_pts)
282 call this%n3%init(n_pts)
283 call this%r1%init(n_pts)
284 call this%r2%init(n_pts)
285 call this%r3%init(n_pts)
286 call this%force1%init(n_pts)
287 call this%force2%init(n_pts)
288 call this%force3%init(n_pts)
289 call this%force4%init(n_pts)
290 call this%force5%init(n_pts)
291 call this%force6%init(n_pts)
292 call this%s11msk%init(n_pts)
293 call this%s22msk%init(n_pts)
294 call this%s33msk%init(n_pts)
295 call this%s12msk%init(n_pts)
296 call this%s13msk%init(n_pts)
297 call this%s23msk%init(n_pts)
298 call this%pmsk%init(n_pts)
299 call this%mu_msk%init(n_pts)
300 end if
301
302 call setup_normals(this%coef, this%bc%msk, this%bc%facet, &
303 this%n1%x, this%n2%x, this%n3%x, n_pts)
304 call masked_gather_copy_0(this%r1%x, this%coef%dof%x, this%bc%msk, &
305 this%u%size(), n_pts)
306 call masked_gather_copy_0(this%r2%x, this%coef%dof%y, this%bc%msk, &
307 this%u%size(), n_pts)
308 call masked_gather_copy_0(this%r3%x, this%coef%dof%z, this%bc%msk, &
309 this%u%size(), n_pts)
310
311 call mpi_allreduce(n_pts, glb_n_pts, 1, &
312 mpi_integer, mpi_sum, neko_comm, ierr)
313 ! Print some information
314 call neko_log%section('Force/torque calculation')
315 write(log_buf, '(A,I4,A,A)') 'Zone ', zone_id, ' ', trim(zone_name)
316 call neko_log%message(log_buf)
317 write(log_buf, '(A,I6, I6)') 'Global number of GLL points in zone: ', &
318 glb_n_pts
319 call neko_log%message(log_buf)
320 write(log_buf, '(A,E15.7,E15.7,E15.7)') 'Average of zone''s coordinates: ',&
321 glsum(this%r1%x, n_pts)/glb_n_pts, &
322 glsum(this%r2%x, n_pts)/glb_n_pts, &
323 glsum(this%r3%x, n_pts)/glb_n_pts
324 call neko_log%message(log_buf)
325 write(log_buf, '(A,E15.7,E15.7,E15.7)') 'Center for torque calculation: ', &
326 center
327 call neko_log%message(log_buf)
328 write(log_buf, '(A,E15.7)') 'Scale: ', scale
329 call neko_log%message(log_buf)
330 call neko_log%end_section()
331
332 call cadd(this%r1%x, -center(1), n_pts)
333 call cadd(this%r2%x, -center(2), n_pts)
334 call cadd(this%r3%x, -center(3), n_pts)
335 if (neko_bcknd_device .eq. 1 .and. n_pts .gt. 0) then
336 call device_memcpy(this%n1%x, this%n1%x_d, n_pts, host_to_device, &
337 .false.)
338 call device_memcpy(this%n2%x, this%n2%x_d, n_pts, host_to_device, &
339 .false.)
340 call device_memcpy(this%n3%x, this%n3%x_d, n_pts, host_to_device, &
341 .true.)
342 call device_memcpy(this%r1%x, this%r1%x_d, n_pts, host_to_device, &
343 .false.)
344 call device_memcpy(this%r2%x, this%r2%x_d, n_pts, host_to_device, &
345 .false.)
346 call device_memcpy(this%r3%x, this%r3%x_d, n_pts, host_to_device, &
347 .true.)
348 end if
349
350 end subroutine force_torque_init_common
351
353 subroutine force_torque_free(this)
354 class(force_torque_t), intent(inout) :: this
355 call this%free_base()
356
357 call this%n1%free()
358 call this%n2%free()
359 call this%n3%free()
360
361 call this%r1%free()
362 call this%r2%free()
363 call this%r3%free()
364
365 call this%force1%free()
366 call this%force2%free()
367 call this%force3%free()
368
369 call this%force4%free()
370 call this%force5%free()
371 call this%force6%free()
372
373 call this%pmsk%free()
374 call this%mu_msk%free()
375 call this%s11msk%free()
376 call this%s22msk%free()
377 call this%s33msk%free()
378 call this%s12msk%free()
379 call this%s13msk%free()
380 call this%s23msk%free()
381
382 nullify(this%u)
383 nullify(this%v)
384 nullify(this%w)
385 nullify(this%p)
386 nullify(this%coef)
387 nullify(this%mu)
388 end subroutine force_torque_free
389
392 subroutine force_torque_compute(this, time)
393 class(force_torque_t), intent(inout) :: this
394 type(time_state_t), intent(in) :: time
395
396 real(kind=rp) :: dgtq(12) = 0.0_rp
397 integer :: n_pts, temp_indices(6)
398 type(field_t), pointer :: s11, s22, s33, s12, s13, s23
399 character(len=1000) :: log_buf
400
401 n_pts = this%bc%msk(0)
402
403 call neko_scratch_registry%request_field(s11, temp_indices(1), .false.)
404 call neko_scratch_registry%request_field(s12, temp_indices(2), .false.)
405 call neko_scratch_registry%request_field(s13, temp_indices(3), .false.)
406 call neko_scratch_registry%request_field(s22, temp_indices(4), .false.)
407 call neko_scratch_registry%request_field(s23, temp_indices(5), .false.)
408 call neko_scratch_registry%request_field(s33, temp_indices(6), .false.)
409
410 call strain_rate(s11%x, s22%x, s33%x, s12%x, &
411 s13%x, s23%x, this%u, this%v, this%w, this%coef)
412
413 ! On the CPU we can actually just use the original subroutines...
414 if (neko_bcknd_device .eq. 0) then
415 call masked_gather_copy_0(this%s11msk%x, s11%x, this%bc%msk, &
416 this%u%size(), n_pts)
417 call masked_gather_copy_0(this%s22msk%x, s22%x, this%bc%msk, &
418 this%u%size(), n_pts)
419 call masked_gather_copy_0(this%s33msk%x, s33%x, this%bc%msk, &
420 this%u%size(), n_pts)
421 call masked_gather_copy_0(this%s12msk%x, s12%x, this%bc%msk, &
422 this%u%size(), n_pts)
423 call masked_gather_copy_0(this%s13msk%x, s13%x, this%bc%msk, &
424 this%u%size(), n_pts)
425 call masked_gather_copy_0(this%s23msk%x, s23%x, this%bc%msk, &
426 this%u%size(), n_pts)
427 call masked_gather_copy_0(this%pmsk%x, this%p%x, this%bc%msk, &
428 this%u%size(), n_pts)
429 call masked_gather_copy_0(this%mu_msk%x, this%mu%x, this%bc%msk, &
430 this%u%size(), n_pts)
431 call calc_force_array(this%force1%x, this%force2%x, this%force3%x, &
432 this%force4%x, this%force5%x, this%force6%x, &
433 this%s11msk%x, &
434 this%s22msk%x, &
435 this%s33msk%x, &
436 this%s12msk%x, &
437 this%s13msk%x, &
438 this%s23msk%x, &
439 this%pmsk%x, &
440 this%n1%x, &
441 this%n2%x, &
442 this%n3%x, &
443 this%mu_msk%x, &
444 n_pts)
445 dgtq(1) = glsum(this%force1%x, n_pts)
446 dgtq(2) = glsum(this%force2%x, n_pts)
447 dgtq(3) = glsum(this%force3%x, n_pts)
448 dgtq(4) = glsum(this%force4%x, n_pts)
449 dgtq(5) = glsum(this%force5%x, n_pts)
450 dgtq(6) = glsum(this%force6%x, n_pts)
451 call vcross(this%s11msk%x, this%s22msk%x, this%s33msk%x, &
452 this%r1%x, this%r2%x, this%r3%x, &
453 this%force1%x, this%force2%x, this%force3%x, n_pts)
454
455 dgtq(7) = glsum(this%s11msk%x, n_pts)
456 dgtq(8) = glsum(this%s22msk%x, n_pts)
457 dgtq(9) = glsum(this%s33msk%x, n_pts)
458 call vcross(this%s11msk%x, this%s22msk%x, this%s33msk%x, &
459 this%r1%x, this%r2%x, this%r3%x, &
460 this%force4%x, this%force5%x, this%force6%x, n_pts)
461 dgtq(10) = glsum(this%s11msk%x, n_pts)
462 dgtq(11) = glsum(this%s22msk%x, n_pts)
463 dgtq(12) = glsum(this%s33msk%x, n_pts)
464 else
465 if (n_pts .gt. 0) then
466 call device_masked_gather_copy_0(this%s11msk%x_d, s11%x_d, &
467 this%bc%msk_d, this%u%size(), n_pts)
468 call device_masked_gather_copy_0(this%s22msk%x_d, s22%x_d, &
469 this%bc%msk_d, this%u%size(), n_pts)
470 call device_masked_gather_copy_0(this%s33msk%x_d, s33%x_d, &
471 this%bc%msk_d, this%u%size(), n_pts)
472 call device_masked_gather_copy_0(this%s12msk%x_d, s12%x_d, &
473 this%bc%msk_d, this%u%size(), n_pts)
474 call device_masked_gather_copy_0(this%s13msk%x_d, s13%x_d, &
475 this%bc%msk_d, this%u%size(), n_pts)
476 call device_masked_gather_copy_0(this%s23msk%x_d, s23%x_d, &
477 this%bc%msk_d, this%u%size(), n_pts)
478 call device_masked_gather_copy_0(this%pmsk%x_d, this%p%x_d, &
479 this%bc%msk_d, this%u%size(), n_pts)
480 call device_masked_gather_copy_0(this%mu_msk%x_d, this%mu%x_d, &
481 this%bc%msk_d, this%u%size(), n_pts)
482
483 call device_calc_force_array(this%force1%x_d, this%force2%x_d, &
484 this%force3%x_d, &
485 this%force4%x_d, &
486 this%force5%x_d, &
487 this%force6%x_d, &
488 this%s11msk%x_d, &
489 this%s22msk%x_d, &
490 this%s33msk%x_d, &
491 this%s12msk%x_d, &
492 this%s13msk%x_d, &
493 this%s23msk%x_d, &
494 this%pmsk%x_d, &
495 this%n1%x_d, &
496 this%n2%x_d, &
497 this%n3%x_d, &
498 this%mu%x_d, &
499 n_pts)
500 !Overwriting masked s11, s22, s33 as they are no longer needed
501 call device_vcross(this%s11msk%x_d, this%s22msk%x_d, &
502 this%s33msk%x_d, &
503 this%r1%x_d, this%r2%x_d, this%r3%x_d, &
504 this%force1%x_d, this%force2%x_d, &
505 this%force3%x_d, n_pts)
506 call device_vcross(this%s12msk%x_d, this%s13msk%x_d, this%s23msk%x_d,&
507 this%r1%x_d, this%r2%x_d, this%r3%x_d, &
508 this%force4%x_d, this%force5%x_d, this%force6%x_d, n_pts)
509 end if
510 dgtq(1) = device_glsum(this%force1%x_d, n_pts)
511 dgtq(2) = device_glsum(this%force2%x_d, n_pts)
512 dgtq(3) = device_glsum(this%force3%x_d, n_pts)
513 dgtq(4) = device_glsum(this%force4%x_d, n_pts)
514 dgtq(5) = device_glsum(this%force5%x_d, n_pts)
515 dgtq(6) = device_glsum(this%force6%x_d, n_pts)
516 dgtq(7) = device_glsum(this%s11msk%x_d, n_pts)
517 dgtq(8) = device_glsum(this%s22msk%x_d, n_pts)
518 dgtq(9) = device_glsum(this%s33msk%x_d, n_pts)
519 dgtq(10) = device_glsum(this%s12msk%x_d, n_pts)
520 dgtq(11) = device_glsum(this%s13msk%x_d, n_pts)
521 dgtq(12) = device_glsum(this%s23msk%x_d, n_pts)
522 end if
523 dgtq = this%scale*dgtq
524 write(log_buf, '(A, I4, A, A)') 'Force and torque on zone ', &
525 this%zone_id, ' ', this%zone_name
526 call neko_log%message(log_buf)
527 write(log_buf, '(A)') &
528 'Time step, time, total force/torque, pressure, viscous, direction'
529 call neko_log%message(log_buf)
530 write(log_buf, this%print_format) &
531 time%tstep, time%t, dgtq(1) + dgtq(4), dgtq(1), dgtq(4), ', forcex'
532 call neko_log%message(log_buf)
533 write(log_buf, this%print_format) &
534 time%tstep, time%t, dgtq(2) + dgtq(5), dgtq(2), dgtq(5), ', forcey'
535 call neko_log%message(log_buf)
536 write(log_buf, this%print_format) &
537 time%tstep, time%t, dgtq(3) + dgtq(6), dgtq(3), dgtq(6), ', forcez'
538 call neko_log%message(log_buf)
539 write(log_buf, this%print_format) &
540 time%tstep, time%t, dgtq(7) + dgtq(10), dgtq(7), dgtq(10), ', torquex'
541 call neko_log%message(log_buf)
542 write(log_buf, this%print_format) &
543 time%tstep, time%t, dgtq(8) + dgtq(11), dgtq(8), dgtq(11), ', torquey'
544 call neko_log%message(log_buf)
545 write(log_buf, this%print_format) &
546 time%tstep, time%t, dgtq(9) + dgtq(12), dgtq(9), dgtq(12), ', torquez'
547 call neko_log%message(log_buf)
548 call neko_scratch_registry%relinquish_field(temp_indices)
549
550 end subroutine force_torque_compute
551
552end 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.
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_init_common(this, name, fluid_name, zone_id, zone_name, center, scale, coef, long_print)
Common part of constructors.
subroutine force_torque_free(this)
Destructor.
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_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:77
integer, parameter, public log_size
Definition log.f90:46
Definition math.f90:60
subroutine, public cadd(a, s, n)
Add a scalar to vector .
Definition math.f90:462
real(kind=rp) function, public glsum(a, n)
Sum a vector of length n.
Definition math.f90:499
subroutine, public masked_gather_copy_0(a, b, mask, n, n_mask)
Gather a masked vector to reduced contigous vector .
Definition math.f90:312
subroutine, public vcross(u1, u2, u3, v1, v2, v3, w1, w2, w3, n)
Compute a cross product assuming vector components etc.
Definition math.f90:653
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:149
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.
Defines a vector.
Definition vector.f90:34
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:56
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.