SYCL: Avoid using no_init read accessor in rocFFT
[alexxy/gromacs.git] / src / gromacs / selection / selectionoption.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2010,2011,2012,2013,2014 by the GROMACS development team.
5  * Copyright (c) 2015,2018,2019,2020,2021, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36 /*! \file
37  * \brief
38  * Declares gmx::SelectionOption and gmx::SelectionOptionInfo.
39  *
40  * \author Teemu Murtola <teemu.murtola@gmail.com>
41  * \inpublicapi
42  * \ingroup module_selection
43  */
44 #ifndef GMX_SELECTION_SELECTIONOPTION_H
45 #define GMX_SELECTION_SELECTIONOPTION_H
46
47 #include "gromacs/options/abstractoption.h"
48 #include "gromacs/selection/selection.h"
49
50 #include "selectionenums.h"
51
52 namespace gmx
53 {
54
55 class SelectionOptionInfo;
56 class SelectionOptionManager;
57 class SelectionOptionStorage;
58
59 /*! \brief
60  * Specifies an option that provides selection(s).
61  *
62  * Public methods in this class do not throw.
63  *
64  * To use options of this type, SelectionOptionManager must first be added to
65  * the Options collection.  For trajectory analysis tools, the framework takes
66  * care of this.
67  *
68  * \todo
69  * Support for specifying that an option accepts, e.g., two to four selections.
70  * Currently, only a fixed count or any number of selections is possible.
71  * \if internal
72  * In addition to allowing this in OptionTemplate, also SelectionOptionManager
73  * needs to be updated.
74  * \endif
75  *
76  * \inpublicapi
77  * \ingroup module_selection
78  */
79 class SelectionOption : public OptionTemplate<Selection, SelectionOption>
80 {
81 public:
82     //! OptionInfo subclass corresponding to this option type.
83     typedef SelectionOptionInfo InfoType;
84
85     //! Initializes an option with the given name.
86     explicit SelectionOption(const char* name) :
87         MyBase(name), defaultText_(""), selectionFlags_(efSelection_DisallowEmpty)
88     {
89     }
90
91     /*! \brief
92      * Request velocity evaluation for output positions.
93      *
94      * \see Selection::setEvaluateVelocities()
95      */
96     MyClass& evaluateVelocities()
97     {
98         selectionFlags_.set(efSelection_EvaluateVelocities);
99         return me();
100     }
101     /*! \brief
102      * Request force evaluation for output positions.
103      *
104      * \see Selection::setEvaluateForces()
105      */
106     MyClass& evaluateForces()
107     {
108         selectionFlags_.set(efSelection_EvaluateForces);
109         return me();
110     }
111     /*! \brief
112      * Only accept selections that evaluate to atom positions.
113      */
114     MyClass& onlyAtoms()
115     {
116         selectionFlags_.set(efSelection_OnlyAtoms);
117         return me();
118     }
119     /*! \brief
120      * Only accept selections that evaluate to atom positions in sorted order.
121      */
122     MyClass& onlySortedAtoms()
123     {
124         selectionFlags_.set(efSelection_OnlyAtoms);
125         selectionFlags_.set(efSelection_OnlySorted);
126         return me();
127     }
128     /*! \brief
129      * Only accept static selections for this option.
130      */
131     MyClass& onlyStatic()
132     {
133         selectionFlags_.set(efSelection_OnlyStatic);
134         return me();
135     }
136     /*! \brief
137      * Handle dynamic selections for this option with position masks.
138      *
139      * \see Selection
140      * \see SelectionPosition::selected()
141      */
142     MyClass& dynamicMask()
143     {
144         selectionFlags_.set(efSelection_DynamicMask);
145         return me();
146     }
147     /*! \brief
148      * Allow specifying an unconditionally empty selection for this option.
149      *
150      * If this option is not set, selections that are unconditionally empty
151      * (i.e., can never match any atoms) result in errors.
152      * Note that even without this option, it is still possible that a
153      * dynamic selection evaluates to zero atoms for some frames.
154      */
155     MyClass& allowEmpty()
156     {
157         selectionFlags_.clear(efSelection_DisallowEmpty);
158         return me();
159     }
160
161     /*! \brief
162      * Sets default selection text for the option.
163      *
164      * If the option is not set by the user, the provided text is parsed as
165      * the value of the selection.
166      */
167     MyClass& defaultSelectionText(const char* text)
168     {
169         defaultText_ = text;
170         return me();
171     }
172
173 private:
174     // Disable possibility to allow multiple occurrences, since it isn't
175     // implemented.
176     using MyBase::allowMultiple;
177     // Disable default value because it is impossible to provide a
178     // Selection object.
179     using MyBase::defaultValue;
180     using MyBase::defaultValueIfSet;
181
182     AbstractOptionStorage* createStorage(const OptionManagerContainer& managers) const override;
183
184     const char*    defaultText_;
185     SelectionFlags selectionFlags_;
186
187     /*! \brief
188      * Needed to initialize SelectionOptionStorage from this class without
189      * otherwise unnecessary accessors.
190      */
191     friend class SelectionOptionStorage;
192 };
193
194 /*! \brief
195  * Wrapper class for accessing and modifying selection option information.
196  *
197  * Allows changes to a selection option after creation.
198  *
199  * This class provides the necessary interface for changing, e.g., the number
200  * of allowed selections for a selection option after the option has been
201  * created with Options::addOption().  This is needed if the number or other
202  * flags are only known after other options have been parsed.  The main
203  * advantage of this class over custom checks is that if used before
204  * interactive selection prompt, the interactive prompt is updated accordingly.
205  *
206  * When using this class, the option should be initially created with the most
207  * permissive flags, and this class should be used to place restrictions where
208  * appropriate.  Otherwise, values that are provided before adjustments will
209  * need to follow the more strict checks.  In most cases in trajectory analysis
210  * (which is the main use case for selection options), the adjustments should
211  * be done in TrajectoryAnalysisModule::optionsFinished() for them to take
212  * place before interactive selection prompts.
213  *
214  * An instance of this class for a selection option can be obtained with
215  * SelectionOption::getAdjuster() when the option is created.
216  *
217  * Example use:
218  * \code
219    SelectionList sel;
220    Options options("example", "Example options");
221    SelectionOptionInfo *info;
222    info = options.addOption(SelectionOption("sel").storeVector(&sel)
223                                 .multiValue());
224    // < ... assign values to options ...>
225    if ( condition )
226    {
227        // Put limitations on the selections based on the condition,
228        // which can depend on other option values.
229        // Throws if input given so far violates the limitations.
230        info->setValueCount(2);
231        info->setOnlyStatic(true);
232    }
233  * \endcode
234  *
235  * \inpublicapi
236  * \ingroup module_selection
237  */
238 class SelectionOptionInfo : public OptionInfo
239 {
240 public:
241     /*! \brief
242      * Creates option info object for given storage object.
243      *
244      * Does not throw.
245      */
246     explicit SelectionOptionInfo(SelectionOptionStorage* option);
247
248     /*! \brief
249      * Sets the number of selections allowed for the option.
250      *
251      * \param[in] count  Number of allowed selections.
252      * \throws    std::bad_alloc if out of memory.
253      * \throws    InvalidInputError if values have already been provided
254      *      and their count does not match.
255      */
256     void setValueCount(int count);
257
258     /*! \brief
259      * Sets whether this option evaluates velocities for positions.
260      *
261      * \param[in] bEnabled  If true, velocities are evaluated.
262      *
263      * Does not throw.
264      *
265      * \see Selection::setEvaluateVelocities()
266      */
267     void setEvaluateVelocities(bool bEnabled);
268     /*! \brief
269      * Sets whether this option evaluates forces for positions.
270      *
271      * \param[in] bEnabled  If true, forces are evaluated.
272      *
273      * Does not throw.
274      *
275      * \see Selection::setEvaluateForces()
276      */
277     void setEvaluateForces(bool bEnabled);
278     /*! \brief
279      * Sets whether this option accepts positions that come from multiple
280      * atoms.
281      *
282      * \param[in] bEnabled  If true, the option accepts only positions that
283      *      evaluate to atom positions.
284      *
285      * \see SelectionOption::onlyAtoms()
286      */
287     void setOnlyAtoms(bool bEnabled);
288     /*! \brief
289      * Sets whether this option accepts dynamic selections.
290      *
291      * \param[in] bEnabled  If true, the option accepts only static
292      *      selections.
293      * \throws    std::bad_alloc if out of memory.
294      * \throws    InvalidInputError if dynamic selections have already been
295      *      provided.
296      *
297      * Strong exception safety guarantee.
298      *
299      * \see SelectionOption::onlyStatic()
300      */
301     void setOnlyStatic(bool bEnabled);
302     /*! \brief
303      * Sets whether this option uses position masks for dynamic selections.
304      *
305      * \param[in] bEnabled  If true, the position masks are used.
306      *
307      * Does not throw.
308      *
309      * \see SelectionOption::dynamicMask()
310      */
311     void setDynamicMask(bool bEnabled);
312
313 private:
314     SelectionOptionStorage&       option();
315     const SelectionOptionStorage& option() const;
316 };
317
318 } // namespace gmx
319
320 #endif