39 use json_module,
only : json_file
60 use mpi_f08,
only : mpi_integer, mpi_sum, mpi_allreduce
88 type(
vector_t) :: s11msk, s22msk, s33msk, s12msk, s13msk, s23msk
89 real(kind=
rp) :: center(3) = 0.0_rp
90 real(kind=
rp) :: scale
92 character(len=20) :: zone_name
93 type(
coef_t),
pointer :: coef => null()
95 character(len=80) :: print_format
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 => &
123 type(json_file),
intent(inout) :: json
124 class(
case_t),
intent(inout),
target :: case
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
133 call this%init_base(json,
case)
141 call this%init_common(name, fluid_name, zone_id, zone_name, center, scale, &
142 case%fluid%c_xh, long_print)
160 preprocess_controller, compute_controller, output_controller, &
161 fluid_name, zone_id, zone_name, center, scale, coef, long_print)
163 character(len=*),
intent(in) :: name
164 class(
case_t),
intent(inout),
target :: case
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
177 call this%init_base_from_components(
case, order, preprocess_controller, &
179 call this%init_common(name, fluid_name, zone_id, zone_name, center, scale, &
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)
207 character(len=*),
intent(in) :: name
208 class(
case_t),
intent(inout),
target :: case
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
224 call this%init_base_from_components(
case, order, preprocess_control, &
225 preprocess_value, compute_control, compute_value, output_control, &
227 call this%init_common(name, fluid_name, zone_id, zone_name, center, scale, &
242 zone_name, center, scale, coef, long_print)
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
257 this%zone_id = zone_id
260 this%zone_name = zone_name
262 if (long_print )
then
263 this%print_format =
'(I7,E20.10,E20.10,E20.10,E20.10,A)'
265 this%print_format =
'(I7,E13.5,E13.5,E13.5,E13.5,A)'
272 this%mu =>
neko_registry%get_field_by_name(fluid_name //
'_mu_tot')
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)
303 this%n1%x, this%n2%x, this%n3%x, n_pts)
305 this%u%size(), n_pts)
307 this%u%size(), n_pts)
309 this%u%size(), n_pts)
311 call mpi_allreduce(n_pts, glb_n_pts, 1, &
314 call neko_log%section(
'Force/torque calculation')
315 write(log_buf,
'(A,I4,A,A)')
'Zone ', zone_id,
' ', trim(zone_name)
317 write(log_buf,
'(A,I6, I6)')
'Global number of GLL points in zone: ', &
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
325 write(log_buf,
'(A,E15.7,E15.7,E15.7)')
'Center for torque calculation: ', &
328 write(log_buf,
'(A,E15.7)')
'Scale: ', scale
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)
355 call this%free_base()
365 call this%force1%free()
366 call this%force2%free()
367 call this%force3%free()
369 call this%force4%free()
370 call this%force5%free()
371 call this%force6%free()
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()
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
401 n_pts = this%bc%msk(0)
411 s13%x, s23%x, this%u, this%v, this%w, this%coef)
416 this%u%size(), n_pts)
418 this%u%size(), n_pts)
420 this%u%size(), n_pts)
422 this%u%size(), n_pts)
424 this%u%size(), n_pts)
426 this%u%size(), n_pts)
428 this%u%size(), n_pts)
430 this%u%size(), n_pts)
432 this%force4%x, this%force5%x, this%force6%x, &
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)
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)
465 if (n_pts .gt. 0)
then
467 this%bc%msk_d, this%u%size(), n_pts)
469 this%bc%msk_d, this%u%size(), n_pts)
471 this%bc%msk_d, this%u%size(), n_pts)
473 this%bc%msk_d, this%u%size(), n_pts)
475 this%bc%msk_d, this%u%size(), n_pts)
477 this%bc%msk_d, this%u%size(), n_pts)
479 this%bc%msk_d, this%u%size(), n_pts)
481 this%bc%msk_d, this%u%size(), n_pts)
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)
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
527 write(log_buf,
'(A)') &
528 'Time step, time, total force/torque, pressure, viscous, direction'
530 write(log_buf, this%print_format) &
531 time%tstep, time%t, dgtq(1) + dgtq(4), dgtq(1), dgtq(4),
', forcex'
533 write(log_buf, this%print_format) &
534 time%tstep, time%t, dgtq(2) + dgtq(5), dgtq(2), dgtq(5),
', forcey'
536 write(log_buf, this%print_format) &
537 time%tstep, time%t, dgtq(3) + dgtq(6), dgtq(3), dgtq(6),
', forcez'
539 write(log_buf, this%print_format) &
540 time%tstep, time%t, dgtq(7) + dgtq(10), dgtq(7), dgtq(10),
', torquex'
542 write(log_buf, this%print_format) &
543 time%tstep, time%t, dgtq(8) + dgtq(11), dgtq(8), dgtq(11),
', torquey'
545 write(log_buf, this%print_format) &
546 time%tstep, time%t, dgtq(9) + dgtq(12), dgtq(9), dgtq(12),
', torquez'
Copy data between host and device (or device and device)
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.
Defines a simulation case.
type(mpi_comm), public neko_comm
MPI communicator.
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.
integer, parameter, public host_to_device
Defines a dirichlet boundary condition.
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.
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.
type(log_t), public neko_log
Global log stream.
integer, parameter, public log_size
subroutine, public cadd(a, s, n)
Add a scalar to vector .
real(kind=rp) function, public glsum(a, n)
Sum a vector of length n.
subroutine, public masked_gather_copy_0(a, b, mask, n, n_mask)
Gather a masked vector to reduced contigous vector .
subroutine, public vcross(u1, u2, u3, v1, v2, v3, w1, w2, w3, n)
Compute a cross product assuming vector components etc.
integer, parameter neko_bcknd_device
integer, parameter, public dp
integer, parameter, public sp
integer, parameter, public rp
Global precision used in computations.
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.
type(registry_t), target, public neko_registry
Global field registry.
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.
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Generic Dirichlet boundary condition on .
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.