PULSAR  2.0.0
Parallel Ultra-Light Systolic Array Runtime
 All Data Structures Files Functions Typedefs Enumerations Macros Groups
prt_request.c
Go to the documentation of this file.
1 
11 #include "prt_request.h"
12 
14 
26  prt_packet_t *packet, size_t size, int peer, int tag)
27 {
28  prt_request_t *request = (prt_request_t*)malloc(sizeof(prt_request_t));
29  prt_assert(request != NULL, "malloc failed");
30  request->packet = packet;
31  request->size = size;
32  request->peer = peer;
33  request->tag = tag;
34  return request;
35 }
36 
38 
44 {
45  free(request);
46 }
47 
49 
55 {
56  int retval;
57  // Post the request.
59  retval = MPI_Isend(
60  request->packet->data,
61  request->size,
62  MPI_BYTE,
63  request->peer,
64  request->tag,
65  MPI_COMM_WORLD,
66  &request->request);
67  svg_trace_stop_cpu(0, -DarkKhaki);
68  prt_assert(retval == MPI_SUCCESS, "error in MPI_Isend");
69 }
70 
72 
78 {
79  int retval;
80  // Post the request.
82  retval = MPI_Irecv(
83  request->packet->data,
84  request->size,
85  MPI_BYTE,
86  request->peer,
87  request->tag,
88  MPI_COMM_WORLD,
89  &request->request);
90  svg_trace_stop_cpu(0, -OliveDrab);
91  prt_assert(retval == MPI_SUCCESS, "error in MPI_Isend");
92 }
93 
95 
105 {
106  int flag;
107  int retval;
108  // Test the request.
109  // Check the status.
110  // Trace if completed.
112  retval = MPI_Test(&request->request, &flag, &request->status);
113  if (flag == 1)
114  svg_trace_stop_cpu(0, -RosyBrown);
115  prt_assert(retval == MPI_SUCCESS, "error in MPI_Test");
116  return flag;
117 }
118 
120 
127 {
128  int retval;
129  // Cancel the request.
130  retval = MPI_Cancel(&request->request);
131  prt_assert(retval == MPI_SUCCESS, "error in MPI_Cancel");
132 
133  // Release the packet.
134  prt_packet_release_host(request->packet);
135 
136  // Free the request.
137  free(request);
138 }