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
Classes
pycompadre.KokkosParser
- class pycompadre.KokkosParser
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__(*args, **kwargs)
Overloaded function.
__init__(self: pycompadre._pycompadre.KokkosParser, args: list[str], print: bool = False) -> None
__init__(self: pycompadre._pycompadre.KokkosParser, print: bool = False) -> None
- static status() str
pycompadre.GMLS
- class pycompadre.GMLS
- __init__(*args, **kwargs)
Overloaded function.
__init__(self: pycompadre._pycompadre.GMLS, 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: pycompadre._pycompadre.GMLS, 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: pycompadre._pycompadre.GMLS, 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(*args, **kwargs)
Overloaded function.
addTargets(self: pycompadre._pycompadre.GMLS, arg0: pycompadre._pycompadre.TargetOperation) -> None
Add a target operation.
addTargets(self: pycompadre._pycompadre.GMLS, arg0: list[pycompadre._pycompadre.TargetOperation]) -> None
Add a list of target operations.
- containsValidAlphas(self: pycompadre._pycompadre.GMLS) bool
- generateAlphas(self: pycompadre._pycompadre.GMLS, number_of_batches: int = 1, keep_coefficients: bool = False, clear_cache: bool = True) None
- getNN(self: int, arg0: int, arg1: pycompadre._pycompadre.ReconstructionSpace) int
Heuristic number of neighbors.
- getNP(self: int, arg0: int, arg1: pycompadre._pycompadre.ReconstructionSpace) int
Get size of basis.
- getSolutionSet(self: pycompadre._pycompadre.GMLS, validity_check: bool = True) pycompadre._pycompadre.SolutionSet
- getWeightingParameter(self: pycompadre._pycompadre.GMLS, parameter index: int = 0) int
Get weighting kernel parameter[index].
- getWeightingType(self: pycompadre._pycompadre.GMLS) pycompadre._pycompadre.WeightingFunctionType
Get the weighting type.
- setSourceExtraData(self: pycompadre._pycompadre.GMLS, arg0: numpy.ndarray[numpy.float64]) None
- setTargetExtraData(self: pycompadre._pycompadre.GMLS, arg0: numpy.ndarray[numpy.float64]) None
- setWeightingParameter(self: pycompadre._pycompadre.GMLS, parameter value: int, parameter index: int = 0) None
Set weighting kernel parameter[index] to parameter value.
- setWeightingPower(self: pycompadre._pycompadre.GMLS, parameter value: int, parameter index: int = 0) None
Set weighting kernel parameter[index] to parameter value. [DEPRECATED]
- setWeightingType(*args, **kwargs)
Overloaded function.
setWeightingType(self: pycompadre._pycompadre.GMLS, arg0: str) -> None
Set the weighting type with a string.
setWeightingType(self: pycompadre._pycompadre.GMLS, arg0: pycompadre._pycompadre.WeightingFunctionType) -> None
Set the weighting type with a WeightingFunctionType.
pycompadre.ParticleHelper
- class pycompadre.ParticleHelper
Class to manage calling PointCloudSearch, moving data to/from Numpy arrays in Kokkos::Views, and applying GMLS solutions to multidimensional data arrays
- __init__(*args, **kwargs)
Overloaded function.
__init__(self: pycompadre._pycompadre.ParticleHelper) -> None
Requires follow-on call to setGMLSObject(gmls_obj) before use
__init__(self: pycompadre._pycompadre.ParticleHelper, gmls_instance: Compadre::GMLS) -> None
- applyStencil(self: pycompadre._pycompadre.ParticleHelper, input_data: numpy.ndarray[numpy.float64], target_operation: pycompadre._pycompadre.TargetOperation = <TargetOperation.ScalarPointEvaluation: 0>, sampling_functional: pycompadre._pycompadre.SamplingFunctional = <pycompadre._pycompadre.SamplingFunctional object at 0x7fb94dcb84b0>, evaluation_site_local_index: int = 0) numpy.ndarray[numpy.float64]
- applyStencilAllTargetsAllAdditionalEvaluationSites(self: pycompadre._pycompadre.ParticleHelper, input_data: numpy.ndarray[numpy.float64], target_operation: pycompadre._pycompadre.TargetOperation = <TargetOperation.ScalarPointEvaluation: 0>, sampling_functional: pycompadre._pycompadre.SamplingFunctional = <pycompadre._pycompadre.SamplingFunctional object at 0x7fb94dcb8230>) numpy.ndarray[numpy.float64]
- applyStencilSingleTarget(self: pycompadre._pycompadre.ParticleHelper, input_data: numpy.ndarray[numpy.float64], target_operation: pycompadre._pycompadre.TargetOperation = <TargetOperation.ScalarPointEvaluation: 0>, sampling_functional: pycompadre._pycompadre.SamplingFunctional = <pycompadre._pycompadre.SamplingFunctional object at 0x7fb94dca3e30>, evaluation_site_local_index: int = 0) float
- generateKDTree(self: pycompadre._pycompadre.ParticleHelper, arg0: numpy.ndarray[numpy.float64]) None
- generateNeighborListsFromKNNSearchAndSet(self: pycompadre._pycompadre.ParticleHelper, target_sites: numpy.ndarray[numpy.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.ParticleHelper) Compadre::NeighborLists<Kokkos::View<int*> >
- getGMLSObject(self: pycompadre._pycompadre.ParticleHelper) Compadre::GMLS
- getNeighborLists(self: pycompadre._pycompadre.ParticleHelper) Compadre::NeighborLists<Kokkos::View<int*> >
- getPolynomialCoefficients(self: pycompadre._pycompadre.ParticleHelper, input_data: numpy.ndarray[numpy.float64]) numpy.ndarray[numpy.float64]
- getReferenceOutwardNormalDirection(self: pycompadre._pycompadre.ParticleHelper) numpy.ndarray[numpy.float64]
- getSourceSites(self: pycompadre._pycompadre.ParticleHelper) numpy.ndarray[numpy.float64]
- getTangentBundle(self: pycompadre._pycompadre.ParticleHelper) numpy.ndarray[numpy.float64]
- getTargetSites(self: pycompadre._pycompadre.ParticleHelper) numpy.ndarray[numpy.float64]
- getWindowSizes(self: pycompadre._pycompadre.ParticleHelper) numpy.ndarray[numpy.float64]
- setAdditionalEvaluationSitesData(self: pycompadre._pycompadre.ParticleHelper, additional_evaluation_sites_neighbor_lists: numpy.ndarray[numpy.int32], additional_evaluation_sites_coordinates: numpy.ndarray[numpy.float64]) None
- setGMLSObject(self: pycompadre._pycompadre.ParticleHelper, gmls_object: Compadre::GMLS) None
- setNeighbors(self: pycompadre._pycompadre.ParticleHelper, neighbor_lists: numpy.ndarray[numpy.int32]) None
Sets neighbor lists from 2D array where first column is number of neighbors for corresponding row’s target site.
- setReferenceOutwardNormalDirection(self: pycompadre._pycompadre.ParticleHelper, reference_normal_directions: numpy.ndarray[numpy.float64], use_to_orient_surface: bool = True) None
- setTangentBundle(self: pycompadre._pycompadre.ParticleHelper, tangent_bundle: numpy.ndarray[numpy.float64]) None
- setTargetSites(self: pycompadre._pycompadre.ParticleHelper, target_coordinates: numpy.ndarray[numpy.float64]) None
- setWindowSizes(self: pycompadre._pycompadre.ParticleHelper, window_sizes: numpy.ndarray[numpy.float64]) None
pycompadre.ReconstructionSpace
- class pycompadre.ReconstructionSpace
Members:
ScalarTaylorPolynomial
VectorTaylorPolynomial
VectorOfScalarClonesTaylorPolynomial
DivergenceFreeVectorTaylorPolynomial
BernsteinPolynomial
- property name
pycompadre.TargetOperation
- class pycompadre.TargetOperation
Members:
ScalarPointEvaluation
VectorPointEvaluation
LaplacianOfScalarPointEvaluation
VectorLaplacianPointEvaluation
GradientOfScalarPointEvaluation
GradientOfVectorPointEvaluation
DivergenceOfVectorPointEvaluation
CurlOfVectorPointEvaluation
CurlCurlOfVectorPointEvaluation
PartialXOfScalarPointEvaluation
PartialYOfScalarPointEvaluation
PartialZOfScalarPointEvaluation
ChainedStaggeredLaplacianOfScalarPointEvaluation
GaussianCurvaturePointEvaluation
CellAverageEvaluation
CellIntegralEvaluation
FaceNormalAverageEvaluation
FaceNormalIntegralEvaluation
EdgeTangentAverageEvaluation
EdgeTangentIntegralEvaluation
- property name
pycompadre.WeightingFunctionType
- class pycompadre.WeightingFunctionType
Members:
Power
Gaussian
CubicSpline
CardinalCubicBSpline
Cosine
Sigmoid
- property name
pycompadre.NeighborLists
- class pycompadre.NeighborLists
- computeMaxNumNeighbors(self: pycompadre._pycompadre.NeighborLists) None
Compute maximum number of neighbors over all neighborhoods.
- computeMinNumNeighbors(self: pycompadre._pycompadre.NeighborLists) None
Compute minimum number of neighbors over all neighborhoods.
- getMaxNumNeighbors(self: pycompadre._pycompadre.NeighborLists) int
Get maximum number of neighbors over all neighborhoods.
- getMinNumNeighbors(self: pycompadre._pycompadre.NeighborLists) int
Get minimum number of neighbors over all neighborhoods.
- getNeighbor(self: pycompadre._pycompadre.NeighborLists, target index: int, local neighbor number: int) int
Get neighbor index from target index and local neighbor number.
- getNumberOfNeighbors(self: pycompadre._pycompadre.NeighborLists, target index: int) int
Get number of neighbors for target.
- getNumberOfTargets(self: pycompadre._pycompadre.NeighborLists) int
Number of targets.
- getTotalNeighborsOverAllLists(self: pycompadre._pycompadre.NeighborLists) int
Get total storage size of all neighbor lists combined.
pycompadre.Quadrature
- class pycompadre.Quadrature
- __init__(self: pycompadre._pycompadre.Quadrature, order_of_quadrature_points: int, dimension_of_quadrature_points: int, quadrature_type: str = 'LINE') None
- getDimensionOfQuadraturePoints(self: pycompadre._pycompadre.Quadrature) int
- getNumberOfQuadraturePoints(self: pycompadre._pycompadre.Quadrature) int
- getOrderOfQuadraturePoints(self: pycompadre._pycompadre.Quadrature) int
- getQuadratureType(self: pycompadre._pycompadre.Quadrature) pycompadre._pycompadre.QuadratureType
- getSites(self: pycompadre._pycompadre.Quadrature) numpy.ndarray[numpy.float64]
- getWeights(self: pycompadre._pycompadre.Quadrature) numpy.ndarray[numpy.float64]
- validQuadrature(self: pycompadre._pycompadre.Quadrature) bool
Has quadrature been generated.
pycompadre.Kokkos
Namespace with a subset of Kokkos functionality. KokkosParser class should be preferred.
- pycompadre.Kokkos.finalize() None
- pycompadre.Kokkos.is_initialized() bool
- pycompadre.Kokkos.status() str