Octopus
linear_response.F90
Go to the documentation of this file.
1!! Copyright (C) 2004 E.S. Kadantsev, M. Marques
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
20#include "global.h"
21
23 use comm_oct_m
24 use debug_oct_m
26 use global_oct_m
27 use, intrinsic :: iso_fortran_env
29 use mesh_oct_m
34 use smear_oct_m
35 use space_oct_m
39
40 implicit none
41
42 private
43
44 public :: &
45 lr_t, &
46 lr_init, &
48 lr_zero, &
49 lr_copy, &
56 lr_alpha_j, &
57 lr_dealloc, &
65
66
67
68 type lr_t
69 ! Components are public by default
70 logical, private :: is_allocated, is_allocated_rho
71
73 real(real64), allocatable :: ddl_rho(:,:)
74 real(real64), allocatable :: ddl_psi(:,:,:,:)
75
77 complex(real64), allocatable :: zdl_rho(:,:)
78 complex(real64), allocatable :: zdl_psi(:,:,:,:)
79
81 complex(real64), allocatable :: dl_j(:,:,:)
82 real(real64), allocatable :: ddl_de(:,:)
83 real(real64), allocatable :: ddl_elf(:,:)
84 complex(real64), allocatable :: zdl_de(:,:)
85 complex(real64), allocatable :: zdl_elf(:,:)
86
87 end type lr_t
88
89contains
90
91 ! ---------------------------------------------------------
92 subroutine lr_init(lr)
93 type(lr_t), intent(out) :: lr
94
95 push_sub(lr_init)
96
97 lr%is_allocated = .false.
98
99 pop_sub(lr_init)
100
101 end subroutine lr_init
102
103
104
105 ! ---------------------------------------------------------
106 subroutine lr_allocate(lr, st, mesh, allocate_rho)
107 type(lr_t), intent(inout) :: lr
108 type(states_elec_t), intent(in) :: st
109 class(mesh_t), intent(in) :: mesh
110 logical, optional, intent(in) :: allocate_rho
111
112 push_sub(lr_allocate)
113
114 lr%is_allocated_rho = .true.
115 if (present(allocate_rho)) lr%is_allocated_rho = allocate_rho
116
117 if (states_are_complex(st)) then
118 safe_allocate(lr%zdl_psi(1:mesh%np_part, 1:st%d%dim, st%st_start:st%st_end, st%d%kpt%start:st%d%kpt%end))
119 if (lr%is_allocated_rho) then
120 safe_allocate(lr%zdl_rho(1:mesh%np, 1:st%d%nspin))
121 end if
122 else
123 safe_allocate(lr%ddl_psi(1:mesh%np_part, 1:st%d%dim, st%st_start:st%st_end, st%d%kpt%start:st%d%kpt%end))
124 if (lr%is_allocated_rho) then
125 safe_allocate(lr%ddl_rho(1:mesh%np, 1:st%d%nspin))
126 end if
127 end if
128
129 lr%is_allocated = .true.
130
131 call lr_zero(lr, st)
132
133 pop_sub(lr_allocate)
134
135 end subroutine lr_allocate
136
137
138
139 ! ---------------------------------------------------------
140 subroutine lr_zero(lr, st)
141 type(lr_t), intent(inout) :: lr
142 type(states_elec_t), intent(in) :: st
143
144 push_sub(lr_zero)
145
146 assert(lr%is_allocated)
147
148 if (states_are_complex(st)) then
149 lr%zdl_psi = m_zero
150 if (lr%is_allocated_rho) lr%zdl_rho = m_zero
151 else
152 lr%ddl_psi = m_zero
153 if (lr%is_allocated_rho) lr%ddl_rho = m_zero
154 end if
155
156 pop_sub(lr_zero)
157
158 end subroutine lr_zero
159
160
161 ! ---------------------------------------------------------
162 subroutine lr_dealloc(lr)
163 type(lr_t), intent(inout) :: lr
164
165 push_sub(lr_dealloc)
167 safe_deallocate_a(lr%ddl_psi)
168 safe_deallocate_a(lr%zdl_psi)
169
170 safe_deallocate_a(lr%ddl_rho)
171 safe_deallocate_a(lr%zdl_rho)
172
173 safe_deallocate_a(lr%dl_j)
174 safe_deallocate_a(lr%ddl_de)
175 safe_deallocate_a(lr%ddl_elf)
176 safe_deallocate_a(lr%zdl_de)
177 safe_deallocate_a(lr%zdl_elf)
179 pop_sub(lr_dealloc)
180 end subroutine lr_dealloc
181
182
183 ! ---------------------------------------------------------
184 subroutine lr_copy(st, mesh, src, dest)
185 type(states_elec_t), intent(in) :: st
186 class(mesh_t), intent(in) :: mesh
187 type(lr_t), intent(in) :: src
188 type(lr_t), intent(inout) :: dest
189
190 integer :: ik, ist
191
192 push_sub(lr_copy)
193
194 if (src%is_allocated_rho .and. dest%is_allocated_rho) then
195 if (states_are_complex(st)) then
196 call lalg_copy(mesh%np, st%d%nspin, src%zdl_rho, dest%zdl_rho)
197 else
198 call lalg_copy(mesh%np, st%d%nspin, src%ddl_rho, dest%ddl_rho)
199 end if
200 else
201 if (dest%is_allocated_rho) then
202 if (states_are_complex(st)) then
203 dest%zdl_rho(:, :) = m_zero
204 else
205 dest%ddl_rho(:, :) = m_zero
206 end if
207 end if
208 end if
209
210 do ik = st%d%kpt%start, st%d%kpt%end
211 do ist = st%st_start, st%st_end
212 if (states_are_complex(st)) then
213 call lalg_copy(mesh%np_part, st%d%dim, src%zdl_psi(:, :, ist, ik), dest%zdl_psi(:, :, ist, ik))
214 else
215 call lalg_copy(mesh%np_part, st%d%dim, src%ddl_psi(:, :, ist, ik), dest%ddl_psi(:, :, ist, ik))
216 end if
217 end do
218 end do
219
220 pop_sub(lr_copy)
221
222 end subroutine lr_copy
223
224
225
226 ! ---------------------------------------------------------
227 logical function lr_is_allocated(this)
228 type(lr_t), intent(in) :: this
229
230 push_sub(lr_is_allocated)
231 lr_is_allocated = this%is_allocated
232
234 end function lr_is_allocated
235
236
237
238 ! ---------------------------------------------------------
239 real(real64) function lr_alpha_j(st, jst, ik)
240 type(states_elec_t), intent(in) :: st
241 integer, intent(in) :: jst
242 integer, intent(in) :: ik
243
244 real(real64) :: dsmear
245
246 push_sub(lr_alpha_j)
247
248 if (st%smear%method == smear_fixed_occ) then
249 lr_alpha_j = st%occ(jst, ik) / st%smear%el_per_state
250 else
251 dsmear = max(1e-14_real64, st%smear%dsmear)
252 lr_alpha_j = max(st%smear%e_fermi + m_three*dsmear - st%eigenval(jst, ik), m_zero)
253 end if
254
255 pop_sub(lr_alpha_j)
256 end function lr_alpha_j
257
258#include "undef.F90"
259#include "real.F90"
260#include "linear_response_inc.F90"
261
262#include "undef.F90"
263#include "complex.F90"
264#include "linear_response_inc.F90"
265
266end module linear_response_oct_m
267
268!! Local Variables:
269!! mode: f90
270!! coding: utf-8
271!! End:
Copies a vector x, to a vector y.
Definition: lalg_basic.F90:185
real(real64), parameter, public m_zero
Definition: global.F90:187
real(real64), parameter, public m_three
Definition: global.F90:190
subroutine, public zlr_orth_response(mesh, st, lr, omega)
subroutine, public dlr_swap_sigma(st, mesh, plus, minus)
subroutine, public dlr_orth_vector(mesh, st, vec, ist, ik, omega, min_proj)
Orthogonalizes vec against all the occupied states. For details on the metallic part,...
subroutine, public lr_copy(st, mesh, src, dest)
subroutine, public lr_zero(lr, st)
subroutine, public zlr_dump_rho(lr, space, mesh, nspin, restart, rho_tag, ierr)
real(real64) function, public lr_alpha_j(st, jst, ik)
subroutine, public dlr_load_rho(dl_rho, space, mesh, nspin, restart, rho_tag, ierr)
subroutine, public dlr_dump_rho(lr, space, mesh, nspin, restart, rho_tag, ierr)
subroutine, public lr_allocate(lr, st, mesh, allocate_rho)
subroutine, public lr_init(lr)
subroutine, public lr_dealloc(lr)
subroutine, public zlr_orth_vector(mesh, st, vec, ist, ik, omega, min_proj)
Orthogonalizes vec against all the occupied states. For details on the metallic part,...
subroutine, public zlr_build_dl_rho(mesh, st, lr, nsigma)
subroutine, public zlr_load_rho(dl_rho, space, mesh, nspin, restart, rho_tag, ierr)
subroutine, public dlr_orth_response(mesh, st, lr, omega)
subroutine, public zlr_swap_sigma(st, mesh, plus, minus)
subroutine, public dlr_build_dl_rho(mesh, st, lr, nsigma)
logical function, public lr_is_allocated(this)
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
integer, parameter, public smear_fixed_occ
Definition: smear.F90:171
pure logical function, public states_are_complex(st)
Describes mesh distribution to nodes.
Definition: mesh.F90:186
The states_elec_t class contains all electronic wave functions.
int true(void)