SYCL: Avoid using no_init read accessor in rocFFT
[alexxy/gromacs.git] / src / gromacs / selection / scanner_internal.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2009-2018, The GROMACS development team.
5  * Copyright (c) 2019,2020, 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 /*! \internal \file
37  * \brief Helper functions for the selection tokenizer.
38  *
39  * This file implements the functions in the headers scanner.h and
40  * scanner_internal.h.
41  *
42  * \author Teemu Murtola <teemu.murtola@gmail.com>
43  * \ingroup module_selection
44  */
45 /*! \cond
46  * \internal \file scanner_flex.h
47  * \brief Generated (from scanner.l) header file by Flex.
48  *
49  * This file contains definitions of functions that are needed in
50  * scanner_internal.cpp.
51  *
52  * \ingroup module_selection
53  * \endcond
54  */
55 #include "gmxpre.h"
56
57 #include "scanner_internal.h"
58
59 #include <cstdlib>
60 #include <cstring>
61
62 #include <string>
63
64 #include "gromacs/utility/cstringutil.h"
65 #include "gromacs/utility/exceptions.h"
66 #include "gromacs/utility/gmxassert.h"
67 #include "gromacs/utility/smalloc.h"
68 #include "gromacs/utility/stringutil.h"
69
70 #include "parser.h"
71 #include "parsetree.h"
72 #include "scanner.h"
73 #include "selectioncollection_impl.h"
74 #include "selelem.h"
75 #include "selmethod.h"
76 #include "symrec.h"
77
78 /* These are defined as macros in the generated scanner_flex.h.
79  * We undefine them here to have them as variable names in the subroutines.
80  * There are other ways of doing this, but this is probably the easiest. */
81 #undef yylval
82 #undef yytext
83 #undef yyleng
84
85 /*! \brief
86  * Handles initialization of method parameter token.
87  */
88 static int init_param_token(YYSTYPE* yylval, gmx_ana_selparam_t* param, bool bBoolNo)
89 {
90     if (bBoolNo)
91     {
92         GMX_RELEASE_ASSERT(param->name != nullptr,
93                            "bBoolNo should only be set for a parameters with a name");
94         snew(yylval->str, strlen(param->name) + 3);
95         yylval->str[0] = 'n';
96         yylval->str[1] = 'o';
97         strcpy(yylval->str + 2, param->name);
98     }
99     else
100     {
101         yylval->str = param->name ? gmx_strdup(param->name) : nullptr;
102     }
103     return PARAM;
104 }
105
106 /*! \brief
107  * Processes a selection method token.
108  */
109 static int init_method_token(YYSTYPE*                          yylval,
110                              YYLTYPE*                          yylloc,
111                              const gmx::SelectionParserSymbol* symbol,
112                              bool                              bPosMod,
113                              gmx_sel_lexer_t*                  state)
114 {
115     gmx_ana_selmethod_t* method = symbol->methodValue();
116     /* If the previous token was not KEYWORD_POS, return EMPTY_POSMOD
117      * before the actual method to work around a limitation in Bison. */
118     if (!bPosMod && method->type != POS_VALUE)
119     {
120         state->nextMethodSymbol = symbol;
121         _gmx_sel_lexer_add_token(yylloc, nullptr, 0, state);
122         return EMPTY_POSMOD;
123     }
124     _gmx_sel_lexer_add_token(yylloc, symbol->name().c_str(), -1, state);
125     yylval->meth = method;
126     if (!(method->flags & SMETH_MODIFIER) && method->nparams == 0)
127     {
128         /* Keyword */
129         switch (method->type)
130         {
131             case INT_VALUE:
132             case REAL_VALUE: state->bMatchOf = true; return KEYWORD_NUMERIC;
133             case STR_VALUE: return KEYWORD_STR;
134             case GROUP_VALUE: return KEYWORD_GROUP;
135             default: GMX_THROW(gmx::InternalError("Unsupported keyword type"));
136         }
137     }
138     else
139     {
140         /* Method with parameters or a modifier */
141         if (method->flags & SMETH_MODIFIER)
142         {
143             /* Remove all methods from the stack */
144             state->msp = -1;
145             if (method->param[1].name == nullptr)
146             {
147                 state->nextparam = &method->param[1];
148             }
149         }
150         else
151         {
152             if (method->param[0].name == nullptr)
153             {
154                 state->nextparam = &method->param[0];
155             }
156         }
157         ++state->msp;
158         if (state->msp >= state->mstack_alloc)
159         {
160             state->mstack_alloc += 10;
161             srenew(state->mstack, state->mstack_alloc);
162         }
163         state->mstack[state->msp] = method;
164         if (method->flags & SMETH_MODIFIER)
165         {
166             return MODIFIER;
167         }
168         switch (method->type)
169         {
170             case INT_VALUE: // Intended fall through
171             case REAL_VALUE: return METHOD_NUMERIC;
172             case POS_VALUE: return METHOD_POS;
173             case GROUP_VALUE: return METHOD_GROUP;
174             default: --state->msp; GMX_THROW(gmx::InternalError("Unsupported method type"));
175         }
176     }
177 }
178
179 int _gmx_sel_lexer_process_pending(YYSTYPE* yylval, YYLTYPE* yylloc, gmx_sel_lexer_t* state)
180 {
181     if (state->nextparam)
182     {
183         gmx_ana_selparam_t* param   = state->nextparam;
184         bool                bBoolNo = state->bBoolNo;
185
186         if (state->neom > 0)
187         {
188             --state->neom;
189             _gmx_sel_lexer_add_token(yylloc, nullptr, 0, state);
190             return END_OF_METHOD;
191         }
192         state->nextparam = nullptr;
193         state->bBoolNo   = false;
194         _gmx_sel_lexer_add_token(yylloc, param->name, -1, state);
195         return init_param_token(yylval, param, bBoolNo);
196     }
197     if (state->prev_pos_kw > 0)
198     {
199         --state->prev_pos_kw;
200     }
201     if (state->nextMethodSymbol)
202     {
203         const gmx::SelectionParserSymbol* symbol = state->nextMethodSymbol;
204         state->nextMethodSymbol                  = nullptr;
205         return init_method_token(yylval, yylloc, symbol, true, state);
206     }
207     return 0;
208 }
209
210 int _gmx_sel_lexer_process_identifier(YYSTYPE* yylval, YYLTYPE* yylloc, char* yytext, size_t yyleng, gmx_sel_lexer_t* state)
211 {
212     /* Check if the identifier matches with a parameter name */
213     if (state->msp >= 0)
214     {
215         gmx_ana_selparam_t* param   = nullptr;
216         bool                bBoolNo = false;
217         int                 sp      = state->msp;
218         while (!param && sp >= 0)
219         {
220             int i;
221             for (i = 0; i < state->mstack[sp]->nparams; ++i)
222             {
223                 /* Skip NULL parameters and too long parameters */
224                 if (state->mstack[sp]->param[i].name == nullptr
225                     || strlen(state->mstack[sp]->param[i].name) > yyleng)
226                 {
227                     continue;
228                 }
229                 if (!strncmp(state->mstack[sp]->param[i].name, yytext, yyleng))
230                 {
231                     param = &state->mstack[sp]->param[i];
232                     break;
233                 }
234                 /* Check separately for a 'no' prefix on boolean parameters */
235                 if (state->mstack[sp]->param[i].val.type == NO_VALUE && yyleng > 2
236                     && yytext[0] == 'n' && yytext[1] == 'o'
237                     && !strncmp(state->mstack[sp]->param[i].name, yytext + 2, yyleng - 2))
238                 {
239                     param   = &state->mstack[sp]->param[i];
240                     bBoolNo = true;
241                     break;
242                 }
243             }
244             if (!param)
245             {
246                 --sp;
247             }
248         }
249         if (param)
250         {
251             if (param->val.type == NO_VALUE && !bBoolNo)
252             {
253                 state->bMatchBool = true;
254             }
255             if (sp < state->msp)
256             {
257                 state->neom      = state->msp - sp - 1;
258                 state->nextparam = param;
259                 state->bBoolNo   = bBoolNo;
260                 return END_OF_METHOD;
261             }
262             _gmx_sel_lexer_add_token(yylloc, param->name, -1, state);
263             return init_param_token(yylval, param, bBoolNo);
264         }
265     }
266
267     /* Check if the identifier matches with a symbol */
268     const gmx::SelectionParserSymbol* symbol = state->sc->symtab->findSymbol(std::string(yytext, yyleng));
269     /* If there is no match, return the token as a string */
270     if (!symbol)
271     {
272         yylval->str = gmx_strndup(yytext, yyleng);
273         _gmx_sel_lexer_add_token(yylloc, yytext, yyleng, state);
274         return IDENTIFIER;
275     }
276     gmx::SelectionParserSymbol::SymbolType symtype = symbol->type();
277     /* For method symbols, we need some extra processing. */
278     if (symtype == gmx::SelectionParserSymbol::MethodSymbol)
279     {
280         return init_method_token(yylval, yylloc, symbol, state->prev_pos_kw > 0, state);
281     }
282     _gmx_sel_lexer_add_token(yylloc, symbol->name().c_str(), -1, state);
283     /* Reserved symbols should have been caught earlier */
284     if (symtype == gmx::SelectionParserSymbol::ReservedSymbol)
285     {
286         GMX_THROW(gmx::InternalError(gmx::formatString(
287                 "Mismatch between tokenizer and reserved symbol table (for '%s')", symbol->name().c_str())));
288     }
289     /* For variable symbols, return the type of the variable value */
290     if (symtype == gmx::SelectionParserSymbol::VariableSymbol)
291     {
292         const gmx::SelectionTreeElementPointer& var = symbol->variableValue();
293         /* Return simple tokens for constant variables */
294         if (var->type == SEL_CONST)
295         {
296             switch (var->v.type)
297             {
298                 case INT_VALUE: yylval->i = var->v.u.i[0]; return TOK_INT;
299                 case REAL_VALUE: yylval->r = var->v.u.r[0]; return TOK_REAL;
300                 case POS_VALUE: break;
301                 default: GMX_THROW(gmx::InternalError("Unsupported variable type"));
302             }
303         }
304         yylval->sel = new gmx::SelectionTreeElementPointer(var);
305         switch (var->v.type)
306         {
307             case INT_VALUE: // Intended fall through
308             case REAL_VALUE: return VARIABLE_NUMERIC;
309             case POS_VALUE: return VARIABLE_POS;
310             case GROUP_VALUE: return VARIABLE_GROUP;
311             default: delete yylval->sel; GMX_THROW(gmx::InternalError("Unsupported variable type"));
312         }
313         /* This position should not be reached. */
314     }
315     /* For position symbols, we need to return KEYWORD_POS, but we also need
316      * some additional handling. */
317     if (symtype == gmx::SelectionParserSymbol::PositionSymbol)
318     {
319         state->bMatchOf    = true;
320         yylval->str        = gmx_strdup(symbol->name().c_str());
321         state->prev_pos_kw = 2;
322         return KEYWORD_POS;
323     }
324     /* Should not be reached */
325     return INVALID;
326 }
327
328 void _gmx_sel_lexer_add_token(YYLTYPE* yylloc, const char* str, int len, gmx_sel_lexer_t* state)
329 {
330     yylloc->startIndex = yylloc->endIndex = state->pselstr.size();
331     /* Do nothing if the string is empty, or if it is a space and there is
332      * no other text yet, or if there already is a space. */
333     if (!str || len == 0 || strlen(str) == 0
334         || (str[0] == ' ' && str[1] == 0 && (state->pselstr.empty() || state->pselstr.back() == ' ')))
335     {
336         return;
337     }
338     if (len < 0)
339     {
340         len = strlen(str);
341     }
342     /* Append the token to the stored string */
343     state->pselstr.append(str, len);
344     yylloc->endIndex = state->pselstr.size();
345 }
346
347 void _gmx_sel_init_lexer(yyscan_t*                       scannerp,
348                          struct gmx_ana_selcollection_t* sc,
349                          gmx::TextWriter*                statusWriter,
350                          int                             maxnr,
351                          bool                            bGroups,
352                          struct gmx_ana_indexgrps_t*     grps)
353 {
354     int rc = _gmx_sel_yylex_init(scannerp);
355     if (rc != 0)
356     {
357         // TODO: Throw a more representative exception.
358         GMX_THROW(gmx::InternalError("Lexer initialization failed"));
359     }
360
361     gmx_sel_lexer_t* state = new gmx_sel_lexer_t;
362
363     state->sc      = sc;
364     state->bGroups = bGroups;
365     state->grps    = grps;
366     state->nexpsel = (maxnr > 0 ? gmx::ssize(sc->sel) + maxnr : -1);
367
368     state->statusWriter = statusWriter;
369
370     state->currentLocation.startIndex = 0;
371     state->currentLocation.endIndex   = 0;
372
373     snew(state->mstack, 20);
374     state->mstack_alloc     = 20;
375     state->msp              = -1;
376     state->neom             = 0;
377     state->nextparam        = nullptr;
378     state->nextMethodSymbol = nullptr;
379     state->prev_pos_kw      = 0;
380     state->bBoolNo          = false;
381     state->bMatchOf         = false;
382     state->bMatchBool       = false;
383     state->bCmdStart        = true;
384     state->bBuffer          = false;
385
386     _gmx_sel_yyset_extra(state, *scannerp);
387 }
388
389 void _gmx_sel_free_lexer(yyscan_t scanner)
390 {
391     gmx_sel_lexer_t* state = _gmx_sel_yyget_extra(scanner);
392
393     sfree(state->mstack);
394     if (state->bBuffer)
395     {
396         _gmx_sel_yy_delete_buffer(state->buffer, scanner);
397     }
398     delete state;
399     _gmx_sel_yylex_destroy(scanner);
400 }
401
402 void _gmx_sel_lexer_set_exception(yyscan_t scanner, const std::exception_ptr& ex)
403 {
404     gmx_sel_lexer_t* state = _gmx_sel_yyget_extra(scanner);
405     state->exception       = ex;
406 }
407
408 void _gmx_sel_lexer_rethrow_exception_if_occurred(yyscan_t scanner)
409 {
410     gmx_sel_lexer_t* state = _gmx_sel_yyget_extra(scanner);
411     if (state->exception)
412     {
413         std::exception_ptr ex = state->exception;
414         state->exception      = std::exception_ptr();
415         std::rethrow_exception(ex);
416     }
417 }
418
419 gmx::TextWriter* _gmx_sel_lexer_get_status_writer(yyscan_t scanner)
420 {
421     gmx_sel_lexer_t* state = _gmx_sel_yyget_extra(scanner);
422     return state->statusWriter;
423 }
424
425 struct gmx_ana_selcollection_t* _gmx_sel_lexer_selcollection(yyscan_t scanner)
426 {
427     gmx_sel_lexer_t* state = _gmx_sel_yyget_extra(scanner);
428     return state->sc;
429 }
430
431 bool _gmx_sel_lexer_has_groups_set(yyscan_t scanner)
432 {
433     gmx_sel_lexer_t* state = _gmx_sel_yyget_extra(scanner);
434     return state->bGroups;
435 }
436
437 struct gmx_ana_indexgrps_t* _gmx_sel_lexer_indexgrps(yyscan_t scanner)
438 {
439     gmx_sel_lexer_t* state = _gmx_sel_yyget_extra(scanner);
440     return state->grps;
441 }
442
443 int _gmx_sel_lexer_exp_selcount(yyscan_t scanner)
444 {
445     gmx_sel_lexer_t* state = _gmx_sel_yyget_extra(scanner);
446     return state->nexpsel;
447 }
448
449 const char* _gmx_sel_lexer_pselstr(yyscan_t scanner)
450 {
451     gmx_sel_lexer_t* state = _gmx_sel_yyget_extra(scanner);
452     return state->pselstr.c_str();
453 }
454
455 void _gmx_sel_lexer_set_current_location(yyscan_t scanner, const gmx::SelectionLocation& location)
456 {
457     gmx_sel_lexer_t* state = _gmx_sel_yyget_extra(scanner);
458     state->currentLocation = location;
459 }
460
461 const gmx::SelectionLocation& _gmx_sel_lexer_get_current_location(yyscan_t scanner)
462 {
463     gmx_sel_lexer_t* state = _gmx_sel_yyget_extra(scanner);
464     return state->currentLocation;
465 }
466
467 std::string _gmx_sel_lexer_get_current_text(yyscan_t scanner)
468 {
469     gmx_sel_lexer_t* state = _gmx_sel_yyget_extra(scanner);
470     return _gmx_sel_lexer_get_text(scanner, state->currentLocation);
471 }
472
473 std::string _gmx_sel_lexer_get_text(yyscan_t scanner, const gmx::SelectionLocation& location)
474 {
475     gmx_sel_lexer_t* state      = _gmx_sel_yyget_extra(scanner);
476     const int        startIndex = location.startIndex;
477     const int        endIndex   = location.endIndex;
478     if (startIndex >= endIndex)
479     {
480         return std::string();
481     }
482     return state->pselstr.substr(startIndex, endIndex - startIndex);
483 }
484
485 void _gmx_sel_lexer_clear_pselstr(yyscan_t scanner)
486 {
487     gmx_sel_lexer_t* state = _gmx_sel_yyget_extra(scanner);
488     state->pselstr.clear();
489 }
490
491 void _gmx_sel_lexer_clear_method_stack(yyscan_t scanner)
492 {
493     gmx_sel_lexer_t* state = _gmx_sel_yyget_extra(scanner);
494
495     state->msp = -1;
496 }
497
498 void _gmx_sel_finish_method(yyscan_t scanner)
499 {
500     gmx_sel_lexer_t* state = _gmx_sel_yyget_extra(scanner);
501
502     if (state->msp >= 0)
503     {
504         --state->msp;
505     }
506 }
507
508 void _gmx_sel_set_lex_input_file(yyscan_t scanner, FILE* fp)
509 {
510     gmx_sel_lexer_t* state = _gmx_sel_yyget_extra(scanner);
511
512     state->bBuffer = true;
513     state->buffer  = _gmx_sel_yy_create_buffer(fp, YY_BUF_SIZE, scanner);
514     _gmx_sel_yy_switch_to_buffer(state->buffer, scanner);
515 }
516
517 void _gmx_sel_set_lex_input_str(yyscan_t scanner, const char* str)
518 {
519     gmx_sel_lexer_t* state = _gmx_sel_yyget_extra(scanner);
520
521     if (state->bBuffer)
522     {
523         _gmx_sel_yy_delete_buffer(state->buffer, scanner);
524     }
525     state->bBuffer = true;
526     state->buffer  = _gmx_sel_yy_scan_string(str, scanner);
527 }