Neko  0.8.99
A portable framework for high-order spectral element flow simulations
probes.F90
Go to the documentation of this file.
1 ! Copyright (c) 2020-2023, The Neko Authors
2 ! All rights reserved.
3 !
4 ! Redistribution and use in source and binary forms, with or without
5 ! modification, are permitted provided that the following conditions
6 ! are met:
7 !
8 ! * Redistributions of source code must retain the above copyright
9 ! notice, this list of conditions and the following disclaimer.
10 !
11 ! * Redistributions in binary form must reproduce the above
12 ! copyright notice, this list of conditions and the following
13 ! disclaimer in the documentation and/or other materials provided
14 ! with the distribution.
15 !
16 ! * Neither the name of the authors nor the names of its
17 ! contributors may be used to endorse or promote products derived
18 ! from this software without specific prior written permission.
19 !
20 ! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21 ! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22 ! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23 ! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24 ! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25 ! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26 ! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27 ! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28 ! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29 ! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30 ! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31 ! POSSIBILITY OF SUCH DAMAGE.
32 !
37 module probes
38  use num_types, only: rp
39  use matrix, only: matrix_t
41  use utils, only: neko_error, nonlinear_index
42  use field_list, only: field_list_t
45  use dofmap, only: dofmap_t
46  use json_module, only : json_file, json_value, json_core
49  use tensor, only: trsp
50  use point_zone, only: point_zone_t
52  use comm
53  use device
54  use file, only : file_t, file_free
55  use csv_file, only : csv_file_t
56  use case, only : case_t
57  use, intrinsic :: iso_c_binding
58  implicit none
59  private
60 
61  type, public, extends(simulation_component_t) :: probes_t
63  integer :: n_fields = 0
64  type(global_interpolation_t) :: global_interp
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
78  type(field_list_t) :: sampled_fields
79  character(len=20), allocatable :: which_fields(:)
81  integer, allocatable :: n_local_probes_tot_offset(:)
82  integer, allocatable :: n_local_probes_tot(:)
84  logical :: seq_io
85  real(kind=rp), allocatable :: global_output_values(:,:)
87  type(file_t) :: fout
88  type(matrix_t) :: mat_out
89  contains
91  procedure, pass(this) :: init => probes_init_from_json
92  ! Actual constructor
93  procedure, pass(this) :: init_from_attributes => &
96  procedure, pass(this) :: free => probes_free
98  procedure, pass(this) :: setup_offset => probes_setup_offset
100  procedure, pass(this) :: compute_ => probes_evaluate_and_write
101 
102  ! ----------------------------------------------------------------------- !
103  ! Private methods
104 
106  procedure, private, pass(this) :: read_file
108  procedure, private, pass(this) :: read_point
110  procedure, private, pass(this) :: read_line
112  procedure, private, pass(this) :: read_circle
114  procedure, private, pass(this) :: read_point_zone
115 
117  procedure, private, pass(this) :: add_points
118  end type probes_t
119 
120 contains
121 
123  subroutine probes_init_from_json(this, json, case)
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
129  integer :: i, ierr
130 
131  ! JSON variables
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
137 
138  ! Initialize the base class
139  call this%free()
140  call this%init_base(json, case)
141 
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)
146 
147  call this%sampled_fields%init(this%n_fields)
148  do i = 1, this%n_fields
149 
150  call this%sampled_fields%assign(i, &
151  & neko_field_registry%get_field(trim(this%which_fields(i))))
152  end do
153 
154  ! Setup the required arrays and initialize variables.
155  this%n_local_probes = 0
156  this%n_global_probes = 0
157 
158  ! Read the legacy point specification from the points file.
159  if (json%valid_path('points_file')) then
160 
161  ! Todo: We should add a deprecation notice here
162  call json_get(json, 'points_file', input_file)
163 
164  ! This is distributed as to make it similar to parallel file
165  ! formats later, Reads all into rank 0
166  call read_probe_locations(this, this%xyz, this%n_local_probes, &
167  this%n_global_probes, input_file)
168  end if
169 
170  ! Go through the points list and construct the probe list
171  call json%get('points', json_point_list)
172  call json%info('points', n_children = n_point_children)
173 
174  do idx = 1, n_point_children
175  call json_extract_item(core, json_point_list, idx, json_point)
176 
177  call json_get_or_default(json_point, 'type', point_type, 'none')
178  select case (point_type)
179 
180  case ('file')
181  call this%read_file(json_point)
182  case ('point')
183  call this%read_point(json_point)
184  case ('line')
185  call this%read_line(json_point)
186  case ('plane')
187  call neko_error('Plane probes not implemented yet.')
188  case ('circle')
189  call this%read_circle(json_point)
190 
191  case ('point_zone')
192  call this%read_point_zone(json_point, case%fluid%dm_Xh)
193 
194  case ('none')
195  call json_point%print()
196  call neko_error('No point type specified.')
197  case default
198  call neko_error('Unknown region type ' // point_type)
199  end select
200  end do
201 
202  call mpi_allreduce(this%n_local_probes, this%n_global_probes, 1, &
203  mpi_integer, mpi_sum, neko_comm, ierr)
204 
205  call probes_show(this)
206  call this%init_from_attributes(case%fluid%dm_Xh, output_file)
207 
208  end subroutine probes_init_from_json
209 
210  ! ========================================================================== !
211  ! Readers for different point types
212 
217  subroutine read_file(this, json)
218  class(probes_t), intent(inout) :: this
219  type(json_file), intent(inout) :: json
220  character(len=:), allocatable :: input_file
221  real(kind=rp), dimension(:,:), allocatable :: point_list
222 
223  integer :: n_local, n_global
224 
225  if (pe_rank .ne. 0) return
226 
227  call json_get(json, 'file_name', input_file)
228 
229  call read_probe_locations(this, point_list, n_local, n_global, input_file)
230 
231  call this%add_points(point_list)
232  end subroutine read_file
233 
238  subroutine read_point(this, json)
239  class(probes_t), intent(inout) :: this
240  type(json_file), intent(inout) :: json
241 
242  real(kind=rp), dimension(:,:), allocatable :: point_list
243  real(kind=rp), dimension(:), allocatable :: rp_list_reader
244  logical :: found
245 
246  ! Ensure only rank 0 reads the coordinates.
247  if (pe_rank .ne. 0) return
248  call json_get(json, 'coordinates', rp_list_reader)
249 
250  ! Check if the coordinates were found and were valid
251  if (.not. found) then
252  call neko_error('No coordinates found.')
253  end if
254 
255  if (mod(size(rp_list_reader), 3) /= 0) then
256  call neko_error('Invalid number of coordinates.')
257  end if
258 
259  ! Allocate list of points and reshape the input array
260  allocate(point_list(3, size(rp_list_reader)/3))
261  point_list = reshape(rp_list_reader, [3, size(rp_list_reader)/3])
262 
263  call this%add_points(point_list)
264  end subroutine read_point
265 
269  subroutine read_line(this, json)
270  class(probes_t), intent(inout) :: this
271  type(json_file), intent(inout) :: json
272 
273  real(kind=rp), dimension(:,:), allocatable :: point_list
274  real(kind=rp), dimension(:), allocatable :: start, end
275  real(kind=rp), dimension(3) :: direction
276  real(kind=rp) :: t
277 
278  integer :: n_points, i
279 
280  ! Ensure only rank 0 reads the coordinates.
281  if (pe_rank .ne. 0) return
282  call json_get(json, "start", start)
283  call json_get(json, "end", end)
284  call json_get(json, "amount", n_points)
285 
286  ! If either start or end is not of length 3, error out
287  if (size(start) /= 3 .or. size(end) /= 3) then
288  call neko_error('Invalid start or end coordinates.')
289  end if
290 
291  ! Calculate the number of points
292  allocate(point_list(3, n_points))
293 
294  ! Calculate the direction vector
295  direction = end - start
296  do i = 1, n_points
297  t = real(i - 1, kind = rp) / real(n_points - 1, kind = rp)
298  point_list(:, i) = start + direction * t
299  end do
300 
301  call this%add_points(point_list)
302  end subroutine read_line
303 
313  subroutine read_circle(this, json)
314  class(probes_t), intent(inout) :: this
315  type(json_file), intent(inout) :: json
316 
317  real(kind=rp), dimension(:,:), allocatable :: point_list
318  real(kind=rp), dimension(:), allocatable :: center, normal
319  real(kind=rp) :: radius
320  real(kind=rp) :: angle
321  integer :: n_points, i
322  character(len=:), allocatable :: axis
323 
324  real(kind=rp), dimension(3) :: zero_line, cross_line, temp
325  real(kind=rp) :: pi
326 
327  ! Ensure only rank 0 reads the coordinates.
328  if (pe_rank .ne. 0) return
329  call json_get(json, "center", center)
330  call json_get(json, "normal", normal)
331  call json_get(json, "radius", radius)
332  call json_get(json, "amount", n_points)
333  call json_get(json, "axis", axis)
334 
335  ! If either center or normal is not of length 3, error out
336  if (size(center) /= 3 .or. size(normal) /= 3) then
337  call neko_error('Invalid center or normal coordinates.')
338  end if
339  if (axis /= 'x' .and. axis /= 'y' .and. axis /= 'z') then
340  call neko_error('Invalid axis.')
341  end if
342  if (radius <= 0) then
343  call neko_error('Invalid radius.')
344  end if
345  if (n_points <= 0) then
346  call neko_error('Invalid number of points.')
347  end if
348 
349  ! Normalize the normal vector
350  normal = normal / norm2(normal)
351 
352  ! Set the zero line
353  if (axis .eq. 'x') zero_line = [1.0, 0.0, 0.0]
354  if (axis .eq. 'y') zero_line = [0.0, 1.0, 0.0]
355  if (axis .eq. 'z') zero_line = [0.0, 0.0, 1.0]
356 
357  if (1.0_rp - dot_product(zero_line, normal) .le. 1e-6) then
358  call neko_error('Invalid axis and normal.')
359  end if
360 
361  zero_line = zero_line - dot_product(zero_line, normal) * normal
362  zero_line = zero_line / norm2(zero_line)
363 
364  cross_line(1) = normal(2) * zero_line(3) - normal(3) * zero_line(2)
365  cross_line(2) = normal(3) * zero_line(1) - normal(1) * zero_line(3)
366  cross_line(3) = normal(1) * zero_line(2) - normal(2) * zero_line(1)
367 
368  ! Calculate the number of points
369  allocate(point_list(3, n_points))
370 
371  pi = 4.0_rp * atan(1.0_rp)
372 
373  ! Calculate the points
374  do i = 1, n_points
375  angle = 2.0_rp * pi * real(i - 1, kind = rp) / real(n_points, kind = rp)
376  temp = cos(angle) * zero_line + sin(angle) * cross_line
377 
378  point_list(:, i) = center + radius * temp
379  end do
380 
381  call this%add_points(point_list)
382  end subroutine read_circle
383 
389  subroutine read_point_zone(this, json, dof)
390  class(probes_t), intent(inout) :: this
391  type(json_file), intent(inout) :: json
392  type(dofmap_t), intent(in) :: dof
393 
394  real(kind=rp), dimension(:,:), allocatable :: point_list
395  character(len=:), allocatable :: point_zone_name
396  class(point_zone_t), pointer :: zone
397  integer :: i, idx, lx, nlindex(4)
398  real(kind=rp) :: x, y, z
399 
400  ! Ensure only rank 0 reads the coordinates.
401  if (pe_rank .ne. 0) return
402 
403  call json_get(json, "name", point_zone_name)
404  zone => neko_point_zone_registry%get_point_zone(point_zone_name)
405 
406  ! Allocate list of points and reshape the input array
407  allocate(point_list(3, zone%size))
408 
409  lx = dof%Xh%lx
410  do i = 1, zone%size
411  idx = zone%mask(i)
412 
413  nlindex = nonlinear_index(idx, lx, lx, lx)
414  x = dof%x(nlindex(1), nlindex(2), nlindex(3), nlindex(4))
415  y = dof%y(nlindex(1), nlindex(2), nlindex(3), nlindex(4))
416  z = dof%z(nlindex(1), nlindex(2), nlindex(3), nlindex(4))
417 
418  point_list(:, i) = [x, y, z]
419  end do
420 
421  call this%add_points(point_list)
422  end subroutine read_point_zone
423 
424  ! ========================================================================== !
425  ! Supporting routines
426 
430  subroutine add_points(this, new_points)
431  class(probes_t), intent(inout) :: this
432  real(kind=rp), dimension(:,:), intent(in) :: new_points
433 
434  real(kind=rp), dimension(:,:), allocatable :: temp
435  integer :: n_old, n_new
436 
437  ! Get the current number of points
438  n_old = this%n_local_probes
439  n_new = size(new_points, 2)
440 
441  ! Move current points to a temporary array
442  if (allocated(this%xyz)) then
443  call move_alloc(this%xyz, temp)
444  end if
445 
446  ! Allocate the new array and copy the full list of points
447  allocate(this%xyz(3, n_old + n_new))
448  if (allocated(temp)) then
449  this%xyz(:, 1:n_old) = temp
450  end if
451  this%xyz(:, n_old+1:n_old+n_new) = new_points
452 
453  this%n_local_probes = this%n_local_probes + n_new
454  end subroutine add_points
455 
456  ! ========================================================================== !
457  ! General initialization routine
458 
462  subroutine probes_init_from_attributes(this, dof, output_file)
463  class(probes_t), intent(inout) :: this
464  type(dofmap_t), intent(in) :: dof
465  character(len=:), allocatable, intent(inout) :: output_file
466  character(len=1024) :: header_line
467  real(kind=rp), allocatable :: global_output_coords(:,:)
468  integer :: i, ierr
469  type(matrix_t) :: mat_coords
470 
472  call this%global_interp%init(dof)
473 
475  call this%global_interp%find_points_and_redist(this%xyz, &
476  this%n_local_probes)
477 
479  allocate(this%out_values(this%n_local_probes, this%n_fields))
480  allocate(this%out_values_d(this%n_fields))
481  allocate(this%out_vals_trsp(this%n_fields, this%n_local_probes))
482 
483  if (neko_bcknd_device .eq. 1) then
484  do i = 1, this%n_fields
485  this%out_values_d(i) = c_null_ptr
486  call device_map(this%out_values(:,i), this%out_values_d(i), &
487  this%n_local_probes)
488  end do
489  end if
490 
492  this%fout = file_t(trim(output_file))
493 
494  select type (ft => this%fout%file_type)
495  type is (csv_file_t)
496 
497  this%seq_io = .true.
498 
499  ! Build the header
500  write(header_line, '(I0,A,I0)') this%n_global_probes, ",", this%n_fields
501  do i = 1, this%n_fields
502  header_line = trim(header_line) // "," // trim(this%which_fields(i))
503  end do
504  call this%fout%set_header(header_line)
505 
509  allocate(this%n_local_probes_tot(pe_size))
510  allocate(this%n_local_probes_tot_offset(pe_size))
511  call this%setup_offset()
512  if (pe_rank .eq. 0) then
513  allocate(global_output_coords(3, this%n_global_probes))
514  call this%mat_out%init(this%n_global_probes, this%n_fields)
515  allocate(this%global_output_values(this%n_fields, &
516  this%n_global_probes))
517  call mat_coords%init(this%n_global_probes,3)
518  end if
519  call mpi_gatherv(this%xyz, 3*this%n_local_probes, &
520  mpi_double_precision, global_output_coords, &
521  3*this%n_local_probes_tot, &
522  3*this%n_local_probes_tot_offset, &
523  mpi_double_precision, 0, neko_comm, ierr)
524  if (pe_rank .eq. 0) then
525  call trsp(mat_coords%x, this%n_global_probes, &
526  global_output_coords, 3)
527  !! Write the data to the file
528  call this%fout%write(mat_coords)
529  end if
530  class default
531  call neko_error("Invalid data. Expected csv_file_t.")
532  end select
533 
534  end subroutine probes_init_from_attributes
535 
537  subroutine probes_free(this)
538  class(probes_t), intent(inout) :: this
539 
540  if (allocated(this%xyz)) then
541  deallocate(this%xyz)
542  end if
543 
544  if (allocated(this%out_values)) then
545  deallocate(this%out_values)
546  end if
547 
548  if (allocated(this%out_vals_trsp)) then
549  deallocate(this%out_vals_trsp)
550  end if
551 
552  call this%sampled_fields%free()
553 
554  if (allocated(this%n_local_probes_tot)) then
555  deallocate(this%n_local_probes_tot)
556  end if
557 
558  if (allocated(this%n_local_probes_tot_offset)) then
559  deallocate(this%n_local_probes_tot_offset)
560  end if
561 
562  if (allocated(this%global_output_values)) then
563  deallocate(this%global_output_values)
564  end if
565 
566  call this%global_interp%free()
567 
568  end subroutine probes_free
569 
571  subroutine probes_show(this)
572  class(probes_t), intent(in) :: this
573  character(len=LOG_SIZE) :: log_buf ! For logging status
574  integer :: i
575 
576  ! Probes summary
577  call neko_log%section('Probes')
578  write(log_buf, '(A,I6)') "Number of probes: ", this%n_global_probes
579  call neko_log%message(log_buf)
580 
581  call neko_log%message("xyz-coordinates:", lvl = neko_log_debug)
582  do i = 1, this%n_local_probes
583  write(log_buf, '("(",F10.6,",",F10.6,",",F10.6,")")') this%xyz(:,i)
584  call neko_log%message(log_buf, lvl = neko_log_debug)
585  end do
586 
587  ! Field summary
588  write(log_buf, '(A,I6)') "Number of fields: ", this%n_fields
589  call neko_log%message(log_buf)
590  do i = 1, this%n_fields
591  write(log_buf, '(A,I6,A,A)') &
592  "Field: ", i, " ", trim(this%which_fields(i))
593  call neko_log%message(log_buf, lvl = neko_log_debug)
594  end do
595  call neko_log%end_section()
596  call neko_log%newline()
597 
598  end subroutine probes_show
599 
601  subroutine probes_debug(this)
602  class(probes_t) :: this
603 
604  character(len=LOG_SIZE) :: log_buf ! For logging status
605  integer :: i
606 
607  do i = 1, this%n_local_probes
608  write (log_buf, *) pe_rank, "/", this%global_interp%proc_owner(i), &
609  "/" , this%global_interp%el_owner(i), &
610  "/", this%global_interp%error_code(i)
611  call neko_log%message(log_buf)
612  write(log_buf, '(A5,"(",F10.6,",",F10.6,",",F10.6,")")') &
613  "rst: ", this%global_interp%rst(:,i)
614  call neko_log%message(log_buf)
615  end do
616  end subroutine probes_debug
617 
619  subroutine probes_setup_offset(this)
620  class(probes_t) :: this
621  integer :: ierr
622  this%n_local_probes_tot = 0
623  this%n_local_probes_tot_offset = 0
624  this%n_probes_offset = 0
625  call mpi_gather(this%n_local_probes, 1, mpi_integer, &
626  this%n_local_probes_tot, 1, mpi_integer, &
627  0, neko_comm, ierr)
628 
629  call mpi_exscan(this%n_local_probes, this%n_probes_offset, 1, &
630  mpi_integer, mpi_sum, neko_comm, ierr)
631  call mpi_gather(this%n_probes_offset, 1, mpi_integer, &
632  this%n_local_probes_tot_offset, 1, mpi_integer, &
633  0, neko_comm, ierr)
634 
635 
636 
637  end subroutine probes_setup_offset
638 
643  subroutine probes_evaluate_and_write(this, t, tstep)
644  class(probes_t), intent(inout) :: this
645  real(kind=rp), intent(in) :: t
646  integer, intent(in) :: tstep
647  integer :: i, ierr
648 
650  do i = 1, this%n_fields
651  call this%global_interp%evaluate(this%out_values(:,i), &
652  this%sampled_fields%items(i)%ptr%x)
653  end do
654 
655  if (neko_bcknd_device .eq. 1) then
656  do i = 1, this%n_fields
657  call device_memcpy(this%out_values(:,i), this%out_values_d(i), &
658  this%n_local_probes, device_to_host, sync = .true.)
659  end do
660  end if
661 
662  if (this%output_controller%check(t, tstep)) then
663  ! Gather all values to rank 0
664  ! If io is only done at root
665  if (this%seq_io) then
666  call trsp(this%out_vals_trsp, this%n_fields, &
667  this%out_values, this%n_local_probes)
668  call mpi_gatherv(this%out_vals_trsp, &
669  this%n_fields*this%n_local_probes, &
670  mpi_double_precision, this%global_output_values, &
671  this%n_fields*this%n_local_probes_tot, &
672  this%n_fields*this%n_local_probes_tot_offset, &
673  mpi_double_precision, 0, neko_comm, ierr)
674  if (pe_rank .eq. 0) then
675  call trsp(this%mat_out%x, this%n_global_probes, &
676  this%global_output_values, this%n_fields)
677  call this%fout%write(this%mat_out, t)
678  end if
679  else
680  call neko_error('probes sim comp, parallel io need implementation')
681  end if
682 
683  !! Register the execution of the activity
684  call this%output_controller%register_execution()
685  end if
686 
687  end subroutine probes_evaluate_and_write
688 
691  subroutine read_probe_locations(this, xyz, n_local_probes, n_global_probes, &
692  points_file)
693  class(probes_t), intent(inout) :: this
694  character(len=:), allocatable :: points_file
695  real(kind=rp), allocatable :: xyz(:,:)
696  integer, intent(inout) :: n_local_probes, n_global_probes
697 
699  type(file_t) :: file_in
700 
701  file_in = file_t(trim(points_file))
703  select type (ft => file_in%file_type)
704  type is (csv_file_t)
705  call read_xyz_from_csv(xyz, n_local_probes, n_global_probes, ft)
706  this%seq_io = .true.
707  class default
708  call neko_error("Invalid data. Expected csv_file_t.")
709  end select
710 
712  call file_free(file_in)
713 
714  end subroutine read_probe_locations
715 
721  subroutine read_xyz_from_csv(xyz, n_local_probes, n_global_probes, f)
722  type(csv_file_t), intent(inout) :: f
723  real(kind=rp), allocatable :: xyz(:,:)
724  integer, intent(inout) :: n_local_probes, n_global_probes
725  type(matrix_t) :: mat_in, mat_in2
726  integer :: n_lines
727 
728  n_lines = f%count_lines()
729 
730  ! Update the number of probes
731  n_global_probes = n_lines
732 
733  ! Initialize the temporal array
734  if (pe_rank .eq. 0) then
735  n_local_probes = n_global_probes
736  allocate(xyz(3, n_local_probes))
737  call mat_in%init(n_global_probes,3)
738  call mat_in2%init(3, n_global_probes)
739  call f%read(mat_in)
740  call trsp(xyz, 3, mat_in%x, n_global_probes)
741  else
742  n_local_probes = 0
743  allocate(xyz(3, n_local_probes))
744  end if
745 
746  end subroutine read_xyz_from_csv
747 end module probes
__device__ void nonlinear_index(const int idx, const int lx, int *index)
double real
Definition: device_config.h:12
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Definition: json_utils.f90:53
Retrieves a parameter by name or throws an error.
Definition: json_utils.f90:44
Defines a simulation case.
Definition: case.f90:34
Definition: comm.F90:1
integer pe_rank
MPI rank.
Definition: comm.F90:26
type(mpi_comm) neko_comm
MPI communicator.
Definition: comm.F90:16
File format for .csv files, used for any read/write operations involving floating point data.
Definition: csv_file.f90:35
Device abstraction, common interface for various accelerators.
Definition: device.F90:34
Defines a mapping of the degrees of freedom.
Definition: dofmap.f90:35
Defines a registry for storing solution fields.
type(field_registry_t), target, public neko_field_registry
Global field registry.
Module for file I/O operations.
Definition: file.f90:34
subroutine file_free(this)
File operation destructor.
Definition: file.f90:134
Implements global_interpolation given a dofmap.
Utilities for retrieving parameters from the case files.
Definition: json_utils.f90:34
Logging routines.
Definition: log.f90:34
integer, parameter, public neko_log_debug
Debug log level.
Definition: log.f90:69
type(log_t), public neko_log
Global log stream.
Definition: log.f90:61
integer, parameter, public log_size
Definition: log.f90:40
Defines a matrix.
Definition: matrix.f90:34
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
type(point_zone_registry_t), target, public neko_point_zone_registry
Global point_zone registry.
Implements probes.
Definition: probes.F90:37
subroutine read_point_zone(this, json, dof)
Construct a list of points from a point zone.
Definition: probes.F90:390
subroutine read_circle(this, json)
Construct a list of points from a circle.
Definition: probes.F90:314
subroutine add_points(this, new_points)
Append a new list of points to the exsiting list.
Definition: probes.F90:431
subroutine probes_init_from_attributes(this, dof, output_file)
Initialize without json things.
Definition: probes.F90:463
subroutine read_line(this, json)
Construct a list of points from a line.
Definition: probes.F90:270
subroutine read_probe_locations(this, xyz, n_local_probes, n_global_probes, points_file)
Initialize the physical coordinates from a csv input file.
Definition: probes.F90:693
subroutine read_point(this, json)
Read a list of points from the json file.
Definition: probes.F90:239
subroutine probes_show(this)
Print current probe status, with number of probes and coordinates.
Definition: probes.F90:572
subroutine probes_free(this)
Destructor.
Definition: probes.F90:538
subroutine probes_debug(this)
Show the status of processor/element owner and error code for each point.
Definition: probes.F90:602
subroutine probes_setup_offset(this)
Setup offset for rank 0.
Definition: probes.F90:620
subroutine probes_evaluate_and_write(this, t, tstep)
Interpolate each probe from its r,s,t coordinates.
Definition: probes.F90:644
subroutine read_file(this, json)
Read a list of points from a csv file.
Definition: probes.F90:218
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.
Definition: probes.F90:722
subroutine probes_init_from_json(this, json, case)
Constructor from json.
Definition: probes.F90:124
Simulation components are objects that encapsulate functionality that can be fit to a particular comp...
Tensor operations.
Definition: tensor.f90:61
subroutine, public trsp(a, lda, b, ldb)
Transpose of a rectangular tensor .
Definition: tensor.f90:121
Utilities.
Definition: utils.f90:35
field_list_t, To be able to group fields together
Definition: field_list.f90:13
A wrapper around a polymorphic generic_file_t that handles its init. This is essentially a factory fo...
Definition: file.f90:54
Implements global interpolation for arbitrary points in the domain.
Base abstract type for point zones.
Definition: point_zone.f90:47
Base abstract class for simulation components.