3 #include "sll_assert.h"
23 real(
wp),
allocatable :: qp_temp(:, :)
24 real(
wp),
allocatable :: l(:, :)
25 real(
wp),
allocatable :: lt(:, :)
44 integer,
intent(in) :: n1
45 integer,
intent(in) :: n2
46 real(wp),
intent(in) :: L(:, :)
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)))
60 sll_assert(
size(l, 1) == 2*n2)
61 sll_assert(
size(l, 2) == 3)
64 self%Lt = transpose(l)
72 real(wp),
intent(in) :: Q(:, :)
73 real(wp),
intent(inout) :: Qp(:, :)
77 associate(n1 => self%n1, n2 => self%n2)
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)
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)
92 qp(1:3, 4:3 + nn) = matmul(self%Lt, q(1:2*n2, 2*n2 + 1:n1*n2))
95 qp(4:3 + nn, 1:3) = matmul(q(2*n2 + 1:n1*n2, 1:2*n2), self%L)
98 qp(4:3 + nn, 4:3 + nn) = q(2*n2 + 1:n1*n2, 2*n2 + 1:n1*n2)
107 real(wp),
intent(in) :: V(:)
108 real(wp),
intent(inout) :: Vp(:)
112 associate(n1 => self%n1, n2 => self%n2)
117 sll_assert(
size(v) == n1*n2)
118 sll_assert(
size(vp) == 3 + nn)
120 vp(1:3) = matmul(self%Lt, v(1:2*n2))
121 vp(4:3 + nn) = v(2*n2 + 1:n1*n2)
130 real(wp),
intent(in) :: Vp(:)
131 real(wp),
intent(inout) :: V(:)
135 associate(n1 => self%n1, n2 => self%n2)
140 sll_assert(
size(vp) == 3 + nn)
141 sll_assert(
size(v) == n1*n2)
143 v(1:2*n2) = matmul(self%L, vp(1:3))
144 v(2*n2 + 1:n1*n2) = vp(4:3 + nn)
154 deallocate (self%Qp_temp)
subroutine s_poisson_2d_fem_sps_dense_projector__change_basis_vector(self, V, Vp)
subroutine s_poisson_2d_fem_sps_dense_projector__change_basis_vector_inv(self, Vp, V)
subroutine s_poisson_2d_fem_sps_dense_projector__free(self)
subroutine s_poisson_2d_fem_sps_dense_projector__init(self, n1, n2, L)
subroutine s_poisson_2d_fem_sps_dense_projector__change_basis_matrix(self, Q, Qp)
Module to select the kind parameter.
integer, parameter, public f64
f64 is the kind type for 64-bit reals (double precision)