Neko  0.9.0
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 ('points')
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  case ('point_zone')
191  call this%read_point_zone(json_point, case%fluid%dm_Xh)
192  case ('none')
193  call json_point%print()
194  call neko_error('No point type specified.')
195  case default
196  call neko_error('Unknown region type ' // point_type)
197  end select
198  end do
199 
200  call mpi_allreduce(this%n_local_probes, this%n_global_probes, 1, &
201  mpi_integer, mpi_sum, neko_comm, ierr)
202 
203  call probes_show(this)
204  call this%init_from_attributes(case%fluid%dm_Xh, output_file)
205 
206  end subroutine probes_init_from_json
207 
208  ! ========================================================================== !
209  ! Readers for different point types
210 
215  subroutine read_file(this, json)
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
220 
221  integer :: n_local, n_global
222 
223  if (pe_rank .ne. 0) return
224 
225  call json_get(json, 'file_name', input_file)
226 
227  call read_probe_locations(this, point_list, n_local, n_global, input_file)
228 
229  call this%add_points(point_list)
230  end subroutine read_file
231 
236  subroutine read_point(this, json)
237  class(probes_t), intent(inout) :: this
238  type(json_file), intent(inout) :: json
239 
240  real(kind=rp), dimension(:,:), allocatable :: point_list
241  real(kind=rp), dimension(:), allocatable :: rp_list_reader
242 
243  ! Ensure only rank 0 reads the coordinates.
244  if (pe_rank .ne. 0) return
245  call json_get(json, 'coordinates', rp_list_reader)
246 
247  if (mod(size(rp_list_reader), 3) /= 0) then
248  call neko_error('Invalid number of coordinates.')
249  end if
250 
251  ! Allocate list of points and reshape the input array
252  allocate(point_list(3, size(rp_list_reader)/3))
253  point_list = reshape(rp_list_reader, [3, size(rp_list_reader)/3])
254 
255  call this%add_points(point_list)
256  end subroutine read_point
257 
261  subroutine read_line(this, json)
262  class(probes_t), intent(inout) :: this
263  type(json_file), intent(inout) :: json
264 
265  real(kind=rp), dimension(:,:), allocatable :: point_list
266  real(kind=rp), dimension(:), allocatable :: start, end
267  real(kind=rp), dimension(3) :: direction
268  real(kind=rp) :: t
269 
270  integer :: n_points, i
271 
272  ! Ensure only rank 0 reads the coordinates.
273  if (pe_rank .ne. 0) return
274  call json_get(json, "start", start)
275  call json_get(json, "end", end)
276  call json_get(json, "amount", n_points)
277 
278  ! If either start or end is not of length 3, error out
279  if (size(start) /= 3 .or. size(end) /= 3) then
280  call neko_error('Invalid start or end coordinates.')
281  end if
282 
283  ! Calculate the number of points
284  allocate(point_list(3, n_points))
285 
286  ! Calculate the direction vector
287  direction = end - start
288  do i = 1, n_points
289  t = real(i - 1, kind = rp) / real(n_points - 1, kind = rp)
290  point_list(:, i) = start + direction * t
291  end do
292 
293  call this%add_points(point_list)
294  end subroutine read_line
295 
305  subroutine read_circle(this, json)
306  class(probes_t), intent(inout) :: this
307  type(json_file), intent(inout) :: json
308 
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
315 
316  real(kind=rp), dimension(3) :: zero_line, cross_line, temp
317  real(kind=rp) :: pi
318 
319  ! Ensure only rank 0 reads the coordinates.
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)
326 
327  ! If either center or normal is not of length 3, error out
328  if (size(center) /= 3 .or. size(normal) /= 3) then
329  call neko_error('Invalid center or normal coordinates.')
330  end if
331  if (axis /= 'x' .and. axis /= 'y' .and. axis /= 'z') then
332  call neko_error('Invalid axis.')
333  end if
334  if (radius <= 0) then
335  call neko_error('Invalid radius.')
336  end if
337  if (n_points <= 0) then
338  call neko_error('Invalid number of points.')
339  end if
340 
341  ! Normalize the normal vector
342  normal = normal / norm2(normal)
343 
344  ! Set the zero line
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]
348 
349  if (1.0_rp - dot_product(zero_line, normal) .le. 1e-6) then
350  call neko_error('Invalid axis and normal.')
351  end if
352 
353  zero_line = zero_line - dot_product(zero_line, normal) * normal
354  zero_line = zero_line / norm2(zero_line)
355 
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)
359 
360  ! Calculate the number of points
361  allocate(point_list(3, n_points))
362 
363  pi = 4.0_rp * atan(1.0_rp)
364 
365  ! Calculate the points
366  do i = 1, n_points
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
369 
370  point_list(:, i) = center + radius * temp
371  end do
372 
373  call this%add_points(point_list)
374  end subroutine read_circle
375 
381  subroutine read_point_zone(this, json, dof)
382  class(probes_t), intent(inout) :: this
383  type(json_file), intent(inout) :: json
384  type(dofmap_t), intent(in) :: dof
385 
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
391 
392  ! Ensure only rank 0 reads the coordinates.
393  if (pe_rank .ne. 0) return
394 
395  call json_get(json, "name", point_zone_name)
396  zone => neko_point_zone_registry%get_point_zone(point_zone_name)
397 
398  ! Allocate list of points and reshape the input array
399  allocate(point_list(3, zone%size))
400 
401  lx = dof%Xh%lx
402  do i = 1, zone%size
403  idx = zone%mask(i)
404 
405  nlindex = nonlinear_index(idx, lx, lx, lx)
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))
409 
410  point_list(:, i) = [x, y, z]
411  end do
412 
413  call this%add_points(point_list)
414  end subroutine read_point_zone
415 
416  ! ========================================================================== !
417  ! Supporting routines
418 
422  subroutine add_points(this, new_points)
423  class(probes_t), intent(inout) :: this
424  real(kind=rp), dimension(:,:), intent(in) :: new_points
425 
426  real(kind=rp), dimension(:,:), allocatable :: temp
427  integer :: n_old, n_new
428 
429  ! Get the current number of points
430  n_old = this%n_local_probes
431  n_new = size(new_points, 2)
432 
433  ! Move current points to a temporary array
434  if (allocated(this%xyz)) then
435  call move_alloc(this%xyz, temp)
436  end if
437 
438  ! Allocate the new array and copy the full list of points
439  allocate(this%xyz(3, n_old + n_new))
440  if (allocated(temp)) then
441  this%xyz(:, 1:n_old) = temp
442  end if
443  this%xyz(:, n_old+1:n_old+n_new) = new_points
444 
445  this%n_local_probes = this%n_local_probes + n_new
446  end subroutine add_points
447 
448  ! ========================================================================== !
449  ! General initialization routine
450 
454  subroutine probes_init_from_attributes(this, dof, output_file)
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(:,:)
460  integer :: i, ierr
461  type(matrix_t) :: mat_coords
462 
464  call this%global_interp%init(dof)
465 
467  call this%global_interp%find_points_and_redist(this%xyz, &
468  this%n_local_probes)
469 
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))
474 
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), &
479  this%n_local_probes)
480  end do
481  end if
482 
484  this%fout = file_t(trim(output_file))
485 
486  select type (ft => this%fout%file_type)
487  type is (csv_file_t)
488 
489  this%seq_io = .true.
490 
491  ! Build the header
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))
495  end do
496  call this%fout%set_header(header_line)
497 
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)
510  end if
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)
519  !! Write the data to the file
520  call this%fout%write(mat_coords)
521  end if
522  class default
523  call neko_error("Invalid data. Expected csv_file_t.")
524  end select
525 
526  end subroutine probes_init_from_attributes
527 
529  subroutine probes_free(this)
530  class(probes_t), intent(inout) :: this
531 
532  if (allocated(this%xyz)) then
533  deallocate(this%xyz)
534  end if
535 
536  if (allocated(this%out_values)) then
537  deallocate(this%out_values)
538  end if
539 
540  if (allocated(this%out_vals_trsp)) then
541  deallocate(this%out_vals_trsp)
542  end if
543 
544  call this%sampled_fields%free()
545 
546  if (allocated(this%n_local_probes_tot)) then
547  deallocate(this%n_local_probes_tot)
548  end if
549 
550  if (allocated(this%n_local_probes_tot_offset)) then
551  deallocate(this%n_local_probes_tot_offset)
552  end if
553 
554  if (allocated(this%global_output_values)) then
555  deallocate(this%global_output_values)
556  end if
557 
558  call this%global_interp%free()
559 
560  end subroutine probes_free
561 
563  subroutine probes_show(this)
564  class(probes_t), intent(in) :: this
565  character(len=LOG_SIZE) :: log_buf ! For logging status
566  integer :: i
567 
568  ! Probes summary
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)
572 
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)
577  end do
578 
579  ! Field summary
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)
586  end do
587  call neko_log%end_section()
588  call neko_log%newline()
589 
590  end subroutine probes_show
591 
593  subroutine probes_debug(this)
594  class(probes_t) :: this
595 
596  character(len=LOG_SIZE) :: log_buf ! For logging status
597  integer :: i
598 
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)
607  end do
608  end subroutine probes_debug
609 
611  subroutine probes_setup_offset(this)
612  class(probes_t) :: this
613  integer :: ierr
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, &
619  0, neko_comm, ierr)
620 
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, &
625  0, neko_comm, ierr)
626 
627 
628 
629  end subroutine probes_setup_offset
630 
635  subroutine probes_evaluate_and_write(this, t, tstep)
636  class(probes_t), intent(inout) :: this
637  real(kind=rp), intent(in) :: t
638  integer, intent(in) :: tstep
639  integer :: i, ierr
640 
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)
645  end do
646 
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.)
651  end do
652  end if
653 
654  if (this%output_controller%check(t, tstep)) then
655  ! Gather all values to rank 0
656  ! If io is only done at root
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)
670  end if
671  else
672  call neko_error('probes sim comp, parallel io need implementation')
673  end if
674 
675  !! Register the execution of the activity
676  call this%output_controller%register_execution()
677  end if
678 
679  end subroutine probes_evaluate_and_write
680 
683  subroutine read_probe_locations(this, xyz, n_local_probes, n_global_probes, &
684  points_file)
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
689 
691  type(file_t) :: file_in
692 
693  file_in = file_t(trim(points_file))
695  select type (ft => file_in%file_type)
696  type is (csv_file_t)
697  call read_xyz_from_csv(xyz, n_local_probes, n_global_probes, ft)
698  this%seq_io = .true.
699  class default
700  call neko_error("Invalid data. Expected csv_file_t.")
701  end select
702 
704  call file_free(file_in)
705 
706  end subroutine read_probe_locations
707 
713  subroutine read_xyz_from_csv(xyz, n_local_probes, n_global_probes, f)
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
718  integer :: n_lines
719 
720  n_lines = f%count_lines()
721 
722  ! Update the number of probes
723  n_global_probes = n_lines
724 
725  ! Initialize the temporal array
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)
731  call f%read(mat_in)
732  call trsp(xyz, 3, mat_in%x, n_global_probes)
733  else
734  n_local_probes = 0
735  allocate(xyz(3, n_local_probes))
736  end if
737 
738  end subroutine read_xyz_from_csv
739 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:54
Retrieves a parameter by name or throws an error.
Definition: json_utils.f90:45
Defines a simulation case.
Definition: case.f90:34
Definition: comm.F90:1
integer pe_rank
MPI rank.
Definition: comm.F90:28
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:73
type(log_t), public neko_log
Global log stream.
Definition: log.f90:65
integer, parameter, public log_size
Definition: log.f90:42
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:382
subroutine read_circle(this, json)
Construct a list of points from a circle.
Definition: probes.F90:306
subroutine add_points(this, new_points)
Append a new list of points to the exsiting list.
Definition: probes.F90:423
subroutine probes_init_from_attributes(this, dof, output_file)
Initialize without json things.
Definition: probes.F90:455
subroutine read_line(this, json)
Construct a list of points from a line.
Definition: probes.F90:262
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:685
subroutine read_point(this, json)
Read a list of points from the json file.
Definition: probes.F90:237
subroutine probes_show(this)
Print current probe status, with number of probes and coordinates.
Definition: probes.F90:564
subroutine probes_free(this)
Destructor.
Definition: probes.F90:530
subroutine probes_debug(this)
Show the status of processor/element owner and error code for each point.
Definition: probes.F90:594
subroutine probes_setup_offset(this)
Setup offset for rank 0.
Definition: probes.F90:612
subroutine probes_evaluate_and_write(this, t, tstep)
Interpolate each probe from its r,s,t coordinates.
Definition: probes.F90:636
subroutine read_file(this, json)
Read a list of points from a csv file.
Definition: probes.F90:216
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:714
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.