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.
37 #include "tngio_for_tools.h"
47 #include "tng/tng_io.h"
50 #include "gromacs/math/units.h"
51 #include "gromacs/utility/common.h"
52 #include "gromacs/utility/fatalerror.h"
53 #include "gromacs/utility/smalloc.h"
55 void gmx_prepare_tng_writing(const char *filename,
57 tng_trajectory_t *input,
58 tng_trajectory_t *output,
60 const gmx_mtop_t *mtop,
62 const char *indexGroupName)
65 /* FIXME after 5.0: Currently only standard block types are read */
66 const int defaultNumIds = 5;
67 static gmx_int64_t fallbackIds[defaultNumIds] =
69 TNG_TRAJ_BOX_SHAPE, TNG_TRAJ_POSITIONS,
70 TNG_TRAJ_VELOCITIES, TNG_TRAJ_FORCES,
73 static char fallbackNames[defaultNumIds][32] =
75 "BOX SHAPE", "POSITIONS", "VELOCITIES",
80 gmx_tng_open(filename, mode, output);
82 /* Do we have an input file in TNG format? If so, then there's
83 more data we can copy over, rather than having to improvise. */
86 /* Set parameters (compression, time per frame, molecule
87 * information, number of frames per frame set and writing
88 * intervals of positions, box shape and lambdas) of the
89 * output tng container based on their respective values int
90 * the input tng container */
91 double time, compression_precision;
92 gmx_int64_t n_frames_per_frame_set, interval = -1;
94 tng_compression_precision_get(*input, &compression_precision);
95 tng_compression_precision_set(*output, compression_precision);
96 // TODO make this configurable in a future version
97 char compression_type = TNG_TNG_COMPRESSION;
99 tng_molecule_system_copy(*input, *output);
101 tng_time_per_frame_get(*input, &time);
102 tng_time_per_frame_set(*output, time);
104 tng_num_frames_per_frame_set_get(*input, &n_frames_per_frame_set);
105 tng_num_frames_per_frame_set_set(*output, n_frames_per_frame_set);
107 for (int i = 0; i < defaultNumIds; i++)
109 if (tng_data_get_stride_length(*input, fallbackIds[i], -1, &interval)
112 switch (fallbackIds[i])
114 case TNG_TRAJ_POSITIONS:
115 case TNG_TRAJ_VELOCITIES:
116 tng_util_generic_write_interval_set(*output, interval, 3, fallbackIds[i],
117 fallbackNames[i], TNG_PARTICLE_BLOCK_DATA,
120 case TNG_TRAJ_FORCES:
121 tng_util_generic_write_interval_set(*output, interval, 3, fallbackIds[i],
122 fallbackNames[i], TNG_PARTICLE_BLOCK_DATA,
123 TNG_GZIP_COMPRESSION);
125 case TNG_TRAJ_BOX_SHAPE:
126 tng_util_generic_write_interval_set(*output, interval, 9, fallbackIds[i],
127 fallbackNames[i], TNG_NON_PARTICLE_BLOCK_DATA,
128 TNG_GZIP_COMPRESSION);
131 tng_util_generic_write_interval_set(*output, interval, 1, fallbackIds[i],
132 fallbackNames[i], TNG_NON_PARTICLE_BLOCK_DATA,
133 TNG_GZIP_COMPRESSION);
143 /* TODO after trjconv is modularized: fix this so the user can
144 change precision when they are doing an operation where
145 this makes sense, and not otherwise.
147 char compression = bUseLossyCompression ? TNG_TNG_COMPRESSION : TNG_GZIP_COMPRESSION;
148 gmx_tng_set_compression_precision(*output, ndec2prec(nDecimalsOfPrecision));
150 gmx_tng_add_mtop(*output, mtop);
151 tng_num_frames_per_frame_set_set(*output, 1);
154 if (index && nAtoms > 0)
156 gmx_tng_setup_atom_subgroup(*output, nAtoms, index, indexGroupName);
159 /* If for some reason there are more requested atoms than there are atoms in the
160 * molecular system create a number of implicit atoms (without atom data) to
161 * compensate for that. */
164 tng_implicit_num_particles_set(*output, nAtoms);
167 GMX_UNUSED_VALUE(filename);
168 GMX_UNUSED_VALUE(mode);
169 GMX_UNUSED_VALUE(input);
170 GMX_UNUSED_VALUE(output);
171 GMX_UNUSED_VALUE(nAtoms);
172 GMX_UNUSED_VALUE(mtop);
173 GMX_UNUSED_VALUE(index);
174 GMX_UNUSED_VALUE(indexGroupName);
178 void gmx_write_tng_from_trxframe(tng_trajectory_t output,
185 double timePerFrame = frame->time * PICO / frame->step;
186 tng_time_per_frame_set(output, timePerFrame);
190 natoms = frame->natoms;
192 gmx_fwrite_tng(output,
197 (const rvec *) frame->box,
199 (const rvec *) frame->x,
200 (const rvec *) frame->v,
201 (const rvec *) frame->f);
203 GMX_UNUSED_VALUE(output);
204 GMX_UNUSED_VALUE(frame);
205 GMX_UNUSED_VALUE(natoms);
211 convert_array_to_real_array(void *from,
223 if (sizeof(real) == sizeof(float))
227 memcpy(to, from, nValues * sizeof(real) * nAtoms);
231 for (i = 0; i < nAtoms; i++)
233 for (j = 0; j < nValues; j++)
235 to[i*nValues+j] = (real)((float *)from)[i*nValues+j] * fact;
242 for (i = 0; i < nAtoms; i++)
244 for (j = 0; j < nValues; j++)
246 to[i*nValues+j] = (real)((float *)from)[i*nValues+j] * fact;
252 for (i = 0; i < nAtoms; i++)
254 for (j = 0; j < nValues; j++)
256 to[i*nValues+j] = (real)((gmx_int64_t *)from)[i*nValues+j] * fact;
260 case TNG_DOUBLE_DATA:
261 if (sizeof(real) == sizeof(double))
265 memcpy(to, from, nValues * sizeof(real) * nAtoms);
269 for (i = 0; i < nAtoms; i++)
271 for (j = 0; j < nValues; j++)
273 to[i*nValues+j] = (real)((double *)from)[i*nValues+j] * fact;
280 for (i = 0; i < nAtoms; i++)
282 for (j = 0; j < nValues; j++)
284 to[i*nValues+j] = (real)((double *)from)[i*nValues+j] * fact;
290 gmx_incons("Illegal datatype when converting values to a real array!");
296 static real getDistanceScaleFactor(tng_trajectory_t in)
298 gmx_int64_t exp = -1;
299 real distanceScaleFactor;
301 // TODO Hopefully, TNG 2.0 will do this kind of thing for us
302 tng_distance_unit_exponential_get(in, &exp);
304 // GROMACS expects distances in nm
308 distanceScaleFactor = NANO/NANO;
311 distanceScaleFactor = NANO/ANGSTROM;
314 distanceScaleFactor = pow(10.0, exp + 9.0);
317 return distanceScaleFactor;
321 void gmx_tng_setup_atom_subgroup(tng_trajectory_t tng,
327 gmx_int64_t nAtoms, cnt, nMols;
328 tng_molecule_t mol, iterMol;
332 tng_function_status stat;
334 tng_num_particles_get(tng, &nAtoms);
341 stat = tng_molecule_find(tng, name, -1, &mol);
342 if (stat == TNG_SUCCESS)
344 tng_molecule_num_atoms_get(tng, mol, &nAtoms);
345 tng_molecule_cnt_get(tng, mol, &cnt);
355 if (stat == TNG_FAILURE)
357 /* The indexed atoms are added to one separate molecule. */
358 tng_molecule_alloc(tng, &mol);
359 tng_molecule_name_set(tng, mol, name);
360 tng_molecule_chain_add(tng, mol, "", &chain);
362 for (int i = 0; i < nind; i++)
364 char temp_name[256], temp_type[256];
366 /* Try to retrieve the residue name of the atom */
367 stat = tng_residue_name_of_particle_nr_get(tng, ind[i], temp_name, 256);
368 if (stat != TNG_SUCCESS)
372 /* Check if the molecule of the selection already contains this residue */
373 if (tng_chain_residue_find(tng, chain, temp_name, -1, &res)
376 tng_chain_residue_add(tng, chain, temp_name, &res);
378 /* Try to find the original name and type of the atom */
379 stat = tng_atom_name_of_particle_nr_get(tng, ind[i], temp_name, 256);
380 if (stat != TNG_SUCCESS)
384 stat = tng_atom_type_of_particle_nr_get(tng, ind[i], temp_type, 256);
385 if (stat != TNG_SUCCESS)
389 tng_residue_atom_w_id_add(tng, res, temp_name, temp_type, ind[i], &atom);
391 tng_molecule_existing_add(tng, &mol);
393 /* Set the count of the molecule containing the selected atoms to 1 and all
394 * other molecules to 0 */
395 tng_molecule_cnt_set(tng, mol, 1);
396 tng_num_molecule_types_get(tng, &nMols);
397 for (gmx_int64_t k = 0; k < nMols; k++)
399 tng_molecule_of_index_get(tng, k, &iterMol);
404 tng_molecule_cnt_set(tng, iterMol, 0);
407 GMX_UNUSED_VALUE(tng);
408 GMX_UNUSED_VALUE(nind);
409 GMX_UNUSED_VALUE(ind);
410 GMX_UNUSED_VALUE(name);
414 /* TODO: If/when TNG acquires the ability to copy data blocks without
415 * uncompressing them, then this implemenation should be reconsidered.
416 * Ideally, gmx trjconv -f a.tng -o b.tng -b 10 -e 20 would be fast
417 * and lose no information. */
418 gmx_bool gmx_read_next_tng_frame(tng_trajectory_t input,
420 gmx_int64_t *requestedIds,
425 tng_function_status stat;
426 gmx_int64_t numberOfAtoms = -1, frameNumber = -1;
427 gmx_int64_t nBlocks, blockId, *blockIds = NULL, codecId;
430 double frameTime = -1.0;
431 int size, blockDependency;
433 const int defaultNumIds = 5;
434 static gmx_int64_t fallbackRequestedIds[defaultNumIds] =
436 TNG_TRAJ_BOX_SHAPE, TNG_TRAJ_POSITIONS,
437 TNG_TRAJ_VELOCITIES, TNG_TRAJ_FORCES,
452 /* If no specific IDs were requested read all block types that can
453 * currently be interpreted */
454 if (!requestedIds || numRequestedIds == 0)
456 numRequestedIds = defaultNumIds;
457 requestedIds = fallbackRequestedIds;
460 stat = tng_num_particles_get(input, &numberOfAtoms);
461 if (stat != TNG_SUCCESS)
463 gmx_file("Cannot determine number of atoms from TNG file.");
465 fr->natoms = numberOfAtoms;
467 if (!gmx_get_tng_data_block_types_of_next_frame(input,
483 for (gmx_int64_t i = 0; i < nBlocks; i++)
485 blockId = blockIds[i];
486 tng_data_block_dependency_get(input, blockId, &blockDependency);
487 if (blockDependency & TNG_PARTICLE_DEPENDENT)
489 stat = tng_util_particle_data_next_frame_read(input,
498 stat = tng_util_non_particle_data_next_frame_read(input,
505 if (stat == TNG_CRITICAL)
507 gmx_file("Cannot read positions from TNG file.");
510 else if (stat == TNG_FAILURE)
516 case TNG_TRAJ_BOX_SHAPE:
520 size = sizeof(gmx_int64_t);
523 size = sizeof(float);
525 case TNG_DOUBLE_DATA:
526 size = sizeof(double);
529 gmx_incons("Illegal datatype of box shape values!");
531 for (int i = 0; i < DIM; i++)
533 convert_array_to_real_array((char *)(values) + size * i * DIM,
535 getDistanceScaleFactor(input),
542 case TNG_TRAJ_POSITIONS:
543 srenew(fr->x, fr->natoms);
544 convert_array_to_real_array(values,
546 getDistanceScaleFactor(input),
551 tng_util_frame_current_compression_get(input, blockId, &codecId, &prec);
552 /* This must be updated if/when more lossy compression methods are added */
553 if (codecId == TNG_TNG_COMPRESSION)
559 case TNG_TRAJ_VELOCITIES:
560 srenew(fr->v, fr->natoms);
561 convert_array_to_real_array(values,
563 getDistanceScaleFactor(input),
568 tng_util_frame_current_compression_get(input, blockId, &codecId, &prec);
569 /* This must be updated if/when more lossy compression methods are added */
570 if (codecId == TNG_TNG_COMPRESSION)
576 case TNG_TRAJ_FORCES:
577 srenew(fr->f, fr->natoms);
578 convert_array_to_real_array(values,
580 getDistanceScaleFactor(input),
590 fr->lambda = (*(float *)values);
592 case TNG_DOUBLE_DATA:
593 fr->lambda = (*(double *)values);
596 gmx_incons("Illegal datatype lambda value!");
601 gmx_warning("Illegal block type! Currently GROMACS tools can only handle certain data types. Skipping block.");
603 /* values does not have to be freed before reading next frame. It will
604 * be reallocated if it is not NULL. */
607 fr->step = (int) frameNumber;
609 // Convert the time to ps
610 fr->time = frameTime / PICO;
613 /* values must be freed before leaving this function */
618 GMX_UNUSED_VALUE(input);
619 GMX_UNUSED_VALUE(fr);
620 GMX_UNUSED_VALUE(requestedIds);
621 GMX_UNUSED_VALUE(numRequestedIds);
626 void gmx_print_tng_molecule_system(tng_trajectory_t input,
630 gmx_int64_t nMolecules, nChains, nResidues, nAtoms, *molCntList;
631 tng_molecule_t molecule;
633 tng_residue_t residue;
635 char str[256], varNAtoms;
637 tng_num_molecule_types_get(input, &nMolecules);
638 tng_molecule_cnt_list_get(input, &molCntList);
639 /* Can the number of particles change in the trajectory or is it constant? */
640 tng_num_particles_variable_get(input, &varNAtoms);
642 for (gmx_int64_t i = 0; i < nMolecules; i++)
644 tng_molecule_of_index_get(input, i, &molecule);
645 tng_molecule_name_get(input, molecule, str, 256);
646 if (varNAtoms == TNG_CONSTANT_N_ATOMS)
648 if ((int)molCntList[i] == 0)
652 fprintf(stream, "Molecule: %s, count: %d\n", str, (int)molCntList[i]);
656 fprintf(stream, "Molecule: %s\n", str);
658 tng_molecule_num_chains_get(input, molecule, &nChains);
661 for (gmx_int64_t j = 0; j < nChains; j++)
663 tng_molecule_chain_of_index_get(input, molecule, j, &chain);
664 tng_chain_name_get(input, chain, str, 256);
665 fprintf(stream, "\tChain: %s\n", str);
666 tng_chain_num_residues_get(input, chain, &nResidues);
667 for (gmx_int64_t k = 0; k < nResidues; k++)
669 tng_chain_residue_of_index_get(input, chain, k, &residue);
670 tng_residue_name_get(input, residue, str, 256);
671 fprintf(stream, "\t\tResidue: %s\n", str);
672 tng_residue_num_atoms_get(input, residue, &nAtoms);
673 for (gmx_int64_t l = 0; l < nAtoms; l++)
675 tng_residue_atom_of_index_get(input, residue, l, &atom);
676 tng_atom_name_get(input, atom, str, 256);
677 fprintf(stream, "\t\t\tAtom: %s", str);
678 tng_atom_type_get(input, atom, str, 256);
679 fprintf(stream, " (%s)\n", str);
684 /* It is possible to have a molecule without chains, in which case
685 * residues in the molecule can be iterated through without going
689 tng_molecule_num_residues_get(input, molecule, &nResidues);
692 for (gmx_int64_t k = 0; k < nResidues; k++)
694 tng_molecule_residue_of_index_get(input, molecule, k, &residue);
695 tng_residue_name_get(input, residue, str, 256);
696 fprintf(stream, "\t\tResidue: %s\n", str);
697 tng_residue_num_atoms_get(input, residue, &nAtoms);
698 for (gmx_int64_t l = 0; l < nAtoms; l++)
700 tng_residue_atom_of_index_get(input, residue, l, &atom);
701 tng_atom_name_get(input, atom, str, 256);
702 fprintf(stream, "\t\t\tAtom: %s", str);
703 tng_atom_type_get(input, atom, str, 256);
704 fprintf(stream, " (%s)\n", str);
710 tng_molecule_num_atoms_get(input, molecule, &nAtoms);
711 for (gmx_int64_t l = 0; l < nAtoms; l++)
713 tng_molecule_atom_of_index_get(input, molecule, l, &atom);
714 tng_atom_name_get(input, atom, str, 256);
715 fprintf(stream, "\t\t\tAtom: %s", str);
716 tng_atom_type_get(input, atom, str, 256);
717 fprintf(stream, " (%s)\n", str);
723 GMX_UNUSED_VALUE(input);
724 GMX_UNUSED_VALUE(stream);
728 gmx_bool gmx_get_tng_data_block_types_of_next_frame(tng_trajectory_t input,
731 gmx_int64_t *requestedIds,
732 gmx_int64_t *nextFrame,
733 gmx_int64_t *nBlocks,
734 gmx_int64_t **blockIds)
737 tng_function_status stat;
739 stat = tng_util_trajectory_next_frame_present_data_blocks_find(input, frame,
740 nRequestedIds, requestedIds,
744 if (stat == TNG_CRITICAL)
746 gmx_file("Cannot read TNG file. Cannot find data blocks of next frame.");
748 else if (stat == TNG_FAILURE)
754 GMX_UNUSED_VALUE(input);
755 GMX_UNUSED_VALUE(frame);
756 GMX_UNUSED_VALUE(nRequestedIds);
757 GMX_UNUSED_VALUE(requestedIds);
758 GMX_UNUSED_VALUE(nextFrame);
759 GMX_UNUSED_VALUE(nBlocks);
760 GMX_UNUSED_VALUE(blockIds);
765 gmx_bool gmx_get_tng_data_next_frame_of_block_type(tng_trajectory_t input,
768 gmx_int64_t *frameNumber,
770 gmx_int64_t *nValuesPerFrame,
778 tng_function_status stat;
785 stat = tng_data_block_name_get(input, blockId, name, maxLen);
786 if (stat != TNG_SUCCESS)
788 gmx_file("Cannot read next frame of TNG file");
790 stat = tng_data_block_dependency_get(input, blockId, &blockDependency);
791 if (stat != TNG_SUCCESS)
793 gmx_file("Cannot read next frame of TNG file");
795 if (blockDependency & TNG_PARTICLE_DEPENDENT)
797 tng_num_particles_get(input, nAtoms);
798 stat = tng_util_particle_data_next_frame_read(input,
807 *nAtoms = 1; /* There are not actually any atoms, but it is used for
809 stat = tng_util_non_particle_data_next_frame_read(input,
816 if (stat == TNG_CRITICAL)
818 gmx_file("Cannot read next frame of TNG file");
820 if (stat == TNG_FAILURE)
826 stat = tng_data_block_num_values_per_frame_get(input, blockId, nValuesPerFrame);
827 if (stat != TNG_SUCCESS)
829 gmx_file("Cannot read next frame of TNG file");
831 snew(*values, sizeof(real) * *nValuesPerFrame * *nAtoms);
832 convert_array_to_real_array(data,
834 getDistanceScaleFactor(input),
839 tng_util_frame_current_compression_get(input, blockId, &codecId, &localPrec);
841 /* This must be updated if/when more lossy compression methods are added */
842 if (codecId != TNG_TNG_COMPRESSION)
854 GMX_UNUSED_VALUE(input);
855 GMX_UNUSED_VALUE(blockId);
856 GMX_UNUSED_VALUE(values);
857 GMX_UNUSED_VALUE(frameNumber);
858 GMX_UNUSED_VALUE(frameTime);
859 GMX_UNUSED_VALUE(nValuesPerFrame);
860 GMX_UNUSED_VALUE(nAtoms);
861 GMX_UNUSED_VALUE(prec);
862 GMX_UNUSED_VALUE(name);
863 GMX_UNUSED_VALUE(maxLen);
864 GMX_UNUSED_VALUE(bOK);