ffc package

Subpackages

Submodules

ffc.analysis module

ffc.analysis.analyze_coordinate_mappings(coordinate_elements, parameters)[source]
ffc.analysis.analyze_elements(elements, parameters)[source]
ffc.analysis.analyze_forms(forms, parameters)[source]

Analyze form(s), returning

form_datas - a tuple of form_data objects unique_elements - a tuple of unique elements across all forms element_numbers - a mapping to unique numbers for all elements
ffc.analysis.analyze_ufl_objects(ufl_objects, kind, parameters)[source]

Analyze ufl object(s), either forms, elements, or coordinate mappings, returning:

form_datas - a tuple of form_data objects unique_elements - a tuple of unique elements across all forms element_numbers - a mapping to unique numbers for all elements

ffc.classname module

This module defines some basics for generating C++ code.

ffc.classname.make_classname(prefix, basename, signature)[source]
ffc.classname.make_integral_classname(prefix, integral_type, form_id, subdomain_id)[source]

ffc.codegeneration module

Compiler stage 4: Code generation

This module implements the generation of C++ code for the body of each UFC function from an (optimized) intermediate representation (OIR).

ffc.codegeneration.generate_code(ir, parameters)[source]

Generate code from intermediate representation.

ffc.compiler module

This is the compiler, acting as the main interface for compilation of forms and breaking the compilation into several sequential stages. The output of each stage is the input of the next stage.

Compiler stage 0: Language, parsing

Input: Python code or .ufl file Output: UFL form

This stage consists of parsing and expressing a form in the UFL form language.

This stage is completely handled by UFL.

Compiler stage 1: Analysis

Input: UFL form Output: Preprocessed UFL form and FormData (metadata)

This stage preprocesses the UFL form and extracts form metadata. It may also perform simplifications on the form.

Compiler stage 2: Code representation

Input: Preprocessed UFL form and FormData (metadata) Output: Intermediate Representation (IR)

This stage examines the input and generates all data needed for code generation. This includes generation of finite element basis functions, extraction of data for mapping of degrees of freedom and possible precomputation of integrals.

Most of the complexity of compilation is handled in this stage.

The IR is stored as a dictionary, mapping names of UFC functions to data needed for generation of the corresponding code.

Compiler stage 3: Optimization

Input: Intermediate Representation (IR) Output: Optimized Intermediate Representation (OIR)

This stage examines the IR and performs optimizations.

Optimization is currently disabled as a separate stage but is implemented as part of the code generation for quadrature representation.

Compiler stage 4: Code generation

Input: Optimized Intermediate Representation (OIR) Output: C++ code

This stage examines the OIR and generates the actual C++ code for the body of each UFC function.

The code is stored as a dictionary, mapping names of UFC functions to strings containing the C++ code of the body of each function.

Compiler stage 5: Code formatting

Input: C++ code Output: C++ code files

This stage examines the generated C++ code and formats it according to the UFC format, generating as output one or more .h/.cpp files conforming to the UFC format.

The main interface is defined by the following two functions:

compile_form compile_element

The compiler stages are implemented by the following functions:

analyze_forms or analyze_elements (stage 1) compute_ir (stage 2) optimize_ir (stage 3) generate_code (stage 4) format_code (stage 5)
ffc.compiler.compile_form(forms, object_names=None, prefix='Form', parameters=None, jit=False)[source]

This function generates UFC code for a given UFL form or list of UFL forms.

ffc.compiler.compile_element(elements, object_names=None, prefix='Element', parameters=None, jit=False)[source]

This function generates UFC code for a given UFL element or list of UFL elements.

ffc.fiatinterface module

class ffc.fiatinterface.SpaceOfReals[source]

Bases: object

Constant over the entire domain, rather than just cellwise.

ffc.fiatinterface.create_element(ufl_element)[source]
ffc.fiatinterface.create_quadrature(shape, degree, scheme='default')[source]

Generate quadrature rule (points, weights) for given shape that will integrate an polynomial of order ‘degree’ exactly.

ffc.fiatinterface.map_facet_points(points, facet, cellname)[source]

Map points from the e (UFC) reference simplex of dimension d - 1 to a given facet on the (UFC) reference simplex of dimension d. This may be used to transform points tabulated for example on the 2D reference triangle to points on a given facet of the reference tetrahedron.

