pymor.discretizers.builtin.cg¶
This module provides some operators for continuous finite element discretizations.
Module Contents¶
- class pymor.discretizers.builtin.cg.AdvectionOperatorP1(grid, boundary_info, advection_function=None, advection_constant=None, dirichlet_clear_columns=False, dirichlet_clear_diag=False, solver=None, name=None)[source]¶
Bases:
pymor.operators.numpy.NumpyMatrixBasedOperatorLinear advection
Operatorfor linear finite elements.The operator is of the form
(Lu)(x) = c ∇ ⋅ [ v(x) u(x) ]
The function
vhas to be vector-valued.- Parameters:
grid – The
Gridfor which to assemble the operator.boundary_info –
BoundaryInfofor the treatment of Dirichlet boundary conditions.advection_function – The
Functionv(x)withshape_range = (grid.dim, ). IfNone, constant one is assumed.advection_constant – The constant
c. IfNone,cis set to one.dirichlet_clear_columns – If
True, set columns of the system matrix corresponding to Dirichlet boundary DOFs to zero to obtain a symmetric system matrix. Otherwise, only the rows will be set to zero.dirichlet_clear_diag – If
True, also set diagonal entries corresponding to Dirichlet boundary DOFs to zero. Otherwise they are set to one.solver – The
Solverfor the operator.name – Name of the operator.
- class pymor.discretizers.builtin.cg.AdvectionOperatorQ1(grid, boundary_info, advection_function=None, advection_constant=None, dirichlet_clear_columns=False, dirichlet_clear_diag=False, solver=None, name=None)[source]¶
Bases:
pymor.operators.numpy.NumpyMatrixBasedOperatorLinear advection
Operatorfor bilinear finite elements.The operator is of the form
(Lu)(x) = c ∇ ⋅ [ v(x) u(x) ]
The function
vhas to be vector-valued.- Parameters:
grid – The
Gridfor which to assemble the operator.boundary_info –
BoundaryInfofor the treatment of Dirichlet boundary conditions.advection_function – The
Functionv(x)withshape_range = (grid.dim, ). IfNone, constant one is assumed.advection_constant – The constant
c. IfNone,cis set to one.dirichlet_clear_columns – If
True, set columns of the system matrix corresponding to Dirichlet boundary DOFs to zero to obtain a symmetric system matrix. Otherwise, only the rows will be set to zero.dirichlet_clear_diag – If
True, also set diagonal entries corresponding to Dirichlet boundary DOFs to zero. Otherwise they are set to one.solver – The
Solverfor the operator.name – Name of the operator.
- class pymor.discretizers.builtin.cg.BoundaryDirichletFunctional(grid, dirichlet_data, boundary_info, solver=None, name=None)[source]¶
Bases:
pymor.operators.numpy.NumpyMatrixBasedOperatorLinear functional for enforcing Dirichlet boundary values.
- Parameters:
grid –
Gridfor which to assemble the functional.dirichlet_data –
Functionproviding the Dirichlet boundary values.boundary_info –
BoundaryInfodetermining the Dirichlet boundaries.solver – The
Solverfor the operator.name – The name of the functional.
- class pymor.discretizers.builtin.cg.BoundaryL2ProductFunctional(grid, function, boundary_type=None, dirichlet_clear_dofs=False, boundary_info=None, solver=None, name=None)[source]¶
Bases:
pymor.operators.numpy.NumpyMatrixBasedOperatorLinear functional representing the inner product with an L2-
Functionon the boundary.- Parameters:
grid –
Gridfor which to assemble the functional.function – The
Functionwith which to take the inner product.boundary_type – The type of domain boundary (e.g. ‘neumann’) on which to assemble the functional. If
Nonethe functional is assembled over the whole boundary.dirichlet_clear_dofs – If
True, set dirichlet boundary DOFs to zero.boundary_info – If
boundary_typeis specified ordirichlet_clear_dofsisTrue, theBoundaryInfodetermining which boundary entity belongs to which physical boundary.solver – The
Solverfor the operator.name – The name of the functional.
- class pymor.discretizers.builtin.cg.DiffusionOperatorP1(grid, boundary_info, diffusion_function=None, diffusion_constant=None, dirichlet_clear_columns=False, dirichlet_clear_diag=False, solver=None, name=None)[source]¶
Bases:
pymor.operators.numpy.NumpyMatrixBasedOperatorDiffusion
Operatorfor linear finite elements.The operator is of the form
(Lu)(x) = c ∇ ⋅ [ d(x) ∇ u(x) ]
The function
dcan be scalar- or matrix-valued.- Parameters:
grid – The
Gridfor which to assemble the operator.boundary_info –
BoundaryInfofor the treatment of Dirichlet boundary conditions.diffusion_function – The
Functiond(x)withshape_range == ()orshape_range = (grid.dim, grid.dim). IfNone, constant one is assumed.diffusion_constant – The constant
c. IfNone,cis set to one.dirichlet_clear_columns – If
True, set columns of the system matrix corresponding to Dirichlet boundary DOFs to zero to obtain a symmetric system matrix. Otherwise, only the rows will be set to zero.dirichlet_clear_diag – If
True, also set diagonal entries corresponding to Dirichlet boundary DOFs to zero. Otherwise they are set to one.solver – The
Solverfor the operator.name – Name of the operator.
- class pymor.discretizers.builtin.cg.DiffusionOperatorQ1(grid, boundary_info, diffusion_function=None, diffusion_constant=None, dirichlet_clear_columns=False, dirichlet_clear_diag=False, solver=None, name=None)[source]¶
Bases:
pymor.operators.numpy.NumpyMatrixBasedOperatorDiffusion
Operatorfor bilinear finite elements.The operator is of the form
(Lu)(x) = c ∇ ⋅ [ d(x) ∇ u(x) ]
The function
dcan be scalar- or matrix-valued.- Parameters:
grid – The
Gridfor which to assemble the operator.boundary_info –
BoundaryInfofor the treatment of Dirichlet boundary conditions.diffusion_function – The
Functiond(x)withshape_range == ()orshape_range = (grid.dim, grid.dim). IfNone, constant one is assumed.diffusion_constant – The constant
c. IfNone,cis set to one.dirichlet_clear_columns – If
True, set columns of the system matrix corresponding to Dirichlet boundary DOFs to zero to obtain a symmetric system matrix. Otherwise, only the rows will be set to zero.dirichlet_clear_diag – If
True, also set diagonal entries corresponding to Dirichlet boundary DOFs to zero. Otherwise they are set to one.solver – The
Solverfor the operator.name – Name of the operator.
- class pymor.discretizers.builtin.cg.InterpolationOperator(grid, function, solver=None, name=None)[source]¶
Bases:
pymor.operators.numpy.NumpyMatrixBasedOperatorVector-like Lagrange interpolation
Operatorfor continuous finite element spaces.- Parameters:
- class pymor.discretizers.builtin.cg.L2ProductFunctionalP1(grid, function, dirichlet_clear_dofs=False, boundary_info=None, solver=None, name=None)[source]¶
Bases:
pymor.operators.numpy.NumpyMatrixBasedOperatorLinear functional representing the inner product with an L2-
Function.- Parameters:
grid –
Gridfor which to assemble the functional.function – The
Functionwith which to take the inner product.dirichlet_clear_dofs – If
True, set dirichlet boundary DOFs to zero.boundary_info –
BoundaryInfodetermining the Dirichlet boundaries in casedirichlet_clear_dofsis set toTrue.solver – The
Solverfor the operator.name – The name of the functional.
- class pymor.discretizers.builtin.cg.L2ProductFunctionalQ1(grid, function, dirichlet_clear_dofs=False, boundary_info=None, solver=None, name=None)[source]¶
Bases:
pymor.operators.numpy.NumpyMatrixBasedOperatorBilinear functional representing the inner product with an L2-
Function.- Parameters:
grid –
Gridfor which to assemble the functional.function – The
Functionwith which to take the inner product.dirichlet_clear_dofs – If
True, set dirichlet boundary DOFs to zero.boundary_info –
BoundaryInfodetermining the Dirichlet boundaries in casedirichlet_clear_dofsis set toTrue.solver – The
Solverfor the operator.name – The name of the functional.
- class pymor.discretizers.builtin.cg.L2ProductP1(grid, boundary_info, dirichlet_clear_rows=True, dirichlet_clear_columns=False, dirichlet_clear_diag=False, coefficient_function=None, solver=None, name=None)[source]¶
Bases:
pymor.operators.numpy.NumpyMatrixBasedOperatorOperatorrepresenting the L2-product between linear finite element functions.- Parameters:
grid – The
Gridfor which to assemble the product.boundary_info –
BoundaryInfofor the treatment of Dirichlet boundary conditions.dirichlet_clear_rows – If
True, set the rows of the system matrix corresponding to Dirichlet boundary DOFs to zero.dirichlet_clear_columns – If
True, set columns of the system matrix corresponding to Dirichlet boundary DOFs to zero.dirichlet_clear_diag – If
True, also set diagonal entries corresponding to Dirichlet boundary DOFs to zero. Otherwise, if eitherdirichlet_clear_rowsordirichlet_clear_columnsisTrue, the diagonal entries are set to one.coefficient_function – Coefficient
Functionfor product withshape_range == (). IfNone, constant one is assumed.solver – The
Solverfor the operator.name – The name of the product.
- class pymor.discretizers.builtin.cg.L2ProductQ1(grid, boundary_info, dirichlet_clear_rows=True, dirichlet_clear_columns=False, dirichlet_clear_diag=False, coefficient_function=None, solver=None, name=None)[source]¶
Bases:
pymor.operators.numpy.NumpyMatrixBasedOperatorOperatorrepresenting the L2-product between bilinear finite element functions.- Parameters:
grid – The
Gridfor which to assemble the product.boundary_info –
BoundaryInfofor the treatment of Dirichlet boundary conditions.dirichlet_clear_rows – If
True, set the rows of the system matrix corresponding to Dirichlet boundary DOFs to zero.dirichlet_clear_columns – If
True, set columns of the system matrix corresponding to Dirichlet boundary DOFs to zero.dirichlet_clear_diag – If
True, also set diagonal entries corresponding to Dirichlet boundary DOFs to zero. Otherwise, if eitherdirichlet_clear_rowsordirichlet_clear_columnsisTrue, the diagonal entries are set to one.coefficient_function – Coefficient
Functionfor product withshape_range == (). IfNone, constant one is assumed.solver – The
Solverfor the operator.name – The name of the product.
- class pymor.discretizers.builtin.cg.RobinBoundaryOperator(grid, boundary_info, robin_data=None, solver=None, name=None)[source]¶
Bases:
pymor.operators.numpy.NumpyMatrixBasedOperatorRobin boundary
Operatorfor linear finite elements.The operator represents the contribution of Robin boundary conditions to the stiffness matrix, where the boundary condition is supposed to be given in the form
-[ d(x) ∇u(x) ] ⋅ n(x) = c(x) (u(x) - g(x))
dandnare the diffusion function (seeDiffusionOperatorP1) and the unit outer normal inx, whilecis the (scalar) Robin parameter function andgis the (also scalar) Robin boundary value function.- Parameters:
grid – The
Gridover which to assemble the operator.boundary_info –
BoundaryInfofor the treatment of Dirichlet boundary conditions.robin_data – Tuple (c(x), g(x)) providing two
Functionsthat represent the Robin parameter and boundary value function. IfNone, the resulting operator is zero.solver – The
Solverfor the operator.name – Name of the operator.
- pymor.discretizers.builtin.cg.discretize_instationary_cg(analytical_problem, diameter=None, domain_discretizer=None, grid_type=None, grid=None, boundary_info=None, num_values=None, time_stepper=None, nt=None, preassemble=True, solver=None)[source]¶
Finite Element discretization of an
InstationaryProblem.Discretizes an
InstationaryProblemwith aStationaryProblemas the stationary part using finite elements.- Parameters:
analytical_problem – The
InstationaryProblemto discretize.diameter – If not
None,diameteris passed as an argument to thedomain_discretizer.domain_discretizer – Discretizer to be used for discretizing the analytical domain. This has to be a function
domain_discretizer(domain_description, diameter, ...). IfNone,discretize_domain_defaultis used.grid_type – If not
None, this parameter is forwarded todomain_discretizerto specify the type of the generatedGrid.grid – Instead of using a domain discretizer, the
Gridcan also be passed directly using this parameter.boundary_info – A
BoundaryInfospecifying the boundary types of the grid boundary entities. Must be provided ifgridis specified.num_values – The number of returned vectors of the solution trajectory. If
None, each intermediate vector that is calculated is returned.time_stepper – The
time-stepperto be used bysolveofInstationaryModel.nt – If
time_stepperis not specified, the number of time steps for implicit Euler time stepping.preassemble – If
True, preassemble all operators in the resultingModel.solver – The
Solverto be used.
- Returns:
m – The
Modelthat has been generated.data – Dictionary with the following entries:
- grid:
The generated
Grid.- boundary_info:
The generated
BoundaryInfo.- unassembled_m:
In case
preassembleisTrue, the generatedModelbefore preassembling operators.
- pymor.discretizers.builtin.cg.discretize_stationary_cg(analytical_problem, diameter=None, domain_discretizer=None, grid_type=None, grid=None, boundary_info=None, preassemble=True, mu_energy_product=None, solver=None)[source]¶
Discretizes a
StationaryProblemusing finite elements.- Parameters:
analytical_problem – The
StationaryProblemto discretize.diameter – If not
None,diameteris passed as an argument to thedomain_discretizer.domain_discretizer – Discretizer to be used for discretizing the analytical domain. This has to be a function
domain_discretizer(domain_description, diameter, ...). IfNone,discretize_domain_defaultis used.grid_type – If not
None, this parameter is forwarded todomain_discretizerto specify the type of the generatedGrid.grid – Instead of using a domain discretizer, the
Gridcan also be passed directly using this parameter.boundary_info – A
BoundaryInfospecifying the boundary types of the grid boundary entities. Must be provided ifgridis specified.preassemble – If
True, preassemble all operators in the resultingModel.mu_energy_product – If not
None,parameter valuesfor which to assemble the symmetric part of theOperatorof the resultingModelfom(ignoring the advection part). Thus, assuming no advection and a symmetric diffusion tensor,fom.products['energy']is equal tofom.operator.assemble(mu), except for the fact that the former has cleared Dirichlet rows and columns, while the latter only has cleared Dirichlet rows).solver – The
Solverto be used.
- Returns:
m – The
Modelthat has been generated.data – Dictionary with the following entries:
- grid:
The generated
Grid.- boundary_info:
The generated
BoundaryInfo.- unassembled_m:
In case
preassembleisTrue, the generatedModelbefore preassembling operators.