Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_binomial_filter.F90
Go to the documentation of this file.
1 
5 
7 
8  !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
9 #include "sll_working_precision.h"
10 
11  use sll_m_filter_base_1d, only: &
13 
14  implicit none
15 
16  public :: &
18 
19  private
20 
21 
22  !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
23 
25 
26  sll_real64, allocatable :: scratch(:)
27 
28  logical :: odd
29  sll_int32 :: double_iters
30 
31  contains
32 
33  procedure :: init => init_binomial
34  procedure :: apply => apply_binomial
35  procedure :: apply_inplace => apply_inplace_binomial
36 
37  end type sll_t_binomial_filter
38 
39 contains
40 
41  subroutine init_binomial( self, iterations, n_dofs, mode )
42  class( sll_t_binomial_filter ), intent( out ) :: self
43  sll_int32, intent( in ) :: iterations
44  sll_int32, intent( in ) :: n_dofs
45  sll_int32, optional, intent( in ) :: mode
46 
47  self%iterations = iterations
48  self%n_dofs = n_dofs
49  allocate( self%scratch(n_dofs) )
50 
51  if (mod(iterations, 2) .eq. 1 ) then
52  self%odd = .true.
53  self%double_iters = iterations/2
54  else
55  self%odd = .false.
56  self%double_iters = iterations/2
57  end if
58 
59  !print*, 'Binomial filter', self%double_iters, self%odd
60 
61  end subroutine init_binomial
62 
63  subroutine apply_binomial( self, field_in, field_out )
64  class( sll_t_binomial_filter ), intent( inout ) :: self
65  sll_real64, intent(in) :: field_in(:)
66  sll_real64, intent(out) :: field_out(:)
67 
68  sll_int32 :: iter
69 
70  if ( self%odd .eqv. .true. ) then
71  call sll_s_binomial_filter( field_in, field_out, self%n_dofs )
72  else
73  field_out = field_in
74  end if
75 
76 
77  do iter=1,self%double_iters
78  call sll_s_binomial_filter( field_out, self%scratch, self%n_dofs )
79  call sll_s_binomial_filter( self%scratch, field_out, self%n_dofs )
80  end do
81  end subroutine apply_binomial
82 
83 
84  subroutine apply_inplace_binomial( self, field )
85  class( sll_t_binomial_filter ), intent( inout ) :: self
86  sll_real64, intent(inout) :: field(:)
87 
88  sll_int32 :: iter
89 
90  do iter=1,self%double_iters
91  call sll_s_binomial_filter( field, self%scratch, self%n_dofs )
92  call sll_s_binomial_filter( self%scratch, field, self%n_dofs )
93  end do
94 
95  if ( self%odd .eqv. .true. ) then
96  call sll_s_binomial_filter( field, self%scratch, self%n_dofs )
97  field = self%scratch
98  end if
99 
100  end subroutine apply_inplace_binomial
101 
102  subroutine sll_s_binomial_filter( field_in, field_out, n_dofs )
103  sll_real64, intent(in) :: field_in(:)
104  sll_real64, intent(out) :: field_out(:)
105  sll_int32, intent(in) :: n_dofs
106 
107  sll_int32 :: j
108 
109  ! Filter
110  field_out(1) = 0.25_f64*( field_in(1)*2.0_f64 + &
111  field_in(n_dofs)+field_in(2))
112  do j=2,n_dofs-1
113  field_out(j) = 0.25_f64*( field_in(j)*2.0_f64 + &
114  field_in(j-1)+field_in(j+1))
115  end do
116  field_out(n_dofs) = 0.25_f64*( field_in(n_dofs)*2.0_f64 + &
117  field_in(n_dofs-1)+field_in(1))
118 
119  end subroutine sll_s_binomial_filter
120 
121 
122 
123 
124 
125 end module sll_m_binomial_filter
Binomial filter for smooting of fields.
subroutine sll_s_binomial_filter(field_in, field_out, n_dofs)
subroutine init_binomial(self, iterations, n_dofs, mode)
subroutine apply_binomial(self, field_in, field_out)
subroutine apply_inplace_binomial(self, field)
    Report Typos and Errors