Octopus
output.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
22module output_oct_m
23 use accel_oct_m
24 use basins_oct_m
25 use box_oct_m
26 use comm_oct_m
27 use cube_oct_m
30 use debug_oct_m
33 use dos_oct_m
35 use elf_oct_m
36#if defined(HAVE_ETSF_IO)
37 use etsf_io
38 use etsf_io_tools
39#endif
42 use fft_oct_m
45 use global_oct_m
46 use grid_oct_m
49 use io_oct_m
51 use ions_oct_m
52 use, intrinsic :: iso_fortran_env
55 use lasers_oct_m
56 use lda_u_oct_m
59 use loct_oct_m
60 use math_oct_m
62 use mesh_oct_m
65 use mpi_oct_m
71 use parser_oct_m
74 use smear_oct_m
75 use space_oct_m
80 use stress_oct_m
81 use string_oct_m
85 use unit_oct_m
87 use utils_oct_m
89 use v_ks_oct_m
90 use vtk_oct_m
91 use xc_oct_m
92 use xc_oep_oct_m
94 use xc_f03_lib_m
95
96 implicit none
97
98 private
99 public :: &
100 output_bgw_t, &
101 output_init, &
104 output_all, &
106 doutput_lr, &
107 zoutput_lr, &
111
112contains
113
114 subroutine output_init(outp, namespace, space, st, gr, nst, ks)
115 type(output_t), intent(out) :: outp
116 type(namespace_t), intent(in) :: namespace
117 class(space_t), intent(in) :: space
118 type(states_elec_t), intent(in) :: st
119 type(grid_t), intent(in) :: gr
120 integer, intent(in) :: nst
121 type(v_ks_t), intent(inout) :: ks
122
123 type(block_t) :: blk
124 real(real64) :: norm
125 character(len=80) :: nst_string, default
126#ifndef HAVE_ETSF_IO
127 integer :: iout
128#endif
129
130 push_sub(output_init)
131 outp%what = .false.
132
133 call io_function_read_what_how_when(namespace, space, outp%what, outp%how, outp%output_interval)
134
135 if (outp%what(option__output__wfs_fourier)) then
136 if (accel_is_enabled()) then
137 message(1) = "Wave functions in Fourier space not supported on GPUs."
138 call messages_fatal(1, namespace=namespace)
139 end if
140 call messages_experimental("Wave-functions in Fourier space", namespace=namespace)
141 end if
142
143 ! cannot calculate the ELF in 1D
144 if (outp%what(option__output__elf) .or. outp%what(option__output__elf_basins)) then
145 if (space%dim /= 2 .and. space%dim /= 3) then
146 outp%what(option__output__elf) = .false.
147 outp%what(option__output__elf_basins) = .false.
148 write(message(1), '(a)') 'Cannot calculate ELF except in 2D and 3D.'
149 call messages_warning(1, namespace=namespace)
150 end if
151 end if
152
153
154 if (outp%what(option__output__mmb_wfs)) then
155 call messages_experimental("Model many-body wfs", namespace=namespace)
156 end if
157
158 if (outp%what(option__output__xc_torque)) then
159 if (st%d%ispin /= spinors) then
160 write(message(1), '(a)') 'The output xc_torque can only be computed for spinors.'
161 call messages_fatal(1, namespace=namespace)
162 end if
163 if (space%dim /= 3) then
164 write(message(1), '(a)') 'The output xc_torque can only be computed in the 3D case.'
165 call messages_fatal(1, namespace=namespace)
166 end if
167 end if
168 if (outp%what(option__output__mmb_den)) then
169 call messages_experimental("Model many-body density matrix", namespace=namespace)
170 ! NOTES:
171 ! could be made into block to be able to specify which dimensions to trace
172 ! in principle all combinations are interesting, but this means we need to
173 ! be able to output density matrices for multiple particles or multiple
174 ! dimensions. The current 1D 1-particle case is simple.
175 end if
176
177 if (outp%what(option__output__energy_density)) call messages_experimental("'Output = energy_density'", namespace=namespace)
178 if (outp%what(option__output__heat_current)) call messages_experimental("'Output = heat_current'", namespace=namespace)
179
180 if (outp%what(option__output__wfs) .or. outp%what(option__output__wfs_sqmod)) then
181
182 !%Variable OutputWfsNumber
183 !%Type string
184 !%Default all states
185 !%Section Output
186 !%Description
187 !% Which wavefunctions to print, in list form: <i>i.e.</i>, "1-5" to print the first
188 !% five states, "2,3" to print the second and the third state, etc.
189 !% If more states are specified than available, extra ones will be ignored.
190 !%End
191
192 write(nst_string,'(i6)') nst
193 write(default,'(a,a)') "1-", trim(adjustl(nst_string))
194 call parse_variable(namespace, 'OutputWfsNumber', default, outp%wfs_list)
195 end if
196
197 if (parse_block(namespace, 'CurrentThroughPlane', blk) == 0) then
198 if (.not. outp%what(option__output__j_flow)) then
199 outp%what(option__output__j_flow) = .true.
200 call parse_variable(namespace, 'OutputInterval', 50, outp%output_interval(option__output__j_flow))
201 end if
202
203 !%Variable CurrentThroughPlane
204 !%Type block
205 !%Section Output
206 !%Description
207 !% The code can calculate current
208 !% traversing a user-defined portion of a plane, as specified by this block.
209 !% A small plain-text file <tt>current-flow</tt> will be written containing this information.
210 !% Only available for 1D, 2D, or 3D.
211 !% In the format below, <tt>origin</tt> is a point in the plane.
212 !% <tt>u</tt> and <tt>v</tt> are the (dimensionless) vectors defining the plane;
213 !% they will be normalized. <tt>spacing</tt> is the fineness of the mesh
214 !% on the plane. Integers <tt>nu</tt> and <tt>mu</tt> are the length and
215 !% width of the portion of the plane, in units of <tt>spacing</tt>.
216 !% Thus, the grid points included in the plane are
217 !% <tt>x_ij = origin + i*spacing*u + j*spacing*v</tt>,
218 !% for <tt>nu <= i <= mu</tt> and <tt>nv <= j <= mv</tt>.
219 !% Analogously, in the 2D case, the current flow is calculated through a line;
220 !% in the 1D case, the current flow is calculated through a point. Note that the spacing
221 !% can differ from the one used in the main calculation; an interpolation will be performed.
222 !%
223 !% Example (3D):
224 !%
225 !% <tt>%CurrentThroughPlane
226 !% <br>&nbsp;&nbsp; 0.0 | 0.0 | 0.0 # origin
227 !% <br>&nbsp;&nbsp; 0.0 | 1.0 | 0.0 # u
228 !% <br>&nbsp;&nbsp; 0.0 | 0.0 | 1.0 # v
229 !% <br>&nbsp;&nbsp; 0.2 # spacing
230 !% <br>&nbsp;&nbsp; 0 | 50 # nu | mu
231 !% <br>&nbsp;&nbsp; -50 | 50 # nv | mv
232 !% <br>%</tt>
233 !%
234 !% Example (2D):
235 !%
236 !% <tt>%CurrentThroughPlane
237 !% <br>&nbsp;&nbsp; 0.0 | 0.0 # origin
238 !% <br>&nbsp;&nbsp; 1.0 | 0.0 # u
239 !% <br>&nbsp;&nbsp; 0.2 # spacing
240 !% <br>&nbsp;&nbsp; 0 | 50 # nu | mu
241 !% <br>%</tt>
242 !%
243 !% Example (1D):
244 !%
245 !% <tt>%CurrentThroughPlane
246 !% <br>&nbsp;&nbsp; 0.0 # origin
247 !% <br>%</tt>
248 !%
249 !%End
250
251 select case (space%dim)
252 case (3)
253
254 call parse_block_float(blk, 0, 0, outp%plane%origin(1), units_inp%length)
255 call parse_block_float(blk, 0, 1, outp%plane%origin(2), units_inp%length)
256 call parse_block_float(blk, 0, 2, outp%plane%origin(3), units_inp%length)
257 call parse_block_float(blk, 1, 0, outp%plane%u(1))
258 call parse_block_float(blk, 1, 1, outp%plane%u(2))
259 call parse_block_float(blk, 1, 2, outp%plane%u(3))
260 call parse_block_float(blk, 2, 0, outp%plane%v(1))
261 call parse_block_float(blk, 2, 1, outp%plane%v(2))
262 call parse_block_float(blk, 2, 2, outp%plane%v(3))
263 call parse_block_float(blk, 3, 0, outp%plane%spacing, units_inp%length)
264 call parse_block_integer(blk, 4, 0, outp%plane%nu)
265 call parse_block_integer(blk, 4, 1, outp%plane%mu)
266 call parse_block_integer(blk, 5, 0, outp%plane%nv)
267 call parse_block_integer(blk, 5, 1, outp%plane%mv)
268
269 norm = norm2(outp%plane%u(1:3))
270 if (norm < m_epsilon) then
271 write(message(1), '(a)') 'u-vector for CurrentThroughPlane cannot have norm zero.'
272 call messages_fatal(1, namespace=namespace)
273 end if
274 outp%plane%u(1:3) = outp%plane%u(1:3) / norm
275
276 norm = norm2(outp%plane%v(1:3))
277 if (norm < m_epsilon) then
278 write(message(1), '(a)') 'v-vector for CurrentThroughPlane cannot have norm zero.'
279 call messages_fatal(1, namespace=namespace)
280 end if
281 outp%plane%v(1:3) = outp%plane%v(1:3) / norm
282
283 outp%plane%n(1:3) = dcross_product(outp%plane%u(1:3), outp%plane%v(1:3))
284
285 case (2)
286
287 call parse_block_float(blk, 0, 0, outp%line%origin(1), units_inp%length)
288 call parse_block_float(blk, 0, 1, outp%line%origin(2), units_inp%length)
289 call parse_block_float(blk, 1, 0, outp%line%u(1))
290 call parse_block_float(blk, 1, 1, outp%line%u(2))
291 call parse_block_float(blk, 2, 0, outp%line%spacing, units_inp%length)
292 call parse_block_integer(blk, 3, 0, outp%line%nu)
293 call parse_block_integer(blk, 3, 1, outp%line%mu)
294
295 norm = norm2(outp%line%u(1:2))
296 if (norm < m_epsilon) then
297 write(message(1), '(a)') 'u-vector for CurrentThroughPlane cannot have norm zero.'
298 call messages_fatal(1, namespace=namespace)
299 end if
300 outp%line%u(1:2) = outp%line%u(1:2) / norm
301
302 outp%line%n(1) = -outp%line%u(2)
303 outp%line%n(2) = outp%line%u(1)
304
305 case (1)
306
307 call parse_block_float(blk, 0, 0, outp%line%origin(1), units_inp%length)
308
309 case default
310
311 call messages_not_implemented("CurrentThroughPlane for 4D or higher", namespace=namespace)
312
313 end select
314 call parse_block_end(blk)
315 end if
316
317 if (outp%what(option__output__matrix_elements)) then
318 call output_me_init(outp%me, namespace, space, st, gr, nst)
319 else
320 outp%me%what = .false.
321 end if
322
323 if (outp%what(option__output__berkeleygw)) then
324 if (accel_is_enabled()) then
325 message(1) = "BerkeleyGW is not compatible with GPUs."
326 call messages_fatal(1, namespace=namespace)
327 end if
328 call output_berkeleygw_init(nst, namespace, outp%bgw, space%periodic_dim)
329 end if
330
331 ! required for output_hamiltonian()
332 if (outp%what(option__output__potential_gradient) .and. .not. outp%what(option__output__potential)) then
333 outp%what(option__output__potential) = .true.
334 outp%output_interval(option__output__potential) = outp%output_interval(option__output__potential_gradient)
335 end if
336
337
338 !%Variable OutputDuringSCF
339 !%Type logical
340 !%Default no
341 !%Section Output
342 !%Description
343 !% During <tt>gs</tt> and <tt>unocc</tt> runs, if this variable is set to yes,
344 !% output will be written after every <tt>OutputInterval</tt> iterations.
345 !%End
346 call parse_variable(namespace, 'OutputDuringSCF', .false., outp%duringscf)
347
348 !%Variable RestartWriteInterval
349 !%Type integer
350 !%Default 50
351 !%Section Execution::IO
352 !%Description
353 !% Restart data is written when the iteration number is a multiple
354 !% of the <tt>RestartWriteInterval</tt> variable. For
355 !% time-dependent runs this includes the update of the output
356 !% controlled by the <tt>TDOutput</tt> variable. (Other output is
357 !% controlled by <tt>OutputInterval</tt>.)
358 !%End
359 call parse_variable(namespace, 'RestartWriteInterval', 50, outp%restart_write_interval)
360 if (outp%restart_write_interval <= 0) then
361 message(1) = "RestartWriteInterval must be > 0."
362 call messages_fatal(1, namespace=namespace)
363 end if
364
365 !%Variable OutputIterDir
366 !%Default "output_iter"
367 !%Type string
368 !%Section Output
369 !%Description
370 !% The name of the directory where <tt>Octopus</tt> stores information
371 !% such as the density, forces, etc. requested by variable <tt>Output</tt>
372 !% in the format specified by <tt>OutputFormat</tt>.
373 !% This information is written while iterating <tt>CalculationMode = gs</tt>, <tt>unocc</tt>, or <tt>td</tt>,
374 !% according to <tt>OutputInterval</tt>, and has nothing to do with the restart information.
375 !%End
376 call parse_variable(namespace, 'OutputIterDir', "output_iter", outp%iter_dir)
377 if (any(outp%what) .and. any(outp%output_interval > 0)) then
378 call io_mkdir(outp%iter_dir, namespace)
379 end if
380 call add_last_slash(outp%iter_dir)
381
382 ! At this point, we don`t know whether the states will be real or complex.
383 ! We therefore pass .false. to states_are_real, and need to check for real states later.
384
385 if (output_needs_current(outp, .false.)) then
386 call v_ks_calculate_current(ks, .true.)
387 else
388 call v_ks_calculate_current(ks, .false.)
389 end if
390
391 if (outp%what(option__output__current_dia)) then
392 message(1) = "The diamagnetic current will be calculated only if CalculateDiamagneticCurrent = yes."
393 call messages_warning(1, namespace=namespace)
394 end if
395
396#ifndef HAVE_ETSF_IO
397 do iout = lbound(outp%how, 1), ubound(outp%how, 1)
398 if (bitand(outp%how(iout), option__outputformat__etsf) /= 0) then
399 message(1) = "ETSF output can only be used if Octopus is compiled with ETSF-IO"
400 call messages_fatal(1, namespace=namespace)
401 end if
402 end do
403#endif
404
405 pop_sub(output_init)
406 end subroutine output_init
407
408 ! ---------------------------------------------------------
409 subroutine output_all(outp, namespace, space, dir, gr, ions, iter, st, hm, ks)
410 type(output_t), intent(in) :: outp
411 type(namespace_t), intent(in) :: namespace
412 class(space_t), intent(in) :: space
413 character(len=*), intent(in) :: dir
414 type(grid_t), intent(in) :: gr
415 type(ions_t), intent(in) :: ions
416 integer, intent(in) :: iter
417 type(states_elec_t), intent(inout) :: st
418 type(hamiltonian_elec_t), intent(inout) :: hm
419 type(v_ks_t), intent(inout) :: ks
420
421 integer :: idir, ierr, iout, iunit
422 character(len=MAX_PATH_LEN) :: fname
423
424 push_sub(output_all)
425 call profiling_in("OUTPUT_ALL")
426
427 if (any(outp%what)) then
428 message(1) = "Info: Writing output to " // trim(dir)
429 call messages_info(1, namespace=namespace)
430 call io_mkdir(dir, namespace)
431 end if
432
433 if (outp%what_now(option__output__mesh_r, iter)) then
434 do idir = 1, space%dim
435 write(fname, '(a,a)') 'mesh_r-', index2axis(idir)
436 call dio_function_output(outp%how(option__output__mesh_r), dir, fname, namespace, space, &
437 gr, gr%x(:,idir), units_out%length, ierr, pos=ions%pos, atoms=ions%atom)
438 end do
439 end if
440
441 call output_states(outp, namespace, space, dir, st, gr, ions, hm, iter)
442 call output_hamiltonian(outp, namespace, space, dir, hm, st, gr%der, ions, gr, iter, st%st_kpt_mpi_grp)
443
444 ! We can only test against the theory level here, as it is not set in the call of output_init().
445 if (outp%what_now(option__output__el_pressure, iter)) then
446 if(ks%theory_level /= kohn_sham_dft) then
447 call messages_not_implemented("el_pressure for TheoryLevel different from kohn_sham", namespace=namespace)
448 end if
449 if (st%d%spin_channels > 1) then
450 call messages_not_implemented("el_pressure for spin-polarized or spinors", namespace=namespace)
451 end if
452 end if
453
454 !hm not initialized yet when calling output_init()
455 if (outp%what(option__output__kanamoriu) .and. hm%lda_u_level /= dft_u_acbn0) then
456 message(1) = "kanamoriU output can only be computed for DFTULevel = dft_u_acbn0"
457 call messages_fatal(1, namespace=namespace)
458 end if
459
460
461 call output_localization_funct(outp, namespace, space, dir, st, hm, gr, ions, iter)
462
463 if (outp%what_now(option__output__j_flow, iter)) then
464 call output_current_flow(outp, namespace, space, dir, gr, st, hm%kpoints)
465 end if
466
467 if (outp%what_now(option__output__geometry, iter)) then
468 if (bitand(outp%how(option__output__geometry), option__outputformat__xcrysden) /= 0) then
469 call write_xsf_geometry_file(dir, "geometry", space, ions%latt, ions%pos, ions%atom, gr, namespace)
470 end if
471 if (bitand(outp%how(option__output__geometry), option__outputformat__xyz) /= 0) then
472 call ions%write_xyz(trim(dir)//'/geometry')
473 if (ions%space%is_periodic()) then
474 call ions%write_crystal(dir)
475 end if
476 end if
477 if (bitand(outp%how(option__output__geometry), option__outputformat__vtk) /= 0) then
478 call ions%write_vtk_geometry(trim(dir)//'/geometry')
479 end if
480 end if
481
482 if (outp%what_now(option__output__forces, iter)) then
483 if (bitand(outp%how(option__output__forces), option__outputformat__bild) /= 0) then
484 call ions%write_bild_forces_file(dir, "forces")
485 else
486 call write_xsf_geometry_file(dir, "forces", space, ions%latt, ions%pos, ions%atom, &
487 gr, namespace, total_forces = ions%tot_force)
488 end if
489 end if
490
491 if (outp%what_now(option__output__matrix_elements, iter)) then
492 call output_me(outp%me, namespace, space, dir, st, gr, ions, hm)
493 end if
494
495 do iout = lbound(outp%how, 1), ubound(outp%how, 1)
496 if (bitand(outp%how(iout), option__outputformat__etsf) /= 0) then
497 call output_etsf(outp, namespace, space, dir, st, gr, hm%kpoints, ions, iter)
498 exit
499 end if
500 end do
501
502 if (outp%what_now(option__output__berkeleygw, iter)) then
503 call output_berkeleygw(outp%bgw, namespace, space, dir, st, gr, ks, hm, ions)
504 end if
505
506 if (outp%what_now(option__output__energy_density, iter)) then
507 call output_energy_density(outp, namespace, space, dir, hm, ks, st, ions, gr)
508 end if
509
510 if (outp%what_now(option__output__stress, iter)) then
511 call io_mkdir(dir, namespace)
512 iunit = io_open(trim(dir)//'/stress', namespace, action='write')
513 call output_stress(iunit, space%periodic_dim, st%stress_tensors)
514 call io_close(iunit)
515 end if
516
517 if (hm%lda_u_level /= dft_u_none) then
518 if (outp%what_now(option__output__occ_matrices, iter))&
519 call lda_u_write_occupation_matrices(dir, hm%lda_u, st, namespace)
520
521 if (outp%what_now(option__output__effectiveu, iter))&
522 call lda_u_write_effectiveu(dir, hm%lda_u, namespace)
523
524 if (outp%what_now(option__output__magnetization, iter))&
525 call lda_u_write_magnetization(dir, hm%lda_u, ions, gr, st, namespace)
526
527 if (outp%what_now(option__output__local_orbitals, iter))&
528 call output_dftu_orbitals(outp, dir, namespace, space, hm%lda_u, st, gr, ions, hm%phase%is_allocated())
529
530 if (outp%what_now(option__output__kanamoriu, iter))&
531 call lda_u_write_kanamoriu(dir, st, hm%lda_u, namespace)
532
533 if (ks%oep_photon%level == oep_level_full) then
534 if (outp%what_now(option__output__photon_correlator, iter)) then
535 write(fname, '(a)') 'photon_correlator'
536 call dio_function_output(outp%how(option__output__photon_correlator), dir, trim(fname), namespace, space, &
537 gr, ks%oep_photon%pt%correlator(:,1), units_out%length, ierr, pos=ions%pos, atoms=ions%atom)
538 end if
539 end if
540 end if
541
542 ! We can only test against the theory level here, as it is not set in the call of output_init().
543 if (outp%what_now(option__output__xc_torque, iter)) then
544 if (ks%theory_level /= kohn_sham_dft .and. ks%theory_level /= generalized_kohn_sham_dft) then
545 write(message(1), '(a)') 'The output xc_torque can only be computed when there is a xc potential.'
546 call messages_fatal(1, namespace=namespace)
547 end if
548 end if
549
550 call output_xc_torque(outp, namespace, dir, gr, hm, st, ions, ions%space)
551
552 call profiling_out("OUTPUT_ALL")
553 pop_sub(output_all)
554 end subroutine output_all
555
556
557 ! ---------------------------------------------------------
558 subroutine output_localization_funct(outp, namespace, space, dir, st, hm, gr, ions, iter)
559 type(output_t), intent(in) :: outp
560 type(namespace_t), intent(in) :: namespace
561 class(space_t), intent(in) :: space
562 character(len=*), intent(in) :: dir
563 type(states_elec_t), intent(inout) :: st
564 type(hamiltonian_elec_t), intent(in) :: hm
565 type(grid_t), intent(in) :: gr
566 type(ions_t), intent(in) :: ions
567 integer, intent(in) :: iter
568
569 real(real64), allocatable :: f_loc(:,:)
570 character(len=MAX_PATH_LEN) :: fname
571 integer :: is, ierr, imax
572 type(mpi_grp_t) :: mpi_grp
573
575
576 mpi_grp = st%dom_st_kpt_mpi_grp
577
578 ! if SPIN_POLARIZED, the ELF contains one extra channel: the total ELF
579 imax = st%d%nspin
580 if (st%d%ispin == spin_polarized) imax = 3
581
582 safe_allocate(f_loc(1:gr%np, 1:imax))
583
584 ! First the ELF in real space
585 if (outp%what_now(option__output__elf, iter) .or. outp%what_now(option__output__elf_basins, iter)) then
586 assert(space%dim /= 1)
587
588 call elf_calc(space, st, gr, hm%kpoints, f_loc)
589
590 ! output ELF in real space
591 if (outp%what_now(option__output__elf, iter)) then
592 write(fname, '(a)') 'elf_rs'
593 call dio_function_output(outp%how(option__output__elf), dir, trim(fname), namespace, space, gr, &
594 f_loc(:,imax), unit_one, ierr, pos=ions%pos, atoms=ions%atom, grp = mpi_grp)
595 ! this quantity is dimensionless
596
597 if (st%d%ispin /= unpolarized) then
598 do is = 1, 2
599 write(fname, '(a,i1)') 'elf_rs-sp', is
600 call dio_function_output(outp%how(option__output__elf), dir, trim(fname), namespace, space, gr, &
601 f_loc(:, is), unit_one, ierr, pos=ions%pos, atoms=ions%atom, grp = mpi_grp)
602 ! this quantity is dimensionless
603 end do
604 end if
605 end if
606
607 if (outp%what_now(option__output__elf_basins, iter)) then
608 call out_basins(f_loc(:,1), "elf_rs_basins", outp%how(option__output__elf_basins))
609 end if
610 end if
611
612 ! Now Bader analysis
613 if (outp%what_now(option__output__bader, iter)) then
614 do is = 1, st%d%nspin
615 call dderivatives_lapl(gr%der, st%rho(:,is), f_loc(:,is))
616
617 fname = get_filename_with_spin('bader', st%d%nspin, is)
618
619 call dio_function_output(outp%how(option__output__bader), dir, trim(fname), namespace, space, gr, &
620 f_loc(:,is), units_out%length**(-2 - space%dim), ierr, &
621 pos=ions%pos, atoms=ions%atom, grp = mpi_grp)
622
623 fname = get_filename_with_spin('bader_basins', st%d%nspin, is)
624 call out_basins(f_loc(:, is), fname, outp%how(option__output__bader))
625 end do
626 end if
627
628 ! Now the pressure
629 if (outp%what_now(option__output__el_pressure, iter)) then
630 call calc_electronic_pressure(st, hm, gr, f_loc(:,1))
631 call dio_function_output(outp%how(option__output__el_pressure), dir, "el_pressure", namespace, space, gr, &
632 f_loc(:,1), unit_one, ierr, pos=ions%pos, atoms=ions%atom, grp = mpi_grp)
633 ! this quantity is dimensionless
634 end if
635
636 safe_deallocate_a(f_loc)
637
639
640 contains
641 ! ---------------------------------------------------------
642 subroutine out_basins(ff, filename, output_how)
643 real(real64), intent(in) :: ff(:)
644 character(len=*), intent(in) :: filename
645 integer(int64), intent(in) :: output_how
646
647 character(len=MAX_PATH_LEN) :: fname
648 type(basins_t) :: basins
649 integer :: iunit
650
652
653 call basins_init(basins, namespace, gr)
654 call basins_analyze(basins, namespace, gr, ff(:), st%rho, 0.01_real64)
655
656 call dio_function_output(output_how, dir, trim(filename), namespace, space, gr, &
657 real(basins%map, real64) , unit_one, ierr, pos=ions%pos, atoms=ions%atom, grp = mpi_grp)
658 ! this quantity is dimensionless
659
660 write(fname,'(4a)') trim(dir), '/', trim(filename), '.info'
661 iunit = io_open(trim(fname), namespace, action = 'write')
662 call basins_write(basins, gr, iunit)
663 call io_close(iunit)
664
665 call basins_end(basins)
666
668 end subroutine out_basins
669
670 end subroutine output_localization_funct
671
672
673 ! ---------------------------------------------------------
674 subroutine calc_electronic_pressure(st, hm, gr, pressure)
675 type(states_elec_t), intent(inout) :: st
676 type(hamiltonian_elec_t), intent(in) :: hm
677 type(grid_t), intent(in) :: gr
678 real(real64), intent(out) :: pressure(:)
679
680 real(real64), allocatable :: rho(:,:), lrho(:), tau(:,:)
681 real(real64) :: p_tf, dens
682 integer :: is, ii
683
685
686 safe_allocate( rho(1:gr%np_part, 1:st%d%nspin))
687 safe_allocate(lrho(1:gr%np))
688 safe_allocate( tau(1:gr%np, 1:st%d%nspin))
689
690 rho = m_zero
691 call density_calc(st, gr, rho)
692 call states_elec_calc_quantities(gr, st, hm%kpoints, .false., kinetic_energy_density = tau)
693
694 pressure = m_zero
695 do is = 1, st%d%spin_channels
696 lrho = m_zero
697 call dderivatives_lapl(gr%der, rho(:, is), lrho)
698
699 pressure(:) = pressure(:) + &
700 tau(:, is)/m_three - lrho(:)/m_four
701 end do
702
703 do ii = 1, gr%np
704 dens = sum(rho(ii,1:st%d%spin_channels))
705
706 p_tf = m_two/m_five*(m_three*m_pi**2)**(m_two/m_three)* &
707 dens**(m_five/m_three)
708
709 ! add XC pressure
710 ! FIXME: Not correct for spinors and spin-polarized here
711 pressure(ii) = pressure(ii) + (dens*hm%ks_pot%vxc(ii,1) - hm%energy%exchange - hm%energy%correlation)
712
713 pressure(ii) = pressure(ii)/p_tf
714 pressure(ii) = m_half*(m_one + pressure(ii)/sqrt(m_one + pressure(ii)**2))
715 end do
716
718 end subroutine calc_electronic_pressure
719
720
721 ! ---------------------------------------------------------
722 subroutine output_energy_density(outp, namespace, space, dir, hm, ks, st, ions, gr)
723 type(output_t), intent(in) :: outp
724 type(namespace_t), intent(in) :: namespace
725 class(space_t), intent(in) :: space
726 character(len=*), intent(in) :: dir
727 type(hamiltonian_elec_t), intent(in) :: hm
728 type(v_ks_t), intent(inout) :: ks
729 type(states_elec_t), intent(in) :: st
730 type(ions_t), intent(in) :: ions
731 type(grid_t), intent(in) :: gr
732
733 integer :: is, ierr, ip
734 character(len=MAX_PATH_LEN) :: fname
735 type(unit_t) :: fn_unit
736 real(real64), allocatable :: energy_density(:, :)
737 real(real64), allocatable :: ex_density(:)
738 real(real64), allocatable :: ec_density(:)
739
740 push_sub(output_energy_density)
741
742 fn_unit = units_out%energy*units_out%length**(-space%dim)
743 safe_allocate(energy_density(1:gr%np, 1:st%d%nspin))
744
745 ! the kinetic energy density
746 call states_elec_calc_quantities(gr, st, hm%kpoints, .true., kinetic_energy_density = energy_density)
747
748 ! the external potential energy density
749 do is = 1, st%d%nspin
750 do ip = 1, gr%np
751 energy_density(ip, is) = energy_density(ip, is) + st%rho(ip, is)*hm%ep%vpsl(ip)
752 end do
753 end do
754
755 ! the hartree energy density
756 do is = 1, st%d%nspin
757 do ip = 1, gr%np
758 energy_density(ip, is) = energy_density(ip, is) + m_half*st%rho(ip, is)*hm%ks_pot%vhartree(ip)
759 end do
760 end do
761
762 ! the XC energy density
763 safe_allocate(ex_density(1:gr%np))
764 safe_allocate(ec_density(1:gr%np))
765
766 call xc_get_vxc(gr, ks%xc, st, hm%kpoints, hm%psolver, namespace, space, st%rho, st%d%ispin, &
767 hm%ions%latt%rcell_volume, ex_density = ex_density, ec_density = ec_density)
768 do is = 1, st%d%nspin
769 do ip = 1, gr%np
770 energy_density(ip, is) = energy_density(ip, is) + ex_density(ip) + ec_density(ip)
771 end do
772 end do
773
774 safe_deallocate_a(ex_density)
775 safe_deallocate_a(ec_density)
776
777 do is = 1, st%d%spin_channels
778 fname = get_filename_with_spin('energy_density', st%d%nspin, is)
779 call dio_function_output(outp%how(option__output__energy_density), dir, trim(fname), namespace, space, gr, &
780 energy_density(:, is), unit_one, ierr, pos=ions%pos, atoms=ions%atom, grp = st%dom_st_kpt_mpi_grp)
781 end do
782 safe_deallocate_a(energy_density)
783
784 pop_sub(output_energy_density)
785 end subroutine output_energy_density
786
787 !--------------------------------------------------------------
788
789 logical function output_need_exchange(outp) result(need_exx)
790 type(output_t), intent(in) :: outp
791
792 need_exx =(outp%what(option__output__berkeleygw) &
793 .or. outp%me%what(option__outputmatrixelements__two_body) &
794 .or. outp%me%what(option__outputmatrixelements__two_body_exc_k))
795 end function output_need_exchange
796
797
798 ! ---------------------------------------------------------
799 subroutine output_dftu_orbitals(outp, dir, namespace, space, this, st, mesh, ions, has_phase)
800 type(output_t), intent(in) :: outp
801 character(len=*), intent(in) :: dir
802 type(namespace_t), intent(in) :: namespace
803 class(space_t), intent(in) :: space
804 type(lda_u_t), intent(in) :: this
805 type(states_elec_t), intent(in) :: st
806 class(mesh_t), intent(in) :: mesh
807 type(ions_t), intent(in) :: ions
808 logical, intent(in) :: has_phase
809
810 integer :: ios, im, ik, idim, ierr
811 complex(real64), allocatable :: tmp(:)
812 real(real64), allocatable :: dtmp(:)
813 type(orbitalset_t), pointer :: os
814 type(unit_t) :: fn_unit
815 character(len=MAX_PATH_LEN) :: fname
816
817 push_sub(output_dftu_orbitals)
818
819 fn_unit = sqrt(units_out%length**(-space%dim))
820
821 if (this%basis%use_submesh) then
822 if (states_are_real(st)) then
823 safe_allocate(dtmp(1:mesh%np))
824 else
825 safe_allocate(tmp(1:mesh%np))
826 end if
827 end if
828
829 do ios = 1, this%norbsets
830 os => this%orbsets(ios)
831 do ik = st%d%kpt%start, st%d%kpt%end
832 do im = 1, this%orbsets(ios)%norbs
833 do idim = 1, min(os%ndim, st%d%dim)
834 if (st%nik > 1) then
835 if (min(os%ndim, st%d%dim) > 1) then
836 write(fname, '(a,i1,a,i3.3,a,i8.8,a,i1)') 'orb', im, '-os', ios, '-k', ik, '-sp', idim
837 else
838 write(fname, '(a,i1,a,i3.3,a,i8.8)') 'orb', im, '-os', ios, '-k', ik
839 end if
840 else
841 if (min(os%ndim, st%d%dim) > 1) then
842 write(fname, '(a,i1,a,i3.3,a,i1)') 'orb', im, '-os', ios, '-sp', idim
843 else
844 write(fname, '(a,i1,a,i3.3)') 'orb', im, '-os', ios
845 end if
846 end if
847 if (has_phase) then
848 if (.not. this%basis%use_submesh) then
849 call zio_function_output(outp%how(option__output__local_orbitals), dir, fname, namespace, space, &
850 mesh, os%eorb_mesh(1:mesh%np,im,idim,ik), fn_unit, ierr, pos=ions%pos, atoms=ions%atom)
851 else
852 tmp = m_z0
853 call submesh_add_to_mesh(os%sphere, os%eorb_submesh(1:os%sphere%np,idim,im,ik), tmp)
854 call zio_function_output(outp%how(option__output__local_orbitals), dir, fname, namespace, space, &
855 mesh, tmp, fn_unit, ierr, pos=ions%pos, atoms=ions%atom)
856 end if
857 else
858 if (.not. this%basis%use_submesh) then
859 if (states_are_real(st)) then
860 call dio_function_output(outp%how(option__output__local_orbitals), dir, fname, namespace, space, mesh, &
861 os%dorb(1:mesh%np,idim,im), fn_unit, ierr, pos=ions%pos, atoms=ions%atom)
862 else
863 call zio_function_output(outp%how(option__output__local_orbitals), dir, fname, namespace, space, mesh, &
864 os%zorb(1:mesh%np,idim,im), fn_unit, ierr, pos=ions%pos, atoms=ions%atom)
865 end if
866 else
867 if (states_are_real(st)) then
868 dtmp = m_zero
869 call submesh_add_to_mesh(os%sphere, os%dorb(1:os%sphere%np,idim,im), dtmp)
870 call dio_function_output(outp%how(option__output__local_orbitals), dir, fname, namespace, space, &
871 mesh, dtmp, fn_unit, ierr, pos=ions%pos, atoms=ions%atom)
872 else
873 tmp = m_z0
874 call submesh_add_to_mesh(os%sphere, os%zorb(1:os%sphere%np,idim,im), tmp)
875 call zio_function_output(outp%how(option__output__local_orbitals), dir, fname, namespace, space, &
876 mesh, tmp, fn_unit, ierr, pos=ions%pos, atoms=ions%atom)
877 end if
878 end if
879 end if
880 end do
881 end do
882 end do
883 end do
884
885 safe_deallocate_a(tmp)
886 safe_deallocate_a(dtmp)
887
888 pop_sub(output_dftu_orbitals)
889 end subroutine output_dftu_orbitals
890
891 ! ---------------------------------------------------------
892 logical function output_needs_current(outp, states_are_real)
893 type(output_t), intent(in) :: outp
894 logical, intent(in) :: states_are_real
895
896 output_needs_current = .false.
897
898 if (outp%what(option__output__current) &
899 .or. outp%what(option__output__current_dia) &
900 .or. outp%what(option__output__heat_current) &
901 .or. outp%what(option__output__current_kpt)) then
902 if (.not. states_are_real) then
904 else
905 message(1) = 'No current density output for real states since it is identically zero.'
906 call messages_warning(1)
907 end if
908 end if
909
910
911 end function
912
913#include "output_etsf_inc.F90"
914
915#include "output_states_inc.F90"
916
917#include "output_h_inc.F90"
918
919#include "undef.F90"
920#include "complex.F90"
921#include "output_linear_response_inc.F90"
922
923#include "undef.F90"
924#include "real.F90"
925#include "output_linear_response_inc.F90"
926
927end module output_oct_m
928
929!! Local Variables:
930!! mode: f90
931!! coding: utf-8
932!! End:
pure logical function, public accel_is_enabled()
Definition: accel.F90:401
subroutine, public basins_write(this, mesh, iunit)
Definition: basins.F90:352
subroutine, public basins_init(this, namespace, mesh)
Definition: basins.F90:151
subroutine, public basins_analyze(this, namespace, mesh, f, rho, threshold)
Definition: basins.F90:188
subroutine, public basins_end(this)
Definition: basins.F90:170
This module implements a calculator for the density and defines related functions.
Definition: density.F90:120
subroutine, public density_calc(st, gr, density, istin)
Computes the density from the orbitals in st.
Definition: density.F90:612
This module calculates the derivatives (gradients, Laplacians, etc.) of a function.
subroutine, public dderivatives_lapl(der, ff, op_ff, ghost_update, set_bc, factor)
apply the Laplacian to a mesh function
Module that handles computing and output of various density of states.
Definition: dos.F90:116
integer, parameter, public unpolarized
Parameters...
integer, parameter, public spinors
integer, parameter, public spin_polarized
subroutine, public elf_calc(space, st, gr, kpoints, elf, de)
(time-dependent) electron localization function, (TD)ELF.
Definition: elf.F90:167
Fast Fourier Transform module. This module provides a single interface that works with different FFT ...
Definition: fft.F90:118
real(real64), parameter, public m_two
Definition: global.F90:190
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public m_four
Definition: global.F90:192
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:186
complex(real64), parameter, public m_z0
Definition: global.F90:198
real(real64), parameter, public m_epsilon
Definition: global.F90:204
real(real64), parameter, public m_half
Definition: global.F90:194
real(real64), parameter, public m_one
Definition: global.F90:189
real(real64), parameter, public m_three
Definition: global.F90:191
real(real64), parameter, public m_five
Definition: global.F90:193
This module implements the underlying real-space grid.
Definition: grid.F90:117
This module defines classes and functions for interaction partners.
subroutine, public zio_function_output(how, dir, fname, namespace, space, mesh, ff, unit, ierr, pos, atoms, grp, root)
Top-level IO routine for functions defined on the mesh.
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)
subroutine, public dio_function_output(how, dir, fname, namespace, space, mesh, ff, unit, ierr, pos, atoms, grp, root)
Top-level IO routine for functions defined on the mesh.
subroutine, public write_xsf_geometry_file(dir, fname, space, latt, pos, atoms, mesh, namespace, total_forces)
Definition: io.F90:114
subroutine, public io_close(iunit, grp)
Definition: io.F90:418
subroutine, public io_mkdir(fname, namespace, parents)
Definition: io.F90:311
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:352
A module to handle KS potential, without the external potential.
integer, parameter, public generalized_kohn_sham_dft
integer, parameter, public kohn_sham_dft
subroutine, public lda_u_write_occupation_matrices(dir, this, st, namespace)
Prints the occupation matrices at the end of the scf calculation.
Definition: lda_u_io.F90:157
subroutine, public lda_u_write_kanamoriu(dir, st, this, namespace)
Definition: lda_u_io.F90:350
subroutine, public lda_u_write_magnetization(dir, this, ions, mesh, st, namespace)
Definition: lda_u_io.F90:478
subroutine, public lda_u_write_effectiveu(dir, this, namespace)
Definition: lda_u_io.F90:265
integer, parameter, public dft_u_none
Definition: lda_u.F90:201
integer, parameter, public dft_u_acbn0
Definition: lda_u.F90:201
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
pure real(real64) function, dimension(1:3), public dcross_product(a, b)
Definition: math.F90:1918
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_not_implemented(feature, namespace)
Definition: messages.F90:1114
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:538
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:161
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:415
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1086
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:617
subroutine, public output_berkeleygw(bgw, namespace, space, dir, st, gr, ks, hm, ions)
subroutine, public output_berkeleygw_init(nst, namespace, bgw, periodic_dim)
this module contains the low-level part of the output system
Definition: output_low.F90:115
character(len=max_path_len) function, public get_filename_with_spin(output, nspin, spin_index)
Returns the filame as output, or output-spX is spin polarized.
Definition: output_low.F90:232
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
this module contains the output system
Definition: output.F90:115
subroutine calc_electronic_pressure(st, hm, gr, pressure)
Definition: output.F90:759
subroutine, public output_states(outp, namespace, space, dir, st, gr, ions, hm, iter)
Definition: output.F90:1899
logical function, public output_needs_current(outp, states_are_real)
Definition: output.F90:977
subroutine, public output_hamiltonian(outp, namespace, space, dir, hm, st, der, ions, gr, iter, grp)
Definition: output.F90:2197
subroutine, public output_all(outp, namespace, space, dir, gr, ions, iter, st, hm, ks)
Definition: output.F90:494
logical function, public output_need_exchange(outp)
Definition: output.F90:874
subroutine, public output_init(outp, namespace, space, st, gr, nst, ks)
Definition: output.F90:208
subroutine output_xc_torque(outp, namespace, dir, mesh, hm, st, ions, space)
Definition: output.F90:2515
subroutine, public output_current_flow(outp, namespace, space, dir, gr, st, kpoints)
Definition: output.F90:2087
subroutine, public zoutput_lr(outp, namespace, space, dir, st, mesh, lr, idir, isigma, ions, pert_unit)
Definition: output.F90:2622
subroutine, public doutput_lr(outp, namespace, space, dir, st, mesh, lr, idir, isigma, ions, pert_unit)
Definition: output.F90:2891
subroutine, public output_scalar_pot(outp, namespace, space, dir, mesh, ions, ext_partners, time)
Definition: output.F90:2478
subroutine output_energy_density(outp, namespace, space, dir, hm, ks, st, ions, gr)
Definition: output.F90:807
subroutine output_dftu_orbitals(outp, dir, namespace, space, this, st, mesh, ions, has_phase)
Definition: output.F90:884
subroutine output_localization_funct(outp, namespace, space, dir, st, hm, gr, ions, iter)
Definition: output.F90:643
subroutine output_etsf(outp, namespace, space, dir, st, gr, kpoints, ions, iter)
To create an etsf file one has to do the following:
Definition: output.F90:1028
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:623
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:552
pure logical function, public states_are_real(st)
This module defines routines to write information about states.
subroutine, public states_elec_calc_quantities(gr, st, kpoints, nlcc, kinetic_energy_density, paramagnetic_current, density_gradient, density_laplacian, gi_kinetic_energy_density, st_end)
calculated selected quantities
This module implements the calculation of the stress tensor.
Definition: stress.F90:118
subroutine, public output_stress(iunit, space_dim, stress_tensors, all_terms)
Definition: stress.F90:1075
subroutine, public add_last_slash(str)
Adds a '/' in the end of the string, only if it missing. Useful for directories.
Definition: string.F90:160
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
This module defines the unit system, used for input and output.
type(unit_system_t), public units_out
type(unit_system_t), public units_inp
the units systems for reading and writing
type(unit_t), public unit_one
some special units required for particular quantities
This module is intended to contain simple general-purpose utility functions and procedures.
Definition: utils.F90:118
character pure function, public index2axis(idir)
Definition: utils.F90:202
subroutine, public v_ks_calculate_current(this, calc_cur)
Definition: v_ks.F90:1454
Definition: xc.F90:114
subroutine, public xc_get_vxc(gr, xcs, st, kpoints, psolver, namespace, space, rho, ispin, rcell_volume, vxc, ex, ec, deltaxc, vtau, ex_density, ec_density, stress_xc, force_orbitalfree)
Definition: xc.F90:757
integer, parameter, public oep_level_full
Definition: xc_oep.F90:171
subroutine out_basins(ff, filename, output_how)
Definition: output.F90:727
Description of the grid, containing information on derivatives, stencil, and symmetries.
Definition: grid.F90:168
Class to describe DFT+U parameters.
Definition: lda_u.F90:214
Describes mesh distribution to nodes.
Definition: mesh.F90:186
This is defined even when running serial.
Definition: mpi.F90:142
Output information for BerkeleyGW.
Definition: output_low.F90:145
output handler class
Definition: output_low.F90:164
The states_elec_t class contains all electronic wave functions.
int true(void)