NetSci
Loading...
Searching...
No Matches
Public Member Functions | Private Attributes | List of all members
Network Class Reference

Public Member Functions

 Network ()
 Default constructor for Network.
 
 ~Network ()
 Destructor for Network.
 
void init (const std::string &trajectoryFile, const std::string &topologyFile, int firstFrame, int lastFrame, int stride=1)
 Initialize the Network with trajectory and topology files.
 
int numNodes () const
 Get the number of nodes in the Network.
 
CuArray< float > * nodeCoordinates ()
 Get the node coordinates as a CuArray.
 
std::vector< Node * > & nodes ()
 Get a reference to the vector of nodes in the Network.
 
int numFrames () const
 Get the number of frames in the Network.
 
NodenodeFromAtomIndex (int atomIndex)
 Get the node corresponding to the Atom with the given index.
 
Atomsatoms () const
 Get the Atoms object associated with the Network.
 
void parsePdb (const std::string &fname)
 Parse a PDB file to populate the Network.
 
void parseDcd (const std::string &nodeCoordinates, int firstFrame, int lastFrame, int stride)
 Parse a DCD file to populate the Network.
 
void save (const std::string &jsonFile)
 Save the Network as a JSON file.
 
void load (const std::string &jsonFile)
 Load a Network from a JSON file.
 
void nodeCoordinates (const std::string &nodeCoordinatesFile)
 Set the node coordinates from a file.
 

Private Attributes

std::vector< Node * > nodeAtomIndexVector_
 
std::vector< Node * > nodes_
 
int numNodes_
 
int numFrames_
 
CuArray< float > * nodeCoordinates_
 
Atomsatoms_
 

Constructor & Destructor Documentation

◆ Network()

Network::Network ( )

Default constructor for Network.

Constructs an empty Network object.

Member Function Documentation

◆ atoms()

Atoms * Network::atoms ( ) const

Get the Atoms object associated with the Network.

Returns a pointer to the Atoms object associated with the Network.

Returns
A pointer to the Atoms object.

Python Example

from netchem import Network, data_files
"""
1000 frame Molecular Dynamics trajectory of a pyrophosphatase system.
The dcd is a binary coordinate file, and the pdb is a single frame
topology file.
"""
dcd, pdb = data_files('pyro')
"""
Construct a network using the entire trajectory.
If you want to use only a subset of the trajectory, you can specify
the first frame, last frame, and stride.
"""
network = Network()
network.init(
trajectoryFile=str(dcd), # Convert the dcd Path object to a string
topologyFile=str(pdb), # Convert the pdb Path object to a string
firstFrame=0,
lastFrame=999,
stride=1
)
"""
The Network atoms method returns an Atoms object, which is a container
class for Atom objects. The reason Atom objects are stored in a custom
class instead of a build in container is to keep minimize cache misses,
as the Atoms class ensures that all Atom objects are stored contiguously.
It also allows for more control over memory management. The Atoms class
is iterable and individual Atom objects can be accessed via indexing. For
more information about the Atoms class, see the Atoms class examples and
documentation.
"""
atoms = network.atoms()
"""
Iterate over atoms and print the index and name of each atom.
"""
for atom in atoms:
print(atom.index(), atom.name())
Definition network.h:12

◆ init()

void Network::init ( const std::string &  trajectoryFile,
const std::string &  topologyFile,
int  firstFrame,
int  lastFrame,
int  stride = 1 
)

Initialize the Network with trajectory and topology files.

Initializes the Network by loading trajectory and topology files.

Parameters
trajectoryFilePath to the trajectory file.
topologyFilePath to the topology file.
firstFrameIndex of the first frame to consider.
lastFrameIndex of the last frame to consider.
strideStride between frames.

Python Example

from netchem import Network, data_files
"""
1000 frame Molecular Dynamics trajectory of a pyrophosphatase system.
The dcd is a binary coordinate file, and the pdb is a single frame
topology file.
"""
dcd, pdb = data_files('pyro')
"""
Construct a network using the entire trajectory.
If you want to use only a subset of the trajectory, you can specify
the first frame, last frame, and stride.
"""
network = Network()
network.init(
trajectoryFile=str(dcd), # Convert the dcd Path object to a string
topologyFile=str(pdb), # Convert the pdb Path object to a string
firstFrame=0,
lastFrame=999,
stride=1
)

◆ load()

void Network::load ( const std::string &  jsonFile)

Load a Network from a JSON file.

Loads a Network from the specified JSON file.

Parameters
jsonFilePath to the JSON file.

Python Example

