YetAnotherCoupler 3.2.0_a
Loading...
Searching...
No Matches
Data Structures | Typedefs | Functions
dist_grid_internal.h File Reference
#include "dist_grid.h"
#include "geometry.h"
Include dependency graph for dist_grid_internal.h:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Data Structures

struct  yac_const_basic_grid_data
 
struct  remote_point_info
 
struct  remote_point_infos
 
struct  remote_point
 
struct  remote_points
 

Typedefs

typedef size_t const (*const const_size_t_2_pointer)[2]
 
typedef int const *const const_int_pointer
 
typedef size_t const *const const_size_t_pointer
 
typedef yac_int const *const const_yac_int_pointer
 
typedef enum yac_edge_type const *const const_yac_edge_type_pointer
 
typedef struct bounding_circle const *const const_bounding_circle_pointer
 
typedef struct remote_point_infos const *const const_remote_point_infos_pointer
 

Functions

MPI_Comm yac_dist_grid_pair_get_MPI_Comm (struct yac_dist_grid_pair *grid_pair)
 
struct yac_dist_gridyac_dist_grid_pair_get_dist_grid (struct yac_dist_grid_pair *grid_pair, char const *grid_name)
 
void yac_dist_grid_pair_do_point_search (struct yac_dist_grid_pair *grid_pair, char const *grid_name, yac_coordinate_pointer search_coords, size_t count, size_t *cells)
 
void yac_dist_grid_pair_do_point_search_gc (struct yac_dist_grid_pair *grid_pair, char const *grid_name, yac_coordinate_pointer search_coords, size_t count, size_t *cells)
 
void yac_dist_grid_pair_do_nnn_search (struct yac_dist_grid_pair *grid_pair, char const *grid_name, yac_coordinate_pointer search_coords, size_t count, size_t *local_ids, size_t n, struct yac_interp_field field)
 
void yac_dist_grid_pair_do_bnd_circle_search (struct yac_dist_grid_pair *grid_pair, char const *grid_name, const_bounding_circle_pointer bnd_circles, size_t count, size_t **cells, size_t *num_results_per_bnd_circle, struct yac_interp_field field)
 
void yac_dist_grid_pair_do_cell_search (struct yac_dist_grid_pair *grid_pair, char const *search_grid_name, char const *result_grid_name, size_t *search_cells, size_t count, size_t **result_cells, size_t *num_results_per_search_cell, struct yac_interp_field result_field)
 
void yac_dist_grid_pair_get_cell_neighbours (struct yac_dist_grid_pair *grid_pair, char const *grid_name, size_t *cells, size_t count, size_t *neighbours)
 
struct yac_const_basic_grid_datayac_dist_grid_get_basic_grid_data (struct yac_dist_grid *dist_grid)
 
size_t yac_dist_grid_get_local_count (struct yac_dist_grid *dist_grid, enum yac_location location)
 
void yac_dist_grid_get_local_unmasked_points (struct yac_dist_grid *grid, struct yac_interp_field field, size_t **indices, size_t *count)
 
MPI_Datatype yac_get_remote_point_info_mpi_datatype (MPI_Comm comm)
 
int yac_remote_point_get_pack_size (struct remote_point *point, MPI_Datatype point_info_dt, MPI_Comm comm)
 
void yac_remote_point_pack (struct remote_point *point, void *buffer, int buffer_size, int *position, MPI_Datatype point_info_dt, MPI_Comm comm)
 
void yac_remote_point_unpack (void *buffer, int buffer_size, int *position, struct remote_point *point, MPI_Datatype point_info_dt, MPI_Comm comm)
 
int yac_remote_points_get_pack_size (struct remote_points *points, MPI_Datatype point_info_dt, MPI_Comm comm)
 
void yac_remote_points_pack (struct remote_points *points, void *buffer, int buffer_size, int *position, MPI_Datatype point_info_dt, MPI_Comm comm)
 
void yac_remote_points_unpack (void *buffer, int buffer_size, int *position, struct remote_points **points, MPI_Datatype point_info_dt, MPI_Comm comm)
 
struct remote_pointyac_dist_grid_get_remote_points (struct yac_dist_grid *dist_grid, enum yac_location location, size_t *points, size_t count)
 
void yac_dist_grid_global_to_local (struct yac_dist_grid *dist_grid, enum yac_location location, yac_int *global_ids, size_t count, size_t *local_ids)
 
