Octopus
x_slater.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!! Copyright (C) 2023 N. Tancogne-Dejean
3!!
4!! This program is free software; you can redistribute it and/or modify
5!! it under the terms of the GNU General Public License as published by
6!! the Free Software Foundation; either version 2, or (at your option)
7!! any later version.
8!!
9!! This program is distributed in the hope that it will be useful,
10!! but WITHOUT ANY WARRANTY; without even the implied warranty of
11!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12!! GNU General Public License for more details.
13!!
14!! You should have received a copy of the GNU General Public License
15!! along with this program; if not, write to the Free Software
16!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17!! 02110-1301, USA.
18!!
19
20#include "global.h"
21
22module x_slater_oct_m
23 use comm_oct_m
24 use debug_oct_m
27 use global_oct_m
28 use grid_oct_m
29 use, intrinsic :: iso_fortran_env
36 use space_oct_m
40
41 implicit none
42
43 private
44 public :: &
46
47contains
48
49 ! ---------------------------------------------------------
51 subroutine x_slater_calc (namespace, gr, space, exxop, st, kpoints, ex, vxc)
52 type(namespace_t), intent(in) :: namespace
53 type(grid_t), intent(in) :: gr
54 class(space_t), intent(in) :: space
55 type(exchange_operator_t), intent(in) :: exxop
56 type(states_elec_t), intent(inout) :: st
57 type(kpoints_t), intent(in) :: kpoints
58 real(real64), intent(inout) :: ex
59 real(real64), optional, intent(inout) :: vxc(:,:)
60
61 push_sub(x_slater_calc)
62
63 if (states_are_real(st)) then
64 call dslater_calc(namespace, gr, space, exxop, st, kpoints, ex, vxc)
65 else
66 call zslater_calc(namespace, gr, space, exxop, st, kpoints, ex, vxc)
67 end if
68
69 pop_sub(x_slater_calc)
70 end subroutine x_slater_calc
71
72 ! ---------------------------------------------------------
85 subroutine get_rotation_matrix(dens, alpha, betar, betai)
86 real(real64), intent(in) :: dens(:)
87 real(real64), intent(out) :: alpha, betar, betai
88
89 real(real64), parameter :: tiny = 1e-12_real64
90 real(real64) :: mz, mm
91
92 mz = dens(1) - dens(2)
93
94 mm = sqrt(mz**2 + m_four*(dens(3)**2 + dens(4)**2))
95
96 !Fully spin unpolarized system
97 if (mm < tiny) then
98 alpha = m_one
99 betar = m_zero
100 betai = m_zero
101 return
102 end if
103
104 alpha = sqrt((mm + abs(mz))/(m_two * mm))
105 !We find the absolute values of real and imaginary parts of beta
106 betar = m_two * dens(3) / sqrt(m_two * mm * (mm + abs(mz)))
107 betai = m_two * dens(4) / sqrt(m_two * mm * (mm + abs(mz)))
108
109 if (mz < m_zero) then
110 betar = -betar
111 betai = -betai
112 end if
113
114 end subroutine get_rotation_matrix
116 ! ---------------------------------------------------------
121 subroutine rotate_to_local(mat, alpha, betar, betai, alpha2, beta2, rot_mat)
122 real(real64), intent(in) :: mat(:)
123 real(real64), intent(in) :: alpha, betar, betai, alpha2, beta2
124 real(real64), intent(out) :: rot_mat(:)
125
126 complex(real64) :: cross
127
128 rot_mat(1) = alpha2 * mat(1) + beta2 * mat(2) + m_two * alpha * (betar * mat(3) + betai * mat(4))
129 rot_mat(2) = alpha2 * mat(2) + beta2 * mat(1) - m_two * alpha * (betar * mat(3) + betai * mat(4))
130 cross = (cmplx(betar, betai, real64) )**2 * cmplx(mat(3), -mat(4), real64)
131 rot_mat(3) = alpha2 * mat(3) + alpha * betar * (mat(2)-mat(1)) - real(cross)
132 rot_mat(4) = alpha2 * mat(4) + alpha * betai * (mat(2)-mat(1)) - aimag(cross)
133
134 end subroutine rotate_to_local
135
136 ! ---------------------------------------------------------
142 subroutine rotate_to_global(mat, alpha, betar, betai, alpha2, beta2, rot_mat)
143 real(real64), intent(in) :: mat(:)
144 real(real64), intent(in) :: alpha, betar, betai, alpha2, beta2
145 real(real64), intent(out) :: rot_mat(:)
146
147 complex(real64) :: cross
148
149 rot_mat(1) = alpha2 * mat(1) + beta2 * mat(2) - m_two * alpha * (betar * mat(3) + betai * mat(4))
150 rot_mat(2) = alpha2 * mat(2) + beta2 * mat(1) + m_two * alpha * (betar * mat(3) + betai * mat(4))
151 cross = (cmplx(betar, betai, real64) )**2 * cmplx(mat(3), -mat(4), real64)
152 rot_mat(3) = alpha2 * mat(3) - alpha * betar * (mat(2)-mat(1)) - real(cross)
153 rot_mat(4) = alpha2 * mat(4) - alpha * betai * (mat(2)-mat(1)) - aimag(cross)
154
155 end subroutine rotate_to_global
156
157
158#include "undef.F90"
159#include "real.F90"
160#include "x_slater_inc.F90"
161
162#include "undef.F90"
163#include "complex.F90"
164#include "x_slater_inc.F90"
165
166end module x_slater_oct_m
167
168!! Local Variables:
169!! mode: f90
170!! coding: utf-8
171!! End:
double sqrt(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_four
Definition: global.F90:191
real(real64), parameter, public m_one
Definition: global.F90:188
This module implements the underlying real-space grid.
Definition: grid.F90:117
This module defines various routines, operating on mesh functions.
pure logical function, public states_are_real(st)
This module handles spin dimensions of the states and the k-point distribution.
subroutine rotate_to_local(mat, alpha, betar, betai, alpha2, beta2, rot_mat)
Rotation to the local frame Given a matrix in spin space, this routine rotates according to the rotat...
Definition: x_slater.F90:215
subroutine rotate_to_global(mat, alpha, betar, betai, alpha2, beta2, rot_mat)
Rotation to the global (Cartesian) frame.
Definition: x_slater.F90:236
subroutine zslater_calc(namespace, gr, space, exxop, st, kpoints, ex, vxc)
Definition: x_slater.F90:580
subroutine, public x_slater_calc(namespace, gr, space, exxop, st, kpoints, ex, vxc)
Interface to X(slater_calc)
Definition: x_slater.F90:145
subroutine dslater_calc(namespace, gr, space, exxop, st, kpoints, ex, vxc)
Definition: x_slater.F90:320
subroutine get_rotation_matrix(dens, alpha, betar, betai)
This routine get the SU(2) matrix that diagonalize the spin-density matrix.
Definition: x_slater.F90:179