implementations.problem_classes.generic_spectral module

class GenericSpectralLinear(bases, components, comm=None, Dirichlet_recombination=True, left_preconditioner=True, solver_type='cached_direct', solver_args=None, preconditioner_args=None, useGPU=False, max_cached_factorizations=12, spectral_space=True, real_spectral_coefficients=False, heterogeneous=False, debug=False)[source]

Bases: Problem

Generic class to solve problems of the form M u_t + L u = y, with mass matrix M, linear operator L and some right hand side y using spectral methods. L may contain algebraic conditions, as long as (M + dt L) is invertible.

Note that the __getattr__ method is overloaded to pass requests on to the spectral helper if they are not attributes of this class itself. For instance, you can add a BC by calling self.spectral.add_BC or equivalently self.add_BC.

You can port problems derived from this more or less seamlessly to GPU by using the numerical libraries that are class attributes of the spectral helper. This class will automatically switch the datatype using the setup_GPU class method.

spectral

Spectral helper

Type:

pySDC.helpers.spectral_helper.SpectralHelper

work_counters

Dictionary for counting work

Type:

dict

cached_factorizations

Dictionary of cached matrix factorizations for solving

Type:

dict

L

Linear operator

Type:

sparse matrix

M

Mass matrix

Type:

sparse matrix

diff_mask

Mask for separating differential and algebraic terms

Type:

list

Pl

Left preconditioner

Type:

sparse matrix

Pr

Right preconditioner

Type:

sparse matrix

getOutputFile(fileName)[source]
heterogeneous_setup()[source]
processSolutionForOutput(u)[source]
setUpFieldsIO()[source]

Set up FieldsIO for MPI with the space decomposition of this problem

setup_GPU()[source]

switch to GPU modules

setup_L(LHS)[source]

Setup the left hand side of the linear operator L and store it in self.L.

The argument is meant to be a dictionary with the line you want to write the equation in as the key and the relationship between components as another dictionary. For instance, you can add an algebraic condition capturing a first derivative relationship between u and ux as follows:

` Dx = self.get_differentiation_matrix(axes=(0,)) I = self.get_Id() LHS = {'ux': {'u': Dx, 'ux': -I}} self.setup_L(LHS) `

If you put zero as right hand side for the solver in the line for ux, ux will contain the x-derivative of u afterwards.

Parameters:

LHS (dict) – Dictionary containing the equations.

setup_M(LHS)[source]

Setup mass matrix, see documentation of GenericSpectralLinear.setup_L.

setup_preconditioner(Dirichlet_recombination=True, left_preconditioner=True)[source]

Get left and right preconditioners.

Parameters:
  • Dirichlet_recombination (bool) – Basis conversion for right preconditioner. Useful for Chebychev and Ultraspherical methods. 10/10 would recommend.

  • left_preconditioner (bool) – If True, it will interleave the variables and reverse the Kronecker product

solve_system(rhs, dt, u0=None, *args, skip_itransform=False, **kwargs)[source]

Do an implicit Euler step to solve M u_t + Lu = rhs, with M the mass matrix and L the linear operator as setup by GenericSpectralLinear.setup_L and GenericSpectralLinear.setup_M.

The implicit Euler step is (M - dt L) u = M rhs. Note that M need not be invertible as long as (M + dt*L) is. This means solving with dt=0 to mimic explicit methods does not work for all problems, in particular simple DAEs.

Note that by putting M rhs on the right hand side, this function can only solve algebraic conditions equal to zero. If you want something else, it should be easy to overload this function.

compute_residual_DAE(self, stage='')[source]

Computation of the residual that does not add u_0 - u_m in algebraic equations.

Parameters:

stage (str) – The current stage of the step the level belongs to

compute_residual_DAE_MPI(self, stage=None)[source]

Computation of the residual using the collocation matrix Q

Parameters:

stage (str) – The current stage of the step the level belongs to

get_extrapolated_error_DAE(self, S, **kwargs)[source]

The extrapolation estimate combines values of u and f from multiple steps to extrapolate and compare to the solution obtained by the time marching scheme. This function can be used in EstimateExtrapolationError.

Parameters:

S (pySDC.Step) – The current step

Returns:

None