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,2014, 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.
47 #include "types/commrec.h"
48 #include "gromacs/utility/smalloc.h"
49 #include "gromacs/utility/fatalerror.h"
54 #include "sighandler.h"
56 #include "gromacs/utility/gmxmpi.h"
59 eCommType_ChargeA, eCommType_ChargeB, eCommType_SQRTC6A, eCommType_SQRTC6B,
60 eCommType_SigmaA, eCommType_SigmaB, eCommType_NR, eCommType_COORD,
64 /* Some parts of the code(gmx_pme_send_q, gmx_pme_recv_q_x) assume
65 * that the six first flags are exactly in this order.
66 * If more PP_PME_...-flags are to be introduced be aware of some of
67 * the PME-specific flags in pme.h. Currently, they are also passed
71 #define PP_PME_CHARGE (1<<0)
72 #define PP_PME_CHARGEB (1<<1)
73 #define PP_PME_SQRTC6 (1<<2)
74 #define PP_PME_SQRTC6B (1<<3)
75 #define PP_PME_SIGMA (1<<4)
76 #define PP_PME_SIGMAB (1<<5)
77 #define PP_PME_COORD (1<<6)
78 #define PP_PME_FEP_Q (1<<7)
79 #define PP_PME_FEP_LJ (1<<8)
80 #define PP_PME_ENER_VIR (1<<9)
81 #define PP_PME_FINISH (1<<10)
82 #define PP_PME_SWITCHGRID (1<<11)
83 #define PP_PME_RESETCOUNTERS (1<<12)
85 #define PME_PP_SIGSTOP (1<<0)
86 #define PME_PP_SIGSTOPNSS (1<<1)
88 typedef struct gmx_pme_pp {
90 MPI_Comm mpi_comm_mysim;
92 int nnode; /* The number of PP node to communicate with */
93 int *node; /* The PP node ranks */
94 int node_peer; /* The peer PP node rank */
95 int *nat; /* The number of atom for each PP node */
96 int flags_charge; /* The flags sent along with the last charges */
112 typedef struct gmx_pme_comm_n_box {
121 ivec grid_size; /* For PME grid tuning */
122 real ewaldcoeff_q; /* For PME grid tuning */
124 } gmx_pme_comm_n_box_t;
134 gmx_stop_cond_t stop_cond;
135 } gmx_pme_comm_vir_ene_t;
140 gmx_pme_pp_t gmx_pme_pp_init(t_commrec gmx_unused *cr)
142 struct gmx_pme_pp *pme_pp;
148 pme_pp->mpi_comm_mysim = cr->mpi_comm_mysim;
149 MPI_Comm_rank(cr->mpi_comm_mygroup, &rank);
150 get_pme_ddnodes(cr, rank, &pme_pp->nnode, &pme_pp->node, &pme_pp->node_peer);
151 snew(pme_pp->nat, pme_pp->nnode);
152 snew(pme_pp->req, eCommType_NR*pme_pp->nnode);
153 snew(pme_pp->stat, eCommType_NR*pme_pp->nnode);
155 pme_pp->flags_charge = 0;
161 /* This should be faster with a real non-blocking MPI implementation */
162 /* #define GMX_PME_DELAYED_WAIT */
164 static void gmx_pme_send_coeffs_coords_wait(gmx_domdec_t gmx_unused *dd)
169 MPI_Waitall(dd->nreq_pme, dd->req_pme, MPI_STATUSES_IGNORE);
175 static void gmx_pme_send_coeffs_coords(t_commrec *cr, int flags,
176 real gmx_unused *chargeA, real gmx_unused *chargeB,
177 real gmx_unused *c6A, real gmx_unused *c6B,
178 real gmx_unused *sigmaA, real gmx_unused *sigmaB,
179 matrix box, rvec gmx_unused *x,
180 real lambda_q, real lambda_lj,
181 int maxshift_x, int maxshift_y,
185 gmx_pme_comm_n_box_t *cnb;
193 fprintf(debug, "PP node %d sending to PME node %d: %d%s%s\n",
194 cr->sim_nodeid, dd->pme_nodeid, n,
195 flags & PP_PME_CHARGE ? " charges" : "",
196 flags & PP_PME_COORD ? " coordinates" : "");
199 #ifdef GMX_PME_DELAYED_WAIT
200 /* When can not use cnb until pending communication has finished */
201 gmx_pme_send_coeffs_coords_wait(dd);
204 if (dd->pme_receive_vir_ener)
206 /* Peer PP node: communicate all data */
215 cnb->maxshift_x = maxshift_x;
216 cnb->maxshift_y = maxshift_y;
217 cnb->lambda_q = lambda_q;
218 cnb->lambda_lj = lambda_lj;
220 if (flags & PP_PME_COORD)
222 copy_mat(box, cnb->box);
225 MPI_Isend(cnb, sizeof(*cnb), MPI_BYTE,
226 dd->pme_nodeid, eCommType_CNB, cr->mpi_comm_mysim,
227 &dd->req_pme[dd->nreq_pme++]);
230 else if (flags & (PP_PME_CHARGE | PP_PME_SQRTC6 | PP_PME_SIGMA))
233 /* Communicate only the number of atoms */
234 MPI_Isend(&n, sizeof(n), MPI_BYTE,
235 dd->pme_nodeid, eCommType_CNB, cr->mpi_comm_mysim,
236 &dd->req_pme[dd->nreq_pme++]);
243 if (flags & PP_PME_CHARGE)
245 MPI_Isend(chargeA, n*sizeof(real), MPI_BYTE,
246 dd->pme_nodeid, eCommType_ChargeA, cr->mpi_comm_mysim,
247 &dd->req_pme[dd->nreq_pme++]);
249 if (flags & PP_PME_CHARGEB)
251 MPI_Isend(chargeB, n*sizeof(real), MPI_BYTE,
252 dd->pme_nodeid, eCommType_ChargeB, cr->mpi_comm_mysim,
253 &dd->req_pme[dd->nreq_pme++]);
255 if (flags & PP_PME_SQRTC6)
257 MPI_Isend(c6A, n*sizeof(real), MPI_BYTE,
258 dd->pme_nodeid, eCommType_SQRTC6A, cr->mpi_comm_mysim,
259 &dd->req_pme[dd->nreq_pme++]);
261 if (flags & PP_PME_SQRTC6B)
263 MPI_Isend(c6B, n*sizeof(real), MPI_BYTE,
264 dd->pme_nodeid, eCommType_SQRTC6B, cr->mpi_comm_mysim,
265 &dd->req_pme[dd->nreq_pme++]);
267 if (flags & PP_PME_SIGMA)
269 MPI_Isend(sigmaA, n*sizeof(real), MPI_BYTE,
270 dd->pme_nodeid, eCommType_SigmaA, cr->mpi_comm_mysim,
271 &dd->req_pme[dd->nreq_pme++]);
273 if (flags & PP_PME_SIGMAB)
275 MPI_Isend(sigmaB, n*sizeof(real), MPI_BYTE,
276 dd->pme_nodeid, eCommType_SigmaB, cr->mpi_comm_mysim,
277 &dd->req_pme[dd->nreq_pme++]);
279 if (flags & PP_PME_COORD)
281 MPI_Isend(x[0], n*sizeof(rvec), MPI_BYTE,
282 dd->pme_nodeid, eCommType_COORD, cr->mpi_comm_mysim,
283 &dd->req_pme[dd->nreq_pme++]);
287 #ifndef GMX_PME_DELAYED_WAIT
288 /* Wait for the data to arrive */
289 /* We can skip this wait as we are sure x and q will not be modified
290 * before the next call to gmx_pme_send_x_q or gmx_pme_receive_f.
292 gmx_pme_send_coeffs_coords_wait(dd);
297 void gmx_pme_send_parameters(t_commrec *cr,
298 gmx_bool bFreeEnergy_q, gmx_bool bFreeEnergy_lj,
299 real *chargeA, real *chargeB,
300 real *sqrt_c6A, real *sqrt_c6B,
301 real *sigmaA, real *sigmaB,
302 int maxshift_x, int maxshift_y)
306 /* We always send the charges, even with only LJ- and no Coulomb-PME */
307 flags = PP_PME_CHARGE;
308 if (sqrt_c6A != NULL)
310 flags |= PP_PME_SQRTC6;
314 flags |= PP_PME_SIGMA;
316 if (bFreeEnergy_q || bFreeEnergy_lj)
318 /* Assumes that the B state flags are in the bits just above
319 * the ones for the A state. */
320 flags |= (flags << 1);
323 gmx_pme_send_coeffs_coords(cr, flags,
325 sqrt_c6A, sqrt_c6B, sigmaA, sigmaB,
326 NULL, NULL, 0, 0, maxshift_x, maxshift_y, -1);
329 void gmx_pme_send_coordinates(t_commrec *cr, matrix box, rvec *x,
330 gmx_bool bFreeEnergy_q, gmx_bool bFreeEnergy_lj,
331 real lambda_q, real lambda_lj,
332 gmx_bool bEnerVir, int pme_flags,
337 flags = pme_flags | PP_PME_COORD;
340 flags |= PP_PME_FEP_Q;
344 flags |= PP_PME_FEP_LJ;
348 flags |= PP_PME_ENER_VIR;
350 gmx_pme_send_coeffs_coords(cr, flags, NULL, NULL, NULL, NULL, NULL, NULL,
351 box, x, lambda_q, lambda_lj, 0, 0, step);
354 void gmx_pme_send_finish(t_commrec *cr)
358 flags = PP_PME_FINISH;
360 gmx_pme_send_coeffs_coords(cr, flags, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, 0, 0, 0, 0, -1);
363 void gmx_pme_send_switchgrid(t_commrec gmx_unused *cr,
364 ivec gmx_unused grid_size,
365 real gmx_unused ewaldcoeff_q,
366 real gmx_unused ewaldcoeff_lj)
369 gmx_pme_comm_n_box_t cnb;
371 /* Only let one PP node signal each PME node */
372 if (cr->dd->pme_receive_vir_ener)
374 cnb.flags = PP_PME_SWITCHGRID;
375 copy_ivec(grid_size, cnb.grid_size);
376 cnb.ewaldcoeff_q = ewaldcoeff_q;
377 cnb.ewaldcoeff_lj = ewaldcoeff_lj;
379 /* We send this, uncommon, message blocking to simplify the code */
380 MPI_Send(&cnb, sizeof(cnb), MPI_BYTE,
381 cr->dd->pme_nodeid, eCommType_CNB, cr->mpi_comm_mysim);
386 void gmx_pme_send_resetcounters(t_commrec gmx_unused *cr, gmx_int64_t gmx_unused step)
389 gmx_pme_comm_n_box_t cnb;
391 /* Only let one PP node signal each PME node */
392 if (cr->dd->pme_receive_vir_ener)
394 cnb.flags = PP_PME_RESETCOUNTERS;
397 /* We send this, uncommon, message blocking to simplify the code */
398 MPI_Send(&cnb, sizeof(cnb), MPI_BYTE,
399 cr->dd->pme_nodeid, eCommType_CNB, cr->mpi_comm_mysim);
404 int gmx_pme_recv_coeffs_coords(struct gmx_pme_pp *pme_pp,
412 matrix gmx_unused box,
415 int gmx_unused *maxshift_x,
416 int gmx_unused *maxshift_y,
417 gmx_bool gmx_unused *bFreeEnergy_q,
418 gmx_bool gmx_unused *bFreeEnergy_lj,
419 real gmx_unused *lambda_q,
420 real gmx_unused *lambda_lj,
421 gmx_bool gmx_unused *bEnerVir,
423 gmx_int64_t gmx_unused *step,
424 ivec gmx_unused grid_size,
425 real gmx_unused *ewaldcoeff_q,
426 real gmx_unused *ewaldcoeff_lj)
428 gmx_pme_comm_n_box_t cnb;
429 int nat = 0, q, messages, sender;
434 /* avoid compiler warning about unused variable without MPI support */
440 /* Receive the send count, box and time step from the peer PP node */
441 MPI_Recv(&cnb, sizeof(cnb), MPI_BYTE,
442 pme_pp->node_peer, eCommType_CNB,
443 pme_pp->mpi_comm_mysim, MPI_STATUS_IGNORE);
447 fprintf(debug, "PME only node receiving:%s%s%s%s%s\n",
448 (cnb.flags & PP_PME_CHARGE) ? " charges" : "",
449 (cnb.flags & PP_PME_COORD ) ? " coordinates" : "",
450 (cnb.flags & PP_PME_FINISH) ? " finish" : "",
451 (cnb.flags & PP_PME_SWITCHGRID) ? " switch grid" : "",
452 (cnb.flags & PP_PME_RESETCOUNTERS) ? " reset counters" : "");
455 if (cnb.flags & PP_PME_SWITCHGRID)
457 /* Special case, receive the new parameters and return */
458 copy_ivec(cnb.grid_size, grid_size);
459 *ewaldcoeff_q = cnb.ewaldcoeff_q;
460 *ewaldcoeff_lj = cnb.ewaldcoeff_lj;
461 return pmerecvqxSWITCHGRID;
464 if (cnb.flags & PP_PME_RESETCOUNTERS)
466 /* Special case, receive the step and return */
469 return pmerecvqxRESETCOUNTERS;
472 if (cnb.flags & (PP_PME_CHARGE | PP_PME_SQRTC6 | PP_PME_SIGMA))
474 /* Receive the send counts from the other PP nodes */
475 for (sender = 0; sender < pme_pp->nnode; sender++)
477 if (pme_pp->node[sender] == pme_pp->node_peer)
479 pme_pp->nat[sender] = cnb.natoms;
483 MPI_Irecv(&(pme_pp->nat[sender]), sizeof(pme_pp->nat[0]),
485 pme_pp->node[sender], eCommType_CNB,
486 pme_pp->mpi_comm_mysim, &pme_pp->req[messages++]);
489 MPI_Waitall(messages, pme_pp->req, pme_pp->stat);
493 for (sender = 0; sender < pme_pp->nnode; sender++)
495 nat += pme_pp->nat[sender];
498 if (nat > pme_pp->nalloc)
500 pme_pp->nalloc = over_alloc_dd(nat);
501 if (cnb.flags & PP_PME_CHARGE)
503 srenew(pme_pp->chargeA, pme_pp->nalloc);
505 if (cnb.flags & PP_PME_CHARGEB)
507 srenew(pme_pp->chargeB, pme_pp->nalloc);
509 if (cnb.flags & PP_PME_SQRTC6)
511 srenew(pme_pp->sqrt_c6A, pme_pp->nalloc);
513 if (cnb.flags & PP_PME_SQRTC6B)
515 srenew(pme_pp->sqrt_c6B, pme_pp->nalloc);
517 if (cnb.flags & PP_PME_SIGMA)
519 srenew(pme_pp->sigmaA, pme_pp->nalloc);
521 if (cnb.flags & PP_PME_SIGMAB)
523 srenew(pme_pp->sigmaB, pme_pp->nalloc);
525 srenew(pme_pp->x, pme_pp->nalloc);
526 srenew(pme_pp->f, pme_pp->nalloc);
529 /* maxshift is sent when the charges are sent */
530 *maxshift_x = cnb.maxshift_x;
531 *maxshift_y = cnb.maxshift_y;
533 /* Receive the charges in place */
534 for (q = 0; q < eCommType_NR; q++)
536 if (!(cnb.flags & (PP_PME_CHARGE<<q)))
542 case eCommType_ChargeA: charge_pp = pme_pp->chargeA; break;
543 case eCommType_ChargeB: charge_pp = pme_pp->chargeB; break;
544 case eCommType_SQRTC6A: charge_pp = pme_pp->sqrt_c6A; break;
545 case eCommType_SQRTC6B: charge_pp = pme_pp->sqrt_c6B; break;
546 case eCommType_SigmaA: charge_pp = pme_pp->sigmaA; break;
547 case eCommType_SigmaB: charge_pp = pme_pp->sigmaB; break;
548 default: gmx_incons("Wrong eCommType");
551 for (sender = 0; sender < pme_pp->nnode; sender++)
553 if (pme_pp->nat[sender] > 0)
555 MPI_Irecv(charge_pp+nat,
556 pme_pp->nat[sender]*sizeof(real),
558 pme_pp->node[sender], q,
559 pme_pp->mpi_comm_mysim,
560 &pme_pp->req[messages++]);
561 nat += pme_pp->nat[sender];
564 fprintf(debug, "Received from PP node %d: %d "
566 pme_pp->node[sender], pme_pp->nat[sender]);
572 pme_pp->flags_charge = cnb.flags;
575 if (cnb.flags & PP_PME_COORD)
577 if (!(pme_pp->flags_charge & (PP_PME_CHARGE | PP_PME_SQRTC6)))
579 gmx_incons("PME-only node received coordinates before charges and/or C6-values"
583 /* The box, FE flag and lambda are sent along with the coordinates
585 copy_mat(cnb.box, box);
586 *bFreeEnergy_q = ((cnb.flags & GMX_PME_DO_COULOMB) &&
587 (cnb.flags & PP_PME_FEP_Q));
588 *bFreeEnergy_lj = ((cnb.flags & GMX_PME_DO_LJ) &&
589 (cnb.flags & PP_PME_FEP_LJ));
590 *lambda_q = cnb.lambda_q;
591 *lambda_lj = cnb.lambda_lj;
592 *bEnerVir = (cnb.flags & PP_PME_ENER_VIR);
593 *pme_flags = cnb.flags;
595 if (*bFreeEnergy_q && !(pme_pp->flags_charge & PP_PME_CHARGEB))
597 gmx_incons("PME-only node received free energy request, but "
598 "did not receive B-state charges");
601 if (*bFreeEnergy_lj && !(pme_pp->flags_charge & PP_PME_SQRTC6B))
603 gmx_incons("PME-only node received free energy request, but "
604 "did not receive B-state C6-values");
607 /* Receive the coordinates in place */
609 for (sender = 0; sender < pme_pp->nnode; sender++)
611 if (pme_pp->nat[sender] > 0)
613 MPI_Irecv(pme_pp->x[nat], pme_pp->nat[sender]*sizeof(rvec),
615 pme_pp->node[sender], eCommType_COORD,
616 pme_pp->mpi_comm_mysim, &pme_pp->req[messages++]);
617 nat += pme_pp->nat[sender];
620 fprintf(debug, "Received from PP node %d: %d "
622 pme_pp->node[sender], pme_pp->nat[sender]);
628 /* Wait for the coordinates and/or charges to arrive */
629 MPI_Waitall(messages, pme_pp->req, pme_pp->stat);
632 while (!(cnb.flags & (PP_PME_COORD | PP_PME_FINISH)));
638 *chargeA = pme_pp->chargeA;
639 *chargeB = pme_pp->chargeB;
640 *sqrt_c6A = pme_pp->sqrt_c6A;
641 *sqrt_c6B = pme_pp->sqrt_c6B;
642 *sigmaA = pme_pp->sigmaA;
643 *sigmaB = pme_pp->sigmaB;
647 return ((cnb.flags & PP_PME_FINISH) ? pmerecvqxFINISH : pmerecvqxX);
650 static void receive_virial_energy(t_commrec *cr,
651 matrix vir_q, real *energy_q,
652 matrix vir_lj, real *energy_lj,
653 real *dvdlambda_q, real *dvdlambda_lj,
656 gmx_pme_comm_vir_ene_t cve;
658 if (cr->dd->pme_receive_vir_ener)
663 "PP node %d receiving from PME node %d: virial and energy\n",
664 cr->sim_nodeid, cr->dd->pme_nodeid);
667 MPI_Recv(&cve, sizeof(cve), MPI_BYTE, cr->dd->pme_nodeid, 1, cr->mpi_comm_mysim,
670 memset(&cve, 0, sizeof(cve));
673 m_add(vir_q, cve.vir_q, vir_q);
674 m_add(vir_lj, cve.vir_lj, vir_lj);
675 *energy_q = cve.energy_q;
676 *energy_lj = cve.energy_lj;
677 *dvdlambda_q += cve.dvdlambda_q;
678 *dvdlambda_lj += cve.dvdlambda_lj;
679 *pme_cycles = cve.cycles;
681 if (cve.stop_cond != gmx_stop_cond_none)
683 gmx_set_stop_condition(cve.stop_cond);
694 void gmx_pme_receive_f(t_commrec *cr,
695 rvec f[], matrix vir_q, real *energy_q,
696 matrix vir_lj, real *energy_lj,
697 real *dvdlambda_q, real *dvdlambda_lj,
702 #ifdef GMX_PME_DELAYED_WAIT
703 /* Wait for the x request to finish */
704 gmx_pme_send_coeffs_coords_wait(cr->dd);
707 natoms = cr->dd->nat_home;
709 if (natoms > cr->dd->pme_recv_f_alloc)
711 cr->dd->pme_recv_f_alloc = over_alloc_dd(natoms);
712 srenew(cr->dd->pme_recv_f_buf, cr->dd->pme_recv_f_alloc);
716 MPI_Recv(cr->dd->pme_recv_f_buf[0],
717 natoms*sizeof(rvec), MPI_BYTE,
718 cr->dd->pme_nodeid, 0, cr->mpi_comm_mysim,
722 for (i = 0; i < natoms; i++)
724 rvec_inc(f[i], cr->dd->pme_recv_f_buf[i]);
728 receive_virial_energy(cr, vir_q, energy_q, vir_lj, energy_lj, dvdlambda_q, dvdlambda_lj, pme_cycles);
731 void gmx_pme_send_force_vir_ener(struct gmx_pme_pp *pme_pp,
733 matrix vir_q, real energy_q,
734 matrix vir_lj, real energy_lj,
735 real dvdlambda_q, real dvdlambda_lj,
738 gmx_pme_comm_vir_ene_t cve;
739 int messages, ind_start, ind_end, receiver;
743 /* Now the evaluated forces have to be transferred to the PP nodes */
746 for (receiver = 0; receiver < pme_pp->nnode; receiver++)
749 ind_end = ind_start + pme_pp->nat[receiver];
751 if (MPI_Isend(f[ind_start], (ind_end-ind_start)*sizeof(rvec), MPI_BYTE,
752 pme_pp->node[receiver], 0,
753 pme_pp->mpi_comm_mysim, &pme_pp->req[messages++]) != 0)
755 gmx_comm("MPI_Isend failed in do_pmeonly");
760 /* send virial and energy to our last PP node */
761 copy_mat(vir_q, cve.vir_q);
762 copy_mat(vir_lj, cve.vir_lj);
763 cve.energy_q = energy_q;
764 cve.energy_lj = energy_lj;
765 cve.dvdlambda_q = dvdlambda_q;
766 cve.dvdlambda_lj = dvdlambda_lj;
767 /* check for the signals to send back to a PP node */
768 cve.stop_cond = gmx_get_stop_condition();
774 fprintf(debug, "PME node sending to PP node %d: virial and energy\n",
778 MPI_Isend(&cve, sizeof(cve), MPI_BYTE,
779 pme_pp->node_peer, 1,
780 pme_pp->mpi_comm_mysim, &pme_pp->req[messages++]);
782 /* Wait for the forces to arrive */
783 MPI_Waitall(messages, pme_pp->req, pme_pp->stat);