Octopus
photon_mode_mf.F90
Go to the documentation of this file.
1!! Copyright (C) 2017 Johannes Flick
2!! Copyright (C) 2021 Davis Welakuh (merged)
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!! v_mxc = - omega*q_alpha*(lambda*r) + (lambda*dipole)*(lambda*r)
20!! see PRL 115, 093001 (2015) and PRL 121, 113002 (2018)
21#include "global.h"
22
24 use debug_oct_m
25 use global_oct_m
26 use grid_oct_m
28 use io_oct_m
29 use ions_oct_m
34 use ps_oct_m
38
39
40 implicit none
41
42 private
43
44 public :: &
45 mf_t, &
46 mf_init, &
47 mf_end, &
48 mf_calc, &
51 type mf_t
52 real(real64), pointer :: vmf(:)
53 real(real64), allocatable :: dipole(:)
54 real(real64), allocatable :: dipole_former(:)
55 complex(real64), allocatable :: integral(:)
56 real(real64), allocatable :: pt_q(:), pt_p(:)
57 real(real64), allocatable :: pt_q_former(:)
58 real(real64) :: time_former
59 real(real64), pointer :: fmf(:) !! meanfield force
60 logical :: has_restart
61
62 end type mf_t
63
64contains
65
66 subroutine mf_init(this, gr, st, ions, pt_mode)
67 type(mf_t), intent(out) :: this
68 type(grid_t), intent(in) :: gr
69 type(states_elec_t), intent(in) :: st
70 type(ions_t), intent(in) :: ions
71 type(photon_mode_t), intent(in) :: pt_mode
72
73 integer :: ions_dim
74
75 push_sub(mf_init)
76
77 ions_dim = gr%box%dim
78
79 safe_allocate(this%vmf(1:gr%np))
80 safe_allocate(this%dipole(1:ions_dim))
81 safe_allocate(this%dipole_former(1:ions_dim))
82
83 safe_allocate(this%integral(1:pt_mode%nmodes))
84 safe_allocate(this%pt_q(1:pt_mode%nmodes))
85 safe_allocate(this%pt_p(1:pt_mode%nmodes))
86 safe_allocate(this%pt_q_former(1:pt_mode%nmodes))
87 safe_allocate(this%fmf(1:ions_dim))
88
89 this%vmf = m_zero
90 this%has_restart = .false.
91
92
93 this%dipole = - ions%dipole() - st%dipole(gr)
94
95 this%dipole_former = m_zero
96 this%integral = m_zero
97
98 if (pt_mode%has_q0_p0) then
99 this%pt_q(1:pt_mode%nmodes) = pt_mode%pt_coord_q0(1:pt_mode%nmodes)
100 this%pt_p(1:pt_mode%nmodes) = pt_mode%pt_momen_p0(1:pt_mode%nmodes)
101 else
102 this%pt_q = m_zero
103 this%pt_p = m_zero
104 end if
105
106 this%pt_q_former = m_zero
107 this%time_former = m_zero
108
109 this%fmf(1:ions_dim) = m_zero
110
111 pop_sub(mf_init)
112 end subroutine mf_init
113
114!------------------------------------------
115
116 subroutine mf_end(this)
117 type(mf_t), intent(inout) :: this
118
119 push_sub(mf_end)
120
121 safe_deallocate_p(this%vmf)
122 safe_deallocate_a(this%dipole_former)
123 safe_deallocate_a(this%dipole)
124
125 safe_deallocate_a(this%integral)
126 safe_deallocate_a(this%pt_p)
127 safe_deallocate_a(this%pt_q)
128 safe_deallocate_a(this%pt_q_former)
129
130 safe_deallocate_p(this%fmf)
131
132 pop_sub(mf_end)
133 end subroutine mf_end
134
135!------------------------------------------
136
137 subroutine mf_calc(this, gr, st, ions, pt_mode, time)
138 type(mf_t), intent(inout) :: this
139 type(grid_t), intent(inout) :: gr
140 type(states_elec_t), intent(inout) :: st
141 type(ions_t), intent(in) :: ions
142 type(photon_mode_t), intent(inout) :: pt_mode
143 real(real64), intent(in) :: time
145 real(real64) :: lambda_pol_dipole, lambda_pol_dipole_former
146 complex(real64) :: integrand
147 real(real64) :: q0, p0
148 integer :: ii, jj, ions_dim
149 logical, save :: first = .true.
151 push_sub(mf_calc)
153 ions_dim = gr%box%dim
154
155 if (.not. (first .and. this%has_restart)) then
156 do ii = 1, pt_mode%nmodes
157 this%pt_q_former(ii) = this%pt_q(ii)
158 end do
160 this%dipole_former = this%dipole
161 this%dipole = - ions%dipole() - st%dipole(gr)
162
163 end if
164
165 this%vmf(1:gr%np) = m_zero
166 this%fmf(1:ions_dim) = m_zero
167
168 do ii = 1, pt_mode%nmodes
169 lambda_pol_dipole = m_zero
170 lambda_pol_dipole_former = m_zero
171 do jj = 1, ions_dim
172 lambda_pol_dipole = lambda_pol_dipole + &
173 pt_mode%lambda(ii)*pt_mode%pol(ii, jj)*this%dipole(jj)
174 lambda_pol_dipole_former = lambda_pol_dipole_former + &
175 pt_mode%lambda(ii)*pt_mode%pol(ii, jj)*this%dipole_former(jj)
176 end do
177
178 if (.not.(first .and. this%has_restart)) then
179 if (abs(time) <= m_epsilon) then
180 this%integral(ii) = m_zero
181 else if (abs(this%integral(ii)) <= m_epsilon) then
182 this%integral(ii) = m_half*lambda_pol_dipole*exp(-m_zi*pt_mode%omega(ii)*(this%time_former))
183 else
184 this%integral(ii) = this%integral(ii) + lambda_pol_dipole*exp(-m_zi*pt_mode%omega(ii)*(this%time_former))
185 end if
186 end if
187
188 integrand = -(this%integral(ii)*exp(m_zi*pt_mode%omega(ii)*(time)*pt_mode%mu))* &
189 (time*pt_mode%mu - this%time_former)
190
191 this%pt_q(ii) = aimag(integrand)
192 this%pt_p(ii) = pt_mode%omega(ii)*real(integrand)
193
194 if (pt_mode%has_q0_p0) then
195 q0 = pt_mode%pt_coord_q0(ii)*cos(pt_mode%omega(ii)*(time)*pt_mode%mu)
196 q0 = q0 + pt_mode%pt_momen_p0(ii)/pt_mode%omega(ii)*sin(pt_mode%omega(ii)*(time)*pt_mode%mu)
197 this%pt_q(ii) = this%pt_q(ii) + q0
198
199 p0 = -pt_mode%pt_coord_q0(ii)*pt_mode%omega(ii)*sin(pt_mode%omega(ii)*(time)*pt_mode%mu)
200 p0 = p0 + pt_mode%pt_momen_p0(ii)*cos(pt_mode%omega(ii)*(time)*pt_mode%mu)
201 this%pt_p(ii) = this%pt_p(ii) + p0
202 end if
203
204 if (.not.allocated(pt_mode%pol_dipole)) then
205 call photon_mode_compute_dipoles(pt_mode, gr)
206 end if
207
208 ! we need the negative sign due to the electric pol_dipole
209 this%vmf(1:gr%np) = this%vmf(1:gr%np) - &
210 m_half*((lambda_pol_dipole + lambda_pol_dipole_former) + pt_mode%omega(ii)* &
211 (this%pt_q(ii) + this%pt_q_former(ii)))*(pt_mode%lambda(ii)*pt_mode%pol_dipole(1:gr%np,ii))
212 do jj = 1, ions_dim
213 this%fmf(jj) = this%fmf(jj) - pt_mode%omega(ii)*pt_mode%lambda(ii)* &
214 pt_mode%pol(ii, jj)*(this%pt_q(ii) + lambda_pol_dipole/pt_mode%omega(ii)) !minus?
215 end do
216 end do
217
218 if (first .and. this%has_restart) then
219 first = .false.
220 end if
221
222 this%time_former = time*pt_mode%mu
223
224 pop_sub(mf_calc)
225 end subroutine mf_calc
226
227! ---------------------------------------------------------
228
229 subroutine mf_photons_dump(restart, this, gr, dt, pt_mode, ierr)
230 type(restart_t), intent(in) :: restart
231 type(mf_t), intent(in) :: this
232 type(grid_t), intent(in) :: gr
233 real(real64), intent(in) :: dt
234 type(photon_mode_t), intent(in) :: pt_mode
235 integer, intent(out) :: ierr
236
237 character(len=80), allocatable :: lines(:)
238 integer :: iunit, err, jj, ions_dim
239
240 push_sub(mf_photons_dump)
241 ions_dim = gr%box%dim
242
243 safe_allocate(lines(1:2*ions_dim + 4))
244
245 ierr = 0
246
247 iunit = restart_open(restart, 'photon_mf')
248 write(lines(1), '(a10,2x,es19.12)') 'pt_integral_real', real(this%integral(1))
249 write(lines(2), '(a10,2x,es19.12)') 'pt_integral_aimag', aimag(this%integral(1))
250 write(lines(3), '(a10,2x,es19.12)') 'pt_q_former', this%pt_q_former(1)
251 write(lines(4), '(a10,2x,es19.12)') 'pt_time_former', this%time_former - dt*pt_mode%mu
252 do jj = 1, ions_dim
253 write(lines(4 + jj), '(a10,2x,es19.12)') 'dipole', this%dipole(jj)
254 end do
255 do jj = 1, ions_dim
256 write(lines(4 + ions_dim + jj), '(a10,2x,es19.12)') 'dipole_former', this%dipole_former(jj)
257 end do
258 call restart_write(restart, iunit, lines, 2*ions_dim + 4, err)
259 if (err /= 0) ierr = ierr + 1
260 call restart_close(restart, iunit)
261
262 safe_deallocate_a(lines)
263
264 pop_sub(mf_photons_dump)
265 end subroutine mf_photons_dump
266
267! ---------------------------------------------------------
268
269 subroutine mf_photons_load(restart, this, gr, ierr)
270 type(restart_t), intent(in) :: restart
271 type(mf_t), intent(inout) :: this
272 type(grid_t), intent(in) :: gr
273 integer, intent(out) :: ierr
274
275 integer :: err, iunit, jj, ions_dim
276 character(len=128), allocatable :: lines(:)
277 character(len=7) :: dummy
278 real(real64), allocatable :: rr(:)
279
280 push_sub(mf_photons_load)
281
282 ierr = 0
283 ions_dim = gr%box%dim
284
285 if (restart_skip(restart)) then
286 ierr = -1
287 pop_sub(mf_photons_load)
288 return
289 end if
290
291 if (debug%info) then
292 message(1) = "Debug: Reading Photons restart."
293 call messages_info(1, namespace=restart%namespace)
294 end if
295
296 safe_allocate(rr(1:2*ions_dim + 4))
297 safe_allocate(lines(1:2*ions_dim + 4))
298 iunit = restart_open(restart, 'photon_mf')
299 call restart_read(restart, iunit, lines, 2*ions_dim + 4, err)
300 if (err /= 0) then
301 ierr = ierr + 1
302 else
303 do jj = 1, 2*ions_dim + 4
304 read(lines(jj),'(a10,2x,es19.12)') dummy, rr(jj)
305 end do
306
307 this%integral(1) = rr(1) + m_zi*rr(2)
308 this%pt_q_former(1) = rr(3)
309 this%time_former = rr(4)
310 do jj = 1, ions_dim
311 this%dipole(jj) = rr(jj + 4)
312 end do
313 do jj = 1, ions_dim
314 this%dipole_former(jj) = rr(jj + ions_dim + 4)
315 end do
316 this%has_restart = .true.
317
318 end if
319 call restart_close(restart, iunit)
320
321 if (debug%info) then
322 message(1) = "Debug: Reading Photons restart done."
323 call messages_info(1, namespace=restart%namespace)
324 end if
325
326 safe_deallocate_a(rr)
327 safe_deallocate_a(lines)
328
329 pop_sub(mf_photons_load)
330 end subroutine mf_photons_load
331
332! ---------------------------------------------------------
333
334end module photon_mode_mf_oct_m
335
336!! Local Variables:
337!! mode: f90
338!! coding: utf-8
339!! End:
double exp(double __x) __attribute__((__nothrow__
double sin(double __x) __attribute__((__nothrow__
double cos(double __x) __attribute__((__nothrow__
type(debug_t), save, public debug
Definition: debug.F90:154
real(real64), parameter, public m_zero
Definition: global.F90:187
complex(real64), parameter, public m_zi
Definition: global.F90:201
real(real64), parameter, public m_epsilon
Definition: global.F90:203
real(real64), parameter, public m_half
Definition: global.F90:193
This module implements the underlying real-space grid.
Definition: grid.F90:117
Definition: io.F90:114
This module defines various routines, operating on mesh functions.
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
Definition: messages.F90:624
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public mf_photons_load(restart, this, gr, ierr)
subroutine, public mf_init(this, gr, st, ions, pt_mode)
subroutine, public mf_photons_dump(restart, this, gr, dt, pt_mode, ierr)
subroutine, public mf_end(this)
subroutine, public mf_calc(this, gr, st, ions, pt_mode, time)
subroutine, public photon_mode_compute_dipoles(this, mesh)
Computes the polarization dipole.
Definition: ps.F90:114
subroutine, public restart_read(restart, iunit, lines, nlines, ierr)
Definition: restart.F90:933
subroutine, public restart_close(restart, iunit)
Close a file previously opened with restart_open.
Definition: restart.F90:950
subroutine, public restart_write(restart, iunit, lines, nlines, ierr)
Definition: restart.F90:906
logical pure function, public restart_skip(restart)
Returns true if the restart information should neither be read nor written. This might happen because...
Definition: restart.F90:967
integer function, public restart_open(restart, filename, status, position, silent)
Open file "filename" found inside the current restart directory. Depending on the type of restart,...
Definition: restart.F90:859
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:168
int true(void)