api.fem#
Functions#
|
This function collect the arguments used in the assembly function |
|
|
|
|
|
|
|
|
|
|
|
|
|
Classes#
|
Discrete bilinear form ready to be assembled into a matrix. |
|
Discrete functional ready to be assembled into a scalar (real or complex). |
|
Discrete linear form ready to be assembled into a vector. |
|
Class that represents the concept of a discrete sesqui-linear form with the antilinearity on the first variable. |
|
Details#
- collect_spaces(space, *args)[source]#
This function collect the arguments used in the assembly function
- Parameters:
- space: <FunctionSpace>
the symbolic space
- args<list>
list of discrete space components like basis values, spans, …
- Returns:
- args<list>
list of discrete space components elements used in the asembly
- 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:
BasicDiscrete
Discrete 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.
- Parameters:
- exprsympde.expr.expr.BilinearForm
The symbolic bilinear form.
- kernel_exprsympde.expr.evaluation.KernelExpression
The atomic representation of the bilinear form.
- domain_hGeometry
The discretized domain.
- spaceslist of 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.
- matrixStencilMatrix or 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, optional
The symbolic mapping which defines the physical domain of the bilinear form.
- property domain#
- property target#
- property spaces#
- property grid#
- property nquads#
- property test_basis#
- property trial_basis#
- property global_matrices#
- property args#
- 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.
- construct_arguments(with_openmp=False)[source]#
Collect the arguments used in the assembly method.
- Parameters:
- with_openmpbool
If set to True we collect some extra arguments used in the assembly method
- Returns:
- args: tuple
The arguments passed to the assembly method.
- threads_args: tuple
Extra arguments used in the assembly method in case with_openmp=True.
- 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.
- class DiscreteFunctional(expr, kernel_expr, domain_h, space, *, nquads, backend=None, symbolic_mapping=None)[source]#
Bases:
BasicDiscrete
Discrete functional ready to be assembled into a scalar (real or complex).
This class represents the concept of a discrete functional in Psydac. Instances of this class generate an appropriate functional assembly kernel, and prepare a list of arguments for the kernel.
- Parameters:
- exprsympde.expr.expr.Functional
The symbolic functional form.
- kernel_exprsympde.expr.evaluation.KernelExpression
The atomic representation of the functional form.
- domain_hGeometry
The discretized domain.
- spaceFemSpace
The discrete space.
- nquadslist or tuple of int
The number of quadrature points used in the assembly kernel along each direction.
- update_ghost_regionsbool, default=True
Accumulate the contributions of the neighbouring processes.
- backenddict
The backend used to accelerate the computing kernels. The backend dictionaries are defined in the file psydac/api/settings.py
- symbolic_mappingSympde.topology.Mapping
The symbolic mapping which defines the physical domain of the functional.
- property domain#
- property target#
- property space#
- property grid#
- property nquads#
- property test_basis#
- construct_arguments()[source]#
Collect the arguments used in the assembly method.
- Returns:
- args: tuple
The arguments passed to the assembly method.
- assemble(**kwargs)[source]#
This method assembles the square of the functional expression with the given arguments and then compute the square root of the absolute value of the result.
Examples
>>> n = SemiNorm(1.0j*v, domain, kind='l2') >>> nh = discretize(n, domain_h, Vh , **kwargs) >>> fh = FemField(Vh) >>> fh.coeffs[:] = 1 >>> n_value = nh.assemble(v=fh)
In n_value we have the value of np.sqrt(abs(sum((1.0jv)**2)))
- class DiscreteLinearForm(expr, kernel_expr, domain_h, space, *, nquads, vector=None, update_ghost_regions=True, backend=None, symbolic_mapping=None)[source]#
Bases:
BasicDiscrete
Discrete linear form ready to be assembled into a vector.
This class represents the concept of a discrete linear form in Psydac. Instances of this class generate an appropriate vector assembly kernel, allocate the vector if not provided, and prepare a list of arguments for the kernel.
- Parameters:
- exprsympde.expr.expr.LinearForm
The symbolic linear form.
- kernel_exprsympde.expr.evaluation.KernelExpression
The atomic representation of the linear form.
- domain_hGeometry
The discretized domain.
- spaceFemSpace
The discrete test space.
- nquadslist or tuple of int
The number of quadrature points used in the assembly kernel along each direction.
- vectorStencilVector or BlockVector, optional
The vector that we assemble into. If not provided, a new vector of the appropriate space is created.
- update_ghost_regionsbool, default=False
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
- symbolic_mappingSympde.topology.Mapping, optional
The symbolic mapping which defines the physical domain of the linear form.
- property domain#
- property target#
- property space#
- property grid#
- property nquads#
- property test_basis#
- property global_matrices#
- property args#
- assemble(*, reset=True, **kwargs)[source]#
This method assembles the right-hand side Vector by calling the private method self._func with proper arguments.
In the complex case, this function returns the vector conjugate. This comes from the fact that the problem a(u,v)=b(v) is discretize as A @ conj(U) = B due to the antilinearity of a in the first variable. Thus, to obtain U, the assemble function for the LinearForm return conj(B).
TODO: remove these lines when the dot product is changed for complex in sympde. For now, since the dot product does not do 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.
- construct_arguments(with_openmp=False)[source]#
Collect the arguments used in the assembly method.
- Parameters:
- with_openmpbool
If set to True we collect some extra arguments used in the assembly method
- Returns:
- args: tuple
The arguments passed to the assembly method.
- threads_args: tuple
Extra arguments used in the assembly method in case with_openmp=True.