Fix GPU X buffer ops with empty domain
[gromacs.git] / src / gromacs / coordinateio / coordinatefile.cpp
blob8d9c81a636b6971959e396ce0fa1e6b3a0fb6b9d
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2019, 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 /*!\internal
36 * \file
37 * \brief
38 * Implements methods from coordinatefile.h
40 * \author Paul Bauer <paul.bauer.q@gmail.com>
41 * \ingroup module_coordinateio
45 #include "gmxpre.h"
47 #include "coordinatefile.h"
49 #include <algorithm>
51 #include "gromacs/coordinateio/outputadapters.h"
52 #include "gromacs/coordinateio/requirements.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/trajectory/trajectoryframe.h"
55 #include "gromacs/utility/exceptions.h"
57 namespace gmx
60 /*!\brief
61 * Get the internal file type from the \p filename.
63 * \param[in] filename Filename of output file.
64 * \throws InvalidInputError When unable to work on an emoty file name.
65 * \returns integer value of file type.
67 static int getFileType(const std::string &filename)
69 int filetype = efNR;
70 if (!filename.empty())
72 filetype = fn2ftp(filename.c_str());
74 else
76 GMX_THROW(InvalidInputError("Can not open file with an empty name"));
78 return filetype;
81 /*!\brief
82 * Get the flag representing the requirements for a given file output.
84 * Also checks if the supplied topology is sufficient through the pointer
85 * to \p mtop.
87 * \param[in] filetype Internal file type used to check requirements.
88 * \throws InvalidInputError When encountering an invalid file type.
89 * \returns Requirements represent by the bitmask in the return type.
91 static unsigned long getSupportedOutputAdapters(int filetype)
93 unsigned long supportedOutputAdapters = 0;
94 supportedOutputAdapters |= convertFlag(CoordinateFileFlags::Base);
95 switch (filetype)
97 case (efTNG):
98 supportedOutputAdapters |= (convertFlag(CoordinateFileFlags::RequireForceOutput) |
99 convertFlag(CoordinateFileFlags::RequireVelocityOutput) |
100 convertFlag(CoordinateFileFlags::RequireAtomConnections) |
101 convertFlag(CoordinateFileFlags::RequireAtomInformation) |
102 convertFlag(CoordinateFileFlags::RequireChangedOutputPrecision));
103 break;
104 case (efPDB):
105 supportedOutputAdapters |= (convertFlag(CoordinateFileFlags::RequireAtomConnections) |
106 convertFlag(CoordinateFileFlags::RequireAtomInformation));
107 break;
108 case (efGRO):
109 supportedOutputAdapters |= (convertFlag(CoordinateFileFlags::RequireAtomInformation) |
110 convertFlag(CoordinateFileFlags::RequireVelocityOutput));
111 break;
112 case (efTRR):
113 supportedOutputAdapters |= (convertFlag(CoordinateFileFlags::RequireForceOutput) |
114 convertFlag(CoordinateFileFlags::RequireVelocityOutput));
115 break;
116 case (efXTC):
117 supportedOutputAdapters |= (convertFlag(CoordinateFileFlags::RequireChangedOutputPrecision));
118 break;
119 case (efG96):
120 break;
121 default:
122 GMX_THROW(InvalidInputError("Invalid file type"));
124 return supportedOutputAdapters;
127 /*! \brief
128 * Creates a new container object with the user requested IOutputAdapter derived
129 * methods attached to it.
131 * \param[in] requirements Specifications for modules to add.
132 * \param[in] atoms Local copy of atom information to use.
133 * \param[in] sel Selection to use for choosing atoms to write out.
134 * \param[in] abilities Specifications for what the output method can do.
135 * \returns New container for IoutputAdapter derived methods.
137 static OutputAdapterContainer
138 addOutputAdapters(const OutputRequirements &requirements,
139 AtomsDataPtr atoms,
140 const Selection &sel,
141 unsigned long abilities)
143 OutputAdapterContainer output(abilities);
145 /* An adapter only gets added if the user has specified a non-default
146 * behaviour. In most cases, this behaviour is passive, meaning that
147 * output gets written if it exists in the input and if the output
148 * type supports it.
150 if (requirements.velocity != ChangeSettingType::PreservedIfPresent)
152 output.addAdapter(
153 std::make_unique<SetVelocities>(requirements.velocity),
154 CoordinateFileFlags::RequireVelocityOutput);
156 if (requirements.force != ChangeSettingType::PreservedIfPresent)
158 output.addAdapter(
159 std::make_unique<SetForces>(requirements.force),
160 CoordinateFileFlags::RequireForceOutput);
162 if (requirements.precision != ChangeFrameInfoType::PreservedIfPresent)
164 output.addAdapter(
165 std::make_unique<SetPrecision>(requirements.prec),
166 CoordinateFileFlags::RequireChangedOutputPrecision);
168 if (requirements.atoms != ChangeAtomsType::PreservedIfPresent)
170 output.addAdapter(
171 std::make_unique<SetAtoms>(requirements.atoms,
172 std::move(atoms)),
173 CoordinateFileFlags::RequireAtomInformation);
175 if (requirements.frameTime != ChangeFrameTimeType::PreservedIfPresent)
177 switch (requirements.frameTime)
179 case (ChangeFrameTimeType::StartTime):
180 output.addAdapter(std::make_unique<SetStartTime>(requirements.startTimeValue),
181 CoordinateFileFlags::RequireNewFrameStartTime);
182 break;
183 case (ChangeFrameTimeType::TimeStep):
184 output.addAdapter(std::make_unique<SetTimeStep>(requirements.timeStepValue),
185 CoordinateFileFlags::RequireNewFrameTimeStep);
186 break;
187 case (ChangeFrameTimeType::Both):
188 output.addAdapter(std::make_unique<SetStartTime>(requirements.startTimeValue),
189 CoordinateFileFlags::RequireNewFrameStartTime);
190 output.addAdapter(std::make_unique<SetTimeStep>(requirements.timeStepValue),
191 CoordinateFileFlags::RequireNewFrameTimeStep);
192 break;
193 default:
194 break;
197 if (requirements.box != ChangeFrameInfoType::PreservedIfPresent)
199 output.addAdapter(
200 std::make_unique<SetBox>(requirements.newBox),
201 CoordinateFileFlags::RequireNewBox);
203 if (sel.isValid())
205 output.addAdapter(
206 std::make_unique<OutputSelector>(sel),
207 CoordinateFileFlags::RequireCoordinateSelection);
209 return output;
212 std::unique_ptr<TrajectoryFrameWriter>
213 createTrajectoryFrameWriter(const gmx_mtop_t *top,
214 const Selection &sel,
215 const std::string &filename,
216 AtomsDataPtr atoms,
217 OutputRequirements requirements)
219 /* TODO
220 * Currently the requirements object is expected to be processed and valid,
221 * meaning that e.g. a new box is specified if requested by the option,
222 * or that time values have been set if the corresponding values are set.
223 * This will need to get revisited when the code that builds this object from
224 * the user options gets merged.
226 int filetype = getFileType(filename);
227 unsigned long abilities = getSupportedOutputAdapters(filetype);
229 // first, check if we have a special output format that needs atoms
230 if ((filetype == efPDB) || (filetype == efGRO))
232 if (requirements.atoms == ChangeAtomsType::Never)
234 GMX_THROW(InconsistentInputError("Can not write to PDB or GRO when"
235 "explicitly turning atom information off"));
237 if (requirements.atoms != ChangeAtomsType::AlwaysFromStructure)
239 requirements.atoms = ChangeAtomsType::Always;
242 OutputAdapterContainer outputAdapters = addOutputAdapters(requirements, std::move(atoms),
243 sel, abilities);
245 TrajectoryFrameWriterPointer trajectoryFrameWriter(
246 new TrajectoryFrameWriter(filename,
247 filetype,
248 sel,
249 top,
250 std::move(outputAdapters)));
251 return trajectoryFrameWriter;
255 /*! \brief
256 * Create a deep copy of a t_trxframe \p input into \p copy
258 * When running the analysis tools and changing values with the
259 * outputadapters, a deep copy of the \p input coordinate frame has to be
260 * created first to ensure that the data is not changed if it is needed for other
261 * tools following with analysis later. Therefore, the data is passed
262 * to \p copy by performing a deep copy first.
264 * The method allocates new storage for coordinates of the x, v, and f arrays
265 * in the new coordinate frame. This means that those arrays need to be free'd
266 * after the frame has been processed and been written to disk.
268 * \param[in] input Reference input coordinate frame.
269 * \param[in,out] copy Pointer to new output frame that will receive the deep copy.
270 * \param[in] xvec Pointer to local coordinate storage vector.
271 * \param[in] vvec Pointer to local velocity storage vector.
272 * \param[in] fvec Pointer to local force storage vector.
273 * \param[in] indexvec Pointer to local index storage vector.
275 static void deepCopy_t_trxframe(const t_trxframe &input,
276 t_trxframe *copy,
277 RVec *xvec,
278 RVec *vvec,
279 RVec *fvec,
280 int *indexvec)
282 copy->not_ok = input.not_ok;
283 copy->bStep = input.bStep;
284 copy->bTime = input.bTime;
285 copy->bLambda = input.bLambda;
286 copy->bFepState = input.bFepState;
287 copy->bAtoms = input.bAtoms;
288 copy->bPrec = input.bPrec;
289 copy->bX = input.bX;
290 copy->bV = input.bV;
291 copy->bF = input.bF;
292 copy->bBox = input.bBox;
293 copy->bDouble = input.bDouble;
294 copy->natoms = input.natoms;
295 copy->step = input.step;
296 copy->time = input.time;
297 copy->lambda = input.lambda;
298 copy->fep_state = input.fep_state;
299 if (input.bAtoms)
301 copy->atoms = input.atoms;
303 copy->prec = input.prec;
304 if (copy->bX)
306 copy->x = as_rvec_array(xvec);
308 if (copy->bV)
310 copy->v = as_rvec_array(vvec);
312 if (copy->bF)
314 copy->f = as_rvec_array(fvec);
316 if (input.index)
318 copy->index = indexvec;
320 else
322 copy->index = nullptr;
324 for (int i = 0; i < copy->natoms; i++)
326 if (copy->bX)
328 copy_rvec(input.x[i], copy->x[i]);
330 if (copy->bV)
332 copy_rvec(input.v[i], copy->v[i]);
334 if (copy->bF)
336 copy_rvec(input.f[i], copy->f[i]);
338 if (input.index)
340 copy->index[i] = input.index[i];
343 copy_mat(input.box, copy->box);
344 copy->bPBC = input.bPBC;
345 copy->ePBC = input.ePBC;
348 /*! \brief
349 * Method to open TNG file.
351 * Only need extra method to open this kind of file as it may need access to
352 * a Selection \p sel, if it is valid. Otherwise atom indices will be taken
353 * from the topology \p mtop.
355 * \param[in] name Name of the output file.
356 * \param[in] sel Reference to selection.
357 * \param[in] mtop Pointer to topology, tested before that it is valid.
358 * \todo Those should be methods in a replacement for t_trxstatus instead.
360 static t_trxstatus *openTNG(const std::string &name, const Selection &sel, const gmx_mtop_t *mtop)
362 const char *filemode = "w";
363 if (sel.isValid())
365 GMX_ASSERT(sel.hasOnlyAtoms(), "Can only work with selections consisting out of atoms");
366 return trjtools_gmx_prepare_tng_writing(name.c_str(),
367 filemode[0],
368 nullptr, //infile_, //how to get the input file here?
369 nullptr,
370 sel.atomCount(),
371 mtop,
372 sel.atomIndices(),
373 sel.name());
375 else
377 return trjtools_gmx_prepare_tng_writing(name.c_str(),
378 filemode[0],
379 nullptr, //infile_, //how to get the input file here?
380 nullptr,
381 mtop->natoms,
382 mtop,
383 get_atom_index(mtop),
384 "System");
388 TrajectoryFileOpener::~TrajectoryFileOpener()
390 close_trx(outputFile_);
393 t_trxstatus *TrajectoryFileOpener::outputFile()
395 if (outputFile_ == nullptr)
397 const char *filemode = "w";
398 switch (filetype_)
400 case (efTNG):
401 outputFile_ = openTNG(outputFileName_,
402 sel_,
403 mtop_);
404 break;
405 case (efPDB):
406 case (efGRO):
407 case (efTRR):
408 case (efXTC):
409 case (efG96):
410 outputFile_ = open_trx(outputFileName_.c_str(), filemode);
411 break;
412 default:
413 GMX_THROW(InvalidInputError("Invalid file type"));
416 return outputFile_;
419 void
420 TrajectoryFrameWriter::prepareAndWriteFrame(const int framenumber, const t_trxframe &input)
422 if (!outputAdapters_.isEmpty())
424 t_trxframe local;
425 clear_trxframe(&local, true);
426 localX_.resize(input.natoms);
427 localIndex_.resize(input.natoms);
428 if (input.bV)
430 localV_.resize(input.natoms);
432 if (input.bF)
434 localF_.resize(input.natoms);
436 deepCopy_t_trxframe(input, &local, localX_.data(), localV_.data(),
437 localF_.data(), localIndex_.data());
438 for (const auto &outputAdapter : outputAdapters_.getAdapters())
440 if (outputAdapter)
442 outputAdapter->processFrame(framenumber, &local);
445 write_trxframe(file_.outputFile(), &local, nullptr);
447 else
449 write_trxframe(file_.outputFile(), const_cast<t_trxframe *>(&input), nullptr);
453 } // namespace gmx