Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_control_variate.F90
Go to the documentation of this file.
1 
5 
7 
8 #include "sll_working_precision.h"
9 #include "sll_memory.h"
10 
11  use sll_m_initial_distribution, only : &
13 
14  implicit none
15  private
16 
17  public :: sll_i_control_variate, &
20 
21 
24 
25  sll_real64, pointer :: control_variate_parameters(:) => null()
26  class(sll_c_distribution_params), pointer :: control_variate_distribution_params => null()
27 
28  procedure(sll_i_control_variate), pointer :: control_variate => null()
29 
30  contains
31 
32  procedure :: update_df_weight
33  procedure :: init => init_control_variate
34  procedure :: free => free_control_variate
35 
36  end type sll_t_control_variate
37 
38 
41 
42  class(sll_t_control_variate), allocatable :: cv(:)
43 
44  end type sll_t_control_variates
45 
46  abstract interface
47 
48  function sll_i_control_variate(self, xi, vi, time)
51  class(sll_t_control_variate) :: self
52  sll_real64, optional, intent( in ) :: xi(:)
53  sll_real64, optional, intent( in ) :: vi(:)
54  sll_real64, optional, intent( in ) :: time
55  sll_real64 :: sll_i_control_variate
56 
57  end function sll_i_control_variate
58  end interface
59 
60 
61  contains
62 
64  function update_df_weight(self, xi, vi, time, weight_ff, g0) result(weight_df)
65  class(sll_t_control_variate) :: self
66  sll_real64, intent( in ) :: xi(:)
67  sll_real64, intent( in ) :: vi(:)
68  sll_real64, intent( in ) :: time
69  sll_real64, intent( in ) :: weight_ff
70  sll_real64, intent( in ) :: g0
71  sll_real64 :: weight_df
72 
73  weight_df = weight_ff - self%control_variate(xi, vi, time)/g0
74 
75  end function update_df_weight
76 
78  subroutine init_control_variate( self, control_function, parameters, distribution_params )
79  class(sll_t_control_variate), intent(out) :: self
80  procedure(sll_i_control_variate) :: control_function
81  sll_real64, target,optional, intent(in) :: parameters(:)
82  class(sll_c_distribution_params), target, optional, intent(in) :: distribution_params
83 
84  self%control_variate => control_function
85 
86  if (present(parameters)) self%control_variate_parameters => parameters
87 
88  if (present(distribution_params)) self%control_variate_distribution_params => distribution_params
89 
90 
91  end subroutine init_control_variate
92 
93 
95  subroutine free_control_variate( self )
96  class(sll_t_control_variate), intent(inout) :: self
97 
98  self%control_variate_parameters => null()
99  self%control_variate_distribution_params => null()
100 
101  end subroutine free_control_variate
102 
103 
104 
105 end module sll_m_control_variate
1d real function, abstract interface for function defining the control variate
real(kind=f64) function update_df_weight(self, xi, vi, time, weight_ff, g0)
Update the delta f weights.
subroutine free_control_variate(self)
Destructor.
subroutine init_control_variate(self, control_function, parameters, distribution_params)
Initialization.
Parameters to define common initial distributions.
Module to select the kind parameter.
real(kind=f64) function, dimension(size(particle, 2)) control_variate(particle)
Abstract data type for parameters of initial distribution.
    Report Typos and Errors