void yac_dist_grid_pair_get_aux_grid (struct yac_dist_grid_pair *grid_pair, char const *grid_name, size_t *cells, size_t count, size_t **vertex_to_cell, size_t **vertex_to_cell_offsets, int **num_cells_per_vertex, struct yac_interp_field field)
 
void yac_dist_grid_pair_relocate_point_pairs (struct yac_dist_grid_pair *grid_pair, int a_is_ref, int to_dist_owner, char const *grid_name_a, size_t **points_a, enum yac_location location_a, char const *grid_name_b, size_t **points_b, enum yac_location location_b, double **weights, size_t *count)
 
int const * yac_dist_grid_get_field_mask (struct yac_dist_grid *dist_grid, struct yac_interp_field field)
 
yac_const_coordinate_pointer yac_dist_grid_get_field_coords (struct yac_dist_grid *dist_grid, struct yac_interp_field field)
 
void yac_const_basic_grid_data_get_grid_cell (struct yac_const_basic_grid_data *grid_data, size_t cell_idx, struct yac_grid_cell *cell)
 
void yac_dist_grid_pair_determine_dist_owner (struct yac_dist_grid_pair *grid_pair, char const *grid_name, size_t *points, size_t count, enum yac_location location, int *ranks)
 
void yac_dist_grid_pair_get_vertex_neighbours (struct yac_dist_grid_pair *grid_pair, char const *grid_name, size_t *vertices, size_t count, size_t **neigh_vertices, int *num_neighs_per_vertex, struct yac_interp_field field)
 
void yac_dist_grid_pair_get_corner_cells (struct yac_dist_grid_pair *grid_pair, char const *grid_name, size_t *vertices, size_t count, size_t **vertex_to_cell, size_t *num_cells_per_vertex)
 

Typedef Documentation

◆ const_bounding_circle_pointer

typedef struct bounding_circle const* const const_bounding_circle_pointer

Definition at line 17 of file dist_grid_internal.h.

◆ const_int_pointer

typedef int const* const const_int_pointer

Definition at line 13 of file dist_grid_internal.h.

◆ const_remote_point_infos_pointer

Definition at line 18 of file dist_grid_internal.h.

◆ const_size_t_2_pointer

typedef size_t const(* const const_size_t_2_pointer)[2]

Definition at line 12 of file dist_grid_internal.h.

◆ const_size_t_pointer

typedef size_t const* const const_size_t_pointer

Definition at line 14 of file dist_grid_internal.h.

◆ const_yac_edge_type_pointer

typedef enum yac_edge_type const* const const_yac_edge_type_pointer

Definition at line 16 of file dist_grid_internal.h.

◆ const_yac_int_pointer

typedef yac_int const* const const_yac_int_pointer

Definition at line 15 of file dist_grid_internal.h.

Function Documentation

◆ yac_const_basic_grid_data_get_grid_cell()

void yac_const_basic_grid_data_get_grid_cell ( struct yac_const_basic_grid_data grid_data,
size_t  cell_idx,
struct yac_grid_cell cell 
)

returns the data for a selected grid cell

Parameters
[in]grid_databasic grid data
[in]cell_idxlocal id of the selected cell
[out]cellcell data
Remarks
the user has to ensure that the cell was probably initialised
the user has to free the memory associated to cell (see yac_free_grid_cell)

Definition at line 2996 of file dist_grid.c.

Here is the caller graph for this function:

◆ yac_dist_grid_get_basic_grid_data()

struct yac_const_basic_grid_data * yac_dist_grid_get_basic_grid_data ( struct yac_dist_grid dist_grid)

gets the basic grid data of a distributed grid

Parameters
[in]dist_griddistributed grid
Returns
basic grid data
Remarks
the contents of the basic grid data may change by some calls (for example by yac_dist_grid_pair_do_bnd_circle_search)

Definition at line 2796 of file dist_grid.c.

Here is the caller graph for this function:

◆ yac_dist_grid_get_field_coords()

yac_const_coordinate_pointer yac_dist_grid_get_field_coords ( struct yac_dist_grid dist_grid,
struct yac_interp_field  field 
)

returns a pointer to a field coordinates

Parameters
[in]dist_griddistributed grid
[in]fieldfield for which the coordinates are to be returned
Returns
field coordinates (NULL if no mask is available)
Remarks
if the field location is YAC_LOC_CORNER and no field coordinates are available/defined, this routine will return a pointer to the vertex coordinates of the grid