from pathlib import Path
from netchem import Network, data_files, Node
"""
If operating on a CuArray, you must import the
appropriate CuArray class. In this case, we operate on
nodeCoordinates, which is a FloatCuArray.
"""
from cuarray import FloatCuArray
"""
1000 frame Molecular Dynamics trajectory of a pyrophosphatase system.
The dcd is a binary coordinate file, and the pdb is a single frame
topology file.
"""
dcd, pdb = data_files('pyro')
"""
Construct a network using the entire trajectory.
If you want to use only a subset of the trajectory, you can specify
the first frame, last frame, and stride.
"""
network = Network()
network.init(
trajectoryFile=str(dcd), # Convert the dcd Path object to a string
topologyFile=str(pdb), # Convert the pdb Path object to a string
firstFrame=0,
lastFrame=999,
stride=1
)
"""
Network objects can be serialized and saved to a JSON file
using the save method. The save method takes a single argument,
which is the absolute path to the JSON file that the Network
object will be saved to. If using the pathlib module,
make sure to convert the Path object to a string prior
to calling the save method. Another import note is that
the save method does not automatically save the node coordinates,
which have to be saved separately using the CuArray
save method.
"""
cwd = Path.cwd()
network_json = str(cwd / 'network.json')
network_node_coordinates_npy = str(cwd / 'network_node_coordinates.npy')
network.save(network_json)
"""Save the network node coordinates to a .npy file."""
network.nodeCoordinates().save(network_node_coordinates_npy)
"""
The load method can be used to deserialize, or load a Network object from
a Network serialization JSON file. Like the save method,
the load method takes a single argument, which is the absolute path
to the JSON file that the Network object will be loaded from.
The network's node coordinates must be loaded separate from the
.npy file they were saved to using the Network nodeCoordinates method.
"""
network_copy = Network()
network_copy.load(str(network_json))
"""Load the network node coordinates from the serialized coordinates .npy file."""
network_copy.nodeCoordinates(network_node_coordinates_npy)
"""
Print the two network Node objects side by side to verify that they are the same.
"""
for node_copy, node in zip(network_copy.nodes(), network.nodes()):
print(node)
print(node_copy)
for frame in range(network_copy.numFrames()):
x_idx = frame
y_idx = frame + network_copy.numFrames()
z_idx = frame + 2*network_copy.numFrames()
x = network_copy.nodeCoordinates()[node.index()][x_idx]
y = network_copy.nodeCoordinates()[node.index()][y_idx]
z = network_copy.nodeCoordinates()[node.index()][z_idx]
print(f'Node {node.index()} in frame {frame} coordinates: ({x}, {y}, {z})', end=' ')
x = network.nodeCoordinates()[node.index()][x_idx]
y = network.nodeCoordinates()[node.index()][y_idx]
z = network.nodeCoordinates()[node.index()][z_idx]
print(f'Node {node.index()} in frame {frame} coordinates: ({x}, {y}, {z})')
print("-"*80)

◆ nodeCoordinates() [1/2]

CuArray< float > * Network::nodeCoordinates ( )

Get the node coordinates as a CuArray.

Returns a pointer to the CuArray object containing the node coordinates.

Returns
A pointer to the CuArray containing the node coordinates.

Python Example

from netchem import Network, data_files
"""
If operating on a CuArray, you must import the
appropriate CuArray class. In this case, we operate on
nodeCoordinates, which is a FloatCuArray.
"""
from cuarray import FloatCuArray
"""
1000 frame Molecular Dynamics trajectory of a pyrophosphatase system.
The dcd is a binary coordinate file, and the pdb is a single frame
topology file.
"""
dcd, pdb = data_files('pyro')
"""
Construct a network using the entire trajectory.
If you want to use only a subset of the trajectory, you can specify
the first frame, last frame, and stride.
"""
network = Network()
network.init(
trajectoryFile=str(dcd), # Convert the dcd Path object to a string
topologyFile=str(pdb), # Convert the pdb Path object to a string
firstFrame=0,
lastFrame=999,
stride=1
)
"""
The nodeCoordinates method returns a Mx3*N float CuArray,
where M is the number of nodes in the network and N is the number
of frames. For each node, the coordinates are stored in a linear,
row-major format, so to access the x, y, and z coordinates of the
i-th node in the j-th frame, you would use the following syntax:
x = nodeCoordinates[i][j]
y = nodeCoordinates[i][j+N]
z = nodeCoordinates[i][j+2*N]
"""
nodeCoordinates = network.nodeCoordinates()
"""
Get the x, y, and z coordinates of node 5 in frame 50.
"""
frame_index = 50
node_index = 4
y_offset = network.numFrames()
z_offset = 2*network.numFrames()
x = nodeCoordinates[node_index][frame_index]
y = nodeCoordinates[node_index][frame_index+y_offset]
z = nodeCoordinates[node_index][frame_index+z_offset]
print(f'Node 5 in frame 50 coordinates: ({x}, {y}, {z})')
print(f'Node 5 in frame 50 coordinates: ({x}, {y}, {z})')

