summaryrefslogtreecommitdiff
diff options
context:
space:
mode:
-rw-r--r--include/tng_io.h156
-rw-r--r--src/lib/tng_io.c137
2 files changed, 221 insertions, 72 deletions
diff --git a/include/tng_io.h b/include/tng_io.h
index c713fbd..16ba849 100644
--- a/include/tng_io.h
+++ b/include/tng_io.h
@@ -1616,66 +1616,102 @@ tng_function_status DECLSPECDLLEXPORT tng_molecule_system_copy(tng_trajectory_t
tng_trajectory_t tng_data_dest);
/**
- * @brief Get the chains of a molecule.
+ * @brief Get the number of chains in a molecule.
* @param tng_data is the trajectory containing the molecule.
- * @param molecule is the molecule of which to get the chains.
- * @param chains is pointing to the list of chains. This is a pointer
- * to the list of the molecule, which means that it should be handled
- * carefully, e.g. not freed.
+ * @param molecule is the molecule of which to get the number of chains.
* @param n is pointing to a value set to the number of chains.
* @pre \code molecule != 0 \endcode The molecule must not be NULL.
- * @pre \code chains != 0 \endcode The pointer to the list of chains must not
- * be a NULL pointer.
* @pre \code n != 0 \endcode The pointer to n must not be a NULL pointer.
* @return TNG_SUCCESS (0) if successful.
*/
-tng_function_status DECLSPECDLLEXPORT tng_molecule_chains_get
+tng_function_status DECLSPECDLLEXPORT tng_molecule_num_chains_get
(tng_trajectory_t tng_data,
tng_molecule_t molecule,
- tng_chain_t *chains,
int64_t *n);
/**
- * @brief Get the residues of a molecule.
+ * @brief Retrieve the chain of a molecule with specified index in the list
+ * of chains.
+ * @param tng_data is the trajectory data container containing the molecule.
+ * @param index is the index (in molecule->chains) of the chain to return
+ * @param molecule is the molecule from which to get the chain.
+ * @param chain is a pointer to the chain if it was found - otherwise 0.
+ * @pre \code molecule != 0 \endcode molecule must not be a NULL pointer.
+ * @pre \code chain != 0 \endcode chain must not be a NULL pointer.
+ * @return TNG_SUCCESS (0) if the chain is found or TNG_FAILURE (1) if the
+ * chain is not found.
+ */
+tng_function_status DECLSPECDLLEXPORT tng_molecule_chain_of_index_get
+ (tng_trajectory_t tng_data,
+ tng_molecule_t molecule,
+ int64_t index,
+ tng_chain_t *chain);
+
+/**
+ * @brief Get the number of residues in a molecule.
* @param tng_data is the trajectory containing the molecule.
- * @param molecule is the molecule of which to get the residues.
- * @param residues is pointing to the list of residues. This is a pointer
- * to the list of the molecule, which means that it should be handled
- * carefully, e.g. not freed.
+ * @param molecule is the molecule of which to get the number residues.
* @param n is pointing to a value set to the number of residues.
* @pre \code molecule != 0 \endcode The molecule must not be NULL.
- * @pre \code residues != 0 \endcode The pointer to the list of residues must not
- * be a NULL pointer.
* @pre \code n != 0 \endcode The pointer to n must not be a NULL pointer.
* @return TNG_SUCCESS (0) if successful.
*/
-tng_function_status DECLSPECDLLEXPORT tng_molecule_residues_get
+tng_function_status DECLSPECDLLEXPORT tng_molecule_num_residues_get
(tng_trajectory_t tng_data,
tng_molecule_t molecule,
- tng_residue_t *residues,
int64_t *n);
/**
- * @brief Get the atoms of a molecule.
+ * @brief Retrieve the residue of a molecule with specified index in the list
+ * of chains.
+ * @param tng_data is the trajectory data container containing the molecule.
+ * @param index is the index (in molecule->residues) of the residue to return
+ * @param molecule is the molecule from which to get the residue.
+ * @param residue is a pointer to the residue if it was found - otherwise 0.
+ * @pre \code molecule != 0 \endcode molecule must not be a NULL pointer.
+ * @pre \code residue != 0 \endcode residue must not be a NULL pointer.
+ * @return TNG_SUCCESS (0) if the residue is found or TNG_FAILURE (1) if the
+ * residue is not found.
+ */
+tng_function_status DECLSPECDLLEXPORT tng_molecule_residue_of_index_get
+ (tng_trajectory_t tng_data,
+ tng_molecule_t molecule,
+ int64_t index,
+ tng_residue_t *residue);
+
+/**
+ * @brief Get the number of atoms in a molecule.
* @param tng_data is the trajectory containing the molecule.
- * @param molecule is the molecule of which to get the atoms.
- * @param atoms is pointing to the list of atoms. This is a pointer
- * to the list of the molecule, which means that it should be handled
- * carefully, e.g. not freed.
+ * @param molecule is the molecule of which to get the number of atoms.
* @param n is pointing to a value set to the number of atoms.
* @pre \code molecule != 0 \endcode The molecule must not be NULL.
- * @pre \code atoms != 0 \endcode The pointer to the list of atoms must not
- * be a NULL pointer.
* @pre \code n != 0 \endcode The pointer to n must not be a NULL pointer.
* @return TNG_SUCCESS (0) if successful.
*/
-tng_function_status DECLSPECDLLEXPORT tng_molecule_atoms_get
+tng_function_status DECLSPECDLLEXPORT tng_molecule_num_atoms_get
(tng_trajectory_t tng_data,
tng_molecule_t molecule,
- tng_atom_t *atoms,
int64_t *n);
/**
+ * @brief Retrieve the atom of a molecule with specified index in the list
+ * of atoms.
+ * @param tng_data is the trajectory data container containing the molecule.
+ * @param index is the index (in molecule->atoms) of the atom to return
+ * @param molecule is the molecule from which to get the atom.
+ * @param atom is a pointer to the atom if it was found - otherwise 0.
+ * @pre \code molecule != 0 \endcode molecule must not be a NULL pointer.
+ * @pre \code atom != 0 \endcode atom must not be a NULL pointer.
+ * @return TNG_SUCCESS (0) if the atom is found or TNG_FAILURE (1) if the
+ * atom is not found.
+ */
+tng_function_status DECLSPECDLLEXPORT tng_molecule_atom_of_index_get
+ (tng_trajectory_t tng_data,
+ tng_molecule_t molecule,
+ int64_t index,
+ tng_atom_t *atom);
+
+/**
* @brief Find a chain in a molecule.
* @param tng_data is the trajectory data container containing the molecule.
* @param molecule is the molecule in which to search for the chain.
@@ -1818,26 +1854,38 @@ tng_function_status DECLSPECDLLEXPORT tng_chain_name_set
const char *new_name);
/**
- * @brief Get the residues of a chain.
+ * @brief Get the number of residues in a molecule chain.
* @param tng_data is the trajectory containing the chain.
- * @param chain is the chain of which to get the residues.
- * @param residues is pointing to the list of residues. This is a pointer
- * to the list of the molecule, which means that it should be handled
- * carefully, e.g. not freed.
+ * @param chain is the chain of which to get the number of residues.
* @param n is pointing to a value set to the number of residues.
* @pre \code chain != 0 \endcode The chain must not be NULL.
- * @pre \code residues != 0 \endcode The pointer to the list of residues must not
- * be a NULL pointer.
* @pre \code n != 0 \endcode The pointer to n must not be a NULL pointer.
* @return TNG_SUCCESS (0) if successful.
*/
-tng_function_status DECLSPECDLLEXPORT tng_chain_residues_get
+tng_function_status DECLSPECDLLEXPORT tng_chain_num_residues_get
(const tng_trajectory_t tng_data,
const tng_chain_t chain,
- tng_residue_t *residues,
int64_t *n);
/**
+ * @brief Retrieve the residue of a chain with specified index in the list
+ * of residues.
+ * @param tng_data is the trajectory data container containing the chain.
+ * @param index is the index (in chain->residues) of the residue to return
+ * @param chain is the chain from which to get the residue.
+ * @param residue is a pointer to the residue if it was found - otherwise 0.
+ * @pre \code chain != 0 \endcode chain must not be a NULL pointer.
+ * @pre \code residue != 0 \endcode residue must not be a NULL pointer.
+ * @return TNG_SUCCESS (0) if the residue is found or TNG_FAILURE (1) if the
+ * residue is not found.
+ */
+tng_function_status DECLSPECDLLEXPORT tng_chain_residue_of_index_get
+ (tng_trajectory_t tng_data,
+ tng_chain_t chain,
+ int64_t index,
+ tng_residue_t *residue);
+
+/**
* @brief Find a residue in a chain.
* @param tng_data is the trajectory data container containing the chain.
* @param chain is the chain in which to search for the residue.
@@ -1938,28 +1986,38 @@ tng_function_status DECLSPECDLLEXPORT tng_residue_name_set
const char *new_name);
/**
- * @brief Get the atoms of a residue.
+ * @brief Get the number of atoms in a residue.
* @param tng_data is the trajectory containing the residue.
- * @param molecule is the molecule containing the residue.
- * @param residue is the residue of which to get the atoms.
- * @param atoms is pointing to the list of atoms. This is a pointer
- * to the list of the molecule, which means that it should be handled
- * carefully, e.g. not freed.
+ * @param residue is the residue of which to get the number atoms.
* @param n is pointing to a value set to the number of atoms.
* @pre \code residue != 0 \endcode The residue must not be NULL.
- * @pre \code atoms != 0 \endcode The pointer to the list of atoms must not
- * be a NULL pointer.
* @pre \code n != 0 \endcode The pointer to n must not be a NULL pointer.
* @return TNG_SUCCESS (0) if successful.
*/
-tng_function_status DECLSPECDLLEXPORT tng_residue_atoms_get
+tng_function_status DECLSPECDLLEXPORT tng_residue_num_atoms_get
(const tng_trajectory_t tng_data,
- const tng_molecule_t molecule,
const tng_residue_t residue,
- tng_atom_t *atoms,
int64_t *n);
/**
+ * @brief Retrieve the atom of a residue with specified index in the list
+ * of atoms.
+ * @param tng_data is the trajectory data container containing the residue.
+ * @param index is the index (in residue->atoms) of the atom to return
+ * @param residue is the residue from which to get the atom.
+ * @param atom is a pointer to the atom if it was found - otherwise 0.
+ * @pre \code residue != 0 \endcode residue must not be a NULL pointer.
+ * @pre \code atom != 0 \endcode atom must not be a NULL pointer.
+ * @return TNG_SUCCESS (0) if the atom is found or TNG_FAILURE (1) if the
+ * atom is not found.
+ */
+tng_function_status DECLSPECDLLEXPORT tng_residue_atom_of_index_get
+ (tng_trajectory_t tng_data,
+ tng_residue_t residue,
+ int64_t index,
+ tng_atom_t *atom);
+
+/**
* @brief Add an atom to a residue.
* @param tng_data is the trajectory containing the residue.
* @param residue is the residue to add an atom to.
@@ -3191,12 +3249,12 @@ tng_function_status DECLSPECDLLEXPORT tng_util_time_of_frame_get
* not be a NULL pointer.
* @return TNG_SUCCESS (0) if successful.
*/
-tng_function_status DECLSPECDLLEXPORT tng_util_trajectory_molecules_get
+/*tng_function_status DECLSPECDLLEXPORT tng_util_trajectory_molecules_get
(tng_trajectory_t tng_data,
int64_t *n_mols,
int64_t **molecule_cnt_list,
tng_molecule_t *mols);
-
+*/
/*
* @brief High-level function for adding a molecule to the mol system.
* @param tng_data is the trajectory containing the mol system.
diff --git a/src/lib/tng_io.c b/src/lib/tng_io.c
index 5a6386d..abb0592 100644
--- a/src/lib/tng_io.c
+++ b/src/lib/tng_io.c
@@ -7590,57 +7590,105 @@ tng_function_status DECLSPECDLLEXPORT tng_molecule_system_copy(tng_trajectory_t
return(TNG_SUCCESS);
}
-tng_function_status DECLSPECDLLEXPORT tng_molecule_chains_get
+tng_function_status DECLSPECDLLEXPORT tng_molecule_num_chains_get
(const tng_trajectory_t tng_data,
const tng_molecule_t molecule,
- tng_chain_t *chains,
int64_t *n)
{
(void) tng_data;
TNG_ASSERT(molecule, "TNG library: molecule must not be NULL");
- TNG_ASSERT(chains, "TNG library: chains must not be a NULL pointer");
TNG_ASSERT(n, "TNG library: n must not be a NULL pointer");
- *chains = molecule->chains;
*n = molecule->n_chains;
return(TNG_SUCCESS);
}
-tng_function_status DECLSPECDLLEXPORT tng_molecule_residues_get
+tng_function_status DECLSPECDLLEXPORT tng_molecule_chain_of_index_get
+ (tng_trajectory_t tng_data,
+ tng_molecule_t molecule,
+ int64_t index,
+ tng_chain_t *chain)
+{
+ (void) tng_data;
+ TNG_ASSERT(molecule, "TNG library: molecule must not be a NULL pointer.");
+ TNG_ASSERT(chain, "TNG library: chain must not be a NULL pointer.");
+
+ if(index >= molecule->n_chains)
+ {
+ *chain = 0;
+ return(TNG_FAILURE);
+ }
+ *chain = &molecule->chains[index];
+ return(TNG_SUCCESS);
+}
+
+tng_function_status DECLSPECDLLEXPORT tng_molecule_num_residues_get
(const tng_trajectory_t tng_data,
const tng_molecule_t molecule,
- tng_residue_t *residues,
int64_t *n)
{
(void) tng_data;
TNG_ASSERT(molecule, "TNG library: molecule must not be NULL");
- TNG_ASSERT(residues, "TNG library: residues must not be a NULL pointer");
TNG_ASSERT(n, "TNG library: n must not be a NULL pointer");
- *residues = molecule->residues;
*n = molecule->n_residues;
return(TNG_SUCCESS);
}
-tng_function_status DECLSPECDLLEXPORT tng_molecule_atoms_get
+tng_function_status DECLSPECDLLEXPORT tng_molecule_residue_of_index_get
+ (tng_trajectory_t tng_data,
+ tng_molecule_t molecule,
+ int64_t index,
+ tng_residue_t *residue)
+{
+ (void) tng_data;
+ TNG_ASSERT(molecule, "TNG library: molecule must not be a NULL pointer.");
+ TNG_ASSERT(residue, "TNG library: residue must not be a NULL pointer.");
+
+ if(index >= molecule->n_residues)
+ {
+ *residue = 0;
+ return(TNG_FAILURE);
+ }
+ *residue = &molecule->residues[index];
+ return(TNG_SUCCESS);
+}
+
+tng_function_status DECLSPECDLLEXPORT tng_molecule_num_atoms_get
(const tng_trajectory_t tng_data,
const tng_molecule_t molecule,
- tng_atom_t *atoms,
int64_t *n)
{
(void) tng_data;
TNG_ASSERT(molecule, "TNG library: molecule must not be NULL");
- TNG_ASSERT(atoms, "TNG library: atoms must not be a NULL pointer");
TNG_ASSERT(n, "TNG library: n must not be a NULL pointer");
- *atoms = molecule->atoms;
*n = molecule->n_atoms;
return(TNG_SUCCESS);
}
+tng_function_status DECLSPECDLLEXPORT tng_molecule_atom_of_index_get
+ (tng_trajectory_t tng_data,
+ tng_molecule_t molecule,
+ int64_t index,
+ tng_atom_t *atom)
+{
+ (void) tng_data;
+ TNG_ASSERT(molecule, "TNG library: molecule must not be a NULL pointer.");
+ TNG_ASSERT(atom, "TNG library: atom must not be a NULL pointer.");
+
+ if(index >= molecule->n_atoms)
+ {
+ *atom = 0;
+ return(TNG_FAILURE);
+ }
+ *atom = &molecule->atoms[index];
+ return(TNG_SUCCESS);
+}
+
tng_function_status DECLSPECDLLEXPORT tng_molecule_chain_find
(tng_trajectory_t tng_data,
tng_molecule_t molecule,
@@ -7875,23 +7923,39 @@ tng_function_status DECLSPECDLLEXPORT tng_chain_name_set
return(TNG_SUCCESS);
}
-tng_function_status DECLSPECDLLEXPORT tng_chain_residues_get
+tng_function_status DECLSPECDLLEXPORT tng_chain_num_residues_get
(const tng_trajectory_t tng_data,
const tng_chain_t chain,
- tng_residue_t *residues,
int64_t *n)
{
(void) tng_data;
TNG_ASSERT(chain, "TNG library: chain must not be NULL");
- TNG_ASSERT(residues, "TNG library: residues must not be a NULL pointer");
TNG_ASSERT(n, "TNG library: n must not be a NULL pointer");
- *residues = chain->residues;
*n = chain->n_residues;
return(TNG_SUCCESS);
}
+tng_function_status DECLSPECDLLEXPORT tng_chain_residue_of_index_get
+ (tng_trajectory_t tng_data,
+ tng_chain_t chain,
+ int64_t index,
+ tng_residue_t *residue)
+{
+ (void) tng_data;
+ TNG_ASSERT(chain, "TNG library: chain must not be a NULL pointer.");
+ TNG_ASSERT(residue, "TNG library: residue must not be a NULL pointer.");
+
+ if(index >= chain->n_residues)
+ {
+ *residue = 0;
+ return(TNG_FAILURE);
+ }
+ *residue = &chain->residues[index];
+ return(TNG_SUCCESS);
+}
+
tng_function_status DECLSPECDLLEXPORT tng_chain_residue_find
(tng_trajectory_t tng_data,
tng_chain_t chain,
@@ -8085,25 +8149,51 @@ tng_function_status DECLSPECDLLEXPORT tng_residue_name_set(tng_trajectory_t tng_
return(TNG_SUCCESS);
}
-tng_function_status DECLSPECDLLEXPORT tng_residue_atoms_get
+tng_function_status DECLSPECDLLEXPORT tng_residue_num_atoms_get
(const tng_trajectory_t tng_data,
- const tng_molecule_t molecule,
const tng_residue_t residue,
- tng_atom_t *atoms,
int64_t *n)
{
(void) tng_data;
TNG_ASSERT(residue, "TNG library: residue must not be NULL");
- TNG_ASSERT(atoms, "TNG library: atoms must not be a NULL pointer");
TNG_ASSERT(n, "TNG library: n must not be a NULL pointer");
- *atoms = molecule->atoms;
- *atoms += residue->atoms_offset;
*n = residue->n_atoms;
return(TNG_SUCCESS);
}
+tng_function_status DECLSPECDLLEXPORT tng_residue_atom_of_index_get
+ (tng_trajectory_t tng_data,
+ tng_residue_t residue,
+ int64_t index,
+ tng_atom_t *atom)
+{
+ tng_chain_t chain;
+ tng_molecule_t molecule;
+
+ (void) tng_data;
+ TNG_ASSERT(residue, "TNG library: residue must not be a NULL pointer.");
+ TNG_ASSERT(atom, "TNG library: atom must not be a NULL pointer.");
+
+ if(index >= residue->n_atoms)
+ {
+ *atom = 0;
+ return(TNG_FAILURE);
+ }
+ chain = residue->chain;
+ molecule = chain->molecule;
+
+ if(index + residue->atoms_offset >= molecule->n_atoms)
+ {
+ *atom = 0;
+ return(TNG_FAILURE);
+ }
+
+ *atom = &molecule->atoms[residue->atoms_offset + index];
+ return(TNG_SUCCESS);
+}
+
tng_function_status DECLSPECDLLEXPORT tng_residue_atom_add
(tng_trajectory_t tng_data,
tng_residue_t residue,
@@ -15708,6 +15798,7 @@ tng_function_status DECLSPECDLLEXPORT tng_util_time_of_frame_get
return(TNG_SUCCESS);
}
+/*
tng_function_status DECLSPECDLLEXPORT tng_util_trajectory_molecules_get
(tng_trajectory_t tng_data,
int64_t *n_mols,
@@ -15735,7 +15826,7 @@ tng_function_status DECLSPECDLLEXPORT tng_util_trajectory_molecules_get
return(TNG_SUCCESS);
}
-
+*/
/*
tng_function_status DECLSPECDLLEXPORT tng_util_trajectory_molecule_add
(tng_trajectory_t tng_data,
contact: Jan Huwald // Impressum