PyCOMPADRE: Compadre Toolkit

Important: Be sure to initialize Kokkos before setting up any objects from this library, by created a scoped KokkosParser object, i.e:

kp = pycompadre.KokkosParser()

When kp goes out of scope, then Kokkos will be finalized.

GMLS type objects are passed to ParticleHelper objects. Be sure to deallocate in reverse order, i.e.:

>> kp = pycompadre.KokkosParser()
>> gmls_obj = GMLS(...)
>> gmls_helper = ParticleHelper(gmls_object)
...
...
>> del gmls_obj
>> del gmls_helper
>> del kp

Project details at: https://github.com/sandialabs/compadre

Implementation details at: https://github.com/sandialabs/compadre/blob/master/pycompadre/pycompadre.cpp

Compadre_DEBUG: False Compadre_BUILD_TYPE: RelWithDebInfo

Classes

pycompadre.KokkosParser

class pycompadre.KokkosParser(*args, **kwargs)

Class wrapping the functionality of Kokkos::ScopeGuard. Multiple instances can can be instantiated and go out of scope with only one instance initiating Kokkos and the last instance finalizing Kokkos when it goes out of scope.

__init__(self, args: collections.abc.Sequence[str], print: bool = False) None
__init__(self, print: bool = False) None
status = <nanobind.nb_func object>

pycompadre.GMLS

class pycompadre.GMLS(*args, **kwargs)
__init__(self, poly_order: int, dimension: int = 3, dense_solver_type: str = 'QR', problem_type: str = 'STANDARD', constraint_type: str = 'NO_CONSTRAINT', curvature_poly_order: int = 2) None
__init__(self, reconstruction_space: pycompadre._pycompadre.ReconstructionSpace, sampling_functional: pycompadre._pycompadre.SamplingFunctional, poly_order: int, dimension: int = 3, dense_solver_type: str = 'QR', problem_type: str = 'STANDARD', constraint_type: str = 'NO_CONSTRAINT', curvature_poly_order: int = 2) None
__init__(self, reconstruction_space: pycompadre._pycompadre.ReconstructionSpace, polynomial_sampling_functional: pycompadre._pycompadre.SamplingFunctional, data_sampling_functional: pycompadre._pycompadre.SamplingFunctional, poly_order: int, dimension: int = 3, dense_solver_type: str = 'QR', problem_type: str = 'STANDARD', constraint_type: str = 'NO_CONSTRAINT', curvature_poly_order: int = 2) None
addTargets(self, arg: pycompadre._pycompadre.TargetOperation, /) None
addTargets(self, arg: collections.abc.Sequence[pycompadre._pycompadre.TargetOperation], /) None

Overloaded function.

  1. addTargets(self, arg: pycompadre._pycompadre.TargetOperation, /) -> None

Add a target operation.

  1. addTargets(self, arg: collections.abc.Sequence[pycompadre._pycompadre.TargetOperation], /) -> None

Add a list of target operations.

containsValidAlphas(self) bool
generateAlphas(self, number_of_batches: int = 1, keep_coefficients: bool = False, clear_cache: bool = True) None
getNN(self, arg0: int, arg1: pycompadre._pycompadre.ReconstructionSpace, /) int

Heuristic number of neighbors.

getNP(self, arg0: int, arg1: pycompadre._pycompadre.ReconstructionSpace, /) int

Get size of basis.

getSolutionSet(self, validity_check: bool = True) pycompadre._pycompadre.SolutionSet
getWeightingParameter(self, parameter index: int = 0) int

Get weighting kernel parameter[index].

getWeightingType(self) pycompadre._pycompadre.WeightingFunctionType

Get the weighting type.

setDimensionOfQuadraturePoints(self, dim: int) None
setOrderOfQuadraturePoints(self, order: int) None
setQuadratureType(self, quadrature_type: str) None
setSourceExtraData(self, arg: numpy.ndarray[dtype=float64], /) None
setTargetExtraData(self, arg: numpy.ndarray[dtype=float64], /) None
setWeightingParameter(self, parameter value: int, parameter index: int = 0) None

Set weighting kernel parameter[index] to parameter value.

setWeightingPower(self, parameter value: int, parameter index: int = 0) None

Set weighting kernel parameter[index] to parameter value. [DEPRECATED]

setWeightingType(self, arg: str, /) None
setWeightingType(self, arg: pycompadre._pycompadre.WeightingFunctionType, /) None

Overloaded function.

  1. setWeightingType(self, arg: str, /) -> None

Set the weighting type with a string.

  1. setWeightingType(self, arg: pycompadre._pycompadre.WeightingFunctionType, /) -> None

Set the weighting type with a WeightingFunctionType.

pycompadre.ParticleHelper

class pycompadre.ParticleHelper(*args, **kwargs)

Class to manage calling PointCloudSearch, moving data to/from Numpy arrays in Kokkos::Views, and applying GMLS solutions to multidimensional data arrays

__init__(self) None
__init__(self, gmls_instance: pycompadre._pycompadre.GMLS) None

Overloaded function.

  1. __init__(self) -> None

    Requires follow-on call to setGMLSObject(gmls_obj) before use

  2. __init__(self, gmls_instance: pycompadre._pycompadre.GMLS) -> None

