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.
46 #include "gromacs/legacyheaders/domdec.h"
47 #include "gromacs/legacyheaders/network.h"
48 #include "gromacs/legacyheaders/pme.h"
49 #include "gromacs/legacyheaders/sighandler.h"
50 #include "gromacs/legacyheaders/typedefs.h"
51 #include "gromacs/legacyheaders/types/commrec.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/gmxmpi.h"
55 #include "gromacs/utility/smalloc.h"
58 eCommType_ChargeA, eCommType_ChargeB, eCommType_SQRTC6A, eCommType_SQRTC6B,
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_SQRTC6 (1<<2)
73 #define PP_PME_SQRTC6B (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 {
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_coeffs_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_coeffs_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,
184 gmx_pme_comm_n_box_t *cnb;
192 fprintf(debug, "PP rank %d sending to PME rank %d: %d%s%s%s%s\n",
193 cr->sim_nodeid, dd->pme_nodeid, n,
194 flags & PP_PME_CHARGE ? " charges" : "",
195 flags & PP_PME_SQRTC6 ? " sqrtC6" : "",
196 flags & PP_PME_SIGMA ? " sigma" : "",
197 flags & PP_PME_COORD ? " coordinates" : "");
200 #ifdef GMX_PME_DELAYED_WAIT
201 /* When can not use cnb until pending communication has finished */
202 gmx_pme_send_coeffs_coords_wait(dd);
205 if (dd->pme_receive_vir_ener)
207 /* Peer PP node: communicate all data */
216 cnb->maxshift_x = maxshift_x;
217 cnb->maxshift_y = maxshift_y;
218 cnb->lambda_q = lambda_q;
219 cnb->lambda_lj = lambda_lj;
221 if (flags & PP_PME_COORD)
223 copy_mat(box, cnb->box);
226 MPI_Isend(cnb, sizeof(*cnb), MPI_BYTE,
227 dd->pme_nodeid, eCommType_CNB, cr->mpi_comm_mysim,
228 &dd->req_pme[dd->nreq_pme++]);
231 else if (flags & (PP_PME_CHARGE | PP_PME_SQRTC6 | PP_PME_SIGMA))
234 /* Communicate only the number of atoms */
235 MPI_Isend(&n, sizeof(n), MPI_BYTE,
236 dd->pme_nodeid, eCommType_CNB, cr->mpi_comm_mysim,
237 &dd->req_pme[dd->nreq_pme++]);
244 if (flags & PP_PME_CHARGE)
246 MPI_Isend(chargeA, n*sizeof(real), MPI_BYTE,
247 dd->pme_nodeid, eCommType_ChargeA, cr->mpi_comm_mysim,
248 &dd->req_pme[dd->nreq_pme++]);
250 if (flags & PP_PME_CHARGEB)
252 MPI_Isend(chargeB, n*sizeof(real), MPI_BYTE,
253 dd->pme_nodeid, eCommType_ChargeB, cr->mpi_comm_mysim,
254 &dd->req_pme[dd->nreq_pme++]);
256 if (flags & PP_PME_SQRTC6)
258 MPI_Isend(c6A, n*sizeof(real), MPI_BYTE,
259 dd->pme_nodeid, eCommType_SQRTC6A, cr->mpi_comm_mysim,
260 &dd->req_pme[dd->nreq_pme++]);
262 if (flags & PP_PME_SQRTC6B)
264 MPI_Isend(c6B, n*sizeof(real), MPI_BYTE,
265 dd->pme_nodeid, eCommType_SQRTC6B, cr->mpi_comm_mysim,
266 &dd->req_pme[dd->nreq_pme++]);
268 if (flags & PP_PME_SIGMA)
270 MPI_Isend(sigmaA, n*sizeof(real), MPI_BYTE,
271 dd->pme_nodeid, eCommType_SigmaA, cr->mpi_comm_mysim,
272 &dd->req_pme[dd->nreq_pme++]);
274 if (flags & PP_PME_SIGMAB)
276 MPI_Isend(sigmaB, n*sizeof(real), MPI_BYTE,
277 dd->pme_nodeid, eCommType_SigmaB, cr->mpi_comm_mysim,
278 &dd->req_pme[dd->nreq_pme++]);
280 if (flags & PP_PME_COORD)
282 MPI_Isend(x[0], n*sizeof(rvec), MPI_BYTE,
283 dd->pme_nodeid, eCommType_COORD, cr->mpi_comm_mysim,
284 &dd->req_pme[dd->nreq_pme++]);
288 #ifndef GMX_PME_DELAYED_WAIT
289 /* Wait for the data to arrive */
290 /* We can skip this wait as we are sure x and q will not be modified
291 * before the next call to gmx_pme_send_x_q or gmx_pme_receive_f.
293 gmx_pme_send_coeffs_coords_wait(dd);
298 void gmx_pme_send_parameters(t_commrec *cr,
299 const interaction_const_t *ic,
300 gmx_bool bFreeEnergy_q, gmx_bool bFreeEnergy_lj,
301 real *chargeA, real *chargeB,
302 real *sqrt_c6A, real *sqrt_c6B,
303 real *sigmaA, real *sigmaB,
304 int maxshift_x, int maxshift_y)
309 if (EEL_PME(ic->eeltype))
311 flags |= PP_PME_CHARGE;
313 if (EVDW_PME(ic->vdwtype))
315 flags |= (PP_PME_SQRTC6 | PP_PME_SIGMA);
317 if (bFreeEnergy_q || bFreeEnergy_lj)
319 /* Assumes that the B state flags are in the bits just above
320 * the ones for the A state. */
321 flags |= (flags << 1);
324 gmx_pme_send_coeffs_coords(cr, flags,
326 sqrt_c6A, sqrt_c6B, sigmaA, sigmaB,
327 NULL, NULL, 0, 0, maxshift_x, maxshift_y, -1);
330 void gmx_pme_send_coordinates(t_commrec *cr, matrix box, rvec *x,
331 gmx_bool bFreeEnergy_q, gmx_bool bFreeEnergy_lj,
332 real lambda_q, real lambda_lj,
333 gmx_bool bEnerVir, int pme_flags,
338 flags = pme_flags | PP_PME_COORD;
341 flags |= PP_PME_FEP_Q;
345 flags |= PP_PME_FEP_LJ;
349 flags |= PP_PME_ENER_VIR;
351 gmx_pme_send_coeffs_coords(cr, flags, NULL, NULL, NULL, NULL, NULL, NULL,
352 box, x, lambda_q, lambda_lj, 0, 0, step);
355 void gmx_pme_send_finish(t_commrec *cr)
359 flags = PP_PME_FINISH;
361 gmx_pme_send_coeffs_coords(cr, flags, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, 0, 0, 0, 0, -1);
364 void gmx_pme_send_switchgrid(t_commrec gmx_unused *cr,
365 ivec gmx_unused grid_size,
366 real gmx_unused ewaldcoeff_q,
367 real gmx_unused ewaldcoeff_lj)
370 gmx_pme_comm_n_box_t cnb;
372 /* Only let one PP node signal each PME node */
373 if (cr->dd->pme_receive_vir_ener)
375 cnb.flags = PP_PME_SWITCHGRID;
376 copy_ivec(grid_size, cnb.grid_size);
377 cnb.ewaldcoeff_q = ewaldcoeff_q;
378 cnb.ewaldcoeff_lj = ewaldcoeff_lj;
380 /* We send this, uncommon, message blocking to simplify the code */
381 MPI_Send(&cnb, sizeof(cnb), MPI_BYTE,
382 cr->dd->pme_nodeid, eCommType_CNB, cr->mpi_comm_mysim);
387 void gmx_pme_send_resetcounters(t_commrec gmx_unused *cr, gmx_int64_t gmx_unused step)
390 gmx_pme_comm_n_box_t cnb;
392 /* Only let one PP node signal each PME node */
393 if (cr->dd->pme_receive_vir_ener)
395 cnb.flags = PP_PME_RESETCOUNTERS;
398 /* We send this, uncommon, message blocking to simplify the code */
399 MPI_Send(&cnb, sizeof(cnb), MPI_BYTE,
400 cr->dd->pme_nodeid, eCommType_CNB, cr->mpi_comm_mysim);
405 int gmx_pme_recv_coeffs_coords(struct gmx_pme_pp *pme_pp,
413 matrix gmx_unused box,
416 int gmx_unused *maxshift_x,
417 int gmx_unused *maxshift_y,
418 gmx_bool gmx_unused *bFreeEnergy_q,
419 gmx_bool gmx_unused *bFreeEnergy_lj,
420 real gmx_unused *lambda_q,
421 real gmx_unused *lambda_lj,
422 gmx_bool gmx_unused *bEnerVir,
424 gmx_int64_t gmx_unused *step,
425 ivec gmx_unused grid_size,
426 real gmx_unused *ewaldcoeff_q,
427 real gmx_unused *ewaldcoeff_lj)
429 gmx_pme_comm_n_box_t cnb;
430 int nat = 0, q, messages, sender;
435 /* avoid compiler warning about unused variable without MPI support */
441 /* Receive the send count, box and time step from the peer PP node */
442 MPI_Recv(&cnb, sizeof(cnb), MPI_BYTE,
443 pme_pp->node_peer, eCommType_CNB,
444 pme_pp->mpi_comm_mysim, MPI_STATUS_IGNORE);
448 fprintf(debug, "PME only rank receiving:%s%s%s%s%s\n",
449 (cnb.flags & PP_PME_CHARGE) ? " charges" : "",
450 (cnb.flags & PP_PME_COORD ) ? " coordinates" : "",
451 (cnb.flags & PP_PME_FINISH) ? " finish" : "",
452 (cnb.flags & PP_PME_SWITCHGRID) ? " switch grid" : "",
453 (cnb.flags & PP_PME_RESETCOUNTERS) ? " reset counters" : "");
456 if (cnb.flags & PP_PME_SWITCHGRID)
458 /* Special case, receive the new parameters and return */
459 copy_ivec(cnb.grid_size, grid_size);
460 *ewaldcoeff_q = cnb.ewaldcoeff_q;
461 *ewaldcoeff_lj = cnb.ewaldcoeff_lj;
462 return pmerecvqxSWITCHGRID;
465 if (cnb.flags & PP_PME_RESETCOUNTERS)
467 /* Special case, receive the step and return */
470 return pmerecvqxRESETCOUNTERS;
473 if (cnb.flags & (PP_PME_CHARGE | PP_PME_SQRTC6 | PP_PME_SIGMA))
475 /* Receive the send counts from the other PP nodes */
476 for (sender = 0; sender < pme_pp->nnode; sender++)
478 if (pme_pp->node[sender] == pme_pp->node_peer)
480 pme_pp->nat[sender] = cnb.natoms;
484 MPI_Irecv(&(pme_pp->nat[sender]), sizeof(pme_pp->nat[0]),
486 pme_pp->node[sender], eCommType_CNB,
487 pme_pp->mpi_comm_mysim, &pme_pp->req[messages++]);
490 MPI_Waitall(messages, pme_pp->req, pme_pp->stat);
494 for (sender = 0; sender < pme_pp->nnode; sender++)
496 nat += pme_pp->nat[sender];
499 if (nat > pme_pp->nalloc)
501 pme_pp->nalloc = over_alloc_dd(nat);
502 if (cnb.flags & PP_PME_CHARGE)
504 srenew(pme_pp->chargeA, pme_pp->nalloc);
506 if (cnb.flags & PP_PME_CHARGEB)
508 srenew(pme_pp->chargeB, pme_pp->nalloc);
510 if (cnb.flags & PP_PME_SQRTC6)
512 srenew(pme_pp->sqrt_c6A, pme_pp->nalloc);
514 if (cnb.flags & PP_PME_SQRTC6B)
516 srenew(pme_pp->sqrt_c6B, pme_pp->nalloc);
518 if (cnb.flags & PP_PME_SIGMA)
520 srenew(pme_pp->sigmaA, pme_pp->nalloc);
522 if (cnb.flags & PP_PME_SIGMAB)
524 srenew(pme_pp->sigmaB, pme_pp->nalloc);
526 srenew(pme_pp->x, pme_pp->nalloc);
527 srenew(pme_pp->f, pme_pp->nalloc);
530 /* maxshift is sent when the charges are sent */
531 *maxshift_x = cnb.maxshift_x;
532 *maxshift_y = cnb.maxshift_y;
534 /* Receive the charges in place */
535 for (q = 0; q < eCommType_NR; q++)
537 if (!(cnb.flags & (PP_PME_CHARGE<<q)))
543 case eCommType_ChargeA: charge_pp = pme_pp->chargeA; break;
544 case eCommType_ChargeB: charge_pp = pme_pp->chargeB; break;
545 case eCommType_SQRTC6A: charge_pp = pme_pp->sqrt_c6A; break;
546 case eCommType_SQRTC6B: charge_pp = pme_pp->sqrt_c6B; break;
547 case eCommType_SigmaA: charge_pp = pme_pp->sigmaA; break;
548 case eCommType_SigmaB: charge_pp = pme_pp->sigmaB; break;
549 default: gmx_incons("Wrong eCommType");
552 for (sender = 0; sender < pme_pp->nnode; sender++)
554 if (pme_pp->nat[sender] > 0)
556 MPI_Irecv(charge_pp+nat,
557 pme_pp->nat[sender]*sizeof(real),
559 pme_pp->node[sender], q,
560 pme_pp->mpi_comm_mysim,
561 &pme_pp->req[messages++]);
562 nat += pme_pp->nat[sender];
565 fprintf(debug, "Received from PP rank %d: %d %s\n",
566 pme_pp->node[sender], pme_pp->nat[sender],
567 (q == eCommType_ChargeA ||
568 q == eCommType_ChargeB) ? "charges" : "params");
574 pme_pp->flags_charge = cnb.flags;
577 if (cnb.flags & PP_PME_COORD)
579 if (!(pme_pp->flags_charge & (PP_PME_CHARGE | PP_PME_SQRTC6)))
581 gmx_incons("PME-only rank received coordinates before charges and/or C6-values"
585 /* The box, FE flag and lambda are sent along with the coordinates
587 copy_mat(cnb.box, box);
588 *bFreeEnergy_q = ((cnb.flags & GMX_PME_DO_COULOMB) &&
589 (cnb.flags & PP_PME_FEP_Q));
590 *bFreeEnergy_lj = ((cnb.flags & GMX_PME_DO_LJ) &&
591 (cnb.flags & PP_PME_FEP_LJ));
592 *lambda_q = cnb.lambda_q;
593 *lambda_lj = cnb.lambda_lj;
594 *bEnerVir = (cnb.flags & PP_PME_ENER_VIR);
595 *pme_flags = cnb.flags;
597 if (*bFreeEnergy_q && !(pme_pp->flags_charge & PP_PME_CHARGEB))
599 gmx_incons("PME-only rank received free energy request, but "
600 "did not receive B-state charges");
603 if (*bFreeEnergy_lj && !(pme_pp->flags_charge & PP_PME_SQRTC6B))
605 gmx_incons("PME-only rank received free energy request, but "
606 "did not receive B-state C6-values");
609 /* Receive the coordinates in place */
611 for (sender = 0; sender < pme_pp->nnode; sender++)
613 if (pme_pp->nat[sender] > 0)
615 MPI_Irecv(pme_pp->x[nat], pme_pp->nat[sender]*sizeof(rvec),
617 pme_pp->node[sender], eCommType_COORD,
618 pme_pp->mpi_comm_mysim, &pme_pp->req[messages++]);
619 nat += pme_pp->nat[sender];
622 fprintf(debug, "Received from PP rank %d: %d "
624 pme_pp->node[sender], pme_pp->nat[sender]);
630 /* Wait for the coordinates and/or charges to arrive */
631 MPI_Waitall(messages, pme_pp->req, pme_pp->stat);
634 while (!(cnb.flags & (PP_PME_COORD | PP_PME_FINISH)));
640 *chargeA = pme_pp->chargeA;
641 *chargeB = pme_pp->chargeB;
642 *sqrt_c6A = pme_pp->sqrt_c6A;
643 *sqrt_c6B = pme_pp->sqrt_c6B;
644 *sigmaA = pme_pp->sigmaA;
645 *sigmaB = pme_pp->sigmaB;
649 return ((cnb.flags & PP_PME_FINISH) ? pmerecvqxFINISH : pmerecvqxX);
652 static void receive_virial_energy(t_commrec *cr,
653 matrix vir_q, real *energy_q,
654 matrix vir_lj, real *energy_lj,
655 real *dvdlambda_q, real *dvdlambda_lj,
658 gmx_pme_comm_vir_ene_t cve;
660 if (cr->dd->pme_receive_vir_ener)
665 "PP rank %d receiving from PME rank %d: virial and energy\n",
666 cr->sim_nodeid, cr->dd->pme_nodeid);
669 MPI_Recv(&cve, sizeof(cve), MPI_BYTE, cr->dd->pme_nodeid, 1, cr->mpi_comm_mysim,
672 memset(&cve, 0, sizeof(cve));
675 m_add(vir_q, cve.vir_q, vir_q);
676 m_add(vir_lj, cve.vir_lj, vir_lj);
677 *energy_q = cve.energy_q;
678 *energy_lj = cve.energy_lj;
679 *dvdlambda_q += cve.dvdlambda_q;
680 *dvdlambda_lj += cve.dvdlambda_lj;
681 *pme_cycles = cve.cycles;
683 if (cve.stop_cond != gmx_stop_cond_none)
685 gmx_set_stop_condition(cve.stop_cond);
696 void gmx_pme_receive_f(t_commrec *cr,
697 rvec f[], matrix vir_q, real *energy_q,
698 matrix vir_lj, real *energy_lj,
699 real *dvdlambda_q, real *dvdlambda_lj,
704 #ifdef GMX_PME_DELAYED_WAIT
705 /* Wait for the x request to finish */
706 gmx_pme_send_coeffs_coords_wait(cr->dd);
709 natoms = cr->dd->nat_home;
711 if (natoms > cr->dd->pme_recv_f_alloc)
713 cr->dd->pme_recv_f_alloc = over_alloc_dd(natoms);
714 srenew(cr->dd->pme_recv_f_buf, cr->dd->pme_recv_f_alloc);
718 MPI_Recv(cr->dd->pme_recv_f_buf[0],
719 natoms*sizeof(rvec), MPI_BYTE,
720 cr->dd->pme_nodeid, 0, cr->mpi_comm_mysim,
724 for (i = 0; i < natoms; i++)
726 rvec_inc(f[i], cr->dd->pme_recv_f_buf[i]);
730 receive_virial_energy(cr, vir_q, energy_q, vir_lj, energy_lj, dvdlambda_q, dvdlambda_lj, pme_cycles);
733 void gmx_pme_send_force_vir_ener(struct gmx_pme_pp *pme_pp,
735 matrix vir_q, real energy_q,
736 matrix vir_lj, real energy_lj,
737 real dvdlambda_q, real dvdlambda_lj,
740 gmx_pme_comm_vir_ene_t cve;
741 int messages, ind_start, ind_end, receiver;
745 /* Now the evaluated forces have to be transferred to the PP nodes */
748 for (receiver = 0; receiver < pme_pp->nnode; receiver++)
751 ind_end = ind_start + pme_pp->nat[receiver];
753 if (MPI_Isend(f[ind_start], (ind_end-ind_start)*sizeof(rvec), MPI_BYTE,
754 pme_pp->node[receiver], 0,
755 pme_pp->mpi_comm_mysim, &pme_pp->req[messages++]) != 0)
757 gmx_comm("MPI_Isend failed in do_pmeonly");
762 /* send virial and energy to our last PP node */
763 copy_mat(vir_q, cve.vir_q);
764 copy_mat(vir_lj, cve.vir_lj);
765 cve.energy_q = energy_q;
766 cve.energy_lj = energy_lj;
767 cve.dvdlambda_q = dvdlambda_q;
768 cve.dvdlambda_lj = dvdlambda_lj;
769 /* check for the signals to send back to a PP node */
770 cve.stop_cond = gmx_get_stop_condition();
776 fprintf(debug, "PME rank sending to PP rank %d: virial and energy\n",
780 MPI_Isend(&cve, sizeof(cve), MPI_BYTE,
781 pme_pp->node_peer, 1,
782 pme_pp->mpi_comm_mysim, &pme_pp->req[messages++]);
784 /* Wait for the forces to arrive */
785 MPI_Waitall(messages, pme_pp->req, pme_pp->stat);