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"
47 #include "tng/tng_io.h"
50 #include "gromacs/legacyheaders/types/atoms.h"
51 #include "gromacs/legacyheaders/physics.h"
52 #include "gromacs/legacyheaders/gmx_fatal.h"
54 #include "gromacs/utility/common.h"
55 #include "gromacs/utility/smalloc.h"
57 void gmx_prepare_tng_writing(const char *filename,
59 tng_trajectory_t *input,
60 tng_trajectory_t *output,
62 const gmx_mtop_t *mtop,
64 const char *indexGroupName)
67 /* FIXME after 5.0: Currently only standard block types are read */
68 const int defaultNumIds = 5;
69 static gmx_int64_t fallbackIds[defaultNumIds] =
71 TNG_TRAJ_BOX_SHAPE, TNG_TRAJ_POSITIONS,
72 TNG_TRAJ_VELOCITIES, TNG_TRAJ_FORCES,
75 static char fallbackNames[defaultNumIds][32] =
77 "BOX SHAPE", "POSITIONS", "VELOCITIES",
81 typedef tng_function_status (*set_writing_interval_func_pointer)(tng_trajectory_t,
89 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_double_set;
91 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_set;
94 gmx_tng_open(filename, mode, output);
96 /* Do we have an input file in TNG format? If so, then there's
97 more data we can copy over, rather than having to improvise. */
100 /* Set parameters (compression, time per frame, molecule
101 * information, number of frames per frame set and writing
102 * intervals of positions, box shape and lambdas) of the
103 * output tng container based on their respective values int
104 * the input tng container */
105 double time, compression_precision;
106 gmx_int64_t n_frames_per_frame_set, interval = -1;
108 tng_compression_precision_get(*input, &compression_precision);
109 tng_compression_precision_set(*output, compression_precision);
110 // TODO make this configurable in a future version
111 char compression_type = TNG_TNG_COMPRESSION;
113 tng_molecule_system_copy(*input, *output);
115 tng_time_per_frame_get(*input, &time);
116 tng_time_per_frame_set(*output, time);
118 tng_num_frames_per_frame_set_get(*input, &n_frames_per_frame_set);
119 tng_num_frames_per_frame_set_set(*output, n_frames_per_frame_set);
121 for (int i = 0; i < defaultNumIds; i++)
123 if (tng_data_get_stride_length(*input, fallbackIds[i], -1, &interval)
126 switch (fallbackIds[i])
128 case TNG_TRAJ_POSITIONS:
129 case TNG_TRAJ_VELOCITIES:
130 set_writing_interval(*output, interval, 3, fallbackIds[i],
131 fallbackNames[i], TNG_PARTICLE_BLOCK_DATA,
134 case TNG_TRAJ_FORCES:
135 set_writing_interval(*output, interval, 3, fallbackIds[i],
136 fallbackNames[i], TNG_PARTICLE_BLOCK_DATA,
137 TNG_GZIP_COMPRESSION);
139 case TNG_TRAJ_BOX_SHAPE:
140 set_writing_interval(*output, interval, 9, fallbackIds[i],
141 fallbackNames[i], TNG_NON_PARTICLE_BLOCK_DATA,
142 TNG_GZIP_COMPRESSION);
145 set_writing_interval(*output, interval, 1, fallbackIds[i],
146 fallbackNames[i], TNG_NON_PARTICLE_BLOCK_DATA,
147 TNG_GZIP_COMPRESSION);
157 /* TODO after trjconv is modularized: fix this so the user can
158 change precision when they are doing an operation where
159 this makes sense, and not otherwise.
161 char compression = bUseLossyCompression ? TNG_TNG_COMPRESSION : TNG_GZIP_COMPRESSION;
162 gmx_tng_set_compression_precision(*output, ndec2prec(nDecimalsOfPrecision));
164 gmx_tng_add_mtop(*output, mtop);
165 tng_num_frames_per_frame_set_set(*output, 1);
168 if (index && nAtoms > 0)
170 gmx_tng_setup_atom_subgroup(*output, nAtoms, index, indexGroupName);
173 /* If for some reason there are more requested atoms than there are atoms in the
174 * molecular system create a number of implicit atoms (without atom data) to
175 * compensate for that. */
178 tng_implicit_num_particles_set(*output, nAtoms);
181 GMX_UNUSED_VALUE(filename);
182 GMX_UNUSED_VALUE(mode);
183 GMX_UNUSED_VALUE(input);
184 GMX_UNUSED_VALUE(output);
185 GMX_UNUSED_VALUE(nAtoms);
186 GMX_UNUSED_VALUE(mtop);
187 GMX_UNUSED_VALUE(index);
188 GMX_UNUSED_VALUE(indexGroupName);
192 void gmx_write_tng_from_trxframe(tng_trajectory_t output,
199 double timePerFrame = frame->time * PICO / frame->step;
200 tng_time_per_frame_set(output, timePerFrame);
204 natoms = frame->natoms;
206 gmx_fwrite_tng(output,
211 (const rvec *) frame->box,
213 (const rvec *) frame->x,
214 (const rvec *) frame->v,
215 (const rvec *) frame->f);
217 GMX_UNUSED_VALUE(output);
218 GMX_UNUSED_VALUE(frame);
219 GMX_UNUSED_VALUE(natoms);
225 convert_array_to_real_array(void *from,
237 if (sizeof(real) == sizeof(float))
241 memcpy(to, from, nValues * sizeof(real) * nAtoms);
245 for (i = 0; i < nAtoms; i++)
247 for (j = 0; j < nValues; j++)
249 to[i*nValues+j] = (real)((float *)from)[i*nValues+j] * fact;
256 for (i = 0; i < nAtoms; i++)
258 for (j = 0; j < nValues; j++)
260 to[i*nValues+j] = (real)((float *)from)[i*nValues+j] * fact;
266 for (i = 0; i < nAtoms; i++)
268 for (j = 0; j < nValues; j++)
270 to[i*nValues+j] = (real)((gmx_int64_t *)from)[i*nValues+j] * fact;
274 case TNG_DOUBLE_DATA:
275 if (sizeof(real) == sizeof(double))
279 memcpy(to, from, nValues * sizeof(real) * nAtoms);
283 for (i = 0; i < nAtoms; i++)
285 for (j = 0; j < nValues; j++)
287 to[i*nValues+j] = (real)((double *)from)[i*nValues+j] * fact;
294 for (i = 0; i < nAtoms; i++)
296 for (j = 0; j < nValues; j++)
298 to[i*nValues+j] = (real)((double *)from)[i*nValues+j] * fact;
304 gmx_incons("Illegal datatype when converting values to a real array!");
310 static real getDistanceScaleFactor(tng_trajectory_t in)
312 gmx_int64_t exp = -1;
313 real distanceScaleFactor;
315 // TODO Hopefully, TNG 2.0 will do this kind of thing for us
316 tng_distance_unit_exponential_get(in, &exp);
318 // GROMACS expects distances in nm
322 distanceScaleFactor = NANO/NANO;
325 distanceScaleFactor = NANO/ANGSTROM;
328 distanceScaleFactor = pow(10.0, exp + 9.0);
331 return distanceScaleFactor;
335 void gmx_tng_setup_atom_subgroup(tng_trajectory_t tng,
341 gmx_int64_t nAtoms, cnt, nMols;
342 tng_molecule_t mol, iterMol;
346 tng_function_status stat;
348 tng_num_particles_get(tng, &nAtoms);
355 stat = tng_molecule_find(tng, name, -1, &mol);
356 if (stat == TNG_SUCCESS)
358 tng_molecule_num_atoms_get(tng, mol, &nAtoms);
359 tng_molecule_cnt_get(tng, mol, &cnt);
369 if (stat == TNG_FAILURE)
371 /* The indexed atoms are added to one separate molecule. */
372 tng_molecule_alloc(tng, &mol);
373 tng_molecule_name_set(tng, mol, name);
374 tng_molecule_chain_add(tng, mol, "", &chain);
376 for (int i = 0; i < nind; i++)
378 char temp_name[256], temp_type[256];
380 /* Try to retrieve the residue name of the atom */
381 stat = tng_residue_name_of_particle_nr_get(tng, ind[i], temp_name, 256);
382 if (stat != TNG_SUCCESS)
386 /* Check if the molecule of the selection already contains this residue */
387 if (tng_chain_residue_find(tng, chain, temp_name, -1, &res)
390 tng_chain_residue_add(tng, chain, temp_name, &res);
392 /* Try to find the original name and type of the atom */
393 stat = tng_atom_name_of_particle_nr_get(tng, ind[i], temp_name, 256);
394 if (stat != TNG_SUCCESS)
398 stat = tng_atom_type_of_particle_nr_get(tng, ind[i], temp_type, 256);
399 if (stat != TNG_SUCCESS)
403 tng_residue_atom_w_id_add(tng, res, temp_name, temp_type, ind[i], &atom);
405 tng_molecule_existing_add(tng, &mol);
407 /* Set the count of the molecule containing the selected atoms to 1 and all
408 * other molecules to 0 */
409 tng_molecule_cnt_set(tng, mol, 1);
410 tng_num_molecule_types_get(tng, &nMols);
411 for (gmx_int64_t k = 0; k < nMols; k++)
413 tng_molecule_of_index_get(tng, k, &iterMol);
418 tng_molecule_cnt_set(tng, iterMol, 0);
421 GMX_UNUSED_VALUE(tng);
422 GMX_UNUSED_VALUE(nind);
423 GMX_UNUSED_VALUE(ind);
424 GMX_UNUSED_VALUE(name);
428 /* TODO: If/when TNG acquires the ability to copy data blocks without
429 * uncompressing them, then this implemenation should be reconsidered.
430 * Ideally, gmx trjconv -f a.tng -o b.tng -b 10 -e 20 would be fast
431 * and lose no information. */
432 gmx_bool gmx_read_next_tng_frame(tng_trajectory_t input,
434 gmx_int64_t *requestedIds,
439 tng_function_status stat;
440 gmx_int64_t numberOfAtoms = -1, frameNumber = -1;
441 gmx_int64_t nBlocks, blockId, *blockIds = NULL, codecId;
444 double frameTime = -1.0;
445 int size, blockDependency;
447 const int defaultNumIds = 5;
448 static gmx_int64_t fallbackRequestedIds[defaultNumIds] =
450 TNG_TRAJ_BOX_SHAPE, TNG_TRAJ_POSITIONS,
451 TNG_TRAJ_VELOCITIES, TNG_TRAJ_FORCES,
466 /* If no specific IDs were requested read all block types that can
467 * currently be interpreted */
468 if (!requestedIds || numRequestedIds == 0)
470 numRequestedIds = defaultNumIds;
471 requestedIds = fallbackRequestedIds;
474 stat = tng_num_particles_get(input, &numberOfAtoms);
475 if (stat != TNG_SUCCESS)
477 gmx_file("Cannot determine number of atoms from TNG file.");
479 fr->natoms = numberOfAtoms;
481 if (!gmx_get_tng_data_block_types_of_next_frame(input,
497 for (gmx_int64_t i = 0; i < nBlocks; i++)
499 blockId = blockIds[i];
500 tng_data_block_dependency_get(input, blockId, &blockDependency);
501 if (blockDependency & TNG_PARTICLE_DEPENDENT)
503 stat = tng_util_particle_data_next_frame_read(input,
512 stat = tng_util_non_particle_data_next_frame_read(input,
519 if (stat == TNG_CRITICAL)
521 gmx_file("Cannot read positions from TNG file.");
524 else if (stat == TNG_FAILURE)
530 case TNG_TRAJ_BOX_SHAPE:
534 size = sizeof(gmx_int64_t);
537 size = sizeof(float);
539 case TNG_DOUBLE_DATA:
540 size = sizeof(double);
543 size = 0; /* Just to make the compiler happy. */
544 gmx_incons("Illegal datatype of box shape values!");
546 for (int i = 0; i < DIM; i++)
548 convert_array_to_real_array((char *)(values) + size * i * DIM,
550 getDistanceScaleFactor(input),
557 case TNG_TRAJ_POSITIONS:
558 srenew(fr->x, 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_VELOCITIES:
575 srenew(fr->v, fr->natoms);
576 convert_array_to_real_array(values,
578 getDistanceScaleFactor(input),
583 tng_util_frame_current_compression_get(input, blockId, &codecId, &prec);
584 /* This must be updated if/when more lossy compression methods are added */
585 if (codecId == TNG_TNG_COMPRESSION)
591 case TNG_TRAJ_FORCES:
592 srenew(fr->f, fr->natoms);
593 convert_array_to_real_array(values,
595 getDistanceScaleFactor(input),
605 fr->lambda = (*(float *)values);
607 case TNG_DOUBLE_DATA:
608 fr->lambda = (*(double *)values);
611 gmx_incons("Illegal datatype lambda value!");
616 gmx_warning("Illegal block type! Currently GROMACS tools can only handle certain data types. Skipping block.");
618 /* values does not have to be freed before reading next frame. It will
619 * be reallocated if it is not NULL. */
622 fr->step = (int) frameNumber;
624 // Convert the time to ps
625 fr->time = frameTime / PICO;
628 /* values must be freed before leaving this function */
633 GMX_UNUSED_VALUE(input);
634 GMX_UNUSED_VALUE(fr);
635 GMX_UNUSED_VALUE(requestedIds);
636 GMX_UNUSED_VALUE(numRequestedIds);
641 void gmx_print_tng_molecule_system(tng_trajectory_t input,
645 gmx_int64_t nMolecules, nChains, nResidues, nAtoms, *molCntList;
646 tng_molecule_t molecule;
648 tng_residue_t residue;
650 char str[256], varNAtoms;
652 tng_num_molecule_types_get(input, &nMolecules);
653 tng_molecule_cnt_list_get(input, &molCntList);
654 /* Can the number of particles change in the trajectory or is it constant? */
655 tng_num_particles_variable_get(input, &varNAtoms);
657 for (gmx_int64_t i = 0; i < nMolecules; i++)
659 tng_molecule_of_index_get(input, i, &molecule);
660 tng_molecule_name_get(input, molecule, str, 256);
661 if (varNAtoms == TNG_CONSTANT_N_ATOMS)
663 if ((int)molCntList[i] == 0)
667 fprintf(stream, "Molecule: %s, count: %d\n", str, (int)molCntList[i]);
671 fprintf(stream, "Molecule: %s\n", str);
673 tng_molecule_num_chains_get(input, molecule, &nChains);
676 for (gmx_int64_t j = 0; j < nChains; j++)
678 tng_molecule_chain_of_index_get(input, molecule, j, &chain);
679 tng_chain_name_get(input, chain, str, 256);
680 fprintf(stream, "\tChain: %s\n", str);
681 tng_chain_num_residues_get(input, chain, &nResidues);
682 for (gmx_int64_t k = 0; k < nResidues; k++)
684 tng_chain_residue_of_index_get(input, chain, k, &residue);
685 tng_residue_name_get(input, residue, str, 256);
686 fprintf(stream, "\t\tResidue: %s\n", str);
687 tng_residue_num_atoms_get(input, residue, &nAtoms);
688 for (gmx_int64_t l = 0; l < nAtoms; l++)
690 tng_residue_atom_of_index_get(input, residue, l, &atom);
691 tng_atom_name_get(input, atom, str, 256);
692 fprintf(stream, "\t\t\tAtom: %s", str);
693 tng_atom_type_get(input, atom, str, 256);
694 fprintf(stream, " (%s)\n", str);
699 /* It is possible to have a molecule without chains, in which case
700 * residues in the molecule can be iterated through without going
704 tng_molecule_num_residues_get(input, molecule, &nResidues);
707 for (gmx_int64_t k = 0; k < nResidues; k++)
709 tng_molecule_residue_of_index_get(input, molecule, k, &residue);
710 tng_residue_name_get(input, residue, str, 256);
711 fprintf(stream, "\t\tResidue: %s\n", str);
712 tng_residue_num_atoms_get(input, residue, &nAtoms);
713 for (gmx_int64_t l = 0; l < nAtoms; l++)
715 tng_residue_atom_of_index_get(input, residue, l, &atom);
716 tng_atom_name_get(input, atom, str, 256);
717 fprintf(stream, "\t\t\tAtom: %s", str);
718 tng_atom_type_get(input, atom, str, 256);
719 fprintf(stream, " (%s)\n", str);
725 tng_molecule_num_atoms_get(input, molecule, &nAtoms);
726 for (gmx_int64_t l = 0; l < nAtoms; l++)
728 tng_molecule_atom_of_index_get(input, molecule, l, &atom);
729 tng_atom_name_get(input, atom, str, 256);
730 fprintf(stream, "\t\t\tAtom: %s", str);
731 tng_atom_type_get(input, atom, str, 256);
732 fprintf(stream, " (%s)\n", str);
738 GMX_UNUSED_VALUE(input);
739 GMX_UNUSED_VALUE(stream);
743 gmx_bool gmx_get_tng_data_block_types_of_next_frame(tng_trajectory_t input,
746 gmx_int64_t *requestedIds,
747 gmx_int64_t *nextFrame,
748 gmx_int64_t *nBlocks,
749 gmx_int64_t **blockIds)
752 tng_function_status stat;
754 stat = tng_util_trajectory_next_frame_present_data_blocks_find(input, frame,
755 nRequestedIds, requestedIds,
759 if (stat == TNG_CRITICAL)
761 gmx_file("Cannot read TNG file. Cannot find data blocks of next frame.");
763 else if (stat == TNG_FAILURE)
769 GMX_UNUSED_VALUE(input);
770 GMX_UNUSED_VALUE(frame);
771 GMX_UNUSED_VALUE(nRequestedIds);
772 GMX_UNUSED_VALUE(requestedIds);
773 GMX_UNUSED_VALUE(nextFrame);
774 GMX_UNUSED_VALUE(nBlocks);
775 GMX_UNUSED_VALUE(blockIds);
780 gmx_bool gmx_get_tng_data_next_frame_of_block_type(tng_trajectory_t input,
783 gmx_int64_t *frameNumber,
785 gmx_int64_t *nValuesPerFrame,
793 tng_function_status stat;
800 stat = tng_data_block_name_get(input, blockId, name, maxLen);
801 if (stat != TNG_SUCCESS)
803 gmx_file("Cannot read next frame of TNG file");
805 stat = tng_data_block_dependency_get(input, blockId, &blockDependency);
806 if (stat != TNG_SUCCESS)
808 gmx_file("Cannot read next frame of TNG file");
810 if (blockDependency & TNG_PARTICLE_DEPENDENT)
812 tng_num_particles_get(input, nAtoms);
813 stat = tng_util_particle_data_next_frame_read(input,
822 *nAtoms = 1; /* There are not actually any atoms, but it is used for
824 stat = tng_util_non_particle_data_next_frame_read(input,
831 if (stat == TNG_CRITICAL)
833 gmx_file("Cannot read next frame of TNG file");
835 if (stat == TNG_FAILURE)
841 stat = tng_data_block_num_values_per_frame_get(input, blockId, nValuesPerFrame);
842 if (stat != TNG_SUCCESS)
844 gmx_file("Cannot read next frame of TNG file");
846 snew(*values, sizeof(real) * *nValuesPerFrame * *nAtoms);
847 convert_array_to_real_array(data,
849 getDistanceScaleFactor(input),
854 tng_util_frame_current_compression_get(input, blockId, &codecId, &localPrec);
856 /* This must be updated if/when more lossy compression methods are added */
857 if (codecId != TNG_TNG_COMPRESSION)
869 GMX_UNUSED_VALUE(input);
870 GMX_UNUSED_VALUE(blockId);
871 GMX_UNUSED_VALUE(values);
872 GMX_UNUSED_VALUE(frameNumber);
873 GMX_UNUSED_VALUE(frameTime);
874 GMX_UNUSED_VALUE(nValuesPerFrame);
875 GMX_UNUSED_VALUE(nAtoms);
876 GMX_UNUSED_VALUE(prec);
877 GMX_UNUSED_VALUE(name);
878 GMX_UNUSED_VALUE(maxLen);
879 GMX_UNUSED_VALUE(bOK);