◆ nodeCoordinates() [2/2]

void Network::nodeCoordinates ( const std::string &  nodeCoordinatesFile)

Set the node coordinates from a file.

Sets the node coordinates from the specified node coordinates file.

Parameters
nodeCoordinatesFilePath to the node coordinates file.

Python Example

from pathlib import Path
from netchem import Network, data_files, Node
"""
If operating on a CuArray, you must import the
appropriate CuArray class. In this case, we operate on
nodeCoordinates, which is a FloatCuArray.
"""
from cuarray import FloatCuArray
"""
1000 frame Molecular Dynamics trajectory of a pyrophosphatase system.
The dcd is a binary coordinate file, and the pdb is a single frame
topology file.
"""
dcd, pdb = data_files('pyro')
"""
Construct a network using the entire trajectory.
If you want to use only a subset of the trajectory, you can specify
the first frame, last frame, and stride.
"""
network = Network()
network.init(
trajectoryFile=str(dcd), # Convert the dcd Path object to a string
topologyFile=str(pdb), # Convert the pdb Path object to a string
firstFrame=0,
lastFrame=999,
stride=1
)
"""Current working directory"""
cwd = Path.cwd()
"""Save the network node coordinates to a .npy file."""
network_node_coordinates_npy = str(cwd / 'network_node_coordinates.npy')
network.nodeCoordinates().save(network_node_coordinates_npy)
"""
Create a new Network object and load the serialized network node
coordinates from the .npy file.
"""
network_copy = Network()
network_copy.nodeCoordinates(network_node_coordinates_npy)
"""
Print the two network Node objects' coordinates side by side to verify that they are the same.
"""
for node in network.nodes():
for frame in range(network.numFrames()):
x_idx = frame
y_idx = frame + network_copy.numFrames()
z_idx = frame + 2*network_copy.numFrames()
x = network_copy.nodeCoordinates()[node.index()][x_idx]
y = network_copy.nodeCoordinates()[node.index()][y_idx]
z = network_copy.nodeCoordinates()[node.index()][z_idx]
print(f'Node {node.index()} in frame {frame} coordinates: ({x}, {y}, {z})', end=' ')
x = network.nodeCoordinates()[node.index()][x_idx]
y = network.nodeCoordinates()[node.index()][y_idx]
z = network.nodeCoordinates()[node.index()][z_idx]
print(f'Node {node.index()} in frame {frame} coordinates: ({x}, {y}, {z})')
print("-"*80)

◆ nodeFromAtomIndex()

Node * Network::nodeFromAtomIndex ( int  atomIndex)

Get the node corresponding to the Atom with the given index.

Returns a pointer to the Node object that the Atom with the specified index is part of

Parameters
atomIndexThe index of the Atom.
Returns
A pointer to the Node corresponding to the Atom index.

Python Example

from netchem import Network, data_files
"""
If operating on a CuArray, you must import the
appropriate CuArray class. In this case, we operate on
nodeCoordinates, which is a FloatCuArray.
"""
from cuarray import FloatCuArray
"""
1000 frame Molecular Dynamics trajectory of a pyrophosphatase system.
The dcd is a binary coordinate file, and the pdb is a single frame
topology file.
"""
dcd, pdb = data_files('pyro')
"""
Construct a network using the entire trajectory.
If you want to use only a subset of the trajectory, you can specify
the first frame, last frame, and stride.
"""
network = Network()
network.init(
trajectoryFile=str(dcd), # Convert the dcd Path object to a string
topologyFile=str(pdb), # Convert the pdb Path object to a string
firstFrame=0,
lastFrame=999,
stride=1
)
"""
Use the nodeFromAtomIndex method to get the node index of the atom.
The below example returns the Node object that the fifth atom (index 4) belongs to.
See the Node class documentation for information on the Node class properties
and methods.
"""
node = network.nodeFromAtomIndex(4)

◆ nodes()

std::vector< Node * > & Network::nodes ( )

Get a reference to the vector of nodes in the Network.

Returns a reference to the vector of nodes in the Network.

Returns
A reference to the vector of nodes.

Python Example

from netchem import Network, data_files
"""
1000 frame Molecular Dynamics trajectory of a pyrophosphatase system.
The dcd is a binary coordinate file, and the pdb is a single frame
topology file.
"""
dcd, pdb = data_files('pyro')
"""
Construct a network using the entire trajectory.
If you want to use only a subset of the trajectory, you can specify
the first frame, last frame, and stride.
"""
network = Network()
network.init(
trajectoryFile=str(dcd), # Convert the dcd Path object to a string
topologyFile=str(pdb), # Convert the pdb Path object to a string
firstFrame=0,
lastFrame=999,
stride=1
)
"""
Use the Network nodes method to get a NodeVector object, which
is essentially a list of Node objects, and can be converted to a
Python list via list(NodeVector). The NodeVector is iterable and
individual Node objects can be accessed via indexing.
"""
nodes = network.nodes()
"""
Iterate over nodes and print the index and number of
atoms in each of each node. See the Node class examples
and documentation for more information about Node methods.
"""
for node in nodes:
print(node.index(), node.numAtoms())

