Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_poisson_2d_fem_sps_weak_form.F90
Go to the documentation of this file.
2 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
3 
5 
7 
8  implicit none
9 
11 
12  private
13 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
14 
15  ! Working precision
16  integer, parameter :: wp = f64
17 
19 
20  contains
21 
22  procedure :: element_mat => s_poisson_2d_fem_sps_weak_form__element_mat
23  procedure :: element_rhs => s_poisson_2d_fem_sps_weak_form__element_rhs
24 
26 
27 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
28 contains
29 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
30 
31  ! Compute elementary contribution to stiffness and mass matrices
33  self, &
34  test_values_and_derivs_eta1, &
35  test_values_and_derivs_eta2, &
36  trial_values_and_derivs_eta1, &
37  trial_values_and_derivs_eta2, &
38  int_volume, &
39  inv_metric, &
40  Aij, &
41  Mij)
42  class(sll_t_poisson_2d_fem_sps_weak_form), intent(in) :: self
43  real(wp), intent(in) :: test_values_and_derivs_eta1(:, :)
44  real(wp), intent(in) :: test_values_and_derivs_eta2(:, :)
45  real(wp), intent(in) :: trial_values_and_derivs_eta1(:, :)
46  real(wp), intent(in) :: trial_values_and_derivs_eta2(:, :)
47  real(wp), intent(in) :: int_volume(:, :)
48  real(wp), intent(in) :: inv_metric(:, :, :, :)
49  real(wp), intent(out) :: Aij
50  real(wp), intent(out) :: Mij
51 
52  integer :: q1, q2, Nq1, Nq2
53 
54  real(wp) :: v(2), w(2), t(2)
55 
56  ! Extract number of quadrature points
57  nq1 = size(test_values_and_derivs_eta1, 1)
58  nq2 = size(test_values_and_derivs_eta2, 1)
59 
60  aij = 0.0_wp
61  mij = 0.0_wp
62 
63  ! Quadrature
64  do q2 = 1, nq2
65  do q1 = 1, nq1
66 
67  ! 2D gradient in logical test space
68  v = [test_values_and_derivs_eta1(q1, 2)*test_values_and_derivs_eta2(q2, 1), &
69  test_values_and_derivs_eta1(q1, 1)*test_values_and_derivs_eta2(q2, 2)]
70 
71  ! 2D gradient in logical trial space
72  w = [trial_values_and_derivs_eta1(q1, 2)*trial_values_and_derivs_eta2(q2, 1), &
73  trial_values_and_derivs_eta1(q1, 1)*trial_values_and_derivs_eta2(q2, 2)]
74 
75  t = matmul(inv_metric(q1, q2, :, :), w)
76 
77  ! Elementary contribution to stiffness matrix
78  aij = aij + dot_product(v, t)*int_volume(q1, q2)
79 
80  ! Elementary contribution to mass matrix
81  mij = mij + &
82  test_values_and_derivs_eta1(q1, 1)*test_values_and_derivs_eta2(q2, 1)* &
83  trial_values_and_derivs_eta1(q1, 1)*trial_values_and_derivs_eta2(q2, 1)* &
84  int_volume(q1, q2)
85 
86  end do
87  end do
88 
90 
91  ! Compute elementary contribution to right hand side
93  self, &
94  test_values_and_derivs_eta1, &
95  test_values_and_derivs_eta2, &
96  data_2d_rhs, &
97  int_volume, &
98  bi)
99  class(sll_t_poisson_2d_fem_sps_weak_form), intent(in) :: self
100  real(wp), intent(in) :: test_values_and_derivs_eta1(:, :)
101  real(wp), intent(in) :: test_values_and_derivs_eta2(:, :)
102  real(wp), intent(in) :: data_2d_rhs(:, :)
103  real(wp), intent(in) :: int_volume(:, :)
104  real(wp), intent(out) :: bi
105 
106  integer :: q1, q2, Nq1, Nq2
107 
108  ! Extract number of quadrature points
109  nq1 = size(test_values_and_derivs_eta1, 1)
110  nq2 = size(test_values_and_derivs_eta2, 1)
111 
112  bi = 0.0_wp
113 
114  ! Quadrature
115  do q2 = 1, nq2
116  do q1 = 1, nq1
117 
118  ! Elementary contribution to mass matrix
119  bi = bi + &
120  test_values_and_derivs_eta1(q1, 1)*test_values_and_derivs_eta2(q2, 1)* &
121  data_2d_rhs(q1, q2)*int_volume(q1, q2)
122 
123  end do
124  end do
125 
127 
subroutine s_poisson_2d_fem_sps_weak_form__element_rhs(self, test_values_and_derivs_eta1, test_values_and_derivs_eta2, data_2d_rhs, int_volume, bi)
subroutine s_poisson_2d_fem_sps_weak_form__element_mat(self, test_values_and_derivs_eta1, test_values_and_derivs_eta2, trial_values_and_derivs_eta1, trial_values_and_derivs_eta2, int_volume, inv_metric, Aij, Mij)
Module to select the kind parameter.
integer, parameter, public f64
f64 is the kind type for 64-bit reals (double precision)
    Report Typos and Errors