48 use json_module,
only : json_file, json_value, json_core
62 use,
intrinsic :: iso_c_binding
66 use mpi_f08,
only : mpi_allreduce, mpi_integer, mpi_sum, &
67 mpi_double_precision, mpi_gatherv, mpi_gather, mpi_exscan
73 real(kind=
rp) :: start_time
75 integer :: n_fields = 0
78 integer :: n_global_probes
80 integer :: n_probes_offset
82 real(kind=
rp),
allocatable :: xyz(:,:)
84 real(kind=
rp),
allocatable :: out_values(:,:)
85 type(c_ptr),
allocatable :: out_values_d(:)
86 real(kind=
rp),
allocatable :: out_vals_trsp(:,:)
88 integer :: n_local_probes
91 character(len=NEKO_VARNAME_LEN),
allocatable :: which_fields(:)
93 integer,
allocatable :: n_local_probes_tot_offset(:)
94 integer,
allocatable :: n_local_probes_tot(:)
97 real(kind=
rp),
allocatable :: global_output_values(:,:)
102 logical :: append_out = .true.
107 procedure, pass(this) :: init_from_components => &
140 class(
probes_t),
intent(inout),
target :: this
141 type(json_file),
intent(inout) :: json
142 class(
case_t),
intent(inout),
target :: case
143 character(len=:),
allocatable :: output_file
144 character(len=:),
allocatable :: input_file
148 character(len=:),
allocatable :: point_type
149 type(json_value),
pointer :: json_point_list
150 type(json_file) :: json_point
151 type(json_core) :: core
152 integer :: idx, n_point_children
154 character(len=:),
allocatable :: name
160 call this%init_base(json,
case)
163 call json%info(
'fields', n_children = this%n_fields)
164 call json_get(json,
'fields', this%which_fields)
165 call json_get(json,
'output_file', output_file)
170 call this%sampled_fields%init(this%n_fields)
171 do i = 1, this%n_fields
173 call this%sampled_fields%assign(i, &
178 this%n_local_probes = 0
179 this%n_global_probes = 0
182 if (json%valid_path(
'points_file'))
then
185 call json_get(json,
'points_file', input_file)
190 this%n_global_probes, input_file)
193 if (json%valid_path(
'points'))
then
196 call json%get(
'points', json_point_list)
197 call json%info(
'points', n_children = n_point_children)
199 do idx = 1, n_point_children
203 select case (point_type)
206 call this%read_file(json_point)
208 call this%read_point(json_point)
210 call this%read_line(json_point)
212 call neko_error(
'Plane probes not implemented yet.')
214 call this%read_circle(json_point)
216 call this%read_point_zone(json_point,
case%fluid%dm_Xh)
218 call json_point%print()
221 call neko_error(
'Unknown region type ' // point_type)
227 call mpi_allreduce(this%n_local_probes, this%n_global_probes, 1, &
234 type(json_file) :: interp_subdict
236 call this%global_interp%init(
case%fluid%dm_Xh, &
237 params_subdict = interp_subdict)
239 call this%init_common(
case%fluid%dm_Xh, output_file, name)
253 class(
probes_t),
intent(inout) :: this
254 type(dofmap_t),
intent(in) :: dof
255 character(len=:),
allocatable,
intent(inout) :: output_file
256 character(len=*),
intent(in) :: name
257 real(kind=dp),
intent(in),
optional :: tolerance, padding
259 character(len=1024) :: header_line
260 real(kind=rp),
allocatable :: global_output_coords(:,:)
262 type(matrix_t) :: mat_coords
264 call this%global_interp%init(dof, tol = tolerance, pad = padding)
266 call this%init_common(dof, output_file, name)
278 class(
probes_t),
intent(inout) :: this
279 type(json_file),
intent(inout) :: json
280 character(len=:),
allocatable :: input_file
281 real(kind=rp),
dimension(:,:),
allocatable :: point_list
283 integer :: n_local, n_global
285 if (pe_rank .ne. 0)
return
287 call json_get(json,
'file_name', input_file)
291 call this%add_points(point_list)
299 class(
probes_t),
intent(inout) :: this
300 type(json_file),
intent(inout) :: json
302 real(kind=rp),
dimension(:,:),
allocatable :: point_list
303 real(kind=rp),
dimension(:),
allocatable :: rp_list_reader
306 if (pe_rank .ne. 0)
return
307 call json_get_or_lookup(json,
'coordinates', rp_list_reader)
309 if (mod(
size(rp_list_reader), 3) /= 0)
then
310 call neko_error(
'Invalid number of coordinates.')
314 allocate(point_list(3,
size(rp_list_reader)/3))
315 point_list = reshape(rp_list_reader, [3,
size(rp_list_reader)/3])
317 call this%add_points(point_list)
324 class(
probes_t),
intent(inout) :: this
325 type(json_file),
intent(inout) :: json
327 real(kind=rp),
dimension(:,:),
allocatable :: point_list
328 real(kind=rp),
dimension(:),
allocatable :: start,
end
329 real(kind=rp),
dimension(3) :: direction
332 integer :: n_points, i
335 if (pe_rank .ne. 0)
return
336 call json_get_or_lookup(json,
"start", start)
337 call json_get_or_lookup(json,
"end",
end)
338 call json_get_or_lookup(json,
"amount", n_points)
341 if (
size(start) /= 3 .or.
size(
end) /= 3) then
342 call neko_error(
'Invalid start or end coordinates.')
346 allocate(point_list(3, n_points))
349 direction =
end - start
351 t =
real(i - 1, kind = rp) /
real(n_points - 1, kind = rp)
352 point_list(:, i) = start + direction * t
355 call this%add_points(point_list)
368 class(
probes_t),
intent(inout) :: this
369 type(json_file),
intent(inout) :: json
371 real(kind=rp),
dimension(:,:),
allocatable :: point_list
372 real(kind=rp),
dimension(:),
allocatable :: center, normal
373 real(kind=rp) :: radius
374 real(kind=rp) :: angle
375 integer :: n_points, i
376 character(len=:),
allocatable :: axis
378 real(kind=rp),
dimension(3) :: zero_line, cross_line, temp
382 if (pe_rank .ne. 0)
return
383 call json_get_or_lookup(json,
"center", center)
384 call json_get_or_lookup(json,
"normal", normal)
385 call json_get_or_lookup(json,
"radius", radius)
386 call json_get_or_lookup(json,
"amount", n_points)
387 call json_get(json,
"axis", axis)
390 if (
size(center) /= 3 .or.
size(normal) /= 3)
then
391 call neko_error(
'Invalid center or normal coordinates.')
393 if (axis /=
'x' .and. axis /=
'y' .and. axis /=
'z')
then
394 call neko_error(
'Invalid axis.')
396 if (radius <= 0)
then
397 call neko_error(
'Invalid radius.')
399 if (n_points <= 0)
then
400 call neko_error(
'Invalid number of points.')
404 normal = normal / norm2(normal)
407 if (axis .eq.
'x') zero_line = [1.0, 0.0, 0.0]
408 if (axis .eq.
'y') zero_line = [0.0, 1.0, 0.0]
409 if (axis .eq.
'z') zero_line = [0.0, 0.0, 1.0]
411 if (1.0_rp - dot_product(zero_line, normal) .le. 1e-6)
then
412 call neko_error(
'Invalid axis and normal.')
415 zero_line = zero_line - dot_product(zero_line, normal) * normal
416 zero_line = zero_line / norm2(zero_line)
418 cross_line(1) = normal(2) * zero_line(3) - normal(3) * zero_line(2)
419 cross_line(2) = normal(3) * zero_line(1) - normal(1) * zero_line(3)
420 cross_line(3) = normal(1) * zero_line(2) - normal(2) * zero_line(1)
423 allocate(point_list(3, n_points))
425 pi = 4.0_rp * atan(1.0_rp)
429 angle = 2.0_rp * pi *
real(i - 1, kind = rp) /
real(n_points, kind = rp)
430 temp = cos(angle) * zero_line + sin(angle) * cross_line
432 point_list(:, i) = center + radius * temp
435 call this%add_points(point_list)
444 class(
probes_t),
intent(inout) :: this
445 type(json_file),
intent(inout) :: json
446 type(dofmap_t),
intent(in) :: dof
448 real(kind=rp),
dimension(:,:),
allocatable :: point_list
449 character(len=:),
allocatable :: point_zone_name
450 class(point_zone_t),
pointer :: zone
451 integer :: i, idx, lx, nlindex(4)
452 real(kind=rp) :: x, y, z
455 if (pe_rank .ne. 0)
return
457 call json_get(json,
"name", point_zone_name)
458 zone => neko_point_zone_registry%get_point_zone(point_zone_name)
461 allocate(point_list(3, zone%size))
465 idx = zone%mask%get(i)
468 x = dof%x(nlindex(1), nlindex(2), nlindex(3), nlindex(4))
469 y = dof%y(nlindex(1), nlindex(2), nlindex(3), nlindex(4))
470 z = dof%z(nlindex(1), nlindex(2), nlindex(3), nlindex(4))
472 point_list(:, i) = [x, y, z]
475 call this%add_points(point_list)
485 class(
probes_t),
intent(inout) :: this
486 real(kind=rp),
dimension(:,:),
intent(in) :: new_points
488 real(kind=rp),
dimension(:,:),
allocatable :: temp
489 integer :: n_old, n_new
492 n_old = this%n_local_probes
493 n_new =
size(new_points, 2)
496 if (
allocated(this%xyz))
then
497 call move_alloc(this%xyz, temp)
501 allocate(this%xyz(3, n_old + n_new))
502 if (
allocated(temp))
then
503 this%xyz(:, 1:n_old) = temp
505 this%xyz(:, n_old+1:n_old+n_new) = new_points
507 this%n_local_probes = this%n_local_probes + n_new
518 class(
probes_t),
intent(inout) :: this
519 type(dofmap_t),
intent(in) :: dof
520 character(len=:),
allocatable,
intent(inout) :: output_file
521 character(len=*),
intent(in) :: name
523 character(len=1024) :: header_line
524 real(kind=rp),
allocatable :: global_output_coords(:,:)
525 integer :: i, ierr, out_int
526 type(matrix_t) :: mat_coords
527 logical :: attr_exist = .false.
532 call this%global_interp%find_points_and_redist(this%xyz, &
536 allocate(this%out_values(this%n_local_probes, this%n_fields))
537 allocate(this%out_values_d(this%n_fields))
538 allocate(this%out_vals_trsp(this%n_fields, this%n_local_probes))
540 if (neko_bcknd_device .eq. 1)
then
541 do i = 1, this%n_fields
542 this%out_values_d(i) = c_null_ptr
543 call device_map(this%out_values(:,i), this%out_values_d(i), &
549 call this%fout%init(this%case%output_directory // trim(output_file))
551 select type (ft => this%fout%file_type)
557 write(header_line,
'(I0,A,I0)') this%n_global_probes,
",", this%n_fields
558 do i = 1, this%n_fields
559 header_line = trim(header_line) //
"," // trim(this%which_fields(i))
561 call this%fout%set_header(header_line)
566 allocate(this%n_local_probes_tot(pe_size))
567 allocate(this%n_local_probes_tot_offset(pe_size))
568 call this%setup_offset()
569 if (pe_rank .eq. 0)
then
570 allocate(global_output_coords(3, this%n_global_probes))
571 call this%mat_out%init(this%n_global_probes, this%n_fields)
572 allocate(this%global_output_values(this%n_fields, &
573 this%n_global_probes))
574 call mat_coords%init(this%n_global_probes,3)
576 call mpi_gatherv(this%xyz, 3*this%n_local_probes, &
577 mpi_real_precision, global_output_coords, &
578 3*this%n_local_probes_tot, &
579 3*this%n_local_probes_tot_offset, &
580 mpi_real_precision, 0, neko_comm, ierr)
581 if (pe_rank .eq. 0)
then
582 call trsp(mat_coords%x, this%n_global_probes, &
583 global_output_coords, 3)
585 call this%fout%write(mat_coords)
586 call mat_coords%free()
588 class is (hdf5_file_t)
591 if (this%append_out)
then
592 call ft%set_overwrite(.false.)
594 call ft%set_overwrite(.true.)
596 call mat_coords%init(3, this%n_local_probes,
"coordinates")
597 call copy(mat_coords%x, this%xyz, 3*this%n_local_probes)
601 call ft%set_active_group(
"probes")
604 call ft%read_attribute(
"NSteps", out_int, attr_exist)
608 this%output_controller%nexecutions = out_int
611 call ft%write_dataset(mat_coords)
612 out_int = this%n_global_probes
613 call ft%write_attribute(
"NProbes", out_int)
618 this%seq_io = .false.
619 call this%vec_out%init(this%n_local_probes,
"interpolated_fields_trsp")
622 call mat_coords%free()
625 call neko_error(
"Invalid data. Expected csv_file_t or hdf5_file_t.")
632 class(
probes_t),
intent(inout) :: this
635 if (
allocated(this%xyz))
then
639 if (
allocated(this%out_values))
then
640 deallocate(this%out_values)
643 if (
allocated(this%out_vals_trsp))
then
644 deallocate(this%out_vals_trsp)
647 call this%sampled_fields%free()
649 if (
allocated(this%n_local_probes_tot))
then
650 deallocate(this%n_local_probes_tot)
653 if (
allocated(this%n_local_probes_tot_offset))
then
654 deallocate(this%n_local_probes_tot_offset)
657 if (
allocated(this%global_output_values))
then
658 deallocate(this%global_output_values)
661 if (
allocated(this%which_fields))
then
662 deallocate(this%which_fields)
665 if (
allocated(this%out_values_d))
then
666 do i = 1,
size(this%out_values_d)
667 if (c_associated(this%out_values_d(i)))
then
668 call device_free(this%out_values_d(i))
671 deallocate(this%out_values_d)
674 call this%global_interp%free()
675 call this%mat_out%free()
676 call this%vec_out%free()
683 character(len=LOG_SIZE) :: log_buf
687 call neko_log%section(
'Probes')
688 write(log_buf,
'(A,I6)')
"Number of probes: ", this%n_global_probes
689 call neko_log%message(log_buf)
691 call neko_log%message(
"xyz-coordinates:", lvl = neko_log_debug)
692 do i = 1, this%n_local_probes
693 write(log_buf,
'("(",F10.6,",",F10.6,",",F10.6,")")') this%xyz(:,i)
694 call neko_log%message(log_buf, lvl = neko_log_debug)
698 write(log_buf,
'(A,I6)')
"Number of fields: ", this%n_fields
699 call neko_log%message(log_buf)
700 do i = 1, this%n_fields
701 write(log_buf,
'(A,I6,A,A)') &
702 "Field: ", i,
" ", trim(this%which_fields(i))
703 call neko_log%message(log_buf, lvl = neko_log_debug)
705 call neko_log%end_section()
706 call neko_log%newline()
714 character(len=LOG_SIZE) :: log_buf
717 do i = 1, this%n_local_probes
718 write (log_buf, *) pe_rank,
"/", this%global_interp%pe_owner(i), &
719 "/" , this%global_interp%el_owner0(i)
720 call neko_log%message(log_buf)
721 write(log_buf,
'(A5,"(",F10.6,",",F10.6,",",F10.6,")")') &
722 "rst: ", this%global_interp%rst(:,i)
723 call neko_log%message(log_buf)
731 this%n_local_probes_tot = 0
732 this%n_local_probes_tot_offset = 0
733 this%n_probes_offset = 0
734 call mpi_gather(this%n_local_probes, 1, mpi_integer, &
735 this%n_local_probes_tot, 1, mpi_integer, &
738 call mpi_exscan(this%n_local_probes, this%n_probes_offset, 1, &
739 mpi_integer, mpi_sum, neko_comm, ierr)
740 call mpi_gather(this%n_probes_offset, 1, mpi_integer, &
741 this%n_local_probes_tot_offset, 1, mpi_integer, &
753 class(
probes_t),
intent(inout) :: this
754 type(time_state_t),
intent(in) :: time
756 logical :: do_interp_on_host = .false.
757 character(len=1000) :: group_name
758 real(kind=rp) :: time_
759 type(vector_t) :: vec_time
763 if (time%t .lt. this%start_time)
return
766 do i = 1, this%n_fields
767 call this%global_interp%evaluate(this%out_values(:,i), &
768 this%sampled_fields%items(i)%ptr%x, &
772 if (neko_bcknd_device .eq. 1)
then
773 do i = 1, this%n_fields
774 call device_memcpy(this%out_values(:,i), this%out_values_d(i), &
775 this%n_local_probes, device_to_host, sync = .true.)
779 if (this%output_controller%check(time))
then
780 select type (ft => this%fout%file_type)
784 if (this%seq_io)
then
785 call trsp(this%out_vals_trsp, this%n_fields, &
786 this%out_values, this%n_local_probes)
787 call mpi_gatherv(this%out_vals_trsp, &
788 this%n_fields*this%n_local_probes, &
789 mpi_real_precision, this%global_output_values, &
790 this%n_fields*this%n_local_probes_tot, &
791 this%n_fields*this%n_local_probes_tot_offset, &
792 mpi_real_precision, 0, neko_comm, ierr)
793 if (pe_rank .eq. 0)
then
794 call trsp(this%mat_out%x, this%n_global_probes, &
795 this%global_output_values, this%n_fields)
796 call this%fout%write(this%mat_out, time%t)
799 call neko_error(
"CSV outputs only works sequentially")
802 type is (hdf5_file_t)
804 if (this%seq_io)
then
805 call neko_error(
"HDF5 outputs only works in parallel currently.")
810 if (this%append_out)
then
813 call ft%set_active_group(
"probes")
815 out_int = this%output_controller%nexecutions + 1
816 call ft%write_attribute(
"NSteps", out_int)
818 do i = 1, this%n_fields
819 call copy(this%vec_out%x, this%out_values(:,i), &
821 this%vec_out%name = trim(this%which_fields(i))
822 call ft%write_dataset(this%vec_out)
826 if (pe_rank .eq. 0)
then
827 call vec_time%init(1,
"time")
828 vec_time%x(1) = time%t
830 call vec_time%init(0,
"time")
832 call ft%write_dataset(vec_time)
839 out_int = this%output_controller%nexecutions + 1
841 write(group_name,
'(A,I0)')
"probes/Step_", out_int
843 call ft%set_active_group(
"probes")
845 call ft%write_attribute(
"NSteps", out_int)
847 call ft%set_active_group(trim(group_name))
848 do i = 1, this%n_fields
849 call copy(this%vec_out%x, this%out_values(:,i), &
851 this%vec_out%name = trim(this%which_fields(i))
852 call ft%write_dataset(this%vec_out)
856 call ft%write_attribute(
"time", time_)
862 call neko_error(
"Invalid data. Expected csv_file_t or hdf5_file_t.")
866 call this%output_controller%register_execution()
875 class(
probes_t),
intent(inout) :: this
876 character(len=:),
allocatable :: points_file
877 real(kind=rp),
allocatable :: xyz(:,:)
878 integer,
intent(inout) :: n_local_probes, n_global_probes
879 type(matrix_t) :: mat_in
883 type(file_t) :: file_in
885 call file_in%init(trim(points_file))
887 select type (ft => file_in%file_type)
891 type is (hdf5_file_t)
893 call ft%read_dataset(
"xyz", mat_in,
"rank_0")
897 n_local_probes = mat_in%get_ncols()
898 call mpi_allreduce(n_local_probes, n_global_probes, 1, mpi_integer, &
899 mpi_sum, neko_comm, ierr)
901 allocate(xyz(3, n_local_probes))
902 call copy(xyz, mat_in%x, 3*n_local_probes)
905 call neko_error(
"Invalid data. Expected csv_file_t.")
909 call file_free(file_in)
919 type(csv_file_t),
intent(inout) :: f
920 real(kind=rp),
allocatable :: xyz(:,:)
921 integer,
intent(inout) :: n_local_probes, n_global_probes
922 type(matrix_t) :: mat_in, mat_in2
925 n_lines = f%count_lines()
928 n_global_probes = n_lines
931 if (pe_rank .eq. 0)
then
932 n_local_probes = n_global_probes
933 allocate(xyz(3, n_local_probes))
934 call mat_in%init(n_global_probes,3)
935 call mat_in2%init(3, n_global_probes)
937 call trsp(xyz, 3, mat_in%x, n_global_probes)
940 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.
subroutine, public device_copy(a_d, b_d, n, strm)
Copy a vector .
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.
subroutine, public json_get_subdict_or_empty(json, key, output)
Extract a sub-object from a json object and returns an empty object if the key is missing.
integer, parameter, public neko_log_debug
Debug log level.
type(log_t), public neko_log
Global log stream.
integer, parameter, public log_size
subroutine, public copy(a, b, n)
Copy a vector .
integer, parameter neko_bcknd_device
integer, parameter, public dp
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 probes_init_from_components(this, dof, output_file, name, tolerance, padding)
Initialize based on individual parameters.
subroutine add_points(this, new_points)
Append a new list of points to the exsiting list.
subroutine probes_init_common(this, dof, output_file, name)
Common constructor.
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_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.
integer, parameter, public neko_varname_len
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.
Interface for HDF5 files.
Base abstract type for point zones.
Base abstract class for simulation components.
A struct that contains all info about the time, expand as needed.