Definition at line 2895 of file dist_grid.c.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ yac_dist_grid_get_field_mask()

int const * yac_dist_grid_get_field_mask ( struct yac_dist_grid dist_grid,
struct yac_interp_field  field 
)

returns a pointer to a field mask

Parameters
[in]dist_griddistributed grid
[in]fieldfield for which the mask is to be returned
Returns
field mask (NULL if no mask is available/defined)

Definition at line 2883 of file dist_grid.c.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ yac_dist_grid_get_local_count()

size_t yac_dist_grid_get_local_count ( struct yac_dist_grid dist_grid,
enum yac_location  location 
)

gets the number of points in the local part of the distributed grid

Parameters
[in]dist_griddistributed grid
[in]locationlocation of the points
Remarks
each global point is only counted once by a single process

Definition at line 2802 of file dist_grid.c.

Here is the caller graph for this function:

◆ yac_dist_grid_get_local_unmasked_points()

void yac_dist_grid_get_local_unmasked_points ( struct yac_dist_grid grid,
struct yac_interp_field  field,
size_t **  indices,
size_t *  count 
)

gets all unmasked points available in the local part of the distributed grid

Parameters
[in]griddistributed grid
[in]fieldfield description
[out]indiceslocal indices of all unmasked local cells
[out]countnumber of entries in indices
Remarks
each unmask point of the global grid is returned only by a single process

Definition at line 2850 of file dist_grid.c.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ yac_dist_grid_get_remote_points()

struct remote_point * yac_dist_grid_get_remote_points ( struct yac_dist_grid dist_grid,
enum yac_location  location,
size_t *  points,
size_t  count 
)

allocates and returns remote_point informations for a list of selected points

Parameters
[in]dist_griddistributed grid
[in]locationlocation of the requested points
[in]pointslocal ids of selected points
[in]countnumber of entries in points
Returns
remote_point information for the selected points
Remarks
the user is responsible for freeing the memory allocated for the remote_point information

Definition at line 5757 of file dist_grid.c.

Here is the caller graph for this function:

◆ yac_dist_grid_global_to_local()

void yac_dist_grid_global_to_local ( struct yac_dist_grid dist_grid,
enum yac_location  location,
yac_int global_ids,
size_t  count,
size_t *  local_ids 
)

returns the local ids for the provided global ids

Parameters
[in]dist_griddistributed grid
[in]locationlocation of the selected points
[in]global_idsglobal ids of the selected points
[in]countnumber of entries in global_ids
[out]local_idslocal ids of the selected points
Remarks
the user has to ensure that the array associated to local_ids is big enough to hold enough elements
in case one or more selected global ids are not available in the local data, the local data will be extended accordingly

Definition at line 5843 of file dist_grid.c.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ yac_dist_grid_pair_determine_dist_owner()

void yac_dist_grid_pair_determine_dist_owner ( struct yac_dist_grid_pair grid_pair,
char const *  grid_name,
size_t *  points,
size_t  count,
enum yac_location  location,
int *  ranks 
)

Determines the ranks of the distributed owners for the selected points

Parameters
[in]grid_pairdistributed grid pair
[in]grid_namegrid name of the grid on which this routine is supposed to work on
[in]pointslocal ids of selected points
[in]countnumber of entries in points
[in]locationlocation of the requested points
[out]ranksdistributed owners of selected points

Definition at line 6714 of file dist_grid.c.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ yac_dist_grid_pair_do_bnd_circle_search()

void yac_dist_grid_pair_do_bnd_circle_search ( struct yac_dist_grid_pair grid_pair,
char const *  grid_name,
const_bounding_circle_pointer  bnd_circles,
size_t  count,
size_t **  cells,
size_t *  num_results_per_bnd_circle,
struct yac_interp_field  field 
)

search for all cells matching the provided bounding circles

Parameters
[in]grid_pairdistributed grid pair
[in]grid_namegrid name of the grid for which the search is to be performed
[in]bnd_circlessearch bounding circles
[in]countnumber of search bounding circles
[out]cellslocal ids of results cells
[out]num_results_per_bnd_circlenumber of results per bounding circle
[in]fieldfield description

Definition at line 5035 of file dist_grid.c.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ yac_dist_grid_pair_do_cell_search()

