47 use json_module,
only : json_file, json_value, json_core
57 use,
intrinsic :: iso_c_binding
61 use mpi_f08,
only : mpi_allreduce, mpi_integer, mpi_sum, &
62 mpi_double_precision, mpi_gatherv, mpi_gather, mpi_exscan
68 real(kind=
rp) :: start_time
70 integer :: n_fields = 0
73 integer :: n_global_probes
75 integer :: n_probes_offset
77 real(kind=
rp),
allocatable :: xyz(:,:)
79 real(kind=
rp),
allocatable :: out_values(:,:)
80 type(c_ptr),
allocatable :: out_values_d(:)
81 real(kind=
rp),
allocatable :: out_vals_trsp(:,:)
83 integer :: n_local_probes
86 character(len=20),
allocatable :: which_fields(:)
88 integer,
allocatable :: n_local_probes_tot_offset(:)
89 integer,
allocatable :: n_local_probes_tot(:)
92 real(kind=
rp),
allocatable :: global_output_values(:,:)
100 procedure, pass(this) :: init_from_components => &
131 class(
probes_t),
intent(inout),
target :: this
132 type(json_file),
intent(inout) :: json
133 class(
case_t),
intent(inout),
target :: case
134 character(len=:),
allocatable :: output_file
135 character(len=:),
allocatable :: input_file
139 character(len=:),
allocatable :: point_type
140 type(json_value),
pointer :: json_point_list
141 type(json_file) :: json_point
142 type(json_core) :: core
143 integer :: idx, n_point_children
145 character(len=:),
allocatable :: name
151 call this%init_base(json,
case)
154 call json%info(
'fields', n_children = this%n_fields)
155 call json_get(json,
'fields', this%which_fields)
156 call json_get(json,
'output_file', output_file)
160 call this%sampled_fields%init(this%n_fields)
161 do i = 1, this%n_fields
163 call this%sampled_fields%assign(i, &
168 this%n_local_probes = 0
169 this%n_global_probes = 0
172 if (json%valid_path(
'points_file'))
then
175 call json_get(json,
'points_file', input_file)
180 this%n_global_probes, input_file)
184 call json%get(
'points', json_point_list)
185 call json%info(
'points', n_children = n_point_children)
187 do idx = 1, n_point_children
191 select case (point_type)
194 call this%read_file(json_point)
196 call this%read_point(json_point)
198 call this%read_line(json_point)
200 call neko_error(
'Plane probes not implemented yet.')
202 call this%read_circle(json_point)
204 call this%read_point_zone(json_point,
case%fluid%dm_Xh)
206 call json_point%print()
209 call neko_error(
'Unknown region type ' // point_type)
213 call mpi_allreduce(this%n_local_probes, this%n_global_probes, 1, &
217 call this%init_from_components(
case%fluid%dm_Xh, output_file, name)
229 class(
probes_t),
intent(inout) :: this
230 type(json_file),
intent(inout) :: json
231 character(len=:),
allocatable :: input_file
232 real(kind=
rp),
dimension(:,:),
allocatable :: point_list
234 integer :: n_local, n_global
238 call json_get(json,
'file_name', input_file)
242 call this%add_points(point_list)
250 class(
probes_t),
intent(inout) :: this
251 type(json_file),
intent(inout) :: json
253 real(kind=
rp),
dimension(:,:),
allocatable :: point_list
254 real(kind=
rp),
dimension(:),
allocatable :: rp_list_reader
260 if (mod(
size(rp_list_reader), 3) /= 0)
then
261 call neko_error(
'Invalid number of coordinates.')
265 allocate(point_list(3,
size(rp_list_reader)/3))
266 point_list = reshape(rp_list_reader, [3,
size(rp_list_reader)/3])
268 call this%add_points(point_list)
275 class(
probes_t),
intent(inout) :: this
276 type(json_file),
intent(inout) :: json
278 real(kind=
rp),
dimension(:,:),
allocatable :: point_list
279 real(kind=
rp),
dimension(:),
allocatable :: start,
end
280 real(kind=
rp),
dimension(3) :: direction
283 integer :: n_points, i
289 call json_get_or_lookup(json,
"amount", n_points)
292 if (
size(start) /= 3 .or.
size(
end) /= 3) then
293 call neko_error(
'Invalid start or end coordinates.')
297 allocate(point_list(3, n_points))
300 direction =
end - start
302 t =
real(i - 1, kind = rp) /
real(n_points - 1, kind = rp)
303 point_list(:, i) = start + direction * t
306 call this%add_points(point_list)
319 class(
probes_t),
intent(inout) :: this
320 type(json_file),
intent(inout) :: json
322 real(kind=rp),
dimension(:,:),
allocatable :: point_list
323 real(kind=rp),
dimension(:),
allocatable :: center, normal
324 real(kind=rp) :: radius
325 real(kind=rp) :: angle
326 integer :: n_points, i
327 character(len=:),
allocatable :: axis
329 real(kind=rp),
dimension(3) :: zero_line, cross_line, temp
333 if (pe_rank .ne. 0)
return
334 call json_get_or_lookup(json,
"center", center)
335 call json_get_or_lookup(json,
"normal", normal)
336 call json_get_or_lookup(json,
"radius", radius)
337 call json_get_or_lookup(json,
"amount", n_points)
338 call json_get(json,
"axis", axis)
341 if (
size(center) /= 3 .or.
size(normal) /= 3)
then
342 call neko_error(
'Invalid center or normal coordinates.')
344 if (axis /=
'x' .and. axis /=
'y' .and. axis /=
'z')
then
345 call neko_error(
'Invalid axis.')
347 if (radius <= 0)
then
348 call neko_error(
'Invalid radius.')
350 if (n_points <= 0)
then
351 call neko_error(
'Invalid number of points.')
355 normal = normal / norm2(normal)
358 if (axis .eq.
'x') zero_line = [1.0, 0.0, 0.0]
359 if (axis .eq.
'y') zero_line = [0.0, 1.0, 0.0]
360 if (axis .eq.
'z') zero_line = [0.0, 0.0, 1.0]
362 if (1.0_rp - dot_product(zero_line, normal) .le. 1e-6)
then
363 call neko_error(
'Invalid axis and normal.')
366 zero_line = zero_line - dot_product(zero_line, normal) * normal
367 zero_line = zero_line / norm2(zero_line)
369 cross_line(1) = normal(2) * zero_line(3) - normal(3) * zero_line(2)
370 cross_line(2) = normal(3) * zero_line(1) - normal(1) * zero_line(3)
371 cross_line(3) = normal(1) * zero_line(2) - normal(2) * zero_line(1)
374 allocate(point_list(3, n_points))
376 pi = 4.0_rp * atan(1.0_rp)
380 angle = 2.0_rp * pi *
real(i - 1, kind = rp) /
real(n_points, kind = rp)
381 temp = cos(angle) * zero_line + sin(angle) * cross_line
383 point_list(:, i) = center + radius * temp
386 call this%add_points(point_list)
395 class(
probes_t),
intent(inout) :: this
396 type(json_file),
intent(inout) :: json
397 type(dofmap_t),
intent(in) :: dof
399 real(kind=rp),
dimension(:,:),
allocatable :: point_list
400 character(len=:),
allocatable :: point_zone_name
401 class(point_zone_t),
pointer :: zone
402 integer :: i, idx, lx, nlindex(4)
403 real(kind=rp) :: x, y, z
406 if (pe_rank .ne. 0)
return
408 call json_get(json,
"name", point_zone_name)
409 zone => neko_point_zone_registry%get_point_zone(point_zone_name)
412 allocate(point_list(3, zone%size))
416 idx = zone%mask%get(i)
419 x = dof%x(nlindex(1), nlindex(2), nlindex(3), nlindex(4))
420 y = dof%y(nlindex(1), nlindex(2), nlindex(3), nlindex(4))
421 z = dof%z(nlindex(1), nlindex(2), nlindex(3), nlindex(4))
423 point_list(:, i) = [x, y, z]
426 call this%add_points(point_list)
436 class(
probes_t),
intent(inout) :: this
437 real(kind=rp),
dimension(:,:),
intent(in) :: new_points
439 real(kind=rp),
dimension(:,:),
allocatable :: temp
440 integer :: n_old, n_new
443 n_old = this%n_local_probes
444 n_new =
size(new_points, 2)
447 if (
allocated(this%xyz))
then
448 call move_alloc(this%xyz, temp)
452 allocate(this%xyz(3, n_old + n_new))
453 if (
allocated(temp))
then
454 this%xyz(:, 1:n_old) = temp
456 this%xyz(:, n_old+1:n_old+n_new) = new_points
458 this%n_local_probes = this%n_local_probes + n_new
468 class(
probes_t),
intent(inout) :: this
469 type(dofmap_t),
intent(in) :: dof
470 character(len=:),
allocatable,
intent(inout) :: output_file
471 character(len=*),
intent(in) :: name
473 character(len=1024) :: header_line
474 real(kind=rp),
allocatable :: global_output_coords(:,:)
476 type(matrix_t) :: mat_coords
481 call this%global_interp%init(dof)
484 call this%global_interp%find_points_and_redist(this%xyz, &
488 allocate(this%out_values(this%n_local_probes, this%n_fields))
489 allocate(this%out_values_d(this%n_fields))
490 allocate(this%out_vals_trsp(this%n_fields, this%n_local_probes))
492 if (neko_bcknd_device .eq. 1)
then
493 do i = 1, this%n_fields
494 this%out_values_d(i) = c_null_ptr
495 call device_map(this%out_values(:,i), this%out_values_d(i), &
501 call this%fout%init(trim(output_file))
503 select type (ft => this%fout%file_type)
509 write(header_line,
'(I0,A,I0)') this%n_global_probes,
",", this%n_fields
510 do i = 1, this%n_fields
511 header_line = trim(header_line) //
"," // trim(this%which_fields(i))
513 call this%fout%set_header(header_line)
518 allocate(this%n_local_probes_tot(pe_size))
519 allocate(this%n_local_probes_tot_offset(pe_size))
520 call this%setup_offset()
521 if (pe_rank .eq. 0)
then
522 allocate(global_output_coords(3, this%n_global_probes))
523 call this%mat_out%init(this%n_global_probes, this%n_fields)
524 allocate(this%global_output_values(this%n_fields, &
525 this%n_global_probes))
526 call mat_coords%init(this%n_global_probes,3)
528 call mpi_gatherv(this%xyz, 3*this%n_local_probes, &
529 mpi_real_precision, global_output_coords, &
530 3*this%n_local_probes_tot, &
531 3*this%n_local_probes_tot_offset, &
532 mpi_real_precision, 0, neko_comm, ierr)
533 if (pe_rank .eq. 0)
then
534 call trsp(mat_coords%x, this%n_global_probes, &
535 global_output_coords, 3)
537 call this%fout%write(mat_coords)
540 call neko_error(
"Invalid data. Expected csv_file_t.")
547 class(
probes_t),
intent(inout) :: this
550 if (
allocated(this%xyz))
then
554 if (
allocated(this%out_values))
then
555 deallocate(this%out_values)
558 if (
allocated(this%out_vals_trsp))
then
559 deallocate(this%out_vals_trsp)
562 call this%sampled_fields%free()
564 if (
allocated(this%n_local_probes_tot))
then
565 deallocate(this%n_local_probes_tot)
568 if (
allocated(this%n_local_probes_tot_offset))
then
569 deallocate(this%n_local_probes_tot_offset)
572 if (
allocated(this%global_output_values))
then
573 deallocate(this%global_output_values)
576 if (
allocated(this%which_fields))
then
577 deallocate(this%which_fields)
580 if (
allocated(this%out_values_d))
then
581 do i = 1,
size(this%out_values_d)
582 if (c_associated(this%out_values_d(i)))
then
583 call device_free(this%out_values_d(i))
586 deallocate(this%out_values_d)
589 call this%global_interp%free()
590 call this%mat_out%free()
597 character(len=LOG_SIZE) :: log_buf
601 call neko_log%section(
'Probes')
602 write(log_buf,
'(A,I6)')
"Number of probes: ", this%n_global_probes
603 call neko_log%message(log_buf)
605 call neko_log%message(
"xyz-coordinates:", lvl = neko_log_debug)
606 do i = 1, this%n_local_probes
607 write(log_buf,
'("(",F10.6,",",F10.6,",",F10.6,")")') this%xyz(:,i)
608 call neko_log%message(log_buf, lvl = neko_log_debug)
612 write(log_buf,
'(A,I6)')
"Number of fields: ", this%n_fields
613 call neko_log%message(log_buf)
614 do i = 1, this%n_fields
615 write(log_buf,
'(A,I6,A,A)') &
616 "Field: ", i,
" ", trim(this%which_fields(i))
617 call neko_log%message(log_buf, lvl = neko_log_debug)
619 call neko_log%end_section()
620 call neko_log%newline()
628 character(len=LOG_SIZE) :: log_buf
631 do i = 1, this%n_local_probes
632 write (log_buf, *) pe_rank,
"/", this%global_interp%pe_owner(i), &
633 "/" , this%global_interp%el_owner0(i)
634 call neko_log%message(log_buf)
635 write(log_buf,
'(A5,"(",F10.6,",",F10.6,",",F10.6,")")') &
636 "rst: ", this%global_interp%rst(:,i)
637 call neko_log%message(log_buf)
645 this%n_local_probes_tot = 0
646 this%n_local_probes_tot_offset = 0
647 this%n_probes_offset = 0
648 call mpi_gather(this%n_local_probes, 1, mpi_integer, &
649 this%n_local_probes_tot, 1, mpi_integer, &
652 call mpi_exscan(this%n_local_probes, this%n_probes_offset, 1, &
653 mpi_integer, mpi_sum, neko_comm, ierr)
654 call mpi_gather(this%n_probes_offset, 1, mpi_integer, &
655 this%n_local_probes_tot_offset, 1, mpi_integer, &
667 class(
probes_t),
intent(inout) :: this
668 type(time_state_t),
intent(in) :: time
670 logical :: do_interp_on_host = .false.
673 if (time%t .lt. this%start_time)
return
676 do i = 1, this%n_fields
677 call this%global_interp%evaluate(this%out_values(:,i), &
678 this%sampled_fields%items(i)%ptr%x, &
682 if (neko_bcknd_device .eq. 1)
then
683 do i = 1, this%n_fields
684 call device_memcpy(this%out_values(:,i), this%out_values_d(i), &
685 this%n_local_probes, device_to_host, sync = .true.)
689 if (this%output_controller%check(time))
then
692 if (this%seq_io)
then
693 call trsp(this%out_vals_trsp, this%n_fields, &
694 this%out_values, this%n_local_probes)
695 call mpi_gatherv(this%out_vals_trsp, &
696 this%n_fields*this%n_local_probes, &
697 mpi_real_precision, this%global_output_values, &
698 this%n_fields*this%n_local_probes_tot, &
699 this%n_fields*this%n_local_probes_tot_offset, &
700 mpi_real_precision, 0, neko_comm, ierr)
701 if (pe_rank .eq. 0)
then
702 call trsp(this%mat_out%x, this%n_global_probes, &
703 this%global_output_values, this%n_fields)
704 call this%fout%write(this%mat_out, time%t)
707 call neko_error(
'probes sim comp, parallel io need implementation')
711 call this%output_controller%register_execution()
720 class(
probes_t),
intent(inout) :: this
721 character(len=:),
allocatable :: points_file
722 real(kind=rp),
allocatable :: xyz(:,:)
723 integer,
intent(inout) :: n_local_probes, n_global_probes
726 type(file_t) :: file_in
728 call file_in%init(trim(points_file))
730 select type (ft => file_in%file_type)
735 call neko_error(
"Invalid data. Expected csv_file_t.")
739 call file_free(file_in)
749 type(csv_file_t),
intent(inout) :: f
750 real(kind=rp),
allocatable :: xyz(:,:)
751 integer,
intent(inout) :: n_local_probes, n_global_probes
752 type(matrix_t) :: mat_in, mat_in2
755 n_lines = f%count_lines()
758 n_global_probes = n_lines
761 if (pe_rank .eq. 0)
then
762 n_local_probes = n_global_probes
763 allocate(xyz(3, n_local_probes))
764 call mat_in%init(n_global_probes,3)
765 call mat_in2%init(3, n_global_probes)
767 call trsp(xyz, 3, mat_in%x, n_global_probes)
770 allocate(xyz(3, n_local_probes))
__inline__ __device__ void nonlinear_index(const int idx, const int lx, int *index)
Map a Fortran array to a device (allocate and associate)
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 simulation case.
type(mpi_datatype), public mpi_real_precision
MPI type for working precision of REAL types.
integer, public pe_size
MPI size of communicator.
integer, public pe_rank
MPI rank.
type(mpi_comm), public neko_comm
MPI communicator.
File format for .csv files, used for any read/write operations involving floating point data.
Device abstraction, common interface for various accelerators.
subroutine, public device_free(x_d)
Deallocate memory on the device.
integer, parameter, public device_to_host
Defines a mapping of the degrees of freedom.
Module for file I/O operations.
subroutine file_free(this)
File operation destructor.
Implements global_interpolation given a dofmap.
Utilities for retrieving parameters from the case files.
integer, parameter, public neko_log_debug
Debug log level.
type(log_t), public neko_log
Global log stream.
integer, parameter, public log_size
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
type(point_zone_registry_t), target, public neko_point_zone_registry
Global point_zone registry.
subroutine read_point_zone(this, json, dof)
Construct a list of points from a point zone.
subroutine read_circle(this, json)
Construct a list of points from a circle.
subroutine add_points(this, new_points)
Append a new list of points to the exsiting list.
subroutine probes_evaluate_and_write(this, time)
Interpolate each probe from its r,s,t coordinates.
subroutine read_line(this, json)
Construct a list of points from a line.
subroutine read_probe_locations(this, xyz, n_local_probes, n_global_probes, points_file)
Initialize the physical coordinates from a csv input file.
subroutine read_point(this, json)
Read a list of points from the json file.
subroutine probes_show(this)
Print current probe status, with number of probes and coordinates.
subroutine probes_free(this)
Destructor.
subroutine probes_init_from_components(this, dof, output_file, name)
Initialize without json things.
subroutine probes_debug(this)
Show the status of processor/element owner and error code for each point.
subroutine probes_setup_offset(this)
Setup offset for rank 0.
subroutine read_file(this, json)
Read a list of points from a csv file.
subroutine read_xyz_from_csv(xyz, n_local_probes, n_global_probes, f)
Read and initialize the number of probes from a csv input file.
subroutine probes_init_from_json(this, json, case)
Constructor from json.
Defines a registry for storing solution fields.
type(registry_t), target, public neko_registry
Global field registry.
Simulation components are objects that encapsulate functionality that can be fit to a particular comp...
subroutine compute_(this, time)
Dummy compute function.
subroutine, public trsp(a, lda, b, ldb)
Transpose of a rectangular tensor .
Module with things related to the simulation time.
field_list_t, To be able to group fields together
A wrapper around a polymorphic generic_file_t that handles its init. This is essentially a factory fo...
Implements global interpolation for arbitrary points in the domain.
Base abstract type for point zones.
Base abstract class for simulation components.
A struct that contains all info about the time, expand as needed.