pynitride.core.fem module

pynitride.core.fem.assemble_load_matrix(w, dzn, n, dirichelet1=False, dirichelet2=False)

Assemble a load matrix from the coefficients of the differential equation.

Cythonized assembly method, call as assemble_load_matrix(w, dzn, n).

For the mathematics/definitions of the \(w\) term, see Finite Element Method. All arguments should be defined on the mid mesh, and the \(w\) matrix should be one-dimensional floats

Parameters:
  • w – load-side coefficient

  • dzn – the spacings between mesh nodes

  • n – the size of the \(C\) matrices (eg 1 for scalar equation, 6 for 6x6 kp equation)

  • dirichelet1 (bool) – whether to treat boundary 1, 2 as Dirichelet (True, default) or Neumann.

  • dirichelet2 (bool) – whether to treat boundary 1, 2 as Dirichelet (True, default) or Neumann.

Returns:

A load matrix in sparse CSC form

pynitride.core.fem.assemble_stiffness_matrix(C0, Cl, Cr, C2, dzn, dirichelet1=True, dirichelet2=True)

Assemble a stiffness matrix from the coefficients of the differential equation.

Cythonized assembly method, call as assemble_stiffness_matrix(C0, Cl, Cr, C2, dzn, dirichelet1, dirichelet2).

For the mathematics/definitions of the \(C\) terms, see Finite Element Method. All arguments should be defined on the mid mesh except \(U\), and all the \(C\) matrices should be three-dimensional, even if they are just scalars.

If Cl and Cr are supplied as None, they will be ignored and the matrix will be of dtype double. If Cl and Cr are supplied, the matrix will be of dtype complex.

Parameters:
  • C0 – material-dependent differential equation coefficients

  • Cl – material-dependent differential equation coefficients

  • Cr – material-dependent differential equation coefficients

  • C2 – material-dependent differential equation coefficients

  • dzn – the spacings between mesh nodes

  • dirichelet1 (bool) – whether to treat boundary 1, 2 as Dirichelet (True, default) or Neumann.

  • dirichelet2 (bool) – whether to treat boundary 1, 2 as Dirichelet (True, default) or Neumann.

Returns:

A stiffness matrix in sparse CSC form

pynitride.core.fem.fem_eigsh()

Solve the eigenvalue problem.

For the mathematics/definitions of terms, see Finite Element Method.

Parameters:
  • stiffness_matrix – a stiffness matrix \(A\) from pynitride.fem.assemble_stiffness_matrix()

  • load_matrix – a load matrix \(M\) from pynitride.fem.assemble_load_matrix()

  • eigval_out – an array into which to fill the eigenvalues, should be shape (number of eigenvalues)

  • eigvec_out – an array into which to fill the eigenvectors, should be shape (number of eigenvalues x n x len(zn))

  • n – dimension of the differential equation

  • dirichelet1 (bool) – whether to treat boundary 1, 2 as Dirichelet (True, default) or Neumann.

  • dirichelet2 (bool) – whether to treat boundary 1, 2 as Dirichelet (True, default) or Neumann.

  • pairwise_GS – in the case of near degeneracy, the eigensolver often returns eigenvectors which are not totally orthogonal. For instance in a 6x6 wurtzite k.p problem with small CR coupling term, the states form pairs that are very close in energy. This option breaks the eigenvectors into those pairs [note, it just assumes the pairing 0-1, 2-3, 4-5, etc, it doesn’t actually check the energies], and in each pair, performs a single Gram-Schmidt to ensure the nearly-degenerate eigenvectors are orthogonal.

  • ascending (bool) – by default, energies are sorted from min to max, but for valence band problems, you want the highest absolute energies first, so set ascending=False

  • *args – passed onto scipy.sparse.eigsh along with the main arguments A (stiffness) and M (load)

  • **kwargs

    passed onto scipy.sparse.eigsh along with the main arguments A (stiffness) and M (load)

Returns:

None

pynitride.core.fem.fem_get_error()

Finds the signed error vector of the test solution.

The error is \(b-M^{-1}Ax\), where \(x\) is the supplied test solution. For the mathematics/definitions of terms, see Finite Element Method.

Parameters:
  • stiffness_matrix – a stiffness matrix \(A\) from pynitride.fem.assemble_stiffness_matrix()

  • load_matrix – a load matrix \(M\) from pynitride.fem.assemble_load_matrix()

  • load_vec – the load vector \(b\), should be shape (len(zn))

  • test – the test solution vector \(x\), should be shape (len(zn))

  • err_out – an array into which to fill the error, should be shape (len(zn))

  • n – dimension of the differential equation

  • dirichelet1 (bool) – whether to treat boundary 1, 2 as Dirichelet (True, default) or Neumann.

  • dirichelet2 (bool) – whether to treat boundary 1, 2 as Dirichelet (True, default) or Neumann.

  • *args – passed onto scipy.sparse.spsolve along with the main arguments M and A@x

  • **kwargs

    passed onto scipy.sparse.spsolve along with the main arguments M and A@x

Returns:

the error (same shape as test)

pynitride.core.fem.fem_solve()

Solve the linear matrix equation

For the mathematics/definitions of terms, see Finite Element Method.

Parameters:
  • stiffness_matrix – a stiffness matrix \(A\) from pynitride.fem.assemble_stiffness_matrix()

  • load_matrix – a load matrix \(M\) from pynitride.fem.assemble_load_matrix()

  • load_vec – the load vector \(b\), should be shape (len(zn))

  • val_out – an array into which to fill the solution, should be shape (len(zn))

  • n – dimension of the differential equation

  • dirichelet1 (bool) – whether to treat boundary 1, 2 as Dirichelet (True, default) or Neumann.

  • dirichelet2 (bool) – whether to treat boundary 1, 2 as Dirichelet (True, default) or Neumann.

  • *args

    passed onto scipy.sparse.spsolve along with the main arguments A (stiffness) and M@b (load)

  • **kwargs

    passed onto scipy.sparse.spsolve along with the main arguments A (stiffness) and M@b (load)

Returns:

the solution (same shape as load_vec)