30 #include "sll_assert.h"
31 #include "sll_memory.h"
32 #include "sll_working_precision.h"
57 sll_real64 :: eta1_min
58 sll_real64 :: eta1_max
59 sll_real64,
pointer :: tmp(:)
60 sll_comp64,
pointer :: rhok(:)
79 sll_int32,
intent(in) :: nc_eta1
80 sll_int32,
intent(out) :: error
81 sll_real64,
intent(in) :: eta1_min
82 sll_real64,
intent(in) :: eta1_max
86 self%nc_eta1 = nc_eta1
87 self%eta1_min = eta1_min
88 self%eta1_max = eta1_max
90 sll_allocate(self%rhok(1:nc_eta1), error)
91 sll_allocate(self%tmp(1:nc_eta1), error)
102 sll_real64,
dimension(:),
intent(out) :: field
103 sll_real64,
dimension(:),
intent(in) :: rhs
105 sll_assert(
size(rhs) == self%nc_eta1 + 1)
107 self%tmp = rhs(1:self%nc_eta1)
116 field(1:self%nc_eta1) = self%tmp/real(self%nc_eta1, f64)
118 field(self%nc_eta1 + 1) = field(1)
125 sll_real64,
dimension(:),
intent(in) :: rho
126 sll_real64,
dimension(:),
intent(out) :: phi
128 stop
' compute phi from rho is not implemented '
135 sll_real64,
dimension(:),
intent(in) :: rho
136 sll_real64,
dimension(:),
intent(out) :: e
146 sll_int32,
intent(inout) :: error
148 deallocate (self%rhok, stat=error)
149 deallocate (self%tmp, stat=error)
Fortran module where set some physical and mathematical constants.
real(kind=f64), parameter, public sll_p_pi
complex(kind=f64), parameter, public sll_p_i1
complex(kind=f64), parameter, public sll_p_i0
subroutine, public sll_s_fft_exec_r2c_1d(plan, array_in, array_out)
Compute fast Fourier transform in real to complex mode.
subroutine, public sll_s_fft_free(plan)
Delete a plan.
subroutine, public sll_s_fft_init_r2c_1d(plan, nx, array_in, array_out, normalized, aligned, optimization)
Create new 1d real to complex plan for forward FFT.
subroutine, public sll_s_fft_exec_c2r_1d(plan, array_in, array_out)
Compute fast Fourier transform in complex to real mode.
subroutine, public sll_s_fft_init_c2r_1d(plan, nx, array_in, array_out, normalized, aligned, optimization)
Create new 1d complex to real plan for backward FFT.
Module interface to solve Poisson equation in 1D.
Module to solve Poisson equation for the HMF model.
subroutine, public sll_s_solve_poisson_1d_hmf(self, field, rhs)
Solve the 1d equation for the Vlasov-HMF model.
subroutine, public sll_s_poisson_1d_hmf_init(self, eta1_min, eta1_max, nc_eta1, error)
Initialize the poisson 1d hmf solver.
subroutine, public sll_s_free_poisson_1d_hmf(self, error)
Free the poisson 1d hmf solver.
subroutine compute_e_from_rho_1d_hmf(poisson, E, rho)
solves
subroutine compute_phi_from_rho_1d_hmf(poisson, phi, rho)
solves
Type for FFT plan in SeLaLib.
Implementation of the poisson 1d solver for the Vlasov-HMF model.