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",
79 typedef tng_function_status (*set_writing_interval_func_pointer)(tng_trajectory_t,
87 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_double_set;
89 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_set;
92 gmx_tng_open(filename, mode, output);
94 /* Do we have an input file in TNG format? If so, then there's
95 more data we can copy over, rather than having to improvise. */
98 /* Set parameters (compression, time per frame, molecule
99 * information, number of frames per frame set and writing
100 * intervals of positions, box shape and lambdas) of the
101 * output tng container based on their respective values int
102 * the input tng container */
103 double time, compression_precision;
104 gmx_int64_t n_frames_per_frame_set, interval = -1;
106 tng_compression_precision_get(*input, &compression_precision);
107 tng_compression_precision_set(*output, compression_precision);
108 // TODO make this configurable in a future version
109 char compression_type = TNG_TNG_COMPRESSION;
111 tng_molecule_system_copy(*input, *output);
113 tng_time_per_frame_get(*input, &time);
114 tng_time_per_frame_set(*output, time);
116 tng_num_frames_per_frame_set_get(*input, &n_frames_per_frame_set);
117 tng_num_frames_per_frame_set_set(*output, n_frames_per_frame_set);
119 for (int i = 0; i < defaultNumIds; i++)
121 if (tng_data_get_stride_length(*input, fallbackIds[i], -1, &interval)
124 switch (fallbackIds[i])
126 case TNG_TRAJ_POSITIONS:
127 case TNG_TRAJ_VELOCITIES:
128 set_writing_interval(*output, interval, 3, fallbackIds[i],
129 fallbackNames[i], TNG_PARTICLE_BLOCK_DATA,
132 case TNG_TRAJ_FORCES:
133 set_writing_interval(*output, interval, 3, fallbackIds[i],
134 fallbackNames[i], TNG_PARTICLE_BLOCK_DATA,
135 TNG_GZIP_COMPRESSION);
137 case TNG_TRAJ_BOX_SHAPE:
138 set_writing_interval(*output, interval, 9, fallbackIds[i],
139 fallbackNames[i], TNG_NON_PARTICLE_BLOCK_DATA,
140 TNG_GZIP_COMPRESSION);
143 set_writing_interval(*output, interval, 1, fallbackIds[i],
144 fallbackNames[i], TNG_NON_PARTICLE_BLOCK_DATA,
145 TNG_GZIP_COMPRESSION);
155 /* TODO after trjconv is modularized: fix this so the user can
156 change precision when they are doing an operation where
157 this makes sense, and not otherwise.
159 char compression = bUseLossyCompression ? TNG_TNG_COMPRESSION : TNG_GZIP_COMPRESSION;
160 gmx_tng_set_compression_precision(*output, ndec2prec(nDecimalsOfPrecision));
162 gmx_tng_add_mtop(*output, mtop);
163 tng_num_frames_per_frame_set_set(*output, 1);
166 if (index && nAtoms > 0)
168 gmx_tng_setup_atom_subgroup(*output, nAtoms, index, indexGroupName);
171 /* If for some reason there are more requested atoms than there are atoms in the
172 * molecular system create a number of implicit atoms (without atom data) to
173 * compensate for that. */
176 tng_implicit_num_particles_set(*output, nAtoms);
179 GMX_UNUSED_VALUE(filename);
180 GMX_UNUSED_VALUE(mode);
181 GMX_UNUSED_VALUE(input);
182 GMX_UNUSED_VALUE(output);
183 GMX_UNUSED_VALUE(nAtoms);
184 GMX_UNUSED_VALUE(mtop);
185 GMX_UNUSED_VALUE(index);
186 GMX_UNUSED_VALUE(indexGroupName);
190 void gmx_write_tng_from_trxframe(tng_trajectory_t output,
197 double timePerFrame = frame->time * PICO / frame->step;
198 tng_time_per_frame_set(output, timePerFrame);
202 natoms = frame->natoms;
204 gmx_fwrite_tng(output,
209 (const rvec *) frame->box,
211 (const rvec *) frame->x,
212 (const rvec *) frame->v,
213 (const rvec *) frame->f);
215 GMX_UNUSED_VALUE(output);
216 GMX_UNUSED_VALUE(frame);
217 GMX_UNUSED_VALUE(natoms);
223 convert_array_to_real_array(void *from,
235 if (sizeof(real) == sizeof(float))
239 memcpy(to, from, nValues * sizeof(real) * nAtoms);
243 for (i = 0; i < nAtoms; i++)
245 for (j = 0; j < nValues; j++)
247 to[i*nValues+j] = (real)((float *)from)[i*nValues+j] * fact;
254 for (i = 0; i < nAtoms; i++)
256 for (j = 0; j < nValues; j++)
258 to[i*nValues+j] = (real)((float *)from)[i*nValues+j] * fact;
264 for (i = 0; i < nAtoms; i++)
266 for (j = 0; j < nValues; j++)
268 to[i*nValues+j] = (real)((gmx_int64_t *)from)[i*nValues+j] * fact;
272 case TNG_DOUBLE_DATA:
273 if (sizeof(real) == sizeof(double))
277 memcpy(to, from, nValues * sizeof(real) * nAtoms);
281 for (i = 0; i < nAtoms; i++)
283 for (j = 0; j < nValues; j++)
285 to[i*nValues+j] = (real)((double *)from)[i*nValues+j] * fact;
292 for (i = 0; i < nAtoms; i++)
294 for (j = 0; j < nValues; j++)
296 to[i*nValues+j] = (real)((double *)from)[i*nValues+j] * fact;
302 gmx_incons("Illegal datatype when converting values to a real array!");
308 static real getDistanceScaleFactor(tng_trajectory_t in)
310 gmx_int64_t exp = -1;
311 real distanceScaleFactor;
313 // TODO Hopefully, TNG 2.0 will do this kind of thing for us
314 tng_distance_unit_exponential_get(in, &exp);
316 // GROMACS expects distances in nm
320 distanceScaleFactor = NANO/NANO;
323 distanceScaleFactor = NANO/ANGSTROM;
326 distanceScaleFactor = pow(10.0, exp + 9.0);
329 return distanceScaleFactor;
333 void gmx_tng_setup_atom_subgroup(tng_trajectory_t tng,
339 gmx_int64_t nAtoms, cnt, nMols;
340 tng_molecule_t mol, iterMol;
344 tng_function_status stat;
346 tng_num_particles_get(tng, &nAtoms);
353 stat = tng_molecule_find(tng, name, -1, &mol);
354 if (stat == TNG_SUCCESS)
356 tng_molecule_num_atoms_get(tng, mol, &nAtoms);
357 tng_molecule_cnt_get(tng, mol, &cnt);
367 if (stat == TNG_FAILURE)
369 /* The indexed atoms are added to one separate molecule. */
370 tng_molecule_alloc(tng, &mol);
371 tng_molecule_name_set(tng, mol, name);
372 tng_molecule_chain_add(tng, mol, "", &chain);
374 for (int i = 0; i < nind; i++)
376 char temp_name[256], temp_type[256];
378 /* Try to retrieve the residue name of the atom */
379 stat = tng_residue_name_of_particle_nr_get(tng, ind[i], temp_name, 256);
380 if (stat != TNG_SUCCESS)
384 /* Check if the molecule of the selection already contains this residue */
385 if (tng_chain_residue_find(tng, chain, temp_name, -1, &res)
388 tng_chain_residue_add(tng, chain, temp_name, &res);
390 /* Try to find the original name and type of the atom */
391 stat = tng_atom_name_of_particle_nr_get(tng, ind[i], temp_name, 256);
392 if (stat != TNG_SUCCESS)
396 stat = tng_atom_type_of_particle_nr_get(tng, ind[i], temp_type, 256);
397 if (stat != TNG_SUCCESS)
401 tng_residue_atom_w_id_add(tng, res, temp_name, temp_type, ind[i], &atom);
403 tng_molecule_existing_add(tng, &mol);
405 /* Set the count of the molecule containing the selected atoms to 1 and all
406 * other molecules to 0 */
407 tng_molecule_cnt_set(tng, mol, 1);
408 tng_num_molecule_types_get(tng, &nMols);
409 for (gmx_int64_t k = 0; k < nMols; k++)
411 tng_molecule_of_index_get(tng, k, &iterMol);
416 tng_molecule_cnt_set(tng, iterMol, 0);
419 GMX_UNUSED_VALUE(tng);
420 GMX_UNUSED_VALUE(nind);
421 GMX_UNUSED_VALUE(ind);
422 GMX_UNUSED_VALUE(name);
426 /* TODO: If/when TNG acquires the ability to copy data blocks without
427 * uncompressing them, then this implemenation should be reconsidered.
428 * Ideally, gmx trjconv -f a.tng -o b.tng -b 10 -e 20 would be fast
429 * and lose no information. */
430 gmx_bool gmx_read_next_tng_frame(tng_trajectory_t input,
432 gmx_int64_t *requestedIds,
437 tng_function_status stat;
438 gmx_int64_t numberOfAtoms = -1, frameNumber = -1;
439 gmx_int64_t nBlocks, blockId, *blockIds = NULL, codecId;
442 double frameTime = -1.0;
443 int size, blockDependency;
445 const int defaultNumIds = 5;
446 static gmx_int64_t fallbackRequestedIds[defaultNumIds] =
448 TNG_TRAJ_BOX_SHAPE, TNG_TRAJ_POSITIONS,
449 TNG_TRAJ_VELOCITIES, TNG_TRAJ_FORCES,
464 /* If no specific IDs were requested read all block types that can
465 * currently be interpreted */
466 if (!requestedIds || numRequestedIds == 0)
468 numRequestedIds = defaultNumIds;
469 requestedIds = fallbackRequestedIds;
472 stat = tng_num_particles_get(input, &numberOfAtoms);
473 if (stat != TNG_SUCCESS)
475 gmx_file("Cannot determine number of atoms from TNG file.");
477 fr->natoms = numberOfAtoms;
479 if (!gmx_get_tng_data_block_types_of_next_frame(input,
495 for (gmx_int64_t i = 0; i < nBlocks; i++)
497 blockId = blockIds[i];
498 tng_data_block_dependency_get(input, blockId, &blockDependency);
499 if (blockDependency & TNG_PARTICLE_DEPENDENT)
501 stat = tng_util_particle_data_next_frame_read(input,
510 stat = tng_util_non_particle_data_next_frame_read(input,
517 if (stat == TNG_CRITICAL)
519 gmx_file("Cannot read positions from TNG file.");
522 else if (stat == TNG_FAILURE)
528 case TNG_TRAJ_BOX_SHAPE:
532 size = sizeof(gmx_int64_t);
535 size = sizeof(float);
537 case TNG_DOUBLE_DATA:
538 size = sizeof(double);
541 gmx_incons("Illegal datatype of box shape values!");
543 for (int i = 0; i < DIM; i++)
545 convert_array_to_real_array((char *)(values) + size * i * DIM,
547 getDistanceScaleFactor(input),
554 case TNG_TRAJ_POSITIONS:
555 srenew(fr->x, fr->natoms);
556 convert_array_to_real_array(values,
558 getDistanceScaleFactor(input),
563 tng_util_frame_current_compression_get(input, blockId, &codecId, &prec);
564 /* This must be updated if/when more lossy compression methods are added */
565 if (codecId == TNG_TNG_COMPRESSION)
571 case TNG_TRAJ_VELOCITIES:
572 srenew(fr->v, fr->natoms);
573 convert_array_to_real_array(values,
575 getDistanceScaleFactor(input),
580 tng_util_frame_current_compression_get(input, blockId, &codecId, &prec);
581 /* This must be updated if/when more lossy compression methods are added */
582 if (codecId == TNG_TNG_COMPRESSION)
588 case TNG_TRAJ_FORCES:
589 srenew(fr->f, fr->natoms);
590 convert_array_to_real_array(values,
592 getDistanceScaleFactor(input),
602 fr->lambda = (*(float *)values);
604 case TNG_DOUBLE_DATA:
605 fr->lambda = (*(double *)values);
608 gmx_incons("Illegal datatype lambda value!");
613 gmx_warning("Illegal block type! Currently GROMACS tools can only handle certain data types. Skipping block.");
615 /* values does not have to be freed before reading next frame. It will
616 * be reallocated if it is not NULL. */
619 fr->step = (int) frameNumber;
621 // Convert the time to ps
622 fr->time = frameTime / PICO;
625 /* values must be freed before leaving this function */
630 GMX_UNUSED_VALUE(input);
631 GMX_UNUSED_VALUE(fr);
632 GMX_UNUSED_VALUE(requestedIds);
633 GMX_UNUSED_VALUE(numRequestedIds);
638 void gmx_print_tng_molecule_system(tng_trajectory_t input,
642 gmx_int64_t nMolecules, nChains, nResidues, nAtoms, *molCntList;
643 tng_molecule_t molecule;
645 tng_residue_t residue;
647 char str[256], varNAtoms;
649 tng_num_molecule_types_get(input, &nMolecules);
650 tng_molecule_cnt_list_get(input, &molCntList);
651 /* Can the number of particles change in the trajectory or is it constant? */
652 tng_num_particles_variable_get(input, &varNAtoms);
654 for (gmx_int64_t i = 0; i < nMolecules; i++)
656 tng_molecule_of_index_get(input, i, &molecule);
657 tng_molecule_name_get(input, molecule, str, 256);
658 if (varNAtoms == TNG_CONSTANT_N_ATOMS)
660 if ((int)molCntList[i] == 0)
664 fprintf(stream, "Molecule: %s, count: %d\n", str, (int)molCntList[i]);
668 fprintf(stream, "Molecule: %s\n", str);
670 tng_molecule_num_chains_get(input, molecule, &nChains);
673 for (gmx_int64_t j = 0; j < nChains; j++)
675 tng_molecule_chain_of_index_get(input, molecule, j, &chain);
676 tng_chain_name_get(input, chain, str, 256);
677 fprintf(stream, "\tChain: %s\n", str);
678 tng_chain_num_residues_get(input, chain, &nResidues);
679 for (gmx_int64_t k = 0; k < nResidues; k++)
681 tng_chain_residue_of_index_get(input, chain, k, &residue);
682 tng_residue_name_get(input, residue, str, 256);
683 fprintf(stream, "\t\tResidue: %s\n", str);
684 tng_residue_num_atoms_get(input, residue, &nAtoms);
685 for (gmx_int64_t l = 0; l < nAtoms; l++)
687 tng_residue_atom_of_index_get(input, residue, l, &atom);
688 tng_atom_name_get(input, atom, str, 256);
689 fprintf(stream, "\t\t\tAtom: %s", str);
690 tng_atom_type_get(input, atom, str, 256);
691 fprintf(stream, " (%s)\n", str);
696 /* It is possible to have a molecule without chains, in which case
697 * residues in the molecule can be iterated through without going
701 tng_molecule_num_residues_get(input, molecule, &nResidues);
704 for (gmx_int64_t k = 0; k < nResidues; k++)
706 tng_molecule_residue_of_index_get(input, molecule, k, &residue);
707 tng_residue_name_get(input, residue, str, 256);
708 fprintf(stream, "\t\tResidue: %s\n", str);
709 tng_residue_num_atoms_get(input, residue, &nAtoms);
710 for (gmx_int64_t l = 0; l < nAtoms; l++)
712 tng_residue_atom_of_index_get(input, residue, l, &atom);
713 tng_atom_name_get(input, atom, str, 256);
714 fprintf(stream, "\t\t\tAtom: %s", str);
715 tng_atom_type_get(input, atom, str, 256);
716 fprintf(stream, " (%s)\n", str);
722 tng_molecule_num_atoms_get(input, molecule, &nAtoms);
723 for (gmx_int64_t l = 0; l < nAtoms; l++)
725 tng_molecule_atom_of_index_get(input, molecule, l, &atom);
726 tng_atom_name_get(input, atom, str, 256);
727 fprintf(stream, "\t\t\tAtom: %s", str);
728 tng_atom_type_get(input, atom, str, 256);
729 fprintf(stream, " (%s)\n", str);
735 GMX_UNUSED_VALUE(input);
736 GMX_UNUSED_VALUE(stream);
740 gmx_bool gmx_get_tng_data_block_types_of_next_frame(tng_trajectory_t input,
743 gmx_int64_t *requestedIds,
744 gmx_int64_t *nextFrame,
745 gmx_int64_t *nBlocks,
746 gmx_int64_t **blockIds)
749 tng_function_status stat;
751 stat = tng_util_trajectory_next_frame_present_data_blocks_find(input, frame,
752 nRequestedIds, requestedIds,
756 if (stat == TNG_CRITICAL)
758 gmx_file("Cannot read TNG file. Cannot find data blocks of next frame.");
760 else if (stat == TNG_FAILURE)
766 GMX_UNUSED_VALUE(input);
767 GMX_UNUSED_VALUE(frame);
768 GMX_UNUSED_VALUE(nRequestedIds);
769 GMX_UNUSED_VALUE(requestedIds);
770 GMX_UNUSED_VALUE(nextFrame);
771 GMX_UNUSED_VALUE(nBlocks);
772 GMX_UNUSED_VALUE(blockIds);
777 gmx_bool gmx_get_tng_data_next_frame_of_block_type(tng_trajectory_t input,
780 gmx_int64_t *frameNumber,
782 gmx_int64_t *nValuesPerFrame,
790 tng_function_status stat;
797 stat = tng_data_block_name_get(input, blockId, name, maxLen);
798 if (stat != TNG_SUCCESS)
800 gmx_file("Cannot read next frame of TNG file");
802 stat = tng_data_block_dependency_get(input, blockId, &blockDependency);
803 if (stat != TNG_SUCCESS)
805 gmx_file("Cannot read next frame of TNG file");
807 if (blockDependency & TNG_PARTICLE_DEPENDENT)
809 tng_num_particles_get(input, nAtoms);
810 stat = tng_util_particle_data_next_frame_read(input,
819 *nAtoms = 1; /* There are not actually any atoms, but it is used for
821 stat = tng_util_non_particle_data_next_frame_read(input,
828 if (stat == TNG_CRITICAL)
830 gmx_file("Cannot read next frame of TNG file");
832 if (stat == TNG_FAILURE)
838 stat = tng_data_block_num_values_per_frame_get(input, blockId, nValuesPerFrame);
839 if (stat != TNG_SUCCESS)
841 gmx_file("Cannot read next frame of TNG file");
843 snew(*values, sizeof(real) * *nValuesPerFrame * *nAtoms);
844 convert_array_to_real_array(data,
846 getDistanceScaleFactor(input),
851 tng_util_frame_current_compression_get(input, blockId, &codecId, &localPrec);
853 /* This must be updated if/when more lossy compression methods are added */
854 if (codecId != TNG_TNG_COMPRESSION)
866 GMX_UNUSED_VALUE(input);
867 GMX_UNUSED_VALUE(blockId);
868 GMX_UNUSED_VALUE(values);
869 GMX_UNUSED_VALUE(frameNumber);
870 GMX_UNUSED_VALUE(frameTime);
871 GMX_UNUSED_VALUE(nValuesPerFrame);
872 GMX_UNUSED_VALUE(nAtoms);
873 GMX_UNUSED_VALUE(prec);
874 GMX_UNUSED_VALUE(name);
875 GMX_UNUSED_VALUE(maxLen);
876 GMX_UNUSED_VALUE(bOK);