void yac_dist_grid_pair_do_cell_search ( struct yac_dist_grid_pair grid_pair,
char const *  search_grid_name,
char const *  result_grid_name,
size_t *  search_cells,
size_t  count,
size_t **  result_cells,
size_t *  num_results_per_search_cell,
struct yac_interp_field  result_field 
)

search for matches between cells of the two grids within the grid pair

Parameters
[in]grid_pairdistributed grid pair
[in]search_grid_namegrid name of the grid from which to take the search cells
[in]result_grid_namegrid name of the grid from which to take the matching cells
[in]search_cellslocal ids of search cells
[in]countnumber of search cells
[out]result_cellslocal ids of results cells
[out]num_results_per_search_cellnumber of results per search cell
[in]result_fieldfield description

Definition at line 5323 of file dist_grid.c.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ yac_dist_grid_pair_do_nnn_search()

void yac_dist_grid_pair_do_nnn_search ( struct yac_dist_grid_pair grid_pair,
char const *  grid_name,
yac_coordinate_pointer  search_coords,
size_t  count,
size_t *  local_ids,
size_t  n,
struct yac_interp_field  field 
)

does a n-nearest-neighbour search

Parameters
[in]grid_pairdistributed grid pair
[in]grid_namegrid name of the grid for which the search is to be performed
[in]search_coordssearch coordinates
[in]countnumber of search coordinates
[out]local_idslocal ids of results points
[in]nnumber of points per search coordinate
[in]fieldfield description

Definition at line 4702 of file dist_grid.c.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ yac_dist_grid_pair_do_point_search()

void yac_dist_grid_pair_do_point_search ( struct yac_dist_grid_pair grid_pair,
char const *  grid_name,
yac_coordinate_pointer  search_coords,
size_t  count,
size_t *  cells 
)

search for cells that map to the provided search coordinates

Parameters
[in]grid_pairdistributed grid pair
[in]grid_namegrid name of the grid for which the search is to be performed
[in]search_coordssearch coordinates
[in]countnumber of search coordinates
[out]cellslocal ids of matching cells
Remarks
in case for a search coordinate no matching cell was found, the respective entry in cells is SIZE_MAX

Definition at line 4553 of file dist_grid.c.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ yac_dist_grid_pair_do_point_search_gc()

void yac_dist_grid_pair_do_point_search_gc ( struct yac_dist_grid_pair grid_pair,
char const *  grid_name,
yac_coordinate_pointer  search_coords,
size_t  count,
size_t *  cells 
)

search for cells that map to the provided search coordinates

Parameters
[in]grid_pairdistributed grid pair
[in]grid_namegrid name of the grid for which the search is to be performed
[in]search_coordssearch coordinates
[in]countnumber of search coordinates
[out]cellslocal ids of matching cells
Remarks
SIZE_MAX is the returned result in case no cell was found
assumes that all grid edges are on great circle edges

Definition at line 4561 of file dist_grid.c.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ yac_dist_grid_pair_get_aux_grid()

void yac_dist_grid_pair_get_aux_grid ( struct yac_dist_grid_pair grid_pair,
char const *  grid_name,
size_t *  cells,
size_t  count,
size_t **  vertex_to_cell,
size_t **  vertex_to_cell_offsets,
int **  num_cells_per_vertex,
struct yac_interp_field  field 
)

generates for all vertices of the grid the list of all cells surrounding the vertices
vertices not belonging to the selected cells will have a zero entry in num_cells_per_vertex

Parameters
[in]grid_pairdistributed grid pair
[in]grid_namegrid name of the grid on which this routine is supposed to work on
[in]cellsselected cells
[in]countnumber of entries in cells
[out]vertex_to_celllist of result cells for all vertices associated to the selected cells
[out]vertex_to_cell_offsetsthe results cells for vertex i can be found at:
vertex_to_cell + vertex_to_cell_offsets[i]
[out]num_cells_per_vertexnumber of result cells for each vertex
[in]fieldif the provided field contains a mask for cells, it will be used
Remarks
the result cells for a vertex are sorted in clock- or counter clockwise order
in case a vertex is not fully surrounded by unmasked cells, the respective entry in num_cells_per_vertex is 0
the user is responsible to free the memory returned through the out arrays

Definition at line 6430 of file dist_grid.c.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ yac_dist_grid_pair_get_cell_neighbours()

void yac_dist_grid_pair_get_cell_neighbours ( struct yac_dist_grid_pair grid_pair,
char const *  grid_name,
size_t *  cells,
size_t  count,
size_t *  neighbours 
)

