Octopus
lattice_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
24 use, intrinsic :: iso_fortran_env
25 use math_oct_m
28 use parser_oct_m
30 use space_oct_m
31 use unit_oct_m
33
34 implicit none
35
36 private
37
38 public :: &
41
43 ! Components are public by default
44 type(space_t), private :: space
45 real(real64), allocatable :: rlattice_primitive(:,:)
46 real(real64), allocatable :: rlattice (:,:)
47 real(real64), allocatable :: klattice_primitive(:,:)
48 real(real64), allocatable :: klattice (:,:)
49 real(real64) :: alpha, beta, gamma
50 real(real64) :: rcell_volume
51 logical :: nonorthogonal = .false.
52 ! Some notes:
53 ! rlattice_primitive is the A matrix from Chelikowski PRB 78 075109 (2008)
54 ! klattice_primitive is the transpose (!) of the B matrix, with no 2 pi factor included
55 ! klattice is the proper reciprocal lattice vectors, with 2 pi factor, and in units of 1/bohr
56 ! The F matrix of Chelikowski is matmul(transpose(latt%klattice_primitive), latt%klattice_primitive)
57 contains
58 procedure :: copy => lattice_vectors_copy
59 generic :: assignment(=) => copy
60 procedure :: scale => lattice_vectors_scale
61 procedure :: write_info => lattice_vectors_write_info
62 procedure :: short_info => lattice_vectors_short_info
63 procedure :: cart_to_red => lattice_vectors_cart_to_red
64 procedure :: red_to_cart => lattice_vectors_red_to_cart
65 procedure :: fold_into_cell => lattice_vectors_fold_into_cell
66 procedure :: max_length => lattice_vectors_max_length
67 procedure :: update => lattice_vectors_update
69 end type lattice_vectors_t
70
71 interface lattice_vectors_t
73 end interface lattice_vectors_t
74
79 private
80 integer, public :: n_cells = 0
81 integer, allocatable :: icell(:,:)
82 type(lattice_vectors_t), pointer :: latt => null()
83 contains
84 procedure :: copy => lattice_iterator_copy
85 generic :: assignment(=) => copy
86 procedure :: get => lattice_iterator_get
88 end type lattice_iterator_t
89
90 interface lattice_iterator_t
91 module procedure lattice_iterator_constructor
92 end interface lattice_iterator_t
93
94contains
95
96 !--------------------------------------------------------------
97 type(lattice_vectors_t) function lattice_vectors_constructor_from_rlattice(namespace, space, rlattice) result(latt)
98 type(namespace_t), intent(in) :: namespace
99 class(space_t), intent(in) :: space
100 real(real64), intent(in) :: rlattice(:, :)
101
102 integer :: idir, idir2
103 real(real64) :: volume_element
104
105 push_sub(lattice_vectors_constructor_from_rlattice)
106
107 latt%space = space
108
109 safe_allocate(latt%rlattice_primitive(1:space%dim, 1:space%dim))
110 safe_allocate(latt%rlattice(1:space%dim, 1:space%dim))
111 safe_allocate(latt%klattice_primitive(1:space%dim, 1:space%dim))
112 safe_allocate(latt%klattice(1:space%dim, 1:space%dim))
113
114 latt%rlattice(1:space%dim, 1:space%dim) = rlattice(1:space%dim, 1:space%dim)
115 do idir = 1, space%dim
116 latt%rlattice_primitive(:, idir) = latt%rlattice(:, idir) / norm2(latt%rlattice(:, idir))
117 end do
118
119 latt%nonorthogonal = .false.
120 do idir = 1, space%dim
121 do idir2 = 1, space%dim
122 if (idir /= idir2 .and. abs(latt%rlattice_primitive(idir, idir2)) > m_epsilon) then
123 latt%nonorthogonal = .true.
124 end if
125 end do
126 end do
127
128 call reciprocal_lattice(latt%rlattice, latt%klattice, latt%rcell_volume, space%dim, namespace)
129 latt%klattice = latt%klattice * m_two*m_pi
130
131 call reciprocal_lattice(latt%rlattice_primitive, latt%klattice_primitive, volume_element, space%dim, namespace)
132
133 if (space%dim == 3) then
134 call angles_from_rlattice_primitive(latt%rlattice_primitive, latt%alpha, latt%beta, latt%gamma)
135 else
136 ! Angles should not be used, so set them to zero
137 latt%alpha = m_zero
138 latt%beta = m_zero
139 latt%gamma = m_zero
140 end if
142 pop_sub(lattice_vectors_constructor_from_rlattice)
145 !--------------------------------------------------------------
146 type(lattice_vectors_t) function lattice_vectors_constructor_from_input(namespace, space, variable_prefix) result(latt)
147 type(namespace_t), intent(in) :: namespace
148 class(space_t), intent(in) :: space
149 character(len=*), optional, intent(in) :: variable_prefix
150
151 type(block_t) :: blk
152 real(real64) :: norm, lparams(space%dim), volume_element
153 integer :: idim, jdim, ncols
154 logical :: has_angles
155 real(real64) :: angles(1:space%dim)
156 character(len=:), allocatable :: prefix
160 if (present(variable_prefix)) then
161 prefix = variable_prefix
162 else
163 prefix = ''
164 end if
165
166 latt%space = space
167
168 safe_allocate(latt%rlattice_primitive(1:space%dim, 1:space%dim))
169 safe_allocate(latt%rlattice(1:space%dim, 1:space%dim))
170 safe_allocate(latt%klattice_primitive(1:space%dim, 1:space%dim))
171 safe_allocate(latt%klattice(1:space%dim, 1:space%dim))
172
173 latt%alpha = 90.0_real64
174 latt%beta = 90.0_real64
175 latt%gamma = 90.0_real64
176
177 has_angles = .false.
178 angles = 90.0_real64
180 if (space%is_periodic()) then
181
182 !%Variable LatticeParameters
183 !%Type block
184 !%Section Mesh::Simulation Box
185 !%Description
186 !% The lattice parameters (a, b, c).
187 !% This variable is mandatory for periodic systems and is ignored otherwise.
188 !% When PeriodicDimensions = 3, a second optional line can be used to
189 !% define the angles between the lattice vectors. If the angles are not
190 !% provided, then the variable LatticeVectors must be set.
191 !% The number of parameters specified in the block must be at least equal
192 !% to the number of periodic dimensions, but it is not mandatory to
193 !% specify parameters for the non-periodic dimensions (in that case they
194 !% are set to 1).
195 !%End
196 if (parse_block(namespace, prefix//'LatticeParameters', blk) == 0) then
197 ncols = parse_block_cols(blk, 0)
198 if (ncols < space%periodic_dim) then
199 call messages_input_error(namespace, prefix//'LatticeParameters', &
200 'The number of columns must be at least PeriodicDimensions')
201 end if
202 do idim = 1, ncols
203 call parse_block_float(blk, 0, idim - 1, lparams(idim))
204 end do
205
206 ! If some parameters for non-periodic dimensions are not set in the input file, with set them to 1.
207 do idim = ncols + 1, space%dim
208 lparams(idim) = m_one
209 end do
210
211 ! Parse angles, if available
212 if (parse_block_n(blk) > 1) then
213 if (space%dim /= 3) then
214 call messages_input_error(namespace, prefix//'LatticeParameters', 'Angles can only be specified when Dimensions = 3')
215 end if
216
217 ncols = parse_block_cols(blk, 1)
218 if (ncols /= space%dim) then
219 call messages_input_error(namespace, prefix//'LatticeParameters', 'You must specify three angles')
220 end if
221 do idim = 1, space%dim
222 call parse_block_float(blk, 1, idim - 1, angles(idim))
223 end do
224 has_angles = .true.
225 end if
226 call parse_block_end(blk)
227
228 else
229 call messages_input_error(namespace, prefix//'LatticeParameters', 'Variable is mandatory for periodic systems')
230 end if
231
232
233 if (has_angles) then
234 latt%alpha = angles(1)
235 latt%beta = angles(2)
236 latt%gamma = angles(3)
237
238 if (parse_is_defined(namespace, prefix//'LatticeVectors')) then
239 message(1) = 'LatticeParameters with angles is incompatible with LatticeVectors'
240 call messages_print_var_info(prefix//"LatticeParameters", namespace=namespace)
241 call messages_fatal(1, namespace=namespace)
242 end if
243
244 call build_metric_from_angles(latt, angles)
245
246 else
247
248 !%Variable LatticeVectors
249 !%Type block
250 !%Default simple cubic
251 !%Section Mesh::Simulation Box
252 !%Description
253 !% Primitive lattice vectors. Vectors are stored in rows.
254 !% LatticeVectors are normalized and then used, if both LatticeParameters and
255 !% LatticeVectors are defined, norm of LatticeVectors are
256 !% used as multipliers of lattice parameters.
257 !% Default:
258 !% <br><br><tt>%LatticeVectors
259 !% <br>&nbsp;&nbsp;1.0 | 0.0 | 0.0
260 !% <br>&nbsp;&nbsp;0.0 | 1.0 | 0.0
261 !% <br>&nbsp;&nbsp;0.0 | 0.0 | 1.0
262 !% <br>%<br></tt>
263 !%End
264 latt%rlattice_primitive = diagonal_matrix(space%dim, m_one)
265 latt%nonorthogonal = .false.
266
267 if (parse_block(namespace, prefix//'LatticeVectors', blk) == 0) then
268 do idim = 1, space%dim
269 do jdim = 1, space%dim
270 call parse_block_float(blk, idim - 1, jdim - 1, latt%rlattice_primitive(jdim, idim))
271 if (idim /= jdim .and. abs(latt%rlattice_primitive(jdim, idim)) > m_epsilon) then
272 latt%nonorthogonal = .true.
273 end if
274 end do
275 end do
276 call parse_block_end(blk)
277
278 end if
279 end if
280
281 else
282 ! Non-periodic
283 lparams = m_one
284 latt%rlattice_primitive = diagonal_matrix(space%dim, m_one)
285 end if
286
287 latt%rlattice = m_zero
288 do idim = 1, space%dim
289 norm = norm2(latt%rlattice_primitive(1:space%dim, idim))
290 lparams(idim) = lparams(idim)*norm
291 do jdim = 1, space%dim
292 latt%rlattice_primitive(jdim, idim) = latt%rlattice_primitive(jdim, idim) / norm
293 latt%rlattice(jdim, idim) = latt%rlattice_primitive(jdim, idim) * lparams(idim)
294 end do
295 end do
296
297 call reciprocal_lattice(latt%rlattice, latt%klattice, latt%rcell_volume, space%dim, namespace)
298 latt%klattice = latt%klattice * m_two*m_pi
299
300 call reciprocal_lattice(latt%rlattice_primitive, latt%klattice_primitive, volume_element, space%dim, namespace)
301
302 if (.not. has_angles .and. space%dim == 3) then
303 !We compute the angles from the lattice vectors
304 call angles_from_rlattice_primitive(latt%rlattice_primitive, latt%alpha, latt%beta, latt%gamma)
305 end if
306
309
310 !--------------------------------------------------------------
311 subroutine angles_from_rlattice_primitive(rlattice_primitive, alpha, beta, gamma)
312 real(real64), intent(in) :: rlattice_primitive(1:3, 1:3)
313 real(real64), intent(out) :: alpha
314 real(real64), intent(out) :: beta
315 real(real64), intent(out) :: gamma
316
317 real(real64) :: rlatt(1:3, 1:3)
318
320
321 rlatt = matmul(transpose(rlattice_primitive), rlattice_primitive)
322 alpha = acos(rlatt(2, 3)/sqrt(rlatt(2, 2)*rlatt(3, 3)))/m_pi*180.0_real64
323 beta = acos(rlatt(1, 3)/sqrt(rlatt(1, 1)*rlatt(3, 3)))/m_pi*180.0_real64
324 gamma = acos(rlatt(1, 2)/sqrt(rlatt(1, 1)*rlatt(2, 2)))/m_pi*180.0_real64
325
327 end subroutine angles_from_rlattice_primitive
328
329 !--------------------------------------------------------------
330 subroutine lattice_vectors_copy(this, source)
331 class(lattice_vectors_t), intent(out) :: this
332 class(lattice_vectors_t), intent(in) :: source
333
334 push_sub(lattice_vectors_copy)
335
336 this%space = source%space
337 safe_allocate_source_a(this%rlattice_primitive, source%rlattice_primitive)
338 safe_allocate_source_a(this%rlattice, source%rlattice)
339 safe_allocate_source_a(this%klattice_primitive, source%klattice_primitive)
340 safe_allocate_source_a(this%klattice, source%klattice)
341 this%alpha = source%alpha
342 this%beta = source%beta
343 this%gamma = source%gamma
344 this%rcell_volume = source%rcell_volume
345 this%nonorthogonal = source%nonorthogonal
346
347 pop_sub(lattice_vectors_copy)
348 end subroutine lattice_vectors_copy
349
350 !--------------------------------------------------------------
351 subroutine lattice_vectors_finalize(this)
352 type(lattice_vectors_t), intent(inout) :: this
353
355
356 safe_deallocate_a(this%rlattice_primitive)
357 safe_deallocate_a(this%rlattice)
358 safe_deallocate_a(this%klattice_primitive)
359 safe_deallocate_a(this%klattice)
360 this%nonorthogonal = .false.
361
363 end subroutine lattice_vectors_finalize
364
365 !--------------------------------------------------------------
366 subroutine lattice_vectors_scale(this, factor)
367 class(lattice_vectors_t), intent(inout) :: this
368 real(real64), intent(in) :: factor(this%space%dim)
369
370 integer :: idir
371
372 push_sub(lattice_vectors_scale)
373
374 ! Scale the lattice in real space
375 do idir = 1, this%space%dim
376 this%rlattice(1:this%space%dim, idir) = this%rlattice(1:this%space%dim, idir)*factor(idir)
377 end do
378
379 ! Regenerate the lattice in reciprocal space
380 call reciprocal_lattice(this%rlattice, this%klattice, this%rcell_volume, this%space%dim)
381 this%klattice = this%klattice * m_two*m_pi
382
383 pop_sub(lattice_vectors_scale)
384 end subroutine lattice_vectors_scale
385
386 !--------------------------------------------------------------
388 subroutine lattice_vectors_update(this, rlattice)
389 class(lattice_vectors_t), intent(inout) :: this
390 real(real64), intent(in) :: rlattice(this%space%dim, this%space%dim)
391
392 integer :: idir
393 real(real64) :: volume_element, maxrlatt
394
395 push_sub(lattice_vectors_update)
396
397 this%rlattice = rlattice
398 maxrlatt = maxval(abs(this%rlattice))
399 where(abs(this%rlattice) < 1e-12_real64*maxrlatt) this%rlattice = m_zero
400
401 ! Update the primitive lattice vectors
402 do idir = 1, this%space%dim
403 this%rlattice_primitive(:, idir) = this%rlattice(:, idir) / norm2(this%rlattice(:, idir))
404 end do
405
406 ! Regenerate the lattice in reciprocal space
407 call reciprocal_lattice(this%rlattice, this%klattice, this%rcell_volume, this%space%dim)
408 this%klattice = this%klattice * m_two*m_pi
409
410 call reciprocal_lattice(this%rlattice_primitive, this%klattice_primitive, volume_element, this%space%dim)
411
412 ! Recompute the angles
413 if (this%space%dim == 3) then
414 call angles_from_rlattice_primitive(this%rlattice_primitive, this%alpha, this%beta, this%gamma)
415 else
416 ! Angles should not be used, so set them to zero
417 this%alpha = m_zero
418 this%beta = m_zero
419 this%gamma = m_zero
420 end if
421
424
425 !--------------------------------------------------------------
426 pure function lattice_vectors_cart_to_red(this, xx_cart) result(xx_red)
427 class(lattice_vectors_t), intent(in) :: this
428 real(real64), intent(in) :: xx_cart(this%space%dim)
429 real(real64) :: xx_red(this%space%dim)
430
431 xx_red = matmul(xx_cart, this%klattice)/(m_two*m_pi)
432
433 end function lattice_vectors_cart_to_red
434
435 !--------------------------------------------------------------
436 pure function lattice_vectors_red_to_cart(this, xx_red) result(xx_cart)
437 class(lattice_vectors_t), intent(in) :: this
438 real(real64), intent(in) :: xx_red(this%space%dim)
439 real(real64) :: xx_cart(this%space%dim)
440
441 xx_cart = matmul(this%rlattice, xx_red)
442
443 end function lattice_vectors_red_to_cart
445 !--------------------------------------------------------------
446 function lattice_vectors_fold_into_cell(this, xx) result(new_xx)
447 class(lattice_vectors_t), intent(in) :: this
448 real(real64), intent(in) :: xx(this%space%dim)
449 real(real64) :: new_xx(this%space%dim)
450
451 integer :: idir
452
453 if (this%space%is_periodic()) then
454 ! Convert the position to reduced coordinates
455 new_xx = this%cart_to_red(xx)
456
457 do idir = 1, this%space%periodic_dim
458 ! Change of origin
459 new_xx(idir) = new_xx(idir) + m_half
460
461 ! Fold into cell
462 new_xx(idir) = new_xx(idir) - anint(new_xx(idir))
463 if (new_xx(idir) < -1.0e-6_real64) then
464 new_xx(idir) = new_xx(idir) + m_one
465 end if
466
467 ! Sanity checks
468 assert(new_xx(idir) >= -1.0e-6_real64)
469 assert(new_xx(idir) < m_one)
470
471 ! Change origin back
472 new_xx(idir) = new_xx(idir) - m_half
473 end do
474
475 ! Convert back to Cartesian coordinates
476 new_xx = this%red_to_cart(new_xx)
477 else
478 new_xx = xx
479 end if
480
482
483 !--------------------------------------------------------------
484 subroutine lattice_vectors_write_info(this, namespace)
485 class(lattice_vectors_t), intent(in) :: this
486 type(namespace_t), intent(in) :: namespace
487
488 integer :: idir, idir2
489
491
492 assert(this%space%is_periodic())
493
494 call messages_print_with_emphasis(msg="Lattice", namespace=namespace)
495
496 write(message(1),'(a,3a,a)') ' Lattice Vectors [', trim(units_abbrev(units_out%length)), ']'
497 do idir = 1, this%space%periodic_dim
498 write(message(1+idir),'(9f12.6)') (units_from_atomic(units_out%length, this%rlattice(idir2, idir)), &
499 idir2 = 1, this%space%dim)
500 end do
501 call messages_info(1+this%space%periodic_dim, namespace=namespace)
502
503 select case (this%space%periodic_dim)
504 case (1)
505 write(message(1),'(a)') ' Cell length ='
506 case (2)
507 write(message(1),'(a)') ' Cell area ='
508 case default
509 write(message(1),'(a)') ' Cell volume ='
510 end select
511 write(message(1),'(a,1x,f18.4,3a,i1.1,a)') &
512 trim(message(1)), units_from_atomic(units_out%length**this%space%periodic_dim, this%rcell_volume), &
513 ' [', trim(units_abbrev(units_out%length**this%space%periodic_dim)), ']'
514 call messages_info(1, namespace=namespace)
515 write(message(1),'(a,3a,a)') ' Reciprocal-Lattice Vectors [', trim(units_abbrev(units_out%length**(-1))), ']'
516 do idir = 1, this%space%periodic_dim
517 write(message(1+idir),'(3f12.6)') (units_from_atomic(unit_one / units_out%length, this%klattice(idir2, idir)), &
518 idir2 = 1, this%space%dim)
519 end do
520 call messages_info(1+this%space%periodic_dim, namespace=namespace)
521
522 if (this%space%dim == 3) then
523 write(message(1),'(a)') ' Cell angles [degree]'
524 write(message(2),'(a, f8.3)') ' alpha = ', this%alpha
525 write(message(3),'(a, f8.3)') ' beta = ', this%beta
526 write(message(4),'(a, f8.3)') ' gamma = ', this%gamma
527 call messages_info(4, namespace=namespace)
528 end if
530 call messages_print_with_emphasis(namespace=namespace)
531
533 end subroutine lattice_vectors_write_info
534
535 !--------------------------------------------------------------
536 character(len=140) function lattice_vectors_short_info(this, unit_length) result(info)
537 class(lattice_vectors_t), intent(in) :: this
538 type(unit_t), intent(in) :: unit_length
540 integer :: idir1, idir2
541
543
544 write(info, '(a,a,a)') 'LatticeVectors [', trim(units_abbrev(unit_length)), '] = '
545
546 do idir1 = 1, this%space%dim
547 write(info, '(a,a)') trim(info), '['
548 do idir2 = 1, this%space%dim
549 write(info, '(a,x,f11.6)') trim(info), units_from_atomic(unit_length, this%rlattice(idir2, idir1))
550 end do
551 write(info, '(a,a)') trim(info), ']'
552 if (idir1 < this%space%dim) then
553 write(info, '(a,a)') trim(info), ', '
554 end if
555 end do
556
558 end function lattice_vectors_short_info
559
560 !--------------------------------------------------------------
561 real(real64) function lattice_vectors_max_length(this) result(length)
562 class(lattice_vectors_t), intent(in) :: this
563
565
566 if (this%space%is_periodic()) then
567 length = maxval(norm2(this%rlattice(:,1:this%space%periodic_dim), dim=1))
568 else
569 length = m_zero
570 end if
571
573 end function lattice_vectors_max_length
574
575 !--------------------------------------------------------------
576 subroutine build_metric_from_angles(this, angles)
577 type(lattice_vectors_t), intent(inout) :: this
578 real(real64), intent(in) :: angles(3)
579
580 real(real64) :: cosang, a2, aa, cc
581 real(real64), parameter :: tol_angle = 1.0e-6_real64
582
584
585 !Converting the angles to LatticeVectors
586 !See 57_iovars/ingeo.F90 in Abinit for details
587 if (abs(angles(1) - angles(2)) < tol_angle .and. abs(angles(2) - angles(3)) < tol_angle .and. &
588 (abs(angles(1) - 90.0_real64) + abs(angles(2) - 90.0_real64) + abs(angles(3) - 90.0_real64)) > tol_angle) then
589
590 cosang = cos(m_pi*angles(1)/180.0_real64)
591 a2 = m_two/m_three*(m_one - cosang)
592 aa = sqrt(a2)
593 cc = sqrt(m_one - a2)
594 this%rlattice_primitive(1, 1) = aa
595 this%rlattice_primitive(2, 1) = m_zero
596 this%rlattice_primitive(3, 1) = cc
597 this%rlattice_primitive(1, 2) = -m_half*aa
598 this%rlattice_primitive(2, 2) = m_half*sqrt(m_three)*aa
599 this%rlattice_primitive(3, 2) = cc
600 this%rlattice_primitive(1, 3) = -m_half*aa
601 this%rlattice_primitive(2, 3) = -m_half*sqrt(m_three)*aa
602 this%rlattice_primitive(3, 3) = cc
603 else
604 this%rlattice_primitive(1, 1) = m_one
605 this%rlattice_primitive(2, 1) = m_zero
606 this%rlattice_primitive(3, 1) = m_zero
607 this%rlattice_primitive(1, 2) = cos(m_pi*angles(3)/180.0_real64)
608 this%rlattice_primitive(2, 2) = sin(m_pi*angles(3)/180.0_real64)
609 this%rlattice_primitive(3, 2) = m_zero
610 this%rlattice_primitive(1, 3) = cos(m_pi*angles(2)/180.0_real64)
611 this%rlattice_primitive(2, 3) = (cos(m_pi*angles(1)/180.0_real64) - &
612 this%rlattice_primitive(1, 2)*this%rlattice_primitive(1, 3))/this%rlattice_primitive(2,2)
613 this%rlattice_primitive(3, 3) = sqrt(m_one - this%rlattice_primitive(1, 3)**2 - this%rlattice_primitive(2, 3)**2)
614 end if
615
616 this%nonorthogonal = any(abs(angles - 90.0_real64) > m_epsilon)
617
619 end subroutine build_metric_from_angles
620
621 !--------------------------------------------------------------
622 subroutine reciprocal_lattice(rv, kv, volume, dim, namespace)
623 real(real64), intent(in) :: rv(:,:)
624 real(real64), intent(out) :: kv(:,:)
625 real(real64), intent(out) :: volume
626 integer, intent(in) :: dim
627 type(namespace_t), optional, intent(in) :: namespace
628
629 integer :: ii
630 real(real64) :: cross(1:3)
631
632 push_sub(reciprocal_lattice)
633
634 select case (dim)
635 case (3)
636 cross(1:3) = dcross_product(rv(1:3, 2), rv(1:3, 3))
637 volume = dot_product(rv(1:3, 1), cross(1:3))
638
639 kv(1:3, 1) = dcross_product(rv(:, 2), rv(:, 3))/volume
640 kv(1:3, 2) = dcross_product(rv(:, 3), rv(:, 1))/volume
641 kv(1:3, 3) = dcross_product(rv(:, 1), rv(:, 2))/volume
642 case (2)
643 volume = rv(1, 1)*rv(2, 2) - rv(2, 1)*rv(1, 2)
644 kv(1, 1) = rv(2, 2)/volume
645 kv(2, 1) = -rv(1, 2)/volume
646 kv(1, 2) = -rv(2, 1)/volume
647 kv(2, 2) = rv(1, 1)/volume
648 case (1)
649 volume = rv(1, 1)
650 kv(1, 1) = m_one / rv(1, 1)
651 case default ! dim > 3
652 message(1) = "Reciprocal lattice for dim > 3 assumes no periodicity."
653 call messages_warning(1, namespace=namespace)
654 volume = m_one
655 kv(:,:) = m_zero
656 do ii = 1, dim
657 kv(ii, ii) = m_one/rv(ii,ii)
658 ! At least initialize the thing
659 volume = volume * norm2(rv(:, ii))
660 end do
661 end select
662
663 if (volume < m_zero) then
664 message(1) = "Your lattice vectors form a left-handed system."
665 call messages_fatal(1, namespace=namespace)
666 end if
667
668 pop_sub(reciprocal_lattice)
669 end subroutine reciprocal_lattice
670
671 ! ---------------------------------------------------------
672 !! The solution is based on the replies found in
673 !! https://stackoverflow.com/questions/29538135/best-way-to-find-all-points-of-lattice-in-sphere
674 !! https://math.stackexchange.com/questions/1230112/maximizing-a-coordinate-of-xt-at-a-x-leq-r2/1230160#1230160
675 function lattice_iterator_constructor(latt, range) result(iter)
676 type(lattice_vectors_t), target, intent(in) :: latt
677 real(real64), intent(in) :: range
678 type(lattice_iterator_t) :: iter
679
680 integer :: ii, jj, idir
681 integer :: n_size(latt%space%periodic_dim)
682 real(real64) :: temp(latt%space%dim)
683
685
686 iter%latt => latt
687
688 ! Determine number of cells
689 iter%n_cells = 1
690 do idir = 1, latt%space%periodic_dim
691 temp = range * matmul(transpose(latt%klattice),latt%klattice(:, idir)) / norm2(latt%klattice(:, idir)) / (m_two * m_pi)
692 n_size(idir) = ceiling(temp(idir))
693 iter%n_cells = iter%n_cells*(2*n_size(idir) + 1)
694 end do
695
696 ! Indexes of the cells
697 safe_allocate(iter%icell(1:latt%space%dim, 1:iter%n_cells))
698
699 do ii = 1, iter%n_cells
700 jj = ii - 1
701 do idir = latt%space%periodic_dim, 1, -1
702 iter%icell(idir, ii) = mod(jj, 2*n_size(idir) + 1) - n_size(idir)
703 if (idir > 1) jj = jj/(2*n_size(idir) + 1)
704 end do
705 iter%icell(latt%space%periodic_dim + 1:latt%space%dim, ii) = 0
706 end do
707
710
711 !--------------------------------------------------------------
712 subroutine lattice_iterator_copy(this, source)
713 class(lattice_iterator_t), intent(out) :: this
714 class(lattice_iterator_t), intent(in) :: source
716 push_sub(lattice_iterator_copy)
717
718 this%n_cells = source%n_cells
719 safe_allocate_source_a(this%icell, source%icell)
720 this%latt => source%latt
721
722 pop_sub(lattice_iterator_copy)
723 end subroutine lattice_iterator_copy
724
727 function lattice_iterator_get(this, ii) result(coord)
728 class(lattice_iterator_t), intent(in) :: this
729 integer, intent(in) :: ii
730 real(real64) :: coord(1:this%latt%space%dim)
731
732 ! No push/pop, this is called too often
733
734 assert(ii <= this%n_cells)
735
736 coord = matmul(this%latt%rlattice, this%icell(:, ii))
737
738 end function lattice_iterator_get
739
740 ! ---------------------------------------------------------
741 subroutine lattice_iterator_finalize(this)
742 type(lattice_iterator_t), intent(inout) :: this
743
745
746 safe_deallocate_a(this%icell)
747 this%n_cells = 0
748 nullify(this%latt)
749
751 end subroutine lattice_iterator_finalize
752
753end module lattice_vectors_oct_m
754
755!! Local Variables:
756!! mode: f90
757!! coding: utf-8
758!! End:
subroutine info()
Definition: em_resp.F90:1096
double acos(double __x) __attribute__((__nothrow__
double sin(double __x) __attribute__((__nothrow__
double sqrt(double __x) __attribute__((__nothrow__
double cos(double __x) __attribute__((__nothrow__
real(real64), parameter, public m_two
Definition: global.F90:189
real(real64), parameter, public m_zero
Definition: global.F90:187
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:185
real(real64), parameter, public m_epsilon
Definition: global.F90:203
real(real64) function, dimension(this%space%dim) lattice_vectors_fold_into_cell(this, xx)
pure real(real64) function, dimension(this%space%dim) lattice_vectors_red_to_cart(this, xx_red)
type(lattice_iterator_t) function lattice_iterator_constructor(latt, range)
subroutine reciprocal_lattice(rv, kv, volume, dim, namespace)
subroutine lattice_vectors_update(this, rlattice)
Updates the lattice vector object from a new set of Cartesian lattice vectors.
pure real(real64) function, dimension(this%space%dim) lattice_vectors_cart_to_red(this, xx_cart)
character(len=140) function lattice_vectors_short_info(this, unit_length)
subroutine angles_from_rlattice_primitive(rlattice_primitive, alpha, beta, gamma)
subroutine lattice_vectors_scale(this, factor)
type(lattice_vectors_t) function lattice_vectors_constructor_from_input(namespace, space, variable_prefix)
real(real64) function, dimension(1:this%latt%space%dim) lattice_iterator_get(this, ii)
real(real64) function lattice_vectors_max_length(this)
subroutine lattice_vectors_copy(this, source)
subroutine lattice_iterator_copy(this, source)
subroutine build_metric_from_angles(this, angles)
subroutine lattice_vectors_write_info(this, namespace)
subroutine lattice_vectors_finalize(this)
type(lattice_vectors_t) function lattice_vectors_constructor_from_rlattice(namespace, space, rlattice)
subroutine lattice_iterator_finalize(this)
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
This module defines the unit system, used for input and output.
The following class implements a lattice iterator. It allows one to loop over all cells that are with...
int true(void)