From: Carsten Kutzner Date: Wed, 24 Aug 2011 15:20:38 +0000 (+0200) Subject: Instead of a segv, mdrun now gives an error msg if npme>np, fixes #795 X-Git-Url: http://biod.pnpi.spb.ru/gitweb/?a=commitdiff_plain;h=e230816fcb44aff25cb35b7c67ddeb67a7ca42c6;p=alexxy%2Fgromacs.git Instead of a segv, mdrun now gives an error msg if npme>np, fixes #795 Instead of a segv, mdrun now gives an error msg if npme>np, fixes #795 Change-Id: I2e93bccf45ace215d6237354559b8750628a13d6 --- diff --git a/src/mdlib/domdec.c b/src/mdlib/domdec.c index 82989fdf23..a16e548e6f 100644 --- a/src/mdlib/domdec.c +++ b/src/mdlib/domdec.c @@ -6459,7 +6459,7 @@ gmx_domdec_t *init_domain_decomposition(FILE *fplog,t_commrec *cr, if (cr->npmenodes > dd->nnodes) { gmx_fatal_collective(FARGS,cr,NULL, - "The number of separate PME node (%d) is larger than the number of PP nodes (%d), this is not supported.",cr->npmenodes,dd->nnodes); + "The number of separate PME nodes (%d) is larger than the number of PP nodes (%d), this is not supported.",cr->npmenodes,dd->nnodes); } if (cr->npmenodes > 0) { diff --git a/src/mdlib/domdec_setup.c b/src/mdlib/domdec_setup.c index 8048e72710..ed095f52e3 100644 --- a/src/mdlib/domdec_setup.c +++ b/src/mdlib/domdec_setup.c @@ -37,6 +37,11 @@ static int factorize(int n,int **fac,int **mfac) { int d,ndiv; + if (n <= 0) + { + gmx_fatal(FARGS, "Can only factorize positive integers."); + } + /* Decompose n in factors */ snew(*fac,n/2); snew(*mfac,n/2); @@ -671,6 +676,12 @@ real dd_choose_grid(FILE *fplog, gmx_fatal(FARGS, "Can not have separate PME nodes with 2 or less nodes"); } + if (cr->npmenodes >= cr->nnodes) + { + gmx_fatal(FARGS, + "Can not have %d separate PME nodes with just %d total nodes", + cr->npmenodes, cr->nnodes); + } /* If the user purposely selected the number of PME nodes, * only check for large primes in the PP node count.