Octopus
basis_vectors.F90
Go to the documentation of this file.
1!! Copyright (C) 2021 N. Tancogne-Dejean, M. Oliveira
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
22 use debug_oct_m
23 use global_oct_m
26 use math_oct_m
30
31 implicit none
32
33 private
34
35 public :: &
37
57 ! Components are public by default
58 integer, private :: dim
59 real(real64), allocatable :: vectors(:, :)
60 real(real64), allocatable :: change_of_basis_matrix(:, :)
61 ! !! a Cartesian basis to this basis
62 logical :: orthogonal
63 contains
64 procedure :: copy => basis_vectors_copy
65 generic :: assignment(=) => copy
66 procedure :: scale => basis_vectors_scale
67 generic :: from_cartesian => from_cartesian1, from_cartesian2
68 generic :: to_cartesian => to_cartesian1, to_cartesian2
69 procedure :: from_cartesian1 => basis_vectors_from_cartesian1
70 procedure :: to_cartesian1 => basis_vectors_to_cartesian1
71 procedure :: from_cartesian2 => basis_vectors_from_cartesian2
72 procedure :: to_cartesian2 => basis_vectors_to_cartesian2
74 end type basis_vectors_t
75
76 interface basis_vectors_t
77 module procedure basis_vectors_constructor
78 end interface basis_vectors_t
79
80contains
81
82 !--------------------------------------------------------------
83 type(basis_vectors_t) function basis_vectors_constructor(namespace, dim, vectors) result(basis)
84 type(namespace_t), intent(in) :: namespace
85 integer, intent(in) :: dim
86 real(real64), intent(in) :: vectors(dim, dim)
87
88 integer :: idir1, idir2
89
90 push_sub(basis_vectors_constructor)
91
92 basis%dim = dim
93
94 safe_allocate(basis%vectors(1:dim, 1:dim))
95 safe_allocate(basis%change_of_basis_matrix(1:dim, 1:dim))
96
97 basis%vectors(:, :) = vectors(:, :)
98
99 basis%orthogonal = .true.
100 do idir1 = 1, dim
101 do idir2 = idir1 + 1, dim
102 if (abs(dot_product(basis%vectors(:, idir1), basis%vectors(:, idir2))) > m_epsilon) then
103 basis%orthogonal = .false.
104 exit
105 end if
106 end do
107 end do
108
109 if (abs(lalg_determinant(dim, basis%vectors, preserve_mat = .true.)) < m_epsilon) then
110 message(1) = "Basis vectors are not linearly independent and therefore the change-of-basis matrix cannot be calculated."
111 call messages_fatal(1, namespace=namespace)
112 end if
113 call calculate_change_of_basis_matrix(dim, basis%vectors, basis%change_of_basis_matrix)
115 pop_sub(basis_vectors_constructor)
116 end function basis_vectors_constructor
117
118 !--------------------------------------------------------------
119 subroutine basis_vectors_copy(this, source)
120 class(basis_vectors_t), intent(out) :: this
121 class(basis_vectors_t), intent(in) :: source
122
123 push_sub(basis_vectors_copy)
124
125 this%dim = source%dim
126 safe_allocate_source_a(this%vectors, source%vectors)
127 safe_allocate_source_a(this%change_of_basis_matrix, source%change_of_basis_matrix)
128 this%orthogonal = source%orthogonal
129
130 pop_sub(basis_vectors_copy)
131 end subroutine basis_vectors_copy
132
133 !--------------------------------------------------------------
134 subroutine basis_vectors_finalize(this)
135 type(basis_vectors_t), intent(inout) :: this
136
137 push_sub(basis_vectors_finalize)
138
139 safe_deallocate_a(this%vectors)
140 safe_deallocate_a(this%change_of_basis_matrix)
141
143 end subroutine basis_vectors_finalize
144
145 !--------------------------------------------------------------
146 subroutine basis_vectors_scale(this, factor)
147 class(basis_vectors_t), intent(inout) :: this
148 real(real64), intent(in) :: factor(1:this%dim)
150 integer :: idir
154 ! Scale the basis in real space
155 do idir = 1, this%dim
156 this%vectors(:, idir) = this%vectors(:, idir)*factor(idir)
157 end do
159 ! Calculate the new change-of-basis matrix
160 call calculate_change_of_basis_matrix(this%dim, this%vectors, this%change_of_basis_matrix)
163 end subroutine basis_vectors_scale
165 !--------------------------------------------------------------
166 pure function basis_vectors_from_cartesian1(this, xx_cart) result(xx)
167 class(basis_vectors_t), intent(in) :: this
168 real(real64), intent(in) :: xx_cart(1:this%dim)
169 real(real64) :: xx(1:this%dim)
170
171 ! no PUSH_SUB, called too often
172
173 xx = matmul(this%change_of_basis_matrix, xx_cart)
174
177 !--------------------------------------------------------------
178 pure function basis_vectors_to_cartesian1(this, xx) result(xx_cart)
179 class(basis_vectors_t), intent(in) :: this
180 real(real64), intent(in) :: xx(1:this%dim)
181 real(real64) :: xx_cart(1:this%dim)
182
183 ! no PUSH_SUB, called too often
184
185 xx_cart = matmul(this%vectors, xx)
186
187 end function basis_vectors_to_cartesian1
188
189 !--------------------------------------------------------------
190 function basis_vectors_from_cartesian2(this, np, xx_cart) result(xx)
191 class(basis_vectors_t), intent(in) :: this
192 integer, intent(in) :: np
193 real(real64), contiguous, intent(in) :: xx_cart(:, :)
194 real(real64) :: xx(np, this%dim)
195
196 ! no PUSH_SUB, called too often
197
198 call lalg_gemm_nc(np, this%dim, this%dim, m_one, xx_cart, this%change_of_basis_matrix, m_zero, xx)
199
201
202 !--------------------------------------------------------------
203 function basis_vectors_to_cartesian2(this, np, xx) result(xx_cart)
204 class(basis_vectors_t), intent(in) :: this
205 integer, intent(in) :: np
206 real(real64), contiguous, intent(in) :: xx(:, :)
207 real(real64) :: xx_cart(np, this%dim)
208
209 ! no PUSH_SUB, called too often
210 call lalg_gemm_nc(np, this%dim, this%dim, m_one, xx, this%vectors, m_zero, xx_cart)
211
213
214 !--------------------------------------------------------------
215 subroutine calculate_change_of_basis_matrix(dim, vectors, matrix)
216 integer, intent(in) :: dim
217 real(real64), intent(in) :: vectors(1:dim, 1:dim)
218 real(real64), intent(out) :: matrix(1:dim, 1:dim)
219
220 real(real64) :: volume, cross(1:3)
221
223
224 select case (dim)
225 case (3)
226 cross(1:3) = dcross_product(vectors(1:3, 2), vectors(1:3, 3))
227 volume = dot_product(vectors(1:3, 1), cross(1:3))
228
229 matrix(1, 1:3) = dcross_product(vectors(:, 2), vectors(:, 3))/volume
230 matrix(2, 1:3) = dcross_product(vectors(:, 3), vectors(:, 1))/volume
231 matrix(3, 1:3) = dcross_product(vectors(:, 1), vectors(:, 2))/volume
232 case (2)
233 volume = vectors(1, 1)*vectors(2, 2) - vectors(2, 1)*vectors(1, 2)
234 matrix(1, 1) = vectors(2, 2)/volume
235 matrix(1, 2) = -vectors(1, 2)/volume
236 matrix(2, 1) = -vectors(2, 1)/volume
237 matrix(2, 2) = vectors(1, 1)/volume
238 case (1)
239 volume = vectors(1, 1)
240 matrix(1, 1) = m_one / vectors(1, 1)
241 case default ! dim > 3
242 matrix = vectors
243 call lalg_inverse(dim, matrix, 'dir')
244 end select
245
248
249end module basis_vectors_oct_m
250
251!! Local Variables:
252!! mode: f90
253!! coding: utf-8
254!! End:
Note that lalg_determinant and lalg_inverse are just wrappers over the same routine.
Definition: lalg_adv.F90:194
subroutine basis_vectors_copy(this, source)
real(real64) function, dimension(np, this%dim) basis_vectors_to_cartesian2(this, np, xx)
real(real64) function, dimension(np, this%dim) basis_vectors_from_cartesian2(this, np, xx_cart)
subroutine calculate_change_of_basis_matrix(dim, vectors, matrix)
pure real(real64) function, dimension(1:this%dim) basis_vectors_from_cartesian1(this, xx_cart)
subroutine basis_vectors_scale(this, factor)
subroutine basis_vectors_finalize(this)
pure real(real64) function, dimension(1:this%dim) basis_vectors_to_cartesian1(this, xx)
type(basis_vectors_t) function basis_vectors_constructor(namespace, dim, vectors)
real(real64), parameter, public m_epsilon
Definition: global.F90:204
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:414
Vectors defining a basis in a vector space. This class provides methods to convert vector coordinates...
int true(void)