ffc.fiatinterface.reference_cell(cellname)[source]

Return FIAT reference cell

ffc.fiatinterface.reference_cell_vertices(cellname)[source]

Return dict of coordinates of reference cell vertices for this ‘cellname’.

ffc.formatting module

Compiler stage 5: Code formatting

This module implements the formatting of UFC code from a given dictionary of generated C++ code for the body of each UFC function.

It relies on templates for UFC code available as part of the module ufc_utils.

ffc.formatting.format_code(code, wrapper_code, prefix, parameters, jit=False)[source]

Format given code in UFC format. Returns two strings with header and source file contents.

ffc.formatting.generate_factory_functions(prefix, kind, classname)[source]
ffc.formatting.generate_jit_factory_functions(code, prefix)[source]
ffc.formatting.write_code(code_h, code_c, prefix, parameters)[source]

ffc.git_commit_hash module

ffc.git_commit_hash.git_commit_hash()[source]

Return git changeset hash (returns “unknown” if changeset is not known)

ffc.jitcompiler module

This module provides a just-in-time (JIT) form compiler. It uses dijitso to wrap the generated code into a Python module.

exception ffc.jitcompiler.FFCError[source]

Bases: exceptions.Exception

exception ffc.jitcompiler.FFCJitError[source]

Bases: ffc.jitcompiler.FFCError

ffc.jitcompiler.compute_jit_prefix(ufl_object, parameters, kind=None)[source]

Compute the prefix (module name) for jit modules.

ffc.jitcompiler.jit(ufl_object, parameters=None, indirect=False)[source]

Just-in-time compile the given form or element

Parameters:

ufl_object : The UFL object to be compiled parameters : A set of parameters
ffc.jitcompiler.jit_build(ufl_object, module_name, parameters)[source]

Wraps dijitso jit with some parameter conversion etc.

ffc.jitcompiler.jit_generate(ufl_object, module_name, signature, parameters)[source]

Callback function passed to dijitso.jit: generate code and return as strings.

ffc.log module

This module provides functions used by the FFC implementation to output messages. These may be redirected by the user of FFC.

This module reuses the corresponding log.py module from UFL which is a wrapper for the standard Python logging module.

ffc.log.add_indent(*message)
ffc.log.add_logfile(*message)
ffc.log.begin(*message)
ffc.log.debug(*message)
ffc.log.debug_code(code, name='')[source]

Debug generated code.

ffc.log.debug_dict(d, title='')[source]

Pretty-print dictionary.

ffc.log.debug_ir(ir, name='')[source]

Debug intermediate representation.

ffc.log.deprecate(*message)
ffc.log.end(*message)
ffc.log.error(*message)
ffc.log.ffc_assert(condition, *message)[source]

Assert that condition is true and otherwise issue an error with given message.

ffc.log.get_handler(*message)
ffc.log.get_logger(*message)
ffc.log.info(*message)
ffc.log.info_blue(*message)
ffc.log.info_green(*message)
ffc.log.info_red(*message)
ffc.log.log(*message)
ffc.log.pop_level(*message)
ffc.log.push_level(*message)
ffc.log.set_handler(*message)
ffc.log.set_indent(*message)
ffc.log.set_level(*message)
ffc.log.set_prefix(*message)
ffc.log.warning(*message)
ffc.log.warning_blue(*message)
ffc.log.warning_green(*message)
ffc.log.warning_red(*message)

ffc.main module

This script is the command-line interface to FFC.

It parses command-line arguments and generates code from input UFL form files.

ffc.main.compile_ufl_data(ufd, prefix, parameters)[source]
ffc.main.info_usage()[source]

Print usage information.

ffc.main.info_version()[source]

Print version number.

ffc.main.main(args=None)[source]

This is the commandline tool for the python module ffc.

ffc.main.print_error(msg)[source]

Print error message (cannot use log system at top level).

ffc.optimization module

Compiler stage 5: optimization

This module implements the optimization of an intermediate code representation.

ffc.optimization.optimize_ir(ir, parameters)[source]

Optimize intermediate form representation.

ffc.parameters module

ffc.parameters.compilation_relevant_parameters(parameters)[source]
ffc.parameters.compute_jit_parameters_signature(parameters)[source]

Return parameters signature (some parameters must be ignored).

ffc.parameters.default_jit_parameters()[source]
ffc.parameters.default_parameters()[source]

