Octopus
pes.F90
Go to the documentation of this file.
1!! Copyright (C) 2006-2011 M. Marques, U. De Giovannini
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
21module pes_oct_m
23 use box_oct_m
24 use debug_oct_m
26 use global_oct_m
27 use grid_oct_m
29 use ions_oct_m
31 use mesh_oct_m
33 use mpi_oct_m
36 use parser_oct_m
41 use space_oct_m
44
45 implicit none
46
47 private
48
49 public :: &
50 pes_t, &
51 pes_end, &
52 pes_init, &
54 pes_calc, &
55 pes_output, &
56 pes_load, &
58
59 type pes_t
60 private
61 logical, public :: calc_spm
62 type(pes_spm_t) :: spm
63
64 logical, public :: calc_mask
65 type(pes_mask_t) :: mask
66
67 logical, public :: calc_flux
68 type(pes_flux_t) :: flux
69
70 end type pes_t
71
72 integer, parameter :: &
73 PHOTOELECTRON_NONE = 0, &
77
78contains
79
80 ! ---------------------------------------------------------
81 subroutine pes_init(pes, namespace, space, mesh, box, st, save_iter, kpoints, abs_boundaries, ext_partners, max_iter, dt)
82 type(pes_t), intent(out) :: pes
83 type(namespace_t), intent(in) :: namespace
84 class(space_t), intent(in) :: space
85 class(mesh_t), intent(in) :: mesh
86 class(box_t), intent(in) :: box
87 type(states_elec_t), intent(in) :: st
88 integer, intent(in) :: save_iter
89 type(kpoints_t), intent(in) :: kpoints
90 type(absorbing_boundaries_t), intent(in) :: abs_boundaries
91 type(partner_list_t), intent(in) :: ext_partners
92 integer, intent(in) :: max_iter
93 real(real64), intent(in) :: dt
94
95 integer :: photoelectron_flags
96
97 push_sub(pes_init)
98
99 !%Variable PhotoElectronSpectrum
100 !%Type integer
101 !%Default none
102 !%Section Time-Dependent::PhotoElectronSpectrum
103 !%Description
104 !% This variable controls the method used for the calculation of
105 !% the photoelectron spectrum. You can specify more than one value
106 !% by giving them as a sum, for example:
107 !% <tt>PhotoElectronSpectrum = pes_spm + pes_mask</tt>
108 !%Option none 0
109 !% The photoelectron spectrum is not calculated. This is the default.
110 !%Option pes_spm 2
111 !% Store the wavefunctions at specific points in order to
112 !% calculate the photoelectron spectrum at a point far in the box as proposed in
113 !% A. Pohl, P.-G. Reinhard, and E. Suraud, <i>Phys. Rev. Lett.</i> <b>84</b>, 5090 (2000).
114 !%Option pes_mask 4
115 !% Calculate the photo-electron spectrum using the mask method.
116 !% U. De Giovannini, D. Varsano, M. A. L. Marques, H. Appel, E. K. U. Gross, and A. Rubio,
117 !% <i>Phys. Rev. A</i> <b>85</b>, 062515 (2012).
118 !%Option pes_flux 8
119 !% Calculate the photo-electron spectrum using the t-surff technique, <i>i.e.</i>,
120 !% spectra are computed from the electron flux through a surface close to the absorbing
121 !% boundaries of the box. (Experimental.)
122 !% L. Tao and A. Scrinzi, <i>New Journal of Physics</i> <b>14</b>, 013021 (2012).
123 !%End
124
125 call parse_variable(namespace, 'PhotoElectronSpectrum', photoelectron_none, photoelectron_flags)
126 if (.not. varinfo_valid_option('PhotoElectronSpectrum', photoelectron_flags, is_flag = .true.)) then
127 call messages_input_error(namespace, 'PhotoElectronSpectrum')
128 end if
129
130 pes%calc_spm = bitand(photoelectron_flags, photoelectron_spm) /= 0
131 pes%calc_mask = bitand(photoelectron_flags, photoelectron_mask) /= 0
132 pes%calc_flux = bitand(photoelectron_flags, photoelectron_flux) /= 0
133
134 !Header Photoelectron info
135 if (pes%calc_spm .or. pes%calc_mask .or. pes%calc_flux) then
136 call messages_print_with_emphasis(msg='Photoelectron', namespace=namespace)
137 end if
138
139
140 if (pes%calc_spm) call pes_spm_init(pes%spm, namespace, mesh, st, save_iter)
141 if (pes%calc_mask) then
142 call pes_mask_init(pes%mask, namespace, space, mesh, box, st, ext_partners, &
143 abs_boundaries, max_iter,dt)
144 end if
145 if (pes%calc_flux) then
146 call pes_flux_init(pes%flux, namespace, space, mesh, st, &
147 ext_partners, kpoints, abs_boundaries, save_iter, max_iter)
148 end if
149
150
151 !Footer Photoelectron info
152 if (pes%calc_spm .or. pes%calc_mask .or. pes%calc_flux) then
153 call messages_print_with_emphasis(namespace=namespace)
154 end if
156 pop_sub(pes_init)
157 end subroutine pes_init
159
160 ! ---------------------------------------------------------
161 subroutine pes_end(pes)
162 type(pes_t), intent(inout) :: pes
163
164 push_sub(pes_end)
166 if (pes%calc_spm) call pes_spm_end (pes%spm)
167 if (pes%calc_mask) call pes_mask_end(pes%mask)
168 if (pes%calc_flux) call pes_flux_end(pes%flux)
169
170 pop_sub(pes_end)
171 end subroutine pes_end
172
173
174 ! ---------------------------------------------------------
175 subroutine pes_calc(pes, namespace, space, mesh, st, dt, iter, der, kpoints, ext_partners, stopping)
176 type(pes_t), intent(inout) :: pes
177 type(namespace_t), intent(in) :: namespace
178 class(space_t), intent(in) :: space
179 class(mesh_t), intent(in) :: mesh
180 type(states_elec_t), intent(inout) :: st
181 type(derivatives_t), intent(in) :: der
182 real(real64), intent(in) :: dt
183 integer, intent(in) :: iter
184 type(kpoints_t), intent(in) :: kpoints
185 type(partner_list_t), intent(in) :: ext_partners
186 logical, intent(in) :: stopping
187
188 push_sub(pes_calc)
189
190 if (pes%calc_spm) call pes_spm_calc(pes%spm, st, mesh, dt, iter, ext_partners)
191 if (pes%calc_mask) call pes_mask_calc(pes%mask, namespace, space, mesh, st, kpoints, dt, iter)
192 if (pes%calc_flux) then
193 call pes_flux_calc(pes%flux, space, mesh, st, der, ext_partners, kpoints, iter, dt, stopping)
194 end if
195
196 pop_sub(pes_calc)
197 end subroutine pes_calc
198
199
200 ! ---------------------------------------------------------
201 subroutine pes_output(pes, namespace, space, gr, st, iter, outp, dt, ions)
202 type(pes_t), intent(inout) :: pes
203 type(namespace_t), intent(in) :: namespace
204 class(space_t), intent(in) :: space
205 type(grid_t), intent(in) :: gr
206 type(states_elec_t), intent(in) :: st
207 integer, intent(in) :: iter
208 type(output_t), intent(in) :: outp
209 real(real64), intent(in) :: dt
210 type(ions_t), intent(in) :: ions
211
212 push_sub(pes_output)
213
214 if (pes%calc_spm) call pes_spm_output(pes%spm, gr, st, namespace, iter, dt)
215
216 if (pes%calc_mask) call pes_mask_output (pes%mask, gr, st, outp, namespace, space, "td.general/PESM", ions,iter)
217
218 if (pes%calc_flux) call pes_flux_output(pes%flux, gr%box, st, namespace)
219
220 pop_sub(pes_output)
221 end subroutine pes_output
222
223
224 ! ---------------------------------------------------------
225 subroutine pes_dump(pes, namespace, restart, st, mesh, ierr)
226 type(pes_t), intent(in) :: pes
227 type(namespace_t), intent(in) :: namespace
228 type(restart_t), intent(in) :: restart
229 type(states_elec_t), intent(in) :: st
230 class(mesh_t), intent(in) :: mesh
231 integer, intent(out) :: ierr
232
233 push_sub(pes_dump)
234
235 ierr = 0
236
237 if (restart_skip(restart)) then
238 pop_sub(pes_dump)
239 return
240 end if
241
242 if (debug%info) then
243 message(1) = "Debug: Writing PES restart."
244 call messages_info(1, namespace=namespace)
245 end if
246
247 if (pes%calc_mask) then
248 call pes_mask_dump(pes%mask, namespace, restart, st, ierr)
249 end if
250
251 if (pes%calc_flux) then
252 call pes_flux_dump(restart, pes%flux, mesh, st, ierr)
253 end if
255 if (pes%calc_spm) then
256 call pes_spm_dump(restart, pes%spm, st, ierr)
257 end if
258
259 if (debug%info) then
260 message(1) = "Debug: Writing PES restart done."
261 call messages_info(1, namespace=namespace)
262 end if
263
264 pop_sub(pes_dump)
265 end subroutine pes_dump
266
267
268 ! ---------------------------------------------------------
269 subroutine pes_load(pes, namespace, restart, st, ierr)
270 type(pes_t), intent(inout) :: pes
271 type(namespace_t), intent(in) :: namespace
272 type(restart_t), intent(in) :: restart
273 type(states_elec_t), intent(inout) :: st
274 integer, intent(out) :: ierr
275
276 push_sub(pes_load)
277
278 ierr = 0
279
280 if (restart_skip(restart)) then
281 ierr = -1
282 pop_sub(pes_load)
283 return
284 end if
285
286 if (debug%info) then
287 message(1) = "Debug: Reading PES restart."
288 call messages_info(1, namespace=namespace)
289 end if
290
291 if (pes%calc_mask) then
292 call pes_mask_load(pes%mask, namespace, restart, st, ierr)
293 end if
295 if (pes%calc_flux) then
296 call pes_flux_load(restart, pes%flux, st, ierr)
297 end if
298
299 if (pes%calc_spm) then
300 call pes_spm_load(restart, pes%spm, st, ierr)
301 end if
302
303 if (debug%info) then
304 message(1) = "Debug: Reading PES restart done."
305 call messages_info(1, namespace=namespace)
306 end if
307
308 pop_sub(pes_load)
309 end subroutine pes_load
310
311
312 ! ---------------------------------------------------------
313 subroutine pes_init_write(pes, mesh, st, namespace)
314 type(pes_t), intent(in) :: pes
315 class(mesh_t), intent(in) :: mesh
316 type(states_elec_t), intent(in) :: st
317 type(namespace_t), intent(in) :: namespace
319
320 push_sub(pes_init_write)
321
322 if (mpi_grp_is_root(mpi_world)) then
323
324 if (pes%calc_spm) call pes_spm_init_write (pes%spm, mesh, st, namespace)
325
326 end if
327
328 pop_sub(pes_init_write)
329 end subroutine pes_init_write
330
331
332
333end module pes_oct_m
334
335
336!! Local Variables:
337!! mode: f90
338!! coding: utf-8
339!! End:
type(debug_t), save, public debug
Definition: debug.F90:154
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
This module implements the underlying real-space grid.
Definition: grid.F90:117
This module defines classes and functions for interaction partners.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:930
character(len=512), private msg
Definition: messages.F90:165
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 messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:723
logical function mpi_grp_is_root(grp)
Is the current MPI process of grpcomm, root.
Definition: mpi.F90:430
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:266
this module contains the output system
Definition: output_low.F90:115
subroutine, public pes_flux_calc(this, space, mesh, st, der, ext_partners, kpoints, iter, dt, stopping)
Definition: pes_flux.F90:1683
subroutine, public pes_flux_output(this, box, st, namespace)
Definition: pes_flux.F90:3849
subroutine, public pes_flux_load(restart, this, st, ierr)
Definition: pes_flux.F90:4475
subroutine, public pes_flux_dump(restart, this, mesh, st, ierr)
Definition: pes_flux.F90:4379
subroutine, public pes_flux_end(this)
Definition: pes_flux.F90:767
subroutine, public pes_flux_init(this, namespace, space, mesh, st, ext_partners, kpoints, abb, save_iter, max_iter)
Definition: pes_flux.F90:270
subroutine, public pes_mask_calc(mask, namespace, space, mesh, st, kpoints, dt, iter)
Definition: pes_mask.F90:1355
subroutine, public pes_mask_dump(mask, namespace, restart, st, ierr)
Definition: pes_mask.F90:3603
subroutine, public pes_mask_output(mask, gr, st, outp, namespace, space, file, ions, iter)
This routine is the main routine dedicated to the output of PES data.
Definition: pes_mask.F90:3369
subroutine, public pes_mask_init(mask, namespace, space, mesh, box, st, ext_partners, abs_boundaries, max_iter, dt)
Definition: pes_mask.F90:277
subroutine, public pes_mask_load(mask, namespace, restart, st, ierr)
Definition: pes_mask.F90:3698
subroutine, public pes_mask_end(mask)
Definition: pes_mask.F90:869
subroutine, public pes_calc(pes, namespace, space, mesh, st, dt, iter, der, kpoints, ext_partners, stopping)
Definition: pes.F90:269
subroutine, public pes_init(pes, namespace, space, mesh, box, st, save_iter, kpoints, abs_boundaries, ext_partners, max_iter, dt)
Definition: pes.F90:175
subroutine, public pes_output(pes, namespace, space, gr, st, iter, outp, dt, ions)
Definition: pes.F90:295
subroutine, public pes_init_write(pes, mesh, st, namespace)
Definition: pes.F90:407
subroutine, public pes_end(pes)
Definition: pes.F90:255
subroutine, public pes_load(pes, namespace, restart, st, ierr)
Definition: pes.F90:363
integer, parameter photoelectron_spm
Definition: pes.F90:165
integer, parameter photoelectron_flux
Definition: pes.F90:165
integer, parameter photoelectron_mask
Definition: pes.F90:165
subroutine, public pes_dump(pes, namespace, restart, st, mesh, ierr)
Definition: pes.F90:319
subroutine, public pes_spm_dump(restart, this, st, ierr)
Definition: pes_spm.F90:906
subroutine, public pes_spm_init_write(this, mesh, st, namespace)
Definition: pes_spm.F90:853
subroutine, public pes_spm_output(this, mesh, st, namespace, iter, dt)
Definition: pes_spm.F90:544
subroutine, public pes_spm_calc(this, st, mesh, dt, iter, ext_partners)
Definition: pes_spm.F90:449
subroutine, public pes_spm_init(this, namespace, mesh, st, save_iter)
Definition: pes_spm.F90:181
subroutine, public pes_spm_load(restart, this, st, ierr)
Definition: pes_spm.F90:951
subroutine, public pes_spm_end(this)
Definition: pes_spm.F90:430
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
class representing derivatives
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:168
Describes mesh distribution to nodes.
Definition: mesh.F90:186
output handler class
Definition: output_low.F90:163
The states_elec_t class contains all electronic wave functions.
int true(void)