1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
//! Nonblocking collective operations: ibroadcast, iallreduce, ireduce, igather, etc.
use crate::comm::Communicator;
use crate::datatype::MpiDatatype;
use crate::error::{Error, Result};
use crate::ffi;
use crate::request::Request;
use crate::ReduceOp;
impl Communicator {
// ========================================================================
// Generic Nonblocking Collectives
// ========================================================================
/// Nonblocking broadcast.
///
/// Returns a request handle that must be waited on before accessing the buffer.
///
/// # Safety Note
///
/// The buffer must remain valid until the request is completed.
///
/// # Example
///
/// ```no_run
/// # use ferrompi::Mpi;
/// # let mpi = Mpi::init().unwrap();
/// # let world = mpi.world();
/// let mut data = vec![0.0f64; 100];
/// let req = world.ibroadcast(&mut data, 0).unwrap();
/// // ... do other work ...
/// req.wait().unwrap();
/// ```
pub fn ibroadcast<T: MpiDatatype>(&self, data: &mut [T], root: i32) -> Result<Request> {
let mut request_handle: i64 = 0;
let ret = unsafe {
ffi::ferrompi_ibcast(
data.as_mut_ptr().cast::<std::ffi::c_void>(),
data.len() as i64,
T::TAG as i32,
root,
self.handle,
&mut request_handle,
)
};
Error::check_with_op(ret, "ibcast")?;
Ok(Request::new(request_handle))
}
/// Nonblocking all-reduce.
///
/// Returns a request handle that must be waited on before accessing the buffer.
///
/// # Example
///
/// ```no_run
/// # use ferrompi::{Mpi, ReduceOp};
/// # let mpi = Mpi::init().unwrap();
/// # let world = mpi.world();
/// let send = vec![1.0f64; 10];
/// let mut recv = vec![0.0f64; 10];
/// let req = world.iallreduce(&send, &mut recv, ReduceOp::Sum).unwrap();
/// req.wait().unwrap();
/// ```
pub fn iallreduce<T: MpiDatatype>(
&self,
send: &[T],
recv: &mut [T],
op: ReduceOp,
) -> Result<Request> {
if send.len() != recv.len() {
return Err(Error::InvalidBuffer);
}
let mut request_handle: i64 = 0;
let ret = unsafe {
ffi::ferrompi_iallreduce(
send.as_ptr().cast::<std::ffi::c_void>(),
recv.as_mut_ptr().cast::<std::ffi::c_void>(),
send.len() as i64,
T::TAG as i32,
op as i32,
self.handle,
&mut request_handle,
)
};
Error::check_with_op(ret, "iallreduce")?;
Ok(Request::new(request_handle))
}
/// Nonblocking reduce to root.
///
/// Initiates a reduction operation and returns immediately with a [`Request`]
/// handle. The buffers must remain valid until the request is completed.
///
/// # Arguments
///
/// * `send` - Data to send from this process
/// * `recv` - Buffer for result (only significant at root)
/// * `op` - Reduction operation
/// * `root` - Rank of the root process
///
/// # Example
///
/// ```no_run
/// # use ferrompi::{Mpi, ReduceOp};
/// # let mpi = Mpi::init().unwrap();
/// # let world = mpi.world();
/// let send = vec![1.0f64; 10];
/// let mut recv = vec![0.0f64; 10];
/// let req = world.ireduce(&send, &mut recv, ReduceOp::Sum, 0).unwrap();
/// req.wait().unwrap();
/// ```
pub fn ireduce<T: MpiDatatype>(
&self,
send: &[T],
recv: &mut [T],
op: ReduceOp,
root: i32,
) -> Result<Request> {
if send.len() != recv.len() {
return Err(Error::InvalidBuffer);
}
let mut request_handle: i64 = 0;
let ret = unsafe {
ffi::ferrompi_ireduce(
send.as_ptr().cast::<std::ffi::c_void>(),
recv.as_mut_ptr().cast::<std::ffi::c_void>(),
send.len() as i64,
T::TAG as i32,
op as i32,
root,
self.handle,
&mut request_handle,
)
};
Error::check_with_op(ret, "ireduce")?;
Ok(Request::new(request_handle))
}
/// Nonblocking gather to root.
///
/// Initiates a gather operation and returns immediately with a [`Request`]
/// handle. Each process sends `send.len()` elements. Root receives
/// `send.len() * size` elements total.
///
/// # Arguments
///
/// * `send` - Data to send from this process
/// * `recv` - Buffer for received data (only significant at root)
/// * `root` - Rank of the root process
///
/// # Example
///
/// ```no_run
/// # use ferrompi::Mpi;
/// # let mpi = Mpi::init().unwrap();
/// # let world = mpi.world();
/// let send = vec![world.rank() as f64; 5];
/// let mut recv = vec![0.0f64; 5 * world.size() as usize];
/// let req = world.igather(&send, &mut recv, 0).unwrap();
/// req.wait().unwrap();
/// ```
pub fn igather<T: MpiDatatype>(
&self,
send: &[T],
recv: &mut [T],
root: i32,
) -> Result<Request> {
let mut request_handle: i64 = 0;
let ret = unsafe {
ffi::ferrompi_igather(
send.as_ptr().cast::<std::ffi::c_void>(),
send.len() as i64,
recv.as_mut_ptr().cast::<std::ffi::c_void>(),
send.len() as i64,
T::TAG as i32,
root,
self.handle,
&mut request_handle,
)
};
Error::check_with_op(ret, "igather")?;
Ok(Request::new(request_handle))
}
/// Nonblocking all-gather.
///
/// Initiates an all-gather operation and returns immediately with a
/// [`Request`] handle. Each process sends `send.len()` elements and
/// receives from all.
///
/// # Example
///
/// ```no_run
/// # use ferrompi::Mpi;
/// # let mpi = Mpi::init().unwrap();
/// # let world = mpi.world();
/// let send = vec![world.rank() as i32; 3];
/// let mut recv = vec![0i32; 3 * world.size() as usize];
/// let req = world.iallgather(&send, &mut recv).unwrap();
/// req.wait().unwrap();
/// ```
pub fn iallgather<T: MpiDatatype>(&self, send: &[T], recv: &mut [T]) -> Result<Request> {
let mut request_handle: i64 = 0;
let ret = unsafe {
ffi::ferrompi_iallgather(
send.as_ptr().cast::<std::ffi::c_void>(),
send.len() as i64,
recv.as_mut_ptr().cast::<std::ffi::c_void>(),
send.len() as i64,
T::TAG as i32,
self.handle,
&mut request_handle,
)
};
Error::check_with_op(ret, "iallgather")?;
Ok(Request::new(request_handle))
}
/// Nonblocking scatter from root.
///
/// Initiates a scatter operation and returns immediately with a [`Request`]
/// handle. Root sends `recv.len() * size` elements total, each process
/// receives `recv.len()` elements.
///
/// # Example
///
/// ```no_run
/// # use ferrompi::Mpi;
/// # let mpi = Mpi::init().unwrap();
/// # let world = mpi.world();
/// let send = vec![0.0f64; 5 * world.size() as usize];
/// let mut recv = vec![0.0f64; 5];
/// let req = world.iscatter(&send, &mut recv, 0).unwrap();
/// req.wait().unwrap();
/// ```
pub fn iscatter<T: MpiDatatype>(
&self,
send: &[T],
recv: &mut [T],
root: i32,
) -> Result<Request> {
let mut request_handle: i64 = 0;
let ret = unsafe {
ffi::ferrompi_iscatter(
send.as_ptr().cast::<std::ffi::c_void>(),
recv.len() as i64,
recv.as_mut_ptr().cast::<std::ffi::c_void>(),
recv.len() as i64,
T::TAG as i32,
root,
self.handle,
&mut request_handle,
)
};
Error::check_with_op(ret, "iscatter")?;
Ok(Request::new(request_handle))
}
/// Nonblocking barrier.
///
/// Initiates a barrier synchronization and returns immediately with a
/// [`Request`] handle. The barrier is complete when the request is waited on.
///
/// # Example
///
/// ```no_run
/// # use ferrompi::Mpi;
/// # let mpi = Mpi::init().unwrap();
/// # let world = mpi.world();
/// let req = world.ibarrier().unwrap();
/// // ... do other work ...
/// req.wait().unwrap();
/// ```
pub fn ibarrier(&self) -> Result<Request> {
let mut request_handle: i64 = 0;
let ret = unsafe { ffi::ferrompi_ibarrier(self.handle, &mut request_handle) };
Error::check_with_op(ret, "ibarrier")?;
Ok(Request::new(request_handle))
}
/// Nonblocking inclusive prefix reduction (scan).
///
/// Initiates an inclusive scan and returns immediately with a [`Request`]
/// handle. On rank `i`, `recv` will contain the reduction of `send` values
/// from ranks `0..=i` once the request completes.
///
/// # Errors
///
/// Returns [`Error::InvalidBuffer`] if `send.len() != recv.len()`.
///
/// # Example
///
/// ```no_run
/// # use ferrompi::{Mpi, ReduceOp};
/// # let mpi = Mpi::init().unwrap();
/// # let world = mpi.world();
/// let send = vec![1.0f64; 10];
/// let mut recv = vec![0.0f64; 10];
/// let req = world.iscan(&send, &mut recv, ReduceOp::Sum).unwrap();
/// req.wait().unwrap();
/// ```
pub fn iscan<T: MpiDatatype>(
&self,
send: &[T],
recv: &mut [T],
op: ReduceOp,
) -> Result<Request> {
if send.len() != recv.len() {
return Err(Error::InvalidBuffer);
}
let mut request_handle: i64 = 0;
let ret = unsafe {
ffi::ferrompi_iscan(
send.as_ptr().cast::<std::ffi::c_void>(),
recv.as_mut_ptr().cast::<std::ffi::c_void>(),
send.len() as i64,
T::TAG as i32,
op as i32,
self.handle,
&mut request_handle,
)
};
Error::check_with_op(ret, "iscan")?;
Ok(Request::new(request_handle))
}
/// Nonblocking exclusive prefix reduction (exscan).
///
/// Initiates an exclusive scan and returns immediately with a [`Request`]
/// handle. On rank `i`, `recv` will contain the reduction of `send` values
/// from ranks `0..i` once the request completes.
///
/// # Rank 0 Behavior
///
/// **Per the MPI standard, the contents of `recv` on rank 0 are undefined.**
///
/// # Errors
///
/// Returns [`Error::InvalidBuffer`] if `send.len() != recv.len()`.
///
/// # Example
///
/// ```no_run
/// # use ferrompi::{Mpi, ReduceOp};
/// # let mpi = Mpi::init().unwrap();
/// # let world = mpi.world();
/// let send = vec![1.0f64; 10];
/// let mut recv = vec![0.0f64; 10];
/// let req = world.iexscan(&send, &mut recv, ReduceOp::Sum).unwrap();
/// req.wait().unwrap();
/// ```
pub fn iexscan<T: MpiDatatype>(
&self,
send: &[T],
recv: &mut [T],
op: ReduceOp,
) -> Result<Request> {
if send.len() != recv.len() {
return Err(Error::InvalidBuffer);
}
let mut request_handle: i64 = 0;
let ret = unsafe {
ffi::ferrompi_iexscan(
send.as_ptr().cast::<std::ffi::c_void>(),
recv.as_mut_ptr().cast::<std::ffi::c_void>(),
send.len() as i64,
T::TAG as i32,
op as i32,
self.handle,
&mut request_handle,
)
};
Error::check_with_op(ret, "iexscan")?;
Ok(Request::new(request_handle))
}
/// Nonblocking all-to-all personalized communication.
///
/// Initiates an all-to-all operation and returns immediately with a
/// [`Request`] handle.
///
/// # Errors
///
/// Returns [`Error::InvalidBuffer`] if `send.len() != recv.len()` or
/// `send.len()` is not evenly divisible by the communicator size.
///
/// # Example
///
/// ```no_run
/// # use ferrompi::Mpi;
/// # let mpi = Mpi::init().unwrap();
/// # let world = mpi.world();
/// let size = world.size() as usize;
/// let send = vec![world.rank() as f64; size * 3];
/// let mut recv = vec![0.0f64; size * 3];
/// let req = world.ialltoall(&send, &mut recv).unwrap();
/// req.wait().unwrap();
/// ```
pub fn ialltoall<T: MpiDatatype>(&self, send: &[T], recv: &mut [T]) -> Result<Request> {
let size = self.size() as usize;
if send.len() != recv.len() || send.len() % size != 0 {
return Err(Error::InvalidBuffer);
}
let count = (send.len() / size) as i64;
let mut request_handle: i64 = 0;
let ret = unsafe {
ffi::ferrompi_ialltoall(
send.as_ptr().cast::<std::ffi::c_void>(),
count,
recv.as_mut_ptr().cast::<std::ffi::c_void>(),
count,
T::TAG as i32,
self.handle,
&mut request_handle,
)
};
Error::check_with_op(ret, "ialltoall")?;
Ok(Request::new(request_handle))
}
/// Nonblocking reduce-scatter with uniform block size.
///
/// Initiates a reduce-scatter operation and returns immediately with a
/// [`Request`] handle. Performs an element-wise reduction across all
/// processes, then scatters the result so that each process receives
/// `recv.len()` elements.
///
/// # Errors
///
/// Returns [`Error::InvalidBuffer`] if `send.len() != recv.len() * size`.
///
/// # Example
///
/// ```no_run
/// # use ferrompi::{Mpi, ReduceOp};
/// # let mpi = Mpi::init().unwrap();
/// # let world = mpi.world();
/// let size = world.size() as usize;
/// let send = vec![1.0f64; size * 5];
/// let mut recv = vec![0.0f64; 5];
/// let req = world.ireduce_scatter_block(&send, &mut recv, ReduceOp::Sum).unwrap();
/// req.wait().unwrap();
/// ```
pub fn ireduce_scatter_block<T: MpiDatatype>(
&self,
send: &[T],
recv: &mut [T],
op: ReduceOp,
) -> Result<Request> {
let size = self.size() as usize;
if send.len() != recv.len() * size {
return Err(Error::InvalidBuffer);
}
let mut request_handle: i64 = 0;
let ret = unsafe {
ffi::ferrompi_ireduce_scatter_block(
send.as_ptr().cast::<std::ffi::c_void>(),
recv.as_mut_ptr().cast::<std::ffi::c_void>(),
recv.len() as i64,
T::TAG as i32,
op as i32,
self.handle,
&mut request_handle,
)
};
Error::check_with_op(ret, "ireduce_scatter_block")?;
Ok(Request::new(request_handle))
}
/// Nonblocking in-place gather at root. Non-root ranks must use
/// `igather` — this method returns `Error::InvalidOp` on non-root.
///
/// # Buffer Layout (root)
///
/// `data` must have length `recvcount * size()` where `recvcount` is
/// the per-rank count. Rank `r`'s contribution lives at offset
/// `r * recvcount`. Root's own contribution must be pre-written into
/// `data[rank() * recvcount .. (rank()+1) * recvcount]` before the call.
///
/// # Errors
///
/// - `Error::InvalidOp` if this rank is not `root`.
/// - `Error::InvalidBuffer` if `data.len()` is not divisible by `size()`.
///
/// # Example
///
/// ```no_run
/// # use ferrompi::Mpi;
/// # let mpi = Mpi::init().unwrap();
/// # let world = mpi.world();
/// if world.rank() == 0 {
/// let mut data = vec![0i32; 4 * world.size() as usize];
/// let req = world.igather_inplace(&mut data, 0).unwrap();
/// req.wait().unwrap();
/// }
/// ```
pub fn igather_inplace<T: MpiDatatype>(&self, data: &mut [T], root: i32) -> Result<Request> {
if self.rank() != root {
return Err(Error::InvalidOp);
}
let size = self.size() as usize;
if size == 0 || data.len() % size != 0 {
return Err(Error::InvalidBuffer);
}
let recvcount = (data.len() / size) as i64;
let mut request_handle: i64 = 0;
let ret = unsafe {
// SAFETY: data is a valid, exclusively-owned mutable slice. We cast to *mut c_void
// as required by the C FFI. The guard above guarantees self.rank() == root, so
// is_root is hardcoded to 1. The buffer outlives the returned Request handle;
// the caller must call wait() before accessing or dropping the buffer.
ffi::ferrompi_igather_inplace(
data.as_mut_ptr().cast::<std::ffi::c_void>(),
recvcount,
T::TAG as i32,
root,
1, // is_root = true by the rank guard above
self.handle,
&mut request_handle,
)
};
Error::check_with_op(ret, "igather_inplace")?;
Ok(Request::new(request_handle))
}
/// Nonblocking in-place all-gather. Every rank's `data` is both send
/// contribution and receive buffer.
///
/// # Buffer Layout
///
/// `data` must have length `recvcount * size()`. Rank `r`'s contribution
/// lives at offset `r * recvcount` and must be pre-written before the call.
///
/// # Errors
///
/// - `Error::InvalidBuffer` if `data.len()` is not divisible by `size()`.
///
/// # Example
///
/// ```no_run
/// # use ferrompi::Mpi;
/// # let mpi = Mpi::init().unwrap();
/// # let world = mpi.world();
/// let rank = world.rank() as usize;
/// let size = world.size() as usize;
/// let mut data = vec![0i32; size];
/// data[rank] = rank as i32 * 10;
/// let req = world.iallgather_inplace(&mut data).unwrap();
/// req.wait().unwrap();
/// ```
pub fn iallgather_inplace<T: MpiDatatype>(&self, data: &mut [T]) -> Result<Request> {
let size = self.size() as usize;
if size == 0 || data.len() % size != 0 {
return Err(Error::InvalidBuffer);
}
let recvcount = (data.len() / size) as i64;
let mut request_handle: i64 = 0;
let ret = unsafe {
// SAFETY: data is a valid, exclusively-owned mutable slice. We cast to *mut c_void
// as required by the C FFI. Each rank's contribution (at offset rank*recvcount)
// must be pre-written by the caller. The buffer outlives the returned Request;
// the caller must call wait() before accessing or dropping the buffer.
ffi::ferrompi_iallgather_inplace(
data.as_mut_ptr().cast::<std::ffi::c_void>(),
recvcount,
T::TAG as i32,
self.handle,
&mut request_handle,
)
};
Error::check_with_op(ret, "iallgather_inplace")?;
Ok(Request::new(request_handle))
}
/// Nonblocking in-place scatter. At root, `data` is the `sendcount * size()`
/// send buffer; root's own slot is retained in place. At non-root, `data` is
/// the `recvcount`-element receive buffer.
///
/// # Buffer Layout (root)
///
/// `data` must have length `sendcount * size()`. Rank `r`'s slot is
/// `data[r*sendcount .. (r+1)*sendcount]`. After the wait, only root's own
/// slot is guaranteed to remain intact.
///
/// # Errors
///
/// - `Error::InvalidBuffer` at root if `data.len()` is not divisible by `size()`.
///
/// # Example
///
/// ```no_run
/// # use ferrompi::Mpi;
/// # let mpi = Mpi::init().unwrap();
/// # let world = mpi.world();
/// if world.rank() == 0 {
/// let mut data = vec![0i32, 10, 20, 30];
/// let req = world.iscatter_inplace(&mut data, 0).unwrap();
/// req.wait().unwrap();
/// } else {
/// let mut data = vec![0i32; 1];
/// let req = world.iscatter_inplace(&mut data, 0).unwrap();
/// req.wait().unwrap();
/// }
/// ```
pub fn iscatter_inplace<T: MpiDatatype>(&self, data: &mut [T], root: i32) -> Result<Request> {
let is_root = self.rank() == root;
let size = self.size() as usize;
let (sendbuf, sendcount, recvbuf, recvcount, is_root_flag) = if is_root {
if size == 0 || data.len() % size != 0 {
return Err(Error::InvalidBuffer);
}
let per = (data.len() / size) as i64;
(
data.as_ptr().cast::<std::ffi::c_void>(),
per,
std::ptr::null_mut::<std::ffi::c_void>(),
0i64,
1i32,
)
} else {
(
std::ptr::null::<std::ffi::c_void>(),
0i64,
data.as_mut_ptr().cast::<std::ffi::c_void>(),
data.len() as i64,
0i32,
)
};
let mut request_handle: i64 = 0;
let ret = unsafe {
// SAFETY: At root, sendbuf points to valid data of length sendcount*size elements
// (guaranteed by the divisibility check above); recvbuf is null (MPI_IN_PLACE path).
// At non-root, recvbuf points to a valid mutable slice of length recvcount elements;
// sendbuf is null (MPI standard ignores sendbuf on non-root scatter). Both pointers
// are cast to *const/*mut c_void as required by the C FFI. The buffer outlives the
// returned Request; the caller must call wait() before accessing or dropping it.
ffi::ferrompi_iscatter_inplace(
sendbuf,
sendcount,
recvbuf,
recvcount,
T::TAG as i32,
root,
is_root_flag,
self.handle,
&mut request_handle,
)
};
Error::check_with_op(ret, "iscatter_inplace")?;
Ok(Request::new(request_handle))
}
/// Nonblocking in-place all-to-all personalized communication. `data` is
/// both send and receive buffer on every rank.
///
/// Before the call, rank `r` must pre-write into slot `s` (at offset
/// `s * count`) the payload it wishes to send to rank `s`. After the wait,
/// the same slot contains the data received FROM rank `s`.
///
/// # Buffer Layout
///
/// `data` must have length `count * size()`. Slot `s` at
/// `data[s*count..(s+1)*count]` holds data sent to (and later received from)
/// rank `s`.
///
/// # Errors
///
/// - `Error::InvalidBuffer` if `data.len()` is not divisible by `size()`.
///
/// # Example
///
/// ```no_run
/// # use ferrompi::Mpi;
/// # let mpi = Mpi::init().unwrap();
/// # let world = mpi.world();
/// let r = world.rank() as i32;
/// let size = world.size() as usize;
/// let mut data: Vec<i32> = (0..size as i32).map(|s| r * 10 + s).collect();
/// let req = world.ialltoall_inplace(&mut data).unwrap();
/// req.wait().unwrap();
/// ```
pub fn ialltoall_inplace<T: MpiDatatype>(&self, data: &mut [T]) -> Result<Request> {
let size = self.size() as usize;
if size == 0 || data.len() % size != 0 {
return Err(Error::InvalidBuffer);
}
let recvcount = (data.len() / size) as i64;
let mut request_handle: i64 = 0;
let ret = unsafe {
// SAFETY: data is a valid, exclusively-owned mutable slice of length recvcount*size
// elements (guaranteed by the divisibility check above). We cast to *mut c_void as
// required by the C FFI. MPI_IN_PLACE is passed as sendbuf in the C wrapper; the
// caller must pre-write each slot before calling this method. The buffer outlives
// the returned Request; the caller must call wait() before accessing or dropping it.
ffi::ferrompi_ialltoall_inplace(
data.as_mut_ptr().cast::<std::ffi::c_void>(),
recvcount,
T::TAG as i32,
self.handle,
&mut request_handle,
)
};
Error::check_with_op(ret, "ialltoall_inplace")?;
Ok(Request::new(request_handle))
}
}
#[cfg(test)]
mod tests {
use crate::comm::Communicator;
use crate::error::Error;
use crate::ReduceOp;
fn dummy_comm() -> Communicator {
Communicator {
handle: 0,
rank: 0,
size: 1,
}
}
#[test]
fn iallreduce_mismatched_buffers_returns_invalid_buffer() {
let comm = dummy_comm();
let send = vec![1.0f64; 10];
let mut recv = vec![0.0f64; 5];
let result = comm.iallreduce(&send, &mut recv, ReduceOp::Sum);
assert!(matches!(result, Err(Error::InvalidBuffer)));
}
#[test]
fn ireduce_mismatched_buffers_returns_invalid_buffer() {
let comm = dummy_comm();
let send = vec![1.0f64; 10];
let mut recv = vec![0.0f64; 5];
let result = comm.ireduce(&send, &mut recv, ReduceOp::Sum, 0);
assert!(matches!(result, Err(Error::InvalidBuffer)));
}
#[test]
fn iscan_mismatched_buffers_returns_invalid_buffer() {
let comm = dummy_comm();
let send = vec![1.0f64; 10];
let mut recv = vec![0.0f64; 5];
let result = comm.iscan(&send, &mut recv, ReduceOp::Sum);
assert!(matches!(result, Err(Error::InvalidBuffer)));
}
#[test]
fn iexscan_mismatched_buffers_returns_invalid_buffer() {
let comm = dummy_comm();
let send = vec![1.0f64; 10];
let mut recv = vec![0.0f64; 5];
let result = comm.iexscan(&send, &mut recv, ReduceOp::Sum);
assert!(matches!(result, Err(Error::InvalidBuffer)));
}
#[test]
fn igather_inplace_nonroot_returns_invalid_op() {
let comm = Communicator {
handle: 0,
rank: 1,
size: 4,
};
let mut data = vec![0u32; 4];
let result = comm.igather_inplace(&mut data, 0);
assert!(matches!(result, Err(Error::InvalidOp)));
}
#[test]
fn iallgather_inplace_mismatched_len_returns_invalid_buffer() {
let comm = Communicator {
handle: 0,
rank: 0,
size: 4,
};
let mut data = vec![0u32; 7];
let result = comm.iallgather_inplace(&mut data);
assert!(matches!(result, Err(Error::InvalidBuffer)));
}
#[test]
fn ialltoall_inplace_mismatched_len_returns_invalid_buffer() {
let comm = Communicator {
handle: 0,
rank: 0,
size: 4,
};
let mut data = vec![0u32; 7];
let result = comm.ialltoall_inplace(&mut data);
assert!(matches!(result, Err(Error::InvalidBuffer)));
}
}