39 use json_module,
only : json_file
59 use mpi_f08,
only : mpi_integer, mpi_sum, mpi_allreduce
87 type(
vector_t) :: s11msk, s22msk, s33msk, s12msk, s13msk, s23msk
88 real(kind=
rp) :: center(3) = 0.0_rp
89 real(kind=
rp) :: scale
91 character(len=20) :: zone_name
92 type(
coef_t),
pointer :: coef => null()
94 character(len=80) :: print_format
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 => &
122 type(json_file),
intent(inout) :: json
123 class(
case_t),
intent(inout),
target :: case
125 real(kind=
rp),
allocatable :: center(:)
126 character(len=:),
allocatable :: zone_name, fluid_name
127 real(kind=
rp) :: scale
128 logical :: long_print
130 call this%init_base(json,
case)
133 call json_get(json,
'zone_id', zone_id)
137 call json_get(json,
'center', center)
139 call this%init_common(fluid_name, zone_id, zone_name, center, scale, &
140 case%fluid%c_xh, long_print)
157 preprocess_controller, compute_controller, output_controller, &
158 fluid_name, zone_id, zone_name, center, scale, coef, long_print)
160 class(
case_t),
intent(inout),
target :: case
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
173 call this%init_base_from_components(
case, order, preprocess_controller, &
175 call this%init_common(fluid_name, zone_id, zone_name, center, scale, &
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)
202 class(
case_t),
intent(inout),
target :: case
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
218 call this%init_base_from_components(
case, order, preprocess_control, &
219 preprocess_value, compute_control, compute_value, output_control, &
221 call this%init_common(fluid_name, zone_id, zone_name, center, scale, &
235 zone_name, center, scale, coef, long_print)
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
248 this%zone_id = zone_id
251 this%zone_name = zone_name
253 if (long_print )
then
254 this%print_format =
'(I7,E20.10,E20.10,E20.10,E20.10,A)'
256 this%print_format =
'(I7,E13.5,E13.5,E13.5,E13.5,A)'
263 this%mu =>
neko_registry%get_field_by_name(fluid_name //
'_mu_tot')
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)
294 this%n1%x, this%n2%x, this%n3%x, n_pts)
296 this%u%size(), n_pts)
298 this%u%size(), n_pts)
300 this%u%size(), n_pts)
302 call mpi_allreduce(n_pts, glb_n_pts, 1, &
305 call neko_log%section(
'Force/torque calculation')
306 write(log_buf,
'(A,I4,A,A)')
'Zone ', zone_id,
' ', trim(zone_name)
308 write(log_buf,
'(A,I6, I6)')
'Global number of GLL points in zone: ', &
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
316 write(log_buf,
'(A,E15.7,E15.7,E15.7)')
'Center for torque calculation: ', &
319 write(log_buf,
'(A,E15.7)')
'Scale: ', scale
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)
346 call this%free_base()
362 real(kind=
rp) :: dgtq(12) = 0.0_rp
363 integer :: n_pts, temp_indices(6)
364 type(
field_t),
pointer :: s11, s22, s33, s12, s13, s23
365 character(len=1000) :: log_buf
367 n_pts = this%bc%msk(0)
377 s13%x, s23%x, this%u, this%v, this%w, this%coef)
382 this%u%size(), n_pts)
384 this%u%size(), n_pts)
386 this%u%size(), n_pts)
388 this%u%size(), n_pts)
390 this%u%size(), n_pts)
392 this%u%size(), n_pts)
394 this%u%size(), n_pts)
396 this%u%size(), n_pts)
398 this%force4%x, this%force5%x, this%force6%x, &
411 dgtq(1) =
glsum(this%force1%x, n_pts)
412 dgtq(2) =
glsum(this%force2%x, n_pts)
413 dgtq(3) =
glsum(this%force3%x, n_pts)
414 dgtq(4) =
glsum(this%force4%x, n_pts)
415 dgtq(5) =
glsum(this%force5%x, n_pts)
416 dgtq(6) =
glsum(this%force6%x, n_pts)
417 call vcross(this%s11msk%x, this%s22msk%x, this%s33msk%x, &
418 this%r1%x, this%r2%x, this%r3%x, &
419 this%force1%x, this%force2%x, this%force3%x, n_pts)
421 dgtq(7) =
glsum(this%s11msk%x, n_pts)
422 dgtq(8) =
glsum(this%s22msk%x, n_pts)
423 dgtq(9) =
glsum(this%s33msk%x, n_pts)
424 call vcross(this%s11msk%x, this%s22msk%x, this%s33msk%x, &
425 this%r1%x, this%r2%x, this%r3%x, &
426 this%force4%x, this%force5%x, this%force6%x, n_pts)
427 dgtq(10) =
glsum(this%s11msk%x, n_pts)
428 dgtq(11) =
glsum(this%s22msk%x, n_pts)
429 dgtq(12) =
glsum(this%s33msk%x, n_pts)
431 if (n_pts .gt. 0)
then
433 this%bc%msk_d, this%u%size(), n_pts)
435 this%bc%msk_d, this%u%size(), n_pts)
437 this%bc%msk_d, this%u%size(), n_pts)
439 this%bc%msk_d, this%u%size(), n_pts)
441 this%bc%msk_d, this%u%size(), n_pts)
443 this%bc%msk_d, this%u%size(), n_pts)
445 this%bc%msk_d, this%u%size(), n_pts)
447 this%bc%msk_d, this%u%size(), n_pts)
469 this%r1%x_d, this%r2%x_d, this%r3%x_d, &
470 this%force1%x_d, this%force2%x_d, &
471 this%force3%x_d, n_pts)
472 call device_vcross(this%s12msk%x_d, this%s13msk%x_d, this%s23msk%x_d,&
473 this%r1%x_d, this%r2%x_d, this%r3%x_d, &
474 this%force4%x_d, this%force5%x_d, this%force6%x_d, n_pts)
489 dgtq = this%scale*dgtq
490 write(log_buf,
'(A, I4, A, A)')
'Force and torque on zone ', &
491 this%zone_id,
' ', this%zone_name
493 write(log_buf,
'(A)') &
494 'Time step, time, total force/torque, pressure, viscous, direction'
496 write(log_buf, this%print_format) &
497 time%tstep, time%t, dgtq(1) + dgtq(4), dgtq(1), dgtq(4),
', forcex'
499 write(log_buf, this%print_format) &
500 time%tstep, time%t, dgtq(2) + dgtq(5), dgtq(2), dgtq(5),
', forcey'
502 write(log_buf, this%print_format) &
503 time%tstep, time%t, dgtq(3) + dgtq(6), dgtq(3), dgtq(6),
', forcez'
505 write(log_buf, this%print_format) &
506 time%tstep, time%t, dgtq(7) + dgtq(10), dgtq(7), dgtq(10),
', torquex'
508 write(log_buf, this%print_format) &
509 time%tstep, time%t, dgtq(8) + dgtq(11), dgtq(8), dgtq(11),
', torquey'
511 write(log_buf, this%print_format) &
512 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_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.
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.