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
//! MPI-specific implementation for distributed integration.
use crate::integrand::{Integrand, SimdIntegrand};
use crate::vegas::VegasResult;
use crate::vegasplus::VegasPlus;
use mpi::collective::SystemOperation;
use mpi::topology::Communicator;
use mpi::traits::*;
impl VegasPlus {
/// Integrates the given function in a distributed manner using MPI.
///
/// This method should be called by all processes in the MPI communicator.
/// The final result is returned on the root process (rank 0), while other
/// processes will return a default `VegasResult`.
///
/// # Arguments
///
/// * `integrand`: The function to integrate. Must be `Sync`.
/// * `world`: The MPI communicator.
///
/// # Examples
///
/// ```no_run
/// use mchep::vegasplus::VegasPlus;
/// use mchep::integrand::Integrand;
/// use mpi::traits::*;
///
/// struct MyIntegrand;
///
/// impl Integrand for MyIntegrand {
/// fn dim(&self) -> usize { 2 }
/// fn eval(&self, x: &[f64]) -> f64 {
/// (-(x[0].powi(2)) - x[1].powi(2)).exp()
/// }
/// }
///
/// // Initialize MPI
/// let universe = mpi::initialize().unwrap();
/// let world = universe.world();
/// let rank = world.rank();
///
/// let integrand = MyIntegrand;
/// let boundaries = &[(-1.0, 1.0), (-1.0, 1.0)];
/// let mut vegas_plus = VegasPlus::new(10, 20_000, 50, 0.5, 4, 0.75, boundaries);
/// vegas_plus.set_seed(1234);
///
/// let result = vegas_plus.integrate_mpi(&integrand, &world, None);
///
/// if rank == 0 {
/// println!("Result: {} +/- {}", result.value, result.error);
/// assert!((result.value - 2.230985).abs() < 3. * result.error);
/// }
/// ```
pub fn integrate_mpi<F: Integrand + Sync, C: Communicator>(
&mut self,
integrand: &F,
world: &C,
target_accuracy: Option<f64>,
) -> VegasResult {
assert_eq!(integrand.dim(), self.dim());
let rank = world.rank();
let size = world.size();
let n_eval_per_rank = self.n_eval / size as usize;
let n_eval_orig = self.n_eval;
self.n_eval = n_eval_per_rank;
let mut iter_results = Vec::new();
let mut iter_errors = Vec::new();
for iter in 0..self.n_iter {
let (iter_val, iter_err) = self.run_iteration(integrand);
let local_data = self.serialize_adaptation_data();
let mut global_data = vec![0.0; local_data.len()];
world.all_reduce_into(
&local_data[..],
&mut global_data[..],
SystemOperation::sum(),
);
self.deserialize_and_update_state(&global_data);
let mut stop_signal = 0i32;
if rank == 0 {
let mut global_iter_val = 0.0;
world.process_at_rank(0).reduce_into_root(&iter_val, &mut global_iter_val, SystemOperation::sum());
let mut global_iter_err_sq: f64 = 0.0;
let iter_err_sq = iter_err.powi(2);
world.process_at_rank(0).reduce_into_root(
&iter_err_sq,
&mut global_iter_err_sq,
SystemOperation::sum(),
);
if iter > 0 {
iter_results.push(global_iter_val);
iter_errors.push(global_iter_err_sq.sqrt());
}
if let Some(acc_req) = target_accuracy {
if !iter_results.is_empty() {
let current_result = self.combine_results(&iter_results, &iter_errors);
if current_result.value != 0.0 {
let current_acc =
(current_result.error / current_result.value.abs()) * 100.0;
if current_acc < acc_req {
stop_signal = 1;
}
}
}
}
} else {
world.process_at_rank(0).reduce_into(&iter_val, SystemOperation::sum());
let iter_err_sq = iter_err.powi(2);
world.process_at_rank(0).reduce_into(&iter_err_sq, SystemOperation::sum());
}
world.process_at_rank(0).broadcast_into(&mut stop_signal);
if stop_signal == 1 {
break;
}
self.reallocate_samples();
for grid in &mut self.grids {
grid.refine();
}
}
self.n_eval = n_eval_orig;
if rank == 0 {
self.combine_results(&iter_results, &iter_errors)
} else {
VegasResult {
value: 0.0,
error: 0.0,
chi2_dof: 0.0,
}
}
}
/// Integrates the given SIMD function in a distributed manner using MPI.
///
/// This method combines MPI distribution across processes with SIMD vectorization
/// within each process, providing multiplicative speedup for heavy integrands.
///
/// This method should be called by all processes in the MPI communicator.
/// The final result is returned on the root process (rank 0), while other
/// processes will return a default `VegasResult`.
///
/// # Arguments
///
/// * `integrand`: The SIMD function to integrate. Must implement `SimdIntegrand + Sync`.
/// * `world`: The MPI communicator.
/// * `target_accuracy`: Optional target accuracy in percent. If provided, integration
/// stops early when this accuracy is reached.
///
/// # Examples
///
/// ```no_run
/// use mchep::vegasplus::VegasPlus;
/// use mchep::integrand::SimdIntegrand;
/// use mpi::traits::*;
/// use wide::f64x4;
///
/// struct MySimdIntegrand;
///
/// impl SimdIntegrand for MySimdIntegrand {
/// fn dim(&self) -> usize { 2 }
/// fn eval_simd(&self, points: &[f64x4]) -> f64x4 {
/// let x = points[0];
/// let y = points[1];
/// (-(x * x) - (y * y)).exp()
/// }
/// }
///
/// // Initialize MPI
/// let universe = mpi::initialize().unwrap();
/// let world = universe.world();
/// let rank = world.rank();
///
/// let integrand = MySimdIntegrand;
/// let boundaries = &[(-1.0, 1.0), (-1.0, 1.0)];
/// let mut vegas_plus = VegasPlus::new(10, 20_000, 50, 0.5, 4, 0.75, boundaries);
/// vegas_plus.set_seed(1234);
///
/// let result = vegas_plus.integrate_mpi_simd(&integrand, &world, None);
///
/// if rank == 0 {
/// println!("Result: {} +/- {}", result.value, result.error);
/// assert!((result.value - 2.230985).abs() < 3. * result.error);
/// }
/// ```
pub fn integrate_mpi_simd<F: SimdIntegrand + Sync, C: Communicator>(
&mut self,
integrand: &F,
world: &C,
target_accuracy: Option<f64>,
) -> VegasResult {
assert_eq!(integrand.dim(), self.dim());
let rank = world.rank();
let size = world.size();
let n_eval_per_rank = self.n_eval / size as usize;
let n_eval_orig = self.n_eval;
self.n_eval = n_eval_per_rank;
let mut iter_results = Vec::new();
let mut iter_errors = Vec::new();
for iter in 0..self.n_iter {
let (iter_val, iter_err) = self.run_iteration_simd(integrand);
let local_data = self.serialize_adaptation_data();
let mut global_data = vec![0.0; local_data.len()];
world.all_reduce_into(
&local_data[..],
&mut global_data[..],
SystemOperation::sum(),
);
self.deserialize_and_update_state(&global_data);
let mut stop_signal = 0i32;
if rank == 0 {
let mut global_iter_val = 0.0;
world.process_at_rank(0).reduce_into_root(&iter_val, &mut global_iter_val, SystemOperation::sum());
let mut global_iter_err_sq: f64 = 0.0;
let iter_err_sq = iter_err.powi(2);
world.process_at_rank(0).reduce_into_root(
&iter_err_sq,
&mut global_iter_err_sq,
SystemOperation::sum(),
);
if iter > 0 {
iter_results.push(global_iter_val);
iter_errors.push(global_iter_err_sq.sqrt());
}
if let Some(acc_req) = target_accuracy {
if !iter_results.is_empty() {
let current_result = self.combine_results(&iter_results, &iter_errors);
if current_result.value != 0.0 {
let current_acc =
(current_result.error / current_result.value.abs()) * 100.0;
if current_acc < acc_req {
stop_signal = 1;
}
}
}
}
} else {
world.process_at_rank(0).reduce_into(&iter_val, SystemOperation::sum());
let iter_err_sq = iter_err.powi(2);
world.process_at_rank(0).reduce_into(&iter_err_sq, SystemOperation::sum());
}
world.process_at_rank(0).broadcast_into(&mut stop_signal);
if stop_signal == 1 {
break;
}
self.reallocate_samples();
for grid in &mut self.grids {
grid.refine();
}
}
self.n_eval = n_eval_orig;
if rank == 0 {
self.combine_results(&iter_results, &iter_errors)
} else {
VegasResult {
value: 0.0,
error: 0.0,
chi2_dof: 0.0,
}
}
}
/// Serializes the grid `d` vectors and hypercube variances into a flat Vec.
fn serialize_adaptation_data(&self) -> Vec<f64> {
let n_grids = self.grids.len();
let n_bins = self.grids[0].n_bins();
let n_hypercubes = self.hypercubes.len();
let mut data = vec![0.0; n_grids * n_bins + n_hypercubes];
for d in 0..n_grids {
data[d * n_bins..(d + 1) * n_bins].copy_from_slice(&self.grids[d].d);
}
let offset = n_grids * n_bins;
for i in 0..n_hypercubes {
data[offset + i] = self.hypercubes[i].variance;
}
data
}
/// Updates the state from a flat Vec of globally reduced adaptation data.
fn deserialize_and_update_state(&mut self, global_data: &[f64]) {
let n_grids = self.grids.len();
let n_bins = self.grids[0].n_bins();
let n_hypercubes = self.hypercubes.len();
for d in 0..n_grids {
self.grids[d]
.d
.copy_from_slice(&global_data[d * n_bins..(d + 1) * n_bins]);
}
let offset = n_grids * n_bins;
for i in 0..n_hypercubes {
self.hypercubes[i].variance = global_data[offset + i];
}
}
}