Skip to content

Commit

Permalink
Merge pull request #151 from precice/python-bindings-v2.5.0.1
Browse files Browse the repository at this point in the history
Release v2.5.0.1
  • Loading branch information
IshaanDesai authored Sep 5, 2022
2 parents cff454d + f48b04d commit fcee5ff
Show file tree
Hide file tree
Showing 6 changed files with 447 additions and 1 deletion.
1 change: 1 addition & 0 deletions .github/workflows/pythonpublish.yml
Original file line number Diff line number Diff line change
Expand Up @@ -24,6 +24,7 @@ jobs:
cython
packaging
numpy
pkgconfig
- name: Build and publish
env:
TWINE_USERNAME: ${{ secrets.PYPI_USERNAME }}
Expand Down
4 changes: 4 additions & 0 deletions CHANGELOG.md
Original file line number Diff line number Diff line change
Expand Up @@ -4,6 +4,10 @@ All notable changes to this project will be documented in this file.

## latest

## 2.5.0.1

* Add pkgconfig as dependency to the pythonpublish workflow https://github.com/precice/python-bindings/commit/200dc2aba160e18a7d1dae44ef3493d546e69eb9

## 2.5.0.0

* Bindings now use pkgconfig to determine flags and link to preCICE. https://github.com/precice/python-bindings/pull/149
Expand Down
12 changes: 12 additions & 0 deletions cyprecice/SolverInterface.pxd
Original file line number Diff line number Diff line change
Expand Up @@ -101,6 +101,18 @@ cdef extern from "precice/SolverInterface.hpp" namespace "precice":

void readScalarData (const int dataID, const int valueIndex, double& value) const

# Gradient related API

bool isGradientDataRequired(int dataID) const;

void writeBlockVectorGradientData(int dataID, int size, const int* valueIndices, const double* gradientValues);

void writeScalarGradientData(int dataID, int valueIndex, const double* gradientValues);

void writeVectorGradientData(int dataID, int valueIndex, const double* gradientValues);

void writeBlockScalarGradientData(int dataID, int size, const int* valueIndices, const double* gradientValues);

# direct mesh access

void setMeshAccessRegion (const int meshID, const double* boundingBox) const
Expand Down
239 changes: 239 additions & 0 deletions cyprecice/cyprecice.pyx
Original file line number Diff line number Diff line change
Expand Up @@ -1211,6 +1211,245 @@ cdef class Interface:
self.thisptr.readScalarData (data_id, vertex_id, _value)
return _value

def write_block_vector_gradient_data (self, data_id, vertex_ids, gradientValues):
"""
Writes vector gradient data given as block. This function writes gradient values of specified vertices to a dataID.
Values are provided as a block of continuous memory. Values are stored in a numpy array [N x D] where N = number
of vertices and D = number of gradient components.
Parameters
----------
data_id : int
Data ID to write to.
vertex_ids : array_like
Indices of the vertices.
gradientValues : array_like
Gradient values differentiated in the spacial direction (dx, dy) for 2D space, (dx, dy, dz) for 3D space
Notes
-----
Previous calls:
Count of available elements at values matches the configured dimension
Count of available elements at vertex_ids matches the given size
Initialize() has been called
Data with dataID has attribute hasGradient = true
Examples
--------
Write block gradient vector data for a 2D problem with 2 vertices:
>>> data_id = 1
>>> vertex_ids = [1, 2]
>>> gradientValues = np.array([[v1x_dx, v1y_dx, v1x_dy, v1y_dy], [v2x_dx, v2y_dx, v2x_dy, v2y_dy]])
>>> interface.write_block_vector_gradient_data(data_id, vertex_ids, gradientValues)
Write block vector data for a 3D (D=3) problem with 2 (N=2) vertices:
>>> data_id = 1
>>> vertex_ids = [1, 2]
>>> gradientValues = np.array([[v1x_dx, v1y_dx, v1z_dx, v1x_dy, v1y_dy, v1z_dy, v1x_dz, v1y_dz, v1z_dz], [v2x_dx, v2y_dx, v2z_dx, v2x_dy, v2y_dy, v2z_dy, v2x_dz, v2y_dz, v2z_dz]])
>>> interface.write_block_vector_gradient_data(data_id, vertex_ids, gradientValues)
"""
check_array_like(vertex_ids, "vertex_ids", "write_block_vector_gradient_data")
check_array_like(gradientValues, "gradientValues", "write_block_vector_gradient_data")

if not isinstance(gradientValues, np.ndarray):
gradientValues = np.asarray(gradientValues)

if len(gradientValues) > 0:
size, dimensions = gradientValues.shape
assert dimensions == self.get_dimensions() * self.get_dimensions(), "Dimensions of vector data in write_block_vector_gradient_data does not match with dimensions in problem definition. Provided dimensions: {}, expected dimensions: {}".format(dimensions, self.get_dimensions() * self.get_dimensions())
if len(gradientValues) == 0:
size = 0

cdef np.ndarray[int, ndim=1] _vertex_ids = np.ascontiguousarray(vertex_ids, dtype=np.int32)
cdef np.ndarray[double, ndim=1] _gradientValues = np.ascontiguousarray(gradientValues.flatten(), dtype=np.double)

