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