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
40 use vector, only : vector_t
43 use field_list, only : field_list_t
44 use time_state, only : time_state_t
46 use registry, only : neko_registry
47 use dofmap, only : dofmap_t
48 use json_module, only : json_file, json_value, json_core
53 use tensor, only : trsp
54 use point_zone, only : point_zone_t
56 use file, only : file_t, file_free
57 use csv_file, only : csv_file_t
58 use hdf5_file, only : hdf5_file_t
59 use math, only : copy
60 use device_math, only : device_copy
61 use case, only : case_t
62 use, intrinsic :: iso_c_binding
66 use mpi_f08, only : mpi_allreduce, mpi_integer, mpi_sum, &
67 mpi_double_precision, mpi_gatherv, mpi_gather, mpi_exscan
68 implicit none
69 private
70
71 type, public, extends(simulation_component_t) :: probes_t
73 real(kind=rp) :: start_time
75 integer :: n_fields = 0
76 type(global_interpolation_t) :: global_interp
78 integer :: n_global_probes
80 integer :: n_probes_offset
82 real(kind=rp), allocatable :: xyz(:,:)
84 real(kind=rp), allocatable :: out_values(:,:)
85 type(c_ptr), allocatable :: out_values_d(:)
86 real(kind=rp), allocatable :: out_vals_trsp(:,:)
88 integer :: n_local_probes
90 type(field_list_t) :: sampled_fields
91 character(len=NEKO_VARNAME_LEN), allocatable :: which_fields(:)
93 integer, allocatable :: n_local_probes_tot_offset(:)
94 integer, allocatable :: n_local_probes_tot(:)
96 logical :: seq_io
97 real(kind=rp), allocatable :: global_output_values(:,:)
99 type(file_t) :: fout
100 type(matrix_t) :: mat_out
101 type(vector_t) :: vec_out
102 logical :: append_out = .true.
103 contains
105 procedure, pass(this) :: init => probes_init_from_json
107 procedure, pass(this) :: init_from_components => &
110 procedure, private, pass(this) :: init_common => probes_init_common
112 procedure, pass(this) :: free => probes_free
114 procedure, pass(this) :: setup_offset => probes_setup_offset
116 procedure, pass(this) :: compute_ => probes_evaluate_and_write
117
118 ! ----------------------------------------------------------------------- !
119 ! Private methods
120
122 procedure, private, pass(this) :: read_file
124 procedure, private, pass(this) :: read_point
126 procedure, private, pass(this) :: read_line
128 procedure, private, pass(this) :: read_circle
130 procedure, private, pass(this) :: read_point_zone
131
133 procedure, private, pass(this) :: add_points
134 end type probes_t
135
136contains
137
139 subroutine probes_init_from_json(this, json, case)
140 class(probes_t), intent(inout), target :: this
141 type(json_file), intent(inout) :: json
142 class(case_t), intent(inout), target :: case
143 character(len=:), allocatable :: output_file
144 character(len=:), allocatable :: input_file
145 integer :: i, ierr
146
147 ! JSON variables
148 character(len=:), allocatable :: point_type
149 type(json_value), pointer :: json_point_list
150 type(json_file) :: json_point
151 type(json_core) :: core
152 integer :: idx, n_point_children
153
154 character(len=:), allocatable :: name
155
156 ! Initialize the base class
157 call this%free()
158
159 call json_get_or_default(json, "name", name, "probes")
160 call this%init_base(json, case)
161
163 call json%info('fields', n_children = this%n_fields)
164 call json_get(json, 'fields', this%which_fields)
165 call json_get(json, 'output_file', output_file)
166 call json_get_or_lookup_or_default(json, 'start_time', this%start_time, &
167 -1.0_rp)
168 call json_get_or_default(json, 'append_output', this%append_out, .true.)
169
170 call this%sampled_fields%init(this%n_fields)
171 do i = 1, this%n_fields
172
173 call this%sampled_fields%assign(i, &
174 & neko_registry%get_field(trim(this%which_fields(i))))
175 end do
176
177 ! Setup the required arrays and initialize variables.
178 this%n_local_probes = 0
179 this%n_global_probes = 0
180
181 ! Read the legacy point specification from the points file.
182 if (json%valid_path('points_file')) then
183
184 ! Todo: We should add a deprecation notice here
185 call json_get(json, 'points_file', input_file)
186
187 ! This is distributed as to make it similar to parallel file
188 ! formats later, Reads all into rank 0
189 call read_probe_locations(this, this%xyz, this%n_local_probes, &
190 this%n_global_probes, input_file)
191 end if
192
193 if (json%valid_path('points')) then
194
195 ! Go through the points list and construct the probe list
196 call json%get('points', json_point_list)
197 call json%info('points', n_children = n_point_children)
198
199 do idx = 1, n_point_children
200 call json_extract_item(core, json_point_list, idx, json_point)
201
202 call json_get_or_default(json_point, 'type', point_type, 'none')
203 select case (point_type)
204
205 case ('file')
206 call this%read_file(json_point)
207 case ('points')
208 call this%read_point(json_point)
209 case ('line')
210 call this%read_line(json_point)
211 case ('plane')
212 call neko_error('Plane probes not implemented yet.')
213 case ('circle')
214 call this%read_circle(json_point)
215 case ('point_zone')
216 call this%read_point_zone(json_point, case%fluid%dm_Xh)
217 case ('none')
218 call json_point%print()
219 call neko_error('No point type specified.')
220 case default
221 call neko_error('Unknown region type ' // point_type)
222 end select
223 end do
224
225 end if
226
227 call mpi_allreduce(this%n_local_probes, this%n_global_probes, 1, &
228 mpi_integer, mpi_sum, neko_comm, ierr)
229
230 call probes_show(this)
231
232 ! Get interpolation parameters from json
233 block
234 type(json_file) :: interp_subdict
235 call json_get_subdict_or_empty(json, "interpolation", interp_subdict)
236 call this%global_interp%init(case%fluid%dm_Xh, &
237 params_subdict = interp_subdict)
238
239 call this%init_common(case%fluid%dm_Xh, output_file, name)
240
241 end block
242
243 end subroutine probes_init_from_json
244
251 subroutine probes_init_from_components(this, dof, output_file, name, &
252 tolerance, padding)
253 class(probes_t), intent(inout) :: this
254 type(dofmap_t), intent(in) :: dof
255 character(len=:), allocatable, intent(inout) :: output_file
256 character(len=*), intent(in) :: name
257 real(kind=dp), intent(in), optional :: tolerance, padding
258
259 character(len=1024) :: header_line
260 real(kind=rp), allocatable :: global_output_coords(:,:)
261 integer :: i, ierr
262 type(matrix_t) :: mat_coords
263
264 call this%global_interp%init(dof, tol = tolerance, pad = padding)
265
266 call this%init_common(dof, output_file, name)
267
268 end subroutine probes_init_from_components
269
270 ! ========================================================================== !
271 ! Readers for different point types
272
277 subroutine read_file(this, json)
278 class(probes_t), intent(inout) :: this
279 type(json_file), intent(inout) :: json
280 character(len=:), allocatable :: input_file
281 real(kind=rp), dimension(:,:), allocatable :: point_list
282
283 integer :: n_local, n_global
284
285 if (pe_rank .ne. 0) return
286
287 call json_get(json, 'file_name', input_file)
288
289 call read_probe_locations(this, point_list, n_local, n_global, input_file)
290
291 call this%add_points(point_list)
292 end subroutine read_file
293
298 subroutine read_point(this, json)
299 class(probes_t), intent(inout) :: this
300 type(json_file), intent(inout) :: json
301
302 real(kind=rp), dimension(:,:), allocatable :: point_list
303 real(kind=rp), dimension(:), allocatable :: rp_list_reader
304
305 ! Ensure only rank 0 reads the coordinates.
306 if (pe_rank .ne. 0) return
307 call json_get_or_lookup(json, 'coordinates', rp_list_reader)
308
309 if (mod(size(rp_list_reader), 3) /= 0) then
310 call neko_error('Invalid number of coordinates.')
311 end if
312
313 ! Allocate list of points and reshape the input array
314 allocate(point_list(3, size(rp_list_reader)/3))
315 point_list = reshape(rp_list_reader, [3, size(rp_list_reader)/3])
316
317 call this%add_points(point_list)
318 end subroutine read_point
319
323 subroutine read_line(this, json)
324 class(probes_t), intent(inout) :: this
325 type(json_file), intent(inout) :: json
326
327 real(kind=rp), dimension(:,:), allocatable :: point_list
328 real(kind=rp), dimension(:), allocatable :: start, end
329 real(kind=rp), dimension(3) :: direction
330 real(kind=rp) :: t
331
332 integer :: n_points, i
333
334 ! Ensure only rank 0 reads the coordinates.
335 if (pe_rank .ne. 0) return
336 call json_get_or_lookup(json, "start", start)
337 call json_get_or_lookup(json, "end", end)
338 call json_get_or_lookup(json, "amount", n_points)
339
340 ! If either start or end is not of length 3, error out
341 if (size(start) /= 3 .or. size(end) /= 3) then
342 call neko_error('Invalid start or end coordinates.')
343 end if
344
345 ! Calculate the number of points
346 allocate(point_list(3, n_points))
347
348 ! Calculate the direction vector
349 direction = end - start
350 do i = 1, n_points
351 t = real(i - 1, kind = rp) / real(n_points - 1, kind = rp)
352 point_list(:, i) = start + direction * t
353 end do
354
355 call this%add_points(point_list)
356 end subroutine read_line
357
367 subroutine read_circle(this, json)
368 class(probes_t), intent(inout) :: this
369 type(json_file), intent(inout) :: json
370
371 real(kind=rp), dimension(:,:), allocatable :: point_list
372 real(kind=rp), dimension(:), allocatable :: center, normal
373 real(kind=rp) :: radius
374 real(kind=rp) :: angle
375 integer :: n_points, i
376 character(len=:), allocatable :: axis
377
378 real(kind=rp), dimension(3) :: zero_line, cross_line, temp
379 real(kind=rp) :: pi
380
381 ! Ensure only rank 0 reads the coordinates.
382 if (pe_rank .ne. 0) return
383 call json_get_or_lookup(json, "center", center)
384 call json_get_or_lookup(json, "normal", normal)
385 call json_get_or_lookup(json, "radius", radius)
386 call json_get_or_lookup(json, "amount", n_points)
387 call json_get(json, "axis", axis)
388
389 ! If either center or normal is not of length 3, error out
390 if (size(center) /= 3 .or. size(normal) /= 3) then
391 call neko_error('Invalid center or normal coordinates.')
392 end if
393 if (axis /= 'x' .and. axis /= 'y' .and. axis /= 'z') then
394 call neko_error('Invalid axis.')
395 end if
396 if (radius <= 0) then
397 call neko_error('Invalid radius.')
398 end if
399 if (n_points <= 0) then
400 call neko_error('Invalid number of points.')
401 end if
402
403 ! Normalize the normal vector
404 normal = normal / norm2(normal)
405
406 ! Set the zero line
407 if (axis .eq. 'x') zero_line = [1.0, 0.0, 0.0]
408 if (axis .eq. 'y') zero_line = [0.0, 1.0, 0.0]
409 if (axis .eq. 'z') zero_line = [0.0, 0.0, 1.0]
410
411 if (1.0_rp - dot_product(zero_line, normal) .le. 1e-6) then
412 call neko_error('Invalid axis and normal.')
413 end if
414
415 zero_line = zero_line - dot_product(zero_line, normal) * normal
416 zero_line = zero_line / norm2(zero_line)
417
418 cross_line(1) = normal(2) * zero_line(3) - normal(3) * zero_line(2)
419 cross_line(2) = normal(3) * zero_line(1) - normal(1) * zero_line(3)
420 cross_line(3) = normal(1) * zero_line(2) - normal(2) * zero_line(1)
421
422 ! Calculate the number of points
423 allocate(point_list(3, n_points))
424
425 pi = 4.0_rp * atan(1.0_rp)
426
427 ! Calculate the points
428 do i = 1, n_points
429 angle = 2.0_rp * pi * real(i - 1, kind = rp) / real(n_points, kind = rp)
430 temp = cos(angle) * zero_line + sin(angle) * cross_line
431
432 point_list(:, i) = center + radius * temp
433 end do
434
435 call this%add_points(point_list)
436 end subroutine read_circle
437
443 subroutine read_point_zone(this, json, dof)
444 class(probes_t), intent(inout) :: this
445 type(json_file), intent(inout) :: json
446 type(dofmap_t), intent(in) :: dof
447
448 real(kind=rp), dimension(:,:), allocatable :: point_list
449 character(len=:), allocatable :: point_zone_name
450 class(point_zone_t), pointer :: zone
451 integer :: i, idx, lx, nlindex(4)
452 real(kind=rp) :: x, y, z
453
454 ! Ensure only rank 0 reads the coordinates.
455 if (pe_rank .ne. 0) return
456
457 call json_get(json, "name", point_zone_name)
458 zone => neko_point_zone_registry%get_point_zone(point_zone_name)
459
460 ! Allocate list of points and reshape the input array
461 allocate(point_list(3, zone%size))
462
463 lx = dof%Xh%lx
464 do i = 1, zone%size
465 idx = zone%mask%get(i)
466
467 nlindex = nonlinear_index(idx, lx, lx, lx)
468 x = dof%x(nlindex(1), nlindex(2), nlindex(3), nlindex(4))
469 y = dof%y(nlindex(1), nlindex(2), nlindex(3), nlindex(4))
470 z = dof%z(nlindex(1), nlindex(2), nlindex(3), nlindex(4))
471
472 point_list(:, i) = [x, y, z]
473 end do
474
475 call this%add_points(point_list)
476 end subroutine read_point_zone
477
478 ! ========================================================================== !
479 ! Supporting routines
480
484 subroutine add_points(this, new_points)
485 class(probes_t), intent(inout) :: this
486 real(kind=rp), dimension(:,:), intent(in) :: new_points
487
488 real(kind=rp), dimension(:,:), allocatable :: temp
489 integer :: n_old, n_new
490
491 ! Get the current number of points
492 n_old = this%n_local_probes
493 n_new = size(new_points, 2)
494
495 ! Move current points to a temporary array
496 if (allocated(this%xyz)) then
497 call move_alloc(this%xyz, temp)
498 end if
499
500 ! Allocate the new array and copy the full list of points
501 allocate(this%xyz(3, n_old + n_new))
502 if (allocated(temp)) then
503 this%xyz(:, 1:n_old) = temp
504 end if
505 this%xyz(:, n_old+1:n_old+n_new) = new_points
506
507 this%n_local_probes = this%n_local_probes + n_new
508 end subroutine add_points
509
510 ! ========================================================================== !
511 ! General initialization routine
512
517 subroutine probes_init_common(this, dof, output_file, name)
518 class(probes_t), intent(inout) :: this
519 type(dofmap_t), intent(in) :: dof
520 character(len=:), allocatable, intent(inout) :: output_file
521 character(len=*), intent(in) :: name
522
523 character(len=1024) :: header_line
524 real(kind=rp), allocatable :: global_output_coords(:,:)
525 integer :: i, ierr, out_int
526 type(matrix_t) :: mat_coords
527 logical :: attr_exist = .false.
528
529 this%name = name
530
532 call this%global_interp%find_points_and_redist(this%xyz, &
533 this%n_local_probes)
534
536 allocate(this%out_values(this%n_local_probes, this%n_fields))
537 allocate(this%out_values_d(this%n_fields))
538 allocate(this%out_vals_trsp(this%n_fields, this%n_local_probes))
539
540 if (neko_bcknd_device .eq. 1) then
541 do i = 1, this%n_fields
542 this%out_values_d(i) = c_null_ptr
543 call device_map(this%out_values(:,i), this%out_values_d(i), &
544 this%n_local_probes)
545 end do
546 end if
547
549 call this%fout%init(this%case%output_directory // trim(output_file))
550
551 select type (ft => this%fout%file_type)
552 type is (csv_file_t)
553
554 this%seq_io = .true.
555
556 ! Build the header
557 write(header_line, '(I0,A,I0)') this%n_global_probes, ",", this%n_fields
558 do i = 1, this%n_fields
559 header_line = trim(header_line) // "," // trim(this%which_fields(i))
560 end do
561 call this%fout%set_header(header_line)
562
566 allocate(this%n_local_probes_tot(pe_size))
567 allocate(this%n_local_probes_tot_offset(pe_size))
568 call this%setup_offset()
569 if (pe_rank .eq. 0) then
570 allocate(global_output_coords(3, this%n_global_probes))
571 call this%mat_out%init(this%n_global_probes, this%n_fields)
572 allocate(this%global_output_values(this%n_fields, &
573 this%n_global_probes))
574 call mat_coords%init(this%n_global_probes,3)
575 end if
576 call mpi_gatherv(this%xyz, 3*this%n_local_probes, &
577 mpi_real_precision, global_output_coords, &
578 3*this%n_local_probes_tot, &
579 3*this%n_local_probes_tot_offset, &
580 mpi_real_precision, 0, neko_comm, ierr)
581 if (pe_rank .eq. 0) then
582 call trsp(mat_coords%x, this%n_global_probes, &
583 global_output_coords, 3)
584 !! Write the data to the file
585 call this%fout%write(mat_coords)
586 call mat_coords%free()
587 end if
588 class is (hdf5_file_t)
589
591 if (this%append_out) then
592 call ft%set_overwrite(.false.)
593 else
594 call ft%set_overwrite(.true.)
595 end if
596 call mat_coords%init(3, this%n_local_probes, "coordinates")
597 call copy(mat_coords%x, this%xyz, 3*this%n_local_probes)
598
600 call ft%open("w")
601 call ft%set_active_group("probes")
602
603 ! Check if the NSteps attribute already exists
604 call ft%read_attribute("NSteps", out_int, attr_exist)
605 if (attr_exist) then
606 ! If the attribute exists,
607 ! do not write the coordinates but register the executions
608 this%output_controller%nexecutions = out_int
609 else
610 ! Write out the mesh
611 call ft%write_dataset(mat_coords)
612 out_int = this%n_global_probes
613 call ft%write_attribute("NProbes", out_int)
614 end if
615 call ft%close()
616
618 this%seq_io = .false.
619 call this%vec_out%init(this%n_local_probes, "interpolated_fields_trsp")
620
621 ! Free temporaries
622 call mat_coords%free()
623
624 class default
625 call neko_error("Invalid data. Expected csv_file_t or hdf5_file_t.")
626 end select
627
628 end subroutine probes_init_common
629
631 subroutine probes_free(this)
632 class(probes_t), intent(inout) :: this
633 integer :: i
634
635 if (allocated(this%xyz)) then
636 deallocate(this%xyz)
637 end if
638
639 if (allocated(this%out_values)) then
640 deallocate(this%out_values)
641 end if
642
643 if (allocated(this%out_vals_trsp)) then
644 deallocate(this%out_vals_trsp)
645 end if
646
647 call this%sampled_fields%free()
648
649 if (allocated(this%n_local_probes_tot)) then
650 deallocate(this%n_local_probes_tot)
651 end if
652
653 if (allocated(this%n_local_probes_tot_offset)) then
654 deallocate(this%n_local_probes_tot_offset)
655 end if
656
657 if (allocated(this%global_output_values)) then
658 deallocate(this%global_output_values)
659 end if
660
661 if (allocated(this%which_fields)) then
662 deallocate(this%which_fields)
663 end if
664
665 if (allocated(this%out_values_d)) then
666 do i = 1, size(this%out_values_d)
667 if (c_associated(this%out_values_d(i))) then
668 call device_free(this%out_values_d(i))
669 end if
670 end do
671 deallocate(this%out_values_d)
672 end if
673
674 call this%global_interp%free()
675 call this%mat_out%free()
676 call this%vec_out%free()
677
678 end subroutine probes_free
679
681 subroutine probes_show(this)
682 class(probes_t), intent(in) :: this
683 character(len=LOG_SIZE) :: log_buf ! For logging status
684 integer :: i
685
686 ! Probes summary
687 call neko_log%section('Probes')
688 write(log_buf, '(A,I6)') "Number of probes: ", this%n_global_probes
689 call neko_log%message(log_buf)
690
691 call neko_log%message("xyz-coordinates:", lvl = neko_log_debug)
692 do i = 1, this%n_local_probes
693 write(log_buf, '("(",F10.6,",",F10.6,",",F10.6,")")') this%xyz(:,i)
694 call neko_log%message(log_buf, lvl = neko_log_debug)
695 end do
696
697 ! Field summary
698 write(log_buf, '(A,I6)') "Number of fields: ", this%n_fields
699 call neko_log%message(log_buf)
700 do i = 1, this%n_fields
701 write(log_buf, '(A,I6,A,A)') &
702 "Field: ", i, " ", trim(this%which_fields(i))
703 call neko_log%message(log_buf, lvl = neko_log_debug)
704 end do
705 call neko_log%end_section()
706 call neko_log%newline()
707
708 end subroutine probes_show
709
711 subroutine probes_debug(this)
712 class(probes_t) :: this
713
714 character(len=LOG_SIZE) :: log_buf ! For logging status
715 integer :: i
716
717 do i = 1, this%n_local_probes
718 write (log_buf, *) pe_rank, "/", this%global_interp%pe_owner(i), &
719 "/" , this%global_interp%el_owner0(i)
720 call neko_log%message(log_buf)
721 write(log_buf, '(A5,"(",F10.6,",",F10.6,",",F10.6,")")') &
722 "rst: ", this%global_interp%rst(:,i)
723 call neko_log%message(log_buf)
724 end do
725 end subroutine probes_debug
726
728 subroutine probes_setup_offset(this)
729 class(probes_t) :: this
730 integer :: ierr
731 this%n_local_probes_tot = 0
732 this%n_local_probes_tot_offset = 0
733 this%n_probes_offset = 0
734 call mpi_gather(this%n_local_probes, 1, mpi_integer, &
735 this%n_local_probes_tot, 1, mpi_integer, &
736 0, neko_comm, ierr)
737
738 call mpi_exscan(this%n_local_probes, this%n_probes_offset, 1, &
739 mpi_integer, mpi_sum, neko_comm, ierr)
740 call mpi_gather(this%n_probes_offset, 1, mpi_integer, &
741 this%n_local_probes_tot_offset, 1, mpi_integer, &
742 0, neko_comm, ierr)
743
744
745
746 end subroutine probes_setup_offset
747
752 subroutine probes_evaluate_and_write(this, time)
753 class(probes_t), intent(inout) :: this
754 type(time_state_t), intent(in) :: time
755 integer :: i, ierr
756 logical :: do_interp_on_host = .false.
757 character(len=1000) :: group_name
758 real(kind=rp) :: time_
759 type(vector_t) :: vec_time
760 integer :: out_int
761
763 if (time%t .lt. this%start_time) return
764
766 do i = 1, this%n_fields
767 call this%global_interp%evaluate(this%out_values(:,i), &
768 this%sampled_fields%items(i)%ptr%x, &
769 do_interp_on_host)
770 end do
771
772 if (neko_bcknd_device .eq. 1) then
773 do i = 1, this%n_fields
774 call device_memcpy(this%out_values(:,i), this%out_values_d(i), &
775 this%n_local_probes, device_to_host, sync = .true.)
776 end do
777 end if
778
779 if (this%output_controller%check(time)) then
780 select type (ft => this%fout%file_type)
781 type is (csv_file_t)
782 ! Gather all values to rank 0
783 ! If io is only done at root
784 if (this%seq_io) then
785 call trsp(this%out_vals_trsp, this%n_fields, &
786 this%out_values, this%n_local_probes)
787 call mpi_gatherv(this%out_vals_trsp, &
788 this%n_fields*this%n_local_probes, &
789 mpi_real_precision, this%global_output_values, &
790 this%n_fields*this%n_local_probes_tot, &
791 this%n_fields*this%n_local_probes_tot_offset, &
792 mpi_real_precision, 0, neko_comm, ierr)
793 if (pe_rank .eq. 0) then
794 call trsp(this%mat_out%x, this%n_global_probes, &
795 this%global_output_values, this%n_fields)
796 call this%fout%write(this%mat_out, time%t)
797 end if
798 else
799 call neko_error("CSV outputs only works sequentially")
800 end if
801
802 type is (hdf5_file_t)
803
804 if (this%seq_io) then
805 call neko_error("HDF5 outputs only works in parallel currently.")
806
807 else
808
809 ! Append output format
810 if (this%append_out) then
811
812 call ft%open("w")
813 call ft%set_active_group("probes")
814 ! Write Nsteps in root
815 out_int = this%output_controller%nexecutions + 1
816 call ft%write_attribute("NSteps", out_int)
817 ! Write out the data
818 do i = 1, this%n_fields
819 call copy(this%vec_out%x, this%out_values(:,i), &
820 this%vec_out%size())
821 this%vec_out%name = trim(this%which_fields(i))
822 call ft%write_dataset(this%vec_out)
823 end do
824
825 ! Write the time by hacking the vector write
826 if (pe_rank .eq. 0) then
827 call vec_time%init(1, "time")
828 vec_time%x(1) = time%t
829 else
830 call vec_time%init(0, "time")
831 end if
832 call ft%write_dataset(vec_time)
833 call vec_time%free()
834 call ft%close()
835
836 ! Write data in different steps
837 else
838
839 out_int = this%output_controller%nexecutions + 1
840 ! Set up the name
841 write(group_name, '(A,I0)') "probes/Step_", out_int
842 call ft%open("w")
843 call ft%set_active_group("probes")
844 ! Write Nsteps in root
845 call ft%write_attribute("NSteps", out_int)
846 ! Write out the data
847 call ft%set_active_group(trim(group_name))
848 do i = 1, this%n_fields
849 call copy(this%vec_out%x, this%out_values(:,i), &
850 this%vec_out%size())
851 this%vec_out%name = trim(this%which_fields(i))
852 call ft%write_dataset(this%vec_out)
853 end do
854 ! Write the time as an attribute
855 time_ = time%t
856 call ft%write_attribute("time", time_)
857 call ft%close()
858 end if
859
860 end if
861 class default
862 call neko_error("Invalid data. Expected csv_file_t or hdf5_file_t.")
863 end select
864
865 !! Register the execution of the activity
866 call this%output_controller%register_execution()
867 end if
868
869 end subroutine probes_evaluate_and_write
870
873 subroutine read_probe_locations(this, xyz, n_local_probes, n_global_probes, &
874 points_file)
875 class(probes_t), intent(inout) :: this
876 character(len=:), allocatable :: points_file
877 real(kind=rp), allocatable :: xyz(:,:)
878 integer, intent(inout) :: n_local_probes, n_global_probes
879 type(matrix_t) :: mat_in
880 integer :: ierr
881
883 type(file_t) :: file_in
884
885 call file_in%init(trim(points_file))
887 select type (ft => file_in%file_type)
888 type is (csv_file_t)
889 call read_xyz_from_csv(xyz, n_local_probes, n_global_probes, ft)
890 this%seq_io = .true.
891 type is (hdf5_file_t)
892 call ft%open("r")
893 call ft%read_dataset("xyz", mat_in, "rank_0")
894 call ft%close()
895
896 ! Copy the data to the xyz location
897 n_local_probes = mat_in%get_ncols()
898 call mpi_allreduce(n_local_probes, n_global_probes, 1, mpi_integer, &
899 mpi_sum, neko_comm, ierr)
900
901 allocate(xyz(3, n_local_probes)) ! We asume that axis 1 has 3 entries
902 call copy(xyz, mat_in%x, 3*n_local_probes)
903 call mat_in%free()
904 class default
905 call neko_error("Invalid data. Expected csv_file_t.")
906 end select
907
909 call file_free(file_in)
910
911 end subroutine read_probe_locations
912
918 subroutine read_xyz_from_csv(xyz, n_local_probes, n_global_probes, f)
919 type(csv_file_t), intent(inout) :: f
920 real(kind=rp), allocatable :: xyz(:,:)
921 integer, intent(inout) :: n_local_probes, n_global_probes
922 type(matrix_t) :: mat_in, mat_in2
923 integer :: n_lines
924
925 n_lines = f%count_lines()
926
927 ! Update the number of probes
928 n_global_probes = n_lines
929
930 ! Initialize the temporal array
931 if (pe_rank .eq. 0) then
932 n_local_probes = n_global_probes
933 allocate(xyz(3, n_local_probes))
934 call mat_in%init(n_global_probes,3)
935 call mat_in2%init(3, n_global_probes)
936 call f%read(mat_in)
937 call trsp(xyz, 3, mat_in%x, n_global_probes)
938 else
939 n_local_probes = 0
940 allocate(xyz(3, n_local_probes))
941 end if
942
943 end subroutine read_xyz_from_csv
944end 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
subroutine, public device_copy(a_d, b_d, n, strm)
Copy a vector .
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:160
Implements global_interpolation given a dofmap.
HDF5 file format.
Definition hdf5_file.F90:34
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
Definition math.f90:60
subroutine, public copy(a, b, n)
Copy a vector .
Definition math.f90:252
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:444
subroutine read_circle(this, json)
Construct a list of points from a circle.
Definition probes.F90:368
subroutine probes_init_from_components(this, dof, output_file, name, tolerance, padding)
Initialize based on individual parameters.
Definition probes.F90:253
subroutine add_points(this, new_points)
Append a new list of points to the exsiting list.
Definition probes.F90:485
subroutine probes_init_common(this, dof, output_file, name)
Common constructor.
Definition probes.F90:518
subroutine probes_evaluate_and_write(this, time)
Interpolate each probe from its r,s,t coordinates.
Definition probes.F90:753
subroutine read_line(this, json)
Construct a list of points from a line.
Definition probes.F90:324
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:875
subroutine read_point(this, json)
Read a list of points from the json file.
Definition probes.F90:299
subroutine probes_show(this)
Print current probe status, with number of probes and coordinates.
Definition probes.F90:682
subroutine probes_free(this)
Destructor.
Definition probes.F90:632
subroutine probes_debug(this)
Show the status of processor/element owner and error code for each point.
Definition probes.F90:712
subroutine probes_setup_offset(this)
Setup offset for rank 0.
Definition probes.F90:729
subroutine read_file(this, json)
Read a list of points from a csv file.
Definition probes.F90:278
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:919
subroutine probes_init_from_json(this, json, case)
Constructor from json.
Definition probes.F90:140
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
integer, parameter, public neko_varname_len
Definition utils.f90:42
Defines a vector.
Definition vector.f90:34
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.
Interface for HDF5 files.
Definition hdf5_file.F90:60
Base abstract type for point zones.
Base abstract class for simulation components.
A struct that contains all info about the time, expand as needed.