Probes
Descriptions of the contents of the Probes class.
- class pysemtools.interpolation.probes.Probes(comm, output_fname: str = './interpolated_fields.csv', probes: ndarray | str | None = None, msh=typing.Union[pysemtools.datatypes.msh.Mesh, str, list], write_coords: bool = True, progress_bar: bool = False, point_interpolator_type: str = 'single_point_legendre', max_pts: int = 128, find_points_iterative: list = [False, 5000], find_points_comm_pattern: str = 'point_to_point', elem_percent_expansion: float = 0.01, global_tree_type: str = 'rank_bbox', global_tree_nbins: int = 1024, use_autograd: bool = False, find_points_tol: float = np.float64(2.220446049250313e-15), find_points_max_iter: int = 50)
Class to interpolate fields at probes from a SEM mesh.
Main interpolation class. This works in parallel.
If the points to interpolate are available only at rank 0, make sure to only pass them at that rank and set the others to None. In that case, the points will be scattered to all ranks.
If the points are passed in all ranks, then they will be considered as different points to be interpolated and the work to interpolate will be multiplied. If you are doing this, make sure that the points that each rank pass are different. Otherwise simply pass None in all ranks but 0.
See example below to observe how to avoid unnecessary replication of data in all ranks when passing probes as argument if the points are the same in all ranks.
If reading probes from file, they will be read on rank 0 and scattered unless parallel hdf5 is used, in which case the probes will be read in all ranks. (In development)
- Parameters:
- commMPI communicator
MPI communicator.
- output_fnamestr
Output file name. Default is “./interpolated_fields.csv”. Note that you can change the file extension to .hdf5 to write in hdf5 format. by using the name “interpolated_fields.hdf5”.
- probesUnion[np.ndarray, str]
Probes coordinates. If a string, it is assumed to be a file name.
- mshUnion[Mesh, str, list]
Mesh data. If a string, it is assumed to be a file name. If it is a list the first entry is the file name and the second is the dtype of the data. if it is a Mesh object, the x, y, z coordinates are taken from the object.
- write_coordsbool
If True, the coordinates of the probes are written to a file. Default is True.
- progress_barbool
If True, a progress bar is shown. Default is False.
- point_interpolator_typestr
Type of point interpolator. Default is single_point_legendre. options are: single_point_legendre, single_point_lagrange, multiple_point_legendre_numpy, multiple_point_legendre_torch.
- max_ptsint, optional
Maximum number of points to interpolate. Default is 128. Used if multiple point interpolator is selected.
- find_points_iterativelist
List with two elements. First element is a boolean that indicates if the search is iterative. Second element is the maximum number of candidate ranks to send the data. This affects memory. Default is [False, 5000].
- find_points_comm_patternstr
Communication pattern for finding points. Default is point_to_point. options are: point_to_point, collective.
- elem_percent_expansionfloat
Percentage expansion of the element bounding box. Default is 0.01.
- global_tree_typestr
How is the global tree constructed to determine rank candidates for the probes. Only really used if using tree structures to determine candidates. Default is rank_bbox. options are: rank_bbox, domain_binning.
- global_tree_nbinsint
Number of bins in the global tree. Only used if the global tree is domain_binning. Default is 1024.
- use_autogradbool
If True, autograd is used. Default is False.
- find_points_tolfloat
The tolerance to use when finding points. Default is np.finfo(np.double).eps * 10.
- find_points_max_iterint
The maximum number of iterations to use when finding points. Default is 50.
- Attributes:
- probesndarray
2D array of probe coordinates. shape = (n_probes, 3).
- interpolated_fieldsndarray
2D array of interpolated fields at probes. shape = (n_probes, n_fields + 1). The first column is always time, the rest are the interpolated fields.
Methods
interpolate_from_field_list
(t, field_list, comm)Interpolate the probes from a list of fields.
Notes
A sample input file can be found in the examples folder of the main repository, However, the file is not used in said example.
Examples
Initialize from file:
>>> from mpi4py import MPI >>> from pysemtools.interpolation.probes import Probes >>> comm = MPI.COMM_WORLD >>> probes = Probes(comm, filename="path/to/params.json")
2. Initialize from code, passing everything as arguments. Assume msh is created. One must then create the probe data in rank 0. A dummy probe_data must be created in all other ranks
>>> from mpi4py import MPI >>> from pysemtools.interpolation.probes import Probes >>> comm = MPI.COMM_WORLD >>> if comm.Get_rank() == 0: >>> probes_data = np.array([[0.0, 0.0, 0.0], [1.0, 1.0, 1.0], [2.0, 2.0, 2.0]]) >>> else: >>> probes_data = None >>> probes = Probes(comm, probes=probes_data, msh=msh)
Note that probes is initialized in all ranks, but the probe_data containing the coordinates are only relevant in rank 0. They are scattered internally.
- interpolate_from_field_list(t, field_list, comm, write_data=True, field_names: list[str] | None = None)
Interpolate the probes from a list of fields.
This method interpolates from a list of fields (ndarrays of shape (nelv, lz, ly, lx)).
- Parameters:
- tfloat
Time of the field data.
- field_listlist
List of fields to interpolate. Each field is an ndarray of shape (nelv, lz, ly, lx).
- commComm
MPI communicator.
- write_databool
If True, the interpolated data is written to a file. Default is True.
- field_nameslist
List of names of the interpolated fields. Useful when writing to file. Default is None.
Examples
This method can be used to interpolate fields from a list of fields. If you have previosly obtained a set of fields u,v,w as ndarrays of shape (nelv, lz, ly, lx), you can:
>>> probes.interpolate_from_field_list(t, [u,v,w], comm)
The results are stored in probes.interpolated_fields attribute. Remember: the first column of this attribute is always the time t given.