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