assert _gradientValues.size == size * self.get_dimensions() * self.get_dimensions(), "Dimension of vector gradient data provided in write_block_vector_gradient_data does not match problem definition. Check length of input data provided. Provided size: {}, expected size: {}".format(_gradientValues.size, size * self.get_dimensions() * self.get_dimensions())
assert _vertex_ids.size == size, "Vertex IDs are of incorrect length in write_block_vector_gradient_data. Check length of vertex ids input. Provided size: {}, expected size: {}".format(_vertex_ids.size, size)

self.thisptr.writeBlockVectorGradientData (data_id, size, <const int*>_vertex_ids.data, <const double*>_gradientValues.data)

def write_scalar_gradient_data (self, data_id, vertex_id, gradientValues):
"""
Writes scalar gradient data to a vertex
This function writes the corresponding gradient matrix value of a specified vertex to a dataID.
The gradients need to be provided in the following format:
The 2D-format of gradientValues is (v_dx, v_dy) vector corresponding to the data block v = (v)
differentiated respectively in x-direction dx and y-direction dy
The 3D-format of gradientValues is (v_dx, v_dy, v_dz) vector
corresponding to the data block v = (v) differentiated respectively in spatial directions x-direction dx and y-direction dy and z-direction dz
Parameters
----------
data_id : int
ID to write to.
vertex_id : int
Index of the vertex.
gradientValue : array_like
A vector of the gradient values.
Notes
-----
Count of available elements at value matches the configured dimension
Vertex with dataID exists and contains data
Data with dataID has attribute hasGradient = true
Previous calls:
initialize() has been called
Examples
--------
Write scalar data for a 2D problem:
>>> data_id = 1
>>> vertex_id = 5
>>> gradientValue = [v5_dx, v5_dy]
>>> interface.write_scalar_gradient_data(data_id, vertex_id, gradientValue)
"""

check_array_like(gradientValues, "gradientValues", "write_scalar_gradient_data")

if not isinstance(gradientValues, np.ndarray):
gradientValues = np.asarray(gradientValues)

cdef np.ndarray[double, ndim=1] _gradientValues = np.ascontiguousarray(gradientValues.flatten(), dtype=np.double)

assert _gradientValues.size == self.get_dimensions(), "Vector data provided for vertex {} in write_scalar_gradient_data does not match problem definition. Check length of input data provided. Provided size: {}, expected size: {}".format(_gradientValues.size, self.get_dimensions())

self.thisptr.writeScalarGradientData(data_id, vertex_id, <const double*>_gradientValues.data)

def write_vector_gradient_data (self, data_id, vertex_id, gradientValues):
"""
Writes vector gradient data to a vertex
This function writes the corresponding gradient matrix value of a specified vertex to a dataID.
The gradients need to be provided in the following format:
The 2D-format of \p gradientValues is (vx_dx, vy_dx, vx_dy, vy_dy) vector corresponding to the data block v = (vx, vy)
differentiated respectively in x-direction dx and y-direction dy
The 3D-format of \p gradientValues is (vx_dx, vy_dx, vz_dx, vx_dy, vy_dy, vz_dy, vx_dz, vy_dz, vz_dz) vector
corresponding to the data block v = (vx, vy, vz) differentiated respectively in spatial directions x-direction dx and y-direction dy and z-direction dz
Parameters
----------
data_id : int
ID to write to.
vertex_id : int
Index of the vertex.
gradientValue : array_like
A vector of the gradient values.
Notes
-----
Count of available elements at value matches the configured dimension
Vertex with dataID exists and contains data
Data with dataID has attribute hasGradient = true
Previous calls:
initialize() has been called
Examples
--------
Write scalar data for a 2D problem:
>>> data_id = 1
>>> vertex_id = 5
>>> gradientValue = [v5x_dx, v5y_dx, v5x_dy,v5y_dy]
>>> interface.write_vector_gradient_data(data_id, vertex_id, gradientValue)
"""

check_array_like(gradientValues, "gradientValues", "write_vector_gradient_data")

if not isinstance(gradientValues, np.ndarray):
gradientValues = np.asarray(gradientValues)

cdef np.ndarray[double, ndim=1] _gradientValues = np.ascontiguousarray(gradientValues.flatten(), dtype=np.double)

assert _gradientValues.size == self.get_dimensions() * self.get_dimensions(), "Dimensions of vector gradient data provided for vertex {} in write_vector_gradient_data does not match problem definition. Check length of input data provided. Provided size: {}, expected size: {}".format(_gradientValues.size, self.get_dimensions() * self.get_dimensions())

self.thisptr.writeVectorGradientData(data_id, vertex_id, <const double*>_gradientValues.data)

