Octopus
|
Data Types | |
type | poisson_t |
Functions/Subroutines | |
subroutine, public | poisson_init (this, namespace, space, der, mc, stencil, qtot, label, solver, verbose, force_serial, force_cmplx) |
subroutine, public | poisson_end (this) |
subroutine | zpoisson_solve_real_and_imag_separately (this, namespace, pot, rho, all_nodes, kernel) |
subroutine, public | zpoisson_solve (this, namespace, pot, rho, all_nodes, kernel) |
subroutine, public | poisson_solve_batch (this, namespace, potb, rhob, all_nodes, kernel) |
subroutine, public | dpoisson_solve (this, namespace, pot, rho, all_nodes, kernel) |
Calculates the Poisson equation. Given the density returns the corresponding potential. More... | |
subroutine, public | poisson_init_sm (this, namespace, space, main, der, sm, method, force_cmplx) |
subroutine, public | poisson_test (this, space, mesh, latt, namespace, repetitions) |
This routine checks the Hartree solver selected in the input file by calculating numerically and analytically the Hartree potential originated by a Gaussian distribution of charge. For periodic systems, the periodic copies of the Gaussian are taken into account up to to a certain threshold that can be specified in the input file. More... | |
logical pure function, public | poisson_solver_is_iterative (this) |
logical pure function, public | poisson_is_multigrid (this) |
logical pure function, public | poisson_solver_has_free_bc (this) |
integer pure function, public | poisson_get_solver (this) |
subroutine, public | poisson_async_init (this, mc) |
subroutine, public | poisson_async_end (this, mc) |
subroutine, public | poisson_slave_work (this, namespace) |
logical pure function, public | poisson_is_async (this) |
subroutine, public | poisson_build_kernel (this, namespace, space, coulb, qq, alpha, beta, mu, singul) |
real(real64) function, public | poisson_get_full_range_weight (this, alpha, beta, mu) |
Most Poisson solvers do not implement Coulomb attenuated potentials, and can only be used for global hybrids like PBE0. For these ones, this routine returns the weight of the full-range hybrids. More... | |
subroutine | poisson_kernel_init (this, namespace, space, mc, stencil) |
subroutine | poisson_solve_direct (this, namespace, pot, rho) |
subroutine | dpoisson_solve_direct_sm (this, namespace, sm, pot, rho) |
subroutine | zpoisson_solve_direct_sm (this, namespace, sm, pot, rho) |
subroutine, public | dpoisson_solve_start (this, rho) |
subroutine, public | dpoisson_solve_finish (this, pot) |
subroutine, public | dpoisson_solve_sm (this, namespace, sm, pot, rho, all_nodes) |
Calculates the Poisson equation. Given the density returns the corresponding potential. More... | |
subroutine, public | zpoisson_solve_start (this, rho) |
subroutine, public | zpoisson_solve_finish (this, pot) |
subroutine, public | zpoisson_solve_sm (this, namespace, sm, pot, rho, all_nodes) |
Calculates the Poisson equation. Given the density returns the corresponding potential. More... | |
Variables | |
integer, parameter, public | poisson_direct_sum = -1 |
integer, parameter, public | poisson_fft = 0 |
integer, parameter, public | poisson_cg = 5 |
integer, parameter, public | poisson_cg_corrected = 6 |
integer, parameter, public | poisson_multigrid = 7 |
integer, parameter, public | poisson_isf = 8 |
integer, parameter, public | poisson_psolver = 10 |
integer, parameter, public | poisson_no = -99 |
integer, parameter, public | poisson_null = -999 |
integer, parameter | cmd_finish = 1 |
integer, parameter | cmd_poisson_solve = 2 |
subroutine, public poisson_oct_m::poisson_init | ( | type(poisson_t), intent(inout) | this, |
type(namespace_t), intent(in) | namespace, | ||
class(space_t), intent(in) | space, | ||
type(derivatives_t), intent(in), target | der, | ||
type(multicomm_t), intent(in) | mc, | ||
type(stencil_t), intent(in) | stencil, | ||
real(real64), intent(in), optional | qtot, | ||
character(len=*), intent(in), optional | label, | ||
integer, intent(in), optional | solver, | ||
logical, intent(in), optional | verbose, | ||
logical, intent(in), optional | force_serial, | ||
logical, intent(in), optional | force_cmplx | ||
) |
[in] | qtot | total charge |
Definition at line 237 of file poisson.F90.
subroutine, public poisson_oct_m::poisson_end | ( | type(poisson_t), intent(inout) | this | ) |
Definition at line 713 of file poisson.F90.
|
private |
[in,out] | pot | pot(meshnp) |
[in] | rho | rho(meshnp) |
Definition at line 766 of file poisson.F90.
subroutine, public poisson_oct_m::zpoisson_solve | ( | type(poisson_t), intent(in) | this, |
type(namespace_t), intent(in) | namespace, | ||
complex(real64), dimension(:), intent(inout), contiguous | pot, | ||
complex(real64), dimension(:), intent(in), contiguous | rho, | ||
logical, intent(in), optional | all_nodes, | ||
type(fourier_space_op_t), intent(in), optional | kernel | ||
) |
[in,out] | pot | pot(meshnp) |
[in] | rho | rho(meshnp) |
Definition at line 815 of file poisson.F90.
subroutine, public poisson_oct_m::poisson_solve_batch | ( | type(poisson_t), intent(inout) | this, |
type(namespace_t), intent(in) | namespace, | ||
type(batch_t), intent(inout) | potb, | ||
type(batch_t), intent(inout) | rhob, | ||
logical, intent(in), optional | all_nodes, | ||
type(fourier_space_op_t), intent(in), optional | kernel | ||
) |
Definition at line 856 of file poisson.F90.
subroutine, public poisson_oct_m::dpoisson_solve | ( | type(poisson_t), intent(in) | this, |
type(namespace_t), intent(in) | namespace, | ||
real(real64), dimension(:), intent(inout), contiguous | pot, | ||
real(real64), dimension(:), intent(in), contiguous | rho, | ||
logical, intent(in), optional | all_nodes, | ||
type(fourier_space_op_t), intent(in), optional | kernel | ||
) |
Calculates the Poisson equation. Given the density returns the corresponding potential.
Different solvers are available that can be chosen in the input file with the "PoissonSolver" parameter
[in,out] | pot | Local size of the potential vector. |
[in] | rho | Local size of the density (rho) vector. |
[in] | all_nodes | Is the Poisson solver allowed to utilise all nodes or only the domain nodes for its calculations? (Defaults to .true.) |
Definition at line 891 of file poisson.F90.
subroutine, public poisson_oct_m::poisson_init_sm | ( | type(poisson_t), intent(inout) | this, |
type(namespace_t), intent(in) | namespace, | ||
class(space_t), intent(in) | space, | ||
type(poisson_t), intent(in) | main, | ||
type(derivatives_t), intent(in), target | der, | ||
type(submesh_t), intent(inout) | sm, | ||
integer, intent(in), optional | method, | ||
logical, intent(in), optional | force_cmplx | ||
) |
Definition at line 993 of file poisson.F90.
subroutine, public poisson_oct_m::poisson_test | ( | type(poisson_t), intent(in) | this, |
class(space_t), intent(in) | space, | ||
class(mesh_t), intent(in) | mesh, | ||
type(lattice_vectors_t), intent(in) | latt, | ||
type(namespace_t), intent(in) | namespace, | ||
integer, intent(in) | repetitions | ||
) |
This routine checks the Hartree solver selected in the input file by calculating numerically and analytically the Hartree potential originated by a Gaussian distribution of charge. For periodic systems, the periodic copies of the Gaussian are taken into account up to to a certain threshold that can be specified in the input file.
Definition at line 1103 of file poisson.F90.
logical pure function, public poisson_oct_m::poisson_solver_is_iterative | ( | type(poisson_t), intent(in) | this | ) |
Definition at line 1299 of file poisson.F90.
logical pure function, public poisson_oct_m::poisson_is_multigrid | ( | type(poisson_t), intent(in) | this | ) |
Definition at line 1307 of file poisson.F90.
logical pure function, public poisson_oct_m::poisson_solver_has_free_bc | ( | type(poisson_t), intent(in) | this | ) |
Definition at line 1316 of file poisson.F90.
integer pure function, public poisson_oct_m::poisson_get_solver | ( | type(poisson_t), intent(in) | this | ) |
Definition at line 1330 of file poisson.F90.
subroutine, public poisson_oct_m::poisson_async_init | ( | type(poisson_t), intent(inout) | this, |
type(multicomm_t), intent(in) | mc | ||
) |
Definition at line 1338 of file poisson.F90.
subroutine, public poisson_oct_m::poisson_async_end | ( | type(poisson_t), intent(inout) | this, |
type(multicomm_t), intent(in) | mc | ||
) |
Definition at line 1363 of file poisson.F90.
subroutine, public poisson_oct_m::poisson_slave_work | ( | type(poisson_t), intent(inout) | this, |
type(namespace_t), intent(in) | namespace | ||
) |
Definition at line 1390 of file poisson.F90.
logical pure function, public poisson_oct_m::poisson_is_async | ( | type(poisson_t), intent(in) | this | ) |
Definition at line 1442 of file poisson.F90.
subroutine, public poisson_oct_m::poisson_build_kernel | ( | type(poisson_t), intent(in) | this, |
type(namespace_t), intent(in) | namespace, | ||
class(space_t), intent(in) | space, | ||
type(fourier_space_op_t), intent(inout) | coulb, | ||
real(real64), dimension(:), intent(in) | qq, | ||
real(real64), intent(in) | alpha, | ||
real(real64), intent(in) | beta, | ||
real(real64), intent(in) | mu, | ||
real(real64), intent(in), optional | singul | ||
) |
Definition at line 1451 of file poisson.F90.
real(real64) function, public poisson_oct_m::poisson_get_full_range_weight | ( | type(poisson_t), intent(in) | this, |
real(real64), intent(in) | alpha, | ||
real(real64), intent(in) | beta, | ||
real(real64), intent(in) | mu | ||
) |
Most Poisson solvers do not implement Coulomb attenuated potentials, and can only be used for global hybrids like PBE0. For these ones, this routine returns the weight of the full-range hybrids.
If the Poisson solver supports CAM hybrids, the returned weight is 1.
For purely short-range, or purely long-range, the corresponding weight is also returned, such that Poisson solvers that support these cases (like ISF) can also be used correctly.
Definition at line 1524 of file poisson.F90.
|
private |
Definition at line 1564 of file poisson.F90.
|
private |
Definition at line 1768 of file poisson.F90.
|
private |
Definition at line 2017 of file poisson.F90.
|
private |
[in,out] | pot | pot(meshnp) |
[in] | rho | rho(meshnp) |
Definition at line 2199 of file poisson.F90.
subroutine, public poisson_oct_m::dpoisson_solve_start | ( | type(poisson_t), intent(in) | this, |
real(real64), dimension(:), intent(in) | rho | ||
) |
Definition at line 2310 of file poisson.F90.
subroutine, public poisson_oct_m::dpoisson_solve_finish | ( | type(poisson_t), intent(in) | this, |
real(real64), dimension(:), intent(inout) | pot | ||
) |
Definition at line 2333 of file poisson.F90.
subroutine, public poisson_oct_m::dpoisson_solve_sm | ( | type(poisson_t), intent(in) | this, |
type(namespace_t), intent(in) | namespace, | ||
type(submesh_t), intent(in) | sm, | ||
real(real64), dimension(:), intent(inout), contiguous | pot, | ||
real(real64), dimension(:), intent(inout), contiguous | rho, | ||
logical, intent(in), optional | all_nodes | ||
) |
Calculates the Poisson equation. Given the density returns the corresponding potential.
Different solvers are available that can be chosen in the input file with the "PoissonSolver" parameter
[in,out] | pot | Local size of the potential vector. |
[in,out] | rho | Local size of the density (rho) vector. |
[in] | all_nodes | Is the Poisson solver allowed to utilise all nodes or only the domain nodes for its calculations? (Defaults to .true.) |
Definition at line 2356 of file poisson.F90.
subroutine, public poisson_oct_m::zpoisson_solve_start | ( | type(poisson_t), intent(in) | this, |
complex(real64), dimension(:), intent(in) | rho | ||
) |
Definition at line 2495 of file poisson.F90.
subroutine, public poisson_oct_m::zpoisson_solve_finish | ( | type(poisson_t), intent(in) | this, |
complex(real64), dimension(:), intent(inout) | pot | ||
) |
Definition at line 2518 of file poisson.F90.
subroutine, public poisson_oct_m::zpoisson_solve_sm | ( | type(poisson_t), intent(in) | this, |
type(namespace_t), intent(in) | namespace, | ||
type(submesh_t), intent(in) | sm, | ||
complex(real64), dimension(:), intent(inout), contiguous | pot, | ||
complex(real64), dimension(:), intent(inout), contiguous | rho, | ||
logical, intent(in), optional | all_nodes | ||
) |
Calculates the Poisson equation. Given the density returns the corresponding potential.
Different solvers are available that can be chosen in the input file with the "PoissonSolver" parameter
[in,out] | pot | Local size of the potential vector. |
[in,out] | rho | Local size of the density (rho) vector. |
[in] | all_nodes | Is the Poisson solver allowed to utilise all nodes or only the domain nodes for its calculations? (Defaults to .true.) |
Definition at line 2541 of file poisson.F90.
integer, parameter, public poisson_oct_m::poisson_direct_sum = -1 |
Definition at line 194 of file poisson.F90.
integer, parameter, public poisson_oct_m::poisson_fft = 0 |
Definition at line 194 of file poisson.F90.
integer, parameter, public poisson_oct_m::poisson_cg = 5 |
Definition at line 194 of file poisson.F90.
integer, parameter, public poisson_oct_m::poisson_cg_corrected = 6 |
Definition at line 194 of file poisson.F90.
integer, parameter, public poisson_oct_m::poisson_multigrid = 7 |
Definition at line 194 of file poisson.F90.
integer, parameter, public poisson_oct_m::poisson_isf = 8 |
Definition at line 194 of file poisson.F90.
integer, parameter, public poisson_oct_m::poisson_psolver = 10 |
Definition at line 194 of file poisson.F90.
integer, parameter, public poisson_oct_m::poisson_no = -99 |
Definition at line 194 of file poisson.F90.
integer, parameter, public poisson_oct_m::poisson_null = -999 |
Definition at line 194 of file poisson.F90.
|
private |
Definition at line 230 of file poisson.F90.
|
private |
Definition at line 230 of file poisson.F90.