2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
48 #include "gmx_fatal.h"
53 #include "sighandler.h"
55 #include "gromacs/utility/gmxmpi.h"
58 eCommType_ChargeA, eCommType_ChargeB, eCommType_C6A, eCommType_C6B,
59 eCommType_SigmaA, eCommType_SigmaB, eCommType_NR, eCommType_COORD,
63 /* Some parts of the code(gmx_pme_send_q, gmx_pme_recv_q_x) assume
64 * that the six first flags are exactly in this order.
65 * If more PP_PME_...-flags are to be introduced be aware of some of
66 * the PME-specific flags in pme.h. Currently, they are also passed
70 #define PP_PME_CHARGE (1<<0)
71 #define PP_PME_CHARGEB (1<<1)
72 #define PP_PME_C6 (1<<2)
73 #define PP_PME_C6B (1<<3)
74 #define PP_PME_SIGMA (1<<4)
75 #define PP_PME_SIGMAB (1<<5)
76 #define PP_PME_COORD (1<<6)
77 #define PP_PME_FEP_Q (1<<7)
78 #define PP_PME_FEP_LJ (1<<8)
79 #define PP_PME_ENER_VIR (1<<9)
80 #define PP_PME_FINISH (1<<10)
81 #define PP_PME_SWITCHGRID (1<<11)
82 #define PP_PME_RESETCOUNTERS (1<<12)
84 #define PME_PP_SIGSTOP (1<<0)
85 #define PME_PP_SIGSTOPNSS (1<<1)
87 typedef struct gmx_pme_pp {
89 MPI_Comm mpi_comm_mysim;
91 int nnode; /* The number of PP node to communicate with */
92 int *node; /* The PP node ranks */
93 int node_peer; /* The peer PP node rank */
94 int *nat; /* The number of atom for each PP node */
95 int flags_charge; /* The flags sent along with the last charges */
111 typedef struct gmx_pme_comm_n_box {
119 gmx_large_int_t step;
120 ivec grid_size; /* For PME grid tuning */
121 real ewaldcoeff_q; /* For PME grid tuning */
123 } gmx_pme_comm_n_box_t;
133 gmx_stop_cond_t stop_cond;
134 } gmx_pme_comm_vir_ene_t;
139 gmx_pme_pp_t gmx_pme_pp_init(t_commrec gmx_unused *cr)
141 struct gmx_pme_pp *pme_pp;
147 pme_pp->mpi_comm_mysim = cr->mpi_comm_mysim;
148 MPI_Comm_rank(cr->mpi_comm_mygroup, &rank);
149 get_pme_ddnodes(cr, rank, &pme_pp->nnode, &pme_pp->node, &pme_pp->node_peer);
150 snew(pme_pp->nat, pme_pp->nnode);
151 snew(pme_pp->req, eCommType_NR*pme_pp->nnode);
152 snew(pme_pp->stat, eCommType_NR*pme_pp->nnode);
154 pme_pp->flags_charge = 0;
160 /* This should be faster with a real non-blocking MPI implementation */
161 /* #define GMX_PME_DELAYED_WAIT */
163 static void gmx_pme_send_params_coords_wait(gmx_domdec_t gmx_unused *dd)
168 MPI_Waitall(dd->nreq_pme, dd->req_pme, MPI_STATUSES_IGNORE);
174 static void gmx_pme_send_params_coords(t_commrec *cr, int flags,
175 real gmx_unused *chargeA, real gmx_unused *chargeB,
176 real gmx_unused *c6A, real gmx_unused *c6B,
177 real gmx_unused *sigmaA, real gmx_unused *sigmaB,
178 matrix box, rvec gmx_unused *x,
179 real lambda_q, real lambda_lj,
180 int maxshift_x, int maxshift_y,
181 gmx_large_int_t step)
184 gmx_pme_comm_n_box_t *cnb;
192 fprintf(debug, "PP node %d sending to PME node %d: %d%s%s\n",
193 cr->sim_nodeid, dd->pme_nodeid, n,
194 flags & PP_PME_CHARGE ? " charges" : "",
195 flags & PP_PME_COORD ? " coordinates" : "");
198 #ifdef GMX_PME_DELAYED_WAIT
199 /* When can not use cnb until pending communication has finished */
200 gmx_pme_send_params_coords_wait(dd);
203 if (dd->pme_receive_vir_ener)
205 /* Peer PP node: communicate all data */
214 cnb->maxshift_x = maxshift_x;
215 cnb->maxshift_y = maxshift_y;
216 cnb->lambda_q = lambda_q;
217 cnb->lambda_lj = lambda_lj;
219 if (flags & PP_PME_COORD)
221 copy_mat(box, cnb->box);
224 MPI_Isend(cnb, sizeof(*cnb), MPI_BYTE,
225 dd->pme_nodeid, eCommType_CNB, cr->mpi_comm_mysim,
226 &dd->req_pme[dd->nreq_pme++]);
229 else if (flags & (PP_PME_CHARGE | PP_PME_C6 | PP_PME_SIGMA))
232 /* Communicate only the number of atoms */
233 MPI_Isend(&n, sizeof(n), MPI_BYTE,
234 dd->pme_nodeid, eCommType_CNB, cr->mpi_comm_mysim,
235 &dd->req_pme[dd->nreq_pme++]);
242 if (flags & PP_PME_CHARGE)
244 MPI_Isend(chargeA, n*sizeof(real), MPI_BYTE,
245 dd->pme_nodeid, eCommType_ChargeA, cr->mpi_comm_mysim,
246 &dd->req_pme[dd->nreq_pme++]);
248 if (flags & PP_PME_CHARGEB)
250 MPI_Isend(chargeB, n*sizeof(real), MPI_BYTE,
251 dd->pme_nodeid, eCommType_ChargeB, cr->mpi_comm_mysim,
252 &dd->req_pme[dd->nreq_pme++]);
254 if (flags & PP_PME_C6)
256 MPI_Isend(c6A, n*sizeof(real), MPI_BYTE,
257 dd->pme_nodeid, eCommType_C6A, cr->mpi_comm_mysim,
258 &dd->req_pme[dd->nreq_pme++]);
260 if (flags & PP_PME_C6B)
262 MPI_Isend(c6B, n*sizeof(real), MPI_BYTE,
263 dd->pme_nodeid, eCommType_C6B, cr->mpi_comm_mysim,
264 &dd->req_pme[dd->nreq_pme++]);
266 if (flags & PP_PME_SIGMA)
268 MPI_Isend(sigmaA, n*sizeof(real), MPI_BYTE,
269 dd->pme_nodeid, eCommType_SigmaA, cr->mpi_comm_mysim,
270 &dd->req_pme[dd->nreq_pme++]);
272 if (flags & PP_PME_SIGMAB)
274 MPI_Isend(sigmaB, n*sizeof(real), MPI_BYTE,
275 dd->pme_nodeid, eCommType_SigmaB, cr->mpi_comm_mysim,
276 &dd->req_pme[dd->nreq_pme++]);
278 if (flags & PP_PME_COORD)
280 MPI_Isend(x[0], n*sizeof(rvec), MPI_BYTE,
281 dd->pme_nodeid, eCommType_COORD, cr->mpi_comm_mysim,
282 &dd->req_pme[dd->nreq_pme++]);
286 #ifndef GMX_PME_DELAYED_WAIT
287 /* Wait for the data to arrive */
288 /* We can skip this wait as we are sure x and q will not be modified
289 * before the next call to gmx_pme_send_x_q or gmx_pme_receive_f.
291 gmx_pme_send_params_coords_wait(dd);
296 void gmx_pme_send_parameters(t_commrec *cr,
297 gmx_bool bFreeEnergy_q, gmx_bool bFreeEnergy_lj,
298 real *chargeA, real *chargeB,
299 real *c6A, real *c6B, real *sigmaA, real *sigmaB,
300 int maxshift_x, int maxshift_y)
304 flags = PP_PME_CHARGE | PP_PME_C6 | PP_PME_SIGMA;
305 if (bFreeEnergy_q || bFreeEnergy_lj)
307 /* Assumes that the B state flags are in the bits just above
308 * the ones for the A state. */
309 flags |= (flags << 1);
312 gmx_pme_send_params_coords(cr, flags, chargeA, chargeB, c6A, c6B, sigmaA, sigmaB,
313 NULL, NULL, 0, 0, maxshift_x, maxshift_y, -1);
316 void gmx_pme_send_coordinates(t_commrec *cr, matrix box, rvec *x,
317 gmx_bool bFreeEnergy_q, gmx_bool bFreeEnergy_lj,
318 real lambda_q, real lambda_lj,
319 gmx_bool bEnerVir, int pme_flags,
320 gmx_large_int_t step)
324 flags = pme_flags | PP_PME_COORD;
327 flags |= PP_PME_FEP_Q;
331 flags |= PP_PME_FEP_LJ;
335 flags |= PP_PME_ENER_VIR;
337 gmx_pme_send_params_coords(cr, flags, NULL, NULL, NULL, NULL, NULL, NULL,
338 box, x, lambda_q, lambda_lj, 0, 0, step);
341 void gmx_pme_send_finish(t_commrec *cr)
345 flags = PP_PME_FINISH;
347 gmx_pme_send_params_coords(cr, flags, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, 0, 0, 0, 0, -1);
350 void gmx_pme_send_switchgrid(t_commrec gmx_unused *cr,
351 ivec gmx_unused grid_size,
352 real gmx_unused ewaldcoeff_q,
353 real gmx_unused ewaldcoeff_lj)
356 gmx_pme_comm_n_box_t cnb;
358 /* Only let one PP node signal each PME node */
359 if (cr->dd->pme_receive_vir_ener)
361 cnb.flags = PP_PME_SWITCHGRID;
362 copy_ivec(grid_size, cnb.grid_size);
363 cnb.ewaldcoeff_q = ewaldcoeff_q;
364 cnb.ewaldcoeff_lj = ewaldcoeff_lj;
366 /* We send this, uncommon, message blocking to simplify the code */
367 MPI_Send(&cnb, sizeof(cnb), MPI_BYTE,
368 cr->dd->pme_nodeid, eCommType_CNB, cr->mpi_comm_mysim);
373 void gmx_pme_send_resetcounters(t_commrec gmx_unused *cr, gmx_large_int_t gmx_unused step)
376 gmx_pme_comm_n_box_t cnb;
378 /* Only let one PP node signal each PME node */
379 if (cr->dd->pme_receive_vir_ener)
381 cnb.flags = PP_PME_RESETCOUNTERS;
384 /* We send this, uncommon, message blocking to simplify the code */
385 MPI_Send(&cnb, sizeof(cnb), MPI_BYTE,
386 cr->dd->pme_nodeid, eCommType_CNB, cr->mpi_comm_mysim);
391 int gmx_pme_recv_params_coords(struct gmx_pme_pp *pme_pp,
399 matrix gmx_unused box,
402 int gmx_unused *maxshift_x,
403 int gmx_unused *maxshift_y,
404 gmx_bool gmx_unused *bFreeEnergy_q,
405 gmx_bool gmx_unused *bFreeEnergy_lj,
406 real gmx_unused *lambda_q,
407 real gmx_unused *lambda_lj,
408 gmx_bool gmx_unused *bEnerVir,
410 gmx_large_int_t gmx_unused *step,
411 ivec gmx_unused grid_size,
412 real gmx_unused *ewaldcoeff_q,
413 real gmx_unused *ewaldcoeff_lj)
415 gmx_pme_comm_n_box_t cnb;
416 int nat = 0, q, messages, sender;
421 /* avoid compiler warning about unused variable without MPI support */
427 /* Receive the send count, box and time step from the peer PP node */
428 MPI_Recv(&cnb, sizeof(cnb), MPI_BYTE,
429 pme_pp->node_peer, eCommType_CNB,
430 pme_pp->mpi_comm_mysim, MPI_STATUS_IGNORE);
434 fprintf(debug, "PME only node receiving:%s%s%s%s%s\n",
435 (cnb.flags & PP_PME_CHARGE) ? " charges" : "",
436 (cnb.flags & PP_PME_COORD ) ? " coordinates" : "",
437 (cnb.flags & PP_PME_FINISH) ? " finish" : "",
438 (cnb.flags & PP_PME_SWITCHGRID) ? " switch grid" : "",
439 (cnb.flags & PP_PME_RESETCOUNTERS) ? " reset counters" : "");
442 if (cnb.flags & PP_PME_SWITCHGRID)
444 /* Special case, receive the new parameters and return */
445 copy_ivec(cnb.grid_size, grid_size);
446 *ewaldcoeff_q = cnb.ewaldcoeff_q;
447 *ewaldcoeff_lj = cnb.ewaldcoeff_lj;
448 return pmerecvqxSWITCHGRID;
451 if (cnb.flags & PP_PME_RESETCOUNTERS)
453 /* Special case, receive the step and return */
456 return pmerecvqxRESETCOUNTERS;
459 if (cnb.flags & (PP_PME_CHARGE | PP_PME_C6 | PP_PME_SIGMA))
461 /* Receive the send counts from the other PP nodes */
462 for (sender = 0; sender < pme_pp->nnode; sender++)
464 if (pme_pp->node[sender] == pme_pp->node_peer)
466 pme_pp->nat[sender] = cnb.natoms;
470 MPI_Irecv(&(pme_pp->nat[sender]), sizeof(pme_pp->nat[0]),
472 pme_pp->node[sender], eCommType_CNB,
473 pme_pp->mpi_comm_mysim, &pme_pp->req[messages++]);
476 MPI_Waitall(messages, pme_pp->req, pme_pp->stat);
480 for (sender = 0; sender < pme_pp->nnode; sender++)
482 nat += pme_pp->nat[sender];
485 if (nat > pme_pp->nalloc)
487 pme_pp->nalloc = over_alloc_dd(nat);
488 if (cnb.flags & PP_PME_CHARGE)
490 srenew(pme_pp->chargeA, pme_pp->nalloc);
492 if (cnb.flags & PP_PME_CHARGEB)
494 srenew(pme_pp->chargeB, pme_pp->nalloc);
496 if (cnb.flags & PP_PME_C6)
498 srenew(pme_pp->c6A, pme_pp->nalloc);
500 if (cnb.flags & PP_PME_C6B)
502 srenew(pme_pp->c6B, pme_pp->nalloc);
504 if (cnb.flags & PP_PME_SIGMA)
506 srenew(pme_pp->sigmaA, pme_pp->nalloc);
508 if (cnb.flags & PP_PME_SIGMAB)
510 srenew(pme_pp->sigmaB, pme_pp->nalloc);
512 srenew(pme_pp->x, pme_pp->nalloc);
513 srenew(pme_pp->f, pme_pp->nalloc);
516 /* maxshift is sent when the charges are sent */
517 *maxshift_x = cnb.maxshift_x;
518 *maxshift_y = cnb.maxshift_y;
520 /* Receive the charges in place */
521 for (q = 0; q < eCommType_NR; q++)
523 if (!(cnb.flags & (PP_PME_CHARGE<<q)))
529 case eCommType_ChargeA: charge_pp = pme_pp->chargeA; break;
530 case eCommType_ChargeB: charge_pp = pme_pp->chargeB; break;
531 case eCommType_C6A: charge_pp = pme_pp->c6A; break;
532 case eCommType_C6B: charge_pp = pme_pp->c6B; break;
533 case eCommType_SigmaA: charge_pp = pme_pp->sigmaA; break;
534 case eCommType_SigmaB: charge_pp = pme_pp->sigmaB; break;
535 default: gmx_incons("Wrong eCommType");
538 for (sender = 0; sender < pme_pp->nnode; sender++)
540 if (pme_pp->nat[sender] > 0)
542 MPI_Irecv(charge_pp+nat,
543 pme_pp->nat[sender]*sizeof(real),
545 pme_pp->node[sender], q,
546 pme_pp->mpi_comm_mysim,
547 &pme_pp->req[messages++]);
548 nat += pme_pp->nat[sender];
551 fprintf(debug, "Received from PP node %d: %d "
553 pme_pp->node[sender], pme_pp->nat[sender]);
559 pme_pp->flags_charge = cnb.flags;
562 if (cnb.flags & PP_PME_COORD)
564 if (!(pme_pp->flags_charge & (PP_PME_CHARGE | PP_PME_C6)))
566 gmx_incons("PME-only node received coordinates before charges and/or C6-values"
570 /* The box, FE flag and lambda are sent along with the coordinates
572 copy_mat(cnb.box, box);
573 *bFreeEnergy_q = (cnb.flags & PP_PME_FEP_Q);
574 *bFreeEnergy_lj = (cnb.flags & PP_PME_FEP_LJ);
575 *lambda_q = cnb.lambda_q;
576 *lambda_lj = cnb.lambda_lj;
577 *bEnerVir = (cnb.flags & PP_PME_ENER_VIR);
578 *pme_flags = cnb.flags;
580 if (*bFreeEnergy_q && !(pme_pp->flags_charge & PP_PME_CHARGEB))
582 gmx_incons("PME-only node received free energy request, but "
583 "did not receive B-state charges");
586 if (*bFreeEnergy_lj && !(pme_pp->flags_charge & PP_PME_C6B))
588 gmx_incons("PME-only node received free energy request, but "
589 "did not receive B-state C6-values");
592 /* Receive the coordinates in place */
594 for (sender = 0; sender < pme_pp->nnode; sender++)
596 if (pme_pp->nat[sender] > 0)
598 MPI_Irecv(pme_pp->x[nat], pme_pp->nat[sender]*sizeof(rvec),
600 pme_pp->node[sender], eCommType_COORD,
601 pme_pp->mpi_comm_mysim, &pme_pp->req[messages++]);
602 nat += pme_pp->nat[sender];
605 fprintf(debug, "Received from PP node %d: %d "
607 pme_pp->node[sender], pme_pp->nat[sender]);
613 /* Wait for the coordinates and/or charges to arrive */
614 MPI_Waitall(messages, pme_pp->req, pme_pp->stat);
617 while (!(cnb.flags & (PP_PME_COORD | PP_PME_FINISH)));
623 *chargeA = pme_pp->chargeA;
624 *chargeB = pme_pp->chargeB;
627 *sigmaA = pme_pp->sigmaA;
628 *sigmaB = pme_pp->sigmaB;
632 return ((cnb.flags & PP_PME_FINISH) ? pmerecvqxFINISH : pmerecvqxX);
635 static void receive_virial_energy(t_commrec *cr,
636 matrix vir_q, real *energy_q,
637 matrix vir_lj, real *energy_lj,
638 real *dvdlambda_q, real *dvdlambda_lj,
641 gmx_pme_comm_vir_ene_t cve;
643 if (cr->dd->pme_receive_vir_ener)
648 "PP node %d receiving from PME node %d: virial and energy\n",
649 cr->sim_nodeid, cr->dd->pme_nodeid);
652 MPI_Recv(&cve, sizeof(cve), MPI_BYTE, cr->dd->pme_nodeid, 1, cr->mpi_comm_mysim,
655 memset(&cve, 0, sizeof(cve));
658 m_add(vir_q, cve.vir_q, vir_q);
659 m_add(vir_lj, cve.vir_lj, vir_lj);
660 *energy_q = cve.energy_q;
661 *energy_lj = cve.energy_lj;
662 *dvdlambda_q += cve.dvdlambda_q;
663 *dvdlambda_lj += cve.dvdlambda_lj;
664 *pme_cycles = cve.cycles;
666 if (cve.stop_cond != gmx_stop_cond_none)
668 gmx_set_stop_condition(cve.stop_cond);
679 void gmx_pme_receive_f(t_commrec *cr,
680 rvec f[], matrix vir_q, real *energy_q,
681 matrix vir_lj, real *energy_lj,
682 real *dvdlambda_q, real *dvdlambda_lj,
687 #ifdef GMX_PME_DELAYED_WAIT
688 /* Wait for the x request to finish */
689 gmx_pme_send_params_coords_wait(cr->dd);
692 natoms = cr->dd->nat_home;
694 if (natoms > cr->dd->pme_recv_f_alloc)
696 cr->dd->pme_recv_f_alloc = over_alloc_dd(natoms);
697 srenew(cr->dd->pme_recv_f_buf, cr->dd->pme_recv_f_alloc);
701 MPI_Recv(cr->dd->pme_recv_f_buf[0],
702 natoms*sizeof(rvec), MPI_BYTE,
703 cr->dd->pme_nodeid, 0, cr->mpi_comm_mysim,
707 for (i = 0; i < natoms; i++)
709 rvec_inc(f[i], cr->dd->pme_recv_f_buf[i]);
713 receive_virial_energy(cr, vir_q, energy_q, vir_lj, energy_lj, dvdlambda_q, dvdlambda_lj, pme_cycles);
716 void gmx_pme_send_force_vir_ener(struct gmx_pme_pp *pme_pp,
718 matrix vir_q, real energy_q,
719 matrix vir_lj, real energy_lj,
720 real dvdlambda_q, real dvdlambda_lj,
723 gmx_pme_comm_vir_ene_t cve;
724 int messages, ind_start, ind_end, receiver;
728 /* Now the evaluated forces have to be transferred to the PP nodes */
731 for (receiver = 0; receiver < pme_pp->nnode; receiver++)
734 ind_end = ind_start + pme_pp->nat[receiver];
736 if (MPI_Isend(f[ind_start], (ind_end-ind_start)*sizeof(rvec), MPI_BYTE,
737 pme_pp->node[receiver], 0,
738 pme_pp->mpi_comm_mysim, &pme_pp->req[messages++]) != 0)
740 gmx_comm("MPI_Isend failed in do_pmeonly");
745 /* send virial and energy to our last PP node */
746 copy_mat(vir_q, cve.vir_q);
747 copy_mat(vir_lj, cve.vir_lj);
748 cve.energy_q = energy_q;
749 cve.energy_lj = energy_lj;
750 cve.dvdlambda_q = dvdlambda_q;
751 cve.dvdlambda_lj = dvdlambda_lj;
752 /* check for the signals to send back to a PP node */
753 cve.stop_cond = gmx_get_stop_condition();
759 fprintf(debug, "PME node sending to PP node %d: virial and energy\n",
763 MPI_Isend(&cve, sizeof(cve), MPI_BYTE,
764 pme_pp->node_peer, 1,
765 pme_pp->mpi_comm_mysim, &pme_pp->req[messages++]);
767 /* Wait for the forces to arrive */
768 MPI_Waitall(messages, pme_pp->req, pme_pp->stat);