Skip to content
GitLab
Menu
Projects
Groups
Snippets
Loading...
Help
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Contribute to GitLab
Sign in
Toggle navigation
Menu
Open sidebar
tsoc
openmm
Commits
c29e203b
Commit
c29e203b
authored
Sep 08, 2022
by
peastman
Browse files
Ensure valid atom order after loading a checkpoint
parent
9d9ffb01
Changes
4
Show whitespace changes
Inline
Side-by-side
Showing
4 changed files
with
31 additions
and
3 deletions
+31
-3
platforms/common/include/openmm/common/ComputeContext.h
platforms/common/include/openmm/common/ComputeContext.h
+6
-0
platforms/common/src/ComputeContext.cpp
platforms/common/src/ComputeContext.cpp
+23
-3
platforms/cuda/src/CudaKernels.cpp
platforms/cuda/src/CudaKernels.cpp
+1
-0
platforms/opencl/src/OpenCLKernels.cpp
platforms/opencl/src/OpenCLKernels.cpp
+1
-0
No files found.
platforms/common/include/openmm/common/ComputeContext.h
View file @
c29e203b
...
...
@@ -496,6 +496,11 @@ public:
* definitions and order to be revalidated.
*/
bool
invalidateMolecules
(
ComputeForceInfo
*
force
);
/**
* Make sure the current atom order is valid, based on the forces. If not, perform reordering
* to generate a new valid order. This method is only needed in very unusual situations.
*/
void
validateAtomOrder
();
/**
* Wait until all work that has been queued (kernel executions, asynchronous data transfers, etc.)
* has been submitted to the device. This does not mean it has necessarily been completed.
...
...
@@ -508,6 +513,7 @@ protected:
struct
MoleculeGroup
;
class
VirtualSiteInfo
;
void
findMoleculeGroups
();
void
resetAtomOrder
();
/**
* This is the internal implementation of reorderAtoms(), templatized by the numerical precision in use.
*/
...
...
platforms/common/src/ComputeContext.cpp
View file @
c29e203b
...
...
@@ -376,6 +376,13 @@ bool ComputeContext::invalidateMolecules(ComputeForceInfo* force) {
// atoms to their original order, rebuild the list of identical molecules, and sort them
// again.
resetAtomOrder
();
findMoleculeGroups
();
reorderAtoms
();
return
true
;
}
void
ComputeContext
::
resetAtomOrder
()
{
ContextSelector
selector
(
*
this
);
vector
<
mm_int4
>
newCellOffsets
(
numAtoms
);
if
(
getUseDoublePrecision
())
{
...
...
@@ -436,12 +443,25 @@ bool ComputeContext::invalidateMolecules(ComputeForceInfo* force) {
posCellOffsets
[
i
]
=
newCellOffsets
[
i
];
}
getAtomIndexArray
().
upload
(
atomIndex
);
findMoleculeGroups
();
for
(
auto
listener
:
reorderListeners
)
listener
->
execute
();
forceNextReorder
=
true
;
}
void
ComputeContext
::
validateAtomOrder
()
{
for
(
auto
&
mol
:
moleculeGroups
)
{
for
(
int
atom
:
mol
.
atoms
)
{
set
<
int
>
identical
;
for
(
int
offset
:
mol
.
offsets
)
identical
.
insert
(
atom
+
offset
);
for
(
int
i
:
identical
)
if
(
identical
.
find
(
atomIndex
[
i
])
==
identical
.
end
())
{
resetAtomOrder
();
reorderAtoms
();
return
true
;
return
;
}
}
}
}
void
ComputeContext
::
forceReorder
()
{
...
...
platforms/cuda/src/CudaKernels.cpp
View file @
c29e203b
...
...
@@ -424,6 +424,7 @@ void CudaUpdateStateDataKernel::loadCheckpoint(ContextImpl& context, istream& st
SimTKOpenMMUtilities
::
loadCheckpoint
(
stream
);
for
(
auto
listener
:
cu
.
getReorderListeners
())
listener
->
execute
();
cu
.
validateAtomOrder
();
}
class
CudaCalcNonbondedForceKernel
::
ForceInfo
:
public
CudaForceInfo
{
...
...
platforms/opencl/src/OpenCLKernels.cpp
View file @
c29e203b
...
...
@@ -441,6 +441,7 @@ void OpenCLUpdateStateDataKernel::loadCheckpoint(ContextImpl& context, istream&
SimTKOpenMMUtilities
::
loadCheckpoint
(
stream
);
for
(
auto
listener
:
cl
.
getReorderListeners
())
listener
->
execute
();
cl
.
validateAtomOrder
();
}
class
OpenCLCalcNonbondedForceKernel
::
ForceInfo
:
public
OpenCLForceInfo
{
...
...
Write
Preview
Markdown
is supported
0%
Try again
or
attach a new file
.
Attach a file
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Cancel
Please
register
or
sign in
to comment