Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_preconditioner_singular.F90
Go to the documentation of this file.
1 
9 
11  !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
12 #include "sll_assert.h"
13 #include "sll_memory.h"
14 #include "sll_working_precision.h"
15 
16  use sll_m_linear_solver_block, only : &
18 
19  use sll_m_linear_solver_abstract, only : &
21 
22  implicit none
23 
24  public :: &
26 
27  private
28  !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
29 
31  type(sll_t_linear_solver_block), pointer :: inverse_mass_fft
32  sll_int32 :: n_total
33  sll_real64, allocatable :: lumped_mass(:)
34 
35 
36  contains
37 
38  procedure :: create => create_3d_trafo
39  procedure :: free => free_3d_trafo
40  procedure :: solve_real => solve_3d_trafo
41  procedure :: set_verbose
42  procedure :: print_info
43  procedure :: read_from_file
44 
46 
47 contains
48 
49 
50  subroutine create_3d_trafo( self, inverse_mass_fft, lumped_mass, n_total)
51  class(sll_t_preconditioner_singular), intent( inout ) :: self
52  type(sll_t_linear_solver_block), target :: inverse_mass_fft
53  sll_real64 :: lumped_mass(:)
54  sll_int32, intent( in ) :: n_total
55 
56  self%inverse_mass_fft => inverse_mass_fft
57  self%n_total = n_total
58 
59  allocate( self%lumped_mass(n_total) )
60  self%lumped_mass = lumped_mass
61 
62  end subroutine create_3d_trafo
63 
64  subroutine solve_3d_trafo(self, rhs, unknown)
65  class(sll_t_preconditioner_singular), intent( inout ) :: self
66  sll_real64, intent(in ) :: rhs(:)
67  sll_real64, intent( out) :: unknown(:)
68  !local variables
69  sll_real64 :: scratch(self%n_total)
70 
71  scratch = rhs*self%lumped_mass
72  call self%inverse_mass_fft%solve( scratch, unknown )
73  unknown = unknown*self%lumped_mass
74 
75  end subroutine solve_3d_trafo
76 
77  subroutine set_verbose(self, verbose)
78  class(sll_t_preconditioner_singular), intent( inout ) :: self
79  logical, intent( in ) :: verbose
80 
81  self%verbose = verbose
82  end subroutine set_verbose
83 
84  subroutine print_info(self)
85  class(sll_t_preconditioner_singular), intent( in ) :: self
86 
87  end subroutine print_info
88 
89  subroutine read_from_file(self, filename)
90  class(sll_t_preconditioner_singular), intent( inout ) :: self
91  character(len=*), intent( in ) :: filename
92 
93  end subroutine read_from_file
94 
95 
96  subroutine free_3d_trafo( self )
97  class(sll_t_preconditioner_singular), intent(inout) :: self
98 
99  self%inverse_mass_fft => null()
100 
101  end subroutine free_3d_trafo
102 
103 
module for abstract linear solver
module for a block linear solver
Module interface to solve Maxwell's equations with coordinate transformation in 3D The linear systems...
subroutine solve_3d_trafo(self, rhs, unknown)
subroutine create_3d_trafo(self, inverse_mass_fft, lumped_mass, n_total)
    Report Typos and Errors