Remove UB in ArrayRef from array
authorPaul Bauer <paul.bauer.q@gmail.com>
Wed, 5 May 2021 15:16:05 +0000 (17:16 +0200)
committerPaul Bauer <paul.bauer.q@gmail.com>
Thu, 6 May 2021 07:44:10 +0000 (09:44 +0200)
Several issues with ArrayRefs from nullptr have been found by UBSAN.

src/gromacs/domdec/mdsetup.cpp
src/gromacs/domdec/partition.cpp
src/gromacs/listed_forces/listed_forces.cpp
src/gromacs/mdlib/force.cpp
src/gromacs/mdlib/sim_util.cpp
src/gromacs/mdlib/tests/leapfrogtestrunners.cpp

index 53236dd3d8c948767a4c1f1eb5c95cbeebde29ee..ead9a61fd1c67f6c38f66624a9253ce24e5216e7 100644 (file)
@@ -148,8 +148,10 @@ void mdAlgorithmsSetupAtomData(const t_commrec*     cr,
         const int numPmeAtoms = numHomeAtoms - fr->n_tpi;
         gmx_pme_reinit_atoms(fr->pmedata,
                              numPmeAtoms,
-                             gmx::arrayRefFromArray(mdatoms->chargeA, mdatoms->nr),
-                             gmx::arrayRefFromArray(mdatoms->chargeB, mdatoms->nr));
+                             mdatoms->chargeA ? gmx::arrayRefFromArray(mdatoms->chargeA, mdatoms->nr)
+                                              : gmx::ArrayRef<real>{},
+                             mdatoms->chargeB ? gmx::arrayRefFromArray(mdatoms->chargeB, mdatoms->nr)
+                                              : gmx::ArrayRef<real>{});
     }
 
     if (constr)