◆ numFrames()

int Network::numFrames ( ) const

Get the number of frames in the Network.

Returns the number of frames in the Network.

Returns
The number of frames.

Python Example

from netchem import Network, data_files
"""
1000 frame Molecular Dynamics trajectory of a pyrophosphatase system.
The dcd is a binary coordinate file, and the pdb is a single frame
topology file.
"""
dcd, pdb = data_files('pyro')
"""
Construct a network using the entire trajectory.
If you want to use only a subset of the trajectory, you can specify
the first frame, last frame, and stride.
"""
network = Network()
network.init(
trajectoryFile=str(dcd), # Convert the dcd Path object to a string
topologyFile=str(pdb), # Convert the pdb Path object to a string
firstFrame=0,
lastFrame=999,
stride=1
)
"""
Use the numFrames method to get the number of frames in the network's
Molecular Dynamics trajectory.
"""
numFrames = network.numFrames()
print('Number of frames in network trajectory: {}'.format(numFrames))

◆ numNodes()

int Network::numNodes ( ) const

Get the number of nodes in the Network.

Returns the number of nodes in the Network.

Returns
The number of nodes.

Python Example

from netchem import Network, data_files
"""
1000 frame Molecular Dynamics trajectory of a pyrophosphatase system.
The dcd is a binary coordinate file, and the pdb is a single frame
topology file.
"""
dcd, pdb = data_files('pyro')
"""
Construct a network using the entire trajectory.
If you want to use only a subset of the trajectory, you can specify
the first frame, last frame, and stride.
"""
network = Network()
network.init(
trajectoryFile=str(dcd), # Convert the dcd Path object to a string
topologyFile=str(pdb), # Convert the pdb Path object to a string
firstFrame=0,
lastFrame=999,
stride=1
)
"""
Use the numNodes method to get the number of nodes in the network.
"""
numNodes = network.numNodes()
print('Number of nodes: {}'.format(numNodes))

◆ parseDcd()

void Network::parseDcd ( const std::string &  nodeCoordinates,
int  firstFrame,
int  lastFrame,
int  stride 
)

Parse a DCD file to populate the Network.

Parses the specified DCD file to populate the Network with node coordinates.

Parameters
fnamePath to the node coordinates file.
firstFrameIndex of the first frame to consider.
lastFrameIndex of the last frame to consider.
strideStride between frames.

◆ parsePdb()

void Network::parsePdb ( const std::string &  fname)

Parse a PDB file to populate the Network.

Parses the specified PDB file to populate the Network with Atom and Node objects.

Parameters
fnamePath to the PDB file.

◆ save()

void Network::save ( const std::string &  jsonFile)

Save the Network as a JSON file.

Saves the Network as a JSON file.

Parameters
jsonFilePath to the JSON file.

Python Example

from pathlib import Path
from netchem import Network, data_files, Node
"""
If operating on a CuArray, you must import the
appropriate CuArray class. In this case, we operate on
nodeCoordinates, which is a FloatCuArray.
"""
from cuarray import FloatCuArray
"""
1000 frame Molecular Dynamics trajectory of a pyrophosphatase system.
The dcd is a binary coordinate file, and the pdb is a single frame
topology file.
"""
dcd, pdb = data_files('pyro')
"""
Construct a network using the entire trajectory.
If you want to use only a subset of the trajectory, you can specify
the first frame, last frame, and stride.
"""
network = Network()
network.init(
trajectoryFile=str(dcd), # Convert the dcd Path object to a string
topologyFile=str(pdb), # Convert the pdb Path object to a string
firstFrame=0,
lastFrame=999,
stride=1
)
"""
Network objects can be serialized and saved to a JSON file
using the save method. The save method takes a single argument,
which is the absolute path to the JSON file that the Network
object will be saved to. If using the pathlib module,
make sure to convert the Path object to a string prior
to calling the save method. Another import note is that
the save method does not automatically save the node coordinates,
which have to be saved separately using the CuArray
save method.
"""
cwd = Path.cwd()
network_json = str(cwd / 'network.json')
network_node_coordinates_npy = str(cwd / 'network_node_coordinates.npy')
network.save(network_json)
"""Save the network node coordinates to a .npy file."""
network.nodeCoordinates().save(network_node_coordinates_npy)

The documentation for this class was generated from the following file: