46 use json_module,
only : json_file, json_value, json_core
57 use,
intrinsic :: iso_c_binding
66 integer :: n_global_probes
68 integer :: n_probes_offset
70 real(kind=
rp),
allocatable :: xyz(:,:)
72 real(kind=
rp),
allocatable :: out_values(:,:)
73 type(c_ptr),
allocatable :: out_values_d(:)
74 real(kind=
rp),
allocatable :: out_vals_trsp(:,:)
76 integer :: n_local_probes
79 character(len=20),
allocatable :: which_fields(:)
81 integer,
allocatable :: n_local_probes_tot_offset(:)
82 integer,
allocatable :: n_local_probes_tot(:)
85 real(kind=
rp),
allocatable :: global_output_values(:,:)
93 procedure, pass(this) :: init_from_attributes => &
124 class(
probes_t),
intent(inout) :: this
125 type(json_file),
intent(inout) :: json
126 class(
case_t),
intent(inout),
target :: case
127 character(len=:),
allocatable :: output_file
128 character(len=:),
allocatable :: input_file
132 character(len=:),
allocatable :: point_type
133 type(json_value),
pointer :: json_point_list
134 type(json_file) :: json_point
135 type(json_core) :: core
136 integer :: idx, n_point_children
140 call this%init_base(json,
case)
143 call json%info(
'fields', n_children = this%n_fields)
144 call json_get(json,
'fields', this%which_fields)
145 call json_get(json,
'output_file', output_file)
147 call this%sampled_fields%init(this%n_fields)
148 do i = 1, this%n_fields
150 call this%sampled_fields%assign(i, &
155 this%n_local_probes = 0
156 this%n_global_probes = 0
159 if (json%valid_path(
'points_file'))
then
162 call json_get(json,
'points_file', input_file)
167 this%n_global_probes, input_file)
171 call json%get(
'points', json_point_list)
172 call json%info(
'points', n_children = n_point_children)
174 do idx = 1, n_point_children
178 select case (point_type)
181 call this%read_file(json_point)
183 call this%read_point(json_point)
185 call this%read_line(json_point)
187 call neko_error(
'Plane probes not implemented yet.')
189 call this%read_circle(json_point)
191 call this%read_point_zone(json_point,
case%fluid%dm_Xh)
193 call json_point%print()
196 call neko_error(
'Unknown region type ' // point_type)
200 call mpi_allreduce(this%n_local_probes, this%n_global_probes, 1, &
204 call this%init_from_attributes(
case%fluid%dm_Xh, output_file)
216 class(
probes_t),
intent(inout) :: this
217 type(json_file),
intent(inout) :: json
218 character(len=:),
allocatable :: input_file
219 real(kind=
rp),
dimension(:,:),
allocatable :: point_list
221 integer :: n_local, n_global
225 call json_get(json,
'file_name', input_file)
229 call this%add_points(point_list)
237 class(
probes_t),
intent(inout) :: this
238 type(json_file),
intent(inout) :: json
240 real(kind=
rp),
dimension(:,:),
allocatable :: point_list
241 real(kind=
rp),
dimension(:),
allocatable :: rp_list_reader
245 call json_get(json,
'coordinates', rp_list_reader)
247 if (mod(
size(rp_list_reader), 3) /= 0)
then
248 call neko_error(
'Invalid number of coordinates.')
252 allocate(point_list(3,
size(rp_list_reader)/3))
253 point_list = reshape(rp_list_reader, [3,
size(rp_list_reader)/3])
255 call this%add_points(point_list)
262 class(
probes_t),
intent(inout) :: this
263 type(json_file),
intent(inout) :: json
265 real(kind=
rp),
dimension(:,:),
allocatable :: point_list
266 real(kind=
rp),
dimension(:),
allocatable :: start,
end
267 real(kind=
rp),
dimension(3) :: direction
270 integer :: n_points, i
276 call json_get(json,
"amount", n_points)
279 if (
size(start) /= 3 .or.
size(
end) /= 3) then
280 call neko_error(
'Invalid start or end coordinates.')
284 allocate(point_list(3, n_points))
287 direction =
end - start
289 t =
real(i - 1, kind = rp) /
real(n_points - 1, kind = rp)
290 point_list(:, i) = start + direction * t
293 call this%add_points(point_list)
306 class(
probes_t),
intent(inout) :: this
307 type(json_file),
intent(inout) :: json
309 real(kind=rp),
dimension(:,:),
allocatable :: point_list
310 real(kind=rp),
dimension(:),
allocatable :: center, normal
311 real(kind=rp) :: radius
312 real(kind=rp) :: angle
313 integer :: n_points, i
314 character(len=:),
allocatable :: axis
316 real(kind=rp),
dimension(3) :: zero_line, cross_line, temp
320 if (pe_rank .ne. 0)
return
321 call json_get(json,
"center", center)
322 call json_get(json,
"normal", normal)
323 call json_get(json,
"radius", radius)
324 call json_get(json,
"amount", n_points)
325 call json_get(json,
"axis", axis)
328 if (
size(center) /= 3 .or.
size(normal) /= 3)
then
329 call neko_error(
'Invalid center or normal coordinates.')
331 if (axis /=
'x' .and. axis /=
'y' .and. axis /=
'z')
then
332 call neko_error(
'Invalid axis.')
334 if (radius <= 0)
then
335 call neko_error(
'Invalid radius.')
337 if (n_points <= 0)
then
338 call neko_error(
'Invalid number of points.')
342 normal = normal / norm2(normal)
345 if (axis .eq.
'x') zero_line = [1.0, 0.0, 0.0]
346 if (axis .eq.
'y') zero_line = [0.0, 1.0, 0.0]
347 if (axis .eq.
'z') zero_line = [0.0, 0.0, 1.0]
349 if (1.0_rp - dot_product(zero_line, normal) .le. 1e-6)
then
350 call neko_error(
'Invalid axis and normal.')
353 zero_line = zero_line - dot_product(zero_line, normal) * normal
354 zero_line = zero_line / norm2(zero_line)
356 cross_line(1) = normal(2) * zero_line(3) - normal(3) * zero_line(2)
357 cross_line(2) = normal(3) * zero_line(1) - normal(1) * zero_line(3)
358 cross_line(3) = normal(1) * zero_line(2) - normal(2) * zero_line(1)
361 allocate(point_list(3, n_points))
363 pi = 4.0_rp * atan(1.0_rp)
367 angle = 2.0_rp * pi *
real(i - 1, kind = rp) /
real(n_points, kind = rp)
368 temp = cos(angle) * zero_line + sin(angle) * cross_line
370 point_list(:, i) = center + radius * temp
373 call this%add_points(point_list)
382 class(
probes_t),
intent(inout) :: this
383 type(json_file),
intent(inout) :: json
384 type(dofmap_t),
intent(in) :: dof
386 real(kind=rp),
dimension(:,:),
allocatable :: point_list
387 character(len=:),
allocatable :: point_zone_name
388 class(point_zone_t),
pointer :: zone
389 integer :: i, idx, lx, nlindex(4)
390 real(kind=rp) :: x, y, z
393 if (pe_rank .ne. 0)
return
395 call json_get(json,
"name", point_zone_name)
396 zone => neko_point_zone_registry%get_point_zone(point_zone_name)
399 allocate(point_list(3, zone%size))
406 x = dof%x(nlindex(1), nlindex(2), nlindex(3), nlindex(4))
407 y = dof%y(nlindex(1), nlindex(2), nlindex(3), nlindex(4))
408 z = dof%z(nlindex(1), nlindex(2), nlindex(3), nlindex(4))
410 point_list(:, i) = [x, y, z]
413 call this%add_points(point_list)
423 class(
probes_t),
intent(inout) :: this
424 real(kind=rp),
dimension(:,:),
intent(in) :: new_points
426 real(kind=rp),
dimension(:,:),
allocatable :: temp
427 integer :: n_old, n_new
430 n_old = this%n_local_probes
431 n_new =
size(new_points, 2)
434 if (
allocated(this%xyz))
then
435 call move_alloc(this%xyz, temp)
439 allocate(this%xyz(3, n_old + n_new))
440 if (
allocated(temp))
then
441 this%xyz(:, 1:n_old) = temp
443 this%xyz(:, n_old+1:n_old+n_new) = new_points
445 this%n_local_probes = this%n_local_probes + n_new
455 class(
probes_t),
intent(inout) :: this
456 type(dofmap_t),
intent(in) :: dof
457 character(len=:),
allocatable,
intent(inout) :: output_file
458 character(len=1024) :: header_line
459 real(kind=rp),
allocatable :: global_output_coords(:,:)
461 type(matrix_t) :: mat_coords
464 call this%global_interp%init(dof)
467 call this%global_interp%find_points_and_redist(this%xyz, &
471 allocate(this%out_values(this%n_local_probes, this%n_fields))
472 allocate(this%out_values_d(this%n_fields))
473 allocate(this%out_vals_trsp(this%n_fields, this%n_local_probes))
475 if (neko_bcknd_device .eq. 1)
then
476 do i = 1, this%n_fields
477 this%out_values_d(i) = c_null_ptr
478 call device_map(this%out_values(:,i), this%out_values_d(i), &
484 this%fout = file_t(trim(output_file))
486 select type (ft => this%fout%file_type)
492 write(header_line,
'(I0,A,I0)') this%n_global_probes,
",", this%n_fields
493 do i = 1, this%n_fields
494 header_line = trim(header_line) //
"," // trim(this%which_fields(i))
496 call this%fout%set_header(header_line)
501 allocate(this%n_local_probes_tot(pe_size))
502 allocate(this%n_local_probes_tot_offset(pe_size))
503 call this%setup_offset()
504 if (pe_rank .eq. 0)
then
505 allocate(global_output_coords(3, this%n_global_probes))
506 call this%mat_out%init(this%n_global_probes, this%n_fields)
507 allocate(this%global_output_values(this%n_fields, &
508 this%n_global_probes))
509 call mat_coords%init(this%n_global_probes,3)
511 call mpi_gatherv(this%xyz, 3*this%n_local_probes, &
512 mpi_double_precision, global_output_coords, &
513 3*this%n_local_probes_tot, &
514 3*this%n_local_probes_tot_offset, &
515 mpi_double_precision, 0, neko_comm, ierr)
516 if (pe_rank .eq. 0)
then
517 call trsp(mat_coords%x, this%n_global_probes, &
518 global_output_coords, 3)
520 call this%fout%write(mat_coords)
523 call neko_error(
"Invalid data. Expected csv_file_t.")
530 class(
probes_t),
intent(inout) :: this
532 if (
allocated(this%xyz))
then
536 if (
allocated(this%out_values))
then
537 deallocate(this%out_values)
540 if (
allocated(this%out_vals_trsp))
then
541 deallocate(this%out_vals_trsp)
544 call this%sampled_fields%free()
546 if (
allocated(this%n_local_probes_tot))
then
547 deallocate(this%n_local_probes_tot)
550 if (
allocated(this%n_local_probes_tot_offset))
then
551 deallocate(this%n_local_probes_tot_offset)
554 if (
allocated(this%global_output_values))
then
555 deallocate(this%global_output_values)
558 call this%global_interp%free()
565 character(len=LOG_SIZE) :: log_buf
569 call neko_log%section(
'Probes')
570 write(log_buf,
'(A,I6)')
"Number of probes: ", this%n_global_probes
571 call neko_log%message(log_buf)
573 call neko_log%message(
"xyz-coordinates:", lvl = neko_log_debug)
574 do i = 1, this%n_local_probes
575 write(log_buf,
'("(",F10.6,",",F10.6,",",F10.6,")")') this%xyz(:,i)
576 call neko_log%message(log_buf, lvl = neko_log_debug)
580 write(log_buf,
'(A,I6)')
"Number of fields: ", this%n_fields
581 call neko_log%message(log_buf)
582 do i = 1, this%n_fields
583 write(log_buf,
'(A,I6,A,A)') &
584 "Field: ", i,
" ", trim(this%which_fields(i))
585 call neko_log%message(log_buf, lvl = neko_log_debug)
587 call neko_log%end_section()
588 call neko_log%newline()
596 character(len=LOG_SIZE) :: log_buf
599 do i = 1, this%n_local_probes
600 write (log_buf, *) pe_rank,
"/", this%global_interp%proc_owner(i), &
601 "/" , this%global_interp%el_owner(i), &
602 "/", this%global_interp%error_code(i)
603 call neko_log%message(log_buf)
604 write(log_buf,
'(A5,"(",F10.6,",",F10.6,",",F10.6,")")') &
605 "rst: ", this%global_interp%rst(:,i)
606 call neko_log%message(log_buf)
614 this%n_local_probes_tot = 0
615 this%n_local_probes_tot_offset = 0
616 this%n_probes_offset = 0
617 call mpi_gather(this%n_local_probes, 1, mpi_integer, &
618 this%n_local_probes_tot, 1, mpi_integer, &
621 call mpi_exscan(this%n_local_probes, this%n_probes_offset, 1, &
622 mpi_integer, mpi_sum, neko_comm, ierr)
623 call mpi_gather(this%n_probes_offset, 1, mpi_integer, &
624 this%n_local_probes_tot_offset, 1, mpi_integer, &
636 class(
probes_t),
intent(inout) :: this
637 real(kind=rp),
intent(in) :: t
638 integer,
intent(in) :: tstep
642 do i = 1, this%n_fields
643 call this%global_interp%evaluate(this%out_values(:,i), &
644 this%sampled_fields%items(i)%ptr%x)
647 if (neko_bcknd_device .eq. 1)
then
648 do i = 1, this%n_fields
649 call device_memcpy(this%out_values(:,i), this%out_values_d(i), &
650 this%n_local_probes, device_to_host, sync = .true.)
654 if (this%output_controller%check(t, tstep))
then
657 if (this%seq_io)
then
658 call trsp(this%out_vals_trsp, this%n_fields, &
659 this%out_values, this%n_local_probes)
660 call mpi_gatherv(this%out_vals_trsp, &
661 this%n_fields*this%n_local_probes, &
662 mpi_double_precision, this%global_output_values, &
663 this%n_fields*this%n_local_probes_tot, &
664 this%n_fields*this%n_local_probes_tot_offset, &
665 mpi_double_precision, 0, neko_comm, ierr)
666 if (pe_rank .eq. 0)
then
667 call trsp(this%mat_out%x, this%n_global_probes, &
668 this%global_output_values, this%n_fields)
669 call this%fout%write(this%mat_out, t)
672 call neko_error(
'probes sim comp, parallel io need implementation')
676 call this%output_controller%register_execution()
685 class(
probes_t),
intent(inout) :: this
686 character(len=:),
allocatable :: points_file
687 real(kind=rp),
allocatable :: xyz(:,:)
688 integer,
intent(inout) :: n_local_probes, n_global_probes
691 type(file_t) :: file_in
693 file_in = file_t(trim(points_file))
695 select type (ft => file_in%file_type)
700 call neko_error(
"Invalid data. Expected csv_file_t.")
704 call file_free(file_in)
714 type(csv_file_t),
intent(inout) :: f
715 real(kind=rp),
allocatable :: xyz(:,:)
716 integer,
intent(inout) :: n_local_probes, n_global_probes
717 type(matrix_t) :: mat_in, mat_in2
720 n_lines = f%count_lines()
723 n_global_probes = n_lines
726 if (pe_rank .eq. 0)
then
727 n_local_probes = n_global_probes
728 allocate(xyz(3, n_local_probes))
729 call mat_in%init(n_global_probes,3)
730 call mat_in2%init(3, n_global_probes)
732 call trsp(xyz, 3, mat_in%x, n_global_probes)
735 allocate(xyz(3, n_local_probes))
__device__ void nonlinear_index(const int idx, const int lx, int *index)
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_comm) 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.
Defines a mapping of the degrees of freedom.
Defines a registry for storing solution fields.
type(field_registry_t), target, public neko_field_registry
Global field registry.
pure integer function n_fields(this)
Get the number of fields stored in the registry.
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, 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_init_from_attributes(this, dof, output_file)
Initialize without json things.
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 probes_evaluate_and_write(this, t, tstep)
Interpolate each probe from its r,s,t coordinates.
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.
Simulation components are objects that encapsulate functionality that can be fit to a particular comp...
subroutine compute_(this, t, tstep)
Dummy compute function.
subroutine, public trsp(a, lda, b, ldb)
Transpose of a rectangular tensor .
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.