Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_maxwell_1d_base.F90
Go to the documentation of this file.
1 
9  !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
10 #include "sll_working_precision.h"
11 
12  use sll_m_utilities, only: &
14 
15  implicit none
16 
17  public :: &
21 
22  private
23  !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
24 
25  type, abstract :: sll_c_maxwell_1d_base
26  sll_real64 :: lx
27  sll_int32 :: n_cells
28  sll_int32 :: n_dofs0
29  sll_int32 :: n_dofs1
30  sll_real64 :: delta_x
31  sll_int32 :: s_deg_0
32  sll_int32 :: s_deg_1
33  sll_real64 :: mass_solver_tolerance = 1d-12
34  sll_real64 :: poisson_solver_tolerance = 1d-12
35  logical :: strong_ampere = .false.
36 
37  sll_real64 :: solver_tolerance = 1d-12
38 
39  contains
40  procedure(compute_field1_from_field2), deferred :: &
41  compute_e_from_b
42  procedure(compute_field1_from_field2), deferred :: &
43  compute_b_from_e
44  procedure(signature_compute_field_from_field), deferred :: &
45  compute_e_from_rho
46  procedure(signature_compute_field_from_field), deferred :: &
47  compute_rho_from_e
48  procedure(signature_compute_e_from_j_1d), deferred :: &
49  compute_e_from_j
50  procedure(signature_compute_phi_from_field), deferred :: &
51  compute_phi_from_rho
52  procedure(signature_compute_phi_from_field), deferred :: &
53  compute_phi_from_j
54  !procedure(signature_solve), deferred :: &
55  ! solve !< Solve Amperes law and Faraday equation
56  procedure(update_dofs_function), deferred :: &
57  compute_rhs_from_function
58  procedure(update_dofs_function), deferred :: &
59  l2projection
60  procedure(norm_squared), deferred :: &
61  l2norm_squared
62  procedure(signature_inner_product), deferred :: &
63  inner_product
64  procedure(empty), deferred :: &
65  free
66  procedure(signature_mass), deferred :: &
67  multiply_mass
68  procedure(signature_mass), deferred :: &
69  invert_mass
70  procedure(multiply_derivative), deferred :: &
71  multiply_g
72  procedure(multiply_derivative), deferred :: &
73  multiply_gt
74  procedure :: compute_curl_part
75 
76  procedure :: transform_dofs
77 
78  end type sll_c_maxwell_1d_base
79 
80  !---------------------------------------------------------------------------!
81  abstract interface
82  subroutine empty(self)
84  class( sll_c_maxwell_1d_base) :: self
85 
86  end subroutine empty
87  end interface
88 
89  !---------------------------------------------------------------------------!
90  abstract interface
91  function norm_squared(self, coefs_dofs, degree) result( r )
94  class( sll_c_maxwell_1d_base) :: self
95  sll_real64 :: coefs_dofs(:)
96  sll_int32 :: degree
97  sll_real64 :: r
98  end function norm_squared
99  end interface
100 
101 
102  !---------------------------------------------------------------------------!
103  abstract interface
104  subroutine signature_mass(self, in, out, degree)
106  import sll_c_maxwell_1d_base
107  class( sll_c_maxwell_1d_base), intent(inout) :: self
108  sll_real64, intent(in) :: in(:)
109  sll_real64, intent(out) :: out(:)
110  sll_int32, intent(in) :: degree
111  end subroutine signature_mass
112  end interface
113 
114 
115 
116 
117  !---------------------------------------------------------------------------!
118  abstract interface
119  function signature_inner_product(self, coefs1_dofs, coefs2_dofs, degree, degree2) result( r )
121  import sll_c_maxwell_1d_base
122  class( sll_c_maxwell_1d_base) :: self
123  sll_real64 :: coefs1_dofs(:)
124  sll_real64 :: coefs2_dofs(:)
125  sll_int32 :: degree
126  sll_int32, optional :: degree2
127  sll_real64 :: r
128  end function signature_inner_product
129  end interface
130 
131  !---------------------------------------------------------------------------!
132  abstract interface
133 
135  use sll_m_working_precision ! can't pass a header file because the
136  ! preprocessor prevents double inclusion.
137  ! It is very rare.
138  sll_real64 :: sll_i_function_1d_real64
139  sll_real64, intent(in) :: x
140  end function sll_i_function_1d_real64
141  end interface
142  !---------------------------------------------------------------------------!
143  abstract interface
144  subroutine update_dofs_function(self, func, degree, coefs_dofs)
146  import sll_c_maxwell_1d_base
148  class( sll_c_maxwell_1d_base) :: self
149  procedure(sll_i_function_1d_real64) :: func
150  sll_int32, intent(in) :: degree
151  sll_real64, intent(out) :: coefs_dofs(:)
152  end subroutine update_dofs_function
153  end interface
154 
155  !---------------------------------------------------------------------------!
156  abstract interface
157  subroutine compute_field1_from_field2(self, delta_t, field_in, field_out)
159  import sll_c_maxwell_1d_base
160  class(sll_c_maxwell_1d_base) :: self
161  sll_real64, intent(in) :: delta_t
162  sll_real64, intent(in) :: field_in(:)
163  sll_real64, intent(inout) :: field_out(:)
164  end subroutine compute_field1_from_field2
165  end interface
166 
167  !---------------------------------------------------------------------------!
168  abstract interface
169  subroutine signature_compute_field_from_field(self, field_in, field_out)
171  import sll_c_maxwell_1d_base
172  class(sll_c_maxwell_1d_base) :: self
173  sll_real64, intent(in) :: field_in(:)
174  sll_real64, intent(out) :: field_out(:)
176  end interface
177 
178  !---------------------------------------------------------------------------!
179  abstract interface
180  subroutine signature_compute_e_from_j_1d(self, current, component, E)
182  import sll_c_maxwell_1d_base
183  class(sll_c_maxwell_1d_base) :: self
184  sll_real64,dimension(:),intent(in) :: current
185  sll_int32, intent(in) :: component
186  sll_real64,dimension(:),intent(inout) :: e
187  end subroutine signature_compute_e_from_j_1d
188  end interface
189 
190  !---------------------------------------------------------------------------!
191  abstract interface
192  subroutine signature_compute_phi_from_field(self, in, phi, efield)
194  import sll_c_maxwell_1d_base
195  class(sll_c_maxwell_1d_base) :: self
196  sll_real64, intent(in) :: in(:)
197  sll_real64, intent(out) :: phi(:)
198  sll_real64, intent(out) :: efield(:)
199  end subroutine signature_compute_phi_from_field
200  end interface
201 
202  !---------------------------------------------------------------------------!
203  abstract interface
204  subroutine multiply_derivative(self, in, out)
206  import sll_c_maxwell_1d_base
207  class(sll_c_maxwell_1d_base), intent(in) :: self
208  sll_real64, intent(in) :: in(:)
209  sll_real64, intent(out) :: out(:)
210  end subroutine multiply_derivative
211  end interface
212 
213 
214 contains
216  subroutine sll_s_plot_two_fields_1d(fname, n1, f1, f2, iplot, time )
217  character(len=*), intent(in) :: fname
218  sll_int32, intent(in) :: n1
219  sll_real64, dimension(n1), intent(in) :: f1
220  sll_real64, dimension(n1), intent(in) :: f2
221  sll_int32, intent(in) :: iplot
222  sll_real64, intent(in) :: time
223 
224  integer :: i
225  character(len=4) :: cplot
226 
227  call sll_s_int2string(iplot, cplot)
228 
229  !write domains
230  open( 80, file = fname//cplot//".dat" )
231  do i=1,n1
232  write(80,*) i, sngl(f1(i)), sngl(f2(i))
233  end do
234  close(80)
235 
236 !!$ open( 90, file = fname//'plots.gnu', position="append" )
237 !!$ if ( iplot == 1 ) then
238 !!$ rewind(90)
239 !!$ !write(90,*)"set xr[-0.1:1.1]"
240 !!$ !write(90,*)"set yr[-0.1:1.1]"
241 !!$ !write(90,*)"set zr[-1.1:1.1]"
242 !!$ !write(90,*)"set cbrange[-1:1]"
243 !!$ !write(90,*)"set pm3d"
244 !!$ !write(90,*)"set surf"
245 !!$ write(90,*)"set term x11"
246 !!$ end if
247 !!$ write(90,*)"set title 'Time = ",time,"'"
248 !!$ write(90,"(a)",advance='no')"plot '"//fname//cplot//".dat' w lines"
249 !!$ write(90,"(a)")",'"//fname//cplot//".dat' u 1:3 w lines"
250 !!$ write(90,*)"pause -1"
251 !!$ close(90)
252 
253  end subroutine sll_s_plot_two_fields_1d
254 
256  subroutine transform_dofs(self, in, out, degree)
257  class(sll_c_maxwell_1d_base), intent(inout) :: self
258  sll_int32, intent(in) :: degree
259  sll_real64, intent(in) :: in(:)
260  sll_real64, intent(out) :: out(:)
261  ! local variables
262 
263  out = in
264 
265  end subroutine transform_dofs
266 
268  subroutine compute_curl_part( self, delta_t, efield, bfield, betar )
269  class(sll_c_maxwell_1d_base) :: self
270  sll_real64, intent(in) :: delta_t
271  sll_real64, intent(inout) :: efield(:)
272  sll_real64, intent(inout) :: bfield(:)
273  sll_real64, optional :: betar
274  !local variables
275  sll_real64 :: factor
276 
277  if( present(betar) ) then
278  factor = betar
279  else
280  factor = 1._f64
281  end if
282 
283 
284  call self%compute_B_from_E( &
285  delta_t, efield, bfield )
286 
287  call self%compute_E_from_B(&
288  delta_t, bfield*factor , efield )
289 
290  end subroutine compute_curl_part
291 
292 end module sll_m_maxwell_1d_base
Module interface to solve Maxwell's equations in 1D.
subroutine compute_curl_part(self, delta_t, efield, bfield, betar)
Solve source-free Maxwell's equations: Default implementation by solving first B_from_E then E_from_B...
subroutine transform_dofs(self, in, out, degree)
Transformation of the dofs (identity by default)
subroutine, public sll_s_plot_two_fields_1d(fname, n1, f1, f2, iplot, time)
write files to visualize 1d fields with gnuplot
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