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