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, achar(44), &
558 this%n_fields
559 do i = 1, this%n_fields
560 header_line = trim(header_line) // achar(44) // &
561 trim(this%which_fields(i))
562 end do
563 call this%fout%set_header(header_line)
564
568 allocate(this%n_local_probes_tot(pe_size))
569 allocate(this%n_local_probes_tot_offset(pe_size))
570 call this%setup_offset()
571 if (pe_rank .eq. 0) then
572 allocate(global_output_coords(3, this%n_global_probes))
573 call this%mat_out%init(this%n_global_probes, this%n_fields)
574 allocate(this%global_output_values(this%n_fields, &
575 this%n_global_probes))
576 call mat_coords%init(this%n_global_probes,3)
577 end if
578 call mpi_gatherv(this%xyz, 3*this%n_local_probes, &
579 mpi_real_precision, global_output_coords, &
580 3*this%n_local_probes_tot, &
581 3*this%n_local_probes_tot_offset, &
582 mpi_real_precision, 0, neko_comm, ierr)
583 if (pe_rank .eq. 0) then
584 call trsp(mat_coords%x, this%n_global_probes, &
585 global_output_coords, 3)
586 !! Write the data to the file
587 call this%fout%write(mat_coords)
588 call mat_coords%free()
589 end if
590 class is (hdf5_file_t)
591
593 if (this%append_out) then
594 call ft%set_overwrite(.false.)
595 else
596 call ft%set_overwrite(.true.)
597 end if
598 call mat_coords%init(3, this%n_local_probes, "coordinates")
599 call copy(mat_coords%x, this%xyz, 3*this%n_local_probes)
600
602 call ft%open("w")
603 call ft%set_active_group("probes")
604
605 ! Check if the NSteps attribute already exists
606 call ft%read_attribute("NSteps", out_int, attr_exist)
607 if (attr_exist) then
608 ! If the attribute exists,
609 ! do not write the coordinates but register the executions
610 this%output_controller%nexecutions = out_int
611 else
612 ! Write out the mesh
613 call ft%write_dataset(mat_coords)
614 out_int = this%n_global_probes
615 call ft%write_attribute("NProbes", out_int)
616 end if
617 call ft%close()
618
620 this%seq_io = .false.
621 call this%vec_out%init(this%n_local_probes, "interpolated_fields_trsp")
622
623 ! Free temporaries
624 call mat_coords%free()
625
626 class default
627 call neko_error("Invalid data. Expected csv_file_t or hdf5_file_t.")
628 end select
629
630 end subroutine probes_init_common
631
633 subroutine probes_free(this)
634 class(probes_t), intent(inout) :: this
635 integer :: i
636
637 if (allocated(this%xyz)) then
638 deallocate(this%xyz)
639 end if
640
641 if (allocated(this%out_vals_trsp)) then
642 deallocate(this%out_vals_trsp)
643 end if
644
645 call this%sampled_fields%free()
646
647 if (allocated(this%n_local_probes_tot)) then
648 deallocate(this%n_local_probes_tot)
649 end if
650
651 if (allocated(this%n_local_probes_tot_offset)) then
652 deallocate(this%n_local_probes_tot_offset)
653 end if
654
655 if (allocated(this%global_output_values)) then
656 deallocate(this%global_output_values)
657 end if
658
659 if (allocated(this%which_fields)) then
660 deallocate(this%which_fields)
661 end if
662
663 if (allocated(this%out_values)) then
664 if (neko_bcknd_device .eq. 1) then
665 do i = 1, this%n_fields
666 call device_unmap(this%out_values(:,i), this%out_values_d(i))
667 end do
668 end if
669 if (allocated(this%out_values_d)) deallocate(this%out_values_d)
670 deallocate(this%out_values)
671 end if
672
673 call this%global_interp%free()
674 call this%mat_out%free()
675 call this%vec_out%free()
676
677 end subroutine probes_free
678
680 subroutine probes_show(this)
681 class(probes_t), intent(in) :: this
682 character(len=LOG_SIZE) :: log_buf ! For logging status
683 integer :: i
684
685 ! Probes summary
686 call neko_log%section('Probes')
687 write(log_buf, '(A,I6)') "Number of probes: ", this%n_global_probes
688 call neko_log%message(log_buf)
689
690 call neko_log%message("xyz-coordinates:", lvl = neko_log_debug)
691 do i = 1, this%n_local_probes
692 write(log_buf, '("(",F10.6,",",F10.6,",",F10.6,")")') this%xyz(:,i)
693 call neko_log%message(log_buf, lvl = neko_log_debug)
694 end do
695
696 ! Field summary
697 write(log_buf, '(A,I6)') "Number of fields: ", this%n_fields
698 call neko_log%message(log_buf)
699 do i = 1, this%n_fields
700 write(log_buf, '(A,I6,A,A)') &
701 "Field: ", i, " ", trim(this%which_fields(i))
702 call neko_log%message(log_buf, lvl = neko_log_debug)
703 end do
704 call neko_log%end_section()
705 call neko_log%newline()
706
707 end subroutine probes_show
708
710 subroutine probes_debug(this)
711 class(probes_t) :: this
712
713 character(len=LOG_SIZE) :: log_buf ! For logging status
714 integer :: i
715
716 do i = 1, this%n_local_probes
717 write (log_buf, *) pe_rank, "/", this%global_interp%pe_owner(i), &
718 "/" , this%global_interp%el_owner0(i)
719 call neko_log%message(log_buf)
720 write(log_buf, '(A5, "(", F10.6, ",", F10.6, ",", F10.6, ")")') &
721 "rst: ", this%global_interp%rst(:, i)
722 call neko_log%message(log_buf)
723 end do
724 end subroutine probes_debug
725
727 subroutine probes_setup_offset(this)
728 class(probes_t) :: this
729 integer :: ierr
730 this%n_local_probes_tot = 0
731 this%n_local_probes_tot_offset = 0
732 this%n_probes_offset = 0
733 call mpi_gather(this%n_local_probes, 1, mpi_integer, &
734 this%n_local_probes_tot, 1, mpi_integer, &
735 0, neko_comm, ierr)
736
737 call mpi_exscan(this%n_local_probes, this%n_probes_offset, 1, &
738 mpi_integer, mpi_sum, neko_comm, ierr)
739 call mpi_gather(this%n_probes_offset, 1, mpi_integer, &
740 this%n_local_probes_tot_offset, 1, mpi_integer, &
741 0, neko_comm, ierr)
742
743
744
745 end subroutine probes_setup_offset
746
751 subroutine probes_evaluate_and_write(this, time)
752 class(probes_t), intent(inout) :: this
753 type(time_state_t), intent(in) :: time
754 integer :: i, ierr
755 logical :: do_interp_on_host = .false.
756 character(len=1000) :: group_name
757 real(kind=rp) :: time_
758 type(vector_t) :: vec_time
759 integer :: out_int
760
762 if (time%t .lt. this%start_time) return
763
765 do i = 1, this%n_fields
766 call this%global_interp%evaluate(this%out_values(:,i), &
767 this%sampled_fields%items(i)%ptr%x, &
768 do_interp_on_host)
769 end do
770
771 if (neko_bcknd_device .eq. 1) then
772 do i = 1, this%n_fields
773 call device_memcpy(this%out_values(:,i), this%out_values_d(i), &
774 this%n_local_probes, device_to_host, sync = .true.)
775 end do
776 end if
777
778 if (this%output_controller%check(time)) then
779 select type (ft => this%fout%file_type)
780 type is (csv_file_t)
781 ! Gather all values to rank 0
782 ! If io is only done at root
783 if (this%seq_io) then
784 call trsp(this%out_vals_trsp, this%n_fields, &
785 this%out_values, this%n_local_probes)
786 call mpi_gatherv(this%out_vals_trsp, &
787 this%n_fields*this%n_local_probes, &
788 mpi_real_precision, this%global_output_values, &
789 this%n_fields*this%n_local_probes_tot, &
790 this%n_fields*this%n_local_probes_tot_offset, &
791 mpi_real_precision, 0, neko_comm, ierr)
792 if (pe_rank .eq. 0) then
793 call trsp(this%mat_out%x, this%n_global_probes, &
794 this%global_output_values, this%n_fields)
795 call this%fout%write(this%mat_out, time%t)
796 end if
797 else
798 call neko_error("CSV outputs only works sequentially")
799 end if
800
801 type is (hdf5_file_t)
802
803 if (this%seq_io) then
804 call neko_error("HDF5 outputs only works in parallel currently.")
805
806 else
807
808 ! Append output format
809 if (this%append_out) then
810
811 call ft%open("w")
812 call ft%set_active_group("probes")
813 ! Write Nsteps in root
814 out_int = this%output_controller%nexecutions + 1
815 call ft%write_attribute("NSteps", out_int)
816 ! Write out the data
817 do i = 1, this%n_fields
818 call copy(this%vec_out%x, this%out_values(:,i), &
819 this%vec_out%size())
820 this%vec_out%name = trim(this%which_fields(i))
821 call ft%write_dataset(this%vec_out)
822 end do
823
824 ! Write the time by hacking the vector write
825 if (pe_rank .eq. 0) then
826 call vec_time%init(1, "time")
827 vec_time%x(1) = time%t
828 else
829 call vec_time%init(0, "time")
830 end if
831 call ft%write_dataset(vec_time)
832 call vec_time%free()
833 call ft%close()
834
835 ! Write data in different steps
836 else
837
838 out_int = this%output_controller%nexecutions + 1
839 ! Set up the name
840 write(group_name, '(A,I0)') "probes/Step_", out_int
841 call ft%open("w")
842 call ft%set_active_group("probes")
843 ! Write Nsteps in root
844 call ft%write_attribute("NSteps", out_int)
845 ! Write out the data
846 call ft%set_active_group(trim(group_name))
847 do i = 1, this%n_fields
848 call copy(this%vec_out%x, this%out_values(:,i), &
849 this%vec_out%size())
850 this%vec_out%name = trim(this%which_fields(i))
851 call ft%write_dataset(this%vec_out)
852 end do
853 ! Write the time as an attribute
854 time_ = time%t
855 call ft%write_attribute("time", time_)
856 call ft%close()
857 end if
858
859 end if
860 class default
861 call neko_error("Invalid data. Expected csv_file_t or hdf5_file_t.")
862 end select
863
864 !! Register the execution of the activity
865 call this%output_controller%register_execution()
866 end if
867
868 end subroutine probes_evaluate_and_write
869
872 subroutine read_probe_locations(this, xyz, n_local_probes, n_global_probes, &
873 points_file)
874 class(probes_t), intent(inout) :: this
875 character(len=:), allocatable :: points_file
876 real(kind=rp), allocatable :: xyz(:,:)
877 integer, intent(inout) :: n_local_probes, n_global_probes
878 type(matrix_t) :: mat_in
879 integer :: ierr
880
882 type(file_t) :: file_in
883
884 call file_in%init(trim(points_file))
886 select type (ft => file_in%file_type)
887 type is (csv_file_t)
888 call read_xyz_from_csv(xyz, n_local_probes, n_global_probes, ft)
889 this%seq_io = .true.
890 type is (hdf5_file_t)
891 call ft%open("r")
892 call ft%read_dataset("xyz", mat_in, "rank_0")
893 call ft%close()
894
895 ! Copy the data to the xyz location
896 n_local_probes = mat_in%get_ncols()
897 call mpi_allreduce(n_local_probes, n_global_probes, 1, mpi_integer, &
898 mpi_sum, neko_comm, ierr)
899
900 allocate(xyz(3, n_local_probes)) ! We asume that axis 1 has 3 entries
901 call copy(xyz, mat_in%x, 3*n_local_probes)
902 call mat_in%free()
903 class default
904 call neko_error("Invalid data. Expected csv_file_t.")
905 end select
906
908 call file_free(file_in)
909
910 end subroutine read_probe_locations
911
917 subroutine read_xyz_from_csv(xyz, n_local_probes, n_global_probes, f)
918 type(csv_file_t), intent(inout) :: f
919 real(kind=rp), allocatable :: xyz(:,:)
920 integer, intent(inout) :: n_local_probes, n_global_probes
921 type(matrix_t) :: mat_in, mat_in2
922 integer :: n_lines
923
924 n_lines = f%count_lines()
925
926 ! Update the number of probes
927 n_global_probes = n_lines
928
929 ! Initialize the temporal array
930 if (pe_rank .eq. 0) then
931 n_local_probes = n_global_probes
932 allocate(xyz(3, n_local_probes))
933 call mat_in%init(n_global_probes,3)
934 call mat_in2%init(3, n_global_probes)
935 call f%read(mat_in)
936 call trsp(xyz, 3, mat_in%x, n_global_probes)
937 else
938 n_local_probes = 0
939 allocate(xyz(3, n_local_probes))
940 end if
941
942 end subroutine read_xyz_from_csv
943end 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:78
Copy data between host and device (or device and device)
Definition device.F90:72
Unmap a Fortran array from a device (deassociate and free)
Definition device.F90:84
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:53
integer, public pe_size
MPI size of communicator.
Definition comm.F90:61
integer, public pe_rank
MPI rank.
Definition comm.F90:58
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:45
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
integer, parameter, public device_to_host
Definition device.F90:48
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:90
type(log_t), public neko_log
Global log stream.
Definition log.f90:80
integer, parameter, public log_size
Definition log.f90:46
Definition math.f90:60
subroutine, public copy(a, b, n)
Copy a vector .
Definition math.f90:289
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:752
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:874
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:681
subroutine probes_free(this)
Destructor.
Definition probes.F90:634
subroutine probes_debug(this)
Show the status of processor/element owner and error code for each point.
Definition probes.F90:711
subroutine probes_setup_offset(this)
Setup offset for rank 0.
Definition probes.F90:728
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:918
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:43
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.