Return (a copy of) the default parameter values for FFC.

ffc.parameters.split_parameters(parameters)[source]

Split a parameters dict into groups based on what parameters are used for.

ffc.parameters.validate_jit_parameters(parameters)[source]

Check parameters and add any missing parameters

ffc.parameters.validate_parameters(parameters)[source]

Initial check of parameters.

ffc.plot module

This module provides functionality for plotting finite elements.

ffc.plot.plot(element, rotate=True)[source]

Plot finite element.

ffc.representation module

Compiler stage 2: Code representation

This module computes intermediate representations of forms, elements and dofmaps. For each UFC function, we extract the data needed for code generation at a later stage.

The representation should conform strictly to the naming and order of functions in UFC. Thus, for code generation of the function “foo”, one should only need to use the data stored in the intermediate representation under the key “foo”.

ffc.representation.all_elements(fiat_element)[source]
ffc.representation.cell_midpoint(cell)[source]
ffc.representation.compute_ir(analysis, prefix, parameters, jit=False)[source]

Compute intermediate representation.

ffc.representation.make_all_element_classnames(prefix, elements, coordinate_elements, element_numbers, parameters, jit)[source]
ffc.representation.make_coordinate_mapping_jit_classname(ufl_mesh, parameters)[source]
ffc.representation.make_dofmap_jit_classname(ufl_element, parameters)[source]
ffc.representation.make_finite_element_jit_classname(ufl_element, parameters)[source]
ffc.representation.needs_oriented_jacobian(fiat_element)[source]
ffc.representation.pick_representation(representation)[source]

Return one of the specialized code generation modules from a representation string.

ffc.representation.uses_integral_moments(fiat_element)[source]

True if element uses integral moments for its degrees of freedom.

ffc.representationutils module

This module contains utility functions for some code shared between quadrature and tensor representation.

ffc.representationutils.create_quadrature_points_and_weights(integral_type, cell, degree, rule)[source]

Create quadrature rule and return points and weights.

ffc.representationutils.entity_type_from_integral_type(integral_type)[source]
ffc.representationutils.generate_enabled_coefficients(enabled_coefficients)[source]
ffc.representationutils.initialize_integral_code(ir, prefix, parameters)[source]

Representation independent default initialization of code dict for integral from intermediate representation.

ffc.representationutils.initialize_integral_ir(representation, itg_data, form_data, form_id)[source]

Initialize a representation dict with common information that is expected independently of which representation is chosen.

ffc.representationutils.integral_type_to_entity_dim(integral_type, tdim)[source]

Given integral_type and domain tdim, return the tdim of the integration entity.

ffc.representationutils.map_integral_points(points, integral_type, cell, entity)[source]

Map points from reference entity to its parent reference cell.

ffc.representationutils.needs_oriented_jacobian(form_data)[source]
ffc.representationutils.transform_component(component, offset, ufl_element)[source]

This function accounts for the fact that if the geometrical and topological dimension does not match, then for native vector elements, in particular the Piola-mapped ones, the physical value dimensions and the reference value dimensions are not the same. This has certain consequences for mixed elements, aka ‘fun with offsets’.

ffc.utils module

ffc.utils.all_equal(sequence)[source]

Check that all items in list are equal.

ffc.utils.compute_permutations(k, n, skip=[])[source]

Compute all permutations of k elements from (0, n) in rising order. Any elements that are contained in the list skip are not included.

ffc.utils.insert_nested_dict(root, keys, value)[source]

Set root[keys[0]][...][keys[-1]] = value, creating subdicts on the way if missing.

ffc.utils.listcopy(sequence)[source]

Create a copy of the list, calling the copy constructor on each object in the list (problems when using copy.deepcopy).

ffc.utils.pick_first(sequence)[source]

Check that all values are equal and return the value.

ffc.wrappers module

ffc.wrappers.generate_wrapper_code(analysis, prefix, object_names, parameters)[source]

Generate code for additional wrappers.

Module contents

FEniCS Form Compiler (FFC)

FFC compiles finite element variational forms into C++ code.

The interface consists of the following functions:

compile_form - Compilation of forms compile_element - Compilation of finite elements jit - Just-In-Time compilation of forms and elements default_parameters - Default parameter values for FFC ufc_signature - Signature of UFC interface (SHA-1 hash of ufc.h)