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