Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_maxwell_3d_base.F90
Go to the documentation of this file.
1 
6 
8  !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
9 #include "sll_working_precision.h"
10 
11  use sll_m_utilities, only: &
13 
14  use sll_m_profile_functions, only: &
16 
17  implicit none
18 
19  public :: &
22  ! sll_s_plot_two_fields_3d
23 
24  private
25  !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
26 
27  type, abstract :: sll_c_maxwell_3d_base
28 
29  sll_real64 :: delta_x(3)
30  sll_real64 :: volume
31  sll_int32 :: n_cells(3)
32  sll_int32 :: n_dofs(3)
33  sll_int32 :: n_total
34  sll_int32 :: n_total0
35  sll_int32 :: n_total1
36  sll_int32 :: s_deg_0(3)
37  sll_int32 :: s_deg_1(3)
38  sll_real64 :: lx(3)
39  sll_real64 :: mass_solver_tolerance = 1d-12
40  sll_real64 :: poisson_solver_tolerance = 1d-12
41  sll_real64 :: solver_tolerance = 1d-12
42  type(sll_t_profile_functions) :: profile
43 
44  contains
45 
46  procedure(compute_field1_from_field2), deferred :: &
47  compute_e_from_b
48  procedure(compute_field1_from_field2), deferred :: &
49  compute_b_from_e
50  procedure(compute_field_from_field), deferred :: &
51  compute_e_from_rho
52  procedure(compute_field_from_field), deferred :: &
53  compute_rho_from_e
54  procedure(compute_e_from_j_3d), deferred :: &
55  compute_e_from_j
56  procedure(compute_phi_e_from_field), deferred :: &
57  compute_phi_from_rho
58  procedure(compute_phi_e_from_field), deferred :: &
59  compute_phi_from_j
60  procedure(update_dofs_function), deferred :: &
61  compute_rhs_from_function
62  procedure(update_dofs_function), deferred :: &
63  l2projection
64  procedure(norm_squared), deferred :: &
65  l2norm_squared
66  procedure(inner_product), deferred :: &
68  procedure(empty), deferred :: &
69  free
70 
71  procedure(compute_field_from_field), deferred :: &
72  multiply_c
73  procedure(compute_field_from_field), deferred :: &
74  multiply_ct
75  procedure(compute_field_from_field), deferred :: &
76  multiply_g
77  procedure(compute_field_from_field), deferred :: &
78  multiply_gt
79  procedure(multiply_mass), deferred :: &
81  procedure(multiply_mass_inverse), deferred :: &
83  procedure :: compute_curl_part
84 
85  end type sll_c_maxwell_3d_base
86 
87 
88 
89  !---------------------------------------------------------------------------!
90  abstract interface
91  subroutine compute_field1_from_field2( self, delta_t, field_in, field_out )
93  import sll_c_maxwell_3d_base
94  class(sll_c_maxwell_3d_base) :: self
95  sll_real64, intent( in ) :: delta_t
96  sll_real64, intent( in ) :: field_in(:)
97  sll_real64, intent( inout ) :: field_out(:)
98 
99  end subroutine compute_field1_from_field2
100  end interface
101 
102 
103 
104  !---------------------------------------------------------------------------!
105  abstract interface
106  subroutine compute_field_from_field( self, field_in, field_out )
108  import sll_c_maxwell_3d_base
109  class(sll_c_maxwell_3d_base) :: self
110  sll_real64, intent( in ) :: field_in(:)
111  sll_real64, intent( out ) :: field_out(:)
112 
113  end subroutine compute_field_from_field
114  end interface
115 
116  !---------------------------------------------------------------------------!
117  abstract interface
118  subroutine compute_phi_e_from_field( self, field_in, field_out, efield_dofs )
120  import sll_c_maxwell_3d_base
121  class(sll_c_maxwell_3d_base) :: self
122  sll_real64, intent( in ) :: field_in(:)
123  sll_real64, intent( inout ) :: field_out(:)
124  sll_real64, intent( out ) :: efield_dofs(:)
125 
126  end subroutine compute_phi_e_from_field
127  end interface
128 
129  !---------------------------------------------------------------------------!
130  abstract interface
131  subroutine compute_e_from_j_3d( self, current, E, component )
133  import sll_c_maxwell_3d_base
134  class(sll_c_maxwell_3d_base) :: self
135  sll_real64, intent( in ) :: current(:)
136  sll_real64, intent( inout ) :: e(:)
137  sll_int32, intent( in ), optional :: component
138 
139  end subroutine compute_e_from_j_3d
140  end interface
141 
142 
143 
144  !---------------------------------------------------------------------------!
145  abstract interface
146 
148  use sll_m_working_precision ! can't pass a header file because the
149  ! preprocessor prevents double inclusion.
150  ! It is very rare.
151  sll_real64 :: sll_i_function_3d_real64
152  sll_real64, intent(in) :: x(3)
153  end function sll_i_function_3d_real64
154  end interface
155 
156  !---------------------------------------------------------------------------!
157  abstract interface
158  subroutine update_dofs_function( self, form, component, coefs_dofs, func1, func2, func3 )
160  import sll_c_maxwell_3d_base
162  class( sll_c_maxwell_3d_base) :: self
163  sll_int32, intent( in ) :: form
164  sll_int32, intent( in ) :: component
165  sll_real64, intent( out ) :: coefs_dofs(:)
166  procedure(sll_i_function_3d_real64) :: func1
167  procedure(sll_i_function_3d_real64), optional :: func2
168  procedure(sll_i_function_3d_real64), optional :: func3
169 
170  end subroutine update_dofs_function
171  end interface
172 
173 
174 
175  !---------------------------------------------------------------------------!
176  abstract interface
177  function norm_squared( self, coefs, form, component ) result( r )
179  import sll_c_maxwell_3d_base
180  class( sll_c_maxwell_3d_base) :: self
181  sll_real64 :: coefs(:)
182  sll_int32 :: form
183  sll_int32 :: component
184  sll_real64 :: r
185  end function norm_squared
186  end interface
187 
188 
189 
190  !---------------------------------------------------------------------------!
191  abstract interface
192  function inner_product( self, coefs1, coefs2, form, component ) result( r )
194  import sll_c_maxwell_3d_base
195  class( sll_c_maxwell_3d_base) :: self
196  sll_real64 :: coefs1(:)
197  sll_real64 :: coefs2(:)
198  sll_int32 :: form
199  sll_int32, optional :: component
200  sll_real64 :: r
201 
202  end function inner_product
203  end interface
204 
205 
206 
207  !---------------------------------------------------------------------------!
208  abstract interface
209  subroutine empty( self )
210  import sll_c_maxwell_3d_base
211  class( sll_c_maxwell_3d_base) :: self
212 
213  end subroutine empty
214  end interface
215 
216 
217 
218  !---------------------------------------------------------------------------!
219  abstract interface
220  subroutine multiply_mass( self, deg, coefs_in, coefs_out )
222  import sll_c_maxwell_3d_base
223  class(sll_c_maxwell_3d_base) :: self
224  sll_int32, intent( in ) :: deg(:)
225  sll_real64, intent( in ) :: coefs_in(:)
226  sll_real64, intent( out ) :: coefs_out(:)
227 
228  end subroutine multiply_mass
229  end interface
230  !---------------------------------------------------------------------------!
231  abstract interface
232  subroutine multiply_mass_inverse( self, form, coefs_in, coefs_out )
234  import sll_c_maxwell_3d_base
235  class(sll_c_maxwell_3d_base) :: self
236  sll_int32, intent( in ) :: form
237  sll_real64, intent( in ) :: coefs_in(:)
238  sll_real64, intent( out ) :: coefs_out(:)
239 
240  end subroutine multiply_mass_inverse
241  end interface
242 
243 contains
244 
245  subroutine compute_curl_part( self, delta_t, efield, bfield, betar )
246  class(sll_c_maxwell_3d_base) :: self
247  sll_real64, intent(in) :: delta_t
248  sll_real64, intent(inout) :: efield(:)
249  sll_real64, intent(inout) :: bfield(:)
250  sll_real64, optional :: betar
251 
252 
253  call self%compute_B_from_E( &
254  delta_t, efield, bfield )
255 
256  call self%compute_E_from_B(&
257  delta_t, bfield*betar , efield )
258 
259  end subroutine compute_curl_part
260 
261 
262 end module sll_m_maxwell_3d_base
Module interface to solve Maxwell's equations in 3D.
subroutine compute_curl_part(self, delta_t, efield, bfield, betar)
functions for initial profile of the particle distribution function
Some common numerical utilities.
subroutine, public sll_s_int2string(istep, cstep)
Convert an integer < 9999 to a 4 characters string.
Module to select the kind parameter.
    Report Typos and Errors