Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
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!
37module probes
38 use num_types, only : rp, dp
39 use matrix, only : matrix_t
42 use field_list, only : field_list_t
43 use time_state, only : time_state_t
45 use registry, only : neko_registry
46 use dofmap, only : dofmap_t
47 use json_module, only : json_file, json_value, json_core
52 use tensor, only : trsp
53 use point_zone, only : point_zone_t
55 use file, only : file_t, file_free
56 use csv_file, only : csv_file_t
57 use case, only : case_t
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
64 implicit none
65 private
66
67 type, public, extends(simulation_component_t) :: probes_t
69 real(kind=rp) :: start_time
71 integer :: n_fields = 0
72 type(global_interpolation_t) :: global_interp
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
86 type(field_list_t) :: sampled_fields
87 character(len=20), allocatable :: which_fields(:)
89 integer, allocatable :: n_local_probes_tot_offset(:)
90 integer, allocatable :: n_local_probes_tot(:)
92 logical :: seq_io
93 real(kind=rp), allocatable :: global_output_values(:,:)
95 type(file_t) :: fout
96 type(matrix_t) :: mat_out
97 contains
99 procedure, pass(this) :: init => probes_init_from_json
101 procedure, pass(this) :: init_from_components => &
104 procedure, private, pass(this) :: init_common => probes_init_common
106 procedure, pass(this) :: free => probes_free
108 procedure, pass(this) :: setup_offset => probes_setup_offset
110 procedure, pass(this) :: compute_ => probes_evaluate_and_write
111
112 ! ----------------------------------------------------------------------- !
113 ! Private methods
114
116 procedure, private, pass(this) :: read_file
118 procedure, private, pass(this) :: read_point
120 procedure, private, pass(this) :: read_line
122 procedure, private, pass(this) :: read_circle
124 procedure, private, pass(this) :: read_point_zone
125
127 procedure, private, pass(this) :: add_points
128 end type probes_t
129
130contains
131
133 subroutine probes_init_from_json(this, json, case)
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
139 integer :: i, ierr
140
141 ! JSON variables
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
147
148 character(len=:), allocatable :: name
149
150 ! Initialize the base class
151 call this%free()
152
153 call json_get_or_default(json, "name", name, "probes")
154 call this%init_base(json, case)
155
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)
160 call json_get_or_lookup_or_default(json, 'start_time', this%start_time, &
161 -1.0_rp)
162
163 call this%sampled_fields%init(this%n_fields)
164 do i = 1, this%n_fields
165
166 call this%sampled_fields%assign(i, &
167 & neko_registry%get_field(trim(this%which_fields(i))))
168 end do
169
170 ! Setup the required arrays and initialize variables.
171 this%n_local_probes = 0
172 this%n_global_probes = 0
173
174 ! Read the legacy point specification from the points file.
175 if (json%valid_path('points_file')) then
176
177 ! Todo: We should add a deprecation notice here
178 call json_get(json, 'points_file', input_file)
179
180 ! This is distributed as to make it similar to parallel file
181 ! formats later, Reads all into rank 0
182 call read_probe_locations(this, this%xyz, this%n_local_probes, &
183 this%n_global_probes, input_file)
184 end if
185
186 if (json%valid_path('points')) then
187
188 ! Go through the points list and construct the probe list
189 call json%get('points', json_point_list)
190 call json%info('points', n_children = n_point_children)
191
192 do idx = 1, n_point_children
193 call json_extract_item(core, json_point_list, idx, json_point)
194
195 call json_get_or_default(json_point, 'type', point_type, 'none')
196 select case (point_type)
197
198 case ('file')
199 call this%read_file(json_point)
200 case ('points')
201 call this%read_point(json_point)
202 case ('line')
203 call this%read_line(json_point)
204 case ('plane')
205 call neko_error('Plane probes not implemented yet.')
206 case ('circle')
207 call this%read_circle(json_point)
208 case ('point_zone')
209 call this%read_point_zone(json_point, case%fluid%dm_Xh)
210 case ('none')
211 call json_point%print()
212 call neko_error('No point type specified.')
213 case default
214 call neko_error('Unknown region type ' // point_type)
215 end select
216 end do
217
218 end if
219
220 call mpi_allreduce(this%n_local_probes, this%n_global_probes, 1, &
221 mpi_integer, mpi_sum, neko_comm, ierr)
222
223 call probes_show(this)
224
225 ! Get interpolation parameters from json
226 block
227 type(json_file) :: interp_subdict
228 call json_get_subdict_or_empty(json, "interpolation", interp_subdict)
229 call this%global_interp%init(case%fluid%dm_Xh, &
230 params_subdict = interp_subdict)
231
232 call this%init_common(case%fluid%dm_Xh, output_file, name)
233
234 end block
235
236 end subroutine probes_init_from_json
237
244 subroutine probes_init_from_components(this, dof, output_file, name, &
245 tolerance, padding)
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
251
252 character(len=1024) :: header_line
253 real(kind=rp), allocatable :: global_output_coords(:,:)
254 integer :: i, ierr
255 type(matrix_t) :: mat_coords
256
257 call this%global_interp%init(dof, tol = tolerance, pad = padding)
258
259 call this%init_common(dof, output_file, name)
260
261 end subroutine probes_init_from_components
262
263 ! ========================================================================== !
264 ! Readers for different point types
265
270 subroutine read_file(this, json)
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
275
276 integer :: n_local, n_global
277
278 if (pe_rank .ne. 0) return
279
280 call json_get(json, 'file_name', input_file)
281
282 call read_probe_locations(this, point_list, n_local, n_global, input_file)
283
284 call this%add_points(point_list)
285 end subroutine read_file
286
291 subroutine read_point(this, json)
292 class(probes_t), intent(inout) :: this
293 type(json_file), intent(inout) :: json
294
295 real(kind=rp), dimension(:,:), allocatable :: point_list
296 real(kind=rp), dimension(:), allocatable :: rp_list_reader
297
298 ! Ensure only rank 0 reads the coordinates.
299 if (pe_rank .ne. 0) return
300 call json_get_or_lookup(json, 'coordinates', rp_list_reader)
301
302 if (mod(size(rp_list_reader), 3) /= 0) then
303 call neko_error('Invalid number of coordinates.')
304 end if
305
306 ! Allocate list of points and reshape the input array
307 allocate(point_list(3, size(rp_list_reader)/3))
308 point_list = reshape(rp_list_reader, [3, size(rp_list_reader)/3])
309
310 call this%add_points(point_list)
311 end subroutine read_point
312
316 subroutine read_line(this, json)
317 class(probes_t), intent(inout) :: this
318 type(json_file), intent(inout) :: json
319
320 real(kind=rp), dimension(:,:), allocatable :: point_list
321 real(kind=rp), dimension(:), allocatable :: start, end
322 real(kind=rp), dimension(3) :: direction
323 real(kind=rp) :: t
324
325 integer :: n_points, i
326
327 ! Ensure only rank 0 reads the coordinates.
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)
332
333 ! If either start or end is not of length 3, error out
334 if (size(start) /= 3 .or. size(end) /= 3) then
335 call neko_error('Invalid start or end coordinates.')
336 end if
337
338 ! Calculate the number of points
339 allocate(point_list(3, n_points))
340
341 ! Calculate the direction vector
342 direction = end - start
343 do i = 1, n_points
344 t = real(i - 1, kind = rp) / real(n_points - 1, kind = rp)
345 point_list(:, i) = start + direction * t
346 end do
347
348 call this%add_points(point_list)
349 end subroutine read_line
350
360 subroutine read_circle(this, json)
361 class(probes_t), intent(inout) :: this
362 type(json_file), intent(inout) :: json
363
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
370
371 real(kind=rp), dimension(3) :: zero_line, cross_line, temp
372 real(kind=rp) :: pi
373
374 ! Ensure only rank 0 reads the coordinates.
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)
381
382 ! If either center or normal is not of length 3, error out
383 if (size(center) /= 3 .or. size(normal) /= 3) then
384 call neko_error('Invalid center or normal coordinates.')
385 end if
386 if (axis /= 'x' .and. axis /= 'y' .and. axis /= 'z') then
387 call neko_error('Invalid axis.')
388 end if
389 if (radius <= 0) then
390 call neko_error('Invalid radius.')
391 end if
392 if (n_points <= 0) then
393 call neko_error('Invalid number of points.')
394 end if
395
396 ! Normalize the normal vector
397 normal = normal / norm2(normal)
398
399 ! Set the zero line
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]
403
404 if (1.0_rp - dot_product(zero_line, normal) .le. 1e-6) then
405 call neko_error('Invalid axis and normal.')
406 end if
407
408 zero_line = zero_line - dot_product(zero_line, normal) * normal
409 zero_line = zero_line / norm2(zero_line)
410
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)
414
415 ! Calculate the number of points
416 allocate(point_list(3, n_points))
417
418 pi = 4.0_rp * atan(1.0_rp)
419
420 ! Calculate the points
421 do i = 1, n_points
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
424
425 point_list(:, i) = center + radius * temp
426 end do
427
428 call this%add_points(point_list)
429 end subroutine read_circle
430
436 subroutine read_point_zone(this, json, dof)
437 class(probes_t), intent(inout) :: this
438 type(json_file), intent(inout) :: json
439 type(dofmap_t), intent(in) :: dof
440
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
446
447 ! Ensure only rank 0 reads the coordinates.
448 if (pe_rank .ne. 0) return
449
450 call json_get(json, "name", point_zone_name)
451 zone => neko_point_zone_registry%get_point_zone(point_zone_name)
452
453 ! Allocate list of points and reshape the input array
454 allocate(point_list(3, zone%size))
455
456 lx = dof%Xh%lx
457 do i = 1, zone%size
458 idx = zone%mask%get(i)
459
460 nlindex = nonlinear_index(idx, lx, lx, lx)
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))
464
465 point_list(:, i) = [x, y, z]
466 end do
467
468 call this%add_points(point_list)
469 end subroutine read_point_zone
470
471 ! ========================================================================== !
472 ! Supporting routines
473
477 subroutine add_points(this, new_points)
478 class(probes_t), intent(inout) :: this
479 real(kind=rp), dimension(:,:), intent(in) :: new_points
480
481 real(kind=rp), dimension(:,:), allocatable :: temp
482 integer :: n_old, n_new
483
484 ! Get the current number of points
485 n_old = this%n_local_probes
486 n_new = size(new_points, 2)
487
488 ! Move current points to a temporary array
489 if (allocated(this%xyz)) then
490 call move_alloc(this%xyz, temp)
491 end if
492
493 ! Allocate the new array and copy the full list of points
494 allocate(this%xyz(3, n_old + n_new))
495 if (allocated(temp)) then
496 this%xyz(:, 1:n_old) = temp
497 end if
498 this%xyz(:, n_old+1:n_old+n_new) = new_points
499
500 this%n_local_probes = this%n_local_probes + n_new
501 end subroutine add_points
502
503 ! ========================================================================== !
504 ! General initialization routine
505
510 subroutine probes_init_common(this, dof, output_file, name)
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
515
516 character(len=1024) :: header_line
517 real(kind=rp), allocatable :: global_output_coords(:,:)
518 integer :: i, ierr
519 type(matrix_t) :: mat_coords
520
521 this%name = name
522
524 call this%global_interp%find_points_and_redist(this%xyz, &
525 this%n_local_probes)
526
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))
531
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), &
536 this%n_local_probes)
537 end do
538 end if
539
541 call this%fout%init(trim(output_file))
542
543 select type (ft => this%fout%file_type)
544 type is (csv_file_t)
545
546 this%seq_io = .true.
547
548 ! Build the header
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))
552 end do
553 call this%fout%set_header(header_line)
554
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)
567 end if
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)
576 !! Write the data to the file
577 call this%fout%write(mat_coords)
578 end if
579 class default
580 call neko_error("Invalid data. Expected csv_file_t.")
581 end select
582
583 end subroutine probes_init_common
584
586 subroutine probes_free(this)
587 class(probes_t), intent(inout) :: this
588 integer :: i
589
590 if (allocated(this%xyz)) then
591 deallocate(this%xyz)
592 end if
593
594 if (allocated(this%out_values)) then
595 deallocate(this%out_values)
596 end if
597
598 if (allocated(this%out_vals_trsp)) then
599 deallocate(this%out_vals_trsp)
600 end if
601
602 call this%sampled_fields%free()
603
604 if (allocated(this%n_local_probes_tot)) then
605 deallocate(this%n_local_probes_tot)
606 end if
607
608 if (allocated(this%n_local_probes_tot_offset)) then
609 deallocate(this%n_local_probes_tot_offset)
610 end if
611
612 if (allocated(this%global_output_values)) then
613 deallocate(this%global_output_values)
614 end if
615
616 if (allocated(this%which_fields)) then
617 deallocate(this%which_fields)
618 end if
619
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))
624 end if
625 end do
626 deallocate(this%out_values_d)
627 end if
628
629 call this%global_interp%free()
630 call this%mat_out%free()
631
632 end subroutine probes_free
633
635 subroutine probes_show(this)
636 class(probes_t), intent(in) :: this
637 character(len=LOG_SIZE) :: log_buf ! For logging status
638 integer :: i
639
640 ! Probes summary
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)
644
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)
649 end do
650
651 ! Field summary
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)
658 end do
659 call neko_log%end_section()
660 call neko_log%newline()
661
662 end subroutine probes_show
663
665 subroutine probes_debug(this)
666 class(probes_t) :: this
667
668 character(len=LOG_SIZE) :: log_buf ! For logging status
669 integer :: i
670
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)
678 end do
679 end subroutine probes_debug
680
682 subroutine probes_setup_offset(this)
683 class(probes_t) :: this
684 integer :: ierr
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, &
690 0, neko_comm, ierr)
691
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, &
696 0, neko_comm, ierr)
697
698
699
700 end subroutine probes_setup_offset
701
706 subroutine probes_evaluate_and_write(this, time)
707 class(probes_t), intent(inout) :: this
708 type(time_state_t), intent(in) :: time
709 integer :: i, ierr
710 logical :: do_interp_on_host = .false.
711
713 if (time%t .lt. this%start_time) return
714
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, &
719 do_interp_on_host)
720 end do
721
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.)
726 end do
727 end if
728
729 if (this%output_controller%check(time)) then
730 ! Gather all values to rank 0
731 ! If io is only done at root
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)
745 end if
746 else
747 call neko_error('probes sim comp, parallel io need implementation')
748 end if
749
750 !! Register the execution of the activity
751 call this%output_controller%register_execution()
752 end if
753
754 end subroutine probes_evaluate_and_write
755
758 subroutine read_probe_locations(this, xyz, n_local_probes, n_global_probes, &
759 points_file)
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
764
766 type(file_t) :: file_in
767
768 call file_in%init(trim(points_file))
770 select type (ft => file_in%file_type)
771 type is (csv_file_t)
772 call read_xyz_from_csv(xyz, n_local_probes, n_global_probes, ft)
773 this%seq_io = .true.
774 class default
775 call neko_error("Invalid data. Expected csv_file_t.")
776 end select
777
779 call file_free(file_in)
780
781 end subroutine read_probe_locations
782
788 subroutine read_xyz_from_csv(xyz, n_local_probes, n_global_probes, f)
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
793 integer :: n_lines
794
795 n_lines = f%count_lines()
796
797 ! Update the number of probes
798 n_global_probes = n_lines
799
800 ! Initialize the temporal array
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)
806 call f%read(mat_in)
807 call trsp(xyz, 3, mat_in%x, n_global_probes)
808 else
809 n_local_probes = 0
810 allocate(xyz(3, n_local_probes))
811 end if
812
813 end subroutine read_xyz_from_csv
814end module probes
__inline__ __device__ void nonlinear_index(const int idx, const int lx, int *index)
Definition bc_utils.h:44
double real
Map a Fortran array to a device (allocate and associate)
Definition device.F90:77
Copy data between host and device (or device and device)
Definition device.F90:71
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.
Definition case.f90:34
Definition comm.F90:1
type(mpi_datatype), public mpi_real_precision
MPI type for working precision of REAL types.
Definition comm.F90:51
integer, public pe_size
MPI size of communicator.
Definition comm.F90:59
integer, public pe_rank
MPI rank.
Definition comm.F90:56
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:43
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
subroutine, public device_free(x_d)
Deallocate memory on the device.
Definition device.F90:225
integer, parameter, public device_to_host
Definition device.F90:47
Defines a mapping of the degrees of freedom.
Definition dofmap.f90:35
Module for file I/O operations.
Definition file.f90:34
subroutine file_free(this)
File operation destructor.
Definition file.f90:161
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.
Logging routines.
Definition log.f90:34
integer, parameter, public neko_log_debug
Debug log level.
Definition log.f90:89
type(log_t), public neko_log
Global log stream.
Definition log.f90:79
integer, parameter, public log_size
Definition log.f90:45
Defines a matrix.
Definition matrix.f90:34
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public dp
Definition num_types.f90:9
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:437
subroutine read_circle(this, json)
Construct a list of points from a circle.
Definition probes.F90:361
subroutine probes_init_from_components(this, dof, output_file, name, tolerance, padding)
Initialize based on individual parameters.
Definition probes.F90:246
subroutine add_points(this, new_points)
Append a new list of points to the exsiting list.
Definition probes.F90:478
subroutine probes_init_common(this, dof, output_file, name)
Common constructor.
Definition probes.F90:511
subroutine probes_evaluate_and_write(this, time)
Interpolate each probe from its r,s,t coordinates.
Definition probes.F90:707
subroutine read_line(this, json)
Construct a list of points from a line.
Definition probes.F90:317
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:760
subroutine read_point(this, json)
Read a list of points from the json file.
Definition probes.F90:292
subroutine probes_show(this)
Print current probe status, with number of probes and coordinates.
Definition probes.F90:636
subroutine probes_free(this)
Destructor.
Definition probes.F90:587
subroutine probes_debug(this)
Show the status of processor/element owner and error code for each point.
Definition probes.F90:666
subroutine probes_setup_offset(this)
Setup offset for rank 0.
Definition probes.F90:683
subroutine read_file(this, json)
Read a list of points from a csv file.
Definition probes.F90:271
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:789
subroutine probes_init_from_json(this, json, case)
Constructor from json.
Definition probes.F90:134
Defines a registry for storing solution fields.
Definition registry.f90:34
type(registry_t), target, public neko_registry
Global field registry.
Definition registry.f90:144
Simulation components are objects that encapsulate functionality that can be fit to a particular comp...
subroutine compute_(this, time)
Dummy compute function.
Tensor operations.
Definition tensor.f90:61
subroutine, public trsp(a, lda, b, ldb)
Transpose of a rectangular tensor .
Definition tensor.f90:124
Module with things related to the simulation time.
Utilities.
Definition utils.f90:35
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...
Definition file.f90:56
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.