gets neighbouring grid cells for all provided cells

Parameters
[in]grid_pairdistributed grid pair
[in]grid_namegrid name of the grid for which the search is to be performed
[in]cellslocal ids of cells for which the neighbours are to be determined
[in]countnumber for entries in cells
[out]neighboursneighbours for all edges of all provided cells
Remarks
the user has to ensure that the array neighbours has the appropriate size, which is the sum of the number of edges for each provided cell
for each provided cell, this routine sets the neighbours for all edges of the cells
in case there is no neighbour of an edge, the respective entry contains SIZE_MAX
the result cells are sorted in clock- or counter clockwise order

Definition at line 5748 of file dist_grid.c.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ yac_dist_grid_pair_get_corner_cells()

void yac_dist_grid_pair_get_corner_cells ( struct yac_dist_grid_pair grid_pair,
char const *  grid_name,
size_t *  vertices,
size_t  count,
size_t **  vertex_to_cell,
size_t *  num_cells_per_vertex 
)

Determines for all provided vertices the cells surrouding them.

Parameters
[in]grid_pairdistributed grid pair
[in]grid_namegrid name of the grid on which this routine is supposed to work on
[in]verticeslocal ids of selected vertices
[in]countnumber of entries in vertices
[out]vertex_to_celllist of result cells for all vertices associated to the selected cells
[out]num_cells_per_vertexnumber of result cells for each vertex
Remarks
the result cells for a vertex is sorted by their global id
the user is responsible for free the memory of the vertex_to_cell array
the user has to make sure that the array num_cells_per_vertex is at least of size "count"

Definition at line 6380 of file dist_grid.c.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ yac_dist_grid_pair_get_dist_grid()

struct yac_dist_grid * yac_dist_grid_pair_get_dist_grid ( struct yac_dist_grid_pair grid_pair,
char const *  grid_name 
)

gets one of the two grids from the distributed grid pairs

Parameters
[in]grid_pairdistributed grid pair
[in]grid_namegrid name of the selected grid
Returns
distributed grid

Definition at line 2784 of file dist_grid.c.

Here is the caller graph for this function:

◆ yac_dist_grid_pair_get_MPI_Comm()

MPI_Comm yac_dist_grid_pair_get_MPI_Comm ( struct yac_dist_grid_pair grid_pair)

gets a communicator containing all ranks of the distributed grid pair

Parameters
[in]grid_pairdistributed grid pair
Returns
comm communicator

Definition at line 2779 of file dist_grid.c.

Here is the caller graph for this function:

◆ yac_dist_grid_pair_get_vertex_neighbours()

void yac_dist_grid_pair_get_vertex_neighbours ( struct yac_dist_grid_pair grid_pair,
char const *  grid_name,
size_t *  vertices,
size_t  count,
size_t **  neigh_vertices,
int *  num_neighs_per_vertex,
struct yac_interp_field  field 
)

Determines all non-masked neighour vertices for the selected vertices

Parameters
[in]grid_pairdistributed grid pair
[in]grid_namegrid name of the grid on which this routine is supposed to work on
[in]verticeslocal ids of selected vertices
[in]countnumber of entries in vertices
[out]neigh_verticesneighbour vertices
[out]num_neighs_per_vertexnumber of neighbours per vertex
[in]fieldif the provided field contains a mask for corners, it will be used
Remarks
the user has to ensure that num_neighs_per_vertex is allocated big enough
the user is responsible for freeing memory associated to neigh_vertices

Definition at line 6294 of file dist_grid.c.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ yac_dist_grid_pair_relocate_point_pairs()

void yac_dist_grid_pair_relocate_point_pairs ( struct yac_dist_grid_pair grid_pair,
int  a_is_ref,
int  to_dist_owner,
char const *  grid_name_a,
size_t **  points_a,
enum yac_location  location_a,
char const *  grid_name_b,
size_t **  points_b,
enum yac_location  location_b,
double **  weights,
size_t *  count 
)

Each point in a distributed grid has a unique owner rank.
This routine relocates point pairs of two grids. If (a_is_ref != 0), this routine will relocate the point pairs such that each process only has pairs with points_a being owned by the respective process.

