api.fem_bilinear_form#
Classes#

|
Discrete bilinear form ready to be assembled into a matrix. |
Details#
- class DiscreteBilinearForm(expr, kernel_expr, domain_h, spaces, *, nquads, matrix=None, update_ghost_regions=True, backend=None, linalg_backend=None, assembly_backend=None, symbolic_mapping=None)[source]#
Bases:
objectDiscrete bilinear form ready to be assembled into a matrix.
This class represents the concept of a discrete bilinear form in PSYDAC. Instances of this class generate an appropriate matrix assembly kernel, allocate the matrix if not provided, and prepare a list of arguments for the kernel.
An implementation of the sum factorization algorithm is used to assemble the matrix.
- Parameters:
- exprsympde.expr.expr.BilinearForm
The symbolic bilinear form.
- kernel_exprlist or tuple of sympde.expr.evaluation.KernelExpression
The atomic representation of the bilinear form.
- domain_hpsydac.cad.geometry.Geometry
The discretized domain.
- spaceslist of psydac.fem.basic.FemSpace
The discrete trial and test spaces.
- nquadslist or tuple of int
The number of quadrature points used in the assembly kernel along each direction.
- matrixpsydac.linalg.stencil.StencilMatrix or psydac.linalg.block.BlockLinearOperator, optional
The matrix that we assemble into. If not provided, a new matrix is created with the appropriate domain and codomain (default: None).
- update_ghost_regionsbool, default=True
Accumulate the contributions of the neighbouring processes.
- backenddict, optional
The backend used to accelerate the computing kernels. The backend dictionaries are defined in the file psydac/api/settings.py
- assembly_backenddict, optional
The backend used to accelerate the assembly kernel. The backend dictionaries are defined in the file psydac/api/settings.py
- linalg_backenddict, optional
The backend used to accelerate the computing kernels of the linear operator. The backend dictionaries are defined in the file psydac/api/settings.py
- symbolic_mappingsympde.topology.mapping.Mapping, optional
The symbolic mapping which defines the physical domain of the bilinear form.
See also
DiscreteLinearFormDiscreteFunctionalDiscreteSumForm
- property comm#
- property expr#
- property kernel_expr#
- property domain#
- property mapping#
- property is_rational_mapping#
- property target#
- property spaces#
- property test_basis#
- property trial_basis#
- property grid#
- property nquads#
- property free_args#
- property max_nderiv#
- property backend#
- property args#
- property global_matrices#
- allocate_matrices(backend=None)[source]#
Allocate the global matrices used in the assembly method. In this method we allocate only the matrices that are computed in the self._target domain, we also avoid double allocation if we have many DiscreteLinearForm that are defined on the same self._target domain.
- Parameters:
- backenddict
The backend used to accelerate the computing kernels.
- assemble(*, reset=True, **kwargs)[source]#
This method assembles the left hand side Matrix by calling the private method self._func with proper arguments.
In the complex case, this function returns the matrix conjugate. This comes from the fact that the problem a(u,v)=b(v) is discretized as A @ conj(U) = B due to the antilinearity of a in the first variable. Thus, to obtain U, the assemble function returns conj(A).
TODO: remove these lines when the dot product is changed for complex. For now, since the dot product does not compute the conjugate in the complex case. We do not use the conjugate in the assemble function. It should work if the complex only comes from the rhs in the linear form.
- make_file(temps, ordered_stmts, field_derivatives, max_logical_derivative, test_mult, trial_mult, test_v_p, trial_u_p, keys_1, keys_2, keys_3, mapping_option)[source]#
Part of the sum factorization algorithm implementation. Generates the correct assembly file. Used at the end of construct_arguments_generate_assembly_file, before eventually pyccelizing that file.
- Parameters:
- tempstuple
Tuple of Assign statements defining temporary values. Arithmetic combinations of these make up the coupling terms.
- ordered_stmtsdict
Dictionary defining the coupling terms. Keys are combinations of test and trial function components, values are Assign statements in terms of temporaries appearing in temps.
- field_derivativesdict
Dictionary containing information on the derivatives of free FemFields. Keys are components of free FemFields. Values are dictionaries again. Their keys are names, as appearing in the assembly file, of partial derivatives of the corresponding FemField component, and their values are dictionaries again. Example: {F1_0_x3 : {‘x1’: 0, ‘x2’: 0, ‘x3’: 1}, F1_0_x2 : …} Meaning: There exists a free FemField named F1. Among other, the partial derivative w.r.t. x3 of its first component F1_0 appears.
- max_logical_derivativeint
The largest appearing derivative order.
- test_multlist
List of length 3(scalar test function) or 9(vector test function) including multiplicity information.
- trial_multlist
List of length 3(scalar trial function) or 9(vector trial function) including multiplicity information.
- test_v_pdict
Dictionary of length 1(scalar test function) or length 3(vector test function). Each key corresponds to a component of the funciton (space), and each corresponding value is a list of Bspline degrees of length 3. Example: Discretizing a de de Rham sequence using a degree vector [2, 3, 4] means that test_v_p for a test function belonging to H(curl) will be {0: [1, 3, 4], 1: [2, 2, 4], 2: [2, 3, 3]}
- trial_u_pdict
Dictionary of length 1(scalar trial function) or length 3(vector trial function). Each key corresponds to a component of the funciton (space), and each corresponding value is a list of Bspline degrees of length 3. Example: Discretizing a de de Rham sequence using a degree vector [2, 3, 4] means that trial_u_p for a trial function belonging to H^1 will be {0: [2, 3, 4]}
- keys_1dict
Dictionary relating subexpressions to x1-derivative combinations. Keys are combinations of test and trial function components. Values are lists, each entry corresponding to one appearing partial derivative combination of these components. Example: keys_1[(u[0], v[1])][3] = [1,0] means that the fourth ([3]) sub-expression (partial derivative combination) corresponding to the trial-test-function-component-product u[0] * v[1] involves a first derivative in x1 direction of the trial function and no derivative in x1 direction of the test function. Information on appearing partial derivatives in x2 and x3 direction is stored in keys_2 and keys_3.
- keys_2dict
See keys_1.
- keys_3dict
See keys_1.
- mapping_optionNone | ‘Bspline’
None in case of no mapping or an analytical mapping, ‘Bspline’ in case of a Bspline mapping.
- Returns:
- file_idstr
random string of length 8, corresponding to the assembly file name located in __psydac__/
- read_BilinearForm()[source]#
Part of the sum factorization algorithm implementation. Used at the beginning of construct_arguments_generate_assembly_file(). It’s output determines both the design of the assembly function, and the arguments passed to it.
- Returns:
- tempstuple
tuple of Assign objects. Often times usable building blocks of complicated coupling terms.
- ordered_stmtsdict
assigns each block (trial&test component combination) a list of coupling term assignment
- ordered_sub_exprs_keysdict
relates each coupling term assignment of ordered_stmts a partial derivative combination
- mapping_optionstr | None
‘Bspline’ if a spline mapping is involved, None if an analytical or no mapping is involved
- field_derivativesdict
contains information regarding appearing free FemFields and appearing partial derivatives of those
- g_mat_information_falselist
possibly wrong list of non-zero blocks
- g_mat_information_truelist
correct list of non-zero blocks
- max_logical_derivativeint
maximum appearing partial derivative (in any fixed direction)
- construct_arguments_generate_assembly_file()[source]#
Collect the arguments used in the assembly method, and generate and possibly pyccelize the assembly function.
Used only when sum factorization is enabled, else the method construct_arguments is called.
- Returns:
- args: tuple
The arguments passed to the assembly method.
- threads_args: None
None as openMP parallelization is not supported by this implementation.
- property test_trial_template#