Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_preconditioner_poisson_fft.F90
Go to the documentation of this file.
1 
7 
9  !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
10 #include "sll_assert.h"
11 #include "sll_memory.h"
12 #include "sll_working_precision.h"
13 
14  use sll_m_linear_solver_abstract, only : &
16 
17  use sll_m_poisson_3d_fem_fft, only : &
19 
20  implicit none
21 
22  public :: &
24 
25  private
26  !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
27 
29  type(sll_t_poisson_3d_fem_fft) :: poisson_solver
30  sll_real64, allocatable :: jacobi_precond(:)
31 
32  sll_int32 :: n_dofs
33 
34  contains
35 
36  procedure :: create => create_poisson_fft
37  procedure :: free => free_poisson_fft
38  procedure :: solve_real => solve_poisson_fft
39  procedure :: set_verbose => set_verbose_poisson_fft
40  procedure :: print_info => print_info_poisson_fft
41  procedure :: read_from_file => read_from_file_poisson_fft
42 
44 
45 contains
46 
47  subroutine create_poisson_fft( self, n_dofs, degree, delta_x, jacobi_in )
48  class(sll_t_preconditioner_poisson_fft), intent( inout ) :: self
49  sll_int32, intent( in ) :: n_dofs(3)
50  sll_int32, intent( in ) :: degree(3)
51  sll_real64, intent( in ) :: delta_x(3)
52  sll_real64, intent( in ) :: jacobi_in(:)
53 
54  call self%poisson_solver%init( n_dofs, degree, delta_x )
55 
56  allocate( self%jacobi_precond( product(n_dofs)) )
57  self%jacobi_precond = jacobi_in
58  self%n_dofs = product(n_dofs)
59 
60  end subroutine create_poisson_fft
61 
62  subroutine free_poisson_fft( self )
63  class(sll_t_preconditioner_poisson_fft), intent( inout ) :: self
64 
65  call self%poisson_solver%free()
66 
67  end subroutine free_poisson_fft
68 
69  subroutine solve_poisson_fft(self, rhs, unknown )
70  class(sll_t_preconditioner_poisson_fft), intent( inout ) :: self
71  sll_real64, intent(in ) :: rhs(:)
72  sll_real64, intent( out) :: unknown(:)
73 
74  sll_real64, allocatable :: tmp(:)
75 
76  allocate( tmp(self%n_dofs) )
77  unknown = self%jacobi_precond*rhs
78  ! tmp = unknown
79  call self%poisson_solver%compute_phi_from_rho( unknown, tmp )
80  unknown = self%jacobi_precond*tmp
81 
82 
83  end subroutine solve_poisson_fft
84 
85  subroutine read_from_file_poisson_fft(self, filename)
86  class(sll_t_preconditioner_poisson_fft), intent( inout ) :: self
87  character(len=*), intent( in ) :: filename
88 
89  end subroutine read_from_file_poisson_fft
90 
91  subroutine set_verbose_poisson_fft(self, verbose)
92  class(sll_t_preconditioner_poisson_fft), intent( inout ) :: self
93  logical, intent( in ) :: verbose
94 
95  self%verbose = verbose
96  end subroutine set_verbose_poisson_fft
97 
98  subroutine print_info_poisson_fft(self)
99  class(sll_t_preconditioner_poisson_fft), intent( in ) :: self
100 
101  end subroutine print_info_poisson_fft
102 
103 
module for abstract linear solver
This module is a wrapper around the spline FEM Poisson solver for the uniform grid with periodic boun...
subroutine solve_poisson_fft(self, rhs, unknown)
subroutine create_poisson_fft(self, n_dofs, degree, delta_x, jacobi_in)
    Report Typos and Errors