1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * GROwing Monsters And Cloning Shrimps
47 #include "gmx_fatal.h"
52 #include "sighandler.h"
61 #define PP_PME_CHARGE (1<<0)
62 #define PP_PME_CHARGEB (1<<1)
63 #define PP_PME_COORD (1<<2)
64 #define PP_PME_FEP (1<<3)
65 #define PP_PME_ENER_VIR (1<<4)
66 #define PP_PME_FINISH (1<<5)
67 #define PP_PME_SWITCHGRID (1<<6)
68 #define PP_PME_RESETCOUNTERS (1<<7)
71 #define PME_PP_SIGSTOP (1<<0)
72 #define PME_PP_SIGSTOPNSS (1<<1)
74 typedef struct gmx_pme_pp {
76 MPI_Comm mpi_comm_mysim;
78 int nnode; /* The number of PP node to communicate with */
79 int *node; /* The PP node ranks */
80 int node_peer; /* The peer PP node rank */
81 int *nat; /* The number of atom for each PP node */
82 int flags_charge; /* The flags sent along with the last charges */
94 typedef struct gmx_pme_comm_n_box {
101 gmx_large_int_t step;
102 ivec grid_size; /* For PME grid tuning */
103 real ewaldcoeff; /* For PME grid tuning */
104 } gmx_pme_comm_n_box_t;
111 gmx_stop_cond_t stop_cond;
112 } gmx_pme_comm_vir_ene_t;
117 gmx_pme_pp_t gmx_pme_pp_init(t_commrec *cr)
119 struct gmx_pme_pp *pme_pp;
125 pme_pp->mpi_comm_mysim = cr->mpi_comm_mysim;
126 MPI_Comm_rank(cr->mpi_comm_mygroup, &rank);
127 get_pme_ddnodes(cr, rank, &pme_pp->nnode, &pme_pp->node, &pme_pp->node_peer);
128 snew(pme_pp->nat, pme_pp->nnode);
129 snew(pme_pp->req, 2*pme_pp->nnode);
130 snew(pme_pp->stat, 2*pme_pp->nnode);
132 pme_pp->flags_charge = 0;
138 /* This should be faster with a real non-blocking MPI implementation */
139 /* #define GMX_PME_DELAYED_WAIT */
141 static void gmx_pme_send_q_x_wait(gmx_domdec_t *dd)
146 MPI_Waitall(dd->nreq_pme, dd->req_pme, MPI_STATUSES_IGNORE);
152 static void gmx_pme_send_q_x(t_commrec *cr, int flags,
153 real *chargeA, real *chargeB,
156 int maxshift_x, int maxshift_y,
157 gmx_large_int_t step)
160 gmx_pme_comm_n_box_t *cnb;
168 fprintf(debug, "PP node %d sending to PME node %d: %d%s%s\n",
169 cr->sim_nodeid, dd->pme_nodeid, n,
170 flags & PP_PME_CHARGE ? " charges" : "",
171 flags & PP_PME_COORD ? " coordinates" : "");
174 #ifdef GMX_PME_DELAYED_WAIT
175 /* When can not use cnb until pending communication has finished */
176 gmx_pme_send_x_q_wait(dd);
179 if (dd->pme_receive_vir_ener)
181 /* Peer PP node: communicate all data */
190 cnb->maxshift_x = maxshift_x;
191 cnb->maxshift_y = maxshift_y;
192 cnb->lambda = lambda;
194 if (flags & PP_PME_COORD)
196 copy_mat(box, cnb->box);
199 MPI_Isend(cnb, sizeof(*cnb), MPI_BYTE,
200 dd->pme_nodeid, 0, cr->mpi_comm_mysim,
201 &dd->req_pme[dd->nreq_pme++]);
204 else if (flags & PP_PME_CHARGE)
207 /* Communicate only the number of atoms */
208 MPI_Isend(&n, sizeof(n), MPI_BYTE,
209 dd->pme_nodeid, 0, cr->mpi_comm_mysim,
210 &dd->req_pme[dd->nreq_pme++]);
217 if (flags & PP_PME_CHARGE)
219 MPI_Isend(chargeA, n*sizeof(real), MPI_BYTE,
220 dd->pme_nodeid, 1, cr->mpi_comm_mysim,
221 &dd->req_pme[dd->nreq_pme++]);
223 if (flags & PP_PME_CHARGEB)
225 MPI_Isend(chargeB, n*sizeof(real), MPI_BYTE,
226 dd->pme_nodeid, 2, cr->mpi_comm_mysim,
227 &dd->req_pme[dd->nreq_pme++]);
229 if (flags & PP_PME_COORD)
231 MPI_Isend(x[0], n*sizeof(rvec), MPI_BYTE,
232 dd->pme_nodeid, 3, cr->mpi_comm_mysim,
233 &dd->req_pme[dd->nreq_pme++]);
237 #ifndef GMX_PME_DELAYED_WAIT
238 /* Wait for the data to arrive */
239 /* We can skip this wait as we are sure x and q will not be modified
240 * before the next call to gmx_pme_send_x_q or gmx_pme_receive_f.
242 gmx_pme_send_q_x_wait(dd);
247 void gmx_pme_send_q(t_commrec *cr,
248 gmx_bool bFreeEnergy, real *chargeA, real *chargeB,
249 int maxshift_x, int maxshift_y)
253 flags = PP_PME_CHARGE;
256 flags |= PP_PME_CHARGEB;
259 gmx_pme_send_q_x(cr, flags,
260 chargeA, chargeB, NULL, NULL, 0, maxshift_x, maxshift_y, -1);
263 void gmx_pme_send_x(t_commrec *cr, matrix box, rvec *x,
264 gmx_bool bFreeEnergy, real lambda,
266 gmx_large_int_t step)
270 flags = PP_PME_COORD;
277 flags |= PP_PME_ENER_VIR;
280 gmx_pme_send_q_x(cr, flags, NULL, NULL, box, x, lambda, 0, 0, step);
283 void gmx_pme_send_finish(t_commrec *cr)
287 flags = PP_PME_FINISH;
289 gmx_pme_send_q_x(cr, flags, NULL, NULL, NULL, NULL, 0, 0, 0, -1);
292 void gmx_pme_send_switchgrid(t_commrec *cr, ivec grid_size, real ewaldcoeff)
295 gmx_pme_comm_n_box_t cnb;
297 /* Only let one PP node signal each PME node */
298 if (cr->dd->pme_receive_vir_ener)
300 cnb.flags = PP_PME_SWITCHGRID;
301 copy_ivec(grid_size, cnb.grid_size);
302 cnb.ewaldcoeff = ewaldcoeff;
304 /* We send this, uncommon, message blocking to simplify the code */
305 MPI_Send(&cnb, sizeof(cnb), MPI_BYTE,
306 cr->dd->pme_nodeid, 0, cr->mpi_comm_mysim);
311 void gmx_pme_send_resetcounters(t_commrec *cr, gmx_large_int_t step)
314 gmx_pme_comm_n_box_t cnb;
316 /* Only let one PP node signal each PME node */
317 if (cr->dd->pme_receive_vir_ener)
319 cnb.flags = PP_PME_RESETCOUNTERS;
322 /* We send this, uncommon, message blocking to simplify the code */
323 MPI_Send(&cnb, sizeof(cnb), MPI_BYTE,
324 cr->dd->pme_nodeid, 0, cr->mpi_comm_mysim);
329 int gmx_pme_recv_q_x(struct gmx_pme_pp *pme_pp,
331 real **chargeA, real **chargeB,
332 matrix box, rvec **x, rvec **f,
333 int *maxshift_x, int *maxshift_y,
334 gmx_bool *bFreeEnergy, real *lambda,
336 gmx_large_int_t *step,
337 ivec grid_size, real *ewaldcoeff)
339 gmx_pme_comm_n_box_t cnb;
340 int nat = 0, q, messages, sender;
345 /* avoid compiler warning about unused variable without MPI support */
350 /* Receive the send count, box and time step from the peer PP node */
351 MPI_Recv(&cnb, sizeof(cnb), MPI_BYTE,
352 pme_pp->node_peer, 0,
353 pme_pp->mpi_comm_mysim, MPI_STATUS_IGNORE);
357 fprintf(debug, "PME only node receiving:%s%s%s%s%s\n",
358 (cnb.flags & PP_PME_CHARGE) ? " charges" : "",
359 (cnb.flags & PP_PME_COORD ) ? " coordinates" : "",
360 (cnb.flags & PP_PME_FINISH) ? " finish" : "",
361 (cnb.flags & PP_PME_SWITCHGRID) ? " switch grid" : "",
362 (cnb.flags & PP_PME_RESETCOUNTERS) ? " reset counters" : "");
365 if (cnb.flags & PP_PME_SWITCHGRID)
367 /* Special case, receive the new parameters and return */
368 copy_ivec(cnb.grid_size, grid_size);
369 *ewaldcoeff = cnb.ewaldcoeff;
371 return pmerecvqxSWITCHGRID;
374 if (cnb.flags & PP_PME_RESETCOUNTERS)
376 /* Special case, receive the step and return */
379 return pmerecvqxRESETCOUNTERS;
382 if (cnb.flags & PP_PME_CHARGE)
384 /* Receive the send counts from the other PP nodes */
385 for (sender = 0; sender < pme_pp->nnode; sender++)
387 if (pme_pp->node[sender] == pme_pp->node_peer)
389 pme_pp->nat[sender] = cnb.natoms;
393 MPI_Irecv(&(pme_pp->nat[sender]), sizeof(pme_pp->nat[0]),
395 pme_pp->node[sender], 0,
396 pme_pp->mpi_comm_mysim, &pme_pp->req[messages++]);
399 MPI_Waitall(messages, pme_pp->req, pme_pp->stat);
403 for (sender = 0; sender < pme_pp->nnode; sender++)
405 nat += pme_pp->nat[sender];
408 if (nat > pme_pp->nalloc)
410 pme_pp->nalloc = over_alloc_dd(nat);
411 srenew(pme_pp->chargeA, pme_pp->nalloc);
412 if (cnb.flags & PP_PME_CHARGEB)
414 srenew(pme_pp->chargeB, pme_pp->nalloc);
416 srenew(pme_pp->x, pme_pp->nalloc);
417 srenew(pme_pp->f, pme_pp->nalloc);
420 /* maxshift is sent when the charges are sent */
421 *maxshift_x = cnb.maxshift_x;
422 *maxshift_y = cnb.maxshift_y;
424 /* Receive the charges in place */
425 for (q = 0; q < ((cnb.flags & PP_PME_CHARGEB) ? 2 : 1); q++)
429 charge_pp = pme_pp->chargeA;
433 charge_pp = pme_pp->chargeB;
436 for (sender = 0; sender < pme_pp->nnode; sender++)
438 if (pme_pp->nat[sender] > 0)
440 MPI_Irecv(charge_pp+nat,
441 pme_pp->nat[sender]*sizeof(real),
443 pme_pp->node[sender], 1+q,
444 pme_pp->mpi_comm_mysim,
445 &pme_pp->req[messages++]);
446 nat += pme_pp->nat[sender];
449 fprintf(debug, "Received from PP node %d: %d "
451 pme_pp->node[sender], pme_pp->nat[sender]);
457 pme_pp->flags_charge = cnb.flags;
460 if (cnb.flags & PP_PME_COORD)
462 if (!(pme_pp->flags_charge & PP_PME_CHARGE))
464 gmx_incons("PME-only node received coordinates before charges"
468 /* The box, FE flag and lambda are sent along with the coordinates
470 copy_mat(cnb.box, box);
471 *bFreeEnergy = (cnb.flags & PP_PME_FEP);
472 *lambda = cnb.lambda;
473 *bEnerVir = (cnb.flags & PP_PME_ENER_VIR);
475 if (*bFreeEnergy && !(pme_pp->flags_charge & PP_PME_CHARGEB))
477 gmx_incons("PME-only node received free energy request, but "
478 "did not receive B-state charges");
481 /* Receive the coordinates in place */
483 for (sender = 0; sender < pme_pp->nnode; sender++)
485 if (pme_pp->nat[sender] > 0)
487 MPI_Irecv(pme_pp->x[nat], pme_pp->nat[sender]*sizeof(rvec),
489 pme_pp->node[sender], 3,
490 pme_pp->mpi_comm_mysim, &pme_pp->req[messages++]);
491 nat += pme_pp->nat[sender];
494 fprintf(debug, "Received from PP node %d: %d "
496 pme_pp->node[sender], pme_pp->nat[sender]);
502 /* Wait for the coordinates and/or charges to arrive */
503 MPI_Waitall(messages, pme_pp->req, pme_pp->stat);
506 while (!(cnb.flags & (PP_PME_COORD | PP_PME_FINISH)));
512 *chargeA = pme_pp->chargeA;
513 *chargeB = pme_pp->chargeB;
517 return ((cnb.flags & PP_PME_FINISH) ? pmerecvqxFINISH : pmerecvqxX);
520 static void receive_virial_energy(t_commrec *cr,
521 matrix vir, real *energy, real *dvdlambda,
524 gmx_pme_comm_vir_ene_t cve;
526 if (cr->dd->pme_receive_vir_ener)
531 "PP node %d receiving from PME node %d: virial and energy\n",
532 cr->sim_nodeid, cr->dd->pme_nodeid);
535 MPI_Recv(&cve, sizeof(cve), MPI_BYTE, cr->dd->pme_nodeid, 1, cr->mpi_comm_mysim,
538 memset(&cve, 0, sizeof(cve));
541 m_add(vir, cve.vir, vir);
542 *energy = cve.energy;
543 *dvdlambda += cve.dvdlambda;
544 *pme_cycles = cve.cycles;
546 if (cve.stop_cond != gmx_stop_cond_none)
548 gmx_set_stop_condition(cve.stop_cond);
558 void gmx_pme_receive_f(t_commrec *cr,
559 rvec f[], matrix vir,
560 real *energy, real *dvdlambda,
565 #ifdef GMX_PME_DELAYED_WAIT
566 /* Wait for the x request to finish */
567 gmx_pme_send_q_x_wait(cr->dd);
570 natoms = cr->dd->nat_home;
572 if (natoms > cr->dd->pme_recv_f_alloc)
574 cr->dd->pme_recv_f_alloc = over_alloc_dd(natoms);
575 srenew(cr->dd->pme_recv_f_buf, cr->dd->pme_recv_f_alloc);
579 MPI_Recv(cr->dd->pme_recv_f_buf[0],
580 natoms*sizeof(rvec), MPI_BYTE,
581 cr->dd->pme_nodeid, 0, cr->mpi_comm_mysim,
585 for (i = 0; i < natoms; i++)
587 rvec_inc(f[i], cr->dd->pme_recv_f_buf[i]);
591 receive_virial_energy(cr, vir, energy, dvdlambda, pme_cycles);
594 void gmx_pme_send_force_vir_ener(struct gmx_pme_pp *pme_pp,
596 real energy, real dvdlambda,
599 gmx_pme_comm_vir_ene_t cve;
600 int messages, ind_start, ind_end, receiver;
604 /* Now the evaluated forces have to be transferred to the PP nodes */
607 for (receiver = 0; receiver < pme_pp->nnode; receiver++)
610 ind_end = ind_start + pme_pp->nat[receiver];
612 if (MPI_Isend(f[ind_start], (ind_end-ind_start)*sizeof(rvec), MPI_BYTE,
613 pme_pp->node[receiver], 0,
614 pme_pp->mpi_comm_mysim, &pme_pp->req[messages++]) != 0)
616 gmx_comm("MPI_Isend failed in do_pmeonly");
621 /* send virial and energy to our last PP node */
622 copy_mat(vir, cve.vir);
624 cve.dvdlambda = dvdlambda;
625 /* check for the signals to send back to a PP node */
626 cve.stop_cond = gmx_get_stop_condition();
632 fprintf(debug, "PME node sending to PP node %d: virial and energy\n",
636 MPI_Isend(&cve, sizeof(cve), MPI_BYTE,
637 pme_pp->node_peer, 1,
638 pme_pp->mpi_comm_mysim, &pme_pp->req[messages++]);
640 /* Wait for the forces to arrive */
641 MPI_Waitall(messages, pme_pp->req, pme_pp->stat);