projection module#

class projection.phaseShift(domain: dolfinx.mesh.Mesh, degree_psi=1)[source]#

Bases: object

This class implements the phase shift of the wave functions \([\psi_1, \psi_2]\).

Parameters:
  • domain (dolfinx.mesh.Mesh) – Mesh imported.

  • degree_psi (int (Default = 1)) – Degree of the polynomials used for the definition of the functional space.

assemble(hbar: float, direct: bool = False)[source]#

This function assembles the forms of the weak formulations and creates the linear structures for the solution (matrix, rhs and solver).

Parameters:
  • hbar (float) – Value of \(\hbar\).

  • direct (boolean (Default: False)) – Boolean variable to identify if a direct solver must be used.

solve(psiTilde_1: dolfinx.fem.Function, psiTilde_2: dolfinx.fem.Function, pressure: dolfinx.fem.Function)[source]#

This function computes the shifted wave function.

Parameters:
  • psiTilde_1 (dolfinx.fem.Function) – First component of the wave function \(\tilde{\psi}_1\).

  • psiTilde_2 (dolfinx.fem.Function) – Second component of the wave function \(\tilde{\psi}_2\).

  • pressure (dolfinx.fem.Function) – Pressure computed by the projection step.

class projection.poisson(domain: dolfinx.mesh.Mesh, ft: dolfinx.cpp.mesh.MeshTags_int32, boundary_marks: dict)[source]#

Bases: object

This class implements the projection step of the Incompressible Schrodinger Flow, in particular the solution of the Poisson problem.

Parameters:
  • domain (dolfinx.mesh.Mesh) – Mesh imported.

  • ft (dolfinx.cpp.mesh.MeshTags_int32) – Face tags of the boundaries.

  • boundary_marks (dict) – Dictionary with the markers for the boundaries.

assemble(bc_values: dict, gauss: bool = False, direct: bool = False)[source]#

This function assign the boundary condition (if inlet is present), assembles the forms of the weak formulations and creates the linear structures for the solution (matrix, rhs and solver).

The weak formulation is as follows:

\[\int_\Omega \nabla \varphi \cdot \nabla q\, d\Omega = - \int_\Omega q\,\nabla \cdot \tilde{\mathbf{u}}\, d\Omega\]

there is an option to apply the integration by parts on the rhs as

\[\int_\Omega \nabla \varphi \cdot \nabla q\, d\Omega = \int_\Omega \nabla q\cdot \tilde{\mathbf{u}}\, d\Omega - \int_{\partial \Omega} \tilde{\mathbf{u}}\cdot \mathbf{n}\,q\,d\sigma\]
Parameters:
  • bc_values (dict) – Dictionary with the boundary conditions parameters.

  • gauss (boolean (Default: False)) – If True, the integration by parts is applied also to the right-hand side.

  • direct (boolean (Default: False)) – Boolean variable to identify if a direct solver must be used.

solve(uTilde: dolfinx.fem.Function)[source]#

This function solves the Poisson problem.

Parameters:

uTilde (dolfinx.fem.Function) – Velocity field coming from the prediction step.