def write_block_scalar_gradient_data (self, data_id, vertex_ids, gradientValues):
"""
Writes scalar gradient data given as block. This function writes values of specified vertices to a dataID.
Values are provided as a block of continuous memory. Values are stored in a numpy array [N x D] where N = number
of vertices and D = dimensions of geometry.
Parameters
----------
data_id : int
Data ID to write to.
vertex_ids : array_like
Indices of the vertices.
gradientValues : array_like
Gradient values differentiated in the spacial direction (dx, dy) for 2D space, (dx, dy, dz) for 3D space
Notes
-----
Previous calls:
Count of available elements at values matches the configured dimension
Count of available elements at vertex_ids matches the given size
Initialize() has been called
Data with dataID has attribute hasGradient = true
Examples
--------
Write block gradient scalar data for a 2D problem with 2 vertices:
>>> data_id = 1
>>> vertex_ids = [1, 2]
>>> gradientValues = np.array([[v1_dx, v1_dy], [v2_dx, v2_dy]])
>>> interface.write_block_scalar_gradient_data(data_id, vertex_ids, gradientValues)
Write block scalar data for a 3D (D=3) problem with 2 (N=2) vertices:
>>> data_id = 1
>>> vertex_ids = [1, 2]
>>> values = np.array([[v1_dx, v1_dy, v1x_dz], [v2_dx, v2_dy, v2_dz]])
>>> interface.write_block_scalar_gradient_data(data_id, vertex_ids, values)
"""
check_array_like(vertex_ids, "vertex_ids", "write_block_scalar_gradient_data")
check_array_like(gradientValues, "gradientValues", "write_block_sclar_gradient_data")

if not isinstance(gradientValues, np.ndarray):
gradientValues = np.asarray(gradientValues)

if len(gradientValues) > 0:
size, dimensions = gradientValues.shape
assert dimensions == self.get_dimensions() , "Dimensions of scalar gradient data provided in write_block_scalar_gradient_data does not match with dimensions in problem definition. Provided dimensions: {}, expected dimensions: {}".format(dimensions, self.get_dimensions())
if len(gradientValues) == 0:
size = 0

cdef np.ndarray[int, ndim=1] _vertex_ids = np.ascontiguousarray(vertex_ids, dtype=np.int32)
cdef np.ndarray[double, ndim=1] _gradientValues = np.ascontiguousarray(gradientValues.flatten(), dtype=np.double)

assert _gradientValues.size == size * self.get_dimensions(), "Scalar gradient data is not provided for all vertices in write_block_scalar_gradient_data. Check length of input data provided. Provided size: {}, expected size: {}".format(_gradientValues.size, size * self.get_dimensions())
assert _vertex_ids.size == size, "Vertex IDs are of incorrect length in write_block_scalar_gradient_data. Check length of vertex ids input. Provided size: {}, expected size: {}".format(_vertex_ids.size, size)

self.thisptr.writeBlockScalarGradientData (data_id, size, <const int*>_vertex_ids.data, <const double*>_gradientValues.data)

def is_gradient_data_required(self,data_id):
"""
Checks if the given data set requires gradient data. We check if the data object has been intialized with the gradient flag.
Parameters
----------
data_id : int
Data ID to check.
Returns
-------
bool
True if gradient data is required for a dataID.
Examples
--------
Check if gradient data is required for a dataID:
>>> data_id = 1
>>> interface.is_gradient_data_required(data_id)
"""
return self.thisptr.isGradientDataRequired(data_id)


def set_mesh_access_region (self, mesh_id, bounding_box):
"""
This function is required if you don't want to use the mapping schemes in preCICE, but rather
Expand Down
52 changes: 51 additions & 1 deletion test/SolverInterface.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -395,6 +395,56 @@ void SolverInterface:: getMeshVerticesAndIDs
}
}

bool SolverInterface::isGradientDataRequired(int dataID) const
{
return 0;
}

void SolverInterface::writeBlockVectorGradientData(
int dataID,
int size,
const int *valueIndices,
const double *gradientValues)
{
fake_read_write_buffer.clear();
for (int i = 0; i < size * this->getDimensions() * this->getDimensions(); i++) {
fake_read_write_buffer.push_back(gradientValues[i]);
}
}

void SolverInterface::writeScalarGradientData(
int dataID,
int valueIndex,
const double *gradientValues)
{
fake_read_write_buffer.clear();
for (int i = 0; i < this->getDimensions(); i++) {
fake_read_write_buffer.push_back(gradientValues[i]);
}
}
void SolverInterface::writeBlockScalarGradientData(
int dataID,
int size,
const int *valueIndices,
const double *gradientValues)
{
fake_read_write_buffer.clear();
for (int i = 0; i < size * this->getDimensions(); i++) {
fake_read_write_buffer.push_back(gradientValues[i]);
}
}

void SolverInterface::writeVectorGradientData(
int dataID,
int valueIndex,
const double *gradientValues)
{
fake_read_write_buffer.clear();
for (int i = 0; i < this->getDimensions() * this->getDimensions(); i++) {
fake_read_write_buffer.push_back(gradientValues[i]);
}
}

std::string getVersionInformation()
{
std::string dummy ("dummy");
Expand Down Expand Up @@ -423,4 +473,4 @@ const std::string& actionReadIterationCheckpoint()

} // namespace precice, constants

} // namespace precice
} // namespace precice
Loading

0 comments on commit fcee5ff

Please sign in to comment.