Neko 1.99.2
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
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
50 use tensor, only: trsp
51 use point_zone, only: point_zone_t
53 use file, only : file_t, file_free
54 use csv_file, only : csv_file_t
55 use case, only : case_t
56 use, intrinsic :: iso_c_binding
60 use mpi_f08, only : mpi_allreduce, mpi_integer, mpi_sum, &
61 mpi_double_precision, mpi_gatherv, mpi_gather, mpi_exscan
62 implicit none
63 private
64
65 type, public, extends(simulation_component_t) :: probes_t
67 real(kind=rp) :: start_time
69 integer :: n_fields = 0
70 type(global_interpolation_t) :: global_interp
72 integer :: n_global_probes
74 integer :: n_probes_offset
76 real(kind=rp), allocatable :: xyz(:,:)
78 real(kind=rp), allocatable :: out_values(:,:)
79 type(c_ptr), allocatable :: out_values_d(:)
80 real(kind=rp), allocatable :: out_vals_trsp(:,:)
82 integer :: n_local_probes
84 type(field_list_t) :: sampled_fields
85 character(len=20), allocatable :: which_fields(:)
87 integer, allocatable :: n_local_probes_tot_offset(:)
88 integer, allocatable :: n_local_probes_tot(:)
90 logical :: seq_io
91 real(kind=rp), allocatable :: global_output_values(:,:)
93 type(file_t) :: fout
94 type(matrix_t) :: mat_out
95 contains
97 procedure, pass(this) :: init => probes_init_from_json
98 ! Actual constructor
99 procedure, pass(this) :: init_from_components => &
102 procedure, pass(this) :: free => probes_free
104 procedure, pass(this) :: setup_offset => probes_setup_offset
106 procedure, pass(this) :: compute_ => probes_evaluate_and_write
107
108 ! ----------------------------------------------------------------------- !
109 ! Private methods
110
112 procedure, private, pass(this) :: read_file
114 procedure, private, pass(this) :: read_point
116 procedure, private, pass(this) :: read_line
118 procedure, private, pass(this) :: read_circle
120 procedure, private, pass(this) :: read_point_zone
121
123 procedure, private, pass(this) :: add_points
124 end type probes_t
125
126contains
127
129 subroutine probes_init_from_json(this, json, case)
130 class(probes_t), intent(inout), target :: this
131 type(json_file), intent(inout) :: json
132 class(case_t), intent(inout), target :: case
133 character(len=:), allocatable :: output_file
134 character(len=:), allocatable :: input_file
135 integer :: i, ierr
136
137 ! JSON variables
138 character(len=:), allocatable :: point_type
139 type(json_value), pointer :: json_point_list
140 type(json_file) :: json_point
141 type(json_core) :: core
142 integer :: idx, n_point_children
143
144 ! Initialize the base class
145 call this%free()
146 call this%init_base(json, case)
147
149 call json%info('fields', n_children = this%n_fields)
150 call json_get(json, 'fields', this%which_fields)
151 call json_get(json, 'output_file', output_file)
152 call json_get_or_default(json, 'start_time', this%start_time, -1.0_rp)
153
154 call this%sampled_fields%init(this%n_fields)
155 do i = 1, this%n_fields
156
157 call this%sampled_fields%assign(i, &
158 & neko_registry%get_field(trim(this%which_fields(i))))
159 end do
160
161 ! Setup the required arrays and initialize variables.
162 this%n_local_probes = 0
163 this%n_global_probes = 0
164
165 ! Read the legacy point specification from the points file.
166 if (json%valid_path('points_file')) then
167
168 ! Todo: We should add a deprecation notice here
169 call json_get(json, 'points_file', input_file)
170
171 ! This is distributed as to make it similar to parallel file
172 ! formats later, Reads all into rank 0
173 call read_probe_locations(this, this%xyz, this%n_local_probes, &
174 this%n_global_probes, input_file)
175 end if
176
177 ! Go through the points list and construct the probe list
178 call json%get('points', json_point_list)
179 call json%info('points', n_children = n_point_children)
180
181 do idx = 1, n_point_children
182 call json_extract_item(core, json_point_list, idx, json_point)
183
184 call json_get_or_default(json_point, 'type', point_type, 'none')
185 select case (point_type)
186
187 case ('file')
188 call this%read_file(json_point)
189 case ('points')
190 call this%read_point(json_point)
191 case ('line')
192 call this%read_line(json_point)
193 case ('plane')
194 call neko_error('Plane probes not implemented yet.')
195 case ('circle')
196 call this%read_circle(json_point)
197 case ('point_zone')
198 call this%read_point_zone(json_point, case%fluid%dm_Xh)
199 case ('none')
200 call json_point%print()
201 call neko_error('No point type specified.')
202 case default
203 call neko_error('Unknown region type ' // point_type)
204 end select
205 end do
206
207 call mpi_allreduce(this%n_local_probes, this%n_global_probes, 1, &
208 mpi_integer, mpi_sum, neko_comm, ierr)
209
210 call probes_show(this)
211 call this%init_from_components(case%fluid%dm_Xh, output_file)
212
213 end subroutine probes_init_from_json
214
215 ! ========================================================================== !
216 ! Readers for different point types
217
222 subroutine read_file(this, json)
223 class(probes_t), intent(inout) :: this
224 type(json_file), intent(inout) :: json
225 character(len=:), allocatable :: input_file
226 real(kind=rp), dimension(:,:), allocatable :: point_list
227
228 integer :: n_local, n_global
229
230 if (pe_rank .ne. 0) return
231
232 call json_get(json, 'file_name', input_file)
233
234 call read_probe_locations(this, point_list, n_local, n_global, input_file)
235
236 call this%add_points(point_list)
237 end subroutine read_file
238
243 subroutine read_point(this, json)
244 class(probes_t), intent(inout) :: this
245 type(json_file), intent(inout) :: json
246
247 real(kind=rp), dimension(:,:), allocatable :: point_list
248 real(kind=rp), dimension(:), allocatable :: rp_list_reader
249
250 ! Ensure only rank 0 reads the coordinates.
251 if (pe_rank .ne. 0) return
252 call json_get(json, 'coordinates', rp_list_reader)
253
254 if (mod(size(rp_list_reader), 3) /= 0) then
255 call neko_error('Invalid number of coordinates.')
256 end if
257
258 ! Allocate list of points and reshape the input array
259 allocate(point_list(3, size(rp_list_reader)/3))
260 point_list = reshape(rp_list_reader, [3, size(rp_list_reader)/3])
261
262 call this%add_points(point_list)
263 end subroutine read_point
264
268 subroutine read_line(this, json)
269 class(probes_t), intent(inout) :: this
270 type(json_file), intent(inout) :: json
271
272 real(kind=rp), dimension(:,:), allocatable :: point_list
273 real(kind=rp), dimension(:), allocatable :: start, end
274 real(kind=rp), dimension(3) :: direction
275 real(kind=rp) :: t
276
277 integer :: n_points, i
278
279 ! Ensure only rank 0 reads the coordinates.
280 if (pe_rank .ne. 0) return
281 call json_get(json, "start", start)
282 call json_get(json, "end", end)
283 call json_get(json, "amount", n_points)
284
285 ! If either start or end is not of length 3, error out
286 if (size(start) /= 3 .or. size(end) /= 3) then
287 call neko_error('Invalid start or end coordinates.')
288 end if
289
290 ! Calculate the number of points
291 allocate(point_list(3, n_points))
292
293 ! Calculate the direction vector
294 direction = end - start
295 do i = 1, n_points
296 t = real(i - 1, kind = rp) / real(n_points - 1, kind = rp)
297 point_list(:, i) = start + direction * t
298 end do
299
300 call this%add_points(point_list)
301 end subroutine read_line
302
312 subroutine read_circle(this, json)
313 class(probes_t), intent(inout) :: this
314 type(json_file), intent(inout) :: json
315
316 real(kind=rp), dimension(:,:), allocatable :: point_list
317 real(kind=rp), dimension(:), allocatable :: center, normal
318 real(kind=rp) :: radius
319 real(kind=rp) :: angle
320 integer :: n_points, i
321 character(len=:), allocatable :: axis
322
323 real(kind=rp), dimension(3) :: zero_line, cross_line, temp
324 real(kind=rp) :: pi
325
326 ! Ensure only rank 0 reads the coordinates.
327 if (pe_rank .ne. 0) return
328 call json_get(json, "center", center)
329 call json_get(json, "normal", normal)
330 call json_get(json, "radius", radius)
331 call json_get(json, "amount", n_points)
332 call json_get(json, "axis", axis)
333
334 ! If either center or normal is not of length 3, error out
335 if (size(center) /= 3 .or. size(normal) /= 3) then
336 call neko_error('Invalid center or normal coordinates.')
337 end if
338 if (axis /= 'x' .and. axis /= 'y' .and. axis /= 'z') then
339 call neko_error('Invalid axis.')
340 end if
341 if (radius <= 0) then
342 call neko_error('Invalid radius.')
343 end if
344 if (n_points <= 0) then
345 call neko_error('Invalid number of points.')
346 end if
347
348 ! Normalize the normal vector
349 normal = normal / norm2(normal)
350
351 ! Set the zero line
352 if (axis .eq. 'x') zero_line = [1.0, 0.0, 0.0]
353 if (axis .eq. 'y') zero_line = [0.0, 1.0, 0.0]
354 if (axis .eq. 'z') zero_line = [0.0, 0.0, 1.0]
355
356 if (1.0_rp - dot_product(zero_line, normal) .le. 1e-6) then
357 call neko_error('Invalid axis and normal.')
358 end if
359
360 zero_line = zero_line - dot_product(zero_line, normal) * normal
361 zero_line = zero_line / norm2(zero_line)
362
363 cross_line(1) = normal(2) * zero_line(3) - normal(3) * zero_line(2)
364 cross_line(2) = normal(3) * zero_line(1) - normal(1) * zero_line(3)
365 cross_line(3) = normal(1) * zero_line(2) - normal(2) * zero_line(1)
366
367 ! Calculate the number of points
368 allocate(point_list(3, n_points))
369
370 pi = 4.0_rp * atan(1.0_rp)
371
372 ! Calculate the points
373 do i = 1, n_points
374 angle = 2.0_rp * pi * real(i - 1, kind = rp) / real(n_points, kind = rp)
375 temp = cos(angle) * zero_line + sin(angle) * cross_line
376
377 point_list(:, i) = center + radius * temp
378 end do
379
380 call this%add_points(point_list)
381 end subroutine read_circle
382
388 subroutine read_point_zone(this, json, dof)
389 class(probes_t), intent(inout) :: this
390 type(json_file), intent(inout) :: json
391 type(dofmap_t), intent(in) :: dof
392
393 real(kind=rp), dimension(:,:), allocatable :: point_list
394 character(len=:), allocatable :: point_zone_name
395 class(point_zone_t), pointer :: zone
396 integer :: i, idx, lx, nlindex(4)
397 real(kind=rp) :: x, y, z
398
399 ! Ensure only rank 0 reads the coordinates.
400 if (pe_rank .ne. 0) return
401
402 call json_get(json, "name", point_zone_name)
403 zone => neko_point_zone_registry%get_point_zone(point_zone_name)
404
405 ! Allocate list of points and reshape the input array
406 allocate(point_list(3, zone%size))
407
408 lx = dof%Xh%lx
409 do i = 1, zone%size
410 idx = zone%mask%get(i)
411
412 nlindex = nonlinear_index(idx, lx, lx, lx)
413 x = dof%x(nlindex(1), nlindex(2), nlindex(3), nlindex(4))
414 y = dof%y(nlindex(1), nlindex(2), nlindex(3), nlindex(4))
415 z = dof%z(nlindex(1), nlindex(2), nlindex(3), nlindex(4))
416
417 point_list(:, i) = [x, y, z]
418 end do
419
420 call this%add_points(point_list)
421 end subroutine read_point_zone
422
423 ! ========================================================================== !
424 ! Supporting routines
425
429 subroutine add_points(this, new_points)
430 class(probes_t), intent(inout) :: this
431 real(kind=rp), dimension(:,:), intent(in) :: new_points
432
433 real(kind=rp), dimension(:,:), allocatable :: temp
434 integer :: n_old, n_new
435
436 ! Get the current number of points
437 n_old = this%n_local_probes
438 n_new = size(new_points, 2)
439
440 ! Move current points to a temporary array
441 if (allocated(this%xyz)) then
442 call move_alloc(this%xyz, temp)
443 end if
444
445 ! Allocate the new array and copy the full list of points
446 allocate(this%xyz(3, n_old + n_new))
447 if (allocated(temp)) then
448 this%xyz(:, 1:n_old) = temp
449 end if
450 this%xyz(:, n_old+1:n_old+n_new) = new_points
451
452 this%n_local_probes = this%n_local_probes + n_new
453 end subroutine add_points
454
455 ! ========================================================================== !
456 ! General initialization routine
457
461 subroutine probes_init_from_components(this, dof, output_file)
462 class(probes_t), intent(inout) :: this
463 type(dofmap_t), intent(in) :: dof
464 character(len=:), allocatable, intent(inout) :: output_file
465 character(len=1024) :: header_line
466 real(kind=rp), allocatable :: global_output_coords(:,:)
467 integer :: i, ierr
468 type(matrix_t) :: mat_coords
469
471 call this%global_interp%init(dof)
472
474 call this%global_interp%find_points_and_redist(this%xyz, &
475 this%n_local_probes)
476
478 allocate(this%out_values(this%n_local_probes, this%n_fields))
479 allocate(this%out_values_d(this%n_fields))
480 allocate(this%out_vals_trsp(this%n_fields, this%n_local_probes))
481
482 if (neko_bcknd_device .eq. 1) then
483 do i = 1, this%n_fields
484 this%out_values_d(i) = c_null_ptr
485 call device_map(this%out_values(:,i), this%out_values_d(i), &
486 this%n_local_probes)
487 end do
488 end if
489
491 call this%fout%init(trim(output_file))
492
493 select type (ft => this%fout%file_type)
494 type is (csv_file_t)
495
496 this%seq_io = .true.
497
498 ! Build the header
499 write(header_line, '(I0,A,I0)') this%n_global_probes, ",", this%n_fields
500 do i = 1, this%n_fields
501 header_line = trim(header_line) // "," // trim(this%which_fields(i))
502 end do
503 call this%fout%set_header(header_line)
504
508 allocate(this%n_local_probes_tot(pe_size))
509 allocate(this%n_local_probes_tot_offset(pe_size))
510 call this%setup_offset()
511 if (pe_rank .eq. 0) then
512 allocate(global_output_coords(3, this%n_global_probes))
513 call this%mat_out%init(this%n_global_probes, this%n_fields)
514 allocate(this%global_output_values(this%n_fields, &
515 this%n_global_probes))
516 call mat_coords%init(this%n_global_probes,3)
517 end if
518 call mpi_gatherv(this%xyz, 3*this%n_local_probes, &
519 mpi_real_precision, global_output_coords, &
520 3*this%n_local_probes_tot, &
521 3*this%n_local_probes_tot_offset, &
522 mpi_real_precision, 0, neko_comm, ierr)
523 if (pe_rank .eq. 0) then
524 call trsp(mat_coords%x, this%n_global_probes, &
525 global_output_coords, 3)
526 !! Write the data to the file
527 call this%fout%write(mat_coords)
528 end if
529 class default
530 call neko_error("Invalid data. Expected csv_file_t.")
531 end select
532
533 end subroutine probes_init_from_components
534
536 subroutine probes_free(this)
537 class(probes_t), intent(inout) :: this
538 integer :: i
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 if (allocated(this%which_fields)) then
567 deallocate(this%which_fields)
568 end if
569
570 if (allocated(this%out_values_d)) then
571 do i = 1, size(this%out_values_d)
572 if (c_associated(this%out_values_d(i))) then
573 call device_free(this%out_values_d(i))
574 end if
575 end do
576 deallocate(this%out_values_d)
577 end if
578
579 call this%global_interp%free()
580 call this%mat_out%free()
581
582 end subroutine probes_free
583
585 subroutine probes_show(this)
586 class(probes_t), intent(in) :: this
587 character(len=LOG_SIZE) :: log_buf ! For logging status
588 integer :: i
589
590 ! Probes summary
591 call neko_log%section('Probes')
592 write(log_buf, '(A,I6)') "Number of probes: ", this%n_global_probes
593 call neko_log%message(log_buf)
594
595 call neko_log%message("xyz-coordinates:", lvl = neko_log_debug)
596 do i = 1, this%n_local_probes
597 write(log_buf, '("(",F10.6,",",F10.6,",",F10.6,")")') this%xyz(:,i)
598 call neko_log%message(log_buf, lvl = neko_log_debug)
599 end do
600
601 ! Field summary
602 write(log_buf, '(A,I6)') "Number of fields: ", this%n_fields
603 call neko_log%message(log_buf)
604 do i = 1, this%n_fields
605 write(log_buf, '(A,I6,A,A)') &
606 "Field: ", i, " ", trim(this%which_fields(i))
607 call neko_log%message(log_buf, lvl = neko_log_debug)
608 end do
609 call neko_log%end_section()
610 call neko_log%newline()
611
612 end subroutine probes_show
613
615 subroutine probes_debug(this)
616 class(probes_t) :: this
617
618 character(len=LOG_SIZE) :: log_buf ! For logging status
619 integer :: i
620
621 do i = 1, this%n_local_probes
622 write (log_buf, *) pe_rank, "/", this%global_interp%pe_owner(i), &
623 "/" , this%global_interp%el_owner0(i)
624 call neko_log%message(log_buf)
625 write(log_buf, '(A5,"(",F10.6,",",F10.6,",",F10.6,")")') &
626 "rst: ", this%global_interp%rst(:,i)
627 call neko_log%message(log_buf)
628 end do
629 end subroutine probes_debug
630
632 subroutine probes_setup_offset(this)
633 class(probes_t) :: this
634 integer :: ierr
635 this%n_local_probes_tot = 0
636 this%n_local_probes_tot_offset = 0
637 this%n_probes_offset = 0
638 call mpi_gather(this%n_local_probes, 1, mpi_integer, &
639 this%n_local_probes_tot, 1, mpi_integer, &
640 0, neko_comm, ierr)
641
642 call mpi_exscan(this%n_local_probes, this%n_probes_offset, 1, &
643 mpi_integer, mpi_sum, neko_comm, ierr)
644 call mpi_gather(this%n_probes_offset, 1, mpi_integer, &
645 this%n_local_probes_tot_offset, 1, mpi_integer, &
646 0, neko_comm, ierr)
647
648
649
650 end subroutine probes_setup_offset
651
656 subroutine probes_evaluate_and_write(this, time)
657 class(probes_t), intent(inout) :: this
658 type(time_state_t), intent(in) :: time
659 integer :: i, ierr
660 logical :: do_interp_on_host = .false.
661
663 if (time%t .lt. this%start_time) return
664
666 do i = 1, this%n_fields
667 call this%global_interp%evaluate(this%out_values(:,i), &
668 this%sampled_fields%items(i)%ptr%x, &
669 do_interp_on_host)
670 end do
671
672 if (neko_bcknd_device .eq. 1) then
673 do i = 1, this%n_fields
674 call device_memcpy(this%out_values(:,i), this%out_values_d(i), &
675 this%n_local_probes, device_to_host, sync = .true.)
676 end do
677 end if
678
679 if (this%output_controller%check(time)) then
680 ! Gather all values to rank 0
681 ! If io is only done at root
682 if (this%seq_io) then
683 call trsp(this%out_vals_trsp, this%n_fields, &
684 this%out_values, this%n_local_probes)
685 call mpi_gatherv(this%out_vals_trsp, &
686 this%n_fields*this%n_local_probes, &
687 mpi_real_precision, this%global_output_values, &
688 this%n_fields*this%n_local_probes_tot, &
689 this%n_fields*this%n_local_probes_tot_offset, &
690 mpi_real_precision, 0, neko_comm, ierr)
691 if (pe_rank .eq. 0) then
692 call trsp(this%mat_out%x, this%n_global_probes, &
693 this%global_output_values, this%n_fields)
694 call this%fout%write(this%mat_out, time%t)
695 end if
696 else
697 call neko_error('probes sim comp, parallel io need implementation')
698 end if
699
700 !! Register the execution of the activity
701 call this%output_controller%register_execution()
702 end if
703
704 end subroutine probes_evaluate_and_write
705
708 subroutine read_probe_locations(this, xyz, n_local_probes, n_global_probes, &
709 points_file)
710 class(probes_t), intent(inout) :: this
711 character(len=:), allocatable :: points_file
712 real(kind=rp), allocatable :: xyz(:,:)
713 integer, intent(inout) :: n_local_probes, n_global_probes
714
716 type(file_t) :: file_in
717
718 call file_in%init(trim(points_file))
720 select type (ft => file_in%file_type)
721 type is (csv_file_t)
722 call read_xyz_from_csv(xyz, n_local_probes, n_global_probes, ft)
723 this%seq_io = .true.
724 class default
725 call neko_error("Invalid data. Expected csv_file_t.")
726 end select
727
729 call file_free(file_in)
730
731 end subroutine read_probe_locations
732
738 subroutine read_xyz_from_csv(xyz, n_local_probes, n_global_probes, f)
739 type(csv_file_t), intent(inout) :: f
740 real(kind=rp), allocatable :: xyz(:,:)
741 integer, intent(inout) :: n_local_probes, n_global_probes
742 type(matrix_t) :: mat_in, mat_in2
743 integer :: n_lines
744
745 n_lines = f%count_lines()
746
747 ! Update the number of probes
748 n_global_probes = n_lines
749
750 ! Initialize the temporal array
751 if (pe_rank .eq. 0) then
752 n_local_probes = n_global_probes
753 allocate(xyz(3, n_local_probes))
754 call mat_in%init(n_global_probes,3)
755 call mat_in2%init(3, n_global_probes)
756 call f%read(mat_in)
757 call trsp(xyz, 3, mat_in%x, n_global_probes)
758 else
759 n_local_probes = 0
760 allocate(xyz(3, n_local_probes))
761 end if
762
763 end subroutine read_xyz_from_csv
764end 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:219
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:156
Implements global_interpolation given a dofmap.
Utilities for retrieving parameters from the case files.
Logging routines.
Definition log.f90:34
integer, parameter, public neko_log_debug
Debug log level.
Definition log.f90:84
type(log_t), public neko_log
Global log stream.
Definition log.f90:76
integer, parameter, public log_size
Definition log.f90:46
Defines a matrix.
Definition matrix.f90:34
Build configurations.
integer, parameter neko_bcknd_device
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:389
subroutine read_circle(this, json)
Construct a list of points from a circle.
Definition probes.F90:313
subroutine add_points(this, new_points)
Append a new list of points to the exsiting list.
Definition probes.F90:430
subroutine probes_evaluate_and_write(this, time)
Interpolate each probe from its r,s,t coordinates.
Definition probes.F90:657
subroutine read_line(this, json)
Construct a list of points from a line.
Definition probes.F90:269
subroutine probes_init_from_components(this, dof, output_file)
Initialize without json things.
Definition probes.F90:462
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:710
subroutine read_point(this, json)
Read a list of points from the json file.
Definition probes.F90:244
subroutine probes_show(this)
Print current probe status, with number of probes and coordinates.
Definition probes.F90:586
subroutine probes_free(this)
Destructor.
Definition probes.F90:537
subroutine probes_debug(this)
Show the status of processor/element owner and error code for each point.
Definition probes.F90:616
subroutine probes_setup_offset(this)
Setup offset for rank 0.
Definition probes.F90:633
subroutine read_file(this, json)
Read a list of points from a csv file.
Definition probes.F90:223
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:739
subroutine probes_init_from_json(this, json, case)
Constructor from json.
Definition probes.F90:130
Defines a registry for storing solution fields.
Definition registry.f90:34
type(registry_t), target, public neko_registry
Global field registry.
Definition registry.f90:128
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:55
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.