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