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_dense_projector.F90
Go to the documentation of this file.
2 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
3 #include "sll_assert.h"
4 
6 
7  implicit none
8 
10 
11  private
12 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
13 
14  ! Working precision
15  integer, parameter :: wp = f64
16 
18 
19  integer :: n1
20  integer :: n2
21 
22  ! Temporary storage needed for projection
23  real(wp), allocatable :: qp_temp(:, :)
24  real(wp), allocatable :: l(:, :)
25  real(wp), allocatable :: lt(:, :) ! transpose
26 
27  contains
28 
32  procedure :: change_basis_vector_inv => s_poisson_2d_fem_sps_dense_projector__change_basis_vector_inv
34 
36 
37 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
38 contains
39 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
40 
41  ! Initializer
42  subroutine s_poisson_2d_fem_sps_dense_projector__init(self, n1, n2, L)
43  class(sll_t_poisson_2d_fem_sps_dense_projector), intent(inout) :: self
44  integer, intent(in) :: n1
45  integer, intent(in) :: n2
46  real(wp), intent(in) :: L(:, :) ! matrix of barycentric coordinates
47 
48  integer :: nn
49 
50  self%n1 = n1
51  self%n2 = n2
52 
53  nn = (n1 - 2)*n2
54 
55  ! Allocate temporary storage
56  allocate (self%Qp_temp(2*n2, 3))
57  allocate (self%L(size(l, 1), size(l, 2)))
58  allocate (self%Lt(size(l, 2), size(l, 1)))
59 
60  sll_assert(size(l, 1) == 2*n2)
61  sll_assert(size(l, 2) == 3)
62 
63  self%L = l
64  self%Lt = transpose(l)
65 
67 
68  ! Change basis: C1 projection of stiffness and mass matrices
69  ! NOTE: 'self' has intent(inout) because temporary storage has to be assigned
71  class(sll_t_poisson_2d_fem_sps_dense_projector), intent(inout) :: self
72  real(wp), intent(in) :: Q(:, :)
73  real(wp), intent(inout) :: Qp(:, :)
74 
75  integer :: nn
76 
77  associate(n1 => self%n1, n2 => self%n2)
78 
79  nn = (n1 - 2)*n2
80 
81  ! Checks
82  sll_assert(size(q, 1) == n1*n2)
83  sll_assert(size(q, 2) == n1*n2)
84  sll_assert(size(qp, 1) == 3 + nn)
85  sll_assert(size(qp, 2) == 3 + nn)
86 
87  ! Fill block 1: 3 x 3
88  self%Qp_temp = matmul(q(1:2*n2, 1:2*n2), self%L)
89  qp(1:3, 1:3) = matmul(self%Lt, self%Qp_temp)
90 
91  ! Fill block 2: 3 x (n1-2)*n2
92  qp(1:3, 4:3 + nn) = matmul(self%Lt, q(1:2*n2, 2*n2 + 1:n1*n2))
93 
94  ! Fill block 3: (n1-2)*n2 x 3
95  qp(4:3 + nn, 1:3) = matmul(q(2*n2 + 1:n1*n2, 1:2*n2), self%L)
96 
97  ! Fill block 4: (n1-2)*n2 x (n1-2)*n2
98  qp(4:3 + nn, 4:3 + nn) = q(2*n2 + 1:n1*n2, 2*n2 + 1:n1*n2)
99 
100  end associate
101 
103 
104  ! Change basis: C1 projection of vectors
106  class(sll_t_poisson_2d_fem_sps_dense_projector), intent(in) :: self
107  real(wp), intent(in) :: V(:)
108  real(wp), intent(inout) :: Vp(:)
109 
110  integer :: nn
111 
112  associate(n1 => self%n1, n2 => self%n2)
113 
114  nn = (n1 - 2)*n2
115 
116  ! Checks
117  sll_assert(size(v) == n1*n2)
118  sll_assert(size(vp) == 3 + nn)
119 
120  vp(1:3) = matmul(self%Lt, v(1:2*n2))
121  vp(4:3 + nn) = v(2*n2 + 1:n1*n2)
122 
123  end associate
124 
126 
127  ! Change basis: C1 projection of vectors
129  class(sll_t_poisson_2d_fem_sps_dense_projector), intent(in) :: self
130  real(wp), intent(in) :: Vp(:)
131  real(wp), intent(inout) :: V(:)
132 
133  integer :: nn
134 
135  associate(n1 => self%n1, n2 => self%n2)
136 
137  nn = (n1 - 2)*n2
138 
139  ! Checks
140  sll_assert(size(vp) == 3 + nn)
141  sll_assert(size(v) == n1*n2)
142 
143  v(1:2*n2) = matmul(self%L, vp(1:3))
144  v(2*n2 + 1:n1*n2) = vp(4:3 + nn)
145 
146  end associate
147 
149 
150  ! Deallocate allocatables
152  class(sll_t_poisson_2d_fem_sps_dense_projector), intent(inout) :: self
153 
154  deallocate (self%Qp_temp)
155  deallocate (self%L)
156  deallocate (self%Lt)
157 
159 
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