9 #include "sll_memory.h"
10 #include "sll_working_precision.h"
27 type:: sll_t_fft_derivative
28 sll_real64,
dimension(:),
pointer :: kx
29 sll_real64,
dimension(:),
pointer :: ky
30 sll_real64,
dimension(:),
pointer :: kpot
31 sll_real64,
dimension(:),
pointer :: kex
32 sll_real64,
dimension(:),
pointer :: key
33 sll_comp64,
dimension(:),
pointer :: fcoeffs, fcoeffs2
46 type(sll_t_sparse_grid_2d),
intent(in) ::interpolator
47 sll_int32 :: ierr, i, j
48 sll_real64,
dimension(:),
allocatable :: data1d
49 sll_real64 :: size_factor
51 sll_allocate(this%kx(interpolator%size_basis), ierr)
52 sll_allocate(this%ky(interpolator%size_basis), ierr)
53 sll_allocate(this%kpot(interpolator%size_basis), ierr)
54 sll_allocate(this%kex(interpolator%size_basis), ierr)
55 sll_allocate(this%key(interpolator%size_basis), ierr)
57 sll_allocate(this%fcoeffs(interpolator%size_basis), ierr)
58 sll_allocate(this%fcoeffs2(interpolator%size_basis), ierr)
60 sll_allocate(data1d(2**(interpolator%max_level)), ierr)
62 size_factor = 2.0_f64*sll_p_pi/interpolator%length(1)
63 do i = 0, interpolator%levels(2)
64 do j = interpolator%index(0, i), &
65 interpolator%index(0, i) + max(2**(i - 1), 1) - 1
67 min(interpolator%max_level - i, interpolator%levels(1)), &
68 j, size_factor, data1d, this%kx)
72 size_factor = 2.0_f64*sll_p_pi/interpolator%length(2)
73 do i = 0, interpolator%levels(1)
74 do j = interpolator%index(i, 0), &
75 interpolator%index(i, 0) + max(2**(i - 1), 1) - 1
77 interpolator%max_level - i), &
78 j, size_factor, data1d, this%ky)
82 this%kpot(1) = 0.0_f64;
83 this%kex(1) = 0.0_f64;
84 this%key(1) = 0.0_f64;
85 do i = 2, interpolator%size_basis
86 this%kpot(i) = 1.0_f64/(this%kx(i)**2 + this%ky(i)**2);
87 this%kex(i) = this%kpot(i)*this%kx(i);
88 this%key(i) = this%kpot(i)*this%ky(i);
96 type(sll_t_sparse_grid_2d),
intent(inout) :: interpolator
97 sll_real64,
dimension(:),
intent(inout) ::phi
98 sll_real64,
dimension(:),
intent(inout) ::rho
101 call interpolator%SPFFT(rho, this%fcoeffs)
103 do i = 1, interpolator%size_basis
104 this%fcoeffs(i) = this%fcoeffs(i)*real(this%kpot(i), f64)
108 call interpolator%ISPFFT(this%fcoeffs, phi)
115 type(sll_t_sparse_grid_2d),
intent(inout) :: interpolator
116 sll_real64,
dimension(:),
intent(inout) :: ex
117 sll_real64,
dimension(:),
intent(inout) :: ey
118 sll_real64,
dimension(:),
intent(inout) :: rho
122 call interpolator%SPFFT(rho, this%fcoeffs)
124 do i = 1, interpolator%size_basis
125 this%fcoeffs2(i) = cmplx(-aimag(this%fcoeffs(i))*this%key(i), &
126 real(this%fcoeffs(i))*this%key(i), kind=f64)
127 this%fcoeffs(i) = cmplx(-aimag(this%fcoeffs(i))*this%kex(i), &
128 real(this%fcoeffs(i))*this%kex(i), kind=f64)
131 call interpolator%ISPFFT(this%fcoeffs, ex)
132 call interpolator%ISPFFT(this%fcoeffs2, ey)
138 class(sll_t_sparse_grid_2d),
intent(in) :: interpolator
139 sll_real64,
dimension(:),
intent(inout) ::
data, data1d
140 sll_int32,
intent(in) :: dim, max_level, index
141 sll_real64,
intent(in) :: size_factor
148 data1d(size) = real(size/2, f64)*size_factor
150 data1d(2*k) = real(k, f64)*size_factor;
151 data1d(2*k + 1) = -data1d(2*k);
160 type(sll_t_sparse_grid_2d),
intent(in) :: sparsegrid
161 sll_int32,
intent(in) :: dim, max_level, index
162 sll_int32 :: n_points, index_running
163 sll_real64,
dimension(:),
intent(in) :: data_in
164 sll_real64,
dimension(:),
intent(out) :: data_out
166 n_points = 2**(max_level)
167 data_out(index) = data_in(1)
168 if (max_level > 0)
then
169 index_running = sparsegrid%hierarchy(index)%children(dim*2)
170 data_out(index_running) = data_in(2)
171 if (max_level > 1)
then
173 2, max_level, dim, data_in, data_out)
180 type(sll_t_sparse_grid_2d),
intent(in) :: sparsegrid
181 sll_int32,
intent(in) :: level, max_level, index_sg, dim, ind
182 sll_real64,
dimension(:),
intent(in) :: data_in
183 sll_real64,
dimension(:),
intent(inout) :: data_out
185 data_out(sparsegrid%hierarchy(index_sg)%children(dim*2 - 1)) = &
186 data_in(2**(level - 1) + 1 + 2*ind)
187 data_out(sparsegrid%hierarchy(index_sg)%children(dim*2)) = &
188 data_in(2**(level - 1) + 1 + 2*ind + 1)
189 if (level < max_level)
then
191 sparsegrid%hierarchy(index_sg)%children(dim*2 - 1), 2*ind, &
192 level + 1, max_level, dim, data_in, data_out)
194 sparsegrid%hierarchy(index_sg)%children(dim*2), 2*ind + 1, &
195 level + 1, max_level, dim, data_in, data_out)
Fortran module where set some physical and mathematical constants.
real(kind=f64), parameter, public sll_p_pi
Implementation of a 3D pseudospectral Poisson solver on sparse grid.
subroutine solve_potential(this, interpolator, rho, phi)
Solve for potential.
subroutine solve_for_electric_field(this, interpolator, rho, ex, ey)
Compute the electric fields from rho.
subroutine derivative_coeffs_1d(interpolator, dim, max_level, index, size_factor, data1d, data)
Helper function to compute the Fourier coefficients for the derivative operation.
recursive subroutine insert_recursive_fourier_real(sparsegrid, index_sg, ind, level, max_level, dim, data_in, data_out)
Helper function to insert the real Fourier coefficients into the sparse grid data structure (recursiv...
subroutine new_poisson_2d_sparse_grid_fft(this, interpolator)
Create Poisson solver object with Fourier spectral method on 2d sparse grid.
subroutine insert_fourier_real(sparsegrid, dim, max_level, index, data_in, data_out)
Helper function to insert the real Fourier coefficient into the sparse grid data structure.
Implementation of a 2D sparse grid with interpolation routines.
sll_t_fft_derivative is the Poisson solver object to solve Poisson's problem in 2d with pseudospectra...
Sparse grid object for 2d with interpolation routines.