Parameters
[in]grid_pairdistributed grid pair
[in]a_is_refdetermines whether points_a or points_b are used as the sorting reference
[in]to_dist_ownerdetermines whether the sorting is based on the dist or orig owner of the respective points
[in]grid_name_agrid name for the a-part of point pairs
[in,out]points_alocal ids of the a-part of the point pairs
[in]location_alocation of the a-part of the point pairs
[in]grid_name_bgrid name for the b-part of point pairs
[in,out]points_blocal ids of the b-part of the point pairs
[in]location_blocation of the b-part of the point pairs
[in,out]weightsweights for all point pairs
[in,out]countnumber of point pairs before and after this routine

Definition at line 6773 of file dist_grid.c.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ yac_get_remote_point_info_mpi_datatype()

MPI_Datatype yac_get_remote_point_info_mpi_datatype ( MPI_Comm  comm)

generates an MPI Datatype for struct remote_point_info

Parameters
[in]commcommunicator
Returns
MPI Datatype for struct remote_point_info
Remarks
the user has to free the returned MPI Datatype using MPI_Type_free

Definition at line 322 of file dist_grid.c.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ yac_remote_point_get_pack_size()

int yac_remote_point_get_pack_size ( struct remote_point point,
MPI_Datatype  point_info_dt,
MPI_Comm  comm 
)

computes the maximum size required by MPI to pack the provided point of type struct remote_point

Parameters
[in]pointpoint for which the pack size is to be determined
[in]point_info_dtMPI Datatype for packing struct point_info
[in]commcommunicator
Returns
maximum packing size

Definition at line 507 of file dist_grid.c.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ yac_remote_point_pack()

void yac_remote_point_pack ( struct remote_point point,
void *  buffer,
int  buffer_size,
int *  position,
MPI_Datatype  point_info_dt,
MPI_Comm  comm 
)

packs a provided remote_point; this is simlar to MPI_Pack

Parameters
[in]pointremote_point to be packed
[in,out]bufferpacking buffer
[in]buffer_sizesize of packing buffer
[in,out]positionpacking position
[in]point_info_dtMPI Datatype for packing struct point_info
[in]commcommunicator

Definition at line 536 of file dist_grid.c.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ yac_remote_point_unpack()

void yac_remote_point_unpack ( void *  buffer,
int  buffer_size,
int *  position,
struct remote_point point,
MPI_Datatype  point_info_dt,
MPI_Comm  comm 
)

unpacks and allocates a remote_point from a buffer; this is similar to MPI_Unpack

Parameters
[in]bufferpacking buffer
[in]buffer_sizesize of packing buffer
[in,out]positionunpacking position
[out]pointunpacked point
[in]point_info_dtMPI Datatype for unpacking struct point_info
[in]commcommunicator

Definition at line 548 of file dist_grid.c.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ yac_remote_points_get_pack_size()

int yac_remote_points_get_pack_size ( struct remote_points points,
MPI_Datatype  point_info_dt,
MPI_Comm  comm 
)

computes the maximum size required by MPI to pack the provided points of type struct remote_points

Parameters
[in]pointspoints for which the pack size is to be determined
[in]point_info_dtMPI Datatype for unpacking struct point_info
[in]commcommunicator
Returns
maximum packing size

Definition at line 574 of file dist_grid.c.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ yac_remote_points_pack()

void yac_remote_points_pack ( struct remote_points points,
void *  buffer,
int  buffer_size,
int *  position,
MPI_Datatype  point_info_dt,
MPI_Comm  comm 
)

packs provided remote_points; this is simlar to MPI_Pack

Parameters
[in]pointsremote_points to be packed
[in,out]bufferpacking buffer
[in]buffer_sizesize of packing buffer
[in,out]positionpacking position
[in]point_info_dtMPI Datatype for packing struct point_info
[in]commcommunicator

Definition at line 592 of file dist_grid.c.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ yac_remote_points_unpack()

void yac_remote_points_unpack ( void *  buffer,
int  buffer_size,
int *  position,
struct remote_points **  points,
MPI_Datatype  point_info_dt,
MPI_Comm  comm 
)

unpacks and allocates remote_points from a buffer; this is similar to MPI_Unpack

Parameters
[in]bufferpacking buffer
[in]buffer_sizesize of packing buffer
[in,out]positionunpacking position
[out]pointsunpacked points
[in]point_info_dtMPI Datatype for unpacking struct point_info
[in]commcommunicator

Definition at line 611 of file dist_grid.c.

Here is the call graph for this function:
Here is the caller graph for this function: