Removed simple.h from nb_kernel_sse4_1_XX
[gromacs.git] / src / gromacs / selection / scanner_internal.cpp
blobd66d647b58dfb524cc635ae20429cea8ae68ccdc
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012,2014,2015, 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 \file
36 * \brief Helper functions for the selection tokenizer.
38 * This file implements the functions in the headers scanner.h and
39 * scanner_internal.h.
41 * \author Teemu Murtola <teemu.murtola@gmail.com>
42 * \ingroup module_selection
44 /*! \cond
45 * \internal \file scanner_flex.h
46 * \brief Generated (from scanner.l) header file by Flex.
48 * This file contains definitions of functions that are needed in
49 * scanner_internal.cpp.
51 * \ingroup module_selection
52 * \endcond
54 #include "gmxpre.h"
56 #include "scanner_internal.h"
58 #include <stdlib.h>
59 #include <string.h>
61 #include <string>
63 #include "gromacs/utility/cstringutil.h"
64 #include "gromacs/utility/exceptions.h"
65 #include "gromacs/utility/gmxassert.h"
66 #include "gromacs/utility/smalloc.h"
67 #include "gromacs/utility/stringutil.h"
69 #include "parser.h"
70 #include "parsetree.h"
71 #include "scanner.h"
72 #include "selectioncollection-impl.h"
73 #include "selelem.h"
74 #include "selmethod.h"
75 #include "symrec.h"
77 /*! \brief
78 * Step in which the allocated memory for pretty-printed input is incremented.
80 #define STRSTORE_ALLOCSTEP 1000
82 /* These are defined as macros in the generated scanner_flex.h.
83 * We undefine them here to have them as variable names in the subroutines.
84 * There are other ways of doing this, but this is probably the easiest. */
85 #undef yylval
86 #undef yytext
87 #undef yyleng
89 /*! \brief
90 * Handles initialization of method parameter token.
92 static int
93 init_param_token(YYSTYPE *yylval, gmx_ana_selparam_t *param, bool bBoolNo)
95 if (bBoolNo)
97 GMX_RELEASE_ASSERT(param->name != NULL,
98 "bBoolNo should only be set for a parameters with a name");
99 snew(yylval->str, strlen(param->name) + 3);
100 yylval->str[0] = 'n';
101 yylval->str[1] = 'o';
102 strcpy(yylval->str+2, param->name);
104 else
106 yylval->str = param->name ? gmx_strdup(param->name) : NULL;
108 return PARAM;
111 /*! \brief
112 * Processes a selection method token.
114 static int
115 init_method_token(YYSTYPE *yylval, YYLTYPE *yylloc,
116 const gmx::SelectionParserSymbol *symbol,
117 bool bPosMod, gmx_sel_lexer_t *state)
119 gmx_ana_selmethod_t *method = symbol->methodValue();
120 /* If the previous token was not KEYWORD_POS, return EMPTY_POSMOD
121 * before the actual method to work around a limitation in Bison. */
122 if (!bPosMod && method->type != POS_VALUE)
124 state->nextMethodSymbol = symbol;
125 _gmx_sel_lexer_add_token(yylloc, NULL, 0, state);
126 return EMPTY_POSMOD;
128 _gmx_sel_lexer_add_token(yylloc, symbol->name().c_str(), -1, state);
129 yylval->meth = method;
130 if (!(method->flags & SMETH_MODIFIER) && method->nparams == 0)
132 /* Keyword */
133 switch (method->type)
135 case INT_VALUE:
136 case REAL_VALUE:
137 state->bMatchOf = true;
138 return KEYWORD_NUMERIC;
139 case STR_VALUE: return KEYWORD_STR;
140 case GROUP_VALUE: return KEYWORD_GROUP;
141 default:
142 GMX_THROW(gmx::InternalError("Unsupported keyword type"));
145 else
147 /* Method with parameters or a modifier */
148 if (method->flags & SMETH_MODIFIER)
150 /* Remove all methods from the stack */
151 state->msp = -1;
152 if (method->param[1].name == NULL)
154 state->nextparam = &method->param[1];
157 else
159 if (method->param[0].name == NULL)
161 state->nextparam = &method->param[0];
164 ++state->msp;
165 if (state->msp >= state->mstack_alloc)
167 state->mstack_alloc += 10;
168 srenew(state->mstack, state->mstack_alloc);
170 state->mstack[state->msp] = method;
171 if (method->flags & SMETH_MODIFIER)
173 return MODIFIER;
175 switch (method->type)
177 case INT_VALUE: return METHOD_NUMERIC;
178 case REAL_VALUE: return METHOD_NUMERIC;
179 case POS_VALUE: return METHOD_POS;
180 case GROUP_VALUE: return METHOD_GROUP;
181 default:
182 --state->msp;
183 GMX_THROW(gmx::InternalError("Unsupported method type"));
186 return INVALID; /* Should not be reached */
190 _gmx_sel_lexer_process_pending(YYSTYPE *yylval, YYLTYPE *yylloc,
191 gmx_sel_lexer_t *state)
193 if (state->nextparam)
195 gmx_ana_selparam_t *param = state->nextparam;
196 bool bBoolNo = state->bBoolNo;
198 if (state->neom > 0)
200 --state->neom;
201 _gmx_sel_lexer_add_token(yylloc, NULL, 0, state);
202 return END_OF_METHOD;
204 state->nextparam = NULL;
205 state->bBoolNo = false;
206 _gmx_sel_lexer_add_token(yylloc, param->name, -1, state);
207 return init_param_token(yylval, param, bBoolNo);
209 if (state->prev_pos_kw > 0)
211 --state->prev_pos_kw;
213 if (state->nextMethodSymbol)
215 const gmx::SelectionParserSymbol *symbol = state->nextMethodSymbol;
216 state->nextMethodSymbol = NULL;
217 return init_method_token(yylval, yylloc, symbol, true, state);
219 return 0;
223 _gmx_sel_lexer_process_identifier(YYSTYPE *yylval, YYLTYPE *yylloc,
224 char *yytext, size_t yyleng,
225 gmx_sel_lexer_t *state)
227 /* Check if the identifier matches with a parameter name */
228 if (state->msp >= 0)
230 gmx_ana_selparam_t *param = NULL;
231 bool bBoolNo = false;
232 int sp = state->msp;
233 while (!param && sp >= 0)
235 int i;
236 for (i = 0; i < state->mstack[sp]->nparams; ++i)
238 /* Skip NULL parameters and too long parameters */
239 if (state->mstack[sp]->param[i].name == NULL
240 || strlen(state->mstack[sp]->param[i].name) > yyleng)
242 continue;
244 if (!strncmp(state->mstack[sp]->param[i].name, yytext, yyleng))
246 param = &state->mstack[sp]->param[i];
247 break;
249 /* Check separately for a 'no' prefix on boolean parameters */
250 if (state->mstack[sp]->param[i].val.type == NO_VALUE
251 && yyleng > 2 && yytext[0] == 'n' && yytext[1] == 'o'
252 && !strncmp(state->mstack[sp]->param[i].name, yytext+2, yyleng-2))
254 param = &state->mstack[sp]->param[i];
255 bBoolNo = true;
256 break;
259 if (!param)
261 --sp;
264 if (param)
266 if (param->val.type == NO_VALUE && !bBoolNo)
268 state->bMatchBool = true;
270 if (sp < state->msp)
272 state->neom = state->msp - sp - 1;
273 state->nextparam = param;
274 state->bBoolNo = bBoolNo;
275 return END_OF_METHOD;
277 _gmx_sel_lexer_add_token(yylloc, param->name, -1, state);
278 return init_param_token(yylval, param, bBoolNo);
282 /* Check if the identifier matches with a symbol */
283 const gmx::SelectionParserSymbol *symbol
284 = state->sc->symtab->findSymbol(std::string(yytext, yyleng));
285 /* If there is no match, return the token as a string */
286 if (!symbol)
288 yylval->str = gmx_strndup(yytext, yyleng);
289 _gmx_sel_lexer_add_token(yylloc, yytext, yyleng, state);
290 return IDENTIFIER;
292 gmx::SelectionParserSymbol::SymbolType symtype = symbol->type();
293 /* For method symbols, we need some extra processing. */
294 if (symtype == gmx::SelectionParserSymbol::MethodSymbol)
296 return init_method_token(yylval, yylloc, symbol, state->prev_pos_kw > 0, state);
298 _gmx_sel_lexer_add_token(yylloc, symbol->name().c_str(), -1, state);
299 /* Reserved symbols should have been caught earlier */
300 if (symtype == gmx::SelectionParserSymbol::ReservedSymbol)
302 GMX_THROW(gmx::InternalError(gmx::formatString(
303 "Mismatch between tokenizer and reserved symbol table (for '%s')",
304 symbol->name().c_str())));
306 /* For variable symbols, return the type of the variable value */
307 if (symtype == gmx::SelectionParserSymbol::VariableSymbol)
309 gmx::SelectionTreeElementPointer var = symbol->variableValue();
310 /* Return simple tokens for constant variables */
311 if (var->type == SEL_CONST)
313 switch (var->v.type)
315 case INT_VALUE:
316 yylval->i = var->v.u.i[0];
317 return TOK_INT;
318 case REAL_VALUE:
319 yylval->r = var->v.u.r[0];
320 return TOK_REAL;
321 case POS_VALUE:
322 break;
323 default:
324 GMX_THROW(gmx::InternalError("Unsupported variable type"));
327 yylval->sel = new gmx::SelectionTreeElementPointer(var);
328 switch (var->v.type)
330 case INT_VALUE: return VARIABLE_NUMERIC;
331 case REAL_VALUE: return VARIABLE_NUMERIC;
332 case POS_VALUE: return VARIABLE_POS;
333 case GROUP_VALUE: return VARIABLE_GROUP;
334 default:
335 delete yylval->sel;
336 GMX_THROW(gmx::InternalError("Unsupported variable type"));
337 return INVALID;
339 /* This position should not be reached. */
341 /* For position symbols, we need to return KEYWORD_POS, but we also need
342 * some additional handling. */
343 if (symtype == gmx::SelectionParserSymbol::PositionSymbol)
345 state->bMatchOf = true;
346 yylval->str = gmx_strdup(symbol->name().c_str());
347 state->prev_pos_kw = 2;
348 return KEYWORD_POS;
350 /* Should not be reached */
351 return INVALID;
354 void
355 _gmx_sel_lexer_add_token(YYLTYPE *yylloc, const char *str, int len,
356 gmx_sel_lexer_t *state)
358 yylloc->startIndex = yylloc->endIndex = state->pslen;
359 /* Do nothing if the string is empty, or if it is a space and there is
360 * no other text yet, or if there already is a space. */
361 if (!str || len == 0 || strlen(str) == 0
362 || (str[0] == ' ' && str[1] == 0
363 && (state->pslen == 0 || state->pselstr[state->pslen - 1] == ' ')))
365 return;
367 if (len < 0)
369 len = strlen(str);
371 /* Allocate more memory if necessary */
372 if (state->nalloc_psel - state->pslen < len)
374 int incr = STRSTORE_ALLOCSTEP < len ? len : STRSTORE_ALLOCSTEP;
375 state->nalloc_psel += incr;
376 srenew(state->pselstr, state->nalloc_psel);
378 /* Append the token to the stored string */
379 strncpy(state->pselstr + state->pslen, str, len);
380 state->pslen += len;
381 state->pselstr[state->pslen] = 0;
382 yylloc->endIndex = state->pslen;
385 void
386 _gmx_sel_init_lexer(yyscan_t *scannerp, struct gmx_ana_selcollection_t *sc,
387 gmx::TextWriter *statusWriter, int maxnr,
388 bool bGroups, struct gmx_ana_indexgrps_t *grps)
390 int rc = _gmx_sel_yylex_init(scannerp);
391 if (rc != 0)
393 // TODO: Throw a more representative exception.
394 GMX_THROW(gmx::InternalError("Lexer initialization failed"));
397 gmx_sel_lexer_t *state = new gmx_sel_lexer_t;
399 state->sc = sc;
400 state->bGroups = bGroups;
401 state->grps = grps;
402 state->nexpsel = (maxnr > 0 ? static_cast<int>(sc->sel.size()) + maxnr : -1);
404 state->statusWriter = statusWriter;
406 snew(state->pselstr, STRSTORE_ALLOCSTEP);
407 state->pselstr[0] = 0;
408 state->pslen = 0;
409 state->nalloc_psel = STRSTORE_ALLOCSTEP;
410 state->currentLocation.startIndex = 0;
411 state->currentLocation.endIndex = 0;
413 snew(state->mstack, 20);
414 state->mstack_alloc = 20;
415 state->msp = -1;
416 state->neom = 0;
417 state->nextparam = NULL;
418 state->nextMethodSymbol = NULL;
419 state->prev_pos_kw = 0;
420 state->bBoolNo = false;
421 state->bMatchOf = false;
422 state->bMatchBool = false;
423 state->bCmdStart = true;
424 state->bBuffer = false;
426 _gmx_sel_yyset_extra(state, *scannerp);
429 void
430 _gmx_sel_free_lexer(yyscan_t scanner)
432 gmx_sel_lexer_t *state = _gmx_sel_yyget_extra(scanner);
434 sfree(state->pselstr);
435 sfree(state->mstack);
436 if (state->bBuffer)
438 _gmx_sel_yy_delete_buffer(state->buffer, scanner);
440 delete state;
441 _gmx_sel_yylex_destroy(scanner);
444 void
445 _gmx_sel_lexer_set_exception(yyscan_t scanner,
446 const boost::exception_ptr &ex)
448 gmx_sel_lexer_t *state = _gmx_sel_yyget_extra(scanner);
449 state->exception = ex;
452 void
453 _gmx_sel_lexer_rethrow_exception_if_occurred(yyscan_t scanner)
455 gmx_sel_lexer_t *state = _gmx_sel_yyget_extra(scanner);
456 if (state->exception)
458 boost::exception_ptr ex = state->exception;
459 state->exception = boost::exception_ptr();
460 rethrow_exception(ex);
464 gmx::TextWriter *
465 _gmx_sel_lexer_get_status_writer(yyscan_t scanner)
467 gmx_sel_lexer_t *state = _gmx_sel_yyget_extra(scanner);
468 return state->statusWriter;
471 struct gmx_ana_selcollection_t *
472 _gmx_sel_lexer_selcollection(yyscan_t scanner)
474 gmx_sel_lexer_t *state = _gmx_sel_yyget_extra(scanner);
475 return state->sc;
478 bool
479 _gmx_sel_lexer_has_groups_set(yyscan_t scanner)
481 gmx_sel_lexer_t *state = _gmx_sel_yyget_extra(scanner);
482 return state->bGroups;
485 struct gmx_ana_indexgrps_t *
486 _gmx_sel_lexer_indexgrps(yyscan_t scanner)
488 gmx_sel_lexer_t *state = _gmx_sel_yyget_extra(scanner);
489 return state->grps;
493 _gmx_sel_lexer_exp_selcount(yyscan_t scanner)
495 gmx_sel_lexer_t *state = _gmx_sel_yyget_extra(scanner);
496 return state->nexpsel;
499 const char *
500 _gmx_sel_lexer_pselstr(yyscan_t scanner)
502 gmx_sel_lexer_t *state = _gmx_sel_yyget_extra(scanner);
503 return state->pselstr;
506 void
507 _gmx_sel_lexer_set_current_location(yyscan_t scanner,
508 const gmx::SelectionLocation &location)
510 gmx_sel_lexer_t *state = _gmx_sel_yyget_extra(scanner);
511 state->currentLocation = location;
514 const gmx::SelectionLocation &
515 _gmx_sel_lexer_get_current_location(yyscan_t scanner)
517 gmx_sel_lexer_t *state = _gmx_sel_yyget_extra(scanner);
518 return state->currentLocation;
521 std::string
522 _gmx_sel_lexer_get_current_text(yyscan_t scanner)
524 gmx_sel_lexer_t *state = _gmx_sel_yyget_extra(scanner);
525 return _gmx_sel_lexer_get_text(scanner, state->currentLocation);
528 std::string
529 _gmx_sel_lexer_get_text(yyscan_t scanner,
530 const gmx::SelectionLocation &location)
532 gmx_sel_lexer_t *state = _gmx_sel_yyget_extra(scanner);
533 const int startIndex = location.startIndex;
534 const int endIndex = location.endIndex;
535 if (startIndex >= endIndex)
537 return std::string();
539 return std::string(&state->pselstr[startIndex], endIndex - startIndex);
542 void
543 _gmx_sel_lexer_clear_pselstr(yyscan_t scanner)
545 gmx_sel_lexer_t *state = _gmx_sel_yyget_extra(scanner);
546 state->pselstr[0] = 0;
547 state->pslen = 0;
550 void
551 _gmx_sel_lexer_clear_method_stack(yyscan_t scanner)
553 gmx_sel_lexer_t *state = _gmx_sel_yyget_extra(scanner);
555 state->msp = -1;
558 void
559 _gmx_sel_finish_method(yyscan_t scanner)
561 gmx_sel_lexer_t *state = _gmx_sel_yyget_extra(scanner);
563 if (state->msp >= 0)
565 --state->msp;
569 void
570 _gmx_sel_set_lex_input_file(yyscan_t scanner, FILE *fp)
572 gmx_sel_lexer_t *state = _gmx_sel_yyget_extra(scanner);
574 state->bBuffer = true;
575 state->buffer = _gmx_sel_yy_create_buffer(fp, YY_BUF_SIZE, scanner);
576 _gmx_sel_yy_switch_to_buffer(state->buffer, scanner);
579 void
580 _gmx_sel_set_lex_input_str(yyscan_t scanner, const char *str)
582 gmx_sel_lexer_t *state = _gmx_sel_yyget_extra(scanner);
584 if (state->bBuffer)
586 _gmx_sel_yy_delete_buffer(state->buffer, scanner);
588 state->bBuffer = true;
589 state->buffer = _gmx_sel_yy_scan_string(str, scanner);