PNSolver

Standard spherical harmonics moment solver.

Source:

  • include/solvers/pnsolver.hpp

  • src/solvers/pnsolver.cpp

class PNSolver : public SolverBase

Subclassed by CSDPNSolver

Public Functions

PNSolver(Config *settings)

PNSolver constructor.

Parameters:

settings – stores all needed information

inline virtual ~PNSolver()

PNSolver destructor.

Protected Functions

virtual void PrepareVolumeOutput() override

Initializes the output groups and fields of this solver and names the fields.

virtual void WriteVolumeOutput(unsigned idx_iter) override

Function that prepares VTK export and csv export of the current solver iteration.

Parameters:

idx_iter – current (pseudo) time iteration

virtual void FVMUpdate(unsigned idx_iter) override

Computes the finite Volume update step for the current iteration.

Parameters:

idx_iter – current (peudo) time iteration

virtual void FluxUpdate() override

Constructs the flux update for the current iteration and stores it in psiNew.

virtual void IterPreprocessing(unsigned idx_iter) override

Performs preprocessing for the current solver iteration.

Parameters:

idx_iter – current (peudo) time iteration

virtual void ComputeScalarFlux() override

Computes the flux of the solution to check conservation properties.

double AParam(int l, int k) const

parameter functions for setting up system matrix

Parameters:
  • l – degree , it must hold: 0 <= l <=_nq

  • k – order , it must hold: -l <=k <= l

double BParam(int l, int k) const

parameter functions for setting up system matrix

Parameters:
  • l – degree , it must hold: 0 <= l <=_nq

  • k – order , it must hold: -l <=k <= l

double CParam(int l, int k) const

parameter functions for setting up system matrix

Parameters:
  • l – degree , it must hold: 0 <= l <=_nq

  • k – order , it must hold: -l <=k <= l

double DParam(int l, int k) const

parameter functions for setting up system matrix

Parameters:
  • l – degree , it must hold: 0 <= l <=_nq

  • k – order , it must hold: -l <=k <= l

double EParam(int l, int k) const

parameter functions for setting up system matrix

Parameters:
  • l – degree , it must hold: 0 <= l <=_nq

  • k – order , it must hold: -l <=k <= l

double FParam(int l, int k) const

parameter functions for setting up system matrix

Parameters:
  • l – degree , it must hold: 0 <= l <=_nq

  • k – order , it must hold: -l <=k <= l

double CTilde(int l, int k) const

parameter functions for setting up system matrix

Parameters:
  • l – degree , it must hold: 0 <= l <=_nq

  • k – order , it must hold: -l <=k <= l

double DTilde(int l, int k) const

parameter functions for setting up system matrix

Parameters:
  • l – degree , it must hold: 0 <= l <=_nq

  • k – order , it must hold: -l <=k <= l

double ETilde(int l, int k) const

parameter functions for setting up system matrix

Parameters:
  • l – degree , it must hold: 0 <= l <=_nq

  • k – order , it must hold: -l <=k <= l

double FTilde(int l, int k) const

parameter functions for setting up system matrix

Parameters:
  • l – degree , it must hold: 0 <= l <=_nq

  • k – order , it must hold: -l <=k <= l

int Sgn(int k) const

mathematical + index functions. Helper functions for setting up system matrix.

Parameters:

k – arbitrary integer

int kPlus(int k) const

mathematical + index functions. Helper functions for setting up system matrix.

Parameters:

k – arbitrary integer

int kMinus(int k) const

mathematical + index functions. Helper functions for setting up system matrix.

Parameters:

k – arbitrary integer

int GlobalIndex(int l, int k) const

computes the global index of the moment corresponding to basis function (l,k)

Parameters:
  • l – degree l, it must hold: 0 <= l <=_nq

  • k – order k, it must hold: -l <=k <= l

Returns:

global index

bool CheckIndex(int l, int k) const

Checks, if index invariant for global index holds.

Returns:

True, if invariant holds, else false

void ComputeSystemMatrices()

Function for computing and setting up system matrices.

void ComputeFluxComponents()

Function for computing and setting up flux matrices for upwinding.

void ComputeScatterMatrix()

Function for computing and setting up EV matrix for scattering kernel.

double LegendrePoly(double x, int l)

Computes Legedre polinomial of oder l at point x.

void CleanFluxMatrices()

Sets Entries of FluxMatrices to zero, if they are below double precision, to prevent floating point inaccuracies later in the solver

void FluxUpdatePseudo1D()

Flux update version for Pseudo1D.

void FluxUpdatePseudo2D()

Flux update version for Pseudo2D

Protected Attributes

unsigned _nSystem

total number of equations in the system

unsigned _polyDegreeBasis

maximal degree of the spherical harmonics basis

SymMatrix _Ax

Flux Jacbioan in x direction.

SymMatrix _Ay

Flux Jacbioan in x direction.

SymMatrix _Az

Flux Jacbioan in x direction.

Matrix _AxPlus

Flux Jacbioan in x direction, positive part.

Matrix _AxMinus

Flux Jacbioan in x direction, negative part.

Matrix _AxAbs

Flux Jacbioan in x direction, absolute part.

Matrix _AyPlus

Flux Jacbioan in y direction, positive part.

Matrix _AyMinus

Flux Jacbioan in y direction, negative part.

Matrix _AyAbs

Flux Jacbioan in y direction, absolute part.

Matrix _AzPlus

Flux Jacbioan in z direction, positive part.

Matrix _AzMinus

Flux Jacbioan in z direction, negative part.

Matrix _AzAbs

Flux Jacbioan in z direction, absolute part.

Vector _scatterMatDiag

diagonal of the scattering matrix (its a diagonal matrix by construction). Contains eigenvalues of the scattering kernel.