ASCOT5
Loading...
Searching...
No Matches
mpi_interface.c
Go to the documentation of this file.
1
7#ifdef MPI
8#include <mpi.h>
9#endif
10#include <stddef.h>
11#include <stdlib.h>
12#include "ascot5.h"
13#include "diag.h"
14#include "mpi_interface.h"
15#include "particle.h"
16#include "simulate.h"
17
30void mpi_interface_init(int argc, char** argv, int* mpi_rank, int* mpi_size,
31 int* mpi_root) {
32#ifdef MPI
33 /* MPI is used. Init MPI and set rank and size as determined by MPI.
34 * Rank of the root process is zero */
35 int provided;
36 int initialized;
37 MPI_Initialized(&initialized);
38 if(initialized == 0)
39 MPI_Init_thread(&argc, &argv, MPI_THREAD_FUNNELED, &provided);
40 MPI_Comm_rank(MPI_COMM_WORLD, mpi_rank);
41 MPI_Comm_size(MPI_COMM_WORLD, mpi_size);
42 *mpi_root = 0;
43#else
44 /* MPI is not used and no command line arguments are given. Assume that
45 * the entire simulation is a single process. */
46 *mpi_rank = 0;
47 *mpi_size = 1;
48 *mpi_root = 0;
49#endif
50}
51
59#ifdef MPI
60 if(err) {
61 MPI_Abort(MPI_COMM_WORLD, err);
62 }
63 else {
64 MPI_Finalize();
65 }
66#endif
67}
68
75#ifdef MPI
76 MPI_Barrier(MPI_COMM_WORLD);
77#endif
78}
79
92void mpi_my_particles(int* start_index, int* n, int n_tot, int mpi_rank,
93 int mpi_size) {
94 if(mpi_rank == mpi_size-1) {
95 *n = n_tot - mpi_rank * (n_tot / mpi_size);
96 }
97 else {
98 *n = n_tot / mpi_size;
99 }
100
101 *start_index = mpi_rank * (n_tot / mpi_size);
102}
103
123 particle_state* ps, particle_state** ps_gather, int* n_gather, int n_tot,
124 int mpi_rank, int mpi_size, int mpi_root) {
125#ifdef MPI
126 const int n_real = 32;
127 const int n_int = 5;
128 const int n_err = 1;
129
130 particle_state* ps_all = malloc(n_tot * sizeof(particle_state));
131
132 if(mpi_rank == mpi_root) {
133 int start_index, n;
134 mpi_my_particles(&start_index, &n, n_tot, mpi_rank, mpi_size);
135
136 for(int j = 0; j < n; j++) {
137 ps_all[j] = ps[j];
138 }
139
140 for(int i = 1; i < mpi_size; i++) {
141 mpi_my_particles(&start_index, &n, n_tot, i, mpi_size);
142
143 real* realdata = (real*) malloc(n_real * n * sizeof(realdata));
144 integer* intdata = (integer*) malloc(n_int * n * sizeof(intdata));
145 a5err* errdata = (a5err*) malloc(n_err * n * sizeof(intdata));
146
147 MPI_Recv(realdata, n_real*n, mpi_type_real, i, 0,
148 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
149 MPI_Recv(intdata, n_int*n, mpi_type_integer, i, 0,
150 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
151 MPI_Recv(errdata, n_err*n, mpi_type_a5err, i, 0,
152 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
153
154 for(int j = 0; j < n; j++) {
155 ps_all[start_index+j].r = realdata[0*n+j];
156 ps_all[start_index+j].phi = realdata[1*n+j];
157 ps_all[start_index+j].z = realdata[2*n+j];
158 ps_all[start_index+j].ppar = realdata[3*n+j];
159 ps_all[start_index+j].mu = realdata[4*n+j];
160 ps_all[start_index+j].zeta = realdata[5*n+j];
161 ps_all[start_index+j].rprt = realdata[6*n+j];
162 ps_all[start_index+j].phiprt = realdata[7*n+j];
163 ps_all[start_index+j].zprt = realdata[8*n+j];
164 ps_all[start_index+j].p_r = realdata[9*n+j];
165 ps_all[start_index+j].p_phi = realdata[10*n+j];
166 ps_all[start_index+j].p_z = realdata[11*n+j];
167 ps_all[start_index+j].mass = realdata[12*n+j];
168 ps_all[start_index+j].charge = realdata[13*n+j];
169 ps_all[start_index+j].anum = intdata[0*n+j];
170 ps_all[start_index+j].znum = intdata[1*n+j];
171 ps_all[start_index+j].weight = realdata[14*n+j];
172 ps_all[start_index+j].time = realdata[15*n+j];
173 ps_all[start_index+j].cputime = realdata[16*n+j];
174 ps_all[start_index+j].rho = realdata[17*n+j];
175 ps_all[start_index+j].theta = realdata[18*n+j];
176 ps_all[start_index+j].id = intdata[2*n+j];
177 ps_all[start_index+j].endcond = intdata[3*n+j];
178 ps_all[start_index+j].walltile = intdata[4*n+j];
179 ps_all[start_index+j].B_r = realdata[19*n+j];
180 ps_all[start_index+j].B_phi = realdata[20*n+j];
181 ps_all[start_index+j].B_z = realdata[21*n+j];
182 ps_all[start_index+j].B_r_dr = realdata[22*n+j];
183 ps_all[start_index+j].B_phi_dr = realdata[23*n+j];
184 ps_all[start_index+j].B_z_dr = realdata[24*n+j];
185 ps_all[start_index+j].B_r_dphi = realdata[25*n+j];
186 ps_all[start_index+j].B_phi_dphi = realdata[26*n+j];
187 ps_all[start_index+j].B_z_dphi = realdata[27*n+j];
188 ps_all[start_index+j].B_r_dz = realdata[28*n+j];
189 ps_all[start_index+j].B_phi_dz = realdata[29*n+j];
190 ps_all[start_index+j].B_z_dz = realdata[30*n+j];
191 ps_all[start_index+j].mileage = realdata[31*n+j];
192 ps_all[start_index+j].err = errdata[j];
193 }
194
195 free(realdata);
196 free(intdata);
197 free(errdata);
198 }
199 }
200 else {
201
202 int start_index, n;
203 mpi_my_particles(&start_index, &n, n_tot, mpi_rank, mpi_size);
204
205 real* realdata;
206 realdata = malloc(n_real * n * sizeof(realdata));
207 integer* intdata;
208 intdata = malloc(n_int * n * sizeof(intdata));
209 a5err* errdata;
210 errdata = malloc(n_err * n * sizeof(intdata));
211
212 for(int j = 0; j < n; j++) {
213 realdata[0*n+j] = ps[j].r;
214 realdata[1*n+j] = ps[j].phi;
215 realdata[2*n+j] = ps[j].z;
216 realdata[3*n+j] = ps[j].ppar;
217 realdata[4*n+j] = ps[j].mu;
218 realdata[5*n+j] = ps[j].zeta;
219 realdata[6*n+j] = ps[j].rprt;
220 realdata[7*n+j] = ps[j].phiprt;
221 realdata[8*n+j] = ps[j].zprt;
222 realdata[9*n+j] = ps[j].p_r;
223 realdata[10*n+j] = ps[j].p_phi;
224 realdata[11*n+j] = ps[j].p_z;
225 realdata[12*n+j] = ps[j].mass;
226 realdata[13*n+j] = ps[j].charge;
227 intdata[0*n+j] = ps[j].anum;
228 intdata[1*n+j] = ps[j].znum;
229 realdata[14*n+j] = ps[j].weight;
230 realdata[15*n+j] = ps[j].time;
231 realdata[16*n+j] = ps[j].cputime;
232 realdata[17*n+j] = ps[j].rho;
233 realdata[18*n+j] = ps[j].theta;
234 intdata[2*n+j] = ps[j].id;
235 intdata[3*n+j] = ps[j].endcond;
236 intdata[4*n+j] = ps[j].walltile;
237 realdata[19*n+j] = ps[j].B_r;
238 realdata[20*n+j] = ps[j].B_phi;
239 realdata[21*n+j] = ps[j].B_z;
240 realdata[22*n+j] = ps[j].B_r_dr;
241 realdata[23*n+j] = ps[j].B_phi_dr;
242 realdata[24*n+j] = ps[j].B_z_dr;
243 realdata[25*n+j] = ps[j].B_r_dphi;
244 realdata[26*n+j] = ps[j].B_phi_dphi;
245 realdata[27*n+j] = ps[j].B_z_dphi;
246 realdata[28*n+j] = ps[j].B_r_dz;
247 realdata[29*n+j] = ps[j].B_phi_dz;
248 realdata[30*n+j] = ps[j].B_z_dz;
249 realdata[31*n+j] = ps[j].mileage;
250 errdata[j] = ps[j].err;
251 }
252
253 MPI_Send(realdata, n_real*n, mpi_type_real, 0, 0, MPI_COMM_WORLD);
254 MPI_Send(intdata, n_int*n, mpi_type_integer, 0, 0, MPI_COMM_WORLD);
255 MPI_Send(errdata, n_err*n, mpi_type_a5err, 0, 0, MPI_COMM_WORLD);
256
257 free(realdata);
258 free(intdata);
259 free(errdata);
260 }
261
262 *ps_gather = ps_all;
263 *n_gather = n_tot;
264
265#else
266
267 int start_index, n;
268 mpi_my_particles(&start_index, &n, n_tot, mpi_rank, mpi_size);
269
270 /* Only store particles for this process in Condor-style execution */
271 particle_state* ps_all = malloc(n * sizeof(particle_state));
272
273 for(int j = 0; j < n; j++) {
274 ps_all[j] = ps[j];
275 }
276
277 *ps_gather = ps_all;
278 *n_gather = n;
279
280#endif
281}
282
296void mpi_gather_diag(diag_data* data, int ntotal, int mpi_rank, int mpi_size,
297 int mpi_root) {
298#ifdef MPI
299
300 if(data->dist5D_collect) {
301 MPI_Reduce(
302 mpi_rank == mpi_root ? MPI_IN_PLACE : data->dist5D.histogram,
303 data->dist5D.histogram,
304 data->dist5D.step_6 * data->dist5D.n_r,
305 mpi_type_real, MPI_SUM, mpi_root, MPI_COMM_WORLD);
306 }
307 if(data->dist6D_collect) {
308 MPI_Reduce(
309 mpi_rank == mpi_root ? MPI_IN_PLACE : data->dist6D.histogram,
310 mpi_rank == mpi_root ? NULL : data->dist6D.histogram,
311 data->dist6D.step_7 * data->dist6D.n_r,
312 mpi_type_real, MPI_SUM, mpi_root, MPI_COMM_WORLD);
313 }
314 if(data->distrho5D_collect) {
315 MPI_Reduce(
316 mpi_rank == mpi_root ? MPI_IN_PLACE : data->distrho5D.histogram,
317 mpi_rank == mpi_root ? NULL : data->distrho5D.histogram,
318 data->distrho5D.step_6 * data->distrho5D.n_rho,
319 mpi_type_real, MPI_SUM, mpi_root, MPI_COMM_WORLD);
320 }
321 if(data->distrho6D_collect) {
322 MPI_Reduce(
323 mpi_rank == mpi_root ? MPI_IN_PLACE : data->distrho6D.histogram,
324 mpi_rank == mpi_root ? NULL : data->distrho6D.histogram,
325 data->distrho6D.step_7*data->distrho6D.n_rho,
326 mpi_type_real, MPI_SUM, mpi_root, MPI_COMM_WORLD);
327 }
328 if(data->distCOM_collect) {
329 MPI_Reduce(
330 mpi_rank == mpi_root ? MPI_IN_PLACE : data->distCOM.histogram,
331 mpi_rank == mpi_root ? NULL : data->distCOM.histogram,
332 data->distCOM.step_2 * data->distCOM.n_mu,
333 mpi_type_real, MPI_SUM, mpi_root, MPI_COMM_WORLD);
334 }
335
336 if(data->diagorb_collect) {
337 size_t nreal = ntotal*data->diagorb.Npnt * sizeof(real);
338 if(mpi_rank == mpi_root) {
339 data->diagorb.id = realloc(data->diagorb.id, nreal);
340 data->diagorb.mileage = realloc(data->diagorb.mileage, nreal);
341 data->diagorb.r = realloc(data->diagorb.r, nreal);
342 data->diagorb.phi = realloc(data->diagorb.phi, nreal);
343 data->diagorb.z = realloc(data->diagorb.z, nreal);
344 if(data->diagorb.record_mode == DIAG_ORB_FO) {
345 data->diagorb.p_r = realloc(data->diagorb.p_r, nreal);
346 data->diagorb.p_phi = realloc(data->diagorb.p_phi, nreal);
347 data->diagorb.p_z = realloc(data->diagorb.p_z, nreal);
348 }
349 if(data->diagorb.record_mode == DIAG_ORB_GC) {
350 data->diagorb.ppar = realloc(data->diagorb.ppar, nreal);
351 data->diagorb.mu = realloc(data->diagorb.mu, nreal);
352 data->diagorb.zeta = realloc(data->diagorb.zeta, nreal);
353 }
354 if(data->diagorb.record_mode != DIAG_ORB_ML) {
355 data->diagorb.charge = realloc(data->diagorb.charge, nreal);
356 }
357 data->diagorb.weight = realloc(data->diagorb.weight, nreal);
358 data->diagorb.rho = realloc(data->diagorb.rho, nreal);
359 data->diagorb.theta = realloc(data->diagorb.theta, nreal);
360 data->diagorb.B_r = realloc(data->diagorb.B_r, nreal);
361 data->diagorb.B_phi = realloc(data->diagorb.B_phi, nreal);
362 data->diagorb.B_z = realloc(data->diagorb.B_z, nreal);
363 data->diagorb.simmode = realloc(data->diagorb.simmode, nreal);
364 if(data->diagorb.mode == DIAG_ORB_POINCARE) {
365 data->diagorb.pncrid = realloc(data->diagorb.pncrid, nreal);
366 data->diagorb.pncrdi = realloc(data->diagorb.pncrdi, nreal);
367 }
368 for(int i = 1; i < mpi_size; i++) {
369 int start_index, n;
370 mpi_my_particles(&start_index, &n, ntotal, i, mpi_size);
371 MPI_Recv(
372 &data->diagorb.id[start_index*data->diagorb.Npnt],
373 n*data->diagorb.Npnt, mpi_type_real, i, 0,
374 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
375 MPI_Recv(
376 &data->diagorb.mileage[start_index*data->diagorb.Npnt],
377 n*data->diagorb.Npnt, mpi_type_real, i, 0,
378 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
379 MPI_Recv(
380 &data->diagorb.r[start_index*data->diagorb.Npnt],
381 n*data->diagorb.Npnt, mpi_type_real, i, 0,
382 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
383 MPI_Recv(
384 &data->diagorb.phi[start_index*data->diagorb.Npnt],
385 n*data->diagorb.Npnt, mpi_type_real, i, 0,
386 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
387 MPI_Recv(
388 &data->diagorb.z[start_index*data->diagorb.Npnt],
389 n*data->diagorb.Npnt, mpi_type_real, i, 0,
390 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
391 if(data->diagorb.record_mode == DIAG_ORB_FO) {
392 MPI_Recv(
393 &data->diagorb.p_r[start_index*data->diagorb.Npnt],
394 n*data->diagorb.Npnt, mpi_type_real, i, 0,
395 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
396 MPI_Recv(
397 &data->diagorb.p_phi[start_index*data->diagorb.Npnt],
398 n*data->diagorb.Npnt, mpi_type_real, i, 0,
399 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
400 MPI_Recv(
401 &data->diagorb.p_z[start_index*data->diagorb.Npnt],
402 n*data->diagorb.Npnt, mpi_type_real, i, 0,
403 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
404 }
405 if(data->diagorb.record_mode == DIAG_ORB_GC) {
406 MPI_Recv(
407 &data->diagorb.ppar[start_index*data->diagorb.Npnt],
408 n*data->diagorb.Npnt, mpi_type_real, i, 0,
409 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
410 MPI_Recv(
411 &data->diagorb.mu[start_index*data->diagorb.Npnt],
412 n*data->diagorb.Npnt, mpi_type_real, i, 0,
413 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
414 MPI_Recv(
415 &data->diagorb.zeta[start_index*data->diagorb.Npnt],
416 n*data->diagorb.Npnt, mpi_type_real, i, 0,
417 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
418 }
419 if(data->diagorb.record_mode != DIAG_ORB_ML) {
420 MPI_Recv(
421 &data->diagorb.charge[start_index*data->diagorb.Npnt],
422 n*data->diagorb.Npnt, mpi_type_real, i, 0,
423 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
424 }
425 MPI_Recv(
426 &data->diagorb.weight[start_index*data->diagorb.Npnt],
427 n*data->diagorb.Npnt, mpi_type_real, i, 0,
428 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
429 MPI_Recv(
430 &data->diagorb.rho[start_index*data->diagorb.Npnt],
431 n*data->diagorb.Npnt, mpi_type_real, i, 0,
432 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
433 MPI_Recv(
434 &data->diagorb.theta[start_index*data->diagorb.Npnt],
435 n*data->diagorb.Npnt, mpi_type_real, i, 0,
436 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
437 MPI_Recv(
438 &data->diagorb.B_r[start_index*data->diagorb.Npnt],
439 n*data->diagorb.Npnt, mpi_type_real, i, 0,
440 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
441 MPI_Recv(
442 &data->diagorb.B_phi[start_index*data->diagorb.Npnt],
443 n*data->diagorb.Npnt, mpi_type_real, i, 0,
444 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
445 MPI_Recv(
446 &data->diagorb.B_z[start_index*data->diagorb.Npnt],
447 n*data->diagorb.Npnt, mpi_type_real, i, 0,
448 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
449 MPI_Recv(
450 &data->diagorb.simmode[start_index*data->diagorb.Npnt],
451 n*data->diagorb.Npnt, mpi_type_real, i, 0,
452 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
453 if(data->diagorb.mode == DIAG_ORB_POINCARE) {
454 MPI_Recv(
455 &data->diagorb.pncrid[start_index*data->diagorb.Npnt],
456 n*data->diagorb.Npnt, mpi_type_real, i, 0,
457 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
458 MPI_Recv(
459 &data->diagorb.pncrdi[start_index*data->diagorb.Npnt],
460 n*data->diagorb.Npnt, mpi_type_real, i, 0,
461 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
462 }
463 }
464 }
465 else {
466 int start_index, n;
467 mpi_my_particles(&start_index, &n, ntotal, mpi_rank, mpi_size);
468 MPI_Send(data->diagorb.id, n*data->diagorb.Npnt,
469 mpi_type_real, 0, 0, MPI_COMM_WORLD);
470 MPI_Send(data->diagorb.mileage, n*data->diagorb.Npnt,
471 mpi_type_real, 0, 0, MPI_COMM_WORLD);
472 MPI_Send(data->diagorb.r, n*data->diagorb.Npnt,
473 mpi_type_real, 0, 0, MPI_COMM_WORLD);
474 MPI_Send(data->diagorb.phi, n*data->diagorb.Npnt,
475 mpi_type_real, 0, 0, MPI_COMM_WORLD);
476 MPI_Send(data->diagorb.z, n*data->diagorb.Npnt,
477 mpi_type_real, 0, 0, MPI_COMM_WORLD);
478 if(data->diagorb.record_mode == DIAG_ORB_FO) {
479 MPI_Send(data->diagorb.p_r, n*data->diagorb.Npnt,
480 mpi_type_real, 0, 0, MPI_COMM_WORLD);
481 MPI_Send(data->diagorb.p_phi, n*data->diagorb.Npnt,
482 mpi_type_real, 0, 0, MPI_COMM_WORLD);
483 MPI_Send(data->diagorb.p_z, n*data->diagorb.Npnt,
484 mpi_type_real, 0, 0, MPI_COMM_WORLD);
485 }
486 if(data->diagorb.record_mode == DIAG_ORB_GC) {
487 MPI_Send(data->diagorb.ppar, n*data->diagorb.Npnt,
488 mpi_type_real, 0, 0, MPI_COMM_WORLD);
489 MPI_Send(data->diagorb.mu, n*data->diagorb.Npnt,
490 mpi_type_real, 0, 0, MPI_COMM_WORLD);
491 MPI_Send(data->diagorb.zeta, n*data->diagorb.Npnt,
492 mpi_type_real, 0, 0, MPI_COMM_WORLD);
493 }
494 if(data->diagorb.record_mode != DIAG_ORB_ML) {
495 MPI_Send(data->diagorb.charge, n*data->diagorb.Npnt,
496 mpi_type_real, 0, 0, MPI_COMM_WORLD);
497 }
498 MPI_Send(data->diagorb.weight, n*data->diagorb.Npnt,
499 mpi_type_real, 0, 0, MPI_COMM_WORLD);
500 MPI_Send(data->diagorb.rho, n*data->diagorb.Npnt,
501 mpi_type_real, 0, 0, MPI_COMM_WORLD);
502 MPI_Send(data->diagorb.theta, n*data->diagorb.Npnt,
503 mpi_type_real, 0, 0, MPI_COMM_WORLD);
504 MPI_Send(data->diagorb.B_r, n*data->diagorb.Npnt,
505 mpi_type_real, 0, 0, MPI_COMM_WORLD);
506 MPI_Send(data->diagorb.B_phi, n*data->diagorb.Npnt,
507 mpi_type_real, 0, 0, MPI_COMM_WORLD);
508 MPI_Send(data->diagorb.B_z, n*data->diagorb.Npnt,
509 mpi_type_real, 0, 0, MPI_COMM_WORLD);
510 MPI_Send(data->diagorb.simmode, n*data->diagorb.Npnt,
511 mpi_type_real, 0, 0, MPI_COMM_WORLD);
512 if(data->diagorb.mode == DIAG_ORB_POINCARE) {
513 MPI_Send(data->diagorb.pncrid, n*data->diagorb.Npnt,
514 mpi_type_real, 0, 0, MPI_COMM_WORLD);
515 MPI_Send(data->diagorb.pncrdi, n*data->diagorb.Npnt,
516 mpi_type_real, 0, 0, MPI_COMM_WORLD);
517 }
518 }
519 }
520
521 if(data->diagtrcof_collect) {
522
523 if(mpi_rank == mpi_root) {
524 data->diagtrcof.id = realloc(
525 data->diagtrcof.id, ntotal*data->diagtrcof.Nmrk*sizeof(int)
526 );
527 data->diagtrcof.Kcoef = realloc(
528 data->diagtrcof.Kcoef, ntotal*data->diagtrcof.Nmrk*sizeof(real)
529 );
530 data->diagtrcof.Dcoef = realloc(
531 data->diagtrcof.Dcoef,
532 ntotal*data->diagtrcof.Nmrk*sizeof(real)
533 );
534 for(int i = 1; i < mpi_size; i++) {
535 int start_index, n;
536 mpi_my_particles(&start_index, &n, ntotal, i, mpi_size);
537
538 MPI_Recv(&data->diagtrcof.id[start_index], n, mpi_type_integer,
539 i, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
540 MPI_Recv(&data->diagtrcof.Kcoef[start_index], n, mpi_type_real,
541 i, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
542 MPI_Recv(&data->diagtrcof.Dcoef[start_index], n, mpi_type_real,
543 i, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
544 }
545 }
546 else {
547 int start_index, n;
548 mpi_my_particles(&start_index, &n, ntotal, mpi_rank, mpi_size);
549 MPI_Send(data->diagtrcof.id, n, mpi_type_integer,
550 0, 0, MPI_COMM_WORLD);
551 MPI_Send(data->diagtrcof.Kcoef, n, mpi_type_real,
552 0, 0, MPI_COMM_WORLD);
553 MPI_Send(data->diagtrcof.Dcoef, n, mpi_type_real,
554 0, 0, MPI_COMM_WORLD);
555 }
556 }
557
558#endif
559}
Main header file for ASCOT5.
double real
Definition ascot5.h:85
long integer
Definition ascot5.h:84
Header file for diag.c.
#define DIAG_ORB_POINCARE
Definition diag_orb.h:13
#define DIAG_ORB_ML
Definition diag_orb.h:24
#define DIAG_ORB_GC
Definition diag_orb.h:23
#define DIAG_ORB_FO
Definition diag_orb.h:22
unsigned long int a5err
Simulation error flag.
Definition error.h:17
void mpi_gather_diag(diag_data *data, int ntotal, int mpi_rank, int mpi_size, int mpi_root)
Gather all diagnostics to the root process.
void mpi_my_particles(int *start_index, int *n, int n_tot, int mpi_rank, int mpi_size)
Divide markers to mpi processes.
void mpi_interface_init(int argc, char **argv, int *mpi_rank, int *mpi_size, int *mpi_root)
Initialize MPI.
void mpi_interface_finalize(int err)
Finalize MPI.
void mpi_interface_barrier()
Wait until all processes have reached this routine.
void mpi_gather_particlestate(particle_state *ps, particle_state **ps_gather, int *n_gather, int n_tot, int mpi_rank, int mpi_size, int mpi_root)
Gather all particle states to the root process.
Header file for mpi_interface.c.
#define mpi_type_integer
ASCOT integer in MPI standard.
#define mpi_type_a5err
ASCOT error in MPI standard.
#define mpi_type_real
ASCOT real in MPI standard.
Header file for particle.c.
Header file for simulate.c.
Diagnostics data struct.
Definition diag.h:21
dist_COM_data distCOM
Definition diag.h:35
int distrho5D_collect
Definition diag.h:25
int diagtrcof_collect
Definition diag.h:28
diag_transcoef_data diagtrcof
Definition diag.h:36
int distrho6D_collect
Definition diag.h:26
diag_orb_data diagorb
Definition diag.h:30
dist_rho6D_data distrho6D
Definition diag.h:34
dist_5D_data dist5D
Definition diag.h:31
int dist5D_collect
Definition diag.h:23
int dist6D_collect
Definition diag.h:24
int distCOM_collect
Definition diag.h:27
dist_6D_data dist6D
Definition diag.h:32
dist_rho5D_data distrho5D
Definition diag.h:33
int diagorb_collect
Definition diag.h:22
real * mileage
Definition diag_orb.h:44
real * ppar
Definition diag_orb.h:51
real * simmode
Definition diag_orb.h:61
real * B_r
Definition diag_orb.h:58
real * pncrdi
Definition diag_orb.h:63
real * zeta
Definition diag_orb.h:53
real * p_phi
Definition diag_orb.h:49
real * theta
Definition diag_orb.h:57
real * rho
Definition diag_orb.h:56
real * pncrid
Definition diag_orb.h:62
real * weight
Definition diag_orb.h:54
real * B_z
Definition diag_orb.h:60
real * phi
Definition diag_orb.h:46
real * mu
Definition diag_orb.h:52
real * p_z
Definition diag_orb.h:50
int record_mode
Definition diag_orb.h:76
real * charge
Definition diag_orb.h:55
real * id
Definition diag_orb.h:43
real * p_r
Definition diag_orb.h:48
real * B_phi
Definition diag_orb.h:59
size_t step_6
Definition dist_5D.h:49
real * histogram
Definition dist_5D.h:51
size_t step_7
Definition dist_6D.h:54
real * histogram
Definition dist_6D.h:56
real * histogram
Definition dist_com.h:32
size_t step_2
Definition dist_com.h:30
General representation of a marker.
Definition particle.h:40
integer id
Definition particle.h:63
integer walltile
Definition particle.h:65
real B_phi_dphi
Definition particle.h:73
integer endcond
Definition particle.h:64