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 ! Initialize matrix element calculator
180 call this%elec_me%init(gr, space, st, this%st_start, this%st_end)
181
182 pop_sub(output_me_init)
183 end subroutine output_me_init
184
185
186 ! ---------------------------------------------------------
187 subroutine output_me(this, namespace, space, dir, st, gr, ions, hm)
188 type(output_me_t), intent(in) :: this
189 type(namespace_t), intent(in) :: namespace
190 class(space_t), intent(in) :: space
191 character(len=*), intent(in) :: dir
192 type(states_elec_t), intent(inout) :: st
193 type(grid_t), intent(in) :: gr
194 type(ions_t), intent(in) :: ions
195 type(hamiltonian_elec_t), intent(inout) :: hm
196
197 integer :: id, ll, iunit
198 character(len=256) :: fname
199 real(real64), allocatable :: doneint(:), dtwoint(:)
200 complex(real64), allocatable :: zoneint(:), ztwoint(:)
201 integer, allocatable :: iindex(:,:), jindex(:,:), kindex(:,:), lindex(:,:)
202 type(singularity_t) :: singul
203
204 push_sub(output_me)
205
206 if (this%what(option__outputmatrixelements__momentum)) then
207 write(fname,'(2a)') trim(dir), '/ks_me_momentum'
208 call output_me_out_momentum(fname, this%elec_me, namespace, hm%kpoints)
209 end if
210
211 if (this%what(option__outputmatrixelements__ang_momentum)) then
212 write(fname,'(2a)') trim(dir), '/ks_me_angular_momentum'
213 call output_me_out_ang_momentum(fname, this%elec_me, namespace, hm%kpoints)
214 end if
215
216 if (this%what(option__outputmatrixelements__ks_multipoles)) then
217 call output_me_out_ks_multipoles(dir, this, this%elec_me, namespace)
218 end if
219
220 if (this%what(option__outputmatrixelements__dipole)) then
221 assert(.not. st%parallel_in_states)
222 call output_me_out_dipole(dir, this%elec_me, namespace, hm, ions)
223 end if
224
225 if (this%what(option__outputmatrixelements__one_body)) then
226 message(1) = "Computing one-body matrix elements"
227 call messages_info(1, namespace=namespace)
228
229 if (st%parallel_in_states) call messages_not_implemented("OutputMatrixElements=one_body with states parallelization", &
230 namespace=namespace)
231 if (st%d%kpt%parallel) call messages_not_implemented("OutputMatrixElements=one_body with k-points parallelization", &
232 namespace=namespace)
233 if (family_is_mgga_with_exc(hm%xc)) then
234 call messages_not_implemented("OutputMatrixElements=one_body with MGGA", namespace=namespace)
235 end if
236 ! how to do this properly? states_elec_matrix
237 iunit = io_open(trim(dir)//'/output_me_one_body', namespace, action='write')
238
239 id = st%nst*(st%nst+1)/2
240
241 safe_allocate(iindex(1:id,1))
242 safe_allocate(jindex(1:id,1))
243
244 if (states_are_real(st)) then
245 safe_allocate(doneint(1:id))
246 call this%elec_me%one_body_me(namespace, hm, iindex(:,1), jindex(:,1), doneint)
247 if (mpi_grp_is_root(gr%mpi_grp)) then
248 do ll = 1, id
249 write(iunit, *) iindex(ll,1), jindex(ll,1), doneint(ll)
250 end do
251 end if
252 safe_deallocate_a(doneint)
253 else
254 safe_allocate(zoneint(1:id))
255 call this%elec_me%one_body_me(namespace, hm, iindex(:,1), jindex(:,1), zoneint)
256 if (mpi_grp_is_root(gr%mpi_grp)) then
257 do ll = 1, id
258 write(iunit, *) iindex(ll,1), jindex(ll,1), zoneint(ll)
259 end do
260 end if
261 safe_deallocate_a(zoneint)
262 end if
263
264 safe_deallocate_a(iindex)
265 safe_deallocate_a(jindex)
266 call io_close(iunit)
267 end if
268
269 if (this%what(option__outputmatrixelements__two_body) .or. this%what(option__outputmatrixelements__two_body_exc_k)) then
270 message(1) = "Computing two-body matrix elements"
271 call messages_info(1, namespace=namespace)
272
273 assert(.not. st%parallel_in_states)
274 if (this%what(option__outputmatrixelements__two_body)) then
275 if (st%parallel_in_states) then
276 call messages_not_implemented("OutputMatrixElements=two_body with states parallelization", namespace=namespace)
277 end if
278 if (st%d%kpt%parallel) then
279 call messages_not_implemented("OutputMatrixElements=two_body with k-points parallelization", namespace=namespace)
280 end if
281 ! how to do this properly? states_matrix
282 iunit = io_open(trim(dir)//'/output_me_two_body', namespace, action='write')
283 if (mpi_grp_is_root(gr%mpi_grp)) then
284 write(iunit, '(a)') '# (n1, k1) (n2, k2) (n3, k3) (n4, k4) (n1-k1, n2-k2|n3-k3, n4-k4)'
285 end if
286 else
287 if (st%parallel_in_states) then
288 call messages_not_implemented("OutputMatrixElements=two_body_exc_k with states parallelization", namespace=namespace)
289 end if
290 if (st%d%kpt%parallel) then
291 call messages_not_implemented("OutputMatrixElements=two_body_exc_k with k-points parallelization", namespace=namespace)
292 end if
293 ! how to do this properly? states_matrix
294 iunit = io_open(trim(dir)//'/output_me_two_body_density', namespace, action='write')
295 if (mpi_grp_is_root(gr%mpi_grp)) then
296 write(iunit, '(a)') '#(n1, k1) (n2, k2) (n1-k1, n1-k2|n2-k2, n2-k1)'
297 end if
298 end if
299
300 if (this%what(option__outputmatrixelements__two_body)) then
301 if (states_are_real(st)) then
302 id = st%nik*this%nst*(st%nik*this%nst+1)*(st%nik**2*this%nst**2+st%nik*this%nst+2)/8
303 else
304 id = (st%nik*this%nst)**4
305 end if
306 else
307 id = (st%nik*this%nst)**2
308 end if
309
310 if (states_are_complex(st)) then
311 call singularity_init(singul, namespace, space, st, hm%kpoints)
312 end if
313
314 safe_allocate(iindex(1:2, 1:id))
315 safe_allocate(jindex(1:2, 1:id))
316 safe_allocate(kindex(1:2, 1:id))
317 safe_allocate(lindex(1:2, 1:id))
318
319 if (states_are_real(st)) then
320 safe_allocate(dtwoint(1:id))
321 call this%elec_me%two_body_me(namespace, hm%kpoints, hm%exxop%psolver, this%st_start, this%st_end, &
322 iindex, jindex, kindex, lindex, dtwoint)
323 if (mpi_grp_is_root(gr%mpi_grp)) then
324 do ll = 1, id
325 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)
326 end do
327 end if
328 safe_deallocate_a(dtwoint)
329 else
330 safe_allocate(ztwoint(1:id))
331 if (hm%phase%is_allocated()) then
332 !We cannot pass the phase array like that if kpt%start is not 1.
333 assert(.not. st%d%kpt%parallel)
334 call this%elec_me%two_body_me(namespace, hm%kpoints, hm%exxop%psolver, this%st_start, this%st_end, iindex, jindex, &
335 kindex, lindex, ztwoint, phase = hm%phase, singularity = singul, &
336 exc_k = (this%what(option__outputmatrixelements__two_body_exc_k)))
337 else
338 call this%elec_me%two_body_me(namespace, hm%kpoints, hm%exxop%psolver, this%st_start, this%st_end, &
339 iindex, jindex, kindex, lindex, ztwoint, exc_k = (this%what(option__outputmatrixelements__two_body_exc_k)))
340 end if
341
342 if (mpi_grp_is_root(gr%mpi_grp)) then
343 if (this%what(option__outputmatrixelements__two_body)) then
344 do ll = 1, id
345 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)
346 end do
347 else
348 do ll = 1, id
349 write(iunit, '(2(i4,i5),2es22.12)') iindex(1:2,ll), kindex(1:2,ll), ztwoint(ll)
350 end do
351 end if
352 end if
353 safe_deallocate_a(ztwoint)
354 end if
355
356 safe_deallocate_a(iindex)
357 safe_deallocate_a(jindex)
358 safe_deallocate_a(kindex)
359 safe_deallocate_a(lindex)
360 call io_close(iunit)
361
362 if (states_are_complex(st)) then
363 call singularity_end(singul)
364 end if
365 end if
366
367 pop_sub(output_me)
368 end subroutine output_me
369
370
371 ! ---------------------------------------------------------
372 subroutine output_me_out_momentum(fname, elec_me, namespace, kpoints)
373 character(len=*), intent(in) :: fname
374 type(elec_matrix_elements_t), intent(in) :: elec_me
375 type(namespace_t), intent(in) :: namespace
376 type(kpoints_t), intent(in) :: kpoints
377
378 integer :: ik, ist, is, ns, iunit, idir, dim
379 character(len=80) :: cspin, str_tmp
380 real(real64) :: kpoint(elec_me%space%dim)
381 real(real64), allocatable :: momentum(:,:,:)
382
383 push_sub(output_me_out_momentum)
384
385 dim = elec_me%space%dim
386
387 safe_allocate(momentum(1:dim, 1:elec_me%states%nst, 1:elec_me%states%nik))
388
389 ! Compute matrix elements
390 call elec_me%momentum_me(kpoints, momentum)
391
392 iunit = io_open(fname, namespace, action='write')
393
394 ns = 1
395 if (elec_me%states%d%nspin == 2) ns = 2
396
397 write(message(1),'(a)') 'Momentum of the KS states [a.u.]:'
398 call messages_info(1, iunit)
399 if (elec_me%states%nik > ns) then
400 message(1) = 'k-points [' // trim(units_abbrev(unit_one/units_out%length)) // ']'
401 call messages_info(1, iunit)
402 end if
403
404 do ik = 1, elec_me%states%nik, ns
405 kpoint(:) = kpoints%get_point(elec_me%states%d%get_kpoint_index(ik))
406
407 if (elec_me%states%nik > ns) then
408 write(message(1), '(a,i4, a)') '#k =', ik, ', k = ('
409 do idir = 1, dim
410 write(str_tmp, '(f12.6, a)') units_from_atomic(unit_one/units_out%length, kpoint(idir)), ','
411 message(1) = trim(message(1)) // trim(str_tmp)
412 if (idir == dim) then
413 message(1) = trim(message(1)) // ")"
414 else
415 message(1) = trim(message(1)) // ","
416 end if
417 end do
418 call messages_info(1, iunit)
419 end if
420
421 write(message(1), '(a4,1x,a5)') '#st',' Spin'
422 do idir = 1, dim
423 write(str_tmp, '(a,a1,a)') ' <p', index2axis(idir), '>'
424 message(1) = trim(message(1)) // trim(str_tmp)
425 end do
426 write(str_tmp, '(4x,a12,1x)') 'Occupation '
427 message(1) = trim(message(1)) // trim(str_tmp)
428 call messages_info(1, iunit)
429
430 do ist = 1, elec_me%states%nst
431 do is = 0, ns-1
432
433 if (is == 0) cspin = 'up'
434 if (is == 1) cspin = 'dn'
435 if (elec_me%states%d%ispin == unpolarized .or. elec_me%states%d%ispin == spinors) cspin = '--'
436
437 write(message(1), '(i4,3x,a2,1x)') ist, trim(cspin)
438 do idir = 1, dim
439 write(str_tmp, '(f12.6)') momentum(idir, ist, ik+is)
440 message(1) = trim(message(1)) // trim(str_tmp)
441 end do
442 write(str_tmp, '(3x,f12.6)') elec_me%states%occ(ist, ik+is)
443 message(1) = trim(message(1)) // trim(str_tmp)
444 call messages_info(1, iunit)
445
446 end do
447 end do
448
449 write(message(1),'(a)') ''
450 call messages_info(1, iunit)
451
452 end do
453
454 safe_deallocate_a(momentum)
455 call io_close(iunit)
456
458 end subroutine output_me_out_momentum
459
460
461 ! ---------------------------------------------------------
462 subroutine output_me_out_ang_momentum(fname, elec_me, namespace, kpoints)
463 character(len=*), intent(in) :: fname
464 type(elec_matrix_elements_t), intent(in) :: elec_me
465 type(namespace_t), intent(in) :: namespace
466 type(kpoints_t), intent(in) :: kpoints
467
468 integer :: iunit, ik, ist, is, ns, idir, k_start, k_end, k_n, dim, st_start, st_end, nst
469 character(len=80) :: tmp_str(elec_me%space%dim), cspin
470 real(real64) :: angular(3), lsquare, kpoint(elec_me%space%dim)
471 real(real64), allocatable :: ang(:, :, :), ang2(:, :)
472#if defined(HAVE_MPI)
473 integer :: tmp
474#endif
475 real(real64), allocatable :: lang(:, :)
476
478
479 ns = 1
480 if (elec_me%states%d%nspin == 2) ns = 2
481 assert(elec_me%space%dim == 3)
482 ! Define useful quantities
483 dim = elec_me%space%dim ! Space dimension
484 st_start = elec_me%states%st_start ! Starting state
485 st_end = elec_me%states%st_end ! Last state
486 nst = elec_me%states%nst ! Number of states
487
488 iunit = io_open(fname, namespace, action='write')
489
490 if (mpi_grp_is_root(mpi_world)) then
491 write(iunit,'(a)') 'Warning: When non-local pseudopotentials are used '
492 write(iunit,'(a)') ' the numbers below may be meaningless. '
493 write(iunit,'(a)') ' '
494 write(iunit,'(a)') 'Angular Momentum of the KS states [dimensionless]:'
495 ! r x k is dimensionless. we do not include hbar.
496 if (elec_me%states%nik > ns) then
497 message(1) = 'k-points [' // trim(units_abbrev(unit_one/units_out%length)) // ']'
498 call messages_info(1, iunit)
499 end if
500 end if
501
502 safe_allocate(ang(1:nst, 1:elec_me%states%nik, 1:3))
503 safe_allocate(ang2(1:nst, 1:elec_me%states%nik))
504
505 ! Compute matrix elements
506 call elec_me%angular_momentum_me(ang, ang2)
507
508 k_start = elec_me%states%d%kpt%start
509 k_end = elec_me%states%d%kpt%end
510 do idir = 1, 3
511 angular(idir) = states_elec_eigenvalues_sum(elec_me%states, ang(st_start:st_end, k_start:k_end, idir))
512 end do
513 lsquare = states_elec_eigenvalues_sum(elec_me%states, ang2(st_start:st_end, k_start:k_end))
514
515 if (elec_me%states%d%kpt%parallel) then
516 k_n = elec_me%states%d%kpt%nlocal
517
518 assert(.not. elec_me%states%parallel_in_states)
519
520 ! note: could use lmpi_gen_allgatherv here?
521 safe_allocate(lang(1:elec_me%states%lnst, 1:k_n))
522 do idir = 1, 3
523 lang(1:elec_me%states%lnst, 1:k_n) = ang(st_start:st_end, k_start:k_end, idir)
524 call elec_me%states%d%kpt%mpi_grp%allgatherv(lang, nst*k_n, mpi_double_precision, ang(:, :, idir), &
525 elec_me%states%d%kpt%num(:)*nst, (elec_me%states%d%kpt%range(1, :) - 1)*nst, mpi_double_precision)
526 end do
527 lang(1:elec_me%states%lnst, 1:k_n) = ang2(st_start:st_end, k_start:k_end)
528 call elec_me%states%d%kpt%mpi_grp%allgatherv(lang, nst*k_n, mpi_double_precision, ang2, &
529 elec_me%states%d%kpt%num(:)*nst, (elec_me%states%d%kpt%range(1, :) - 1)*nst, mpi_double_precision)
530 safe_deallocate_a(lang)
531 end if
532
533 if (elec_me%states%parallel_in_states) then
534 safe_allocate(lang(1:elec_me%states%lnst, 1))
535 end if
536
537 do ik = 1, elec_me%states%nik, ns
538 if (elec_me%states%nik > ns) then
539
540 kpoint(:) = kpoints%get_point(elec_me%states%d%get_kpoint_index(ik))
541
542 write(message(1), '(a,i4, a)') '#k =', ik, ', k = ('
543 do idir = 1, dim
544 write(tmp_str(1), '(f12.6, a)') units_from_atomic(unit_one/units_out%length, kpoint(idir)), ','
545 message(1) = trim(message(1)) // trim(tmp_str(1))
546 if (idir == dim) then
547 message(1) = trim(message(1)) // ")"
548 else
549 message(1) = trim(message(1)) // ","
550 end if
551 end do
552 call messages_info(1, iunit)
553 end if
554
555 ! Exchange ang and ang2.
556 if (elec_me%states%parallel_in_states) then
557 assert(.not. elec_me%states%d%kpt%parallel)
558#if defined(HAVE_MPI)
559 do is = 1, ns
560 do idir = 1, 3
561 lang(1:elec_me%states%lnst, 1) = ang(st_start:st_end, ik+is-1, idir)
562 call lmpi_gen_allgatherv(elec_me%states%lnst, lang(:, 1), tmp, ang(:, ik+is-1, idir), elec_me%states%mpi_grp)
563 end do
564 lang(1:elec_me%states%lnst, 1) = ang2(st_start:st_end, ik+is-1)
565 call lmpi_gen_allgatherv(elec_me%states%lnst, lang(:, 1), tmp, ang2(:, ik+is-1), elec_me%states%mpi_grp)
566 end do
567#endif
568 end if
569
570 write(message(1), '(a4,1x,a5,4a12,4x,a12,1x)') &
571 '#st',' Spin',' <Lx>', ' <Ly>', ' <Lz>', ' <L2>', 'Occupation '
572 call messages_info(1, iunit)
573
574 if (mpi_grp_is_root(mpi_world)) then
575 do ist = 1, nst
576 do is = 0, ns-1
577
578 if (is == 0) cspin = 'up'
579 if (is == 1) cspin = 'dn'
580 if (elec_me%states%d%ispin == unpolarized .or. elec_me%states%d%ispin == spinors) cspin = '--'
581
582 write(tmp_str(1), '(i4,3x,a2)') ist, trim(cspin)
583 write(tmp_str(2), '(1x,4f12.6,3x,f12.6)') &
584 (ang(ist, ik+is, idir), idir = 1, 3), ang2(ist, ik+is), elec_me%states%occ(ist, ik+is)
585 message(1) = trim(tmp_str(1))//trim(tmp_str(2))
586 call messages_info(1, iunit)
587 end do
588 end do
589 end if
590 write(message(1),'(a)') ''
591 call messages_info(1, iunit)
592
593 end do
594
595 write(message(1),'(a)') 'Total Angular Momentum L [dimensionless]'
596 write(message(2),'(10x,4f12.6)') angular(1:3), lsquare
597 call messages_info(2, iunit)
598
599 call io_close(iunit)
600
601 if (elec_me%states%parallel_in_states) then
602 safe_deallocate_a(lang)
603 end if
604
605 safe_deallocate_a(ang)
606 safe_deallocate_a(ang2)
607
609 end subroutine output_me_out_ang_momentum
610
611 subroutine output_me_out_dipole(dir, elec_me, namespace, hm, ions)
612 character(len=*), intent(in) :: dir
613 type(elec_matrix_elements_t), intent(in) :: elec_me
614 type(namespace_t), intent(in) :: namespace
615 type(hamiltonian_elec_t), intent(in) :: hm
616 type(ions_t), intent(in) :: ions
617
618 integer :: ik, idir, ist, jst, iunit
619 character(len=256) :: fname
620 real(real64), allocatable :: ddip_elements(:,:,:)
621 complex(real64), allocatable :: zdip_elements(:,:,:)
622
623 push_sub(output_me_out_dipole)
624
625 ! The content of each file should be clear from the header of each file.
626 do ik = elec_me%states%d%kpt%start, elec_me%states%d%kpt%end
627 write(fname,'(i8)') ik
628 write(fname,'(a)') trim(dir)//'/ks_me_dipole.k'//trim(adjustl(fname))//'_'
629
630 if (states_are_real(elec_me%states)) then
631 ! Allocate dipole array
632 safe_allocate(ddip_elements(1:elec_me%space%dim, elec_me%st_start:elec_me%st_end, elec_me%st_start:elec_me%st_end))
633 ! Compute dipole elements
634 call elec_me%dipole_me(namespace, hm, ions, ik, ddip_elements)
635 ! Output elements
636 do idir = 1, elec_me%space%dim
637 iunit = io_open(trim(fname)//index2axis(idir), namespace, action = 'write')
638 write(iunit, '(a)') '# Dipole matrix elements file: |<Phi_i | r | Phi_j>|'
639 write(iunit, '(a,i8)') '# ik =', ik
640 write(iunit, '(a)') '# Units = ['//trim(units_abbrev(units_out%length))//']'
641
642 do ist = elec_me%st_start, elec_me%st_end
643 do jst = elec_me%st_start, elec_me%st_end
644 write(iunit, fmt='(f20.12)', advance = 'no') units_from_atomic(units_out%length, abs(ddip_elements(idir, ist, jst)))
645 write(iunit, fmt='(a)', advance = 'no') ' '
646 end do
647 write(iunit, '(a)') '' ! New line
648 end do
649 call io_close(iunit)
650 end do
651 safe_deallocate_a(ddip_elements)
652 else
653 ! Allocate dipole array
654 safe_allocate(zdip_elements(1:elec_me%space%dim, elec_me%st_start:elec_me%st_end, elec_me%st_start:elec_me%st_end))
655 ! Compute dipole elements
656 call elec_me%dipole_me(namespace, hm, ions, ik, zdip_elements)
657 ! Output elements
658 do idir = 1, elec_me%space%dim
659 iunit = io_open(trim(fname)//index2axis(idir), namespace, action = 'write')
660 write(iunit, '(a)') '# Dipole matrix elements file: |<Phi_i | r | Phi_j>|'
661 write(iunit, '(a,i8)') '# ik =', ik
662 write(iunit, '(a)') '# Units = ['//trim(units_abbrev(units_out%length))//']'
663
664 do ist = elec_me%st_start, elec_me%st_end
665 do jst = elec_me%st_start, elec_me%st_end
666 write(iunit, fmt='(f20.12)', advance = 'no') units_from_atomic(units_out%length, abs(zdip_elements(idir, ist, jst)))
667 write(iunit, fmt='(a)', advance = 'no') ' '
668 end do
669 write(iunit, '(a)') '' ! New line
670 end do
671 call io_close(iunit)
672 end do
673 safe_deallocate_a(zdip_elements)
674 end if
675 end do
676
677 pop_sub(output_me_out_dipole)
678 end subroutine output_me_out_dipole
679
680 subroutine output_me_out_ks_multipoles(dir, this, elec_me, namespace)
681 character(len=*), intent(in) :: dir
682 type(output_me_t), intent(in) :: this
683 type(elec_matrix_elements_t), intent(in) :: elec_me
684 type(namespace_t), intent(in) :: namespace
685
686 integer :: id, ll, mm, ik, iunit
687 character(len=256) :: fname
688 real(real64), allocatable :: delements(:,:)
689 complex(real64), allocatable :: zelements(:,:)
690
692
693 ! The content of each file should be clear from the header of each file.
694 id = 1
695 do ik = 1, elec_me%states%nik
696 select case (elec_me%space%dim)
697 case (3)
698 do ll = 1, this%ks_multipoles
699 do mm = -ll, ll
700 ! Define file name
701 write(fname,'(i4)') id
702 write(fname,'(a)') trim(dir) // '/ks_me_multipoles.' // trim(adjustl(fname))
703 iunit = io_open(fname, namespace, action = 'write')
704 ! Print header
705 call print_ks_multipole_header(iunit, ll, mm, ik)
706 ! Compute and print elements
707 if (states_are_real(elec_me%states)) then
708 safe_allocate(delements(1:elec_me%states%nst, 1:elec_me%states%nst))
709 call elec_me%ks_multipoles_3d(ll, mm, ik, delements)
710 call dprint_ks_multipole(iunit, elec_me%states%nst, delements, ll)
711 safe_deallocate_a(delements)
712 else
713 safe_allocate(zelements(1:elec_me%states%nst, 1:elec_me%states%nst))
714 call elec_me%ks_multipoles_3d(ll, mm, ik, zelements)
715 call zprint_ks_multipole(iunit, elec_me%states%nst, zelements, ll)
716 safe_deallocate_a(zelements)
717 end if
718 ! Close file
719 call io_close(iunit)
720 id = id + 1
721 end do
722 end do
723 case (2)
724 do ll = 1, 2
725 ! Define file name
726 write(fname,'(i4)') id
727 write(fname,'(a)') trim(dir) // '/ks_me_multipoles.' // trim(adjustl(fname))
728 iunit = io_open(fname, namespace, action = 'write')
729 ! Print header
730 call print_ks_multipole_header(iunit, ll, mm, ik)
731 ! Compute and print elements
732 if (states_are_real(elec_me%states)) then
733 safe_allocate(delements(1:elec_me%states%nst, 1:elec_me%states%nst))
734 call elec_me%ks_multipoles_2d(ll, ik, delements)
735 call dprint_ks_multipole(iunit, elec_me%states%nst, delements)
736 safe_deallocate_a(delements)
737 else
738 safe_allocate(zelements(1:elec_me%states%nst, 1:elec_me%states%nst))
739 call elec_me%ks_multipoles_2d(ll, ik, zelements)
740 call zprint_ks_multipole(iunit, elec_me%states%nst, zelements)
741 safe_deallocate_a(zelements)
742 end if
743 ! Close file
744 call io_close(iunit)
745 id = id + 1
746 end do
747 case (1)
748 do ll = 1, this%ks_multipoles
749 ! Define file name
750 write(fname,'(i4)') id
751 write(fname,'(a)') trim(dir) // '/ks_me_multipoles.' // trim(adjustl(fname))
752 iunit = io_open(fname, namespace, action = 'write')
753 ! Print header
754 call print_ks_multipole_header(iunit, ll, mm, ik)
755 ! Compute and print elements
756 if (states_are_real(elec_me%states)) then
757 safe_allocate(delements(1:elec_me%states%nst, 1:elec_me%states%nst))
758 call elec_me%ks_multipoles_1d(ll, ik, delements)
759 call dprint_ks_multipole(iunit, elec_me%states%nst, delements, ll)
760 safe_deallocate_a(delements)
761 else
762 safe_allocate(zelements(1:elec_me%states%nst, 1:elec_me%states%nst))
763 call elec_me%ks_multipoles_1d(ll, ik, zelements)
764 call zprint_ks_multipole(iunit, elec_me%states%nst, zelements, ll)
765 safe_deallocate_a(zelements)
766 end if
767 ! Close file
768 call io_close(iunit)
769 id = id + 1
770 end do
771 end select
772 end do
775
776 contains
777
778 subroutine print_ks_multipole_header(iunit, ll, mm, ik)
779 integer, intent(in) :: iunit, ll, mm, ik
780
781 write(iunit, fmt = '(a)') '# Multipole matrix elements file: <Phi_i | r**l * Y_{lm}(theta,phi) | Phi_j>'
782 write(iunit, fmt = '(a,i2,a,i2)') '# l =', ll, '; m =', mm
783 write(iunit, fmt = '(a,i8)') '# ik =', ik
784 if (ll>1) then
785 write(iunit, fmt = '(a,i1)') '# Units = ['//trim(units_abbrev(units_out%length))//']^', ll
786 else
787 write(iunit, fmt = '(a)') '# Units = ['//trim(units_abbrev(units_out%length))//']'
788 end if
789 end subroutine print_ks_multipole_header
790
791 subroutine dprint_ks_multipole(iunit, nst, elements, ll)
792 integer, intent(in) :: iunit, nst
793 real(real64), intent(in) :: elements(:,:)
794 integer, optional, intent(in) :: ll
795
796 integer :: ist, jst
797 real(real64) :: element
798
799 do ist = 1, nst
800 do jst = 1, nst
801 element = units_from_atomic(units_out%length**optional_default(ll, 1), elements(ist, jst))
802 write(iunit, fmt='(f20.12, a)', advance = 'no') element, ' '
803 end do
804 write(iunit, '(a)') ''
805 end do
806 end subroutine dprint_ks_multipole
807
808 subroutine zprint_ks_multipole(iunit, nst, elements, ll)
809 integer, intent(in) :: iunit, nst
810 complex(real64), intent(in) :: elements(:,:)
811 integer, optional, intent(in) :: ll
812
813 integer :: ist, jst
814 complex(real64) :: element
815
816 do ist = 1, nst
817 do jst = 1, nst
818 element = units_from_atomic(units_out%length**optional_default(ll, 1), elements(ist, jst))
819 write(iunit, fmt='(f20.12,a1,f20.12,a)', advance = 'no') real(element),',',aimag(element), ' '
820 end do
821 write(iunit, '(a)') ''
822 end do
823 end subroutine zprint_ks_multipole
824
825 end subroutine output_me_out_ks_multipoles
826
827
828#include "undef.F90"
829#include "real.F90"
830
831#include "undef.F90"
832#include "complex.F90"
833
834end module output_me_oct_m
835
836!! Local Variables:
837!! mode: f90
838!! coding: utf-8
839!! End:
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:468
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:395
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1125
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:624
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 output system
Definition: output_low.F90:115
subroutine output_me_out_momentum(fname, elec_me, namespace, kpoints)
Definition: output_me.F90:466
subroutine output_me_out_ks_multipoles(dir, this, elec_me, namespace)
Definition: output_me.F90:774
subroutine output_me_out_ang_momentum(fname, elec_me, namespace, kpoints)
Definition: output_me.F90:556
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:281
subroutine output_me_out_dipole(dir, elec_me, namespace, hm, ions)
Definition: output_me.F90:705
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
Definition: xc.F90:114
logical pure function, public family_is_mgga_with_exc(xcs)
Is the xc function part of the mGGA family with an energy functional.
Definition: xc.F90:592
subroutine print_ks_multipole_header(iunit, ll, mm, ik)
Definition: output_me.F90:872
subroutine zprint_ks_multipole(iunit, nst, elements, ll)
Definition: output_me.F90:902
subroutine dprint_ks_multipole(iunit, nst, elements, ll)
Definition: output_me.F90:885
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:168
The states_elec_t class contains all electronic wave functions.