Octopus
local_write.F90
Go to the documentation of this file.
1!! Copyright (C) 2014 M. Marques, A. Castro, A. Rubio, G. Bertsch, J.Jornet-Somoza
2!! Copyright (C) 2020 M. Oliveira
3!!
4!! This program is free software; you can redistribute it and/or modify
5!! it under the terms of the GNU General Public License as published by
6!! the Free Software Foundation; either version 2, or (at your option)
7!! any later version.
8!!
9!! This program is distributed in the hope that it will be useful,
10!! but WITHOUT ANY WARRANTY; without even the implied warranty of
11!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12!! GNU General Public License for more details.
13!!
14!! You should have received a copy of the GNU General Public License
15!! along with this program; if not, write to the Free Software
16!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17!! 02110-1301, USA.
18!!
19
20#include "global.h"
21
23 use iso_c_binding
24 use debug_oct_m
26 use global_oct_m
29 use io_oct_m
31 use ions_oct_m
32 use kick_oct_m
34 use mesh_oct_m
36 use mpi_oct_m
38 use parser_oct_m
41 use space_oct_m
43 use unit_oct_m
46 use v_ks_oct_m
48
49 implicit none
50
51 private
52 public :: &
58
59 integer, parameter :: &
60 LOCAL_OUT_MULTIPOLES = 1, &
63 local_out_energy = 4, &
65
67 private
68 type(c_ptr) :: handle
69 logical :: write = .false.
70 end type local_write_prop_t
71
72 type local_write_t
73 private
74 type(local_write_prop_t), allocatable :: out(:)
75 integer(int64) :: how
76 integer :: lmax
77 end type local_write_t
78
79contains
80
81 ! ---------------------------------------------------------
82 logical function local_write_check_hm(writ) result(check)
83 type(local_write_t), intent(in) :: writ
84
85 push_sub(local_write_check_hm)
86
87 check = writ%out(local_out_potential)%write .or. writ%out(local_out_energy)%write
88
90 end function local_write_check_hm
91
92 ! ---------------------------------------------------------
93 subroutine local_write_init(writ, namespace, lab, iter, dt)
94 type(local_write_t), intent(inout) :: writ
95 type(namespace_t), intent(in) :: namespace
96 character(len=15), intent(in) :: lab
97 integer, intent(in) :: iter
98 real(real64), intent(in) :: dt
99
100 integer :: first, flags, iout, default
101
102 push_sub(local_write_init)
103
104 ! FIXME: if and when these routines are called from a normal run, the Section can be Output.
105 ! but then it will need to be in a different folder, since src/util is not linked by the other folders.
106
107 !%Variable LDOutput
108 !%Type flag
109 !%Default multipoles
110 !%Section Utilities::oct-local_multipoles
111 !%Description
112 !% Defines what should be output during the local domains
113 !% simulation. Many of the options can increase the computational
114 !% cost of the simulation, so only use the ones that you need. In
115 !% most cases the default value is enough, as it is adapted to the
116 !% details.
117 !%Option multipoles 1
118 !% Outputs the (electric) multipole moments of the density to the file <tt>ld.general/multipoles</tt>.
119 !% This is required to, <i>e.g.</i>, calculate optical absorption spectra of finite systems. The
120 !% maximum value of <math>l</math> can be set with the variable <tt>LDMultipoleLmax</tt>.
121 !%Option density 2
122 !% If set, <tt>octopus</tt> outputs the densities corresponding to the local domains to
123 !% the folder <tt>ld.general/densities</tt>.
124 !% The output format is set by the <tt>LDOutputFormat</tt> input variable.
125 !%Option local_v 4
126 !% If set, <tt>octopus</tt> outputs the different components of the potential
127 !% to the folder <tt>ld.general/potential</tt>.
128 !% The output format is set by the <tt>LDOutputFormat</tt> input variable.
129 !%Option energy 8
130 !% If set, <tt>octopus</tt> outputs the different components of the energy of the local domains
131 !% to the folder <tt>ld.general/energy</tt>.
132 !%End
133
134 default = 2**(local_out_multipoles - 1)
135
136 call parse_variable(namespace, 'LDOutput', default, flags)
137
138 if (.not. varinfo_valid_option('LDOutput', flags, is_flag = .true.)) then
139 call messages_input_error(namespace, 'LDOutput')
140 end if
141
142 safe_allocate(writ%out(1:local_out_max))
143 do iout = 1, local_out_max
144 writ%out(iout)%write = (bitand(flags, 2**(iout - 1)) /= 0)
145 end do
146
147 call messages_obsolete_variable(namespace, 'LDOutputHow', 'LDOutputFormat')
148
149 !%Variable LDOutputFormat
150 !%Type flag
151 !%Default none
152 !%Section Utilities::oct-local_multipoles
153 !%Description
154 !% Describes the format of the output files (see <tt>LDOutput</tt>).
155 !% It can take the same values as <tt>OutputFormat</tt> flag.
156 !%End
157 call parse_variable(namespace, 'LDOutputFormat', 0, writ%how)
158 if (.not. varinfo_valid_option('OutputFormat', writ%how, is_flag=.true.)) then
159 call messages_input_error(namespace, 'LDOutputFormat')
160 end if
162 !%Variable LDMultipoleLmax
163 !%Type integer
164 !%Default 1
165 !%Section Utilities::oct-local_multipoles
166 !%Description
167 !% Maximum electric multipole of the density output to the file <tt>local.multipoles/<>domain%<>.multipoles</tt>
168 !% during a time-dependent simulation. Must be non-negative.
169 !%End
170 call parse_variable(namespace, 'LDMultipoleLmax', 1, writ%lmax)
171 if (writ%lmax < 0) then
172 write(message(1), '(a,i6,a)') "Input: '", writ%lmax, "' is not a valid LDMultipoleLmax."
173 message(2) = '(Must be LDMultipoleLmax >= 0)'
174 call messages_fatal(2, namespace=namespace)
175 end if
176
177 call io_mkdir('local.general', namespace)
178
179 if (mpi_world%is_root()) then
180 if (writ%out(local_out_multipoles)%write) then
181 call io_mkdir('local.general/multipoles', namespace)
182 call write_iter_init(writ%out(local_out_multipoles)%handle, iter, units_from_atomic(units_out%time, dt), &
183 trim(io_workpath("local.general/multipoles/"//trim(lab)//".multipoles", namespace)))
184 end if
185
186 if (writ%out(local_out_potential)%write) then
187 call io_mkdir('local.general/potential', namespace)
188 call write_iter_init(writ%out(local_out_potential)%handle, first, units_from_atomic(units_out%time, dt), &
189 trim(io_workpath("local.general/potential/"//trim(lab)//".potential", namespace)))
190 end if
191
192 if (writ%out(local_out_density)%write) then
193 call io_mkdir('local.general/densities', namespace)
194 call write_iter_init(writ%out(local_out_density)%handle, first, units_from_atomic(units_out%time, dt), &
195 trim(io_workpath("local.general/densities/"//trim(lab)//".densities", namespace)))
196 end if
197
198 if (writ%out(local_out_energy)%write) then
199 call io_mkdir('local.general/energy', namespace)
200 call write_iter_init(writ%out(local_out_energy)%handle, iter, units_from_atomic(units_out%time, dt), &
201 trim(io_workpath("local.general/energy/"//trim(lab)//".energy", namespace)))
202 end if
203 end if
204
205 pop_sub(local_write_init)
206 end subroutine local_write_init
207
208 ! ---------------------------------------------------------
209 subroutine local_write_end(writ)
210 type(local_write_t), intent(inout) :: writ
211
212 integer :: iout
213
214 push_sub(local_write_end)
215
216 if (mpi_world%is_root()) then
217 do iout = 1, local_out_max
218 if (writ%out(iout)%write) then
219 call write_iter_end(writ%out(iout)%handle)
220 end if
221 end do
222 end if
223 safe_deallocate_a(writ%out)
224
225 pop_sub(local_write_end)
226 end subroutine local_write_end
227
228 ! ---------------------------------------------------------
229 subroutine local_write_iter(writ, namespace, space, lab, ions_mask, mesh_mask, mesh, st, hm, ks, &
230 ions, ext_partners, kick, iter, l_start, ldoverwrite)
231 type(local_write_t), intent(inout) :: writ
232 type(namespace_t), intent(in) :: namespace
233 type(electron_space_t), intent(in) :: space
234 character(len=15), intent(in) :: lab
235 logical, intent(in) :: ions_mask(:)
236 logical, intent(in) :: mesh_mask(:)
237 class(mesh_t), intent(in) :: mesh
238 type(states_elec_t), intent(inout) :: st
239 type(hamiltonian_elec_t), intent(inout) :: hm
240 type(v_ks_t), intent(inout) :: ks
241 type(ions_t), intent(inout) :: ions
242 type(partner_list_t), intent(in) :: ext_partners
243 type(kick_t), intent(inout) :: kick
244 integer, intent(in) :: iter
245 integer, intent(in) :: l_start
246 logical, intent(in) :: ldoverwrite
247
248
249 push_sub(local_write_iter)
250 call profiling_in("LOCAL_WRITE_ITER")
251
252 if (writ%out(local_out_multipoles)%write) then
253 call local_write_multipole(writ%out(local_out_multipoles), namespace, lab, ions_mask, mesh_mask, mesh, ions, st, &
254 writ%lmax, kick, iter, l_start, ldoverwrite, writ%how)
255 if (mpi_world%is_root()) then
256 call write_iter_flush(writ%out(local_out_multipoles)%handle)
257 end if
258 end if
259
260 if (writ%out(local_out_density)%write .or. writ%out(local_out_potential)%write) then
261 call local_write_density(writ%out(local_out_density), namespace, space, writ%out(local_out_potential), lab, mesh_mask, &
262 mesh, ions, ext_partners, st, hm, ks, iter, writ%how)
263 end if
264
265 if (writ%out(local_out_energy)%write) then
266 call local_write_energy(writ%out(local_out_energy), namespace, space, lab, mesh_mask, mesh, ions, &
267 ext_partners, st, hm, ks, iter, l_start, ldoverwrite)
268 if (mpi_world%is_root()) then
269 call write_iter_flush(writ%out(local_out_energy)%handle)
270 end if
271 end if
272
273 call profiling_out("LOCAL_WRITE_ITER")
274 pop_sub(local_write_iter)
275 end subroutine local_write_iter
276
277 ! ---------------------------------------------------------
278 subroutine local_write_density(out_dens, namespace, space, out_pot, lab, mesh_mask, mesh, ions, &
279 ext_partners, st, hm, ks, iter, how)
280 type(local_write_prop_t), intent(inout) :: out_dens
281 type(namespace_t), intent(in) :: namespace
282 type(electron_space_t), intent(in) :: space
283 type(local_write_prop_t), intent(inout) :: out_pot
284 character(len=15), intent(in) :: lab
285 logical, intent(in) :: mesh_mask(:)
286 class(mesh_t), intent(in) :: mesh
287 type(ions_t), intent(inout) :: ions
288 type(partner_list_t), intent(in) :: ext_partners
289 type(states_elec_t), intent(inout) :: st
290 type(hamiltonian_elec_t), intent(inout) :: hm
291 type(v_ks_t), intent(inout) :: ks
292 integer, intent(in) :: iter
293 integer(int64), intent(in) :: how
294
295 integer :: is, ierr
296 character(len=120) :: folder, out_name
297 real(real64), allocatable :: tmp_rho(:), st_rho(:)
298 real(real64), allocatable :: tmp_vh(:)
299
300 push_sub(local_write_density)
301
302 safe_allocate(tmp_rho(1:mesh%np))
303 safe_allocate(st_rho(1:mesh%np_part))
304 safe_allocate(tmp_vh(1:mesh%np))
305
306 ! FIXME: use just v_ks_calc subroutine for either compute vhartree and vxc
307 do is = 1, st%d%nspin
308 st_rho(:) = st%rho(:, is)
309 tmp_rho = m_zero
310 where (mesh_mask)
311 tmp_rho = st_rho(1:mesh%np)
312 end where
313
314 if (out_dens%write) then
315 folder = 'local.general/densities/'//trim(lab)//'.densities/'
316 write(out_name, '(a,a1,i0,a1,i7.7)') trim(lab), '.', is,'.', iter
317 call dio_function_output(how, trim(folder), trim(out_name), namespace, space, mesh, tmp_rho(1:mesh%np), &
318 units_out%length, ierr, pos=ions%pos, atoms=ions%atom)
319 end if
320
321 if (out_pot%write) then
322 !Computes Hartree potential just for n[r], r belongs to id domain.
323 tmp_vh = m_zero
324 call dpoisson_solve(hm%psolver, namespace, tmp_vh, tmp_rho)
325 folder = 'local.general/potential/'//trim(lab)//'.potential/'
326 write(out_name, '(a,i0,a1,i7.7)') 'vh.', is, '.', iter
327 call dio_function_output(how, trim(folder), trim(out_name), namespace, space, mesh, tmp_vh, units_out%length, ierr, &
328 pos=ions%pos, atoms=ions%atom)
329 !Computes XC potential
330 st%rho(:, is) = m_zero
331 where (mesh_mask)
332 st%rho(1:mesh%np, is) = st_rho(1:mesh%np)
333 end where
334 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval = .false. , calc_energy = .false.)
335 folder = 'local.general/potential/'//trim(lab)//'.potential/'
336 write(out_name, '(a,i0,a1,i7.7)') 'vxc.', is, '.', iter
337 call dio_function_output(how, trim(folder), trim(out_name), namespace, space, mesh, &
338 hm%ks_pot%vxc(:,is), units_out%length, ierr, pos=ions%pos, atoms=ions%atom)
339 st%rho(:,is) = st_rho(:)
340 end if
341 end do
342
343 if (out_pot%write) then
344 do is = 1, st%d%nspin
345 !Computes Hartree potential
346 call dpoisson_solve(hm%psolver, namespace, hm%ks_pot%vhartree, st%rho(1:mesh%np, is))
347 folder = 'local.general/potential/'
348 write(out_name, '(a,i0,a1,i7.7)') 'global-vh.', is, '.', iter
349 call dio_function_output(how, trim(folder), trim(out_name), namespace, space, mesh, hm%ks_pot%vhartree, &
350 units_out%length, ierr, pos=ions%pos, atoms=ions%atom)
351 !Computes global XC potential
352 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval = .false. , calc_energy = .false.)
353 folder = 'local.general/potential/'
354 write(out_name, '(a,i0,a1,i7.7)') 'global-vxc.', is, '.', iter
355 call dio_function_output(how, trim(folder), trim(out_name), namespace, space, mesh, hm%ks_pot%vxc(:,is), &
356 units_out%length, ierr, pos=ions%pos, atoms=ions%atom)
357 end do
358 end if
359
360 safe_deallocate_a(tmp_vh)
361 safe_deallocate_a(st_rho)
362 safe_deallocate_a(tmp_rho)
363
364 pop_sub(local_write_density)
365 end subroutine local_write_density
366
367 ! ---------------------------------------------------------
368 subroutine local_write_energy(out_energy, namespace, space, lab, mesh_mask, mesh, ions, ext_partners, &
369 st, hm, ks, iter, l_start, start)
370 type(local_write_prop_t), intent(inout) :: out_energy
371 type(namespace_t), intent(in) :: namespace
372 type(electron_space_t), intent(in) :: space
373 character(len=15), intent(in) :: lab
374 logical, intent(in) :: mesh_mask(:)
375 class(mesh_t), intent(in) :: mesh
376 type(ions_t), intent(inout) :: ions
377 type(partner_list_t), intent(in) :: ext_partners
378 type(states_elec_t), intent(inout) :: st
379 type(hamiltonian_elec_t), intent(inout) :: hm
380 type(v_ks_t), intent(inout) :: ks
381 integer, intent(in) :: iter
382 integer, intent(in) :: l_start
383 logical, intent(in) :: start
384
385 !integer :: ix
386 integer :: ii, is
387 character(len=120) :: aux
388 real(real64) :: geh, gexc, leh, lexc
389 real(real64), allocatable :: tmp_rhoi(:), st_rho(:)
390 real(real64), allocatable :: hm_vxc(:), hm_vh(:)
391
392 push_sub(local_write_energy)
393
394 safe_allocate(tmp_rhoi(1:mesh%np))
395 safe_allocate(st_rho(1:mesh%np_part))
396 safe_allocate(hm_vxc(1:mesh%np))
397 safe_allocate(hm_vh(1:mesh%np_part))
398
399 if (mpi_world%is_root()) then ! only first node outputs
400
401 if (iter == l_start .and. start) then
402 call local_write_print_header_init(out_energy%handle)
403
404 ! first line -> column names
405 write(aux,'(a,2x,a)')'# Domain ID:', trim(lab)
406 call write_iter_header(out_energy%handle, trim(aux))
407 write(aux,'(a)') trim(lab)
408 call write_iter_header(out_energy%handle, trim(aux))
409 call write_iter_nl(out_energy%handle)
410 call write_iter_header_start(out_energy%handle)
411 aux = 'Total Hartree'
412 call write_iter_header(out_energy%handle, trim(aux))
413 aux = 'Total Int[n*vxc]'
414 call write_iter_header(out_energy%handle, trim(aux))
415 write(aux,'(a)') 'Int[n*vh]'
416 call write_iter_header(out_energy%handle, trim(aux))
417 write(aux,'(a)')'Int[n*vxc]'
418 call write_iter_header(out_energy%handle, trim(aux))
419! Cross-terms between different domains are currently disabled
420! do jd = 1, nd
421! if (jd /= id) then
422! write(aux,'(a,i2.2,a,i2.2,a)')'Int[n(',id,')*vh(',jd,')]'
423! call write_iter_header(out_energy(id)%handle, trim(aux))
424! write(aux,'(a,i2.2,a,i2.2,a)')'Int[n(',id,')*vxc(',jd,')]'
425! call write_iter_header(out_energy(id)%handle, trim(aux))
426! end if
427! end do
428 call write_iter_nl(out_energy%handle)
429
430 ! units
431 call write_iter_string(out_energy%handle, '#[Iter n.]')
432 call write_iter_header(out_energy%handle, '[' // trim(units_abbrev(units_out%time)) // ']')
433 do ii = 1, 4
434 call write_iter_header(out_energy%handle, '[' // trim(units_abbrev(units_out%energy)) // ']')
435 end do
436 call write_iter_nl(out_energy%handle)
437
438 call local_write_print_header_end(out_energy%handle)
439 end if
440 end if
441
442 ! TODO: Make new files for each nspin value.
443 do is = 1, st%d%nspin
444 ! Compute Hartree potential
445 call dpoisson_solve(hm%psolver, namespace, hm%ks_pot%vhartree, st%rho(1:mesh%np, is))
446 ! Compute XC potential
447 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval = .false. , calc_energy = .false.)
448
449 st_rho(:) = st%rho(:, is)
450 hm_vxc(:) = hm%ks_pot%vxc(:, is)
451 hm_vh(:) = hm%ks_pot%vhartree(:)
452 ! Compute Global Values
453 geh = dmf_integrate(mesh, st_rho(1:mesh%np)*hm_vh(1:mesh%np))
454 gexc = dmf_integrate(mesh, st_rho(1:mesh%np)*hm_vxc(1:mesh%np))
455 if (mpi_world%is_root()) then
456 call write_iter_set(out_energy%handle, iter)
457 call write_iter_start(out_energy%handle)
458 call write_iter_double(out_energy%handle, units_from_atomic(units_out%energy, geh), 1)
459 call write_iter_double(out_energy%handle, units_from_atomic(units_out%energy, gexc), 1)
460 end if
461
462 !Compute Local Values
463 hm%ks_pot%vxc(:,is) = m_zero
464 hm%ks_pot%vhartree(:) = m_zero
465 st%rho(:, is) = m_zero
466 where (mesh_mask)
467 st%rho(1:mesh%np, is) = st_rho(1:mesh%np)
468 end where
469
470 call v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, &
471 calc_eigenval = .false. , calc_energy = .false.)
472 tmp_rhoi(1:mesh%np) = st%rho(1:mesh%np, is)
473 !eh = Int[n*v_h]
474 leh = dmf_integrate(mesh, tmp_rhoi*hm%ks_pot%vhartree(1:mesh%np))
475 !exc = Int[n*v_xc]
476 lexc = dmf_integrate(mesh, tmp_rhoi*hm%ks_pot%vxc(1:mesh%np,is))
477 if (mpi_world%is_root()) then
478 call write_iter_double(out_energy%handle, units_from_atomic(units_out%energy, leh), 1)
479 call write_iter_double(out_energy%handle, units_from_atomic(units_out%energy, lexc), 1)
480 end if
481
482! Cross-terms between different domains are currently disabled
483! do jd = 1, nd
484! if (jd /= id) then
485 ! TODO: Check for domains overlaping to avoid segmentation fault.
486 !eh = Int[n(id)*v_h(jd)]
487! st%rho(:,is) = M_ZERO
488! hm%vxc(:,is) = M_ZERO
489! hm%vhartree(:) = M_ZERO
490! do ix = 1, mesh%np
491! if (mesh_mask(ix, jd)) st%rho(ix, is) = st_rho(ix)
492! end do
493! call v_ks_calc(ks, namespace, hm, st, ions, calc_eigenval = .false. , calc_energy = .false.)
494! leh = dmf_integrate(mesh, tmp_rhoi*hm%vhartree(1:mesh%np))
495 !exc = Int[n(id)*v_xc(jd)]
496! lexc = dmf_integrate(mesh, tmp_rhoi*hm%vxc(1:mesh%np, is))
497! if (mpi_world%is_root()) then
498! call write_iter_double(out_energy(id)%handle, units_from_atomic(units_out%energy, leh), 1)
499! call write_iter_double(out_energy(id)%handle, units_from_atomic(units_out%energy, lexc), 1)
500! end if
501! end if
502! end do
503 st%rho(:,is) = st_rho(:)
504 hm%ks_pot%vxc(:,is) = hm_vxc(:)
505 hm%ks_pot%vhartree(:) = hm_vh(:)
506 if (mpi_world%is_root()) then
507 call write_iter_nl(out_energy%handle)
508 end if
509 end do
510
511 safe_deallocate_a(hm_vh)
512 safe_deallocate_a(hm_vxc)
513 safe_deallocate_a(st_rho)
514 safe_deallocate_a(tmp_rhoi)
515
516 pop_sub(local_write_energy)
517 end subroutine local_write_energy
518
519 ! ---------------------------------------------------------
520 subroutine local_write_multipole(out_multip, namespace, lab, ions_mask, mesh_mask, mesh, ions, st, lmax, kick, iter, l_start, &
521 start, how)
522 type(local_write_prop_t), intent(inout) :: out_multip
523 type(namespace_t), intent(in) :: namespace
524 character(len=15), intent(in) :: lab
525 logical, intent(in) :: ions_mask(:)
526 logical, intent(in) :: mesh_mask(:)
527 class(mesh_t), intent(in) :: mesh
528 type(ions_t), intent(in) :: ions
529 type(states_elec_t), intent(in) :: st
530 integer, intent(in) :: lmax
531 type(kick_t), intent(in) :: kick
532 integer, intent(in) :: iter
533 integer, intent(in) :: l_start
534 logical, intent(in) :: start
535 integer(int64), intent(in) :: how
536
537 integer :: is, ll, mm, add_lm
538 character(len=120) :: aux
539 real(real64) :: ionic_dipole(ions%space%dim), center(mesh%box%dim)
540 real(real64), allocatable :: multipole(:,:)
541 logical :: use_ionic_dipole
542
543 push_sub(local_write_multipole)
544
545 if (mpi_world%is_root().and. iter == l_start .and. start) then
546 call local_write_print_header_init(out_multip%handle)
547
548 write(aux, '(a15,i2)') '# nspin ', st%d%nspin
549 call write_iter_string(out_multip%handle, aux)
550 call write_iter_nl(out_multip%handle)
551
552 write(aux, '(a15,i2)') '# lmax ', lmax
553 call write_iter_string(out_multip%handle, aux)
554 call write_iter_nl(out_multip%handle)
555
556 call kick_write(kick, out = out_multip%handle)
557
558 call write_iter_header_start(out_multip%handle)
559
560 do is = 1, st%d%nspin
561 write(aux,'(a18,i1,a1)') 'Electronic charge(', is,')'
562 call write_iter_header(out_multip%handle, aux)
563 if (lmax > 0) then
564 write(aux, '(a3,a1,i1,a1)') '<x>', '(', is,')'
565 call write_iter_header(out_multip%handle, aux)
566 write(aux, '(a3,a1,i1,a1)') '<y>', '(', is,')'
567 call write_iter_header(out_multip%handle, aux)
568 write(aux, '(a3,a1,i1,a1)') '<z>', '(', is,')'
569 call write_iter_header(out_multip%handle, aux)
570 end if
571 do ll = 2, lmax
572 do mm = -ll, ll
573 write(aux, '(a2,i2,a4,i2,a2,i1,a1)') 'l=', ll, ', m=', mm, ' (', is,')'
574 call write_iter_header(out_multip%handle, aux)
575 end do
576 end do
577 end do
578 call write_iter_nl(out_multip%handle)
579
580 ! units
581 call write_iter_string(out_multip%handle, '#[Iter n.]')
582 call write_iter_header(out_multip%handle, '[' // trim(units_abbrev(units_out%time)) // ']')
583
584 do is = 1, st%d%nspin
585 do ll = 0, lmax
586 do mm = -ll, ll
587 select case (ll)
588 case (0)
589 call write_iter_header(out_multip%handle, 'Electrons')
590 case (1)
591 call write_iter_header(out_multip%handle, '[' // trim(units_abbrev(units_out%length)) // ']')
592 case default
593 write(aux, '(a,a2,i1)') trim(units_abbrev(units_out%length)), "**", ll
594 call write_iter_header(out_multip%handle, '[' // trim(aux) // ']')
595 end select
596 end do
597 end do
598 end do
599 call write_iter_nl(out_multip%handle)
600
601 call local_write_print_header_end(out_multip%handle)
602 call write_iter_flush(out_multip%handle)
603 end if
604
605 center(1:ions%space%dim) = ions%center_of_mass(mask=ions_mask)
606
607 safe_allocate(multipole(1:(lmax + 1)**2, 1:st%d%nspin))
608 multipole = m_zero
609
610 do is = 1, st%d%nspin
611 call dmf_multipoles(mesh, st%rho(:,is), lmax, multipole(:,is), mask=mesh_mask)
612 end do
613 ! Setting center of mass as reference needed for non-neutral systems.
614 do is = 1, st%d%nspin
615 multipole(2:mesh%box%dim+1, is) = multipole(2:mesh%box%dim+1, is) - center(1:mesh%box%dim)*multipole(1, is)
616 end do
617
618 ! For transition densities or density differences there is no need to add the ionic dipole
619
620 !%Variable LDIonicDipole
621 !%Type logical
622 !%Default yes
623 !%Section Utilities::oct-local_multipoles
624 !%Description
625 !% Describes if the ionic dipole has to be take into account
626 !% when computing the multipoles.
627 !%End
628 call parse_variable(namespace, 'LDIonicDipole', .true., use_ionic_dipole)
629 if (use_ionic_dipole) then
630 ! Calculate ionic dipole, setting the center of mass as reference, as that is needed for non-neutral systems.
631 ionic_dipole(1:ions%space%dim) = ions%dipole(mask=ions_mask) + &
632 p_proton_charge*ions%val_charge(mask=ions_mask)*center(1:ions%space%dim)
633
634 do is = 1, st%d%nspin
635 multipole(2:ions%space%dim+1, is) = -ionic_dipole(1:ions%space%dim)/st%d%nspin - multipole(2:ions%space%dim+1, is)
636 end do
637 end if
638
639 if (mpi_world%is_root()) then
640 call write_iter_set(out_multip%handle, iter)
641 call write_iter_start(out_multip%handle)
642 do is = 1, st%d%nspin
643 add_lm = 1
644 do ll = 0, lmax
645 do mm = -ll, ll
646 call write_iter_double(out_multip%handle, units_from_atomic(units_out%length**ll, multipole(add_lm, is)), 1)
647 add_lm = add_lm + 1
648 end do
649 end do
650 end do
651 call write_iter_nl(out_multip%handle)
652 end if
653
654 ! Write multipoles in BILD format
655 if (bitand(how, option__outputformat__bild) /= 0 .and. mpi_world%is_root()) then
656 assert(mesh%box%dim == 3)
657 !FIXME: to include spin larger than 1.
658 is = 1
659 call out_bld_multipoles(namespace, multipole(2:4, is), center, lab, iter)
660 end if
661
662 safe_deallocate_a(multipole)
663
664 pop_sub(local_write_multipole)
665 end subroutine local_write_multipole
666
667 ! ---------------------------------------------------------
668 subroutine out_bld_multipoles(namespace, multipoles, center, label, iter)
669 type(namespace_t), intent(in) :: namespace
670 real(real64), intent(in) :: multipoles(:)
671 real(real64), intent(in) :: center(:)
672 character(15), intent(in) :: label
673 integer, intent(in) :: iter
674
675 integer :: ll, out_bld
676 character(len=80) :: filename, folder
677 real(real64) :: dipolearrow(3,2)
678
679 push_sub(out_bld_multipoles)
680
681 write(folder,'(a,a)')'local.general/multipoles/',trim(label)
682 call io_mkdir(folder, namespace)
683 write(filename,'(a,a,a,a,i7.7,a)')trim(folder),'/',trim(label),'.',iter,'.bld'
684 out_bld = io_open(trim(filename), namespace, action='write')
685
686 write(out_bld,'(a,a,a,i7)')'.comment ** Arrow for the dipole moment centered at the center of mass for ', &
687 trim(label), ' domain and iteration number: ',iter
688 write(out_bld,'(a)')''
689 write(out_bld,'(a)')'.color red'
690 write(out_bld,'(a,3(f12.6,2x),a)')'.sphere ',(units_from_atomic(units_out%length,center(ll)), ll= 1, 3),' 0.2'
691 do ll = 1, 3
692 dipolearrow(ll,1) = units_from_atomic(units_out%length, center(ll) - multipoles(ll)/m_two)
693 dipolearrow(ll,2) = units_from_atomic(units_out%length, center(ll) + multipoles(ll)/m_two)
694 end do
695 write(out_bld,'(a,6(f12.6,2x),a)')'.arrow ',(dipolearrow(ll,1), ll= 1, 3), &
696 (dipolearrow(ll,2), ll= 1, 3), ' 0.1 0.5 0.90'
697 call io_close(out_bld)
698
699 pop_sub(out_bld_multipoles)
700 end subroutine out_bld_multipoles
701
702 ! ---------------------------------------------------------
703 subroutine local_write_print_header_init(out)
704 type(c_ptr), intent(inout) :: out
705
707 call write_iter_clear(out)
708 call write_iter_string(out,'################################################################################')
709 call write_iter_nl(out)
710 call write_iter_string(out,'# HEADER')
711 call write_iter_nl(out)
712
714 end subroutine local_write_print_header_init
715
716
717 ! ---------------------------------------------------------
718 subroutine local_write_print_header_end(out)
719 type(c_ptr), intent(inout) :: out
720
722
723 call write_iter_string(out,'################################################################################')
724 call write_iter_nl(out)
725
727 end subroutine local_write_print_header_end
728
729end module local_write_oct_m
730
731!! Local Variables:
732!! mode: f90
733!! coding: utf-8
734!! End:
Sets the iteration number to the C object.
Definition: write_iter.F90:170
Writes to the corresponding file and adds one to the iteration. Must be called after write_iter_init(...
Definition: write_iter.F90:163
real(real64), parameter, public m_two
Definition: global.F90:193
real(real64), parameter, public p_proton_charge
Definition: global.F90:236
real(real64), parameter, public m_zero
Definition: global.F90:191
This module defines classes and functions for interaction partners.
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.
Definition: io.F90:116
subroutine, public io_close(iunit, grp)
Definition: io.F90:467
character(len=max_path_len) function, public io_workpath(path, namespace)
construct path name from given name and namespace
Definition: io.F90:318
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
subroutine, public kick_write(kick, iunit, out)
Definition: kick.F90:892
integer, parameter local_out_potential
subroutine local_write_print_header_init(out)
subroutine, public local_write_init(writ, namespace, lab, iter, dt)
subroutine, public local_write_end(writ)
integer, parameter local_out_energy
logical function, public local_write_check_hm(writ)
subroutine local_write_print_header_end(out)
integer, parameter local_out_density
integer, parameter local_out_max
subroutine local_write_density(out_dens, namespace, space, out_pot, lab, mesh_mask, mesh, ions, ext_partners, st, hm, ks, iter, how)
subroutine out_bld_multipoles(namespace, multipoles, center, label, iter)
subroutine local_write_multipole(out_multip, namespace, lab, ions_mask, mesh_mask, mesh, ions, st, lmax, kick, iter, l_start, start, how)
subroutine local_write_energy(out_energy, namespace, space, lab, mesh_mask, mesh, ions, ext_partners, st, hm, ks, iter, l_start, start)
subroutine, public local_write_iter(writ, namespace, space, lab, ions_mask, mesh_mask, mesh, st, hm, ks, ions, ext_partners, kick, iter, l_start, ldoverwrite)
This module defines various routines, operating on mesh functions.
subroutine, public dmf_multipoles(mesh, ff, lmax, multipole, mask)
This routine calculates the multipoles of a function ff.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:999
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_input_error(namespace, var, details, row, column)
Definition: messages.F90:690
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:272
subroutine, public dpoisson_solve(this, namespace, pot, rho, all_nodes, kernel, reset)
Calculates the Poisson equation. Given the density returns the corresponding potential.
Definition: poisson.F90:872
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
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:134
character(len=20) pure function, public units_abbrev(this)
Definition: unit.F90:225
This module defines the unit system, used for input and output.
type(unit_system_t), public units_out
subroutine, public v_ks_calc(ks, namespace, space, hm, st, ions, ext_partners, calc_eigenval, time, calc_energy, calc_current, force_semilocal)
Definition: v_ks.F90:754
Explicit interfaces to C functions, defined in write_iter_low.cc.
Definition: write_iter.F90:116
subroutine, public write_iter_header(out, string)
Definition: write_iter.F90:249
subroutine, public write_iter_string(out, string)
Definition: write_iter.F90:265
subroutine, public write_iter_init(out, iter, factor, file)
Definition: write_iter.F90:229
Extension of space that contains the knowledge of the spin dimension.
Describes mesh distribution to nodes.
Definition: mesh.F90:187
The states_elec_t class contains all electronic wave functions.
int true(void)