Octopus
output_me.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 F. Troisi
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 output_me_oct_m
23 use debug_oct_m
26 use global_oct_m
27 use grid_oct_m
30 use io_oct_m
31 use ions_oct_m
33 use mpi_oct_m
38 use parser_oct_m
41 use space_oct_m
45 use unit_oct_m
47 use utils_oct_m
49 use xc_oct_m
50
51 implicit none
52
53 private
54 public :: &
57
58contains
59
60 ! ---------------------------------------------------------
61 subroutine output_me_init(this, namespace, space, st, gr, nst)
62 type(output_me_t), intent(out) :: this
63 type(namespace_t), intent(in) :: namespace
64 class(space_t), intent(in) :: space
65 type(states_elec_t), intent(in) :: st
66 type(grid_t), intent(in) :: gr
67 integer, intent(in) :: nst
68
69 integer(int64) :: how(0:MAX_OUTPUT_TYPES)
70 integer :: output_interval(0:MAX_OUTPUT_TYPES)
71
72 push_sub(output_me_init)
73
74 !%Variable OutputMatrixElements
75 !%Type block
76 !%Default none
77 !%Section Output
78 !%Description
79 !% Specifies what matrix elements to print.
80 !% Enabled only if <tt>Output</tt> block includes <tt>matrix_elements</tt>.
81 !% The output files go into the <tt>static</tt> directory, except when
82 !% running a time-dependent simulation, when the directory <tt>td.XXXXXXX</tt> is used.
83 !%
84 !% Example:
85 !% <br><br><tt>%OutputMatrixElements
86 !% <br>&nbsp;&nbsp;momentum
87 !% <br>&nbsp;&nbsp;ks_multipoles
88 !% <br>%<br></tt>
89 !%
90 !% It is possible to specify only compute the matrix elements for some of the states
91 !% using the variables <tt>OutptMEStart</tt> and <tt>OutputMEEnd</tt>.
92 !%Option momentum 1
93 !% Momentum. Filename: <tt>ks_me_momentum</tt>.
94 !%Option ang_momentum 2
95 !% Dimensionless angular momentum <math>\vec{r} \times \vec{k}</math>. Filename: <tt>ks_me_angular_momentum</tt>.
96 !%Option one_body 3
97 !% <math>\left< i \left| \hat{T} + V_{ext} \right| j \right></math>. Not available with states parallelization.
98 !%Option two_body 4
99 !% <math>\left< ij \left| \frac{1}{\left|\vec{r}_1-\vec{r}_2\right|} \right| kl \right></math>.
100 !% Not available with states parallelization.
101 !% Not available with states parallelization. For periodic system, this is not available for k-point parallelization neither.
102 !%Option two_body_exc_k 5
103 !% <math>\left< n1-k1, n2-k2 \left| \frac{1}{\left|\vec{r}_1-\vec{r}_2\right|} \right| n2-k1 n1-k2 \right></math>.
104 !% Not available with states parallelization. For periodic system, this is not available for k-point parallelization neither.
105 !%Option ks_multipoles 6
106 !% See <tt>OutputMEMultipoles</tt>. Not available with states parallelization.
107 !%Option dipole 7
108 !% Prints the dipole matrix elements. Not available with states parallelization.
109 !% For periodic systems, the intraband terms (dipole matrix elements between degenerated states)
110 !% are set to zero, and only the absolute value of the dipole matrix element is printed.
111 !% Not yet supported for spinors.
112 !%End
113
114 this%what = .false.
115 call io_function_read_what_how_when(namespace, space, this%what, how, output_interval, &
116 'OutputMatrixElements')
117
118 if (st%parallel_in_states) then
119 if (this%what(option__outputmatrixelements__two_body)) then
120 call messages_not_implemented("OutputMatrixElements=two_body is not implemented in states parallelization.", &
121 namespace=namespace)
122 end if
123 if (this%what(option__outputmatrixelements__dipole)) then
124 call messages_not_implemented("OutputMatrixElements=dipole is not implemented in states parallelization.", &
125 namespace=namespace)
126 end if
127 end if
128
129 if (space%dim /= 2 .and. space%dim /= 3) this%what(option__outputmatrixelements__ang_momentum) = .false.
130
131 if (this%what(option__outputmatrixelements__ks_multipoles)) then
132 !%Variable OutputMEMultipoles
133 !%Type integer
134 !%Default 1
135 !%Section Output
136 !%Description
137 !% This variable decides which multipole moments are printed out for
138 !% <tt>OutputMatrixElements = ks_multipoles</tt>:
139 !%
140 !% In 3D, if, for example, <tt>OutputMEMultipoles = 1</tt>, then the program will print three
141 !% files, <tt>ks_me_multipoles.x</tt> (<tt>x</tt>=1,2,3), containing
142 !% respectively the (1,-1), (1,0) and (1,1) multipole matrix elements
143 !% between Kohn-Sham states.
144 !%
145 !% In 2D, this variable is ignored: it will always print two files,
146 !% <tt>ks_me_multipoles.i</tt> (<tt>i</tt>=1,2), containing the <math>x</math> and
147 !% <math>y</math> dipole matrix elements.
148 !%
149 !% In 1D, if, for example, <tt>OutputMEMultipoles = 2</tt>, the program will print two files, containing the
150 !% <math>x</math> and <math>x^2</math> matrix elements between Kohn-Sham states.
151 !%End
152 call parse_variable(namespace, 'OutputMEMultipoles', 1, this%ks_multipoles)
153 end if
155 !%Variable OutputMEStart
156 !%Type integer
157 !%Default 1
158 !%Section Output
159 !%Description
160 !% Specifies the state/band index for starting to compute the matrix element.
161 !% So far, this is only used for dipole matrix elements.
162 !%End
163 call parse_variable(namespace, 'OutputMEStart', 1, this%st_start)
164 assert(this%st_start > 0 .and. this%st_start <= nst)
165
166 !%Variable OutputMEEnd
167 !%Type integer
168 !%Default 1
169 !%Section Output
170 !%Description
171 !% Specifies the highest state/band index used to compute the matrix element.
172 !% So far, this is only used for dipole matrix elements.
173 !%End
174 call parse_variable(namespace, 'OutputMEEnd', nst, this%st_end)
175 assert(this%st_end > 0 .and. this%st_end <= nst)
176 assert(this%st_start <= this%st_end)
177 this%nst = this%st_end - this%st_start +1
178
179 pop_sub(output_me_init)
180 end subroutine output_me_init
181
182
183 ! ---------------------------------------------------------
184 subroutine output_me(this, namespace, space, dir, st, gr, ions, hm)
185 type(output_me_t), intent(in) :: this
186 type(namespace_t), intent(in) :: namespace
187 class(space_t), intent(in) :: space
188 character(len=*), intent(in) :: dir
189 type(states_elec_t), intent(inout) :: st
190 type(grid_t), intent(in) :: gr
191 type(ions_t), intent(in) :: ions
192 type(hamiltonian_elec_t), intent(inout) :: hm
193
194 integer :: id, ll, iunit
195 character(len=256) :: fname
196 real(real64), allocatable :: doneint(:), dtwoint(:)
197 complex(real64), allocatable :: zoneint(:), ztwoint(:)
198 integer, allocatable :: iindex(:,:), jindex(:,:), kindex(:,:), lindex(:,:)
199 type(singularity_t) :: singul
200
201 push_sub(output_me)
202
203 if (this%what(option__outputmatrixelements__momentum)) then
204 write(fname,'(2a)') trim(dir), '/ks_me_momentum'
205 call output_me_out_momentum(gr, st, space, fname, namespace, hm%kpoints)
206 end if
207
208 if (this%what(option__outputmatrixelements__ang_momentum)) then
209 write(fname,'(2a)') trim(dir), '/ks_me_angular_momentum'
210 call output_me_out_ang_momentum(gr, st, space, fname, namespace, hm%kpoints)
211 end if
212
213 if (this%what(option__outputmatrixelements__ks_multipoles)) then
214 call output_me_out_ks_multipoles(gr, st, space, dir, this, namespace)
215 end if
216
217 if (this%what(option__outputmatrixelements__dipole)) then
218 assert(.not. st%parallel_in_states)
219 call output_me_out_dipole(gr, st, space, dir, namespace, hm, ions, this%st_start, this%st_end)
220 end if
221
222 if (this%what(option__outputmatrixelements__one_body)) then
223 message(1) = "Computing one-body matrix elements"
224 call messages_info(1, namespace=namespace)
225
226 if (st%parallel_in_states) call messages_not_implemented("OutputMatrixElements=one_body with states parallelization", &
227 namespace=namespace)
228 if (st%d%kpt%parallel) call messages_not_implemented("OutputMatrixElements=one_body with k-points parallelization", &
229 namespace=namespace)
230 if (family_is_mgga_with_exc(hm%xc)) then
231 call messages_not_implemented("OutputMatrixElements=one_body with MGGA", namespace=namespace)
232 end if
233 ! how to do this properly? states_elec_matrix
234 iunit = io_open(trim(dir)//'/output_me_one_body', namespace, action='write')
235
236 id = st%nst*(st%nst+1)/2
237
238 safe_allocate(iindex(1:id,1))
239 safe_allocate(jindex(1:id,1))
240
241 if (states_are_real(st)) then
242 safe_allocate(doneint(1:id))
243 call one_body_me(gr, st, namespace, hm, iindex(:,1), jindex(:,1), doneint)
244 if (mpi_grp_is_root(gr%mpi_grp)) then
245 do ll = 1, id
246 write(iunit, *) iindex(ll,1), jindex(ll,1), doneint(ll)
247 end do
248 end if
249 safe_deallocate_a(doneint)
250 else
251 safe_allocate(zoneint(1:id))
252 call one_body_me(gr, st, namespace, hm, iindex(:,1), jindex(:,1), zoneint)
253 if (mpi_grp_is_root(gr%mpi_grp)) then
254 do ll = 1, id
255 write(iunit, *) iindex(ll,1), jindex(ll,1), zoneint(ll)
256 end do
257 end if
258 safe_deallocate_a(zoneint)
259 end if
260
261 safe_deallocate_a(iindex)
262 safe_deallocate_a(jindex)
263 call io_close(iunit)
264 end if
265
266 if (this%what(option__outputmatrixelements__two_body) .or. this%what(option__outputmatrixelements__two_body_exc_k)) then
267 message(1) = "Computing two-body matrix elements"
268 call messages_info(1, namespace=namespace)
269
270 assert(.not. st%parallel_in_states)
271 if (this%what(option__outputmatrixelements__two_body)) then
272 if (st%parallel_in_states) then
273 call messages_not_implemented("OutputMatrixElements=two_body with states parallelization", namespace=namespace)
274 end if
275 if (st%d%kpt%parallel) then
276 call messages_not_implemented("OutputMatrixElements=two_body with k-points parallelization", namespace=namespace)
277 end if
278 ! how to do this properly? states_matrix
279 iunit = io_open(trim(dir)//'/output_me_two_body', namespace, action='write')
280 if (mpi_grp_is_root(gr%mpi_grp)) then
281 write(iunit, '(a)') '# (n1, k1) (n2, k2) (n3, k3) (n4, k4) (n1-k1, n2-k2|n3-k3, n4-k4)'
282 end if
283 else
284 if (st%parallel_in_states) then
285 call messages_not_implemented("OutputMatrixElements=two_body_exc_k with states parallelization", namespace=namespace)
286 end if
287 if (st%d%kpt%parallel) then
288 call messages_not_implemented("OutputMatrixElements=two_body_exc_k with k-points parallelization", namespace=namespace)
289 end if
290 ! how to do this properly? states_matrix
291 iunit = io_open(trim(dir)//'/output_me_two_body_density', namespace, action='write')
292 if (mpi_grp_is_root(gr%mpi_grp)) then
293 write(iunit, '(a)') '#(n1, k1) (n2, k2) (n1-k1, n1-k2|n2-k2, n2-k1)'
294 end if
295 end if
296
297 if (this%what(option__outputmatrixelements__two_body)) then
298 if (states_are_real(st)) then
299 id = st%nik*this%nst*(st%nik*this%nst+1)*(st%nik**2*this%nst**2+st%nik*this%nst+2)/8
300 else
301 id = (st%nik*this%nst)**4
302 end if
303 else
304 id = (st%nik*this%nst)**2
305 end if
306
307 if (states_are_complex(st)) then
308 call singularity_init(singul, namespace, space, st, hm%kpoints)
309 end if
310
311 safe_allocate(iindex(1:2, 1:id))
312 safe_allocate(jindex(1:2, 1:id))
313 safe_allocate(kindex(1:2, 1:id))
314 safe_allocate(lindex(1:2, 1:id))
315
316 if (states_are_real(st)) then
317 safe_allocate(dtwoint(1:id))
318 call two_body_me(gr, st, space, namespace, hm%kpoints, hm%exxop%psolver, this%st_start, this%st_end, &
319 iindex, jindex, kindex, lindex, dtwoint)
320 if (mpi_grp_is_root(gr%mpi_grp)) then
321 do ll = 1, id
322 write(iunit, '(4(i4,i5,1x),es22.12)') iindex(1:2,ll), jindex(1:2,ll), kindex(1:2,ll), lindex(1:2,ll), dtwoint(ll)
323 end do
324 end if
325 safe_deallocate_a(dtwoint)
326 else
327 safe_allocate(ztwoint(1:id))
328 if (hm%phase%is_allocated()) then
329 !We cannot pass the phase array like that if kpt%start is not 1.
330 assert(.not. st%d%kpt%parallel)
331 call two_body_me(gr, st, space, namespace, hm%kpoints, hm%exxop%psolver, this%st_start, this%st_end, iindex, jindex, &
332 kindex, lindex, ztwoint, phase = hm%phase, singularity = singul, &
333 exc_k = (this%what(option__outputmatrixelements__two_body_exc_k)))
334 else
335 call two_body_me(gr, st, space, namespace, hm%kpoints, hm%exxop%psolver, this%st_start, this%st_end, &
336 iindex, jindex, kindex, lindex, ztwoint, exc_k = (this%what(option__outputmatrixelements__two_body_exc_k)))
337 end if
338
339 if (mpi_grp_is_root(gr%mpi_grp)) then
340 if (this%what(option__outputmatrixelements__two_body)) then
341 do ll = 1, id
342 write(iunit, '(4(i4,i5,1x),2es22.12)') iindex(1:2,ll), jindex(1:2,ll), kindex(1:2,ll), lindex(1:2,ll), ztwoint(ll)
343 end do
344 else
345 do ll = 1, id
346 write(iunit, '(2(i4,i5),2es22.12)') iindex(1:2,ll), kindex(1:2,ll), ztwoint(ll)
347 end do
348 end if
349 end if
350 safe_deallocate_a(ztwoint)
351 end if
352
353 safe_deallocate_a(iindex)
354 safe_deallocate_a(jindex)
355 safe_deallocate_a(kindex)
356 safe_deallocate_a(lindex)
357 call io_close(iunit)
358
359 if (states_are_complex(st)) then
360 call singularity_end(singul)
361 end if
362 end if
363
364 pop_sub(output_me)
365 end subroutine output_me
366
367
368 ! ---------------------------------------------------------
369 subroutine output_me_out_momentum(gr, st, space, fname, namespace, kpoints)
370 type(grid_t), intent(in) :: gr
371 type(states_elec_t), intent(in) :: st
372 type(space_t), intent(in) :: space
373 character(len=*), intent(in) :: fname
374 type(namespace_t), intent(in) :: namespace
375 type(kpoints_t), intent(in) :: kpoints
376
377 integer :: ik, ist, is, ns, iunit, idir, dim
378 character(len=80) :: cspin, str_tmp
379 real(real64) :: kpoint(space%dim)
380 real(real64), allocatable :: momentum(:,:,:)
381
382 push_sub(output_me_out_momentum)
383
384 dim = space%dim
385
386 safe_allocate(momentum(1:dim, 1:st%nst, 1:st%nik))
387
388 ! Compute matrix elements
389 call elec_momentum_me(gr, st, space, kpoints, momentum)
390
391 iunit = io_open(fname, namespace, action='write')
392
393 ns = 1
394 if (st%d%nspin == 2) ns = 2
395
396 write(message(1),'(a)') 'Momentum of the KS states [a.u.]:'
397 call messages_info(1, iunit)
398 if (st%nik > ns) then
399 message(1) = 'k-points [' // trim(units_abbrev(unit_one/units_out%length)) // ']'
400 call messages_info(1, iunit)
401 end if
402
403 do ik = 1, st%nik, ns
404 kpoint(:) = kpoints%get_point(st%d%get_kpoint_index(ik))
405
406 if (st%nik > ns) then
407 write(message(1), '(a,i4, a)') '#k =', ik, ', k = ('
408 do idir = 1, dim
409 write(str_tmp, '(f12.6, a)') units_from_atomic(unit_one/units_out%length, kpoint(idir)), ','
410 message(1) = trim(message(1)) // trim(str_tmp)
411 if (idir == dim) then
412 message(1) = trim(message(1)) // ")"
413 else
414 message(1) = trim(message(1)) // ","
415 end if
416 end do
417 call messages_info(1, iunit)
418 end if
419
420 write(message(1), '(a4,1x,a5)') '#st',' Spin'
421 do idir = 1, dim
422 write(str_tmp, '(a,a1,a)') ' <p', index2axis(idir), '>'
423 message(1) = trim(message(1)) // trim(str_tmp)
424 end do
425 write(str_tmp, '(4x,a12,1x)') 'Occupation '
426 message(1) = trim(message(1)) // trim(str_tmp)
427 call messages_info(1, iunit)
428
429 do ist = 1, st%nst
430 do is = 0, ns-1
431
432 if (is == 0) cspin = 'up'
433 if (is == 1) cspin = 'dn'
434 if (st%d%ispin == unpolarized .or. st%d%ispin == spinors) cspin = '--'
435
436 write(message(1), '(i4,3x,a2,1x)') ist, trim(cspin)
437 do idir = 1, dim
438 write(str_tmp, '(f12.6)') momentum(idir, ist, ik+is)
439 message(1) = trim(message(1)) // trim(str_tmp)
440 end do
441 write(str_tmp, '(3x,f12.6)') st%occ(ist, ik+is)
442 message(1) = trim(message(1)) // trim(str_tmp)
443 call messages_info(1, iunit)
444
445 end do
446 end do
447
448 write(message(1),'(a)') ''
449 call messages_info(1, iunit)
450
451 end do
452
453 safe_deallocate_a(momentum)
454 call io_close(iunit)
455
457 end subroutine output_me_out_momentum
458
459
460 ! ---------------------------------------------------------
461 subroutine output_me_out_ang_momentum(gr, st, space, fname, namespace, kpoints)
462 type(grid_t), intent(in) :: gr
463 type(states_elec_t), intent(in) :: st
464 type(space_t), intent(in) :: space
465 character(len=*), intent(in) :: fname
466 type(namespace_t), intent(in) :: namespace
467 type(kpoints_t), intent(in) :: kpoints
468
469 integer :: iunit, ik, ist, is, ns, idir, k_start, k_end, k_n, dim, st_start, st_end, nst
470 character(len=80) :: tmp_str(space%dim), cspin
471 real(real64) :: angular(3), lsquare, kpoint(space%dim)
472 real(real64), allocatable :: ang(:, :, :), ang2(:, :)
473#if defined(HAVE_MPI)
474 integer :: tmp
475#endif
476 real(real64), allocatable :: lang(:, :)
477
479
480 ns = 1
481 if (st%d%nspin == 2) ns = 2
482 assert(space%dim == 3)
483 ! Define useful quantities
484 dim = space%dim ! Space dimension
485 st_start = st%st_start ! Starting state
486 st_end = st%st_end ! Last state
487 nst = st%nst ! Number of states
488
489 iunit = io_open(fname, namespace, action='write')
490
491 if (mpi_grp_is_root(mpi_world)) then
492 write(iunit,'(a)') 'Warning: When non-local pseudopotentials are used '
493 write(iunit,'(a)') ' the numbers below may be meaningless. '
494 write(iunit,'(a)') ' '
495 write(iunit,'(a)') 'Angular Momentum of the KS states [dimensionless]:'
496 ! r x k is dimensionless. we do not include hbar.
497 if (st%nik > ns) then
498 message(1) = 'k-points [' // trim(units_abbrev(unit_one/units_out%length)) // ']'
499 call messages_info(1, iunit)
500 end if
501 end if
502
503 safe_allocate(ang(1:nst, 1:st%nik, 1:3))
504 safe_allocate(ang2(1:nst, 1:st%nik))
505
506 ! Compute matrix elements
507 call elec_angular_momentum_me(gr, st, space, ang, ang2)
508
509 k_start = st%d%kpt%start
510 k_end = st%d%kpt%end
511 do idir = 1, 3
512 angular(idir) = states_elec_eigenvalues_sum(st, ang(st_start:st_end, k_start:k_end, idir))
513 end do
514 lsquare = states_elec_eigenvalues_sum(st, ang2(st_start:st_end, k_start:k_end))
515
516 if (st%d%kpt%parallel) then
517 k_n = st%d%kpt%nlocal
518
519 assert(.not. st%parallel_in_states)
520
521 ! note: could use lmpi_gen_allgatherv here?
522 safe_allocate(lang(1:st%lnst, 1:k_n))
523 do idir = 1, 3
524 lang(1:st%lnst, 1:k_n) = ang(st_start:st_end, k_start:k_end, idir)
525 call st%d%kpt%mpi_grp%allgatherv(lang, nst*k_n, mpi_double_precision, ang(:, :, idir), &
526 st%d%kpt%num(:)*nst, (st%d%kpt%range(1, :) - 1)*nst, mpi_double_precision)
527 end do
528 lang(1:st%lnst, 1:k_n) = ang2(st_start:st_end, k_start:k_end)
529 call st%d%kpt%mpi_grp%allgatherv(lang, nst*k_n, mpi_double_precision, ang2, &
530 st%d%kpt%num(:)*nst, (st%d%kpt%range(1, :) - 1)*nst, mpi_double_precision)
531 safe_deallocate_a(lang)
532 end if
533
534 if (st%parallel_in_states) then
535 safe_allocate(lang(1:st%lnst, 1))
536 end if
537
538 do ik = 1, st%nik, ns
539 if (st%nik > ns) then
540
541 kpoint(:) = kpoints%get_point(st%d%get_kpoint_index(ik))
542
543 write(message(1), '(a,i4, a)') '#k =', ik, ', k = ('
544 do idir = 1, dim
545 write(tmp_str(1), '(f12.6, a)') units_from_atomic(unit_one/units_out%length, kpoint(idir)), ','
546 message(1) = trim(message(1)) // trim(tmp_str(1))
547 if (idir == dim) then
548 message(1) = trim(message(1)) // ")"
549 else
550 message(1) = trim(message(1)) // ","
551 end if
552 end do
553 call messages_info(1, iunit)
554 end if
555
556 ! Exchange ang and ang2.
557 if (st%parallel_in_states) then
558 assert(.not. st%d%kpt%parallel)
559#if defined(HAVE_MPI)
560 do is = 1, ns
561 do idir = 1, 3
562 lang(1:st%lnst, 1) = ang(st_start:st_end, ik+is-1, idir)
563 call lmpi_gen_allgatherv(st%lnst, lang(:, 1), tmp, ang(:, ik+is-1, idir), st%mpi_grp)
564 end do
565 lang(1:st%lnst, 1) = ang2(st_start:st_end, ik+is-1)
566 call lmpi_gen_allgatherv(st%lnst, lang(:, 1), tmp, ang2(:, ik+is-1), st%mpi_grp)
567 end do
568#endif
569 end if
570
571 write(message(1), '(a4,1x,a5,4a12,4x,a12,1x)') &
572 '#st',' Spin',' <Lx>', ' <Ly>', ' <Lz>', ' <L2>', 'Occupation '
573 call messages_info(1, iunit)
574
575 if (mpi_grp_is_root(mpi_world)) then
576 do ist = 1, nst
577 do is = 0, ns-1
578
579 if (is == 0) cspin = 'up'
580 if (is == 1) cspin = 'dn'
581 if (st%d%ispin == unpolarized .or. st%d%ispin == spinors) cspin = '--'
582
583 write(tmp_str(1), '(i4,3x,a2)') ist, trim(cspin)
584 write(tmp_str(2), '(1x,4f12.6,3x,f12.6)') &
585 (ang(ist, ik+is, idir), idir = 1, 3), ang2(ist, ik+is), st%occ(ist, ik+is)
586 message(1) = trim(tmp_str(1))//trim(tmp_str(2))
587 call messages_info(1, iunit)
588 end do
589 end do
590 end if
591 write(message(1),'(a)') ''
592 call messages_info(1, iunit)
593
594 end do
595
596 write(message(1),'(a)') 'Total Angular Momentum L [dimensionless]'
597 write(message(2),'(10x,4f12.6)') angular(1:3), lsquare
598 call messages_info(2, iunit)
599
600 call io_close(iunit)
601
602 if (st%parallel_in_states) then
603 safe_deallocate_a(lang)
604 end if
605
606 safe_deallocate_a(ang)
607 safe_deallocate_a(ang2)
608
610 end subroutine output_me_out_ang_momentum
611
612 subroutine output_me_out_dipole(gr, st, space, dir, namespace, hm, ions, st_start, st_end)
613 type(grid_t), intent(in) :: gr
614 type(states_elec_t), intent(in) :: st
615 type(space_t), intent(in) :: space
616 character(len=*), intent(in) :: dir
617 type(namespace_t), intent(in) :: namespace
618 type(hamiltonian_elec_t), intent(in) :: hm
619 type(ions_t), intent(in) :: ions
620 integer, intent(in) :: st_start, st_end
621
622 integer :: ik, idir, ist, jst, iunit
623 character(len=256) :: fname
624 real(real64), allocatable :: ddip_elements(:,:,:)
625 complex(real64), allocatable :: zdip_elements(:,:,:)
626
627 push_sub(output_me_out_dipole)
628
629 ! The content of each file should be clear from the header of each file.
630 do ik = st%d%kpt%start, st%d%kpt%end
631 write(fname,'(i8)') ik
632 write(fname,'(a)') trim(dir)//'/ks_me_dipole.k'//trim(adjustl(fname))//'_'
633
634 if (states_are_real(st)) then
635 ! Allocate dipole array
636 safe_allocate(ddip_elements(1:space%dim, st_start:st_end, st_start:st_end))
637 ! Compute dipole elements
638 call dipole_me(gr, st, namespace, hm, ions, ik, st_start, st_end, ddip_elements)
639 ! Output elements
640 do idir = 1, space%dim
641 iunit = io_open(trim(fname)//index2axis(idir), namespace, action = 'write')
642 write(iunit, '(a)') '# Dipole matrix elements file: |<Phi_i | r | Phi_j>|'
643 write(iunit, '(a,i8)') '# ik =', ik
644 write(iunit, '(a)') '# Units = ['//trim(units_abbrev(units_out%length))//']'
645
646 do ist = st_start, st_end
647 do jst = st_start, st_end
648 write(iunit, fmt='(f20.12)', advance = 'no') units_from_atomic(units_out%length, abs(ddip_elements(idir, ist, jst)))
649 write(iunit, fmt='(a)', advance = 'no') ' '
650 end do
651 write(iunit, '(a)') '' ! New line
652 end do
653 call io_close(iunit)
654 end do
655 safe_deallocate_a(ddip_elements)
656 else
657 ! Allocate dipole array
658 safe_allocate(zdip_elements(1:space%dim, st_start:st_end, st_start:st_end))
659 ! Compute dipole elements
660 call dipole_me(gr, st, namespace, hm, ions, ik, st_start, st_end, zdip_elements)
661 ! Output elements
662 do idir = 1, space%dim
663 iunit = io_open(trim(fname)//index2axis(idir), namespace, action = 'write')
664 write(iunit, '(a)') '# Dipole matrix elements file: |<Phi_i | r | Phi_j>|'
665 write(iunit, '(a,i8)') '# ik =', ik
666 write(iunit, '(a)') '# Units = ['//trim(units_abbrev(units_out%length))//']'
667
668 do ist = st_start, st_end
669 do jst = st_start, st_end
670 write(iunit, fmt='(f20.12)', advance = 'no') units_from_atomic(units_out%length, abs(zdip_elements(idir, ist, jst)))
671 write(iunit, fmt='(a)', advance = 'no') ' '
672 end do
673 write(iunit, '(a)') '' ! New line
674 end do
675 call io_close(iunit)
676 end do
677 safe_deallocate_a(zdip_elements)
678 end if
679 end do
680
681 pop_sub(output_me_out_dipole)
682 end subroutine output_me_out_dipole
683
684 subroutine output_me_out_ks_multipoles(gr, st, space, dir, this, namespace)
685 type(grid_t), intent(in) :: gr
686 type(states_elec_t), intent(in) :: st
687 type(space_t), intent(in) :: space
688 character(len=*), intent(in) :: dir
689 type(output_me_t), intent(in) :: this
690 type(namespace_t), intent(in) :: namespace
691
692 integer :: id, ll, mm, ik, iunit
693 character(len=256) :: fname
694 real(real64), allocatable :: delements(:,:)
695 complex(real64), allocatable :: zelements(:,:)
696
698
699 ! The content of each file should be clear from the header of each file.
700 id = 1
701 do ik = 1, st%nik
702 select case (space%dim)
703 case (3)
704 do ll = 1, this%ks_multipoles
705 do mm = -ll, ll
706 ! Define file name
707 write(fname,'(i4)') id
708 write(fname,'(a)') trim(dir) // '/ks_me_multipoles.' // trim(adjustl(fname))
709 iunit = io_open(fname, namespace, action = 'write')
710 ! Print header
711 call print_ks_multipole_header(iunit, ll, mm, ik)
712 ! Compute and print elements
713 if (states_are_real(st)) then
714 safe_allocate(delements(1:st%nst, 1:st%nst))
715 call ks_multipoles_3d(gr, st, space, ll, mm, ik, delements)
716 call dprint_ks_multipole(iunit, st%nst, delements, ll)
717 safe_deallocate_a(delements)
718 else
719 safe_allocate(zelements(1:st%nst, 1:st%nst))
720 call ks_multipoles_3d(gr, st, space, ll, mm, ik, zelements)
721 call zprint_ks_multipole(iunit, st%nst, zelements, ll)
722 safe_deallocate_a(zelements)
723 end if
724 ! Close file
725 call io_close(iunit)
726 id = id + 1
727 end do
728 end do
729 case (2)
730 do ll = 1, 2
731 ! Define file name
732 write(fname,'(i4)') id
733 write(fname,'(a)') trim(dir) // '/ks_me_multipoles.' // trim(adjustl(fname))
734 iunit = io_open(fname, namespace, action = 'write')
735 ! Print header
736 call print_ks_multipole_header(iunit, ll, mm, ik)
737 ! Compute and print elements
738 if (states_are_real(st)) then
739 safe_allocate(delements(1:st%nst, 1:st%nst))
740 call ks_multipoles_2d(gr, st, ll, ik, delements)
741 call dprint_ks_multipole(iunit, st%nst, delements)
742 safe_deallocate_a(delements)
743 else
744 safe_allocate(zelements(1:st%nst, 1:st%nst))
745 call ks_multipoles_2d(gr, st, ll, ik, zelements)
746 call zprint_ks_multipole(iunit, st%nst, zelements)
747 safe_deallocate_a(zelements)
748 end if
749 ! Close file
750 call io_close(iunit)
751 id = id + 1
752 end do
753 case (1)
754 do ll = 1, this%ks_multipoles
755 ! Define file name
756 write(fname,'(i4)') id
757 write(fname,'(a)') trim(dir) // '/ks_me_multipoles.' // trim(adjustl(fname))
758 iunit = io_open(fname, namespace, action = 'write')
759 ! Print header
760 call print_ks_multipole_header(iunit, ll, mm, ik)
761 ! Compute and print elements
762 if (states_are_real(st)) then
763 safe_allocate(delements(1:st%nst, 1:st%nst))
764 call ks_multipoles_1d(gr, st, ll, ik, delements)
765 call dprint_ks_multipole(iunit, st%nst, delements, ll)
766 safe_deallocate_a(delements)
767 else
768 safe_allocate(zelements(1:st%nst, 1:st%nst))
769 call ks_multipoles_1d(gr, st, ll, ik, zelements)
770 call zprint_ks_multipole(iunit, st%nst, zelements, ll)
771 safe_deallocate_a(zelements)
772 end if
773 ! Close file
774 call io_close(iunit)
775 id = id + 1
776 end do
777 end select
778 end do
779
781
782 contains
783
784 subroutine print_ks_multipole_header(iunit, ll, mm, ik)
785 integer, intent(in) :: iunit, ll, mm, ik
786
787 write(iunit, fmt = '(a)') '# Multipole matrix elements file: <Phi_i | r**l * Y_{lm}(theta,phi) | Phi_j>'
788 write(iunit, fmt = '(a,i2,a,i2)') '# l =', ll, '; m =', mm
789 write(iunit, fmt = '(a,i8)') '# ik =', ik
790 if (ll>1) then
791 write(iunit, fmt = '(a,i1)') '# Units = ['//trim(units_abbrev(units_out%length))//']^', ll
792 else
793 write(iunit, fmt = '(a)') '# Units = ['//trim(units_abbrev(units_out%length))//']'
794 end if
795 end subroutine print_ks_multipole_header
796
797 subroutine dprint_ks_multipole(iunit, nst, elements, ll)
798 integer, intent(in) :: iunit, nst
799 real(real64), intent(in) :: elements(:,:)
800 integer, optional, intent(in) :: ll
801
802 integer :: ist, jst
803 real(real64) :: element
804
805 do ist = 1, nst
806 do jst = 1, nst
807 element = units_from_atomic(units_out%length**optional_default(ll, 1), elements(ist, jst))
808 write(iunit, fmt='(f20.12, a)', advance = 'no') element, ' '
809 end do
810 write(iunit, '(a)') ''
811 end do
812 end subroutine dprint_ks_multipole
813
814 subroutine zprint_ks_multipole(iunit, nst, elements, ll)
815 integer, intent(in) :: iunit, nst
816 complex(real64), intent(in) :: elements(:,:)
817 integer, optional, intent(in) :: ll
818
819 integer :: ist, jst
820 complex(real64) :: element
821
822 do ist = 1, nst
823 do jst = 1, nst
824 element = units_from_atomic(units_out%length**optional_default(ll, 1), elements(ist, jst))
825 write(iunit, fmt='(f20.12,a1,f20.12,a)', advance = 'no') real(element),',',aimag(element), ' '
826 end do
827 write(iunit, '(a)') ''
828 end do
829 end subroutine zprint_ks_multipole
830
831 end subroutine output_me_out_ks_multipoles
832
833
834#include "undef.F90"
835#include "real.F90"
836
837#include "undef.F90"
838#include "complex.F90"
839
840end module output_me_oct_m
841
842!! Local Variables:
843!! mode: f90
844!! coding: utf-8
845!! End:
subroutine, public elec_angular_momentum_me(gr, st, space, ll, l2)
subroutine, public elec_momentum_me(gr, st, space, kpoints, momentum)
integer, parameter, public unpolarized
Parameters...
integer, parameter, public spinors
This module implements the underlying real-space grid.
Definition: grid.F90:117
subroutine, public io_function_read_what_how_when(namespace, space, what, how, output_interval, what_tag_in, how_tag_in, output_interval_tag_in, ignore_error)
Definition: io.F90:114
subroutine, public io_close(iunit, grp)
Definition: io.F90:418
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:352
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1113
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:616
This module contains some common usage patterns of MPI routines.
Definition: mpi_lib.F90:115
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 low-level part of the output system
Definition: output_low.F90:115
subroutine output_me_out_ang_momentum(gr, st, space, fname, namespace, kpoints)
Definition: output_me.F90:555
subroutine output_me_out_momentum(gr, st, space, fname, namespace, kpoints)
Definition: output_me.F90:463
subroutine output_me_out_dipole(gr, st, space, dir, namespace, hm, ions, st_start, st_end)
Definition: output_me.F90:706
subroutine, public output_me_init(this, namespace, space, st, gr, nst)
Definition: output_me.F90:155
subroutine, public output_me(this, namespace, space, dir, st, gr, ions, hm)
Definition: output_me.F90:278
subroutine output_me_out_ks_multipoles(gr, st, space, dir, this, namespace)
Definition: output_me.F90:778
subroutine, public singularity_end(this)
subroutine, public singularity_init(this, namespace, space, st, kpoints)
pure logical function, public states_are_complex(st)
pure logical function, public states_are_real(st)
This module handles spin dimensions of the states and the k-point distribution.
real(real64) function, public states_elec_eigenvalues_sum(st, alt_eig)
function to calculate the eigenvalues sum using occupations as weights
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
character(len=20) pure function, public units_abbrev(this)
Definition: unit.F90:223
This module defines the unit system, used for input and output.
type(unit_system_t), public units_out
type(unit_t), public unit_one
some special units required for particular quantities
This module is intended to contain simple general-purpose utility functions and procedures.
Definition: utils.F90:118
character pure function, public index2axis(idir)
Definition: utils.F90:202
subroutine print_ks_multipole_header(iunit, ll, mm, ik)
Definition: output_me.F90:878
subroutine zprint_ks_multipole(iunit, nst, elements, ll)
Definition: output_me.F90:908
subroutine dprint_ks_multipole(iunit, nst, elements, ll)
Definition: output_me.F90:891
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:168
Output information for matrix elements.
Definition: output_low.F90:131
The states_elec_t class contains all electronic wave functions.