Octopus
local_multipoles.F90
Go to the documentation of this file.
1!! Copyright (C) 2014 M. Oliveira, J. Jornet-Somoza
2!! Copyright (C) 2020-2021 M. Oliveira, K. Lively, A. Obzhirov, I. Albar
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 atom_oct_m
24 use basins_oct_m
27 use box_oct_m
32 use comm_oct_m
33 use debug_oct_m
34 use global_oct_m
36 use io_oct_m
38 use ions_oct_m
39 use kick_oct_m
40 use, intrinsic :: iso_fortran_env
41 use loct_oct_m
43 use mesh_oct_m
45 use mpi_oct_m
48 use parser_oct_m
52 use space_oct_m
56 use string_oct_m
57 use unit_oct_m
59 use utils_oct_m
61 use xc_oct_m
62
63 implicit none
64
66 class(box_t), pointer :: box
67 character(len=500) :: clist
68 character(len=15) :: lab
69 integer :: dshape
70 logical, allocatable :: mesh_mask(:)
71 logical, allocatable :: ions_mask(:)
72 type(local_write_t) :: writ
73 end type local_domain_t
74
75 type(electrons_t), pointer :: sys
76 integer, parameter :: BADER = 9 ! should be = OPTION__OUTPUT__BADER
77
78 ! Initialize stuff
79 call global_init()
80
81 call parser_init()
82
83 call messages_init()
84
85 call io_init()
87
88 call print_header()
89 call messages_print_with_emphasis(msg="Local Domains mode", namespace=global_namespace)
91
92 call messages_experimental("oct-local_multipoles utility")
93
95
96 call calc_mode_par%set_parallelization(p_strategy_states, default = .false.)
98 call sys%init_parallelization(mpi_world)
99
100 call local_domains()
101
102 safe_deallocate_p(sys)
104 call io_end()
105 call print_date("Calculation ended on ")
106 call messages_end()
107 call parser_end()
108 call global_end()
109
110contains
111 ! -------------
115 subroutine local_domains()
116 integer :: nd
117 type(local_domain_t), allocatable :: loc_domains(:)
118 integer :: err, iter, l_start, l_end, l_step, id
119 integer :: ia, n_spec_def, read_data, iunit, ispec
120 integer :: length, folder_index
121 real(real64) :: default_dt, dt
122 character(len=MAX_PATH_LEN) :: filename, folder, folder_default, radiifile
123 character(len=MAX_PATH_LEN) :: in_folder, restart_folder, basename, frmt, ldrestart_folder
124 character(len=15) :: lab
125 logical :: iterate, ldupdate, ldoverwrite, ldrestart, ld_need_hm
126 type(kick_t) :: kick
127 type(restart_t) :: restart, restart_ld
128
129 push_sub(local_domains)
130
131 message(1) = 'Info: Creating local domains'
132 message(2) = ''
133 call messages_info(2)
134
135 write(folder_default,'(a)') 'restart/gs/'
136
137 !%Variable LDFolder
138 !%Type string
139 !%Section Utilities::oct-local_multipoles
140 !%Description
141 !% The folder name where the density used as input file is.
142 !%End
143 call parse_variable(global_namespace, 'LDFolder', folder_default, folder)
144
145 ! Check if the folder is finished by an /
146 if (index(folder, '/', back = .true.) /= len_trim(folder)) then
147 folder = trim(folder)//'/'
148 end if
149
150 default_dt = m_zero
151 call parse_variable(global_namespace, 'TDTimeStep', default_dt, dt, unit = units_inp%time)
152 if (dt < m_zero) then
153 write(message(1),'(a)') 'Input: TDTimeStep must be positive.'
154 write(message(2),'(a)') 'Input: TDTimeStep reset to 0. Check input file.'
155 call messages_info(2)
156 end if
157
158 !%Variable LDFilename
159 !%Type string
160 !%Default 'density'
161 !%Section Utilities::oct-local_multipoles
162 !%Description
163 !% Input filename. The original filename for the density which is going to be
164 !% fragmented into domains.
165 !%End
166 call parse_variable(global_namespace, 'LDFilename', 'density', basename)
167 if (basename == " ") basename = ""
168 ! Delete the extension if present
169 length = len_trim(basename)
170 if (length > 4) then
171 if (basename(length-3:length) == '.obf') then
172 basename = trim(basename(1:length-4))
173 end if
174 end if
175
176 !%Variable LDUpdate
177 !%Type logical
178 !%Default false
179 !%Section Utilities::oct-local_multipoles
180 !%Description
181 !% Controls if the calculation of the local domains is desired at each iteration.
182 !%End
183 call parse_variable(global_namespace, 'LDUpdate', .false., ldupdate)
184
185 !%Variable LDOverWrite
186 !%Type logical
187 !%Default true
188 !%Section Utilities::oct-local_multipoles
189 !%Description
190 !% Controls whether to over-write existing files.
191 !%End
192 call parse_variable(global_namespace, 'LDOverWrite', .true., ldoverwrite)
193
194 !%Variable LDRadiiFile
195 !%Type string
196 !%Default 'default'
197 !%Section Utilities::oct-local_multipoles
198 !%Description
199 !% Full path for the radii file. If set, def_rsize will be reset to the new values.
200 !% This file should have the same format as share/PP/default.
201 !%End
202 call parse_variable(global_namespace, 'LDRadiiFile', 'default', radiifile)
203
204 if (trim(radiifile) /= "default") then
205 n_spec_def = max(0, loct_number_of_lines(string_f_to_c(radiifile)))
206 if (n_spec_def > 0) n_spec_def = n_spec_def - 1 ! First line is a comment
207
208 ! get parameters from file
209 do ia = 1, sys%ions%nspecies
210 read_data = 0
211 iunit = io_open(radiifile, global_namespace, action='read', status='old', die=.false.)
212 if (iunit /= -1) then
213 if (ia == 1) then
214 write(message(1),'(a,a)') 'Redefining def_rsize from file:', trim(radiifile)
215 call messages_info(1)
216 end if
217 read(iunit,*)
218 default_file: do ispec = 1, n_spec_def
219 read(iunit,*) lab
220 if (trim(lab) == trim(sys%ions%species(ia)%s%get_label())) then
221 select type(spec=>sys%ions%species(ia)%s)
222 class is(pseudopotential_t)
223 call read_from_default_file(iunit, read_data, spec)
224 end select
225 exit default_file
226 end if
227 end do default_file
228 call io_close(iunit)
229 else
230 write(message(1),'(a,a)') 'Octopus could not open then file:', trim(radiifile)
231 call messages_warning(1)
232 end if
233 end do
234 end if
235
236 !TODO: use standart FromSratch and %RestartOptions
237 !%Variable LDRestart
238 !%Type logical
239 !%Default false
240 !%Section Utilities::oct-local_multipoles
241 !%Description
242 !% Restart information will be read from <tt>LDRestartFolder</tt>.
243 !%End
244 call parse_variable(global_namespace, 'LDRestart', .false., ldrestart)
245
246 if (ldrestart) then
247 write(folder_default,'(a)')'ld.general'
248
249 !%Variable LDRestartFolder
250 !%Type string
251 !%Default "ld.general"
252 !%Section Utilities::oct-local_multipoles
253 !%Description
254 !% The folder name where the density used as input file is.
255 !%End
256 call parse_variable(global_namespace, 'LDRestartFolder', folder_default, ldrestart_folder)
257
258 ! Check if the folder is finished by an /
259 if (index(ldrestart_folder, '/', .true.) /= len_trim(ldrestart_folder)) then
260 write(ldrestart_folder,'(a,a1)') trim(ldrestart_folder), '/'
261 end if
262 end if
263
264
265 !%Variable LDIterateFolder
266 !%Type logical
267 !%Default false
268 !%Section Utilities::oct-local_multipoles
269 !%Description
270 !% This variable decides if a folder is going to be iterated.
271 !%End
272 call parse_variable(global_namespace, 'LDIterateFolder', .false., iterate)
273
274 !%Variable LDStart
275 !%Type integer
276 !%Default 0
277 !%Section Utilities::oct-local_multipoles
278 !%Description
279 !% The starting number of the filename or folder.
280 !%End
281 call parse_variable(global_namespace, 'LDStart', 0, l_start)
282
283 !%Variable LDEnd
284 !%Type integer
285 !%Default 0
286 !%Section Utilities::oct-local_multipoles
287 !%Description
288 !% The last number of the filename or folder.
289 !%End
290 call parse_variable(global_namespace, 'LDEnd', 0, l_end)
291
292 !%Variable LDStep
293 !%Type integer
294 !%Default 1
295 !%Section Utilities::oct-local_multipoles
296 !%Description
297 !% The padding between the filenames or folder.
298 !%End
299 call parse_variable(global_namespace, 'LDStep', 1, l_step)
300
301 message(1) = 'Info: Computing local multipoles'
302 message(2) = ''
303 call messages_info(2)
304
305 call local_init(sys%space, sys%gr, sys%ions, nd, loc_domains)
306
307 ! Starting loop over selected densities.
308 if (any(loc_domains(:)%dshape == bader)) then
309 call messages_experimental('Bader volumes in oct-local_multipoles')
310 end if
311
312 call kick_init(kick, global_namespace, sys%space, sys%kpoints, sys%st%d%ispin)
313 ld_need_hm = .false.
314 do id = 1, nd
315 call local_write_init(loc_domains(id)%writ, global_namespace, loc_domains(id)%lab, 0, dt)
316 ld_need_hm = ld_need_hm .or. local_write_check_hm(loc_domains(id)%writ)
317 end do
318
319 if (ld_need_hm) then
320 message(1) = 'Info: Local domains output includes local potential and/or local energy.'
321 call messages_info(1)
322
323 call hamiltonian_elec_epot_generate(sys%hm, global_namespace, sys%space, sys%gr, sys%ions, &
324 sys%ext_partners, sys%st, time=m_zero)
325 call sys%hm%update(sys%gr, global_namespace, sys%space, sys%ext_partners, time=m_zero)
326 end if
327
328 if (ldrestart) then
329 !TODO: check for domains & mesh compatibility
330 call restart_ld%init(global_namespace, restart_undefined, restart_type_load, sys%mc, err, &
331 dir=trim(ldrestart_folder), mesh = sys%gr)
332 call local_restart_read(sys%space, sys%gr, sys%ions, nd, loc_domains, restart_ld)
333 call restart_ld%end()
334 end if
335
336 ! Initialize the restart directory from <tt>LDFolder</tt> value.
337 ! This directory has to have the files 'grid', 'mesh' and 'lxyz.obf'
338 ! and the files that are going to be used, must be inside this folder
339 if (iterate) then
340 ! Delete the last / and find the previous /, if any
341 in_folder = folder(1:len_trim(folder)-1)
342 folder_index = index(in_folder, '/', .true.)
343 restart_folder = in_folder(1:folder_index)
344 else
345 restart_folder = folder
346 end if
347 call restart%init(global_namespace, restart_undefined, restart_type_load, sys%mc, err, &
348 dir=trim(restart_folder), mesh = sys%gr)
349
350!!$ call loct_progress_bar(-1, l_end-l_start)
351 do iter = l_start, l_end, l_step
352 if (iterate) then
353 ! Delete the last / and add the corresponding folder number
354 write(in_folder,'(a,i0.7,a)') folder(folder_index+1:len_trim(folder)-1),iter, "/"
355 write(filename, '(a,a,a)') trim(in_folder), trim(basename)
356 else
357 in_folder = "."
358 if (l_start /= l_end) then
359 ! Here, we are only considering 10 character long filenames.
360 ! Subtract the initial part given at 'ConvertFilename' from the format and pad
361 ! with zeros.
362 write(frmt,'(a,i0,a)')"(a,i0.",10-len_trim(basename),")"
363 write(filename, fmt=trim(frmt)) trim(basename), iter
364 else
365 ! Assuming filename is given complete in the 'LDFilename'
366 write(filename, '(a,a,a,a)') trim(in_folder),"/", trim(basename)
367 filename = basename
368 end if
369 end if
370 ! FIXME: there is a special function for reading the density. Why not use that?
371 ! TODO: Find the function that reads the density. Which one?
372 ! FIXME: why only real functions? Please generalize.
373 ! TODO: up to know the local_multipoles utlity acts over density functions, which are real.
374 if (err == 0) then
375 call restart%read_mesh_function(trim(filename), sys%gr, sys%st%rho(:,1), err)
376 end if
377 if (err /= 0) then
378 write(message(1),*) 'While reading density: "', trim(filename), '", error code:', err
379 call messages_fatal(1)
380 end if
381
382 ! Look for the mesh points inside local domains
383 if ((iter == l_start .and. .not. ldrestart) .or. ldupdate) then
384 call local_inside_domain(sys%space, sys%gr, sys%ions, nd, loc_domains, global_namespace, sys%st%rho(:,1))
385 call local_restart_write(global_namespace, sys%space, sys%gr, sys%mc, nd, loc_domains)
386 end if
387
388 do id = 1, nd
389 call local_write_iter(loc_domains(id)%writ, global_namespace, sys%space, loc_domains(id)%lab, loc_domains(id)%ions_mask, &
390 loc_domains(id)%mesh_mask, sys%gr, sys%st, sys%hm, sys%ks, sys%ions, &
391 sys%ext_partners, kick, iter, l_start, ldoverwrite)
392 end do
393 call loct_progress_bar(iter-l_start, l_end-l_start)
394 end do
395
396 call restart%end()
397 call kick_end(kick)
398
399 do id = 1, nd
400 call local_end(loc_domains(id))
401 end do
402 safe_deallocate_a(loc_domains)
403
404 message(1) = 'Info: Exiting local domains'
405 message(2) = ''
406 call messages_info(2)
407
408 pop_sub(local_domains)
409 end subroutine local_domains
410
411 ! ---------------------------------------------------------
414 ! ---------------------------------------------------------
415 subroutine local_init(space, mesh, ions, nd, loc_domains)
416 class(space_t), intent(in) :: space
417 class(mesh_t), intent(in) :: mesh
418 type(ions_t), intent(in) :: ions
419 integer, intent(out) :: nd
420 type(local_domain_t), allocatable, intent(out) :: loc_domains(:)
421
422 integer :: id
423 type(block_t) :: blk
424
425 push_sub(local_init)
426
427 !TODO: use same strategy used in %Species block to define variables
428 !%Variable LocalDomains
429 !%Type block
430 !%Section Utilities::oct-local_multipoles
431 !%Description
432 !% The LocalDomains are by definition part of the global grid. The domains are defined by
433 !% selecting a type shape. The domain box will be constructed using the given parameters.
434 !% A local domain could be construct by addition of several box centered on the ions.
435 !% The grid points inside this box will belong to the local domain.
436 !%
437 !% The format of this block is the following:<br>
438 !% <tt> 'Label' | Shape | %< | Shape dependencies >% </tt>
439 !% <br>The first field is the label of the domain.
440 !% Label = string with the name of the new local domain.
441 !% The second is the shape type of the box used to define the domain.
442 !% Shape = SPHERE, CYLINDER, PARALLELEPIPED, MINIMUM, BADER.
443 !% Some types may need some parameters given in the remaining fields of the row.
444 !% (the valid options are detailed below).
445 !%
446 !% <tt>%LocalDomains
447 !% <br>case (SPHERE): | rsize | %<dim origin coordinates>
448 !% <br>case (CYLINDER): | rsize | xsize | %<origin coordinates>
449 !% <br>case (PARALLELEPIPED): | %<lsize> | %<origin coordinates>
450 !% <br>case (MINIMUM): | rsize | 'center_list'
451 !% <br>case (BADER): | 'center_list'
452 !% <br>%</tt>
453 !% <br>rsize: Radius in input length units
454 !% <br>xsize: the length of the cylinder in the x-direction
455 !% <br>origin coordinates: in input length units separated by |, where the box is centered.
456 !% <br>lsize: half of the length of the parallelepiped in each direction.
457 !% <br>center_list: string containing the list of atoms in xyz file for each domain in the form "2,16-23"
458 !%End
459
460 ! First, find out if there is a LocalDomains block.
461 nd = 0
462 if (parse_block(global_namespace, 'LocalDomains', blk) == 0) then
463 nd = parse_block_n(blk)
464 end if
465
466 safe_allocate(loc_domains(1:nd))
467
468 block: do id = 1, nd
469 safe_allocate(loc_domains(id)%mesh_mask(1:mesh%np))
470 safe_allocate(loc_domains(id)%ions_mask(1:ions%natoms))
471 call parse_block_string(blk, id-1, 0, loc_domains(id)%lab)
472 call local_read_from_block(loc_domains(id), space, ions, blk, id-1, global_namespace)
473 end do block
474 call parse_block_end(blk)
475 message(1) = ''
476 call messages_info(1)
477
478 pop_sub(local_init)
479 end subroutine local_init
480
481 ! ---------------------------------------------------------
484 ! ---------------------------------------------------------
485 subroutine local_end(domain)
486 type(local_domain_t), intent(inout) :: domain
487
488 class(box_t), pointer :: box
489
490 push_sub(local_end)
491
492 if (domain%dshape /= bader) then
493 box => domain%box
494 safe_deallocate_p(box)
495 end if
496 call local_write_end(domain%writ)
497 safe_deallocate_a(domain%mesh_mask)
498 safe_deallocate_a(domain%ions_mask)
499
500 pop_sub(local_end)
501 end subroutine local_end
502
503 ! ---------------------------------------------------------
504 subroutine local_read_from_block(domain, space, ions, blk, row, namespace)
505 type(local_domain_t), intent(inout) :: domain
506 class(space_t), intent(in) :: space
507 type(ions_t), intent(in) :: ions
508 type(block_t), intent(in) :: blk
509 integer, intent(in) :: row
510 type(namespace_t), intent(in) :: namespace
511
512 integer :: ic, ia
513 real(real64) :: radius, center(space%dim), half_length(space%dim)
514 class(box_union_t), pointer :: minimum_box
515
516 push_sub(local_read_from_block)
517
518 call parse_block_integer(blk, row, 1, domain%dshape)
519
520 select case (domain%dshape)
521 case (minimum)
522 call parse_block_float(blk, row, 2, radius, unit = units_inp%length)
523 if (radius < m_zero) call messages_input_error(namespace, 'radius', row=row, column=2)
524 call parse_block_string(blk, row, 3, domain%clist)
525
526 minimum_box => box_union_t(space%dim)
527 do ia = 1, ions%natoms
528 if (loct_isinstringlist(ia, domain%clist)) then
529 call minimum_box%add_box(box_sphere_t(space%dim, ions%pos(:, ia), radius, namespace))
530 end if
531 end do
532 domain%box => minimum_box
533
534 ic = 0
535 do ia = 1, ions%natoms
536 if (domain%box%contains_point(ions%pos(:, ia)) .and. .not. loct_isinstringlist(ia, domain%clist)) then
537 ic = ic + 1
538 if (ic <= 20) write(message(ic),'(a,a,I0,a,a)')'Atom: ',trim(ions%atom(ia)%species%get_label()), ia, &
539 ' is inside the union box BUT not in list: ', trim(domain%clist)
540 end if
541 end do
542 if (ic > 0) then
543 if (ic > 20) ic = 20
544 call messages_warning(ic)
545 message(1) = ' THIS COULD GIVE INCORRECT RESULTS '
546 call messages_warning(1)
547 if (ic == 20) then
548 message(1) = ' AT LEAST 19 ATOMS ARE NOT PRESENT IN THE LIST'
549 call messages_warning(1)
550 end if
551 end if
552
553 case (sphere)
554 call parse_block_float(blk, row, 2, radius, unit = units_inp%length)
555 if (radius < m_zero) call messages_input_error(namespace, 'radius', row=row, column=2)
556 do ic = 1, space%dim
557 call parse_block_float(blk, row, 2 + ic, center(ic), unit = units_inp%length)
558 end do
559
560 domain%box => box_sphere_t(space%dim, center, radius, namespace)
561
562 case (cylinder)
563 call parse_block_float(blk, row, 2, radius, unit = units_inp%length)
564 if (radius < m_zero) then
565 call messages_input_error(namespace, 'radius', row=row, column=2)
566 end if
567 call parse_block_float(blk, row, 3, half_length(1), unit = units_inp%length)
568 do ic = 1, space%dim
569 call parse_block_float(blk, row, 3 + ic, center(ic), unit = units_inp%length)
570 end do
571
572 domain%box => box_cylinder_t(space%dim, center, ions%latt%rlattice, radius, half_length(1)*m_two, namespace)
573
574 case (parallelepiped)
575 do ic = 1, space%dim
576 call parse_block_float(blk, row, 1 + ic, half_length(ic), unit = units_inp%length)
577 end do
578 do ic = 1, space%dim
579 call parse_block_float(blk, row, 1 + space%dim + ic, center(ic), unit = units_inp%length)
580 end do
581
582 domain%box => box_parallelepiped_t(space%dim, center, ions%latt%rlattice, half_length*m_two, namespace)
583
584 case (bader)
585 call parse_block_string(blk, row, 2, domain%clist)
586 pop_sub(local_read_from_block)
587 return
588 end select
589
590 pop_sub(local_read_from_block)
591 end subroutine local_read_from_block
592
593 ! ---------------------------------------------------------
595 subroutine local_restart_write(namespace, space, mesh, mc, nd, loc_domains)
596 type(namespace_t), intent(in) :: namespace
597 class(space_t), intent(in) :: space
598 class(mesh_t), intent(in) :: mesh
599 type(multicomm_t), intent(in) :: mc
600 integer, intent(in) :: nd
601 type(local_domain_t), intent(in) :: loc_domains(:)
602
603 integer :: id, ip, iunit, ierr
604 character(len=MAX_PATH_LEN) :: dirname
605 character(len=140), allocatable :: lines(:)
606 real(real64), allocatable :: ff(:,:)
607 type(restart_t) :: restart
608
609 push_sub(local_restart_write)
610
611 dirname = "./restart/ld"
612 write(message(1),'(a,a)')'Info: Writing restart info to ', trim(dirname)
613 call messages_info(1)
614
615 call restart%init(namespace, restart_undefined, restart_type_dump, mc, ierr, mesh=mesh, dir=trim(dirname))
616
617 safe_allocate(lines(1:nd+2))
618 write(lines(1),'(a,1x,i5)') 'Number of local domains =', nd
619 write(lines(2),'(a3,1x,a15,1x,a5,1x)') '#id', 'label', 'shape'
620 do id = 1, nd
621 write(lines(id+2), '(i3,1x,a15,1x,i5)') id, trim(loc_domains(id)%lab), loc_domains(id)%dshape
622 end do
623 iunit = restart%open("ldomains.info")
624 call restart%write(iunit, lines, nd+2, ierr)
625 call io_close(iunit)
626 safe_deallocate_a(lines)
627
628 safe_allocate(ff(1:mesh%np, 1))
629 ff = m_zero
630 do id = 1, nd
631 do ip = 1, mesh%np
632 if (loc_domains(id)%mesh_mask(ip)) ff(ip, 1) = ff(ip, 1) + 2**real(id, real64)
633 end do
634 end do
635 call restart%write_mesh_function("ldomains", mesh, ff(1:mesh%np, 1), ierr)
636 safe_deallocate_a(ff)
637
638 call restart%end()
639
640 pop_sub(local_restart_write)
641 end subroutine local_restart_write
642
643 ! ---------------------------------------------------------
644 subroutine local_restart_read(space, mesh, ions, nd, loc_domains, restart)
645 class(space_t), intent(in) :: space
646 class(mesh_t), intent(in) :: mesh
647 type(ions_t), intent(in) :: ions
648 integer, intent(in) :: nd
649 type(local_domain_t), intent(inout) :: loc_domains(:)
650 type(restart_t), intent(in) :: restart
651
652 integer :: id, ip, iunit, ierr
653 character(len=31) :: line(1), tmp
654 real(real64), allocatable :: mask(:)
655
656 push_sub(local_restart_read)
657
658 message(1) = 'Info: Reading mesh points inside each local domain'
659 call messages_info(1)
660
661 safe_allocate(mask(1:mesh%np))
662
663 !Read local domain information from ldomains.info
664 iunit = restart%open("ldomains.info", status='old')
665 call restart%read(iunit, line, 1, ierr)
666 read(line(1),'(a25,1x,i5)') tmp, ierr
667 call restart%close(iunit)
668
669 call restart%read_mesh_function("ldomains", mesh, mask, ierr)
670
671 do id = 1, nd
672 loc_domains(id)%mesh_mask = .false.
673 do ip = 1 , mesh%np
674 if (bitand(int(mask(ip)), 2**id) /= 0) loc_domains(id)%mesh_mask(ip) = .true.
675 end do
676
677 !Check for atom list inside each domain
678 call local_ions_mask(loc_domains(id)%mesh_mask, ions, mesh, loc_domains(id)%ions_mask)
679 end do
680
681 safe_deallocate_a(mask)
682
683 pop_sub(local_restart_read)
684 end subroutine local_restart_read
685
686 ! ---------------------------------------------------------
687 subroutine local_inside_domain(space, mesh, ions, nd, loc_domains, namespace, ff)
688 class(space_t), intent(in) :: space
689 class(mesh_t), intent(in) :: mesh
690 type(ions_t), intent(in) :: ions
691 integer, intent(in) :: nd
692 type(local_domain_t), intent(inout) :: loc_domains(:)
693 type(namespace_t), intent(in) :: namespace
694 real(real64), intent(in) :: ff(:)
695
696 integer(int64) :: how
697 integer :: id, ip, ierr
698 type(basins_t) :: basins
699 character(len=MAX_PATH_LEN) :: filename
700 real(real64), allocatable :: dble_domain_map(:), domain_mesh(:)
701
702 push_sub(local_inside_domain)
703
704 !TODO: find a parallel algorithm to perform the charge-density fragmentation.
705 message(1) = 'Info: Assigning mesh points inside each local domain'
706 call messages_info(1)
707
708 ! If we have Bader domains, then we need to get the basins
709 if (any(loc_domains(:)%dshape == bader)) then
710 if (mesh%parallel_in_domains) then
711 write(message(1),'(a)') 'Bader volumes can only be computed in serial'
712 call messages_fatal(1)
713 end if
714 call create_basins(namespace, space, mesh, ions, ff, basins)
715 end if
716
717 do id = 1, nd
718 ! Create the mask that tells which mesh points are inside the domain
719 select case (loc_domains(id)%dshape)
720 case (bader)
721 call bader_domain_create_mask(loc_domains(id), basins, mesh, ions)
722 case default
723 call box_domain_create_mask(loc_domains(id), mesh)
724 end select
725
726 !Check for atom list inside each domain
727 call local_ions_mask(loc_domains(id)%mesh_mask, ions, mesh, loc_domains(id)%ions_mask)
728 end do
729
730 if (any(loc_domains(:)%dshape == bader)) then
731 call basins_end(basins)
732 end if
733
734 if (debug%info) then
735 call parse_variable(namespace, 'LDOutputFormat', 0, how)
736 if (.not. varinfo_valid_option('OutputFormat', how, is_flag=.true.)) then
737 call messages_input_error(namespace, 'LDOutputFormat')
738 end if
740 ! Write extra information about domains
741 safe_allocate(dble_domain_map(1:mesh%np))
742 safe_allocate(domain_mesh(1:mesh%np))
743
744 domain_mesh = m_zero
745 do id = 1, nd
746 dble_domain_map = m_zero
747 do ip = 1, mesh%np
748 if (loc_domains(id)%mesh_mask(ip)) then
749 dble_domain_map(ip) = real(id, real64)
750 domain_mesh(ip) = domain_mesh(ip) + dble_domain_map(ip)
751 end if
752 end do
753
754 write(filename,'(a,a)') 'domain.', trim(loc_domains(id)%lab)
755 call dio_function_output(how, 'local.general', trim(filename), namespace, space, mesh, dble_domain_map, unit_one, ierr, &
756 pos=ions%pos, atoms=ions%atom)
757 end do
758
759 call dio_function_output(how, 'local.general', 'domain.mesh', namespace, space, mesh, domain_mesh, unit_one, ierr, &
760 pos=ions%pos, atoms=ions%atom)
761
762 safe_deallocate_a(dble_domain_map)
763 safe_deallocate_a(domain_mesh)
764 end if
765
766 pop_sub(local_inside_domain)
767 end subroutine local_inside_domain
768
769 ! ---------------------------------------------------------
770 subroutine create_basins(namespace, space, mesh, ions, ff, basins)
771 type(namespace_t), intent(in) :: namespace
772 class(space_t), intent(in) :: space
773 class(mesh_t), intent(in) :: mesh
774 type(ions_t), intent(in) :: ions
775 real(real64), intent(in) :: ff(:)
776 type(basins_t), intent(inout) :: basins
777
778 logical :: use_atomic_radii
779 integer :: ip, ierr, ia, is, rankmin
780 integer(int64) :: how
781 real(real64) :: bader_threshold, dmin, dd
782 integer, allocatable :: ion_map(:)
783 real(real64), allocatable :: ff2(:,:), ffs(:)
784
785 push_sub(create_basins)
786
787 !%Variable LDBaderThreshold
788 !%Type float
789 !%Default 0.01
790 !%Section Utilities::oct-local_multipoles
791 !%Description
792 !% This variable sets the threshold for the basins calculations. Recommended values:
793 !% 0.01 -> intramolecular volumes; 0.2 -> intermolecular volumes.
794 !%End
795 call parse_variable(namespace, 'LDBaderThreshold', 0.01_real64, bader_threshold)
796
797 ! Copy function used to construct the basins, as we need to modify it
798 safe_allocate(ff2(1:mesh%np,1))
799 ff2(1:mesh%np,1) = ff(1:mesh%np)
800
801 ! Add long range density from atoms
802 safe_allocate(ffs(1:mesh%np))
803 do ia = 1, ions%natoms
804 call species_get_long_range_density(ions%atom(ia)%species, ions%namespace, ions%space, ions%latt, &
805 ions%pos(:, ia), mesh, ffs)
806 do is = 1, ubound(ff2, dim=2)
807 ff2(1:mesh%np, is) = ff2(1:mesh%np, is) - ffs(1:mesh%np)
808 end do
809 end do
810 safe_deallocate_a(ffs)
811
812 call basins_init(basins, namespace, mesh)
813 call basins_analyze(basins, namespace, mesh, ff2(:,1), ff2, bader_threshold)
814
815 if (debug%info) then
816 call parse_variable(namespace, 'LDOutputFormat', 0, how)
817 if (.not. varinfo_valid_option('OutputFormat', how, is_flag=.true.)) then
818 call messages_input_error(namespace, 'LDOutputFormat')
819 end if
820
821 call dio_function_output(how, 'local.general', 'basinsmap', namespace, space, mesh, real(basins%map(1:mesh%np), real64) , &
822 unit_one, ierr, pos=ions%pos, atoms=ions%atom)
823 call dio_function_output(how, 'local.general', 'dens_ff2', namespace, space, mesh, ff2(:,1), &
824 unit_one, ierr, pos=ions%pos, atoms=ions%atom)
825 end if
826 safe_deallocate_a(ff2)
827
828 !%Variable LDUseAtomicRadii
829 !%Type logical
830 !%Default false
831 !%Section Utilities::oct-local_multipoles
832 !%Description
833 !% If set, atomic radii will be used to assign lone pairs to ion.
834 !%End
835 call parse_variable(global_namespace, 'LDUseAtomicRadii', .false., use_atomic_radii)
836
837 if (use_atomic_radii) then
838 safe_allocate(ion_map(1:ions%natoms))
839 ion_map = 0
840 do ia = 1, ions%natoms
841 ion_map(ia) = basins%map(mesh_nearest_point(mesh, ions%pos(:, ia), dmin, rankmin))
842 end do
843
844 do ip = 1, mesh%np
845 ! Check if lonely pair has no ions assigned to it
846 if (all(ion_map(:) /= basins%map(ip))) then
847 ! Assign lonely pair to ion in a atomic radii distance
848 do ia = 1, ions%natoms
849 dd = norm2(mesh%x(:, ip) - ions%pos(:, ia))
850 if (dd <= ions%atom(ia)%species%get_vdw_radius()) basins%map(ip) = ion_map(ia)
851 end do
852 end if
853 end do
854
855 safe_deallocate_a(ion_map)
856 end if
857
858 pop_sub(create_basins)
859 end subroutine create_basins
860
861 ! ---------------------------------------------------------
862 subroutine box_domain_create_mask(domain, mesh)
863 class(local_domain_t), intent(inout) :: domain
864 class(mesh_t), intent(in) :: mesh
866 push_sub(box_domain_create_mask)
867
868 domain%mesh_mask = domain%box%contains_points(mesh%np, mesh%x_t)
869
871 end subroutine box_domain_create_mask
872
873 ! ---------------------------------------------------------
874 subroutine bader_domain_create_mask(domain, basins, mesh, ions)
875 class(local_domain_t), intent(inout) :: domain
876 type(basins_t), intent(in) :: basins
877 class(mesh_t), intent(in) :: mesh
878 type(ions_t), intent(in) :: ions
879
880 integer :: ia, ib, ip, ix, n_basins, rankmin
881 real(real64) :: dmin
882 integer, allocatable :: domain_map(:)
883
885
886 ! Count then number of basins to consider. That is the number of ions specified in the input by the user
887 n_basins = 0
888 do ia = 1, ions%natoms
889 if (loct_isinstringlist(ia, domain%clist)) n_basins = n_basins + 1
890 end do
891
892 safe_allocate(domain_map(1:n_basins))
893
894 ! Assign basins to ions
895 domain_map = 0
896 ib = 0
897 do ia = 1, ions%natoms
898 if (.not. loct_isinstringlist(ia, domain%clist)) cycle
899 ib = ib + 1
900 ix = mesh_nearest_point(mesh, ions%pos(:, ia), dmin, rankmin)
901 domain_map(ib) = basins%map(ix)
902 end do
903
904 ! Now construct the mask: a point belongs to this Bader domain if it is part
905 ! of at least one basin that is part of this domain
906 do ip = 1, mesh%np
907 domain%mesh_mask(ip) = any(domain_map(1:n_basins) == basins%map(ip))
908 end do
909
910 safe_deallocate_a(domain_map)
911
913 end subroutine bader_domain_create_mask
914
915 ! ---------------------------------------------------------
916 !Check for the ions inside each local domain.
917 subroutine local_ions_mask(mesh_mask, ions, mesh, ions_mask)
918 logical, intent(in) :: mesh_mask(:)
919 type(ions_t), intent(in) :: ions
920 class(mesh_t), intent(in) :: mesh
921 logical, intent(out) :: ions_mask(:)
922
923 integer :: ia, ix
924 integer, allocatable :: mask_tmp(:)
925 integer :: rankmin
926 real(real64) :: dmin
927
928 push_sub(local_ions_mask)
929
930 safe_allocate(mask_tmp(1:ions%natoms))
931 mask_tmp = 0
932
933 do ia = 1, ions%natoms
934 ix = mesh_nearest_point(mesh, ions%pos(:, ia), dmin, rankmin)
935 if (rankmin /= mesh%mpi_grp%rank) cycle
936 if (mesh_mask(ix)) mask_tmp(ia) = 1
937 end do
938
939 call mesh%allreduce(mask_tmp, ions%natoms)
940 ions_mask = mask_tmp == 1
941
942 safe_deallocate_a(mask_tmp)
943
944 pop_sub(local_ions_mask)
945 end subroutine local_ions_mask
946
947end program oct_local_multipoles
948
949!! Local Variables:
950!! mode: f90
951!! coding: utf-8
952!! End:
subroutine create_basins(namespace, space, mesh, ions, ff, basins)
program oct_local_multipoles
subroutine local_restart_read(space, mesh, ions, nd, loc_domains, restart)
subroutine local_ions_mask(mesh_mask, ions, mesh, ions_mask)
subroutine local_restart_write(namespace, space, mesh, mc, nd, loc_domains)
Write restart files for local domains.
subroutine local_inside_domain(space, mesh, ions, nd, loc_domains, namespace, ff)
subroutine bader_domain_create_mask(domain, basins, mesh, ions)
subroutine local_read_from_block(domain, space, ions, blk, row, namespace)
subroutine local_domains()
Reads a global grid and density and make local domains This is a high-level interface that reads the ...
subroutine local_init(space, mesh, ions, nd, loc_domains)
Initialize local_domain_t variable, allocating variable and reading parameters from input file.
subroutine local_end(domain)
Ending local_domain_t variable, allocating variable and reading parameters from input file.
subroutine box_domain_create_mask(domain, mesh)
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
Module, implementing a factory for boxes.
integer, parameter, public sphere
integer, parameter, public minimum
integer, parameter, public parallelepiped
integer, parameter, public cylinder
This module handles the calculation mode.
type(calc_mode_par_t), public calc_mode_par
Singleton instance of parallel calculation mode.
integer, parameter, public p_strategy_states
parallelization in states
type(debug_t), save, public debug
Definition: debug.F90:158
real(real64), parameter, public m_two
Definition: global.F90:193
subroutine, public global_end()
Finalise parser varinfo file, and MPI.
Definition: global.F90:458
real(real64), parameter, public m_zero
Definition: global.F90:191
subroutine, public global_init(communicator)
Initialise Octopus.
Definition: global.F90:365
integer, parameter, public length
subroutine, public hamiltonian_elec_epot_generate(this, namespace, space, gr, ions, ext_partners, st, time)
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_init(defaults)
If the argument defaults is present and set to true, then the routine will not try to read anything f...
Definition: io.F90:165
subroutine, public io_close(iunit, grp)
Definition: io.F90:467
subroutine, public io_end()
Definition: io.F90:271
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:402
subroutine, public kick_end(kick)
Definition: kick.F90:796
subroutine, public kick_init(kick, namespace, space, kpoints, nspin)
Definition: kick.F90:225
subroutine, public local_write_init(writ, namespace, lab, iter, dt)
subroutine, public local_write_end(writ)
logical function, public local_write_check_hm(writ)
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)
System information (time, memory, sysname)
Definition: loct.F90:117
logical function, public loct_isinstringlist(a, s)
Definition: loct.F90:288
subroutine, public loct_progress_bar(a, maxcount)
A wrapper around the progress bar, such that it can be silenced without needing to dress the call wit...
Definition: loct.F90:276
character(kind=c_char, len=1) function, dimension(len_trim(f_string)+1), private string_f_to_c(f_string)
convert a Fortran string to a C string
Definition: loct.F90:240
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
integer function, public mesh_nearest_point(mesh, pos, dmin, rankmin)
Returns the index of the point which is nearest to a given vector position pos.
Definition: mesh.F90:386
subroutine, public messages_end()
Definition: messages.F90:272
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:897
character(len=512), private msg
Definition: messages.F90:166
subroutine, public messages_init(output_dir)
Definition: messages.F90:219
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:524
subroutine, public print_date(str)
Definition: messages.F90:982
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
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
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:272
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:147
type(namespace_t), public global_namespace
Definition: namespace.F90:135
subroutine, public parse_block_string(blk, l, c, res, convert_to_c)
Definition: parser.F90:810
subroutine, public parser_init()
Initialise the Octopus parser.
Definition: parser.F90:402
subroutine, public parser_end()
End the Octopus parser.
Definition: parser.F90:434
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:615
subroutine, public profiling_end(namespace)
Definition: profiling.F90:415
subroutine, public profiling_init(namespace)
Create profiling subdirectory.
Definition: profiling.F90:257
subroutine, public read_from_default_file(iunit, read_data, spec)
integer, parameter, public restart_undefined
Definition: restart.F90:156
integer, parameter, public restart_type_dump
Definition: restart.F90:183
integer, parameter, public restart_type_load
Definition: restart.F90:183
subroutine, public species_get_long_range_density(species, namespace, space, latt, pos, mesh, rho, sphere_inout, nlr_x)
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.
subroutine, public unit_system_init(namespace)
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
subroutine, public print_header()
This subroutine prints the logo followed by information about the compilation and the system....
Definition: utils.F90:304
Definition: xc.F90:116
Class implementing a cylinder box. The cylinder axis is always along the first direction defined by t...
class to tell whether a point is inside or outside
Definition: box.F90:143
Class implementing a parallelepiped box. Currently this is restricted to a rectangular cuboid (all th...
Class implementing a spherical box.
Definition: box_sphere.F90:134
Class implementing a box that is an union other boxes.
Definition: box_union.F90:135
Class describing the electron system.
Definition: electrons.F90:220
Describes mesh distribution to nodes.
Definition: mesh.F90:187
Stores all communicators and groups.
Definition: multicomm.F90:208
int true(void)