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