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) = outp%plane%u(2)*outp%plane%v(3) - outp%plane%u(3)*outp%plane%v(2)
285 outp%plane%n(2) = outp%plane%u(3)*outp%plane%v(1) - outp%plane%u(1)*outp%plane%v(3)
286 outp%plane%n(3) = outp%plane%u(1)*outp%plane%v(2) - outp%plane%u(2)*outp%plane%v(1)
287
288 case (2)
289
290 call parse_block_float(blk, 0, 0, outp%line%origin(1), units_inp%length)
291 call parse_block_float(blk, 0, 1, outp%line%origin(2), units_inp%length)
292 call parse_block_float(blk, 1, 0, outp%line%u(1))
293 call parse_block_float(blk, 1, 1, outp%line%u(2))
294 call parse_block_float(blk, 2, 0, outp%line%spacing, units_inp%length)
295 call parse_block_integer(blk, 3, 0, outp%line%nu)
296 call parse_block_integer(blk, 3, 1, outp%line%mu)
297
298 norm = norm2(outp%line%u(1:2))
299 if (norm < m_epsilon) then
300 write(message(1), '(a)') 'u-vector for CurrentThroughPlane cannot have norm zero.'
301 call messages_fatal(1, namespace=namespace)
302 end if
303 outp%line%u(1:2) = outp%line%u(1:2) / norm
304
305 outp%line%n(1) = -outp%line%u(2)
306 outp%line%n(2) = outp%line%u(1)
307
308 case (1)
309
310 call parse_block_float(blk, 0, 0, outp%line%origin(1), units_inp%length)
311
312 case default
313
314 call messages_not_implemented("CurrentThroughPlane for 4D or higher", namespace=namespace)
315
316 end select
317 call parse_block_end(blk)
318 end if
319
320 if (outp%what(option__output__matrix_elements)) then
321 call output_me_init(outp%me, namespace, space, st, gr, nst)
322 else
323 outp%me%what = .false.
324 end if
325
326 if (outp%what(option__output__berkeleygw)) then
327 if (accel_is_enabled()) then
328 message(1) = "BerkeleyGW is not compatible with GPUs."
329 call messages_fatal(1, namespace=namespace)
330 end if
331 call output_berkeleygw_init(nst, namespace, outp%bgw, space%periodic_dim)
332 end if
333
334 ! required for output_hamiltonian()
335 if (outp%what(option__output__potential_gradient) .and. .not. outp%what(option__output__potential)) then
336 outp%what(option__output__potential) = .true.
337 outp%output_interval(option__output__potential) = outp%output_interval(option__output__potential_gradient)
338 end if
339
340
341 !%Variable OutputDuringSCF
342 !%Type logical
343 !%Default no
344 !%Section Output
345 !%Description
346 !% During <tt>gs</tt> and <tt>unocc</tt> runs, if this variable is set to yes,
347 !% output will be written after every <tt>OutputInterval</tt> iterations.
348 !%End
349 call parse_variable(namespace, 'OutputDuringSCF', .false., outp%duringscf)
350
351 if (parse_is_defined(namespace, 'RestartWriteInterval')) then
352 write(message(1), '(a)') 'Input variable RestartWriteInterval is obsolete.'
353 write(message(2), '(a)') 'Restart files are now written in periods of the wallclock time'
354 write(message(3), '(a)') 'given by RestartWallTimePeriod, so you can simply delete this variable.'
355 call messages_fatal(3, only_root_writes=.true., namespace=namespace)
356 end if
357
358 !%Variable OutputIterDir
359 !%Default "output_iter"
360 !%Type string
361 !%Section Output
362 !%Description
363 !% The name of the directory where <tt>Octopus</tt> stores information
364 !% such as the density, forces, etc. requested by variable <tt>Output</tt>
365 !% in the format specified by <tt>OutputFormat</tt>.
366 !% This information is written while iterating <tt>CalculationMode = gs</tt>, <tt>unocc</tt>, or <tt>td</tt>,
367 !% according to <tt>OutputInterval</tt>, and has nothing to do with the restart information.
368 !%End
369 call parse_variable(namespace, 'OutputIterDir', "output_iter", outp%iter_dir)
370 if (any(outp%what) .and. any(outp%output_interval > 0)) then
371 call io_mkdir(outp%iter_dir, namespace)
372 end if
373 call add_last_slash(outp%iter_dir)
374
375 ! At this point, we don`t know whether the states will be real or complex.
376 ! We therefore pass .false. to states_are_real, and need to check for real states later.
377
378 if (output_needs_current(outp, .false.)) then
379 call v_ks_calculate_current(ks, .true.)
380 else
381 call v_ks_calculate_current(ks, .false.)
382 end if
383
384 if (outp%what(option__output__current_dia)) then
385 message(1) = "The diamagnetic current will be calculated only if CalculateDiamagneticCurrent = yes."
386 call messages_warning(1, namespace=namespace)
387 end if
388
389#ifndef HAVE_ETSF_IO
390 do iout = lbound(outp%how, 1), ubound(outp%how, 1)
391 if (bitand(outp%how(iout), option__outputformat__etsf) /= 0) then
392 message(1) = "ETSF output can only be used if Octopus is compiled with ETSF-IO"
393 call messages_fatal(1, namespace=namespace)
394 end if
395 end do
396#endif
397
398 pop_sub(output_init)
399 end subroutine output_init
400
401 ! ---------------------------------------------------------
402 subroutine output_all(outp, namespace, space, dir, gr, ions, iter, st, hm, ks)
403 type(output_t), intent(in) :: outp
404 type(namespace_t), intent(in) :: namespace
405 class(space_t), intent(in) :: space
406 character(len=*), intent(in) :: dir
407 type(grid_t), intent(in) :: gr
408 type(ions_t), intent(in) :: ions
409 integer, intent(in) :: iter
410 type(states_elec_t), intent(inout) :: st
411 type(hamiltonian_elec_t), intent(inout) :: hm
412 type(v_ks_t), intent(inout) :: ks
413
414 integer :: idir, ierr, iout, iunit
415 character(len=MAX_PATH_LEN) :: fname
416
417 push_sub(output_all)
418 call profiling_in("OUTPUT_ALL")
419
420 if (any(outp%what)) then
421 message(1) = "Info: Writing output to " // trim(dir)
422 call messages_info(1, namespace=namespace)
423 call io_mkdir(dir, namespace)
424 end if
425
426 if (outp%what_now(option__output__mesh_r, iter)) then
427 do idir = 1, space%dim
428 write(fname, '(a,a)') 'mesh_r-', index2axis(idir)
429 call dio_function_output(outp%how(option__output__mesh_r), dir, fname, namespace, space, &
430 gr, gr%x_t(:,idir), units_out%length, ierr, pos=ions%pos, atoms=ions%atom)
431 end do
432 end if
433
434 call output_states(outp, namespace, space, dir, st, gr, ions, hm, iter)
435 call output_hamiltonian(outp, namespace, space, dir, hm, st, gr%der, ions, gr, iter, st%st_kpt_mpi_grp)
436
437 ! We can only test against the theory level here, as it is not set in the call of output_init().
438 if (outp%what_now(option__output__el_pressure, iter)) then
439 if(ks%theory_level /= kohn_sham_dft) then
440 call messages_not_implemented("el_pressure for TheoryLevel different from kohn_sham", namespace=namespace)
441 end if
442 if (st%d%spin_channels > 1) then
443 call messages_not_implemented("el_pressure for spin-polarized or spinors", namespace=namespace)
444 end if
445 end if
446
447 !hm not initialized yet when calling output_init()
448 if (outp%what(option__output__kanamoriu) .and. hm%lda_u_level /= dft_u_acbn0) then
449 message(1) = "kanamoriU output can only be computed for DFTULevel = dft_u_acbn0"
450 call messages_fatal(1, namespace=namespace)
451 end if
452
453
454 call output_localization_funct(outp, namespace, space, dir, st, hm, gr, ions, iter)
455
456 if (outp%what_now(option__output__j_flow, iter)) then
457 call output_current_flow(outp, namespace, space, dir, gr, st, hm%kpoints)
458 end if
459
460 if (outp%what_now(option__output__geometry, iter)) then
461 if (bitand(outp%how(option__output__geometry), option__outputformat__xcrysden) /= 0) then
462 call write_xsf_geometry_file(dir, "geometry", space, ions%latt, ions%pos, ions%atom, gr, namespace)
463 end if
464 if (bitand(outp%how(option__output__geometry), option__outputformat__xyz) /= 0) then
465 call ions%write_xyz(trim(dir)//'/geometry')
466 if (ions%space%is_periodic()) then
467 call ions%write_crystal(dir)
468 end if
469 end if
470 if (bitand(outp%how(option__output__geometry), option__outputformat__vtk) /= 0) then
471 call ions%write_vtk_geometry(trim(dir)//'/geometry')
472 end if
473 if (bitand(outp%how(option__output__geometry), option__outputformat__poscar) /= 0) then
474 call ions%write_poscar(trim(dir)//'/POSCAR')
475 end if
476 end if
477
478 if (outp%what_now(option__output__forces, iter)) then
479 if (bitand(outp%how(option__output__forces), option__outputformat__bild) /= 0) then
480 call ions%write_bild_forces_file(dir, "forces")
481 else
482 call write_xsf_geometry_file(dir, "forces", space, ions%latt, ions%pos, ions%atom, &
483 gr, namespace, total_forces = ions%tot_force)
484 end if
485 end if
486
487 if (outp%what_now(option__output__matrix_elements, iter)) then
488 call output_me(outp%me, namespace, space, dir, st, gr, ions, hm)
489 end if
490
491 do iout = lbound(outp%how, 1), ubound(outp%how, 1)
492 if (bitand(outp%how(iout), option__outputformat__etsf) /= 0) then
493 call output_etsf(outp, namespace, space, dir, st, gr, hm%kpoints, ions, iter)
494 exit
495 end if
496 end do
498 if (outp%what_now(option__output__berkeleygw, iter)) then
499 call output_berkeleygw(outp%bgw, namespace, space, dir, st, gr, ks, hm, ions)
500 end if
501
502 if (outp%what_now(option__output__energy_density, iter)) then
503 call output_energy_density(outp, namespace, space, dir, hm, ks, st, ions, gr)
504 end if
505
506 if (outp%what_now(option__output__stress, iter)) then
507 call io_mkdir(dir, namespace)
508 iunit = io_open(trim(dir)//'/stress', namespace, action='write')
509 call output_stress(iunit, space%periodic_dim, st%stress_tensors)
510 call io_close(iunit)
511 end if
512
513 if (hm%lda_u_level /= dft_u_none) then
514 if (outp%what_now(option__output__occ_matrices, iter))&
515 call lda_u_write_occupation_matrices(dir, hm%lda_u, st, namespace)
516
517 if (outp%what_now(option__output__effectiveu, iter))&
518 call lda_u_write_effectiveu(dir, hm%lda_u, st, namespace)
519
520 if (outp%what_now(option__output__magnetization, iter))&
521 call lda_u_write_magnetization(dir, hm%lda_u, ions, gr, st, namespace)
522
523 if (outp%what_now(option__output__local_orbitals, iter))&
524 call output_dftu_orbitals(outp, dir, namespace, space, hm%lda_u, st, gr, ions, hm%phase%is_allocated())
525
526 if (outp%what_now(option__output__kanamoriu, iter))&
527 call lda_u_write_kanamoriu(dir, st, hm%lda_u, namespace)
528
529 if (ks%oep_photon%level == oep_level_full) then
530 if (outp%what_now(option__output__photon_correlator, iter)) then
531 write(fname, '(a)') 'photon_correlator'
532 call dio_function_output(outp%how(option__output__photon_correlator), dir, trim(fname), namespace, space, &
533 gr, ks%oep_photon%pt%correlator(:,1), units_out%length, ierr, pos=ions%pos, atoms=ions%atom)
534 end if
535 end if
536 end if
537
538 ! We can only test against the theory level here, as it is not set in the call of output_init().
539 if (outp%what_now(option__output__xc_torque, iter)) then
540 if (ks%theory_level /= kohn_sham_dft .and. ks%theory_level /= generalized_kohn_sham_dft) then
541 write(message(1), '(a)') 'The output xc_torque can only be computed when there is a xc potential.'
542 call messages_fatal(1, namespace=namespace)
543 end if
544 end if
545
546 call output_xc_torque(outp, namespace, dir, gr, hm, st, ions, ions%space)
547
548 call profiling_out("OUTPUT_ALL")
549 pop_sub(output_all)
550 end subroutine output_all
551
552
553 ! ---------------------------------------------------------
554 subroutine output_localization_funct(outp, namespace, space, dir, st, hm, gr, ions, iter)
555 type(output_t), intent(in) :: outp
556 type(namespace_t), intent(in) :: namespace
557 class(space_t), intent(in) :: space
558 character(len=*), intent(in) :: dir
559 type(states_elec_t), intent(inout) :: st
560 type(hamiltonian_elec_t), intent(in) :: hm
561 type(grid_t), intent(in) :: gr
562 type(ions_t), intent(in) :: ions
563 integer, intent(in) :: iter
564
565 real(real64), allocatable :: f_loc(:,:)
566 character(len=MAX_PATH_LEN) :: fname
567 integer :: is, ierr, imax
568 type(mpi_grp_t) :: mpi_grp
569
571
572 mpi_grp = st%dom_st_kpt_mpi_grp
573
574 ! if SPIN_POLARIZED, the ELF contains one extra channel: the total ELF
575 imax = st%d%nspin
576 if (st%d%ispin == spin_polarized) imax = 3
577
578 safe_allocate(f_loc(1:gr%np, 1:imax))
579
580 ! First the ELF in real space
581 if (outp%what_now(option__output__elf, iter) .or. outp%what_now(option__output__elf_basins, iter)) then
582 assert(space%dim /= 1)
583
584 call elf_calc(space, st, gr, hm%kpoints, f_loc)
585
586 ! output ELF in real space
587 if (outp%what_now(option__output__elf, iter)) then
588 write(fname, '(a)') 'elf_rs'
589 call dio_function_output(outp%how(option__output__elf), dir, trim(fname), namespace, space, gr, &
590 f_loc(:,imax), unit_one, ierr, pos=ions%pos, atoms=ions%atom, grp = mpi_grp)
591 ! this quantity is dimensionless
592
593 if (st%d%ispin /= unpolarized) then
594 do is = 1, 2
595 write(fname, '(a,i1)') 'elf_rs-sp', is
596 call dio_function_output(outp%how(option__output__elf), dir, trim(fname), namespace, space, gr, &
597 f_loc(:, is), unit_one, ierr, pos=ions%pos, atoms=ions%atom, grp = mpi_grp)
598 ! this quantity is dimensionless
599 end do
600 end if
601 end if
602
603 if (outp%what_now(option__output__elf_basins, iter)) then
604 call out_basins(f_loc(:,1), "elf_rs_basins", outp%how(option__output__elf_basins))
605 end if
606 end if
607
608 ! Now Bader analysis
609 if (outp%what_now(option__output__bader, iter)) then
610 do is = 1, st%d%nspin
611 call dderivatives_lapl(gr%der, st%rho(:,is), f_loc(:,is))
612
613 fname = get_filename_with_spin('bader', st%d%nspin, is)
614
615 call dio_function_output(outp%how(option__output__bader), dir, trim(fname), namespace, space, gr, &
616 f_loc(:,is), units_out%length**(-2 - space%dim), ierr, &
617 pos=ions%pos, atoms=ions%atom, grp = mpi_grp)
618
619 fname = get_filename_with_spin('bader_basins', st%d%nspin, is)
620 call out_basins(f_loc(:,1), fname, outp%how(option__output__bader))
621 end do
622 end if
623
624 ! Now the pressure
625 if (outp%what_now(option__output__el_pressure, iter)) then
626 call calc_electronic_pressure(st, hm, gr, f_loc(:,1))
627 call dio_function_output(outp%how(option__output__el_pressure), dir, "el_pressure", namespace, space, gr, &
628 f_loc(:,1), unit_one, ierr, pos=ions%pos, atoms=ions%atom, grp = mpi_grp)
629 ! this quantity is dimensionless
630 end if
631
632 safe_deallocate_a(f_loc)
633
635
636 contains
637 ! ---------------------------------------------------------
638 subroutine out_basins(ff, filename, output_how)
639 real(real64), intent(in) :: ff(:)
640 character(len=*), intent(in) :: filename
641 integer(int64), intent(in) :: output_how
642
643 character(len=MAX_PATH_LEN) :: fname
644 type(basins_t) :: basins
645 integer :: iunit
646
648
649 call basins_init(basins, namespace, gr)
650 call basins_analyze(basins, namespace, gr, ff(:), st%rho, 0.01_real64)
651
652 call dio_function_output(output_how, dir, trim(filename), namespace, space, gr, &
653 real(basins%map, real64) , unit_one, ierr, pos=ions%pos, atoms=ions%atom, grp = mpi_grp)
654 ! this quantity is dimensionless
655
656 write(fname,'(4a)') trim(dir), '/', trim(filename), '.info'
657 iunit = io_open(trim(fname), namespace, action = 'write')
658 call basins_write(basins, gr, iunit)
659 call io_close(iunit)
660
661 call basins_end(basins)
662
664 end subroutine out_basins
665
666 end subroutine output_localization_funct
667
668
669 ! ---------------------------------------------------------
670 subroutine calc_electronic_pressure(st, hm, gr, pressure)
671 type(states_elec_t), intent(inout) :: st
672 type(hamiltonian_elec_t), intent(in) :: hm
673 type(grid_t), intent(in) :: gr
674 real(real64), intent(out) :: pressure(:)
675
676 real(real64), allocatable :: rho(:,:), lrho(:), tau(:,:)
677 real(real64) :: p_tf, dens
678 integer :: is, ii
679
681
682 safe_allocate( rho(1:gr%np_part, 1:st%d%nspin))
683 safe_allocate(lrho(1:gr%np))
684 safe_allocate( tau(1:gr%np, 1:st%d%nspin))
685
686 rho = m_zero
687 call density_calc(st, gr, rho)
688 call states_elec_calc_quantities(gr, st, hm%kpoints, .false., kinetic_energy_density = tau)
689
690 pressure = m_zero
691 do is = 1, st%d%spin_channels
692 lrho = m_zero
693 call dderivatives_lapl(gr%der, rho(:, is), lrho)
694
695 pressure(:) = pressure(:) + &
696 tau(:, is)/m_three - lrho(:)/m_four
697 end do
698
699 do ii = 1, gr%np
700 dens = sum(rho(ii,1:st%d%spin_channels))
701
702 p_tf = m_two/m_five*(m_three*m_pi**2)**(m_two/m_three)* &
703 dens**(m_five/m_three)
704
705 ! add XC pressure
706 ! FIXME: Not correct for spinors and spin-polarized here
707 pressure(ii) = pressure(ii) + (dens*hm%ks_pot%vxc(ii,1) - hm%energy%exchange - hm%energy%correlation)
708
709 pressure(ii) = pressure(ii)/p_tf
710 pressure(ii) = m_half*(m_one + pressure(ii)/sqrt(m_one + pressure(ii)**2))
711 end do
712
714 end subroutine calc_electronic_pressure
715
716
717 ! ---------------------------------------------------------
718 subroutine output_energy_density(outp, namespace, space, dir, hm, ks, st, ions, gr)
719 type(output_t), intent(in) :: outp
720 type(namespace_t), intent(in) :: namespace
721 class(space_t), intent(in) :: space
722 character(len=*), intent(in) :: dir
723 type(hamiltonian_elec_t), intent(in) :: hm
724 type(v_ks_t), intent(inout) :: ks
725 type(states_elec_t), intent(in) :: st
726 type(ions_t), intent(in) :: ions
727 type(grid_t), intent(in) :: gr
728
729 integer :: is, ierr, ip
730 character(len=MAX_PATH_LEN) :: fname
731 type(unit_t) :: fn_unit
732 real(real64), allocatable :: energy_density(:, :)
733 real(real64), allocatable :: ex_density(:)
734 real(real64), allocatable :: ec_density(:)
735
736 push_sub(output_energy_density)
737
738 fn_unit = units_out%energy*units_out%length**(-space%dim)
739 safe_allocate(energy_density(1:gr%np, 1:st%d%nspin))
740
741 ! the kinetic energy density
742 call states_elec_calc_quantities(gr, st, hm%kpoints, .true., kinetic_energy_density = energy_density)
743
744 ! the external potential energy density
745 do is = 1, st%d%nspin
746 do ip = 1, gr%np
747 energy_density(ip, is) = energy_density(ip, is) + st%rho(ip, is)*hm%ep%vpsl(ip)
748 end do
749 end do
750
751 ! the hartree energy density
752 do is = 1, st%d%nspin
753 do ip = 1, gr%np
754 energy_density(ip, is) = energy_density(ip, is) + m_half*st%rho(ip, is)*hm%ks_pot%vhartree(ip)
755 end do
756 end do
757
758 ! the XC energy density
759 safe_allocate(ex_density(1:gr%np))
760 safe_allocate(ec_density(1:gr%np))
761
762 call xc_get_vxc(gr, ks%xc, st, hm%kpoints, hm%psolver, namespace, space, st%rho, st%d%ispin, &
763 hm%ions%latt%rcell_volume, ex_density = ex_density, ec_density = ec_density)
764 do is = 1, st%d%nspin
765 do ip = 1, gr%np
766 energy_density(ip, is) = energy_density(ip, is) + ex_density(ip) + ec_density(ip)
767 end do
768 end do
769
770 safe_deallocate_a(ex_density)
771 safe_deallocate_a(ec_density)
772
773 do is = 1, st%d%spin_channels
774 fname = get_filename_with_spin('energy_density', st%d%nspin, is)
775 call dio_function_output(outp%how(option__output__energy_density), dir, trim(fname), namespace, space, gr, &
776 energy_density(:, is), unit_one, ierr, pos=ions%pos, atoms=ions%atom, grp = st%dom_st_kpt_mpi_grp)
777 end do
778 safe_deallocate_a(energy_density)
779
780 pop_sub(output_energy_density)
781 end subroutine output_energy_density
782
783 !--------------------------------------------------------------
784
785 logical function output_need_exchange(outp) result(need_exx)
786 type(output_t), intent(in) :: outp
787
788 need_exx =(outp%what(option__output__berkeleygw) &
789 .or. outp%me%what(option__outputmatrixelements__two_body) &
790 .or. outp%me%what(option__outputmatrixelements__two_body_exc_k))
791 end function output_need_exchange
792
793
794 ! ---------------------------------------------------------
795 subroutine output_dftu_orbitals(outp, dir, namespace, space, this, st, mesh, ions, has_phase)
796 type(output_t), intent(in) :: outp
797 character(len=*), intent(in) :: dir
798 type(namespace_t), intent(in) :: namespace
799 class(space_t), intent(in) :: space
800 type(lda_u_t), intent(in) :: this
801 type(states_elec_t), intent(in) :: st
802 class(mesh_t), intent(in) :: mesh
803 type(ions_t), intent(in) :: ions
804 logical, intent(in) :: has_phase
805
806 integer :: ios, im, ik, idim, ierr
807 complex(real64), allocatable :: tmp(:)
808 real(real64), allocatable :: dtmp(:)
809 type(orbitalset_t), pointer :: os
810 type(unit_t) :: fn_unit
811 character(len=MAX_PATH_LEN) :: fname
812
814
815 fn_unit = sqrt(units_out%length**(-space%dim))
816
817 if (this%basis%use_submesh) then
818 if (states_are_real(st)) then
819 safe_allocate(dtmp(1:mesh%np))
820 else
821 safe_allocate(tmp(1:mesh%np))
822 end if
823 end if
824
825 do ios = 1, this%norbsets
826 os => this%orbsets(ios)
827 do ik = st%d%kpt%start, st%d%kpt%end
828 do im = 1, this%orbsets(ios)%norbs
829 do idim = 1, min(os%ndim, st%d%dim)
830 if (st%nik > 1) then
831 if (min(os%ndim, st%d%dim) > 1) then
832 write(fname, '(a,i1,a,i3.3,a,i8.8,a,i1)') 'orb', im, '-os', ios, '-k', ik, '-sp', idim
833 else
834 write(fname, '(a,i1,a,i3.3,a,i8.8)') 'orb', im, '-os', ios, '-k', ik
835 end if
836 else
837 if (min(os%ndim, st%d%dim) > 1) then
838 write(fname, '(a,i1,a,i3.3,a,i1)') 'orb', im, '-os', ios, '-sp', idim
839 else
840 write(fname, '(a,i1,a,i3.3)') 'orb', im, '-os', ios
841 end if
842 end if
843 if (has_phase) then
844 if (.not. this%basis%use_submesh) then
845 call zio_function_output(outp%how(option__output__local_orbitals), dir, fname, namespace, space, &
846 mesh, os%eorb_mesh(1:mesh%np,im,idim,ik), fn_unit, ierr, pos=ions%pos, atoms=ions%atom)
847 else
848 tmp = m_z0
849 call submesh_add_to_mesh(os%sphere, os%eorb_submesh(1:os%sphere%np,idim,im,ik), tmp)
850 call zio_function_output(outp%how(option__output__local_orbitals), dir, fname, namespace, space, &
851 mesh, tmp, fn_unit, ierr, pos=ions%pos, atoms=ions%atom)
852 end if
853 else
854 if (.not. this%basis%use_submesh) then
855 if (states_are_real(st)) then
856 call dio_function_output(outp%how(option__output__local_orbitals), dir, fname, namespace, space, mesh, &
857 os%dorb(1:mesh%np,idim,im), fn_unit, ierr, pos=ions%pos, atoms=ions%atom)
858 else
859 call zio_function_output(outp%how(option__output__local_orbitals), dir, fname, namespace, space, mesh, &
860 os%zorb(1:mesh%np,idim,im), fn_unit, ierr, pos=ions%pos, atoms=ions%atom)
861 end if
862 else
863 if (states_are_real(st)) then
864 dtmp = m_z0
865 call submesh_add_to_mesh(os%sphere, os%dorb(1:os%sphere%np,idim,im), dtmp)
866 call dio_function_output(outp%how(option__output__local_orbitals), dir, fname, namespace, space, &
867 mesh, dtmp, fn_unit, ierr, pos=ions%pos, atoms=ions%atom)
868 else
869 tmp = m_z0
870 call submesh_add_to_mesh(os%sphere, os%zorb(1:os%sphere%np,idim,im), tmp)
871 call zio_function_output(outp%how(option__output__local_orbitals), dir, fname, namespace, space, &
872 mesh, tmp, fn_unit, ierr, pos=ions%pos, atoms=ions%atom)
873 end if
874 end if
875 end if
876 end do
877 end do
878 end do
879 end do
881 safe_deallocate_a(tmp)
882 safe_deallocate_a(dtmp)
883
884 pop_sub(output_dftu_orbitals)
885 end subroutine output_dftu_orbitals
886
887 ! ---------------------------------------------------------
888 logical function output_needs_current(outp, states_are_real)
889 type(output_t), intent(in) :: outp
890 logical, intent(in) :: states_are_real
891
892 output_needs_current = .false.
893
894 if (outp%what(option__output__current) &
895 .or. outp%what(option__output__current_dia) &
896 .or. outp%what(option__output__heat_current) &
897 .or. outp%what(option__output__current_kpt)) then
898 if (.not. states_are_real) then
900 else
901 message(1) = 'No current density output for real states since it is identically zero.'
902 call messages_warning(1)
903 end if
904 end if
905
906
907 end function
908
909#include "output_etsf_inc.F90"
910
911#include "output_states_inc.F90"
912
913#include "output_h_inc.F90"
914
915#include "undef.F90"
916#include "complex.F90"
917#include "output_linear_response_inc.F90"
918
919#include "undef.F90"
920#include "real.F90"
921#include "output_linear_response_inc.F90"
922
923end module output_oct_m
924
925!! Local Variables:
926!! mode: f90
927!! coding: utf-8
928!! End:
pure logical function, public accel_is_enabled()
Definition: accel.F90:402
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:616
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:193
real(real64), parameter, public m_zero
Definition: global.F90:191
real(real64), parameter, public m_four
Definition: global.F90:195
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:189
complex(real64), parameter, public m_z0
Definition: global.F90:201
integer, parameter, public generalized_kohn_sham_dft
Definition: global.F90:241
integer, parameter, public kohn_sham_dft
Definition: global.F90:241
real(real64), parameter, public m_epsilon
Definition: global.F90:207
real(real64), parameter, public m_half
Definition: global.F90:197
real(real64), parameter, public m_one
Definition: global.F90:192
real(real64), parameter, public m_three
Definition: global.F90:194
real(real64), parameter, public m_five
Definition: global.F90:196
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:355
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:486
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
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:1067
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:524
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:409
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1039
subroutine, public messages_info(no_lines, iunit, debug_only, stress, all_nodes, namespace)
Definition: messages.F90:593
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:766
subroutine, public output_states(outp, namespace, space, dir, st, gr, ions, hm, iter)
Definition: output.F90:1094
logical function, public output_needs_current(outp, states_are_real)
Definition: output.F90:984
subroutine, public output_hamiltonian(outp, namespace, space, dir, hm, st, der, ions, gr, iter, grp)
Definition: output.F90:1390
subroutine, public output_all(outp, namespace, space, dir, gr, ions, iter, st, hm, ks)
Definition: output.F90:498
logical function, public output_need_exchange(outp)
Definition: output.F90:881
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:1709
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:1816
subroutine, public doutput_lr(outp, namespace, space, dir, st, mesh, lr, idir, isigma, ions, pert_unit)
Definition: output.F90:2085
subroutine, public output_scalar_pot(outp, namespace, space, dir, mesh, ions, ext_partners, time)
Definition: output.F90:1672
subroutine output_energy_density(outp, namespace, space, dir, hm, ks, st, ions, gr)
Definition: output.F90:814
subroutine output_dftu_orbitals(outp, dir, namespace, space, this, st, mesh, ions, has_phase)
Definition: output.F90:891
subroutine output_localization_funct(outp, namespace, space, dir, st, hm, gr, ions, iter)
Definition: output.F90:650
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:1035
logical function, public parse_is_defined(namespace, name)
Definition: parser.F90:455
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:615
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:163
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:1499
Definition: xc.F90:116
integer, parameter, public oep_level_full
Definition: xc_oep.F90:175
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_vxc.F90:185
subroutine out_basins(ff, filename, output_how)
Definition: output.F90:734
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)