pub unsafe trait Source: AsCommunicator {
Show 19 methods
// Required method
fn source_rank(&self) -> Rank;
// Provided methods
fn probe_with_tag(&self, tag: Tag) -> Status { ... }
fn probe(&self) -> Status { ... }
fn matched_probe_with_tag(&self, tag: Tag) -> (Message, Status) { ... }
fn matched_probe(&self) -> (Message, Status) { ... }
fn receive_with_tag<Msg>(&self, tag: Tag) -> (Msg, Status)
where Msg: Equivalence { ... }
fn receive<Msg>(&self) -> (Msg, Status)
where Msg: Equivalence { ... }
fn receive_into_with_tag<Buf>(&self, buf: &mut Buf, tag: Tag) -> Status
where Buf: BufferMut + ?Sized { ... }
fn receive_into<Buf>(&self, buf: &mut Buf) -> Status
where Buf: BufferMut + ?Sized { ... }
fn receive_vec_with_tag<Msg>(&self, tag: Tag) -> (Vec<Msg>, Status)
where Msg: Equivalence { ... }
fn receive_vec<Msg>(&self) -> (Vec<Msg>, Status)
where Msg: Equivalence { ... }
fn immediate_receive_into_with_tag<'a, Sc, Buf>(
&self,
scope: Sc,
buf: &'a mut Buf,
tag: Tag,
) -> Request<'a, Sc>
where Buf: 'a + BufferMut + ?Sized,
Sc: Scope<'a> { ... }
fn immediate_receive_into<'a, Sc, Buf>(
&self,
scope: Sc,
buf: &'a mut Buf,
) -> Request<'a, Sc>
where Buf: 'a + BufferMut + ?Sized,
Sc: Scope<'a> { ... }
fn immediate_receive_with_tag<Msg>(&self, tag: Tag) -> ReceiveFuture<Msg>
where Msg: Equivalence { ... }
fn immediate_receive<Msg>(&self) -> ReceiveFuture<Msg>
where Msg: Equivalence { ... }
fn immediate_probe_with_tag(&self, tag: Tag) -> Option<Status> { ... }
fn immediate_probe(&self) -> Option<Status> { ... }
fn immediate_matched_probe_with_tag(
&self,
tag: Tag,
) -> Option<(Message, Status)> { ... }
fn immediate_matched_probe(&self) -> Option<(Message, Status)> { ... }
}
Expand description
Required Methods§
Sourcefn source_rank(&self) -> Rank
fn source_rank(&self) -> Rank
Rank
that identifies the source
Provided Methods§
Sourcefn probe_with_tag(&self, tag: Tag) -> Status
fn probe_with_tag(&self, tag: Tag) -> Status
Probe a source for incoming messages.
Probe Source
&self
for incoming messages with a certain tag.
An ordinary probe()
returns a Status
which allows inspection of the properties of the
incoming message, but does not guarantee reception by a subsequent receive()
(especially
in a multi-threaded set-up). For a probe operation with stronger guarantees, see
matched_probe()
.
§Standard section(s)
3.8.1
Sourcefn probe(&self) -> Status
fn probe(&self) -> Status
Probe a source for incoming messages.
Probe Source
&self
for incoming messages with any tag.
An ordinary probe()
returns a Status
which allows inspection of the properties of the
incoming message, but does not guarantee reception by a subsequent receive()
(especially
in a multi-threaded set-up). For a probe operation with stronger guarantees, see
matched_probe()
.
§Standard section(s)
3.8.1
Sourcefn matched_probe_with_tag(&self, tag: Tag) -> (Message, Status)
fn matched_probe_with_tag(&self, tag: Tag) -> (Message, Status)
Probe a source for incoming messages with guaranteed reception.
Probe Source
&self
for incoming messages with a certain tag.
A matched_probe()
returns both a Status
that describes the properties of a pending
incoming message and a Message
which can and must subsequently be used in a
matched_receive()
to receive the probed message.
§Standard section(s)
3.8.2
Sourcefn matched_probe(&self) -> (Message, Status)
fn matched_probe(&self) -> (Message, Status)
Probe a source for incoming messages with guaranteed reception.
Probe Source
&self
for incoming messages with any tag.
A matched_probe()
returns both a Status
that describes the properties of a pending
incoming message and a Message
which can and must subsequently be used in a
matched_receive()
to receive the probed message.
§Standard section(s)
3.8.2
Sourcefn receive_with_tag<Msg>(&self, tag: Tag) -> (Msg, Status)where
Msg: Equivalence,
fn receive_with_tag<Msg>(&self, tag: Tag) -> (Msg, Status)where
Msg: Equivalence,
Receive a message containing a single instance of type Msg
.
Receive a message from Source
&self
tagged tag
containing a single instance of type
Msg
.
§Standard section(s)
3.2.4
Sourcefn receive<Msg>(&self) -> (Msg, Status)where
Msg: Equivalence,
fn receive<Msg>(&self) -> (Msg, Status)where
Msg: Equivalence,
Receive a message containing a single instance of type Msg
.
Receive a message from Source
&self
containing a single instance of type Msg
.
§Examples
extern crate mpi_fork_fnsp as mpi;
use mpi::traits::*;
let universe = mpi::initialize().unwrap();
let world = universe.world();
let x = world.any_process().receive::<f64>();
§Standard section(s)
3.2.4
Sourcefn receive_into_with_tag<Buf>(&self, buf: &mut Buf, tag: Tag) -> Status
fn receive_into_with_tag<Buf>(&self, buf: &mut Buf, tag: Tag) -> Status
Receive a message into a Buffer
.
Receive a message from Source
&self
tagged tag
into Buffer
buf
.
§Standard section(s)
3.2.4
Sourcefn receive_into<Buf>(&self, buf: &mut Buf) -> Status
fn receive_into<Buf>(&self, buf: &mut Buf) -> Status
Receive a message into a Buffer
.
Receive a message from Source
&self
into Buffer
buf
.
§Standard section(s)
3.2.4
Examples found in repository?
6fn main() {
7 let universe = mpi::initialize().unwrap();
8 let world = universe.world();
9 let receiver_rank = 0;
10
11 if world.rank() == receiver_rank {
12 let n = (world.size() - 1) as usize;
13 let mut buf = vec![0u64; 2 * n];
14 for x in buf[0..n].iter_mut() {
15 world.any_process().receive_into(x);
16 }
17 world.barrier();
18 for x in buf[n..2 * n].iter_mut() {
19 world.any_process().receive_into(x);
20 }
21 println!("{:?}", buf);
22 assert!(buf[0..n].iter().all(|&x| x == 1));
23 assert!(buf[n..2 * n].iter().all(|&x| x == 2));
24 } else {
25 world.process_at_rank(0).send(&1u64);
26 world.barrier();
27 world.process_at_rank(0).send(&2u64);
28 }
29}
More examples
11fn main() {
12 use mpi::traits::*;
13 let universe = mpi::initialize().unwrap();
14 let world = universe.world();
15 let size = world.size();
16 let receiver_rank = 0;
17
18 if world.rank() == receiver_rank {
19 // receiver process
20 let n = (size - 1) as usize;
21 let mut buf = vec![0u64; 3 * n];
22 // receive first 2 * n messages
23 for x in buf[0..2 * n].iter_mut() {
24 world.any_process().receive_into(x);
25 }
26 // signal the waiting senders that 2 * n messages have been received
27 let breq = world.immediate_barrier();
28 // receive remaining n messages
29 for x in buf[2 * n..3 * n].iter_mut() {
30 world.any_process().receive_into(x);
31 }
32 println!("{:?}", buf);
33 // messages "1" and "2" may be interleaved, but all have to be contained within the first
34 // 2 * n slots of the buffer
35 assert_eq!(buf[0..2 * n].iter().filter(|&&x| x == 1).count(), n);
36 assert_eq!(buf[0..2 * n].iter().filter(|&&x| x == 2).count(), n);
37 // the last n slots in the buffer may only contain message "3"
38 assert!(buf[2 * n..3 * n].iter().all(|&x| x == 3));
39 // clean up the barrier request
40 breq.wait();
41 } else {
42 // sender processes
43 // send message "1"
44 world.process_at_rank(0).send(&1u64);
45 // join barrier, but do not block
46 let breq = world.immediate_barrier();
47 // send message "2"
48 world.process_at_rank(0).send(&2u64);
49 // wait for receiver process to receive the first 2 * n messages
50 breq.wait();
51 // send message "3"
52 world.process_at_rank(0).send(&3u64);
53 }
54}
Sourcefn receive_vec_with_tag<Msg>(&self, tag: Tag) -> (Vec<Msg>, Status)where
Msg: Equivalence,
fn receive_vec_with_tag<Msg>(&self, tag: Tag) -> (Vec<Msg>, Status)where
Msg: Equivalence,
Receive a message containing multiple instances of type Msg
into a Vec
.
Receive a message from Source
&self
tagged tag
containing multiple instances of type
Msg
into a Vec
.
§Standard section(s)
3.2.4
Sourcefn receive_vec<Msg>(&self) -> (Vec<Msg>, Status)where
Msg: Equivalence,
fn receive_vec<Msg>(&self) -> (Vec<Msg>, Status)where
Msg: Equivalence,
Receive a message containing multiple instances of type Msg
into a Vec
.
Receive a message from Source
&self
containing multiple instances of type Msg
into a
Vec
.
§Examples
See examples/send_receive.rs
§Standard section(s)
3.2.4
Examples found in repository?
7fn main() {
8 let universe = mpi::initialize().unwrap();
9 let world = universe.world();
10 let size = world.size();
11 let rank = world.rank();
12
13 let next_rank = if rank + 1 < size { rank + 1 } else { 0 };
14 let previous_rank = if rank > 0 { rank - 1 } else { size - 1 };
15
16 let msg = vec![rank, 2 * rank, 4 * rank];
17 mpi::request::scope(|scope| {
18 let _sreq = WaitGuard::from(
19 world
20 .process_at_rank(next_rank)
21 .immediate_send(scope, &msg[..]),
22 );
23
24 let (msg, status) = world.any_process().receive_vec();
25
26 println!(
27 "Process {} got message {:?}.\nStatus is: {:?}",
28 rank, msg, status
29 );
30 let x = status.source_rank();
31 assert_eq!(x, previous_rank);
32 assert_eq!(vec![x, 2 * x, 4 * x], msg);
33
34 let root_rank = 0;
35 let root_process = world.process_at_rank(root_rank);
36
37 let mut a;
38 if world.rank() == root_rank {
39 a = vec![2, 4, 8, 16];
40 println!("Root broadcasting value: {:?}.", &a[..]);
41 } else {
42 a = vec![0; 4];
43 }
44 root_process.broadcast_into(&mut a[..]);
45 println!("Rank {} received value: {:?}.", world.rank(), &a[..]);
46 assert_eq!(&a[..], &[2, 4, 8, 16]);
47 });
48}
More examples
8fn main() {
9 let universe = mpi::initialize().unwrap();
10 let world = universe.world();
11 let size = world.size();
12 let rank = world.rank();
13
14 let next_rank = if rank + 1 < size { rank + 1 } else { 0 };
15 let next_process = world.process_at_rank(next_rank);
16 let previous_rank = if rank > 0 { rank - 1 } else { size - 1 };
17 let previous_process = world.process_at_rank(previous_rank);
18
19 let (msg, status): (Rank, _) = p2p::send_receive(&rank, &previous_process, &next_process);
20 println!(
21 "Process {} got message {}.\nStatus is: {:?}",
22 rank, msg, status
23 );
24 world.barrier();
25 assert_eq!(msg, next_rank);
26
27 if rank > 0 {
28 let msg = vec![rank, rank + 1, rank - 1];
29 world.process_at_rank(0).send(&msg[..]);
30 } else {
31 for _ in 1..size {
32 let (msg, status) = world.any_process().receive_vec::<Rank>();
33 println!(
34 "Process {} got long message {:?}.\nStatus is: {:?}",
35 rank, msg, status
36 );
37
38 let x = status.source_rank();
39 let v = vec![x, x + 1, x - 1];
40 assert_eq!(v, msg);
41 }
42 }
43 world.barrier();
44
45 let mut x = rank;
46 p2p::send_receive_replace_into(&mut x, &next_process, &previous_process);
47 assert_eq!(x, previous_rank);
48}
Sourcefn immediate_receive_into_with_tag<'a, Sc, Buf>(
&self,
scope: Sc,
buf: &'a mut Buf,
tag: Tag,
) -> Request<'a, Sc>
fn immediate_receive_into_with_tag<'a, Sc, Buf>( &self, scope: Sc, buf: &'a mut Buf, tag: Tag, ) -> Request<'a, Sc>
Initiate an immediate (non-blocking) receive operation.
Initiate receiving a message matching tag
into buf
.
§Standard section(s)
3.7.2
Sourcefn immediate_receive_into<'a, Sc, Buf>(
&self,
scope: Sc,
buf: &'a mut Buf,
) -> Request<'a, Sc>
fn immediate_receive_into<'a, Sc, Buf>( &self, scope: Sc, buf: &'a mut Buf, ) -> Request<'a, Sc>
Initiate an immediate (non-blocking) receive operation.
Initiate receiving a message into buf
.
§Examples
See examples/immediate.rs
§Standard section(s)
3.7.2
Examples found in repository?
9fn main() {
10 let mut universe = mpi::initialize().unwrap();
11 // Try to attach a buffer.
12 universe.set_buffer_size(BUFFER_SIZE);
13 // Check buffer size matches.
14 assert_eq!(universe.buffer_size(), BUFFER_SIZE);
15 // Try to detach the buffer.
16 universe.detach_buffer();
17 // Attach another buffer.
18 universe.set_buffer_size(BUFFER_SIZE);
19
20 let world = universe.world();
21
22 let x = vec![std::f32::consts::PI; 1024];
23 let mut y = vec![0.0; 1024];
24 mpi::request::scope(|scope| {
25 let _rreq = WaitGuard::from(
26 world
27 .any_process()
28 .immediate_receive_into(scope, &mut y[..]),
29 );
30 world.this_process().buffered_send(&x[..]);
31 });
32 assert_eq!(x, y);
33}
More examples
7fn main() {
8 let universe = mpi::initialize().unwrap();
9 let world = universe.world();
10 let size = world.size();
11 let rank = world.rank();
12
13 if rank > 0 {
14 let msg = rank as u8;
15 world.barrier();
16 world.process_at_rank(0).ready_send(&msg);
17 } else {
18 let mut v = vec![0u8; (size - 1) as usize];
19 mpi::request::scope(|scope| {
20 let reqs = v
21 .iter_mut()
22 .zip(1..)
23 .map(|(x, i)| {
24 world
25 .process_at_rank(i as Rank)
26 .immediate_receive_into(scope, x)
27 })
28 .collect::<Vec<_>>();
29 world.barrier();
30 for req in reqs {
31 req.wait();
32 }
33 });
34 println!("Got message: {:?}", v);
35 assert!(v.iter().zip(1..).all(|(x, i)| i == *x as usize));
36 }
37}
8fn main() {
9 let universe = mpi::initialize().unwrap();
10 let world = universe.world();
11
12 let x = std::f32::consts::PI;
13 let mut y: f32 = 0.0;
14
15 mpi::request::scope(|scope| {
16 if world.rank() == 0 {
17 let mut requests = Vec::new();
18 for i in 1..world.size() {
19 requests.push(
20 world
21 .process_at_rank(i)
22 .immediate_synchronous_send(scope, &x),
23 );
24 }
25
26 println!("World size {}", world.size());
27 while let Some((index, _status)) = mpi::request::wait_any(&mut requests) {
28 println!("Request with index {} completed", index);
29 }
30 println!("All requests completed");
31 } else {
32 let secs = time::Duration::from_secs(world.rank() as u64);
33
34 thread::sleep(secs);
35
36 let rreq = world.any_process().immediate_receive_into(scope, &mut y);
37 rreq.wait();
38 println!("Process {} received data", world.rank());
39 }
40 });
41}
8fn main() {
9 let universe = mpi::initialize().unwrap();
10 let world = universe.world();
11
12 let x = std::f32::consts::PI;
13 let mut y: f32 = 0.0;
14
15 mpi::request::scope(|scope| {
16 let mut sreq = world.this_process().immediate_send(scope, &x);
17 let rreq = world.any_process().immediate_receive_into(scope, &mut y);
18 rreq.wait();
19 loop {
20 match sreq.test() {
21 Ok(_) => {
22 break;
23 }
24 Err(req) => {
25 sreq = req;
26 }
27 }
28 }
29 });
30 assert_eq!(x, y);
31
32 y = 0.0;
33 mpi::request::scope(|scope| {
34 let _rreq = WaitGuard::from(world.any_process().immediate_receive_into(scope, &mut y));
35 let _sreq = WaitGuard::from(world.this_process().immediate_ready_send(scope, &x));
36 });
37 assert_eq!(x, y);
38
39 assert!(world.any_process().immediate_probe().is_none());
40 assert!(world.any_process().immediate_matched_probe().is_none());
41
42 y = 0.0;
43 mpi::request::scope(|scope| {
44 let _sreq: WaitGuard<_> = world
45 .this_process()
46 .immediate_synchronous_send(scope, &x)
47 .into();
48 let preq = world.any_process().immediate_matched_probe();
49 assert!(preq.is_some());
50 let (msg, _) = preq.unwrap();
51 let _rreq: WaitGuard<_> = msg.immediate_matched_receive_into(scope, &mut y).into();
52 });
53 assert_eq!(x, y);
54
55 let future = world.any_process().immediate_receive();
56 world.this_process().send(&x);
57 let (msg, _) = future.get();
58 assert_eq!(x, msg);
59
60 let future = world.any_process().immediate_receive();
61 let res = future.r#try();
62 assert!(res.is_err());
63 let mut future = res.err().unwrap();
64 world.this_process().send(&x);
65 loop {
66 match future.r#try() {
67 Ok((msg, _)) => {
68 assert_eq!(x, msg);
69 break;
70 }
71 Err(f) => {
72 future = f;
73 }
74 }
75 }
76
77 mpi::request::scope(|scope| {
78 let sreq = world.this_process().immediate_send(scope, &x);
79 sreq.cancel();
80 sreq.wait();
81
82 let _sreq = CancelGuard::from(world.this_process().immediate_receive_into(scope, &mut y));
83 });
84}
Sourcefn immediate_receive_with_tag<Msg>(&self, tag: Tag) -> ReceiveFuture<Msg>where
Msg: Equivalence,
fn immediate_receive_with_tag<Msg>(&self, tag: Tag) -> ReceiveFuture<Msg>where
Msg: Equivalence,
Sourcefn immediate_receive<Msg>(&self) -> ReceiveFuture<Msg>where
Msg: Equivalence,
fn immediate_receive<Msg>(&self) -> ReceiveFuture<Msg>where
Msg: Equivalence,
Initiate a non-blocking receive operation.
§Examples
See examples/immediate.rs
§Standard section(s)
3.7.2
Examples found in repository?
8fn main() {
9 let universe = mpi::initialize().unwrap();
10 let world = universe.world();
11
12 let x = std::f32::consts::PI;
13 let mut y: f32 = 0.0;
14
15 mpi::request::scope(|scope| {
16 let mut sreq = world.this_process().immediate_send(scope, &x);
17 let rreq = world.any_process().immediate_receive_into(scope, &mut y);
18 rreq.wait();
19 loop {
20 match sreq.test() {
21 Ok(_) => {
22 break;
23 }
24 Err(req) => {
25 sreq = req;
26 }
27 }
28 }
29 });
30 assert_eq!(x, y);
31
32 y = 0.0;
33 mpi::request::scope(|scope| {
34 let _rreq = WaitGuard::from(world.any_process().immediate_receive_into(scope, &mut y));
35 let _sreq = WaitGuard::from(world.this_process().immediate_ready_send(scope, &x));
36 });
37 assert_eq!(x, y);
38
39 assert!(world.any_process().immediate_probe().is_none());
40 assert!(world.any_process().immediate_matched_probe().is_none());
41
42 y = 0.0;
43 mpi::request::scope(|scope| {
44 let _sreq: WaitGuard<_> = world
45 .this_process()
46 .immediate_synchronous_send(scope, &x)
47 .into();
48 let preq = world.any_process().immediate_matched_probe();
49 assert!(preq.is_some());
50 let (msg, _) = preq.unwrap();
51 let _rreq: WaitGuard<_> = msg.immediate_matched_receive_into(scope, &mut y).into();
52 });
53 assert_eq!(x, y);
54
55 let future = world.any_process().immediate_receive();
56 world.this_process().send(&x);
57 let (msg, _) = future.get();
58 assert_eq!(x, msg);
59
60 let future = world.any_process().immediate_receive();
61 let res = future.r#try();
62 assert!(res.is_err());
63 let mut future = res.err().unwrap();
64 world.this_process().send(&x);
65 loop {
66 match future.r#try() {
67 Ok((msg, _)) => {
68 assert_eq!(x, msg);
69 break;
70 }
71 Err(f) => {
72 future = f;
73 }
74 }
75 }
76
77 mpi::request::scope(|scope| {
78 let sreq = world.this_process().immediate_send(scope, &x);
79 sreq.cancel();
80 sreq.wait();
81
82 let _sreq = CancelGuard::from(world.this_process().immediate_receive_into(scope, &mut y));
83 });
84}
Sourcefn immediate_probe_with_tag(&self, tag: Tag) -> Option<Status>
fn immediate_probe_with_tag(&self, tag: Tag) -> Option<Status>
Asynchronously probe a source for incoming messages.
Asynchronously probe Source
&self
for incoming messages with a certain tag.
Like Probe
but returns a None
immediately if there is no incoming message to be probed.
§Standard section(s)
3.8.1
Sourcefn immediate_probe(&self) -> Option<Status>
fn immediate_probe(&self) -> Option<Status>
Asynchronously probe a source for incoming messages.
Asynchronously probe Source
&self
for incoming messages with any tag.
Like Probe
but returns a None
immediately if there is no incoming message to be probed.
§Standard section(s)
3.8.1
Examples found in repository?
8fn main() {
9 let universe = mpi::initialize().unwrap();
10 let world = universe.world();
11
12 let x = std::f32::consts::PI;
13 let mut y: f32 = 0.0;
14
15 mpi::request::scope(|scope| {
16 let mut sreq = world.this_process().immediate_send(scope, &x);
17 let rreq = world.any_process().immediate_receive_into(scope, &mut y);
18 rreq.wait();
19 loop {
20 match sreq.test() {
21 Ok(_) => {
22 break;
23 }
24 Err(req) => {
25 sreq = req;
26 }
27 }
28 }
29 });
30 assert_eq!(x, y);
31
32 y = 0.0;
33 mpi::request::scope(|scope| {
34 let _rreq = WaitGuard::from(world.any_process().immediate_receive_into(scope, &mut y));
35 let _sreq = WaitGuard::from(world.this_process().immediate_ready_send(scope, &x));
36 });
37 assert_eq!(x, y);
38
39 assert!(world.any_process().immediate_probe().is_none());
40 assert!(world.any_process().immediate_matched_probe().is_none());
41
42 y = 0.0;
43 mpi::request::scope(|scope| {
44 let _sreq: WaitGuard<_> = world
45 .this_process()
46 .immediate_synchronous_send(scope, &x)
47 .into();
48 let preq = world.any_process().immediate_matched_probe();
49 assert!(preq.is_some());
50 let (msg, _) = preq.unwrap();
51 let _rreq: WaitGuard<_> = msg.immediate_matched_receive_into(scope, &mut y).into();
52 });
53 assert_eq!(x, y);
54
55 let future = world.any_process().immediate_receive();
56 world.this_process().send(&x);
57 let (msg, _) = future.get();
58 assert_eq!(x, msg);
59
60 let future = world.any_process().immediate_receive();
61 let res = future.r#try();
62 assert!(res.is_err());
63 let mut future = res.err().unwrap();
64 world.this_process().send(&x);
65 loop {
66 match future.r#try() {
67 Ok((msg, _)) => {
68 assert_eq!(x, msg);
69 break;
70 }
71 Err(f) => {
72 future = f;
73 }
74 }
75 }
76
77 mpi::request::scope(|scope| {
78 let sreq = world.this_process().immediate_send(scope, &x);
79 sreq.cancel();
80 sreq.wait();
81
82 let _sreq = CancelGuard::from(world.this_process().immediate_receive_into(scope, &mut y));
83 });
84}
Sourcefn immediate_matched_probe_with_tag(
&self,
tag: Tag,
) -> Option<(Message, Status)>
fn immediate_matched_probe_with_tag( &self, tag: Tag, ) -> Option<(Message, Status)>
Asynchronously probe a source for incoming messages with guaranteed reception.
Asynchronously probe Source
&self
for incoming messages with a certain tag.
Like MatchedProbe
but returns a None
immediately if there is no incoming message to be
probed.
§Standard section(s)
3.8.2
Sourcefn immediate_matched_probe(&self) -> Option<(Message, Status)>
fn immediate_matched_probe(&self) -> Option<(Message, Status)>
Asynchronously probe a source for incoming messages with guaranteed reception.
Asynchronously probe Source
&self
for incoming messages with any tag.
Like MatchedProbe
but returns a None
immediately if there is no incoming message to be
probed.
§Standard section(s)
3.8.2
Examples found in repository?
8fn main() {
9 let universe = mpi::initialize().unwrap();
10 let world = universe.world();
11
12 let x = std::f32::consts::PI;
13 let mut y: f32 = 0.0;
14
15 mpi::request::scope(|scope| {
16 let mut sreq = world.this_process().immediate_send(scope, &x);
17 let rreq = world.any_process().immediate_receive_into(scope, &mut y);
18 rreq.wait();
19 loop {
20 match sreq.test() {
21 Ok(_) => {
22 break;
23 }
24 Err(req) => {
25 sreq = req;
26 }
27 }
28 }
29 });
30 assert_eq!(x, y);
31
32 y = 0.0;
33 mpi::request::scope(|scope| {
34 let _rreq = WaitGuard::from(world.any_process().immediate_receive_into(scope, &mut y));
35 let _sreq = WaitGuard::from(world.this_process().immediate_ready_send(scope, &x));
36 });
37 assert_eq!(x, y);
38
39 assert!(world.any_process().immediate_probe().is_none());
40 assert!(world.any_process().immediate_matched_probe().is_none());
41
42 y = 0.0;
43 mpi::request::scope(|scope| {
44 let _sreq: WaitGuard<_> = world
45 .this_process()
46 .immediate_synchronous_send(scope, &x)
47 .into();
48 let preq = world.any_process().immediate_matched_probe();
49 assert!(preq.is_some());
50 let (msg, _) = preq.unwrap();
51 let _rreq: WaitGuard<_> = msg.immediate_matched_receive_into(scope, &mut y).into();
52 });
53 assert_eq!(x, y);
54
55 let future = world.any_process().immediate_receive();
56 world.this_process().send(&x);
57 let (msg, _) = future.get();
58 assert_eq!(x, msg);
59
60 let future = world.any_process().immediate_receive();
61 let res = future.r#try();
62 assert!(res.is_err());
63 let mut future = res.err().unwrap();
64 world.this_process().send(&x);
65 loop {
66 match future.r#try() {
67 Ok((msg, _)) => {
68 assert_eq!(x, msg);
69 break;
70 }
71 Err(f) => {
72 future = f;
73 }
74 }
75 }
76
77 mpi::request::scope(|scope| {
78 let sreq = world.this_process().immediate_send(scope, &x);
79 sreq.cancel();
80 sreq.wait();
81
82 let _sreq = CancelGuard::from(world.this_process().immediate_receive_into(scope, &mut y));
83 });
84}
Dyn Compatibility§
This trait is not dyn compatible.
In older versions of Rust, dyn compatibility was called "object safety", so this trait is not object safe.