Fixed nullptr derefence in LINCS error
authorBerk Hess <hess@kth.se>
Wed, 12 Dec 2018 09:23:46 +0000 (10:23 +0100)
committerPaul Bauer <paul.bauer.q@gmail.com>
Fri, 14 Dec 2018 13:25:51 +0000 (14:25 +0100)
Somehow a dereferenced nullptr could be passed as a reference
to gmx_multisim_t to the Constraint factory function.
Changed the reference in the Constraint object to a pointer.

Fixed #2803

Change-Id: I4806069973067d27078a1324d18a406c7b3e227d

src/gromacs/mdlib/constr.cpp
src/gromacs/mdlib/constr.h
src/gromacs/mdlib/lincs.cpp
src/gromacs/mdlib/lincs.h
src/gromacs/mdrun/runner.cpp

index f1ec807362582d86c24f4d07275cf225af602575..0a20f8baaf0941ce8a5ab94ffc72e4afcbab6393 100644 (file)
@@ -104,7 +104,7 @@ class Constraints::Impl
              FILE                 *log_p,
              const t_mdatoms      &md_p,
              const t_commrec      *cr_p,
-             const gmx_multisim_t &ms,
+             const gmx_multisim_t *ms,
              t_nrnb               *nrnb,
              gmx_wallcycle        *wcycle_p,
              bool                  pbcHandlingRequired,
@@ -171,7 +171,7 @@ class Constraints::Impl
         //! Communication support.
         const t_commrec      *cr = nullptr;
         //! Multi-sim support.
-        const gmx_multisim_t &ms;
+        const gmx_multisim_t *ms = nullptr;
         /*!\brief Input options.
          *
          * \todo Replace with IMdpOptions */
@@ -958,7 +958,7 @@ Constraints::Constraints(const gmx_mtop_t     &mtop,
                          FILE                 *log,
                          const t_mdatoms      &md,
                          const t_commrec      *cr,
-                         const gmx_multisim_t &ms,
+                         const gmx_multisim_t *ms,
                          t_nrnb               *nrnb,
                          gmx_wallcycle        *wcycle,
                          bool                  pbcHandlingRequired,
@@ -983,7 +983,7 @@ Constraints::Impl::Impl(const gmx_mtop_t     &mtop_p,
                         FILE                 *log_p,
                         const t_mdatoms      &md_p,
                         const t_commrec      *cr_p,
-                        const gmx_multisim_t &ms_p,
+                        const gmx_multisim_t *ms_p,
                         t_nrnb               *nrnb_p,
                         gmx_wallcycle        *wcycle_p,
                         bool                  pbcHandlingRequired,
index 804f51bbdf24e3a9dbcaae4486356d8d44b503bc..ecb713ce14012353677716e763a1dcdeceabb651 100644 (file)
@@ -103,7 +103,7 @@ class Constraints
                     FILE                 *log,
                     const t_mdatoms      &md,
                     const t_commrec      *cr,
-                    const gmx_multisim_t &ms,
+                    const gmx_multisim_t *ms,
                     t_nrnb               *nrnb,
                     gmx_wallcycle        *wcycle,
                     bool                  pbcHandlingRequired,
index 543c1848dd2a2239f1848676de5e3aa04f22f474..e6e6a33b6667604a1395a5ec78a6501ead4ff6a1 100644 (file)
@@ -2414,7 +2414,7 @@ bool constrain_lincs(bool computeRmsd,
                      int64_t step,
                      Lincs *lincsd, const t_mdatoms &md,
                      const t_commrec *cr,
-                     const gmx_multisim_t &ms,
+                     const gmx_multisim_t *ms,
                      const rvec *x, rvec *xprime, rvec *min_proj,
                      matrix box, t_pbc *pbc,
                      real lambda, real *dvdlambda,
@@ -2560,9 +2560,9 @@ bool constrain_lincs(bool computeRmsd,
             {
                 cconerr(lincsd, xprime, pbc,
                         &ncons_loc, &p_ssd, &p_max, &p_imax);
-                if (isMultiSim(&ms))
+                if (isMultiSim(ms))
                 {
-                    sprintf(buf3, " in simulation %d", ms.sim);
+                    sprintf(buf3, " in simulation %d", ms->sim);
                 }
                 else
                 {
index fe172a4746fc4d9a1209c4b5839d499def5e1502..c04741e45fd10c53c205c02f158ecc4ae5e9b103 100644 (file)
@@ -97,7 +97,7 @@ constrain_lincs(bool computeRmsd,
                 int64_t step,
                 Lincs *lincsd, const t_mdatoms &md,
                 const t_commrec *cr,
-                const gmx_multisim_t &ms,
+                const gmx_multisim_t *ms,
                 const rvec *x, rvec *xprime, rvec *min_proj,
                 matrix box, t_pbc *pbc,
                 real lambda, real *dvdlambda,
index c507f477b0db5a3a0f7482d2a924e410d2639891..060c83f33a562d917c606e4fe8c26529d489320c 100644 (file)
@@ -1390,7 +1390,7 @@ int Mdrunner::mdrunner()
                                     || observablesHistory.edsamHistory);
         auto constr              = makeConstraints(mtop, *inputrec, doEssentialDynamics,
                                                    fplog, *mdAtoms->mdatoms(),
-                                                   cr, *ms, nrnb, wcycle, fr->bMolPBC);
+                                                   cr, ms, nrnb, wcycle, fr->bMolPBC);
 
         if (DOMAINDECOMP(cr))
         {