FEENOX

2024-03-19

NAME

FeenoX - a cloud-first free no-X uniX-like finite-element(ish) computational engineering tool

SYNOPSIS

The basic usage is to execute the feenox binary passing a path to an input file that defines the problem, along with other options and command-line replacement arguments which are explained below:

feenox [options ...] input-file [optional_commandline_replacement_arguments ...]

For large problems that do not fit in a single computer, a parallel run using mpirun(1) will be needed:

mpirun -n number_of_threads feenox [options ...] input-file [optional_commandline_replacement_arguments ...]

DESCRIPTION

FeenoX is a computational tool that can solve engineering problems which are usually casted as differential-algebraic equations (DAEs) or partial differential equations (PDEs). It is to finite elements programs and libraries what Markdown is to Word and TeX, respectively. In particular, it can solve

FeenoX reads a plain-text input file which contains the problem definition and writes 100%-user defined results in ASCII (through PRINT or other user-defined output instructions within the input file). For PDE problems, it needs a reference to at least one Gmsh mesh file for the discretization of the domain. It can write post-processing views in either .msh or .vtk formats.

Keep in mind that FeenoX is just a back end reading a set of input files and writing a set of output files following the design philosophy of Unix (separation, composition, representation, economy, extensibility, etc). Think of it as a transfer function (or a filter in computer-science jargon) between input files and output files:

                             +------------+
 mesh (*.msh)  }             |            |             { terminal
 data (*.dat)  } input ----> |   FeenoX   |----> output { data files
 input (*.fee) }             |            |             { post (vtk/msh)
                             +------------+

Following the Unix programming philosophy, there are no graphical interfaces attached to the FeenoX core, although a wide variety of pre and post-processors can be used with FeenoX. To illustrate the transfer-function approach, consider the following input file that solves Laplace’s equation ∇^2^ϕ = 0 on a square with some space-dependent boundary conditions:

ϕ(x,y) =  + y   for x =  − 1 (left)

ϕ(x,y) =  − y   for x =  + 1 (right)

ϕ ⋅  = sin (π/2x)   for y =  − 1 (bottom)

ϕ ⋅  = 0   for y =  + 1 (top)

PROBLEM laplace 2d
READ_MESH square-centered.msh # [-1:+1]x[-1:+1]

# boundary conditions
BC left    phi=+y
BC right   phi=-y
BC bottom  dphidn=sin(pi/2*x)
BC top     dphidn=0

SOLVE_PROBLEM

# same output in .msh and in .vtk formats
WRITE_MESH laplace-square.msh phi VECTOR dphidx dphidy 0
WRITE_MESH laplace-square.vtk phi VECTOR dphidx dphidy 0

Laplace’s equation solved with FeenoX

The .msh file can be post-processed with Gmsh, and the .vtk file can be post-processed with Paraview. See https://www.caeplex.com for a mobile-friendly web-based interface for solving finite elements in the cloud directly from the browser.

OPTIONS

feenox [options] inputfile [replacement arguments] [petsc options]  
-h, --help

display options and detailed explanations of commmand-line usage

-v, --version

display brief version information and exit

-V, --versions

display detailed version information

-c, --check

validates if the input file is sane or not

--pdes

list the types of PROBLEMs that FeenoX can solve, one per line

--elements_info

output a document with information about the supported element types

--linear

force FeenoX to solve the PDE problem as linear

--non-linear

force FeenoX to solve the PDE problem as non-linear

--progress

print ASCII progress bars when solving PDEs

--mumps

ask PETSc to use the direct linear solver MUMPS

Instructions will be read from standard input if “-” is passed as inputfile, i.e.

$ echo 'PRINT 2+2' | feenox -
4

The optional [replacement arguments] part of the command line mean that each argument after the input file that does not start with an hyphen will be expanded verbatim in the input file in each occurrence of $1, $2, etc. For example

$ echo 'PRINT $1+$2' | feenox - 3 4
7

PETSc and SLEPc options can be passed in [petsc options] (or [options]) as well, with the difference that two hyphens have to be used instead of only once. For example, to pass the PETSc option -ksp_view the actual FeenoX invocation should be

$ feenox input.fee --ksp_view

For PETSc options that take values, en equal sign has to be used:

$ feenox input.fee --mg_levels_pc_type=sor

See https://www.seamplex.com/feenox/examples for annotated examples.

EXAMPLES

EXIT STATUS

ENVIRONMENT

FILES

CONFORMING TO

INPUT-FILE KEYWORDS

GENERIC KEYWORDS

ABORT

Catastrophically abort the execution and quit FeenoX.

ABORT  

Whenever the instruction ABORT is executed, FeenoX quits with a non-zero error leve. It does not close files nor unlock shared memory objects. The objective of this instruction is to either debug complex input files by using only parts of them or to conditionally abort the execution using IF clauses.

ALIAS

Define a scalar alias of an already-defined indentifier.

ALIAS { <new_var_name> IS <existing_object> | <existing_object> AS <new_name> }  

The existing object can be a variable, a vector element or a matrix element. In the first case, the name of the variable should be given as the existing object. In the second case, to alias the second element of vector v to the new name new, v(2) should be given as the existing object. In the third case, to alias second element (2,3) of matrix M to the new name new, M(2,3) should be given as the existing object.

CLOSE

Explicitly close a file after input/output.

CLOSE <name>  

The given <name> can be either a fixed-string path or an already-defined FILE.

DEFAULT_ARGUMENT_VALUE

Give a default value for an optional commandline argument.

DEFAULT_ARGUMENT_VALUE <constant> <string>  

If a $n construction is found in the input file but the commandline argument was not given, the default behavior is to fail complaining that an extra argument has to be given in the commandline. With this keyword, a default value can be assigned if no argument is given, thus avoiding the failure and making the argument optional. The <constant> should be 1, 2, 3, etc. and <string> will be expanded character-by-character where the $n construction is. Whether the resulting expression is to be interpreted as a string or as a numerical expression will depend on the context.

FILE

Define a file with a particularly formatted name to be used either as input or as output.

< FILE | OUTPUT_FILE | INPUT_FILE > <name> PATH <format> expr_1 expr_2 ... expr_n [ INPUT | OUTPUT | APPEND | MODE <fopen_mode> ]  

For reading or writing into files with a fixed path, this instruction is usually not needed as the FILE keyword of other instructions (such as PRINT or MESH) can take a fixed-string path as an argument. However, if the file name changes as the execution progresses (say because one file for each step is needed), then an explicit FILE needs to be defined with this keyword and later referenced by the given name.

The path should be given as a printf-like format string followed by the expressions which shuold be evaluated in order to obtain the actual file path. The expressions will always be floating-point expressions, but the particular integer specifier %d is allowed and internally transformed to %.0f. The file can be explicitly defined and INPUT, OUTPUT or a certain fopen() mode can be given (i.e. “a”). If not explicitly given, the nature of the file will be taken from context, i.e. FILEs in PRINT will be OUTPUT and FILEs in FUNCTION will be INPUT. This keyword justs defines the FILE, it does not open it. The file will be actually openened (and eventually closed) automatically. In the rare case where the automated opening and closing does not fit the expected workflow, the file can be explicitly opened or closed with the instructions FILE_OPEN and FILE_CLOSE.

FIT

Find parameters to fit an analytical function to a pointwise-defined function.

FIT <function_to_be_fitted>  TO <function_with_data> VIA <var_1> <var_2> ... <var_n>
 [ GRADIENT <expr_1> <expr_2> ... <expr_n> ]
 [ RANGE_MIN <expr_1> <expr_2> ... <expr_j> ]
 [ RANGE_MAX <expr_1> <expr_2> ... <expr_n> ]
 [ TOL_REL <expr> ] [ TOL_ABS <expr> ] [ MAX_ITER <expr> ]
 [ VERBOSE ]  

The function with the data has to be point-wise defined (i.e. a FUNCTION read from a file, with inline DATA or defined over a mesh). The function to be fitted has to be parametrized with at least one of the variables provided after the USING keyword. For example to fit f(x,y) = ax^2^ + bsqrt(y) to a pointwise-defined function g(x,y) one gives FIT f TO g VIA a b. Only the names of the functions have to be given, not the arguments. Both functions have to have the same number of arguments. The initial guess of the solution is given by the initial value of the variables after the VIA keyword. Analytical expressions for the gradient of the function to be fitted with respect to the parameters to be fitted can be optionally given with the GRADIENT keyword. If none is provided, the gradient will be computed numerically using finite differences. A range over which the residuals are to be minimized can be given with RANGE_MIN and RANGE_MAX. The expressions give the range of the arguments of the functions, not of the parameters. For multidimensional fits, the range is an hypercube. If no range is given, all the definition points of the function with the data are used for the fit. Convergence can be controlled by giving the relative and absolute tolreances with TOL_REL (default DEFAULT_NLIN_FIT_EPSREL) and TOL_ABS (default DEFAULT_NLIN_FIT_EPSABS), and with the maximum number of iterations MAX_ITER (default DEFAULT_NLIN_FIT_MAX_ITER). If the optional keyword VERBOSE is given, some data of the intermediate steps is written in the standard output.

FUNCTION

Define a scalar function of one or more variables.

FUNCTION <function_name>(<var_1>[,var2,...,var_n]) { 
   = <expr> | 
   FILE { <file> } | 
   VECTORS <vector_1> <vector_2> ... <vector_n> <vector_data> | 
   MESH <mesh> | 
   DATA <num_1> <num_2> ... <num_N> 
 } 
 [ COLUMNS <expr_1> <expr_2> ... <expr_n> <expr_n+1> ] 
 [ INTERPOLATION { linear | polynomial | spline | spline_periodic | akima | akima_periodic | steffen | 
 nearest | shepard | shepard_kd | bilinear } ] 
 [ INTERPOLATION_THRESHOLD <expr> ] [ SHEPARD_RADIUS <expr> ] [ SHEPARD_EXPONENT <expr> ]  

The number of variables n is given by the number of arguments given between parenthesis after the function name. The arguments are defined as new variables if they had not been already defined explictly as scalar variables. If the function is given as an algebraic expression, the short-hand operator = (or := for compatiblity with Maxima) can be used. That is to say, FUNCTION f(x) = x^2 is equivalent to f(x) = x^2 (or f(x) := x^2). If a FILE is given, an ASCII file containing at least n + 1 columns is expected. By default, the first n columns are the values of the arguments and the last column is the value of the function at those points. The order of the columns can be changed with the keyword COLUMNS, which expects n + 1 expressions corresponding to the column numbers. If VECTORS is given, a set of n + 1 vectors of the same size is expected. The first n correspond to the arguments and the last one to the function values. If MESH is given, the function is point-wise defined over the mesh topology. That is to say, the independent variables (i.e. the spatial coordinates) coincide with the mesh nodes. The dependent variable (i.e. the function value) is set by “filling” a vector named vec_f (where f has to be replaced with the function name) of size equal to the number of nodes.

The function can be pointwise-defined inline in the input using DATA. This should be the last keyword of the line, followed by N = k ⋅ (n+1) expresions giving k definition points: n arguments and the value of the function. Multiline continuation using brackets { and } can be used for a clean data organization. Interpolation schemes can be given for either one or multi-dimensional functions with INTERPOLATION. Available schemes for n = 1 are:

Default interpolation scheme for one-dimensional functions is DEFAULT_INTERPOLATION.

Available schemes for n > 1 are:

For n > 1, if the euclidean distance between the arguments and the definition points is smaller than INTERPOLATION_THRESHOLD, the definition point is returned and no interpolation is performed. Default value is square root of DEFAULT_MULTIDIM_INTERPOLATION_THRESHOLD.

The initial radius of points to take into account in shepard_kd is given by SHEPARD_RADIUS. If no points are found, the radius is double until at least one definition point is found. The radius is doubled until at least one point is found. Default is DEFAULT_SHEPARD_RADIUS. The exponent of the shepard method is given by SHEPARD_EXPONENT. Default is DEFAULT_SHEPARD_EXPONENT.

IF

Execute a set of instructions if a condition is met.

IF expr 
  <block_of_instructions_if_expr_is_true> 
 [ ELSE  
  <block_of_instructions_if_expr_is_false> ] 
 ENDIF  

IMPLICIT

Define whether implicit definition of variables is allowed or not.

IMPLICIT { NONE | ALLOWED }  

By default, FeenoX allows variables (but not vectors nor matrices) to be implicitly declared. To avoid introducing errors due to typos, explicit declaration of variables can be forced by giving IMPLICIT NONE. Whether implicit declaration is allowed or explicit declaration is required depends on the last IMPLICIT keyword given, which by default is ALLOWED.

INCLUDE

Include another FeenoX input file.

INCLUDE <file_path> [ FROM <num_expr> ] [ TO <num_expr> ]  

Includes the input file located in the string file_path at the current location. The effect is the same as copying and pasting the contents of the included file at the location of the INCLUDE keyword. The path can be relative or absolute. Note, however, that when including files inside IF blocks that instructions are conditionally-executed but all definitions (such as function definitions) are processed at parse-time independently from the evaluation of the conditional. The included file has to be an actual file path (i.e. it cannot be a FeenoX FILE) because it needs to be resolved at parse time. Yet, the name can contain a commandline replacement argument such as $1 so INCLUDE $1.fee will include the file specified after the main input file in the command line. The optional FROM and TO keywords can be used to include only portions of a file.

MATRIX

Define a matrix.

MATRIX <name> ROWS <expr> COLS <expr> [ DATA <expr_1> <expr_2> ... <expr_n> |  

A new matrix of the prescribed size is defined. The number of rows and columns can be an expression which will be evaluated the very first time the matrix is used and then kept at those constant values. All elements will be initialized to zero unless DATA is given (which should be the last keyword of the line), in which case the expressions will be evaluated the very first time the matrix is used and row-major-assigned to each of the elements. If there are less elements than the matrix size, the remaining values will be zero. If there are more elements than the matrix size, the values will be ignored.

OPEN

Explicitly open a file for input/output.

OPEN <name> [ MODE <fopen_mode> ]  

The given <name> can be either a fixed-string path or an already-defined FILE. The mode is only taken into account if the file is not already defined. Default is write w.

PRINT

Write plain-text and/or formatted data to the standard output or into an output file.

PRINT [ <object_1> <object_2> ... <object_n> ] [ TEXT <string_1> ... TEXT <string_n> ] 
 [ FILE { <file_path> | <file_id> } ] [ HEADER ] [ NONEWLINE ] [ SEP <string> ] 
 [ SKIP_STEP <expr> ] [ SKIP_STATIC_STEP <expr> ] [ SKIP_TIME <expr> ] [ SKIP_HEADER_STEP <expr> ] 
  

Each argument object which is not a keyword of the PRINT instruction will be part of the output. Objects can be either a matrix, a vector or any valid scalar algebraic expression. If the given object cannot be solved into a valid matrix, vector or expression, it is treated as a string literal if IMPLICIT is ALLOWED, otherwise a parser error is raised. To explicitly interpret an object as a literal string even if it resolves to a valid numerical expression, it should be prefixed with the TEXT keyword such as PRINT TEXT 1+1 that would print 1+1 instead of 2. Objects and string literals can be mixed and given in any order. Hashes # appearing literal in text strings have to be quoted to prevent the parser to treat them as comments within the FeenoX input file and thus ignoring the rest of the line, like PRINT "\# this is a printed comment". Whenever an argument starts with a porcentage sign %, it is treated as a C printf-compatible format specifier and all the objects that follow it are printed using the given format until a new format definition is found. The objects are treated as double-precision floating point numbers, so only floating point formats should be given. See the printf(3) man page for further details. The default format is DEFAULT_PRINT_FORMAT. Matrices, vectors, scalar expressions, format modifiers and string literals can be given in any desired order, and are processed from left to right. Vectors are printed element-by-element in a single row. See PRINT_VECTOR to print one or more vectors with one element per line (i.e. vertically). Matrices are printed element-by-element in a single line using row-major ordering if mixed with other objects but in the natural row and column fashion if it is the only given object in the PRINT instruction. If the FILE keyword is not provided, default is to write to stdout. If the HEADER keyword is given, a single line containing the literal text given for each object is printed at the very first time the PRINT instruction is processed, starting with a hash # character.

If the NONEWLINE keyword is not provided, default is to write a newline \n character after all the objects are processed. Otherwise, if the last token to be printed is a numerical value, a separator string will be printed but not the newline \n character. If the last token is a string, neither the separator nor the newline will be printed. The SEP keyword expects a string used to separate printed objects. To print objects without any separation in between give an empty string like SEP "". The default is a tabulator character `DEFAULT_PRINT_SEPARATOR' character. To print an empty line write PRINT without arguments. By default the PRINT instruction is evaluated every step. If the SKIP_STEP (SKIP_STATIC_STEP) keyword is given, the instruction is processed only every the number of transient (static) steps that results in evaluating the expression, which may not be constant. The SKIP_HEADER_STEP keyword works similarly for the optional HEADER but by default it is only printed once. The SKIP_TIME keyword use time advancements to choose how to skip printing and may be useful for non-constant time-step problems.

PRINTF

Instruction akin to C’s printf. Instruction akin to C’s printf executed locally from all MPI ranks.

PRINTF PRINTF_ALL format_string [ expr_1 [ expr_2 [ ... ] ] ]  

The format_string should be a printf-like string containing double-precision format specifiers. A matching number of expressions should be given. No newline is written if not explicitly asked for in the format string with \n.

Do not ask for string literals %s.

As always, to get a literal % use %% in the format string.

PRINT_FUNCTION

Print one or more functions as a table of values of dependent and independent variables.

PRINT_FUNCTION <function_1> [ { function | expr } ... { function | expr } ] 
 [ FILE { <file_path> | <file_id> } ] [ HEADER ] 
 [ MIN <expr_1> <expr_2> ... <expr_k> ] [ MAX <expr_1> <expr_2> ... <expr_k> ] 
 [ STEP <expr_1> <expr_2> ... <expr_k> ] [ NSTEPs <expr_1> <expr_2> ... <expr_k> ] 
 [ FORMAT <print_format> ] <vector_1> [ { vector | expr } ... { vector | expr } ]  

Each argument should be either a function or an expression. The output of this instruction consists of n + k columns, where n is the number of arguments of the first function of the list and k is the number of functions and expressions given. The first n columns are the arguments (independent variables) and the last k one has the evaluated functions and expressions. The columns are separated by a tabulator, which is the format that most plotting tools understand. Only function names without arguments are expected. All functions should have the same number of arguments. Expressions can involve the arguments of the first function. If the FILE keyword is not provided, default is to write to stdout. If HEADER is given, the output is prepended with a single line containing the names of the arguments and the names of the functions, separated by tabs. The header starts with a hash # that usually acts as a comment and is ignored by most plotting tools. If there is no explicit range where to evaluate the functions and the first function is point-wise defined, they are evalauted at the points of definition of the first one. The range can be explicitly given as a product of n ranges [x~i, min ~,x~i, max ~] for i = 1, ..., n.

The values x~i, min ~ and x~i, max ~ are given with the MIN and MAX keywords. The discretization steps of the ranges are given by either STEP that gives δx or NSTEPS that gives the number of steps. If the first function is not point-wise defined, the ranges are mandatory.

PRINT_VECTOR

Print the elements of one or more vectors, one element per line.

PRINT_VECTOR 
 [ FILE { <file_path> | <file_id> } ] [ HEADER ] 
 [ SEP <string> ]  

Each argument should be either a vector or an expression of the integer i. If the FILE keyword is not provided, default is to write to stdout. If HEADER is given, the output is prepended with a single line containing the names of the arguments and the names of the functions, separated by tabs. The header starts with a hash # that usually acts as a comment and is ignored by most plotting tools. The SEP keyword expects a string used to separate printed objects. To print objects without any separation in between give an empty string like SEP "". The default is a tabulator character `DEFAULT_PRINT_SEPARATOR' character.

SOLVE

Solve a (small) system of non-linear equations.

SOLVE FOR <n> UNKNOWNS <var_1> <var_2> ... <var_n> [ METHOD { dnewton | hybrid | hybrids | broyden } ]
 [ EPSABS <expr> ] [ EPSREL <expr> ] [ MAX_ITER <expr> ]  

SORT_VECTOR

Sort the elements of a vector, optionally making the same rearrangement in another vector.

SORT_VECTOR <vector> [ ASCENDING | DESCENDING ] [ <other_vector> ]  

This instruction sorts the elements of <vector> into either ascending or descending numerical order. If <other_vector> is given, the same rearrangement is made on it. Default is ascending order.

VAR

Explicitly define one or more scalar variables.

VAR <name_1> [ <name_2> ] ... [ <name_n> ]  

When implicit definition is allowed (see IMPLICIT), scalar variables need not to be defined before being used if from the context FeenoX can tell that an scalar variable is needed. For instance, when defining a function like f(x) = x^2 it is not needed to declare x explictly as a scalar variable. But if one wants to define a function like g(x) = integral(f(x'), x', 0, x) then the variable x' needs to be explicitly defined as VAR x' before the integral.

VECTOR

Define a vector.

VECTOR <name> SIZE <expr> [ FUNCTION_DATA <function> ] [ DATA <expr_1> <expr_2> ... <expr_n> |  

A new vector of the prescribed size is defined. The size can be an expression which will be evaluated the very first time the vector is used and then kept at that constant value. If the keyword FUNCTION_DATA is given, the elements of the vector will be synchronized with the inpedendent values of the function, which should be point-wise defined. The sizes of both the function and the vector should match. All elements will be initialized to zero unless DATA is given (which should be the last keyword of the line), in which case the expressions will be evaluated the very first time the vector is used and assigned to each of the elements. If there are less elements than the vector size, the remaining values will be zero. If there are more elements than the vector size, the values will be ignored.

DAE-RELATED KEYWORDS

INITIAL_CONDITIONS

Define how initial conditions of DAE problems are computed.

INITIAL_CONDITIONS { AS_PROVIDED | FROM_VARIABLES | FROM_DERIVATIVES }  

In DAE problems, initial conditions may be either:

In the first case, it is up to the user to fulfill the DAE system at t = 0. If the residuals are not small enough, a convergence error will occur. The FROM_VARIABLES option means calling IDA’s IDACalcIC routine with the parameter IDA_YA_YDP_INIT. The FROM_DERIVATIVES option means calling IDA’s IDACalcIC routine with the parameter IDA_Y_INIT. Wasora should be able to automatically detect which variables in phase-space are differential and which are purely algebraic. However, the DIFFERENTIAL keyword may be used to explicitly define them. See the SUNDIALS documentation for further information.

PHASE_SPACE

Ask FeenoX to solve a set of algebraic-differntial equations and define the variables, vectors and/or matrices that span the phase space.

PHASE_SPACE PHASE_SPACE <vars> ... <vectors> ... <matrices> ...   

TIME_PATH

Force time-dependent problems to pass through specific instants of time.

TIME_PATH <expr_1> [ <expr_2>  [ ... <expr_n> ] ]  

The time step dt will be reduced whenever the distance between the current time t and the next expression in the list is greater than dt so as to force t to coincide with the expressions given. The list of expresssions should evaluate to a sorted list of values for all times.

PDE-RELATED KEYWORDS

BC

Define a boundary condition to be applied to faces, edges and/or vertices.

BC <name> [ MESH <name> ] [ PHYSICAL_GROUP <name_1>  PHYSICAL_GROUP <name_2> ... ] [ <bc_data1> <bc_data2> ... ]  

If the name of the boundary condition matches a physical group in the mesh, it is automatically linked to that physical group. If there are many meshes, the mesh this keyword refers to has to be given with MESH. If the boundary condition applies to more than one physical group in the mesh, they can be added using as many PHYSICAL_GROUP keywords as needed. If at least one PHYSICAL_GROUP is given explicitly, then the BC name is not used to try to implicitly link it to a physical group in the mesh. Each <bc_data> argument is a single string whose meaning depends on the type of problem being solved. For instance T=150*sin(x/pi) prescribes the temperature to depend on space as the provided expression in a thermal problem and fixed fixes the displacements in all the directions in a mechanical or modal problem. See the particular section on boundary conditions for further details.

COMPUTE_REACTION

Compute the reaction (force, moment, power, etc.) at selected face, edge or vertex.

COMPUTE_REACTION <physical_group> [ MOMENT [ X0 <expr> ] [ Y0 <expr> ] [ Z0 <expr> ] ] RESULT { <variable> | <vector> }  

If the MOMENT keyword is not given, the zero-th order reaction is computed, i.e. force in elasticity and power in thermal. If the MOMENT keyword is given, then the coordinates of the center can be given with X0, Y0 and Z0. If they are not, the moment is computed about the barycenter of the physical group. The resulting reaction will be stored in the variable (thermal) or vector (elasticity) provided. If the variable or vector does not exist, it will be created.

DUMP

Dump raw PETSc objects used to solve PDEs into files.

DUMP [ FORMAT { binary | ascii | octave } ] [ K |   K_bc |   b |   b_bc |   M |   M_bc |  

FIND_EXTREMA

Find and/or compute the absolute extrema of a function or expression over a mesh (or a subset of it).

FIND_EXTREMA { <expression> | <function> } [ OVER <physical_group> ] [ MESH <mesh_identifier> ] [ NODES | CELLS | GAUSS ]
 [ MIN <variable> ] [ MAX <variable> ] [ X_MIN <variable> ] [ X_MAX <variable> ] [ Y_MIN <variable> ] [ Y_MAX <variable> ] [ Z_MIN <variable> ] [ Z_MAX <variable> ] [ I_MIN <variable> ] [ I_MAX <variable> ]  

Either an expression or a function of space x, y and/or z should be given. By default the search is performed over the highest-dimensional elements of the mesh, i.e. over the whole volume, area or length for three, two and one-dimensional meshes, respectively. If the search is to be carried out over just a physical group, it has to be given in OVER. If there are more than one mesh defined, an explicit one has to be given with MESH. If neither NODES, CELLS or GAUSS is given then the search is performed over the three of them. With NODES only the function or expression is evalauted at the mesh nodes. With CELLS only the function or expression is evalauted at the element centers. With GAUSS only the function or expression is evalauted at the Gauss points. The value of the absolute minimum (maximum) is stored in the variable indicated by MIN (MAX). If the variable does not exist, it is created. The value of the x-y-z coordinate of the absolute minimum (maximum) is stored in the variable indicated by X_MIN-Y_MIN-Z_MIN (X_MAX-Y_MAX-Z_MAX). If the variable does not exist, it is created. The index (either node or cell) where the absolute minimum (maximum) is found is stored in the variable indicated by I_MIN (I_MAX).

INTEGRATE

Spatially integrate a function or expression over a mesh (or a subset of it).

INTEGRATE { <expression> | <function> } [ OVER <physical_group> ] [ MESH <mesh_identifier> ] [ GAUSS | CELLS ]
 RESULT <variable>
  

Either an expression or a function of space x, y and/or z should be given. If the integrand is a function, do not include the arguments, i.e. instead of f(x,y,z) just write f. The results should be the same but efficiency will be different (faster for pure functions). By default the integration is performed over the highest-dimensional elements of the mesh, i.e. over the whole volume, area or length for three, two and one-dimensional meshes, respectively. If the integration is to be carried out over just a physical group, it has to be given in OVER. If there are more than one mesh defined, an explicit one has to be given with MESH. Either GAUSS or CELLS define how the integration is to be performed. With GAUSS the integration is performed using the Gauss points and weights associated to each element type. With CELLS the integral is computed as the sum of the product of the integrand at the center of each cell (element) and the cell’s volume. Do expect differences in the results and efficiency between these two approaches depending on the nature of the integrand. The scalar result of the integration is stored in the variable given by the mandatory keyword RESULT. If the variable does not exist, it is created.

MATERIAL

Define a material its and properties to be used in volumes.

MATERIAL <name> [ MESH <name> ] [ LABEL <name_1>  [ LABEL <name_2> [ ... ] ] ] 
 [ <property_name_1>=<expr_1> [ <property_name_2>=<expr_2> [ ... ] ] ]  

If the name of the material matches a physical group in the mesh, it is automatically linked to that physical group. If there are many meshes, the mesh this keyword refers to has to be given with MESH. If the material applies to more than one physical group in the mesh, they can be added using as many LABEL keywords as needed. The names of the properties in principle can be arbitrary, but each problem type needs a minimum set of properties defined with particular names. For example, steady-state thermal problems need at least the conductivity which should be named k. If the problem is transient, it will also need heat capacity rhocp or diffusivity alpha. Mechanical problems need Young modulus E and Poisson’s ratio nu. Modal also needs density rho. Check the particular documentation for each problem type. Besides these mandatory properties, any other one can be defined. For instance, if one mandatory property dependend on the concentration of boron in the material, a new per-material property can be added named boron and then the function boron(x,y,z) can be used in the expression that defines the mandatory property.

PETSC_OPTIONS

Pass verbatim options to PETSc.

PETSC_OPTIONS "command-line options for PETSc"  

Options for PETSc can be passed either in at run time in the command line (run with -h to see how) or they can be set in the input file with PETSC_OPTIONS. This is handy when a particular problem is best suited to be solved using a particular set of options which can the be embedded into the problem definition. @

The string is passed verbatim to PETSc as if the options were set in the command line. Note that in this case, the string is passed verbatim to PETSc. This means that they are non-POSIX options but they have to be in the native PETSc format. That is to say, while in the command line one would give --ksp_view, here one has to give -ksp_view. Conversely, instead of --mg_levels_pc_type=sor one has to give -mg_levels_pc_type sor.

PHYSICAL_GROUP

Explicitly defines a physical group of elements on a mesh.

PHYSICAL_GROUP <name> [ MESH <name> ] [ DIMENSION <expr> ] [ ID <expr> ]
 [ MATERIAL <name> | | BC <name> [ BC ... ] ]
  

This keyword should seldom be needed. Most of the times, a combination of MATERIAL and BC ought to be enough for most purposes. The name of the PHYSICAL_GROUP keyword should match the name of the physical group defined within the input file. If there is no physical group with the provided name in the mesh, this instruction has no effect. If there are many meshes, an explicit mesh can be given with MESH. Otherwise, the physical group is defined on the main mesh. An explicit dimension of the physical group can be provided with DIMENSION. An explicit id can be given with ID. Both dimension and id should match the values in the mesh. For volumetric elements, physical groups can be linked to materials using MATERIAL. Note that if a material is created with the same name as a physical group in the mesh, they will be linked automatically, so there is no need to use PHYSCAL_GROUP for this. The MATERIAL keyword in PHYSICAL_GROUP is used to link a physical group in a mesh file and a material in the feenox input file with different names.

Likewise, for non-volumetric elements, physical groups can be linked to boundary using BC. As in the preceeding case, if a boundary condition is created with the same name as a physical group in the mesh, they will be linked automatically, so there is no need to use PHYSCAL_GROUP for this. The BC keyword in PHYSICAL_GROUP is used to link a physical group in a mesh file and a boundary condition in the feenox input file with different names. Note that while there can be only one MATERIAL associated to a physical group, there can be many BCs associated to a physical group.

PROBLEM

Ask FeenoX to solve a partial differential equation problem.

PROBLEM { laplace | mechanical | modal | neutron_diffusion | neutron_sn | thermal }
 [ 1D | 2D | 3D | DIM <expr> ] [ AXISYMMETRIC { x | y } ] 
 [ MESH <identifier> ] [ PROGRESS ] [ DO_NOT_DETECT_HANGING_NODES | DETECT_HANGING_NODES | HANDLE_HANGING_NODES ]
 [ DETECT_UNRESOLVED_BCS | ALLOW_UNRESOLVED_BCS ]
 [ PREALLOCATE ] [ ALLOW_NEW_NONZEROS ] [ CACHE_J ] [ CACHE_B ] 
 [ TRANSIENT | QUASISTATIC ] [ LINEAR | NON_LINEAR ] 
 [ MODES <expr> ] 
 [ PRECONDITIONER { gamg | mumps | lu | hypre | sor | bjacobi | cholesky | ... } ]
 [ LINEAR_SOLVER { gmres | mumps | bcgs | bicg | richardson | chebyshev | ... } ]
 [ NONLINEAR_SOLVER { newtonls | newtontr | nrichardson | ngmres | qn | ngs | ... } ]
 [ TRANSIENT_SOLVER { bdf | beuler | arkimex | rosw | glee | ... } ]
 [ TIME_ADAPTATION { basic | none | dsp | cfl | glee | ... } ]
 [ EIGEN_SOLVER { krylovschur | lanczos | arnoldi | power | gd | ... } ]
 [ SPECTRAL_TRANSFORMATION { shift | sinvert | cayley | ... } ]
 [ EIGEN_FORMULATION { omega | lambda } ]
 [ DIRICHLET_SCALING { absolute <expr> | relative <expr> } ]
  

Currently, FeenoX can solve the following types of PDE-casted problems:

The number of spatial dimensions of the problem needs to be given either as 1d, 2d, 3d or after the keyword DIM. Alternatively, one can define a MESH with an explicit DIMENSIONS keyword before PROBLEM. Default is 3D. If the AXISYMMETRIC keyword is given, the mesh is expected to be two-dimensional in the x-y plane and the problem is assumed to be axi-symmetric around the given axis. If there are more than one MESHes defined, the one over which the problem is to be solved can be defined by giving the explicit mesh name with MESH. By default, the first mesh to be defined in the input file with READ_MESH (which can be defined after the PROBLEM keyword) is the one over which the problem is solved. If the keyword PROGRESS is given, three ASCII lines will show in the terminal the progress of the ensamble of the stiffness matrix (or matrices), the solution of the system of equations and the computation of gradients (stresses, heat fluxes, etc.), if applicable. If either DETECT_HANGING_NODES or HANDLE_HANGING_NODES are given, an intermediate check for nodes without any associated elements will be performed. For well-behaved meshes this check is redundant so by detault it is not done (DO_NOT_DETEC_HANGING_NODES). With DETECT_HANGING_NODES, FeenoX will report the tag of the hanging nodes and stop. With HANDLE_HANGING_NODES, FeenoX will fix those nodes and try to solve the problem anyway. By default, FeenoX checks that all physical groups referred to in the BC keywors exists (DETECT_UNRESOLVED_BCS). If ALLOW_UNRESOLVED_BCS is given, FeenoX will ignore unresolved boundary conditions instead of complaining. This is handy when using the same input for different meshes which might have different groups, for example solving the same problem using a full geometry or a symmetric geometry. The latter should have at least one symmetry BC whilst the former does not. If the special variable end_time is zero, FeenoX solves a static problem—although the variable static_steps is still honored. If end_time is non-zero, FeenoX solves a transient or quasistatic problem. This can be controlled by TRANSIENT or QUASISTATIC. By default FeenoX tries to detect wheter the computation should be linear or non-linear. An explicit mode can be set with either LINEAR on NON_LINEAR. The number of modes to be computed when solving eigenvalue problems is given by MODES. The default value is problem dependent. The preconditioner (PC), linear (KSP), non-linear (SNES) and time-stepper (TS) solver types be any of those available in PETSc (first option is the default):

If the EIGEN_FORMULATION is omega then  = ω^2^ is solved, and  = λKϕ if it is lambda. Default is lambda, although some particular PDEs might change it (for example free-free modal switches to omega). The EIGEN_DIRICHLET_ZERO keyword controls which of the matrices has a zero and which one has a non-zero in the diagonal when setting Dirichlet boundary conditions. Default is M, i.e. matrix K has a non-zero and matrix M has a zero. This setting, along with EIGEN_FORMULATION determines which spectral transforms can a cannot be used: you cannot invert the matrix with the zero in the diagonal. The DIRICHLET_SCALING keyword controls the way Dirichlet boundary conditions are scaled when computing the residual. Roughly, it defines how to compute the parameter α. If absolute, then α is equal to the given expression. If relative, then α is equal to the given fraction of the average diagonal entries in the stiffness matrix. Default is α = 1.

READ_MESH

Read an unstructured mesh and (optionally) functions of space-time from a file.

READ_MESH { <file_path> | <file_id> } [ DIM <num_expr> ]
 [ SCALE <expr> ] [ OFFSET <expr_x> <expr_y> <expr_z> ]
 [ INTEGRATION { full | reduced } ]
 [ MAIN ] [ UPDATE_EACH_STEP ]
 [ READ_FIELD <name_in_mesh> AS <function_name> ] [ READ_FIELD ... ] 
 [ READ_FUNCTION <function_name> ] [READ_FUNCTION ...] 
  

Either a file identifier (defined previously with a FILE keyword) or a file path should be given. The format is read from the extension, which should be either

Note than only MSH is suitable for defining PDE domains, as it is the only one that provides physical groups (a.k.a labels) which are needed in order to define materials and boundary conditions. The other formats are primarily supported to read function data contained in the file and eventually, to operate over these functions (i.e. take differences with other functions contained in other files to compare results). The file path or file id can be used to refer to a particular mesh when reading more than one, for instance in a WRITE_MESH or INTEGRATE keyword. If a file path is given such as cool_mesh.msh, it can be later referred to as either cool_mesh.msh or just cool_mesh.

The spatial dimensions can be given with DIM. If material properties are uniform and given with variables, the number of dimensions are not needed and will be read from the file at runtime. But if either properties are given by spatial functions or if functions are to be read from the mesh with READ_DATA or READ_FUNCTION, then the number of dimensions ought to be given explicitly because FeenoX needs to know how many arguments these functions take. If either OFFSET and/or SCALE are given, the node locations are first shifted and then scaled by the provided values. When defining several meshes and solving a PDE problem, the mesh used as the PDE domain is the one marked with MAIN. If none of the meshes is explicitly marked as main, the first one is used. If UPDATE_EACH_STEP is given, then the mesh data is re-read from the file at each time step. Default is to read the mesh once, except if the file path changes with time. For each READ_FIELD keyword, a point-wise defined scalar function of space named <function_name> is defined and filled with the scalar data named <name_in_mesh> contained in the mesh file. The READ_FUNCTION keyword is a shortcut when the scalar name and the to-be-defined function are the same. If no mesh is marked as MAIN, the first one is the main one.

SOLVE_PROBLEM

Explicitly solve the PDE problem.

SOLVE_PROBLEM  

Whenever the instruction SOLVE_PROBLEM is executed, FeenoX solves the PDE problem. For static problems, that means solving the equations and filling in the result functions. For transient or quasisstatic problems, that means advancing one time step.

WRITE_MESH

Write a mesh and/or generic functions of space-time to a post-processing file.

WRITE_MESH <file> [ MESH <mesh_identifier> ] [ FILE_FORMAT { gmsh | vtk } ] [ NO_MESH ] [ NO_PHYSICAL_NAMES ]
  [ NODE | CELL ] [ <printf_specification> ]
 [ <scalar_field_1> ] [ <scalar_field_2> ] [...] 
 [ VECTOR [ NAME <name> ] <field_x> <field_y> <field_z> ] [...] 
 [ SYMMETRIC_TENSOR [ NAME <name> ] <field_xx> <field_yy> <field_zz> <field_xy> <field_yz> <field_zx> ] [...] 
  

The format is automatically detected from the extension, which should be either msh (version 2.2 ASCII) or vtk (legacy ASCII). Otherwise, the keyword FILE_FORMAT has to be given to set the format explicitly. If there are several meshes defined by READ_MESH, the mesh used to write the data has be given explicitly with MESH. If the NO_MESH keyword is given, only the results are written into the output file without any mesh data. Depending on the output format, this can be used to avoid repeating data and/or creating partial output files which can the be latter assembled by post-processing scripts. When targetting the .msh output format, if NO_PHYSICAL_NAMES is given then the section that sets the actual names of the physical entities is not written.

This might be needed in some cases to avoid name clashes when dealing with multiple .msh files. The output is node-based by default. This can be controlled with both the NODE and CELL keywords. All fields that come after a NODE (CELL) keyword will be written at the node (cells). These keywords can be used several times and mixed with fields. For example CELL k(x,y,z) NODE T sqrt(x^2+y^2) CELL 1+z will write the conductivity and the expression 1 + z as cell-based and the temperature T(x,y,z) and the expression $\sqrt{x^2+y^2}$ as a node-based fields. If a printf-like format specifier starting with % is given, that format is used for the fields that follow. Make sure the format reads floating-point data, i.e. do not use %d. Default is %g. The data to be written has to be given as a list of fields, i.e. distributions (such as k or E), functions of space (such as T) and/or expressions (such as T(x,y,z)*sqrt(x^2+y^2+z^2)). Each field is written as a scalar, unless either the keywords VECTOR or SYMMETRIC_TENSOR are given. In the first case, the next three fields following the VECTOR keyword are taken as the vector elements. In the latter, the next six fields following the SYMMETRIC_TENSOR keyword are taken as the tensor elements.

WRITE_RESULTS

Write the problem mesh and problem results to a file for post-processing.

WRITE_RESULTS [ FORMAT { gmsh | vtk } ] [ FILE <file> ]
 [ NO_PHYSICAL_NAMES ] [ <printf_specification> ]
  

Default format is gmsh. If no FILE is provided, the output file is the same as the input file replacing the .fee extension with the format extension, i.e. $0.msh. If there are further optional command line arguments, they are added pre-pending a dash, i.e. $0-[$1-[$2...]].msh Otherwise the given FILE is used. If no explicit FORMAT is given, the format is read from the FILE extension. When targetting the .msh output format, if NO_PHYSICAL_NAMES is given then the section that sets the actual names of the physical entities is not written.

This might be needed in some cases to avoid name clashes when dealing with multiple .msh files. If a printf-like format specifier starting with % is given, that format is used for the fields that follow. Make sure the format reads floating-point data, i.e. do not use %d. Default is %g.

SPECIAL VARIABLES

done

Flag that indicates whether the overall calculation is over.

This variable is set to true by FeenoX when the computation finished so it can be checked in an IF block to do something only in the last step. But this variable can also be set to true from the input file, indicating that the current step should also be the last one. For example, one can set end_time = infinite and then finish the computation at t = 10 by setting done = t > 10. This done variable can also come from (and sent to) other sources, like a shared memory object for coupled calculations.

done_static

Flag that indicates whether the static calculation is over or not.

It is set to true (i.e.  ≠ 0) by feenox if step_static ≥ static_steps. If the user sets it to true, the current step is marked as the last static step and the static calculation ends after finishing the step. It can be used in IF blocks to check if the static step is finished or not.

done_transient

Flag that indicates whether the transient calculation is over or not.

It is set to true (i.e.  ≠ 0) by feenox if t ≥ end_time. If the user sets it to true, the current step is marked as the last transient step and the transient calculation ends after finishing the step. It can be used in IF blocks to check if the transient steps are finished or not.

dt

Actual value of the time step for transient calculations.

When solving DAE systems, this variable is set by feenox. It can be written by the user for example by importing it from another transient code by means of shared-memory objects. Care should be taken when solving DAE systems and overwriting t. Default value is DEFAULT_DT, which is a power of two and roundoff errors are thus reduced.

end_time

Final time of the transient calculation, to be set by the user.

The default value is zero, meaning no transient calculation.

i

Dummy index, used mainly in vector and matrix row subindex expressions.

infinite

A very big positive number.

It can be used as end_time = infinite or to define improper integrals with infinite limits. Default is 2^50^ ≈ 1 × 10^15^.

in_static

Flag that indicates if FeenoX is solving the iterative static calculation.

This is a read-only variable that is non zero if the static calculation.

in_static_first

Flag that indicates if feenox is in the first step of the iterative static calculation.

in_static_last

Flag that indicates if feenox is in the last step of the iterative static calculation.

in_transient

Flag that indicates if feenox is solving transient calculation.

in_transient_first

Flag that indicates if feenox is in the first step of the transient calculation.

in_transient_last

Flag that indicates if feenox is in the last step of the transient calculation.

j

Dummy index, used mainly in matrix column subindex expressions.

max_dt

Maximum bound for the time step that feenox should take when solving DAE systems.

min_dt

Minimum bound for the time step that feenox should take when solving DAE systems.

mpi_rank

The current rank in an MPI execution. Mind the PRINTF_ALL instruction.

mpi_size

The number of total ranks in an MPI execution.

on_gsl_error

This should be set to a mask that indicates how to proceed if an error ir raised in any routine of the GNU Scientific Library.

on_ida_error

This should be set to a mask that indicates how to proceed if an error ir raised in any routine of the SUNDIALS Library.

on_nan

This should be set to a mask that indicates how to proceed if Not-A-Number signal (such as a division by zero) is generated when evaluating any expression within feenox.

pi

A double-precision floating point representaion of the number π

It is equal to the M_PI constant in math.h .

pid

The Unix process id of the FeenoX instance.

static_steps

Number of steps that ought to be taken during the static calculation, to be set by the user.

The default value is one, meaning only one static step.

step_static

Indicates the current step number of the iterative static calculation.

This is a read-only variable that contains the current step of the static calculation.

step_transient

Indicates the current step number of the transient static calculation.

This is a read-only variable that contains the current step of the transient calculation.

t

Actual value of the time for transient calculations.

This variable is set by FeenoX, but can be written by the user for example by importing it from another transient code by means of shared-memory objects. Care should be taken when solving DAE systems and overwriting t.

zero

A very small positive number.

It is taken to avoid roundoff errors when comparing floating point numbers such as replacing a ≤ a~max~ with a < a~max~+ zero. Default is (1/2)^−50^ ≈ 9 × 10^−16^ .

MATERIAL PROPERTIES

TBD.

BOUNDARY CONDITIONS

TBD.

RESULTING DISTRIBUTIONS

TBD.

BUILT-IN FUNCTIONS

abs

Returns the absolute value of the argument x.

abs(x)  

acos

Computes the arc in radians whose cosine is equal to the argument x. A NaN error is raised if |x| > 1.

acos(x)  

asin

Computes the arc in radians whose sine is equal to the argument x. A NaN error is raised if |x| > 1.

asin(x)  

atan

Computes, in radians, the arc tangent of the argument x.

atan(x)  

atan2

Computes, in radians, the arc tangent of quotient y/x, using the signs of the two arguments to determine the quadrant of the result, which is in the range [−π,π].

atan2(y,x)  

ceil

Returns the smallest integral value not less than the argument x.

ceil(x)  

clock

Returns the value of a certain clock in seconds measured from a certain (but specific) milestone. The kind of clock and the initial milestone depend on the optional integer argument f. It defaults to one, meaning CLOCK_MONOTONIC. The list and the meanings of the other available values for f can be checked in the clock_gettime (2) system call manual page.

clock([f])  

cos

Computes the cosine of the argument x, where x is in radians. A cosine wave can be generated by passing as the argument x a linear function of time such as ωt + ϕ, where ω controls the frequency of the wave and ϕ controls its phase.

cos(x)  

cosh

Computes the hyperbolic cosine of the argument x, where x is in radians.

cosh(x)  

cpu_time

Returns the CPU time used by the local FeenoX rank, in seconds. If the optional argument f is not provided or it is zero (default), the sum of times for both user-space and kernel-space usage is returned. For f=1 only user time is returned. For f=2 only system time is returned.

cpu_time([f])  

d_dt

Computes the time derivative of the expression given in the argument x during a transient problem using the difference between the value of the signal in the previous time step and the actual value divided by the time step δt stored in dt. The argument x does not neet to be a variable, it can be an expression involving one or more variables that change in time. For t = 0, the return value is zero. Unlike the functional derivative, the full dependence of these variables with time does not need to be known beforehand, i.e. the expression x might involve variables read from a shared-memory object at each time step.

d_dt(x)  

deadband

Filters the first argument x with a deadband centered at zero with an amplitude given by the second argument a.

deadband(x, a)  

equal

Checks if the two first expressions a and b are equal, up to the tolerance given by the third optional argument ϵ. If either |a| > 1 or |b| > 1, the arguments are compared using GSL’s gsl_fcmp, otherwise the absolute value of their difference is compared against ϵ. This function returns zero if the arguments are not equal and one otherwise. Default value for ϵ = 10^−9^.

equal(a, b, [eps])  

exp

Computes the exponential function the argument x, i.e. the base of the natural logarithm e raised to the x-th power.

exp(x)  

expint1

Computes the first exponential integral function of the argument x. If x is zero, a NaN error is issued.

expint1(x)  

expint2

Computes the second exponential integral function of the argument x.

expint2(x)  

expint3

Computes the third exponential integral function of the argument x.

expint3(x)  

expintn

Computes the n-th exponential integral function of the argument x. If n is zero or one and x is zero, a NaN error is issued.

expintn(n,x)  

floor

Returns the largest integral value not greater than the argument x.

floor(x)  

gammaf

Computes the Gamma function Γ(x).

gammaf(x)  

heaviside

Computes the zero-centered Heaviside step function of the argument x. If the optional second argument δ is provided, the discontinuous step at x = 0 is replaced by a ramp starting at x = 0 and finishing at x = δ.

heaviside(x, [delta])  

if

Performs a conditional testing of the first argument a, and returns either the second optional argument b if a is different from zero or the third optional argument c if a evaluates to zero. The comparison of the condition a with zero is performed within the precision given by the optional fourth argument ϵ. If the second argument c is not given and a is not zero, the function returns one. If the third argument c is not given and a is zero, the function returns zero. The default precision is ϵ = 10^−9^. Even though if is a logical operation, all the arguments and the returned value are double-precision floating point numbers.

if(a, [b], [c], [eps])  

integral_dt

Computes the time integral of the expression given in the argument x during a transient problem with the trapezoidal rule using the value of the signal in the previous time step and the current value. At t = 0 the integral is initialized to zero. Unlike the functional integral, the full dependence of these variables with time does not need to be known beforehand, i.e. the expression x might involve variables read from a shared-memory object at each time step.

integral_dt(x)  

integral_euler_dt

Idem as integral_dt but uses the backward Euler rule to update the instantaenous integral value. This function is provided in case this particular way of approximating time integrals is needed, for instance to compare FeenoX solutions with other computer codes. In general, it is recommended to use integral_dt.

integral_euler_dt(x)  

is_even

Returns one if the argument x rounded to the nearest integer is even.

is_even(x)  

is_in_interval

Returns true if the argument x is in the interval [a, b), i.e. including a but excluding b.

is_in_interval(x, a, b)  

is_odd

Returns one if the argument x rounded to the nearest integer is odd.

is_odd(x)  

j0

Computes the regular cylindrical Bessel function of zeroth order evaluated at the argument x.

j0(x)  

lag

Filters the first argument x(t) with a first-order lag of characteristic time τ, i.e. this function applies the transfer function $G(s) = \frac{1}{1 + s\tau}$ to the time-dependent signal x(t) to obtain a filtered signal y(t), by assuming that it is constant during the time interval [tΔt,t] and using the analytical solution of the differential equation for that case at t = Δt with the initial condition y(0) = y(tΔt).

lag(x, tau)  

lag_bilinear

Filters the first argument x(t) with a first-order lag of characteristic time τ to the time-dependent signal x(t) by using the bilinear transformation formula.

lag_bilinear(x, tau)  

lag_euler

Filters the first argument x(t) with a first-order lag of characteristic time τ to the time-dependent signal x(t) by using the Euler forward rule.

lag_euler(x, tau)  

last

Returns the value the variable x had in the previous time step. This function is equivalent to the Z-transform operator “delay” denoted by z^−1^[x]. For t = 0 the function returns the actual value of x. The optional flag p should be set to one if the reference to last is done in an assignment over a variable that already appears inside expression x such as x = last(x). See example number 2.

last(x,[p])  

limit

Limits the first argument x to the interval [a,b]. The second argument a should be less than the third argument b.

limit(x, a, b)  

limit_dt

Limits the value of the first argument x(t) so to that its time derivative is bounded to the interval [a,b]. The second argument a should be less than the third argument b.

limit_dt(x, a, b)  

log

Computes the natural logarithm of the argument x. If x is zero or negative, a NaN error is issued.

log(x)  

mark_max

Returns the integer index i of the maximum of the arguments x~i~ provided. Currently only maximum of ten arguments can be provided.

mark_max(x1, x2, [...], [x10])  

mark_min

Returns the integer index i of the minimum of the arguments x~i~ provided. Currently only maximum of ten arguments can be provided.

mark_max(x1, x2, [...], [x10])  

max

Returns the maximum of the arguments x~i~ provided. Currently only maximum of ten arguments can be given.

max(x1, x2, [...], [x10])  

memory

Returns the maximum memory (resident set size) used by FeenoX, in Gigabytes.

memory()  

min

Returns the minimum of the arguments x~i~ provided. Currently only maximum of ten arguments can be given.

min(x1, x2, [...], [x10])  

mod

Returns the remainder of the division between the first argument a and the second one b. Both arguments may be non-integral.

mod(a, b)  

mpi_memory_local

Returns the memory usage as reported by PETSc in the give rank, in Gigabytes. If no rank is given, each rank returns a local value which should be printed with PRINTF_ALL. Returns the memory global usage as reported by PETSc summing over all ranks, in Gigabytes.

mpi_memory_local([rank]) mpi_memory_global()  

not

Returns one if the first argument x is zero and zero otherwise. The second optional argument ϵ gives the precision of the “zero” evaluation. If not given, default is ϵ = 10^−9^.

not(x, [eps])  

random

Returns a random real number uniformly distributed between the first real argument x~1~ and the second one x~2~. If the third integer argument s is given, it is used as the seed and thus repetitive sequences can be obtained. If no seed is provided, the current time (in seconds) plus the internal address of the expression is used. Therefore, two successive calls to the function without seed (hopefully) do not give the same result. This function uses a second-order multiple recursive generator described by Knuth in Seminumerical Algorithms, 3rd Ed., Section 3.6.

random(x1, x2, [s])  

random_gauss

Returns a random real number with a Gaussian distribution with a mean equal to the first argument x~1~ and a standard deviation equatl to the second one x~2~. If the third integer argument s is given, it is used as the seed and thus repetitive sequences can be obtained. If no seed is provided, the current time (in seconds) plus the internal address of the expression is used. Therefore, two successive calls to the function without seed (hopefully) do not give the same result. This function uses a second-order multiple recursive generator described by Knuth in Seminumerical Algorithms, 3rd Ed., Section 3.6.

random_gauss(x1, x2, [s])  

round

Rounds the argument x to the nearest integer. Halfway cases are rounded away from zero.

round(x)  

sawtooth_wave

Computes a sawtooth wave between zero and one with a period equal to one. As with the sine wave, a sawtooh wave can be generated by passing as the argument x a linear function of time such as ωt + ϕ, where ω controls the frequency of the wave and ϕ controls its phase.

sawtooth_wave(x)  

sech

Computes the hyperbolic secant of the argument x, where x is in radians.

sech(x)  

sgn

Returns minus one, zero or plus one depending on the sign of the first argument x. The second optional argument ϵ gives the precision of the “zero” evaluation. If not given, default is ϵ = 10^−9^.

sgn(x, [eps])  

sin

Computes the sine of the argument x, where x is in radians. A sine wave can be generated by passing as the argument x a linear function of time such as ωt + ϕ, where ω controls the frequency of the wave and ϕ controls its phase.

sin(x)  

sinh

Computes the hyperbolic sine of the argument x, where x is in radians.

sinh(x)  

sqrt

Computes the positive square root of the argument x. If x is negative, a NaN error is issued.

sqrt(x)  

square_wave

Computes a square function between zero and one with a period equal to one. The output is one for 0 < x < 1/2 and zero for 1/2 ≤ x < 1. As with the sine wave, a square wave can be generated by passing as the argument x a linear function of time such as ωt + ϕ, where ω controls the frequency of the wave and ϕ controls its phase.

square_wave(x)  

tan

Computes the tangent of the argument x, where x is in radians.

tan(x)  

tanh

Computes the hyperbolic tangent of the argument x, where x is in radians.

tanh(x)  

threshold_max

Returns one if the first argument x is greater than the threshold given by the second argument a, and zero otherwise. If the optional third argument b is provided, an hysteresis of width b is needed in order to reset the function value. Default is no hysteresis, i.e. b = 0.

threshold_max(x, a, [b])  

threshold_min

Returns one if the first argument x is less than the threshold given by the second argument a, and zero otherwise. If the optional third argument b is provided, an hysteresis of width b is needed in order to reset the function value. Default is no hysteresis, i.e. b = 0.

threshold_min(x, a, [b])  

triangular_wave

Computes a triangular wave between zero and one with a period equal to one. As with the sine wave, a triangular wave can be generated by passing as the argument x a linear function of time such as ωt + ϕ, where ω controls the frequency of the wave and ϕ controls its phase.

triangular_wave(x)  

wall_time

Returns the time ellapsed since the invocation of FeenoX, in seconds.

wall_time()  

BUILT-IN FUNCTIONALS

TBD.

BUILT-IN VECTOR FUNCTIONS

TBD.

NOTES

TBD.

BUGS

Report on GitHub https://github.com/seamplex/feenox or at jeremy@seamplex.com

SEE ALSO

gmsh(1), mpirun(1), paraview(1)

The FeenoX web page contains links to the full source code, binary versions, updates, examples, verification & validation cases and full documentation: https://www.seamplex.com/feenox.

AUTHORS

Jeremy Theler jeremy@seamplex.com.