Octopus
elec_matrix_elements.F90
Go to the documentation of this file.
1!! Copyright (C) 2023 F. Troisi
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 batch_oct_m
24 use debug_oct_m
28 use global_oct_m
29 use grid_oct_m
32 use ions_oct_m
33 use, intrinsic :: iso_fortran_env
35 use lda_u_oct_m
37 use mesh_oct_m
41 use mpi_oct_m
44 use phase_oct_m
51 use space_oct_m
56 use xc_oct_m
57
58 implicit none
59
60 private
62
64 type(grid_t), pointer :: grid
65 type(space_t), pointer :: space
66 type(states_elec_t), pointer :: states
67
68 integer :: st_start
69 integer :: st_end
70 contains
71 procedure :: init => elec_matrix_elements_init
72 ! Momentum
73 procedure :: momentum_me => elec_momentum_me
74 ! Angular momentum
75 procedure :: angular_momentum_me => elec_angular_momentum_me
76 ! Dipole
77 procedure :: delec_dipole_me, zelec_dipole_me
78 generic :: dipole_me => delec_dipole_me, zelec_dipole_me
79 ! KS multipoles
80 procedure :: delec_ks_multipoles_3d_me, zelec_ks_multipoles_3d_me
81 generic :: ks_multipoles_3d => delec_ks_multipoles_3d_me, zelec_ks_multipoles_3d_me
82 procedure :: delec_ks_multipoles_2d_me, zelec_ks_multipoles_2d_me
83 generic :: ks_multipoles_2d => delec_ks_multipoles_2d_me, zelec_ks_multipoles_2d_me
84 procedure :: delec_ks_multipoles_1d_me, zelec_ks_multipoles_1d_me
85 generic :: ks_multipoles_1d => delec_ks_multipoles_1d_me, zelec_ks_multipoles_1d_me
86 ! One Body
87 procedure :: delec_one_body_me, zelec_one_body_me
88 generic :: one_body_me => delec_one_body_me, zelec_one_body_me
89 ! Two Body
90 procedure :: delec_two_body_me, zelec_two_body_me
91 generic :: two_body_me => delec_two_body_me, zelec_two_body_me
94
95contains
96
98 subroutine elec_matrix_elements_init(this, grid, space, states, st_start, st_end)
99 class(elec_matrix_elements_t), intent(inout) :: this
100 type(grid_t), target, intent(in) :: grid
101 class(space_t), target, intent(in) :: space
102 type(states_elec_t), target, intent(in) :: states
103 integer, optional, intent(in) :: st_start, st_end
104
106
107 ! Associate pointers
108 this%grid => grid
109 this%space => space
110 this%states => states
111
112 this%st_start = optional_default(st_start, states%st_start)
113 this%st_end = optional_default(st_end, states%st_end)
116 end subroutine elec_matrix_elements_init
117
118 subroutine elec_matrix_elements_finalize(this)
119 type(elec_matrix_elements_t), intent(inout) :: this
121
123 end subroutine elec_matrix_elements_finalize
124
125 ! -----------------------------------------------------------------------------
126 subroutine elec_momentum_me(this, kpoints, momentum)
127 class(elec_matrix_elements_t), intent(in) :: this
128 type(kpoints_t), intent(in) :: kpoints
129 real(real64), intent(out) :: momentum(:,:,:)
130
131 push_sub(elec_momentum_me)
132
133 if (states_are_real(this%states)) then
134 call delec_momentum_me(this, kpoints, momentum)
135 else
136 call zelec_momentum_me(this, kpoints, momentum)
137 end if
138 pop_sub(elec_momentum_me)
139 end subroutine elec_momentum_me
140
141 ! -----------------------------------------------------------------------------
142 subroutine elec_angular_momentum_me(this, ll, l2)
143 class(elec_matrix_elements_t), intent(in) :: this
144 real(real64), contiguous, intent(out) :: ll(:, :, :)
145 real(real64), contiguous, optional, intent(out) :: l2(:, :)
146
148
149 if (states_are_real(this%states)) then
150 call delec_angular_momentum_me(this, ll, l2)
151 else
152 call zelec_angular_momentum_me(this, ll, l2)
153 end if
155 end subroutine elec_angular_momentum_me
157#include "undef.F90"
158#include "complex.F90"
159#include "elec_matrix_elements_inc.F90"
160
161#include "undef.F90"
162#include "real.F90"
163#include "elec_matrix_elements_inc.F90"
167!! Local Variables:
168!! mode: f90
169!! coding: utf-8
170!! End:
This module implements batches of mesh functions.
Definition: batch.F90:133
Module implementing boundary conditions in Octopus.
Definition: boundaries.F90:122
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
subroutine elec_matrix_elements_finalize(this)
subroutine zelec_momentum_me(this, kpoints, momentum)
The routine calculates the expectation value of the momentum operator.
subroutine delec_momentum_me(this, kpoints, momentum)
The routine calculates the expectation value of the momentum operator.
subroutine delec_angular_momentum_me(this, ll, l2)
It calculates the expectation value of the angular momentum of the states. If l2 is passed,...
subroutine zelec_angular_momentum_me(this, ll, l2)
It calculates the expectation value of the angular momentum of the states. If l2 is passed,...
subroutine elec_angular_momentum_me(this, ll, l2)
subroutine elec_momentum_me(this, kpoints, momentum)
subroutine elec_matrix_elements_init(this, grid, space, states, st_start, st_end)
Initialize the elec_matrix_elements object.
This module implements the underlying real-space grid.
Definition: grid.F90:117
This module defines functions over batches of mesh functions.
Definition: mesh_batch.F90:116
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
This module contains some common usage patterns of MPI routines.
Definition: mpi_lib.F90:115
pure logical function, public states_are_real(st)
This module handles spin dimensions of the states and the k-point distribution.
Definition: xc.F90:114