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 MPI_Finalize();
61#endif
62}
63
70#ifdef MPI
71 MPI_Barrier(MPI_COMM_WORLD);
72#endif
73}
74
87void mpi_my_particles(int* start_index, int* n, int n_tot, int mpi_rank,
88 int mpi_size) {
89 if(mpi_rank == mpi_size-1) {
90 *n = n_tot - mpi_rank * (n_tot / mpi_size);
91 }
92 else {
93 *n = n_tot / mpi_size;
94 }
95
96 *start_index = mpi_rank * (n_tot / mpi_size);
97}
98
118 particle_state* ps, particle_state** ps_gather, int* n_gather, int n_tot,
119 int mpi_rank, int mpi_size, int mpi_root) {
120#ifdef MPI
121 const int n_real = 32;
122 const int n_int = 5;
123 const int n_err = 1;
124
125 particle_state* ps_all = malloc(n_tot * sizeof(particle_state));
126
127 if(mpi_rank == mpi_root) {
128 int start_index, n;
129 mpi_my_particles(&start_index, &n, n_tot, mpi_rank, mpi_size);
130
131 for(int j = 0; j < n; j++) {
132 ps_all[j] = ps[j];
133 }
134
135 for(int i = 1; i < mpi_size; i++) {
136 mpi_my_particles(&start_index, &n, n_tot, i, mpi_size);
137
138 real* realdata;
139 realdata = malloc(n_real * n * sizeof(realdata));
140 integer* intdata;
141 intdata = malloc(n_int * n * sizeof(intdata));
142 a5err* errdata;
143 errdata = malloc(n_err * n * sizeof(intdata));
144
145 MPI_Recv(realdata, n_real*n, mpi_type_real, i, 0,
146 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
147 MPI_Recv(intdata, n_int*n, mpi_type_integer, i, 0,
148 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
149 MPI_Recv(errdata, n_err*n, mpi_type_a5err, i, 0,
150 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
151
152 for(int j = 0; j < n; j++) {
153 ps_all[start_index+j].r = realdata[0*n+j];
154 ps_all[start_index+j].phi = realdata[1*n+j];
155 ps_all[start_index+j].z = realdata[2*n+j];
156 ps_all[start_index+j].ppar = realdata[3*n+j];
157 ps_all[start_index+j].mu = realdata[4*n+j];
158 ps_all[start_index+j].zeta = realdata[5*n+j];
159 ps_all[start_index+j].rprt = realdata[6*n+j];
160 ps_all[start_index+j].phiprt = realdata[7*n+j];
161 ps_all[start_index+j].zprt = realdata[8*n+j];
162 ps_all[start_index+j].p_r = realdata[9*n+j];
163 ps_all[start_index+j].p_phi = realdata[10*n+j];
164 ps_all[start_index+j].p_z = realdata[11*n+j];
165 ps_all[start_index+j].mass = realdata[12*n+j];
166 ps_all[start_index+j].charge = realdata[13*n+j];
167 ps_all[start_index+j].anum = intdata[0*n+j];
168 ps_all[start_index+j].znum = intdata[1*n+j];
169 ps_all[start_index+j].weight = realdata[14*n+j];
170 ps_all[start_index+j].time = realdata[15*n+j];
171 ps_all[start_index+j].cputime = realdata[16*n+j];
172 ps_all[start_index+j].rho = realdata[17*n+j];
173 ps_all[start_index+j].theta = realdata[18*n+j];
174 ps_all[start_index+j].id = intdata[2*n+j];
175 ps_all[start_index+j].endcond = intdata[3*n+j];
176 ps_all[start_index+j].walltile = intdata[4*n+j];
177 ps_all[start_index+j].B_r = realdata[19*n+j];
178 ps_all[start_index+j].B_phi = realdata[20*n+j];
179 ps_all[start_index+j].B_z = realdata[21*n+j];
180 ps_all[start_index+j].B_r_dr = realdata[22*n+j];
181 ps_all[start_index+j].B_phi_dr = realdata[23*n+j];
182 ps_all[start_index+j].B_z_dr = realdata[24*n+j];
183 ps_all[start_index+j].B_r_dphi = realdata[25*n+j];
184 ps_all[start_index+j].B_phi_dphi = realdata[26*n+j];
185 ps_all[start_index+j].B_z_dphi = realdata[27*n+j];
186 ps_all[start_index+j].B_r_dz = realdata[28*n+j];
187 ps_all[start_index+j].B_phi_dz = realdata[29*n+j];
188 ps_all[start_index+j].B_z_dz = realdata[30*n+j];
189 ps_all[start_index+j].mileage = realdata[31*n+j];
190 ps_all[start_index+j].err = errdata[j];
191 }
192
193 free(realdata);
194 free(intdata);
195 free(errdata);
196 }
197 }
198 else {
199
200 int start_index, n;
201 mpi_my_particles(&start_index, &n, n_tot, mpi_rank, mpi_size);
202
203 real* realdata;
204 realdata = malloc(n_real * n * sizeof(realdata));
205 integer* intdata;
206 intdata = malloc(n_int * n * sizeof(intdata));
207 a5err* errdata;
208 errdata = malloc(n_err * n * sizeof(intdata));
209
210 for(int j = 0; j < n; j++) {
211 realdata[0*n+j] = ps[j].r;
212 realdata[1*n+j] = ps[j].phi;
213 realdata[2*n+j] = ps[j].z;
214 realdata[3*n+j] = ps[j].ppar;
215 realdata[4*n+j] = ps[j].mu;
216 realdata[5*n+j] = ps[j].zeta;
217 realdata[6*n+j] = ps[j].rprt;
218 realdata[7*n+j] = ps[j].phiprt;
219 realdata[8*n+j] = ps[j].zprt;
220 realdata[9*n+j] = ps[j].p_r;
221 realdata[10*n+j] = ps[j].p_phi;
222 realdata[11*n+j] = ps[j].p_z;
223 realdata[12*n+j] = ps[j].mass;
224 realdata[13*n+j] = ps[j].charge;
225 intdata[0*n+j] = ps[j].anum;
226 intdata[1*n+j] = ps[j].znum;
227 realdata[14*n+j] = ps[j].weight;
228 realdata[15*n+j] = ps[j].time;
229 realdata[16*n+j] = ps[j].cputime;
230 realdata[17*n+j] = ps[j].rho;
231 realdata[18*n+j] = ps[j].theta;
232 intdata[2*n+j] = ps[j].id;
233 intdata[3*n+j] = ps[j].endcond;
234 intdata[4*n+j] = ps[j].walltile;
235 realdata[19*n+j] = ps[j].B_r;
236 realdata[20*n+j] = ps[j].B_phi;
237 realdata[21*n+j] = ps[j].B_z;
238 realdata[22*n+j] = ps[j].B_r_dr;
239 realdata[23*n+j] = ps[j].B_phi_dr;
240 realdata[24*n+j] = ps[j].B_z_dr;
241 realdata[25*n+j] = ps[j].B_r_dphi;
242 realdata[26*n+j] = ps[j].B_phi_dphi;
243 realdata[27*n+j] = ps[j].B_z_dphi;
244 realdata[28*n+j] = ps[j].B_r_dz;
245 realdata[29*n+j] = ps[j].B_phi_dz;
246 realdata[30*n+j] = ps[j].B_z_dz;
247 realdata[31*n+j] = ps[j].mileage;
248 errdata[j] = ps[j].err;
249 }
250
251 MPI_Send(realdata, n_real*n, mpi_type_real, 0, 0, MPI_COMM_WORLD);
252 MPI_Send(intdata, n_int*n, mpi_type_integer, 0, 0, MPI_COMM_WORLD);
253 MPI_Send(errdata, n_err*n, mpi_type_a5err, 0, 0, MPI_COMM_WORLD);
254
255 free(realdata);
256 free(intdata);
257 free(errdata);
258 }
259
260 *ps_gather = ps_all;
261 *n_gather = n_tot;
262
263#else
264
265 int start_index, n;
266 mpi_my_particles(&start_index, &n, n_tot, mpi_rank, mpi_size);
267
268 /* Only store particles for this process in Condor-style execution */
269 particle_state* ps_all = malloc(n * sizeof(particle_state));
270
271 for(int j = 0; j < n; j++) {
272 ps_all[j] = ps[j];
273 }
274
275 *ps_gather = ps_all;
276 *n_gather = n;
277
278#endif
279}
280
294void mpi_gather_diag(diag_data* data, int ntotal, int mpi_rank, int mpi_size,
295 int mpi_root) {
296#ifdef MPI
297
298 if(data->dist5D_collect) {
299 MPI_Reduce(
300 mpi_rank == mpi_root ? MPI_IN_PLACE : data->dist5D.histogram,
301 data->dist5D.histogram, data->dist5D.step_6 * data->dist5D.n_r,
302 MPI_INT, MPI_SUM, 0, MPI_COMM_WORLD);
303 }
304 if(data->dist6D_collect) {
305 MPI_Reduce(
306 mpi_rank == mpi_root ? MPI_IN_PLACE : data->dist6D.histogram,
307 data->dist6D.histogram, data->dist6D.step_7 * data->dist6D.n_r,
308 MPI_INT, MPI_SUM, 0, MPI_COMM_WORLD);
309 }
310 if(data->distrho5D_collect) {
311 MPI_Reduce(
312 mpi_rank == mpi_root ? MPI_IN_PLACE : data->distrho5D.histogram,
313 data->distrho5D.histogram,
314 data->distrho5D.step_6 * data->distrho5D.n_rho,
315 MPI_INT, MPI_SUM, 0, MPI_COMM_WORLD);
316 }
317 if(data->distrho6D_collect) {
318 MPI_Reduce(
319 mpi_rank == mpi_root ? MPI_IN_PLACE : data->distrho6D.histogram,
320 data->distrho6D.histogram,
321 data->distrho6D.step_7*data->distrho6D.n_rho,
322 MPI_INT, MPI_SUM, 0, MPI_COMM_WORLD);
323 }
324 if(data->distCOM_collect) {
325 MPI_Reduce(
326 mpi_rank == mpi_root ? MPI_IN_PLACE : data->distCOM.histogram,
327 data->distCOM.histogram, data->distCOM.step_2 * data->distCOM.n_mu,
328 MPI_INT, MPI_SUM, 0, MPI_COMM_WORLD);
329 }
330
331 if(data->diagorb_collect) {
332 size_t nreal = ntotal*data->diagorb.Npnt * sizeof(real);
333 if(mpi_rank == mpi_root) {
334 data->diagorb.id = realloc(data->diagorb.id, nreal);
335 data->diagorb.mileage = realloc(data->diagorb.mileage, nreal);
336 data->diagorb.r = realloc(data->diagorb.r, nreal);
337 data->diagorb.phi = realloc(data->diagorb.phi, nreal);
338 data->diagorb.z = realloc(data->diagorb.z, nreal);
339 if(data->diagorb.record_mode == DIAG_ORB_FO) {
340 data->diagorb.p_r = realloc(data->diagorb.p_r, nreal);
341 data->diagorb.p_phi = realloc(data->diagorb.p_phi, nreal);
342 data->diagorb.p_z = realloc(data->diagorb.p_z, nreal);
343 }
344 if(data->diagorb.record_mode == DIAG_ORB_GC) {
345 data->diagorb.ppar = realloc(data->diagorb.ppar, nreal);
346 data->diagorb.mu = realloc(data->diagorb.mu, nreal);
347 data->diagorb.zeta = realloc(data->diagorb.zeta, nreal);
348 }
349 if(data->diagorb.record_mode != DIAG_ORB_ML) {
350 data->diagorb.charge = realloc(data->diagorb.charge, nreal);
351 }
352 data->diagorb.weight = realloc(data->diagorb.weight, nreal);
353 data->diagorb.rho = realloc(data->diagorb.rho, nreal);
354 data->diagorb.theta = realloc(data->diagorb.theta, nreal);
355 data->diagorb.B_r = realloc(data->diagorb.B_r, nreal);
356 data->diagorb.B_phi = realloc(data->diagorb.B_phi, nreal);
357 data->diagorb.B_z = realloc(data->diagorb.B_z, nreal);
358 data->diagorb.simmode = realloc(data->diagorb.simmode, nreal);
359 if(data->diagorb.mode == DIAG_ORB_POINCARE) {
360 data->diagorb.pncrid = realloc(data->diagorb.pncrid, nreal);
361 data->diagorb.pncrdi = realloc(data->diagorb.pncrdi, nreal);
362 }
363 for(int i = 1; i < mpi_size; i++) {
364 int start_index, n;
365 mpi_my_particles(&start_index, &n, ntotal, i, mpi_size);
366 MPI_Recv(
367 &data->diagorb.id[start_index*data->diagorb.Npnt],
368 n*data->diagorb.Npnt, mpi_type_real, i, 0,
369 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
370 MPI_Recv(
371 &data->diagorb.mileage[start_index*data->diagorb.Npnt],
372 n*data->diagorb.Npnt, mpi_type_real, i, 0,
373 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
374 MPI_Recv(
375 &data->diagorb.r[start_index*data->diagorb.Npnt],
376 n*data->diagorb.Npnt, mpi_type_real, i, 0,
377 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
378 MPI_Recv(
379 &data->diagorb.phi[start_index*data->diagorb.Npnt],
380 n*data->diagorb.Npnt, mpi_type_real, i, 0,
381 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
382 MPI_Recv(
383 &data->diagorb.z[start_index*data->diagorb.Npnt],
384 n*data->diagorb.Npnt, mpi_type_real, i, 0,
385 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
386 if(data->diagorb.record_mode == DIAG_ORB_FO) {
387 MPI_Recv(
388 &data->diagorb.p_r[start_index*data->diagorb.Npnt],
389 n*data->diagorb.Npnt, mpi_type_real, i, 0,
390 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
391 MPI_Recv(
392 &data->diagorb.p_phi[start_index*data->diagorb.Npnt],
393 n*data->diagorb.Npnt, mpi_type_real, i, 0,
394 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
395 MPI_Recv(
396 &data->diagorb.p_z[start_index*data->diagorb.Npnt],
397 n*data->diagorb.Npnt, mpi_type_real, i, 0,
398 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
399 }
400 if(data->diagorb.record_mode == DIAG_ORB_GC) {
401 MPI_Recv(
402 &data->diagorb.ppar[start_index*data->diagorb.Npnt],
403 n*data->diagorb.Npnt, mpi_type_real, i, 0,
404 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
405 MPI_Recv(
406 &data->diagorb.mu[start_index*data->diagorb.Npnt],
407 n*data->diagorb.Npnt, mpi_type_real, i, 0,
408 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
409 MPI_Recv(
410 &data->diagorb.zeta[start_index*data->diagorb.Npnt],
411 n*data->diagorb.Npnt, mpi_type_real, i, 0,
412 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
413 }
414 if(data->diagorb.record_mode != DIAG_ORB_ML) {
415 MPI_Recv(
416 &data->diagorb.charge[start_index*data->diagorb.Npnt],
417 n*data->diagorb.Npnt, mpi_type_real, i, 0,
418 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
419 }
420 MPI_Recv(
421 &data->diagorb.weight[start_index*data->diagorb.Npnt],
422 n*data->diagorb.Npnt, mpi_type_real, i, 0,
423 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
424 MPI_Recv(
425 &data->diagorb.rho[start_index*data->diagorb.Npnt],
426 n*data->diagorb.Npnt, mpi_type_real, i, 0,
427 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
428 MPI_Recv(
429 &data->diagorb.theta[start_index*data->diagorb.Npnt],
430 n*data->diagorb.Npnt, mpi_type_real, i, 0,
431 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
432 MPI_Recv(
433 &data->diagorb.B_r[start_index*data->diagorb.Npnt],
434 n*data->diagorb.Npnt, mpi_type_real, i, 0,
435 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
436 MPI_Recv(
437 &data->diagorb.B_phi[start_index*data->diagorb.Npnt],
438 n*data->diagorb.Npnt, mpi_type_real, i, 0,
439 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
440 MPI_Recv(
441 &data->diagorb.B_z[start_index*data->diagorb.Npnt],
442 n*data->diagorb.Npnt, mpi_type_real, i, 0,
443 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
444 MPI_Recv(
445 &data->diagorb.simmode[start_index*data->diagorb.Npnt],
446 n*data->diagorb.Npnt, mpi_type_real, i, 0,
447 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
448 if(data->diagorb.mode == DIAG_ORB_POINCARE) {
449 MPI_Recv(
450 &data->diagorb.pncrid[start_index*data->diagorb.Npnt],
451 n*data->diagorb.Npnt, mpi_type_real, i, 0,
452 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
453 MPI_Recv(
454 &data->diagorb.pncrdi[start_index*data->diagorb.Npnt],
455 n*data->diagorb.Npnt, mpi_type_real, i, 0,
456 MPI_COMM_WORLD, MPI_STATUS_IGNORE);
457 }
458 }
459 }
460 else {
461 int start_index, n;
462 mpi_my_particles(&start_index, &n, ntotal, mpi_rank, mpi_size);
463 MPI_Send(data->diagorb.id, n*data->diagorb.Npnt,
464 mpi_type_real, 0, 0, MPI_COMM_WORLD);
465 MPI_Send(data->diagorb.mileage, n*data->diagorb.Npnt,
466 mpi_type_real, 0, 0, MPI_COMM_WORLD);
467 MPI_Send(data->diagorb.r, n*data->diagorb.Npnt,
468 mpi_type_real, 0, 0, MPI_COMM_WORLD);
469 MPI_Send(data->diagorb.phi, n*data->diagorb.Npnt,
470 mpi_type_real, 0, 0, MPI_COMM_WORLD);
471 MPI_Send(data->diagorb.z, n*data->diagorb.Npnt,
472 mpi_type_real, 0, 0, MPI_COMM_WORLD);
473 if(data->diagorb.record_mode == DIAG_ORB_FO) {
474 MPI_Send(data->diagorb.p_r, n*data->diagorb.Npnt,
475 mpi_type_real, 0, 0, MPI_COMM_WORLD);
476 MPI_Send(data->diagorb.p_phi, n*data->diagorb.Npnt,
477 mpi_type_real, 0, 0, MPI_COMM_WORLD);
478 MPI_Send(data->diagorb.p_z, n*data->diagorb.Npnt,
479 mpi_type_real, 0, 0, MPI_COMM_WORLD);
480 }
481 if(data->diagorb.record_mode == DIAG_ORB_GC) {
482 MPI_Send(data->diagorb.ppar, n*data->diagorb.Npnt,
483 mpi_type_real, 0, 0, MPI_COMM_WORLD);
484 MPI_Send(data->diagorb.mu, n*data->diagorb.Npnt,
485 mpi_type_real, 0, 0, MPI_COMM_WORLD);
486 MPI_Send(data->diagorb.zeta, n*data->diagorb.Npnt,
487 mpi_type_real, 0, 0, MPI_COMM_WORLD);
488 }
489 if(data->diagorb.record_mode != DIAG_ORB_ML) {
490 MPI_Send(data->diagorb.charge, n*data->diagorb.Npnt,
491 mpi_type_real, 0, 0, MPI_COMM_WORLD);
492 }
493 MPI_Send(data->diagorb.weight, n*data->diagorb.Npnt,
494 mpi_type_real, 0, 0, MPI_COMM_WORLD);
495 MPI_Send(data->diagorb.rho, n*data->diagorb.Npnt,
496 mpi_type_real, 0, 0, MPI_COMM_WORLD);
497 MPI_Send(data->diagorb.theta, n*data->diagorb.Npnt,
498 mpi_type_real, 0, 0, MPI_COMM_WORLD);
499 MPI_Send(data->diagorb.B_r, n*data->diagorb.Npnt,
500 mpi_type_real, 0, 0, MPI_COMM_WORLD);
501 MPI_Send(data->diagorb.B_phi, n*data->diagorb.Npnt,
502 mpi_type_real, 0, 0, MPI_COMM_WORLD);
503 MPI_Send(data->diagorb.B_z, n*data->diagorb.Npnt,
504 mpi_type_real, 0, 0, MPI_COMM_WORLD);
505 MPI_Send(data->diagorb.simmode, n*data->diagorb.Npnt,
506 mpi_type_real, 0, 0, MPI_COMM_WORLD);
507 if(data->diagorb.mode == DIAG_ORB_POINCARE) {
508 MPI_Send(data->diagorb.pncrid, n*data->diagorb.Npnt,
509 mpi_type_real, 0, 0, MPI_COMM_WORLD);
510 MPI_Send(data->diagorb.pncrdi, n*data->diagorb.Npnt,
511 mpi_type_real, 0, 0, MPI_COMM_WORLD);
512 }
513 }
514 }
515
516 if(data->diagtrcof_collect) {
517
518 if(mpi_rank == mpi_root) {
519 data->diagtrcof.id = realloc(
520 data->diagtrcof.id, ntotal*data->diagtrcof.Nmrk*sizeof(int)
521 );
522 data->diagtrcof.Kcoef = realloc(
523 data->diagtrcof.Kcoef, ntotal*data->diagtrcof.Nmrk*sizeof(real)
524 );
525 data->diagtrcof.Dcoef = realloc(
526 data->diagtrcof.Dcoef,
527 ntotal*data->diagtrcof.Nmrk*sizeof(real)
528 );
529 for(int i = 1; i < mpi_size; i++) {
530 int start_index, n;
531 mpi_my_particles(&start_index, &n, ntotal, i, mpi_size);
532
533 MPI_Recv(&data->diagtrcof.id[start_index], n, mpi_type_integer,
534 i, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
535 MPI_Recv(&data->diagtrcof.Kcoef[start_index], n, mpi_type_real,
536 i, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
537 MPI_Recv(&data->diagtrcof.Dcoef[start_index], n, mpi_type_real,
538 i, 0, MPI_COMM_WORLD, MPI_STATUS_IGNORE);
539 }
540 }
541 else {
542 int start_index, n;
543 mpi_my_particles(&start_index, &n, ntotal, mpi_rank, mpi_size);
544 MPI_Send(data->diagtrcof.id, n, mpi_type_integer,
545 0, 0, MPI_COMM_WORLD);
546 MPI_Send(data->diagtrcof.Kcoef, n, mpi_type_real,
547 0, 0, MPI_COMM_WORLD);
548 MPI_Send(data->diagtrcof.Dcoef, n, mpi_type_real,
549 0, 0, MPI_COMM_WORLD);
550 }
551 }
552
553#endif
554}
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_finalize()
Finalize MPI.
void mpi_interface_init(int argc, char **argv, int *mpi_rank, int *mpi_size, int *mpi_root)
Initialize 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