47 use json_module,
only : json_file, json_value, json_core
58 use,
intrinsic :: iso_c_binding
62 use mpi_f08,
only : mpi_allreduce, mpi_integer, mpi_sum, &
63 mpi_double_precision, mpi_gatherv, mpi_gather, mpi_exscan
69 real(kind=
rp) :: start_time
71 integer :: n_fields = 0
74 integer :: n_global_probes
76 integer :: n_probes_offset
78 real(kind=
rp),
allocatable :: xyz(:,:)
80 real(kind=
rp),
allocatable :: out_values(:,:)
81 type(c_ptr),
allocatable :: out_values_d(:)
82 real(kind=
rp),
allocatable :: out_vals_trsp(:,:)
84 integer :: n_local_probes
87 character(len=20),
allocatable :: which_fields(:)
89 integer,
allocatable :: n_local_probes_tot_offset(:)
90 integer,
allocatable :: n_local_probes_tot(:)
93 real(kind=
rp),
allocatable :: global_output_values(:,:)
101 procedure, pass(this) :: init_from_components => &
134 class(
probes_t),
intent(inout),
target :: this
135 type(json_file),
intent(inout) :: json
136 class(
case_t),
intent(inout),
target :: case
137 character(len=:),
allocatable :: output_file
138 character(len=:),
allocatable :: input_file
142 character(len=:),
allocatable :: point_type
143 type(json_value),
pointer :: json_point_list
144 type(json_file) :: json_point
145 type(json_core) :: core
146 integer :: idx, n_point_children
148 character(len=:),
allocatable :: name
154 call this%init_base(json,
case)
157 call json%info(
'fields', n_children = this%n_fields)
158 call json_get(json,
'fields', this%which_fields)
159 call json_get(json,
'output_file', output_file)
163 call this%sampled_fields%init(this%n_fields)
164 do i = 1, this%n_fields
166 call this%sampled_fields%assign(i, &
171 this%n_local_probes = 0
172 this%n_global_probes = 0
175 if (json%valid_path(
'points_file'))
then
178 call json_get(json,
'points_file', input_file)
183 this%n_global_probes, input_file)
186 if (json%valid_path(
'points'))
then
189 call json%get(
'points', json_point_list)
190 call json%info(
'points', n_children = n_point_children)
192 do idx = 1, n_point_children
196 select case (point_type)
199 call this%read_file(json_point)
201 call this%read_point(json_point)
203 call this%read_line(json_point)
205 call neko_error(
'Plane probes not implemented yet.')
207 call this%read_circle(json_point)
209 call this%read_point_zone(json_point,
case%fluid%dm_Xh)
211 call json_point%print()
214 call neko_error(
'Unknown region type ' // point_type)
220 call mpi_allreduce(this%n_local_probes, this%n_global_probes, 1, &
227 type(json_file) :: interp_subdict
229 call this%global_interp%init(
case%fluid%dm_Xh, &
230 params_subdict = interp_subdict)
232 call this%init_common(
case%fluid%dm_Xh, output_file, name)
246 class(
probes_t),
intent(inout) :: this
247 type(dofmap_t),
intent(in) :: dof
248 character(len=:),
allocatable,
intent(inout) :: output_file
249 character(len=*),
intent(in) :: name
250 real(kind=dp),
intent(in),
optional :: tolerance, padding
252 character(len=1024) :: header_line
253 real(kind=rp),
allocatable :: global_output_coords(:,:)
255 type(matrix_t) :: mat_coords
257 call this%global_interp%init(dof, tol = tolerance, pad = padding)
259 call this%init_common(dof, output_file, name)
271 class(
probes_t),
intent(inout) :: this
272 type(json_file),
intent(inout) :: json
273 character(len=:),
allocatable :: input_file
274 real(kind=rp),
dimension(:,:),
allocatable :: point_list
276 integer :: n_local, n_global
278 if (pe_rank .ne. 0)
return
280 call json_get(json,
'file_name', input_file)
284 call this%add_points(point_list)
292 class(
probes_t),
intent(inout) :: this
293 type(json_file),
intent(inout) :: json
295 real(kind=rp),
dimension(:,:),
allocatable :: point_list
296 real(kind=rp),
dimension(:),
allocatable :: rp_list_reader
299 if (pe_rank .ne. 0)
return
300 call json_get_or_lookup(json,
'coordinates', rp_list_reader)
302 if (mod(
size(rp_list_reader), 3) /= 0)
then
303 call neko_error(
'Invalid number of coordinates.')
307 allocate(point_list(3,
size(rp_list_reader)/3))
308 point_list = reshape(rp_list_reader, [3,
size(rp_list_reader)/3])
310 call this%add_points(point_list)
317 class(
probes_t),
intent(inout) :: this
318 type(json_file),
intent(inout) :: json
320 real(kind=rp),
dimension(:,:),
allocatable :: point_list
321 real(kind=rp),
dimension(:),
allocatable :: start,
end
322 real(kind=rp),
dimension(3) :: direction
325 integer :: n_points, i
328 if (pe_rank .ne. 0)
return
329 call json_get_or_lookup(json,
"start", start)
330 call json_get_or_lookup(json,
"end",
end)
331 call json_get_or_lookup(json,
"amount", n_points)
334 if (
size(start) /= 3 .or.
size(
end) /= 3) then
335 call neko_error(
'Invalid start or end coordinates.')
339 allocate(point_list(3, n_points))
342 direction =
end - start
344 t =
real(i - 1, kind = rp) /
real(n_points - 1, kind = rp)
345 point_list(:, i) = start + direction * t
348 call this%add_points(point_list)
361 class(
probes_t),
intent(inout) :: this
362 type(json_file),
intent(inout) :: json
364 real(kind=rp),
dimension(:,:),
allocatable :: point_list
365 real(kind=rp),
dimension(:),
allocatable :: center, normal
366 real(kind=rp) :: radius
367 real(kind=rp) :: angle
368 integer :: n_points, i
369 character(len=:),
allocatable :: axis
371 real(kind=rp),
dimension(3) :: zero_line, cross_line, temp
375 if (pe_rank .ne. 0)
return
376 call json_get_or_lookup(json,
"center", center)
377 call json_get_or_lookup(json,
"normal", normal)
378 call json_get_or_lookup(json,
"radius", radius)
379 call json_get_or_lookup(json,
"amount", n_points)
380 call json_get(json,
"axis", axis)
383 if (
size(center) /= 3 .or.
size(normal) /= 3)
then
384 call neko_error(
'Invalid center or normal coordinates.')
386 if (axis /=
'x' .and. axis /=
'y' .and. axis /=
'z')
then
387 call neko_error(
'Invalid axis.')
389 if (radius <= 0)
then
390 call neko_error(
'Invalid radius.')
392 if (n_points <= 0)
then
393 call neko_error(
'Invalid number of points.')
397 normal = normal / norm2(normal)
400 if (axis .eq.
'x') zero_line = [1.0, 0.0, 0.0]
401 if (axis .eq.
'y') zero_line = [0.0, 1.0, 0.0]
402 if (axis .eq.
'z') zero_line = [0.0, 0.0, 1.0]
404 if (1.0_rp - dot_product(zero_line, normal) .le. 1e-6)
then
405 call neko_error(
'Invalid axis and normal.')
408 zero_line = zero_line - dot_product(zero_line, normal) * normal
409 zero_line = zero_line / norm2(zero_line)
411 cross_line(1) = normal(2) * zero_line(3) - normal(3) * zero_line(2)
412 cross_line(2) = normal(3) * zero_line(1) - normal(1) * zero_line(3)
413 cross_line(3) = normal(1) * zero_line(2) - normal(2) * zero_line(1)
416 allocate(point_list(3, n_points))
418 pi = 4.0_rp * atan(1.0_rp)
422 angle = 2.0_rp * pi *
real(i - 1, kind = rp) /
real(n_points, kind = rp)
423 temp = cos(angle) * zero_line + sin(angle) * cross_line
425 point_list(:, i) = center + radius * temp
428 call this%add_points(point_list)
437 class(
probes_t),
intent(inout) :: this
438 type(json_file),
intent(inout) :: json
439 type(dofmap_t),
intent(in) :: dof
441 real(kind=rp),
dimension(:,:),
allocatable :: point_list
442 character(len=:),
allocatable :: point_zone_name
443 class(point_zone_t),
pointer :: zone
444 integer :: i, idx, lx, nlindex(4)
445 real(kind=rp) :: x, y, z
448 if (pe_rank .ne. 0)
return
450 call json_get(json,
"name", point_zone_name)
451 zone => neko_point_zone_registry%get_point_zone(point_zone_name)
454 allocate(point_list(3, zone%size))
458 idx = zone%mask%get(i)
461 x = dof%x(nlindex(1), nlindex(2), nlindex(3), nlindex(4))
462 y = dof%y(nlindex(1), nlindex(2), nlindex(3), nlindex(4))
463 z = dof%z(nlindex(1), nlindex(2), nlindex(3), nlindex(4))
465 point_list(:, i) = [x, y, z]
468 call this%add_points(point_list)
478 class(
probes_t),
intent(inout) :: this
479 real(kind=rp),
dimension(:,:),
intent(in) :: new_points
481 real(kind=rp),
dimension(:,:),
allocatable :: temp
482 integer :: n_old, n_new
485 n_old = this%n_local_probes
486 n_new =
size(new_points, 2)
489 if (
allocated(this%xyz))
then
490 call move_alloc(this%xyz, temp)
494 allocate(this%xyz(3, n_old + n_new))
495 if (
allocated(temp))
then
496 this%xyz(:, 1:n_old) = temp
498 this%xyz(:, n_old+1:n_old+n_new) = new_points
500 this%n_local_probes = this%n_local_probes + n_new
511 class(
probes_t),
intent(inout) :: this
512 type(dofmap_t),
intent(in) :: dof
513 character(len=:),
allocatable,
intent(inout) :: output_file
514 character(len=*),
intent(in) :: name
516 character(len=1024) :: header_line
517 real(kind=rp),
allocatable :: global_output_coords(:,:)
519 type(matrix_t) :: mat_coords
524 call this%global_interp%find_points_and_redist(this%xyz, &
528 allocate(this%out_values(this%n_local_probes, this%n_fields))
529 allocate(this%out_values_d(this%n_fields))
530 allocate(this%out_vals_trsp(this%n_fields, this%n_local_probes))
532 if (neko_bcknd_device .eq. 1)
then
533 do i = 1, this%n_fields
534 this%out_values_d(i) = c_null_ptr
535 call device_map(this%out_values(:,i), this%out_values_d(i), &
541 call this%fout%init(trim(output_file))
543 select type (ft => this%fout%file_type)
549 write(header_line,
'(I0,A,I0)') this%n_global_probes,
",", this%n_fields
550 do i = 1, this%n_fields
551 header_line = trim(header_line) //
"," // trim(this%which_fields(i))
553 call this%fout%set_header(header_line)
558 allocate(this%n_local_probes_tot(pe_size))
559 allocate(this%n_local_probes_tot_offset(pe_size))
560 call this%setup_offset()
561 if (pe_rank .eq. 0)
then
562 allocate(global_output_coords(3, this%n_global_probes))
563 call this%mat_out%init(this%n_global_probes, this%n_fields)
564 allocate(this%global_output_values(this%n_fields, &
565 this%n_global_probes))
566 call mat_coords%init(this%n_global_probes,3)
568 call mpi_gatherv(this%xyz, 3*this%n_local_probes, &
569 mpi_real_precision, global_output_coords, &
570 3*this%n_local_probes_tot, &
571 3*this%n_local_probes_tot_offset, &
572 mpi_real_precision, 0, neko_comm, ierr)
573 if (pe_rank .eq. 0)
then
574 call trsp(mat_coords%x, this%n_global_probes, &
575 global_output_coords, 3)
577 call this%fout%write(mat_coords)
580 call neko_error(
"Invalid data. Expected csv_file_t.")
587 class(
probes_t),
intent(inout) :: this
590 if (
allocated(this%xyz))
then
594 if (
allocated(this%out_values))
then
595 deallocate(this%out_values)
598 if (
allocated(this%out_vals_trsp))
then
599 deallocate(this%out_vals_trsp)
602 call this%sampled_fields%free()
604 if (
allocated(this%n_local_probes_tot))
then
605 deallocate(this%n_local_probes_tot)
608 if (
allocated(this%n_local_probes_tot_offset))
then
609 deallocate(this%n_local_probes_tot_offset)
612 if (
allocated(this%global_output_values))
then
613 deallocate(this%global_output_values)
616 if (
allocated(this%which_fields))
then
617 deallocate(this%which_fields)
620 if (
allocated(this%out_values_d))
then
621 do i = 1,
size(this%out_values_d)
622 if (c_associated(this%out_values_d(i)))
then
623 call device_free(this%out_values_d(i))
626 deallocate(this%out_values_d)
629 call this%global_interp%free()
630 call this%mat_out%free()
637 character(len=LOG_SIZE) :: log_buf
641 call neko_log%section(
'Probes')
642 write(log_buf,
'(A,I6)')
"Number of probes: ", this%n_global_probes
643 call neko_log%message(log_buf)
645 call neko_log%message(
"xyz-coordinates:", lvl = neko_log_debug)
646 do i = 1, this%n_local_probes
647 write(log_buf,
'("(",F10.6,",",F10.6,",",F10.6,")")') this%xyz(:,i)
648 call neko_log%message(log_buf, lvl = neko_log_debug)
652 write(log_buf,
'(A,I6)')
"Number of fields: ", this%n_fields
653 call neko_log%message(log_buf)
654 do i = 1, this%n_fields
655 write(log_buf,
'(A,I6,A,A)') &
656 "Field: ", i,
" ", trim(this%which_fields(i))
657 call neko_log%message(log_buf, lvl = neko_log_debug)
659 call neko_log%end_section()
660 call neko_log%newline()
668 character(len=LOG_SIZE) :: log_buf
671 do i = 1, this%n_local_probes
672 write (log_buf, *) pe_rank,
"/", this%global_interp%pe_owner(i), &
673 "/" , this%global_interp%el_owner0(i)
674 call neko_log%message(log_buf)
675 write(log_buf,
'(A5,"(",F10.6,",",F10.6,",",F10.6,")")') &
676 "rst: ", this%global_interp%rst(:,i)
677 call neko_log%message(log_buf)
685 this%n_local_probes_tot = 0
686 this%n_local_probes_tot_offset = 0
687 this%n_probes_offset = 0
688 call mpi_gather(this%n_local_probes, 1, mpi_integer, &
689 this%n_local_probes_tot, 1, mpi_integer, &
692 call mpi_exscan(this%n_local_probes, this%n_probes_offset, 1, &
693 mpi_integer, mpi_sum, neko_comm, ierr)
694 call mpi_gather(this%n_probes_offset, 1, mpi_integer, &
695 this%n_local_probes_tot_offset, 1, mpi_integer, &
707 class(
probes_t),
intent(inout) :: this
708 type(time_state_t),
intent(in) :: time
710 logical :: do_interp_on_host = .false.
713 if (time%t .lt. this%start_time)
return
716 do i = 1, this%n_fields
717 call this%global_interp%evaluate(this%out_values(:,i), &
718 this%sampled_fields%items(i)%ptr%x, &
722 if (neko_bcknd_device .eq. 1)
then
723 do i = 1, this%n_fields
724 call device_memcpy(this%out_values(:,i), this%out_values_d(i), &
725 this%n_local_probes, device_to_host, sync = .true.)
729 if (this%output_controller%check(time))
then
732 if (this%seq_io)
then
733 call trsp(this%out_vals_trsp, this%n_fields, &
734 this%out_values, this%n_local_probes)
735 call mpi_gatherv(this%out_vals_trsp, &
736 this%n_fields*this%n_local_probes, &
737 mpi_real_precision, this%global_output_values, &
738 this%n_fields*this%n_local_probes_tot, &
739 this%n_fields*this%n_local_probes_tot_offset, &
740 mpi_real_precision, 0, neko_comm, ierr)
741 if (pe_rank .eq. 0)
then
742 call trsp(this%mat_out%x, this%n_global_probes, &
743 this%global_output_values, this%n_fields)
744 call this%fout%write(this%mat_out, time%t)
747 call neko_error(
'probes sim comp, parallel io need implementation')
751 call this%output_controller%register_execution()
760 class(
probes_t),
intent(inout) :: this
761 character(len=:),
allocatable :: points_file
762 real(kind=rp),
allocatable :: xyz(:,:)
763 integer,
intent(inout) :: n_local_probes, n_global_probes
766 type(file_t) :: file_in
768 call file_in%init(trim(points_file))
770 select type (ft => file_in%file_type)
775 call neko_error(
"Invalid data. Expected csv_file_t.")
779 call file_free(file_in)
789 type(csv_file_t),
intent(inout) :: f
790 real(kind=rp),
allocatable :: xyz(:,:)
791 integer,
intent(inout) :: n_local_probes, n_global_probes
792 type(matrix_t) :: mat_in, mat_in2
795 n_lines = f%count_lines()
798 n_global_probes = n_lines
801 if (pe_rank .eq. 0)
then
802 n_local_probes = n_global_probes
803 allocate(xyz(3, n_local_probes))
804 call mat_in%init(n_global_probes,3)
805 call mat_in2%init(3, n_global_probes)
807 call trsp(xyz, 3, mat_in%x, n_global_probes)
810 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.
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
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.
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.