index 6bd203df67f4adadb3428f6a02c3bfc0d78278dd..9a9409422c21accab930841751f0f0a07e4c4291 100644 (file)
@@ -3285,21 +3285,25 @@ void dd_partition_system(FILE*                     fplog,
     if (!thisRankHasDuty(cr, DUTY_PME))
     {
         /* Send the charges and/or c6/sigmas to our PME only node */
-        gmx_pme_send_parameters(cr,
-                                *fr->ic,
-                                mdatoms->nChargePerturbed != 0,
-                                mdatoms->nTypePerturbed != 0,
-                                gmx::arrayRefFromArray(mdatoms->chargeA, mdatoms->nr),
-                                mdatoms->chargeB ? gmx::arrayRefFromArray(mdatoms->chargeB, mdatoms->nr)
-                                                 : gmx::ArrayRef<real>{},
-                                gmx::arrayRefFromArray(mdatoms->sqrt_c6A, mdatoms->nr),
-                                mdatoms->sqrt_c6B ? gmx::arrayRefFromArray(mdatoms->sqrt_c6B, mdatoms->nr)
-                                                  : gmx::ArrayRef<real>{},
-                                gmx::arrayRefFromArray(mdatoms->sigmaA, mdatoms->nr),
-                                mdatoms->sigmaB ? gmx::arrayRefFromArray(mdatoms->sigmaB, mdatoms->nr)
-                                                : gmx::ArrayRef<real>{},
-                                dd_pme_maxshift_x(*dd),
-                                dd_pme_maxshift_y(*dd));
+        gmx_pme_send_parameters(
+                cr,
+                *fr->ic,
+                mdatoms->nChargePerturbed != 0,
+                mdatoms->nTypePerturbed != 0,
+                mdatoms->chargeA ? gmx::arrayRefFromArray(mdatoms->chargeA, mdatoms->nr)
+                                 : gmx::ArrayRef<real>{},
+                mdatoms->chargeB ? gmx::arrayRefFromArray(mdatoms->chargeB, mdatoms->nr)
+                                 : gmx::ArrayRef<real>{},
+                mdatoms->sqrt_c6A ? gmx::arrayRefFromArray(mdatoms->sqrt_c6A, mdatoms->nr)
+                                  : gmx::ArrayRef<real>{},
+                mdatoms->sqrt_c6B ? gmx::arrayRefFromArray(mdatoms->sqrt_c6B, mdatoms->nr)
+                                  : gmx::ArrayRef<real>{},
+                mdatoms->sigmaA ? gmx::arrayRefFromArray(mdatoms->sigmaA, mdatoms->nr)
+                                : gmx::ArrayRef<real>{},
+                mdatoms->sigmaB ? gmx::arrayRefFromArray(mdatoms->sigmaB, mdatoms->nr)
+                                : gmx::ArrayRef<real>{},
+                dd_pme_maxshift_x(*dd),
+                dd_pme_maxshift_y(*dd));
     }
 
     if (dd->atomSets != nullptr)
index be7ba91f0fa0f028e74f12180f992af538d8e3f0..85967f1b7f45caa70472ca266741b086630fc894 100644 (file)
@@ -512,8 +512,8 @@ real calc_one_bond(int                           thread,
                  pbc,
                  lambda.data(),
                  dvdl.data(),
-                 gmx::arrayRefFromArray(md->chargeA, md->nr),
-                 gmx::arrayRefFromArray(md->chargeB, md->nr),
+                 md->chargeA ? gmx::arrayRefFromArray(md->chargeA, md->nr) : gmx::ArrayRef<real>{},
+                 md->chargeB ? gmx::arrayRefFromArray(md->chargeB, md->nr) : gmx::ArrayRef<real>{},
                  md->bPerturbed ? gmx::arrayRefFromArray(md->bPerturbed, md->nr) : gmx::ArrayRef<bool>(),
                  gmx::arrayRefFromArray(md->cENER, md->nr),
                  md->nPerturbed,
index bfb1ca7b1d969ef004edb8f369e25fe25710eb33..a939841283a504e4231d7bf914b51128adb81464 100644 (file)
@@ -167,8 +167,10 @@ void calculateLongRangeNonbondeds(t_forcerec*                    fr,
                                 t,
                                 *fr,
                                 ir,
-                                gmx::constArrayRefFromArray(md->chargeA, md->nr),
-                                gmx::constArrayRefFromArray(md->chargeB, md->nr),
+                                md->chargeA ? gmx::constArrayRefFromArray(md->chargeA, md->nr)
+                                            : gmx::ArrayRef<const real>{},
+                                md->chargeB ? gmx::constArrayRefFromArray(md->chargeB, md->nr)
+                                            : gmx::ArrayRef<const real>{},
                                 (md->nChargePerturbed != 0),
                                 coordinates,
                                 box,
@@ -217,12 +219,18 @@ void calculateLongRangeNonbondeds(t_forcerec*                    fr,
                             fr->pmedata,
                             gmx::constArrayRefFromArray(coordinates.data(), md->homenr - fr->n_tpi),
                             forceWithVirial->force_,
-                            gmx::constArrayRefFromArray(md->chargeA, md->nr),
-                            gmx::constArrayRefFromArray(md->chargeB, md->nr),
-                            gmx::constArrayRefFromArray(md->sqrt_c6A, md->nr),
-                            gmx::constArrayRefFromArray(md->sqrt_c6B, md->nr),
-                            gmx::constArrayRefFromArray(md->sigmaA, md->nr),
-                            gmx::constArrayRefFromArray(md->sigmaB, md->nr),
+                            md->chargeA ? gmx::constArrayRefFromArray(md->chargeA, md->nr)
+                                        : gmx::ArrayRef<const real>{},
+                            md->chargeB ? gmx::constArrayRefFromArray(md->chargeB, md->nr)
+                                        : gmx::ArrayRef<const real>{},
+                            md->sqrt_c6A ? gmx::constArrayRefFromArray(md->sqrt_c6A, md->nr)
+                                         : gmx::ArrayRef<const real>{},
+                            md->sqrt_c6B ? gmx::constArrayRefFromArray(md->sqrt_c6B, md->nr)
+                                         : gmx::ArrayRef<const real>{},
+                            md->sigmaA ? gmx::constArrayRefFromArray(md->sigmaA, md->nr)
+                                       : gmx::ArrayRef<const real>{},
+                            md->sigmaB ? gmx::constArrayRefFromArray(md->sigmaB, md->nr)
+                                       : gmx::ArrayRef<const real>{},
                             box,
                             cr,
                             DOMAINDECOMP(cr) ? dd_pme_maxshift_x(*cr->dd) : 0,
index 342c08bb625d865c9eb287baee7ec84e4279ae4c..2f82ddc6d347f618b519630d6ca72236cbe5596c 100644 (file)
@@ -1670,8 +1670,10 @@ void do_force(FILE*                               fplog,
         calc_mu(start,
                 mdatoms->homenr,
                 xRef,
-                gmx::arrayRefFromArray(mdatoms->chargeA, mdatoms->nr),
-                gmx::arrayRefFromArray(mdatoms->chargeB, mdatoms->nr),
+                mdatoms->chargeA ? gmx::arrayRefFromArray(mdatoms->chargeA, mdatoms->nr)
+                                 : gmx::ArrayRef<real>{},
+                mdatoms->chargeB ? gmx::arrayRefFromArray(mdatoms->chargeB, mdatoms->nr)
+                                 : gmx::ArrayRef<real>{},
                 mdatoms->nChargePerturbed != 0,
                 dipoleData.muStaging[0],
                 dipoleData.muStaging[1]);
@@ -1763,35 +1765,45 @@ void do_force(FILE*                               fplog,
         /* Calculate the local and non-local free energy interactions here.
          * Happens here on the CPU both with and without GPU.
          */
-        nbv->dispatchFreeEnergyKernel(InteractionLocality::Local,
-                                      *fr,
-                                      x.unpaddedArrayRef(),
-                                      &forceOutNonbonded->forceWithShiftForces(),
-                                      gmx::arrayRefFromArray(mdatoms->chargeA, mdatoms->nr),
-                                      gmx::arrayRefFromArray(mdatoms->chargeB, mdatoms->nr),
-                                      gmx::arrayRefFromArray(mdatoms->typeA, mdatoms->nr),
-                                      gmx::arrayRefFromArray(mdatoms->typeB, mdatoms->nr),
-                                      inputrec.fepvals.get(),
-                                      lambda,
-                                      enerd,
-                                      stepWork,
-                                      nrnb);
+        nbv->dispatchFreeEnergyKernel(
+                InteractionLocality::Local,
+                *fr,
+                x.unpaddedArrayRef(),
+                &forceOutNonbonded->forceWithShiftForces(),
+                mdatoms->chargeA ? gmx::arrayRefFromArray(mdatoms->chargeA, mdatoms->nr)
+                                 : gmx::ArrayRef<real>{},
+                mdatoms->chargeB ? gmx::arrayRefFromArray(mdatoms->chargeB, mdatoms->nr)
+                                 : gmx::ArrayRef<real>{},
+                mdatoms->typeA ? gmx::arrayRefFromArray(mdatoms->typeA, mdatoms->nr)
+                               : gmx::ArrayRef<int>{},
+                mdatoms->typeB ? gmx::arrayRefFromArray(mdatoms->typeB, mdatoms->nr)
+                               : gmx::ArrayRef<int>{},
+                inputrec.fepvals.get(),
+                lambda,
+                enerd,
+                stepWork,
+                nrnb);
 
         if (havePPDomainDecomposition(cr))
         {
-            nbv->dispatchFreeEnergyKernel(InteractionLocality::NonLocal,
-                                          *fr,
-                                          x.unpaddedArrayRef(),
-                                          &forceOutNonbonded->forceWithShiftForces(),
-                                          gmx::arrayRefFromArray(mdatoms->chargeA, mdatoms->nr),
-                                          gmx::arrayRefFromArray(mdatoms->chargeB, mdatoms->nr),
-                                          gmx::arrayRefFromArray(mdatoms->typeA, mdatoms->nr),
-                                          gmx::arrayRefFromArray(mdatoms->typeB, mdatoms->nr),
-                                          inputrec.fepvals.get(),
-                                          lambda,
-                                          enerd,
-                                          stepWork,
-                                          nrnb);
+            nbv->dispatchFreeEnergyKernel(
+                    InteractionLocality::NonLocal,
+                    *fr,
+                    x.unpaddedArrayRef(),
+                    &forceOutNonbonded->forceWithShiftForces(),
+                    mdatoms->chargeA ? gmx::arrayRefFromArray(mdatoms->chargeA, mdatoms->nr)
+                                     : gmx::ArrayRef<real>{},
+                    mdatoms->chargeB ? gmx::arrayRefFromArray(mdatoms->chargeB, mdatoms->nr)
+                                     : gmx::ArrayRef<real>{},
+                    mdatoms->typeA ? gmx::arrayRefFromArray(mdatoms->typeA, mdatoms->nr)
+                                   : gmx::ArrayRef<int>{},
+                    mdatoms->typeB ? gmx::arrayRefFromArray(mdatoms->typeB, mdatoms->nr)
+                                   : gmx::ArrayRef<int>{},
+                    inputrec.fepvals.get(),
+                    lambda,
+                    enerd,
+                    stepWork,
+                    nrnb);
         }
     }
 
@@ -1841,9 +1853,12 @@ void do_force(FILE*                               fplog,
         real dvdl_walls = do_walls(inputrec,
                                    *fr,
                                    box,
-                                   gmx::arrayRefFromArray(mdatoms->typeA, mdatoms->nr),
-                                   gmx::arrayRefFromArray(mdatoms->typeB, mdatoms->nr),
-                                   gmx::arrayRefFromArray(mdatoms->cENER, mdatoms->nr),
+                                   mdatoms->typeA ? gmx::arrayRefFromArray(mdatoms->typeA, mdatoms->nr)
+                                                  : gmx::ArrayRef<int>{},
+                                   mdatoms->typeB ? gmx::arrayRefFromArray(mdatoms->typeB, mdatoms->nr)
+                                                  : gmx::ArrayRef<int>{},
+                                   mdatoms->cENER ? gmx::arrayRefFromArray(mdatoms->cENER, mdatoms->nr)
+                                                  : gmx::ArrayRef<unsigned short>{},
                                    mdatoms->homenr,
                                    mdatoms->nPerturbed,
                                    x.unpaddedConstArrayRef(),
index 5badaad2725b988524ab51f93864d8dfa4a01b95..a82d9fc9e797674f32c1b0ec3ee091b9d9ab8285 100644 (file)
@@ -76,9 +76,15 @@ void LeapFrogHostTestRunner::integrate(LeapFrogTestData* testData, int numSteps)
                 step,
                 testData->mdAtoms_.homenr,
                 testData->mdAtoms_.havePartiallyFrozenAtoms,
-                gmx::arrayRefFromArray(testData->mdAtoms_.ptype, testData->mdAtoms_.nr),
-                gmx::arrayRefFromArray(testData->mdAtoms_.invmass, testData->mdAtoms_.nr),
-                gmx::arrayRefFromArray(testData->mdAtoms_.invMassPerDim, testData->mdAtoms_.nr),
+                testData->mdAtoms_.ptype
+                        ? gmx::arrayRefFromArray(testData->mdAtoms_.ptype, testData->mdAtoms_.nr)
+                        : gmx::ArrayRef<ParticleType>{},
+                testData->mdAtoms_.invmass
+                        ? gmx::arrayRefFromArray(testData->mdAtoms_.invmass, testData->mdAtoms_.nr)
+                        : gmx::ArrayRef<real>{},
+                testData->mdAtoms_.invMassPerDim ? gmx::arrayRefFromArray(testData->mdAtoms_.invMassPerDim,
+                                                                          testData->mdAtoms_.nr)
+                                                 : gmx::ArrayRef<rvec>{},
                 &testData->state_,
                 testData->f_,
                 testData->forceCalculationData_,