applyStencil(self, input_data: numpy.ndarray[dtype=float64], target_operation: pycompadre._pycompadre.TargetOperation = TargetOperation.ScalarPointEvaluation, sampling_functional: pycompadre._pycompadre.SamplingFunctional = <pycompadre._pycompadre.SamplingFunctional object at 0x7f133926efa0>, evaluation_site_local_index: int = 0) numpy.ndarray[dtype=float64]
applyStencilAllTargetsAllAdditionalEvaluationSites(self, input_data: numpy.ndarray[dtype=float64], target_operation: pycompadre._pycompadre.TargetOperation = TargetOperation.ScalarPointEvaluation, sampling_functional: pycompadre._pycompadre.SamplingFunctional = <pycompadre._pycompadre.SamplingFunctional object at 0x7f1338011dd0>) numpy.ndarray[dtype=float64]
applyStencilSingleTarget(self, input_data: numpy.ndarray[dtype=float64], target_operation: pycompadre._pycompadre.TargetOperation = TargetOperation.ScalarPointEvaluation, sampling_functional: pycompadre._pycompadre.SamplingFunctional = <pycompadre._pycompadre.SamplingFunctional object at 0x7f1338011da0>, evaluation_site_local_index: int = 0) float
generateKDTree(self, arg: numpy.ndarray[dtype=float64], /) None
generateNeighborListsFromKNNSearchAndSet(self, target_sites: numpy.ndarray[dtype=float64], poly_order: int, dimension: int = 3, epsilon_multiplier: float = 1.6, max_search_radius: float = 0.0, scale_k_neighbor_radius: bool = True, scale_num_neighbors: bool = False) None
getAdditionalEvaluationIndices(self) pycompadre._pycompadre.NeighborLists
getGMLSObject(self) pycompadre._pycompadre.GMLS
getNeighborLists(self) pycompadre._pycompadre.NeighborLists
getPolynomialCoefficients(self, input_data: numpy.ndarray[dtype=float64]) numpy.ndarray[dtype=float64]
getReferenceOutwardNormalDirection(self) numpy.ndarray[dtype=float64]
getSourceSites(self) numpy.ndarray[dtype=float64]
getTangentBundle(self) numpy.ndarray[dtype=float64]
getTargetSites(self) numpy.ndarray[dtype=float64]
getWindowSizes(self) numpy.ndarray[dtype=float64]
setAdditionalEvaluationSitesData(self, additional_evaluation_sites_neighbor_lists: numpy.ndarray[dtype=int32], additional_evaluation_sites_coordinates: numpy.ndarray[dtype=float64]) None
setGMLSObject(self, gmls_object: pycompadre._pycompadre.GMLS) None
setNeighbors(self, neighbor_lists: numpy.ndarray[dtype=int32]) None

Sets neighbor lists from 2D array where first column is number of neighbors for corresponding row’s target site.

setReferenceOutwardNormalDirection(self, reference_normal_directions: numpy.ndarray[dtype=float64], use_to_orient_surface: bool = True) None
setTangentBundle(self, tangent_bundle: numpy.ndarray[dtype=float64]) None
setTargetSites(self, target_coordinates: numpy.ndarray[dtype=float64]) None
setWindowSizes(self, window_sizes: numpy.ndarray[dtype=float64]) None

pycompadre.ReconstructionSpace

class pycompadre.ReconstructionSpace(*values)

pycompadre.TargetOperation

class pycompadre.TargetOperation(*values)

pycompadre.WeightingFunctionType

class pycompadre.WeightingFunctionType(*values)

pycompadre.NeighborLists

class pycompadre.NeighborLists
computeMaxNumNeighbors(self) None

Compute maximum number of neighbors over all neighborhoods.

computeMinNumNeighbors(self) None

Compute minimum number of neighbors over all neighborhoods.

getMaxNumNeighbors(self) int

Get maximum number of neighbors over all neighborhoods.

getMinNumNeighbors(self) int

Get minimum number of neighbors over all neighborhoods.

getNeighbor(self, target index: int, local neighbor number: int) int

Get neighbor index from target index and local neighbor number.

getNumberOfNeighbors(self, target index: int) int

Get number of neighbors for target.

getNumberOfTargets(self) int

Number of targets.

getTotalNeighborsOverAllLists(self) int

Get total storage size of all neighbor lists combined.

pycompadre.Quadrature

class pycompadre.Quadrature(*args, **kwargs)
__init__(self, order_of_quadrature_points: int, dimension_of_quadrature_points: int, quadrature_type: str = 'LINE') None
getDimensionOfQuadraturePoints(self) int
getNumberOfQuadraturePoints(self) int
getOrderOfQuadraturePoints(self) int
getQuadratureType(self) pycompadre._pycompadre.QuadratureType
getSites(self) numpy.ndarray[dtype=float64]
getWeights(self) numpy.ndarray[dtype=float64]
validQuadrature(self) bool

Has quadrature been generated.

pycompadre.Kokkos

Namespace with a subset of Kokkos functionality. KokkosParser class should be preferred.