2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2013,2014, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
35 #include "tngio_for_tools.h"
45 #include "tng/tng_io.h"
48 #include "gromacs/math/units.h"
49 #include "gromacs/utility/common.h"
50 #include "gromacs/utility/fatalerror.h"
51 #include "gromacs/utility/smalloc.h"
53 void gmx_prepare_tng_writing(const char *filename,
55 tng_trajectory_t *input,
56 tng_trajectory_t *output,
58 const gmx_mtop_t *mtop,
60 const char *indexGroupName)
63 /* FIXME after 5.0: Currently only standard block types are read */
64 const int defaultNumIds = 5;
65 static gmx_int64_t fallbackIds[defaultNumIds] =
67 TNG_TRAJ_BOX_SHAPE, TNG_TRAJ_POSITIONS,
68 TNG_TRAJ_VELOCITIES, TNG_TRAJ_FORCES,
71 static char fallbackNames[defaultNumIds][32] =
73 "BOX SHAPE", "POSITIONS", "VELOCITIES",
78 gmx_tng_open(filename, mode, output);
80 /* Do we have an input file in TNG format? If so, then there's
81 more data we can copy over, rather than having to improvise. */
84 /* Set parameters (compression, time per frame, molecule
85 * information, number of frames per frame set and writing
86 * intervals of positions, box shape and lambdas) of the
87 * output tng container based on their respective values int
88 * the input tng container */
89 double time, compression_precision;
90 gmx_int64_t n_frames_per_frame_set, interval = -1;
92 tng_compression_precision_get(*input, &compression_precision);
93 tng_compression_precision_set(*output, compression_precision);
94 // TODO make this configurable in a future version
95 char compression_type = TNG_TNG_COMPRESSION;
97 tng_molecule_system_copy(*input, *output);
99 tng_time_per_frame_get(*input, &time);
100 tng_time_per_frame_set(*output, time);
102 tng_num_frames_per_frame_set_get(*input, &n_frames_per_frame_set);
103 tng_num_frames_per_frame_set_set(*output, n_frames_per_frame_set);
105 for (int i = 0; i < defaultNumIds; i++)
107 if (tng_data_get_stride_length(*input, fallbackIds[i], -1, &interval)
110 switch (fallbackIds[i])
112 case TNG_TRAJ_POSITIONS:
113 case TNG_TRAJ_VELOCITIES:
114 tng_util_generic_write_interval_set(*output, interval, 3, fallbackIds[i],
115 fallbackNames[i], TNG_PARTICLE_BLOCK_DATA,
118 case TNG_TRAJ_FORCES:
119 tng_util_generic_write_interval_set(*output, interval, 3, fallbackIds[i],
120 fallbackNames[i], TNG_PARTICLE_BLOCK_DATA,
121 TNG_GZIP_COMPRESSION);
123 case TNG_TRAJ_BOX_SHAPE:
124 tng_util_generic_write_interval_set(*output, interval, 9, fallbackIds[i],
125 fallbackNames[i], TNG_NON_PARTICLE_BLOCK_DATA,
126 TNG_GZIP_COMPRESSION);
129 tng_util_generic_write_interval_set(*output, interval, 1, fallbackIds[i],
130 fallbackNames[i], TNG_NON_PARTICLE_BLOCK_DATA,
131 TNG_GZIP_COMPRESSION);
141 /* TODO after trjconv is modularized: fix this so the user can
142 change precision when they are doing an operation where
143 this makes sense, and not otherwise.
145 char compression = bUseLossyCompression ? TNG_TNG_COMPRESSION : TNG_GZIP_COMPRESSION;
146 gmx_tng_set_compression_precision(*output, ndec2prec(nDecimalsOfPrecision));
148 gmx_tng_add_mtop(*output, mtop);
149 tng_num_frames_per_frame_set_set(*output, 1);
152 if (index && nAtoms > 0)
154 gmx_tng_setup_atom_subgroup(*output, nAtoms, index, indexGroupName);
157 /* If for some reason there are more requested atoms than there are atoms in the
158 * molecular system create a number of implicit atoms (without atom data) to
159 * compensate for that. */
162 tng_implicit_num_particles_set(*output, nAtoms);
165 GMX_UNUSED_VALUE(filename);
166 GMX_UNUSED_VALUE(mode);
167 GMX_UNUSED_VALUE(input);
168 GMX_UNUSED_VALUE(output);
169 GMX_UNUSED_VALUE(nAtoms);
170 GMX_UNUSED_VALUE(mtop);
171 GMX_UNUSED_VALUE(index);
172 GMX_UNUSED_VALUE(indexGroupName);
176 void gmx_write_tng_from_trxframe(tng_trajectory_t output,
183 double timePerFrame = frame->time * PICO / frame->step;
184 tng_time_per_frame_set(output, timePerFrame);
188 natoms = frame->natoms;
190 gmx_fwrite_tng(output,
195 (const rvec *) frame->box,
197 (const rvec *) frame->x,
198 (const rvec *) frame->v,
199 (const rvec *) frame->f);
201 GMX_UNUSED_VALUE(output);
202 GMX_UNUSED_VALUE(frame);
203 GMX_UNUSED_VALUE(natoms);
209 convert_array_to_real_array(void *from,
221 if (sizeof(real) == sizeof(float))
225 memcpy(to, from, nValues * sizeof(real) * nAtoms);
229 for (i = 0; i < nAtoms; i++)
231 for (j = 0; j < nValues; j++)
233 to[i*nValues+j] = (real)((float *)from)[i*nValues+j] * fact;
240 for (i = 0; i < nAtoms; i++)
242 for (j = 0; j < nValues; j++)
244 to[i*nValues+j] = (real)((float *)from)[i*nValues+j] * fact;
250 for (i = 0; i < nAtoms; i++)
252 for (j = 0; j < nValues; j++)
254 to[i*nValues+j] = (real)((gmx_int64_t *)from)[i*nValues+j] * fact;
258 case TNG_DOUBLE_DATA:
259 if (sizeof(real) == sizeof(double))
263 memcpy(to, from, nValues * sizeof(real) * nAtoms);
267 for (i = 0; i < nAtoms; i++)
269 for (j = 0; j < nValues; j++)
271 to[i*nValues+j] = (real)((double *)from)[i*nValues+j] * fact;
278 for (i = 0; i < nAtoms; i++)
280 for (j = 0; j < nValues; j++)
282 to[i*nValues+j] = (real)((double *)from)[i*nValues+j] * fact;
288 gmx_incons("Illegal datatype when converting values to a real array!");
294 static real getDistanceScaleFactor(tng_trajectory_t in)
296 gmx_int64_t exp = -1;
297 real distanceScaleFactor;
299 // TODO Hopefully, TNG 2.0 will do this kind of thing for us
300 tng_distance_unit_exponential_get(in, &exp);
302 // GROMACS expects distances in nm
306 distanceScaleFactor = NANO/NANO;
309 distanceScaleFactor = NANO/ANGSTROM;
312 distanceScaleFactor = pow(10.0, exp + 9.0);
315 return distanceScaleFactor;
319 void gmx_tng_setup_atom_subgroup(tng_trajectory_t tng,
325 gmx_int64_t nAtoms, cnt, nMols;
326 tng_molecule_t mol, iterMol;
330 tng_function_status stat;
332 tng_num_particles_get(tng, &nAtoms);
339 stat = tng_molecule_find(tng, name, -1, &mol);
340 if (stat == TNG_SUCCESS)
342 tng_molecule_num_atoms_get(tng, mol, &nAtoms);
343 tng_molecule_cnt_get(tng, mol, &cnt);
353 if (stat == TNG_FAILURE)
355 /* The indexed atoms are added to one separate molecule. */
356 tng_molecule_alloc(tng, &mol);
357 tng_molecule_name_set(tng, mol, name);
358 tng_molecule_chain_add(tng, mol, "", &chain);
360 for (int i = 0; i < nind; i++)
362 char temp_name[256], temp_type[256];
364 /* Try to retrieve the residue name of the atom */
365 stat = tng_residue_name_of_particle_nr_get(tng, ind[i], temp_name, 256);
366 if (stat != TNG_SUCCESS)
370 /* Check if the molecule of the selection already contains this residue */
371 if (tng_chain_residue_find(tng, chain, temp_name, -1, &res)
374 tng_chain_residue_add(tng, chain, temp_name, &res);
376 /* Try to find the original name and type of the atom */
377 stat = tng_atom_name_of_particle_nr_get(tng, ind[i], temp_name, 256);
378 if (stat != TNG_SUCCESS)
382 stat = tng_atom_type_of_particle_nr_get(tng, ind[i], temp_type, 256);
383 if (stat != TNG_SUCCESS)
387 tng_residue_atom_w_id_add(tng, res, temp_name, temp_type, ind[i], &atom);
389 tng_molecule_existing_add(tng, &mol);
391 /* Set the count of the molecule containing the selected atoms to 1 and all
392 * other molecules to 0 */
393 tng_molecule_cnt_set(tng, mol, 1);
394 tng_num_molecule_types_get(tng, &nMols);
395 for (gmx_int64_t k = 0; k < nMols; k++)
397 tng_molecule_of_index_get(tng, k, &iterMol);
402 tng_molecule_cnt_set(tng, iterMol, 0);
405 GMX_UNUSED_VALUE(tng);
406 GMX_UNUSED_VALUE(nind);
407 GMX_UNUSED_VALUE(ind);
408 GMX_UNUSED_VALUE(name);
412 /* TODO: If/when TNG acquires the ability to copy data blocks without
413 * uncompressing them, then this implemenation should be reconsidered.
414 * Ideally, gmx trjconv -f a.tng -o b.tng -b 10 -e 20 would be fast
415 * and lose no information. */
416 gmx_bool gmx_read_next_tng_frame(tng_trajectory_t input,
418 gmx_int64_t *requestedIds,
423 tng_function_status stat;
424 gmx_int64_t numberOfAtoms = -1, frameNumber = -1;
425 gmx_int64_t nBlocks, blockId, *blockIds = NULL, codecId;
428 double frameTime = -1.0;
429 int size, blockDependency;
431 const int defaultNumIds = 5;
432 static gmx_int64_t fallbackRequestedIds[defaultNumIds] =
434 TNG_TRAJ_BOX_SHAPE, TNG_TRAJ_POSITIONS,
435 TNG_TRAJ_VELOCITIES, TNG_TRAJ_FORCES,
450 /* If no specific IDs were requested read all block types that can
451 * currently be interpreted */
452 if (!requestedIds || numRequestedIds == 0)
454 numRequestedIds = defaultNumIds;
455 requestedIds = fallbackRequestedIds;
458 stat = tng_num_particles_get(input, &numberOfAtoms);
459 if (stat != TNG_SUCCESS)
461 gmx_file("Cannot determine number of atoms from TNG file.");
463 fr->natoms = numberOfAtoms;
465 if (!gmx_get_tng_data_block_types_of_next_frame(input,
481 for (gmx_int64_t i = 0; i < nBlocks; i++)
483 blockId = blockIds[i];
484 tng_data_block_dependency_get(input, blockId, &blockDependency);
485 if (blockDependency & TNG_PARTICLE_DEPENDENT)
487 stat = tng_util_particle_data_next_frame_read(input,
496 stat = tng_util_non_particle_data_next_frame_read(input,
503 if (stat == TNG_CRITICAL)
505 gmx_file("Cannot read positions from TNG file.");
508 else if (stat == TNG_FAILURE)
514 case TNG_TRAJ_BOX_SHAPE:
518 size = sizeof(gmx_int64_t);
521 size = sizeof(float);
523 case TNG_DOUBLE_DATA:
524 size = sizeof(double);
527 gmx_incons("Illegal datatype of box shape values!");
529 for (int i = 0; i < DIM; i++)
531 convert_array_to_real_array((char *)(values) + size * i * DIM,
533 getDistanceScaleFactor(input),
540 case TNG_TRAJ_POSITIONS:
541 srenew(fr->x, fr->natoms);
542 convert_array_to_real_array(values,
544 getDistanceScaleFactor(input),
549 tng_util_frame_current_compression_get(input, blockId, &codecId, &prec);
550 /* This must be updated if/when more lossy compression methods are added */
551 if (codecId == TNG_TNG_COMPRESSION)
557 case TNG_TRAJ_VELOCITIES:
558 srenew(fr->v, fr->natoms);
559 convert_array_to_real_array(values,
561 getDistanceScaleFactor(input),
566 tng_util_frame_current_compression_get(input, blockId, &codecId, &prec);
567 /* This must be updated if/when more lossy compression methods are added */
568 if (codecId == TNG_TNG_COMPRESSION)
574 case TNG_TRAJ_FORCES:
575 srenew(fr->f, fr->natoms);
576 convert_array_to_real_array(values,
578 getDistanceScaleFactor(input),
588 fr->lambda = (*(float *)values);
590 case TNG_DOUBLE_DATA:
591 fr->lambda = (*(double *)values);
594 gmx_incons("Illegal datatype lambda value!");
599 gmx_warning("Illegal block type! Currently GROMACS tools can only handle certain data types. Skipping block.");
601 /* values does not have to be freed before reading next frame. It will
602 * be reallocated if it is not NULL. */
605 fr->step = (int) frameNumber;
607 // Convert the time to ps
608 fr->time = frameTime / PICO;
611 /* values must be freed before leaving this function */
616 GMX_UNUSED_VALUE(input);
617 GMX_UNUSED_VALUE(fr);
618 GMX_UNUSED_VALUE(requestedIds);
619 GMX_UNUSED_VALUE(numRequestedIds);
624 void gmx_print_tng_molecule_system(tng_trajectory_t input,
628 gmx_int64_t nMolecules, nChains, nResidues, nAtoms, *molCntList;
629 tng_molecule_t molecule;
631 tng_residue_t residue;
633 char str[256], varNAtoms;
635 tng_num_molecule_types_get(input, &nMolecules);
636 tng_molecule_cnt_list_get(input, &molCntList);
637 /* Can the number of particles change in the trajectory or is it constant? */
638 tng_num_particles_variable_get(input, &varNAtoms);
640 for (gmx_int64_t i = 0; i < nMolecules; i++)
642 tng_molecule_of_index_get(input, i, &molecule);
643 tng_molecule_name_get(input, molecule, str, 256);
644 if (varNAtoms == TNG_CONSTANT_N_ATOMS)
646 if ((int)molCntList[i] == 0)
650 fprintf(stream, "Molecule: %s, count: %d\n", str, (int)molCntList[i]);
654 fprintf(stream, "Molecule: %s\n", str);
656 tng_molecule_num_chains_get(input, molecule, &nChains);
659 for (gmx_int64_t j = 0; j < nChains; j++)
661 tng_molecule_chain_of_index_get(input, molecule, j, &chain);
662 tng_chain_name_get(input, chain, str, 256);
663 fprintf(stream, "\tChain: %s\n", str);
664 tng_chain_num_residues_get(input, chain, &nResidues);
665 for (gmx_int64_t k = 0; k < nResidues; k++)
667 tng_chain_residue_of_index_get(input, chain, k, &residue);
668 tng_residue_name_get(input, residue, str, 256);
669 fprintf(stream, "\t\tResidue: %s\n", str);
670 tng_residue_num_atoms_get(input, residue, &nAtoms);
671 for (gmx_int64_t l = 0; l < nAtoms; l++)
673 tng_residue_atom_of_index_get(input, residue, l, &atom);
674 tng_atom_name_get(input, atom, str, 256);
675 fprintf(stream, "\t\t\tAtom: %s", str);
676 tng_atom_type_get(input, atom, str, 256);
677 fprintf(stream, " (%s)\n", str);
682 /* It is possible to have a molecule without chains, in which case
683 * residues in the molecule can be iterated through without going
687 tng_molecule_num_residues_get(input, molecule, &nResidues);
690 for (gmx_int64_t k = 0; k < nResidues; k++)
692 tng_molecule_residue_of_index_get(input, molecule, k, &residue);
693 tng_residue_name_get(input, residue, str, 256);
694 fprintf(stream, "\t\tResidue: %s\n", str);
695 tng_residue_num_atoms_get(input, residue, &nAtoms);
696 for (gmx_int64_t l = 0; l < nAtoms; l++)
698 tng_residue_atom_of_index_get(input, residue, l, &atom);
699 tng_atom_name_get(input, atom, str, 256);
700 fprintf(stream, "\t\t\tAtom: %s", str);
701 tng_atom_type_get(input, atom, str, 256);
702 fprintf(stream, " (%s)\n", str);
708 tng_molecule_num_atoms_get(input, molecule, &nAtoms);
709 for (gmx_int64_t l = 0; l < nAtoms; l++)
711 tng_molecule_atom_of_index_get(input, molecule, l, &atom);
712 tng_atom_name_get(input, atom, str, 256);
713 fprintf(stream, "\t\t\tAtom: %s", str);
714 tng_atom_type_get(input, atom, str, 256);
715 fprintf(stream, " (%s)\n", str);
721 GMX_UNUSED_VALUE(input);
722 GMX_UNUSED_VALUE(stream);
726 gmx_bool gmx_get_tng_data_block_types_of_next_frame(tng_trajectory_t input,
729 gmx_int64_t *requestedIds,
730 gmx_int64_t *nextFrame,
731 gmx_int64_t *nBlocks,
732 gmx_int64_t **blockIds)
735 tng_function_status stat;
737 stat = tng_util_trajectory_next_frame_present_data_blocks_find(input, frame,
738 nRequestedIds, requestedIds,
742 if (stat == TNG_CRITICAL)
744 gmx_file("Cannot read TNG file. Cannot find data blocks of next frame.");
746 else if (stat == TNG_FAILURE)
752 GMX_UNUSED_VALUE(input);
753 GMX_UNUSED_VALUE(frame);
754 GMX_UNUSED_VALUE(nRequestedIds);
755 GMX_UNUSED_VALUE(requestedIds);
756 GMX_UNUSED_VALUE(nextFrame);
757 GMX_UNUSED_VALUE(nBlocks);
758 GMX_UNUSED_VALUE(blockIds);
763 gmx_bool gmx_get_tng_data_next_frame_of_block_type(tng_trajectory_t input,
766 gmx_int64_t *frameNumber,
768 gmx_int64_t *nValuesPerFrame,
776 tng_function_status stat;
783 stat = tng_data_block_name_get(input, blockId, name, maxLen);
784 if (stat != TNG_SUCCESS)
786 gmx_file("Cannot read next frame of TNG file");
788 stat = tng_data_block_dependency_get(input, blockId, &blockDependency);
789 if (stat != TNG_SUCCESS)
791 gmx_file("Cannot read next frame of TNG file");
793 if (blockDependency & TNG_PARTICLE_DEPENDENT)
795 tng_num_particles_get(input, nAtoms);
796 stat = tng_util_particle_data_next_frame_read(input,
805 *nAtoms = 1; /* There are not actually any atoms, but it is used for
807 stat = tng_util_non_particle_data_next_frame_read(input,
814 if (stat == TNG_CRITICAL)
816 gmx_file("Cannot read next frame of TNG file");
818 if (stat == TNG_FAILURE)
824 stat = tng_data_block_num_values_per_frame_get(input, blockId, nValuesPerFrame);
825 if (stat != TNG_SUCCESS)
827 gmx_file("Cannot read next frame of TNG file");
829 snew(*values, sizeof(real) * *nValuesPerFrame * *nAtoms);
830 convert_array_to_real_array(data,
832 getDistanceScaleFactor(input),
837 tng_util_frame_current_compression_get(input, blockId, &codecId, &localPrec);
839 /* This must be updated if/when more lossy compression methods are added */
840 if (codecId != TNG_TNG_COMPRESSION)
852 GMX_UNUSED_VALUE(input);
853 GMX_UNUSED_VALUE(blockId);
854 GMX_UNUSED_VALUE(values);
855 GMX_UNUSED_VALUE(frameNumber);
856 GMX_UNUSED_VALUE(frameTime);
857 GMX_UNUSED_VALUE(nValuesPerFrame);
858 GMX_UNUSED_VALUE(nAtoms);
859 GMX_UNUSED_VALUE(prec);
860 GMX_UNUSED_VALUE(name);
861 GMX_UNUSED_VALUE(maxLen);
862 GMX_UNUSED_VALUE(bOK);