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