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 unit_oct_m
58 use utils_oct_m
60 use xc_oct_m
61
62 implicit none
63
65 class(box_t), pointer :: box
66 character(len=500) :: clist
67 character(len=15) :: lab
68 integer :: dshape
69 logical, allocatable :: mesh_mask(:)
70 logical, allocatable :: ions_mask(:)
71 type(local_write_t) :: writ
72 end type local_domain_t
73
74 type(electrons_t), pointer :: sys
75 integer, parameter :: BADER = 9 ! should be = OPTION__OUTPUT__BADER
76
77 ! Initialize stuff
78 call global_init()
79
80 call parser_init()
81
82 call messages_init()
83
84 call io_init()
86
87 call print_header()
88 call messages_print_with_emphasis(msg="Local Domains mode", namespace=global_namespace)
90
91 call messages_experimental("oct-local_multipoles utility")
92
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
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
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(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 > 0) 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 do id = 1, nd
314 call local_write_init(loc_domains(id)%writ, global_namespace, loc_domains(id)%lab, 0, dt)
315 end do
316
317 !TODO: initialize hamiltonian if needed: check for LDOuput = energy or potential, using local_write_check_hm(local%writ)
318
319 if (ldrestart) then
320 !TODO: check for domains & mesh compatibility
321 call restart_init(restart_ld, global_namespace, restart_undefined, restart_type_load, sys%mc, err, &
322 dir=trim(ldrestart_folder), mesh = sys%gr)
323 call local_restart_read(sys%space, sys%gr, sys%ions, nd, loc_domains, restart_ld)
324 call restart_end(restart_ld)
325 end if
326
327 ! Initialize the restart directory from <tt>LDFolder</tt> value.
328 ! This directory has to have the files 'grid', 'mesh' and 'lxyz.obf'
329 ! and the files that are going to be used, must be inside this folder
330 if (iterate) then
331 ! Delete the last / and find the previous /, if any
332 in_folder = folder(1:len_trim(folder)-1)
333 folder_index = index(in_folder, '/', .true.)
334 restart_folder = in_folder(1:folder_index)
335 else
336 restart_folder = folder
337 end if
339 dir=trim(restart_folder), mesh = sys%gr)
340
341!!$ call loct_progress_bar(-1, l_end-l_start)
342 do iter = l_start, l_end, l_step
343 if (iterate) then
344 ! Delete the last / and add the corresponding folder number
345 write(in_folder,'(a,i0.7,a)') folder(folder_index+1:len_trim(folder)-1),iter, "/"
346 write(filename, '(a,a,a)') trim(in_folder), trim(basename)
347 else
348 in_folder = "."
349 if (l_start /= l_end) then
350 ! Here, we are only considering 10 character long filenames.
351 ! Subtract the initial part given at 'ConvertFilename' from the format and pad
352 ! with zeros.
353 write(frmt,'(a,i0,a)')"(a,i0.",10-len_trim(basename),")"
354 write(filename, fmt=trim(frmt)) trim(basename), iter
355 else
356 ! Assuming filename is given complete in the 'LDFilename'
357 write(filename, '(a,a,a,a)') trim(in_folder),"/", trim(basename)
358 filename = basename
359 end if
360 end if
361 ! FIXME: there is a special function for reading the density. Why not use that?
362 ! TODO: Find the function that reads the density. Which one?
363 ! FIXME: why only real functions? Please generalize.
364 ! TODO: up to know the local_multipoles utlity acts over density functions, which are real.
365 if (err == 0) then
366 call restart_read_mesh_function(restart, sys%space, trim(filename), sys%gr, sys%st%rho(:,1), err)
367 end if
368 if (err /= 0) then
369 write(message(1),*) 'While reading density: "', trim(filename), '", error code:', err
370 call messages_fatal(1)
371 end if
372
373 ! Look for the mesh points inside local domains
374 if ((iter == l_start .and. .not. ldrestart) .or. ldupdate) then
375 call local_inside_domain(sys%space, sys%gr, sys%ions, nd, loc_domains, global_namespace, sys%st%rho(:,1))
376 call local_restart_write(global_namespace, sys%space, sys%gr, sys%mc, nd, loc_domains)
377 end if
378
379 do id = 1, nd
380 call local_write_iter(loc_domains(id)%writ, global_namespace, sys%space, loc_domains(id)%lab, loc_domains(id)%ions_mask, &
381 loc_domains(id)%mesh_mask, sys%gr, sys%st, sys%hm, sys%ks, sys%ions, &
382 sys%ext_partners, kick, iter, l_start, ldoverwrite)
383 end do
384 call loct_progress_bar(iter-l_start, l_end-l_start)
385 end do
386
387 call restart_end(restart)
388 call kick_end(kick)
389
390 do id = 1, nd
391 call local_end(loc_domains(id))
392 end do
393 safe_deallocate_a(loc_domains)
394
395 message(1) = 'Info: Exiting local domains'
396 message(2) = ''
397 call messages_info(2)
398
399 pop_sub(local_domains)
400 end subroutine local_domains
401
402 ! ---------------------------------------------------------
405 ! ---------------------------------------------------------
406 subroutine local_init(space, mesh, ions, nd, loc_domains)
407 class(space_t), intent(in) :: space
408 class(mesh_t), intent(in) :: mesh
409 type(ions_t), intent(in) :: ions
410 integer, intent(out) :: nd
411 type(local_domain_t), allocatable, intent(out) :: loc_domains(:)
412
413 integer :: id
414 type(block_t) :: blk
415
416 push_sub(local_init)
417
418 !TODO: use same strategy used in %Species block to define variables
419 !%Variable LocalDomains
420 !%Type block
421 !%Section Utilities::oct-local_multipoles
422 !%Description
423 !% The LocalDomains are by definition part of the global grid. The domains are defined by
424 !% selecting a type shape. The domain box will be constructed using the given parameters.
425 !% A local domain could be construct by addition of several box centered on the ions.
426 !% The grid points inside this box will belong to the local domain.
427 !%
428 !% The format of this block is the following:<br>
429 !% <tt> 'Label' | Shape | %< | Shape dependencies >% </tt>
430 !% <br>The first field is the label of the domain.
431 !% Label = string with the name of the new local domain.
432 !% The second is the shape type of the box used to define the domain.
433 !% Shape = SPHERE, CYLINDER, PARALLELEPIPED, MINIMUM, BADER.
434 !% Some types may need some parameters given in the remaining fields of the row.
435 !% (the valid options are detailed below).
436 !%
437 !% <tt>%LocalDomains
438 !% <br>case (SPHERE): | rsize | %<dim origin coordinates>
439 !% <br>case (CYLINDER): | rsize | xsize | %<origin coordinates>
440 !% <br>case (PARALLELEPIPED): | %<lsize> | %<origin coordinates>
441 !% <br>case (MINIMUM): | rsize | 'center_list'
442 !% <br>case (BADER): | 'center_list'
443 !% <br>%</tt>
444 !% <br>rsize: Radius in input length units
445 !% <br>xsize: the length of the cylinder in the x-direction
446 !% <br>origin coordinates: in input length units separated by |, where the box is centered.
447 !% <br>lsize: half of the length of the parallelepiped in each direction.
448 !% <br>center_list: string containing the list of atoms in xyz file for each domain in the form "2,16-23"
449 !%End
450
451 ! First, find out if there is a LocalDomains block.
452 nd = 0
453 if (parse_block(global_namespace, 'LocalDomains', blk) == 0) then
454 nd = parse_block_n(blk)
455 end if
456
457 safe_allocate(loc_domains(1:nd))
458
459 block: do id = 1, nd
460 safe_allocate(loc_domains(id)%mesh_mask(1:mesh%np))
461 safe_allocate(loc_domains(id)%ions_mask(1:ions%natoms))
462 call parse_block_string(blk, id-1, 0, loc_domains(id)%lab)
463 call local_read_from_block(loc_domains(id), space, ions, blk, id-1, global_namespace)
464 end do block
465 call parse_block_end(blk)
466 message(1) = ''
467 call messages_info(1)
468
469 pop_sub(local_init)
470 end subroutine local_init
471
472 ! ---------------------------------------------------------
475 ! ---------------------------------------------------------
476 subroutine local_end(domain)
477 type(local_domain_t), intent(inout) :: domain
478
479 class(box_t), pointer :: box
480
481 push_sub(local_end)
482
483 if (domain%dshape /= bader) then
484 box => domain%box
485 safe_deallocate_p(box)
486 end if
487 call local_write_end(domain%writ)
488 safe_deallocate_a(domain%mesh_mask)
489 safe_deallocate_a(domain%ions_mask)
490
491 pop_sub(local_end)
492 end subroutine local_end
493
494 ! ---------------------------------------------------------
495 subroutine local_read_from_block(domain, space, ions, blk, row, namespace)
496 type(local_domain_t), intent(inout) :: domain
497 class(space_t), intent(in) :: space
498 type(ions_t), intent(in) :: ions
499 type(block_t), intent(in) :: blk
500 integer, intent(in) :: row
501 type(namespace_t), intent(in) :: namespace
502
503 integer :: ic, ia
504 real(real64) :: radius, center(space%dim), half_length(space%dim)
505 class(box_union_t), pointer :: minimum_box
506
507 push_sub(local_read_from_block)
508
509 call parse_block_integer(blk, row, 1, domain%dshape)
510
511 select case (domain%dshape)
512 case (minimum)
513 call parse_block_float(blk, row, 2, radius, unit = units_inp%length)
514 if (radius < m_zero) call messages_input_error(namespace, 'radius', row=row, column=2)
515 call parse_block_string(blk, row, 3, domain%clist)
516
517 minimum_box => box_union_t(space%dim)
518 do ia = 1, ions%natoms
519 if (loct_isinstringlist(ia, domain%clist)) then
520 call minimum_box%add_box(box_sphere_t(space%dim, ions%pos(:, ia), radius, namespace))
521 end if
522 end do
523 domain%box => minimum_box
524
525 ic = 0
526 do ia = 1, ions%natoms
527 if (domain%box%contains_point(ions%pos(:, ia)) .and. .not. loct_isinstringlist(ia, domain%clist)) then
528 ic = ic + 1
529 if (ic <= 20) write(message(ic),'(a,a,I0,a,a)')'Atom: ',trim(ions%atom(ia)%species%get_label()), ia, &
530 ' is inside the union box BUT not in list: ', trim(domain%clist)
531 end if
532 end do
533 if (ic > 0) then
534 if (ic > 20) ic = 20
535 call messages_warning(ic)
536 message(1) = ' THIS COULD GIVE INCORRECT RESULTS '
537 call messages_warning(1)
538 if (ic == 20) then
539 message(1) = ' AT LEAST 19 ATOMS ARE NOT PRESENT IN THE LIST'
540 call messages_warning(1)
541 end if
542 end if
543
544 case (sphere)
545 call parse_block_float(blk, row, 2, radius, unit = units_inp%length)
546 if (radius < m_zero) call messages_input_error(namespace, 'radius', row=row, column=2)
547 do ic = 1, space%dim
548 call parse_block_float(blk, row, 2 + ic, center(ic), unit = units_inp%length)
549 end do
550
551 domain%box => box_sphere_t(space%dim, center, radius, namespace)
552
553 case (cylinder)
554 call parse_block_float(blk, row, 2, radius, unit = units_inp%length)
555 if (radius < m_zero) then
556 call messages_input_error(namespace, 'radius', row=row, column=2)
557 end if
558 call parse_block_float(blk, row, 3, half_length(1), unit = units_inp%length)
559 do ic = 1, space%dim
560 call parse_block_float(blk, row, 3 + ic, center(ic), unit = units_inp%length)
561 end do
562
563 domain%box => box_cylinder_t(space%dim, center, ions%latt%rlattice, radius, half_length(1)*m_two, namespace)
564
565 case (parallelepiped)
566 do ic = 1, space%dim
567 call parse_block_float(blk, row, 1 + ic, half_length(ic), unit = units_inp%length)
568 end do
569 do ic = 1, space%dim
570 call parse_block_float(blk, row, 1 + space%dim + ic, center(ic), unit = units_inp%length)
571 end do
572
573 domain%box => box_parallelepiped_t(space%dim, center, ions%latt%rlattice, half_length*m_two, namespace)
574
575 case (bader)
576 call parse_block_string(blk, row, 2, domain%clist)
577 pop_sub(local_read_from_block)
578 return
579 end select
580
581 pop_sub(local_read_from_block)
582 end subroutine local_read_from_block
583
584 ! ---------------------------------------------------------
586 subroutine local_restart_write(namespace, space, mesh, mc, nd, loc_domains)
587 type(namespace_t), intent(in) :: namespace
588 class(space_t), intent(in) :: space
589 class(mesh_t), intent(in) :: mesh
590 type(multicomm_t), intent(in) :: mc
591 integer, intent(in) :: nd
592 type(local_domain_t), intent(in) :: loc_domains(:)
593
594 integer :: id, ip, iunit, ierr
595 character(len=MAX_PATH_LEN) :: dirname
596 character(len=140), allocatable :: lines(:)
597 real(real64), allocatable :: ff(:,:)
598 type(restart_t) :: restart
599
600 push_sub(local_restart_write)
601
602 dirname = "./restart/ld"
603 write(message(1),'(a,a)')'Info: Writing restart info to ', trim(dirname)
604 call messages_info(1)
605
606 call restart_init(restart, namespace, restart_undefined, restart_type_dump, mc, ierr, mesh=mesh, dir=trim(dirname))
607
608 safe_allocate(lines(1:nd+2))
609 write(lines(1),'(a,1x,i5)') 'Number of local domains =', nd
610 write(lines(2),'(a3,1x,a15,1x,a5,1x)') '#id', 'label', 'shape'
611 do id = 1, nd
612 write(lines(id+2), '(i3,1x,a15,1x,i5)') id, trim(loc_domains(id)%lab), loc_domains(id)%dshape
613 end do
614 iunit = restart_open(restart, "ldomains.info")
615 call restart_write(restart, iunit, lines, nd+2, ierr)
616 call io_close(iunit)
617 safe_deallocate_a(lines)
618
619 safe_allocate(ff(1:mesh%np, 1))
620 ff = m_zero
621 do id = 1, nd
622 do ip = 1, mesh%np
623 if (loc_domains(id)%mesh_mask(ip)) ff(ip, 1) = ff(ip, 1) + 2**real(id, real64)
624 end do
625 end do
626 call restart_write_mesh_function(restart, space, "ldomains", mesh, ff(1:mesh%np, 1), ierr)
627 safe_deallocate_a(ff)
628
629 call restart_end(restart)
630
631 pop_sub(local_restart_write)
632 end subroutine local_restart_write
633
634 ! ---------------------------------------------------------
635 subroutine local_restart_read(space, mesh, ions, nd, loc_domains, restart)
636 class(space_t), intent(in) :: space
637 class(mesh_t), intent(in) :: mesh
638 type(ions_t), intent(in) :: ions
639 integer, intent(in) :: nd
640 type(local_domain_t), intent(inout) :: loc_domains(:)
641 type(restart_t), intent(in) :: restart
642
643 integer :: id, ip, iunit, ierr
644 character(len=31) :: line(1), tmp
645 real(real64), allocatable :: mask(:)
646
647 push_sub(local_restart_read)
648
649 message(1) = 'Info: Reading mesh points inside each local domain'
650 call messages_info(1)
651
652 safe_allocate(mask(1:mesh%np))
653
654 !Read local domain information from ldomains.info
655 iunit = restart_open(restart, "ldomains.info", status='old')
656 call restart_read(restart, iunit, line, 1, ierr)
657 read(line(1),'(a25,1x,i5)') tmp, ierr
658 call restart_close(restart, iunit)
659
660 call restart_read_mesh_function(restart, space, "ldomains", mesh, mask, ierr)
661
662 do id = 1, nd
663 loc_domains(id)%mesh_mask = .false.
664 do ip = 1 , mesh%np
665 if (bitand(int(mask(ip)), 2**id) /= 0) loc_domains(id)%mesh_mask(ip) = .true.
666 end do
667
668 !Check for atom list inside each domain
669 call local_ions_mask(loc_domains(id)%mesh_mask, ions, mesh, loc_domains(id)%ions_mask)
670 end do
671
672 safe_deallocate_a(mask)
673
674 pop_sub(local_restart_read)
675 end subroutine local_restart_read
676
677 ! ---------------------------------------------------------
678 subroutine local_inside_domain(space, mesh, ions, nd, loc_domains, namespace, ff)
679 class(space_t), intent(in) :: space
680 class(mesh_t), intent(in) :: mesh
681 type(ions_t), intent(in) :: ions
682 integer, intent(in) :: nd
683 type(local_domain_t), intent(inout) :: loc_domains(:)
684 type(namespace_t), intent(in) :: namespace
685 real(real64), intent(in) :: ff(:)
686
687 integer(int64) :: how
688 integer :: id, ip, ierr
689 type(basins_t) :: basins
690 character(len=MAX_PATH_LEN) :: filename
691 real(real64), allocatable :: dble_domain_map(:), domain_mesh(:)
692
693 push_sub(local_inside_domain)
694
695 !TODO: find a parallel algorithm to perform the charge-density fragmentation.
696 message(1) = 'Info: Assigning mesh points inside each local domain'
697 call messages_info(1)
698
699 ! If we have Bader domains, then we need to get the basins
700 if (any(loc_domains(:)%dshape == bader)) then
701 if (mesh%parallel_in_domains) then
702 write(message(1),'(a)') 'Bader volumes can only be computed in serial'
703 call messages_fatal(1)
704 end if
705 call create_basins(namespace, space, mesh, ions, ff, basins)
706 end if
707
708 do id = 1, nd
709 ! Create the mask that tells which mesh points are inside the domain
710 select case (loc_domains(id)%dshape)
711 case (bader)
712 call bader_domain_create_mask(loc_domains(id), basins, mesh, ions)
713 case default
714 call box_domain_create_mask(loc_domains(id), mesh)
715 end select
716
717 !Check for atom list inside each domain
718 call local_ions_mask(loc_domains(id)%mesh_mask, ions, mesh, loc_domains(id)%ions_mask)
719 end do
720
721 if (any(loc_domains(:)%dshape == bader)) then
722 call basins_end(basins)
723 end if
724
725 if (debug%info) then
726 call parse_variable(namespace, 'LDOutputFormat', 0, how)
727 if (.not. varinfo_valid_option('OutputFormat', how, is_flag=.true.)) then
728 call messages_input_error(namespace, 'LDOutputFormat')
729 end if
730
731 ! Write extra information about domains
732 safe_allocate(dble_domain_map(1:mesh%np))
733 safe_allocate(domain_mesh(1:mesh%np))
734
735 domain_mesh = m_zero
736 do id = 1, nd
737 dble_domain_map = m_zero
738 do ip = 1, mesh%np
739 if (loc_domains(id)%mesh_mask(ip)) then
740 dble_domain_map(ip) = real(id, real64)
741 domain_mesh(ip) = domain_mesh(ip) + dble_domain_map(ip)
742 end if
743 end do
744
745 write(filename,'(a,a)') 'domain.', trim(loc_domains(id)%lab)
746 call dio_function_output(how, 'local.general', trim(filename), namespace, space, mesh, dble_domain_map, unit_one, ierr, &
747 pos=ions%pos, atoms=ions%atom)
748 end do
749
750 call dio_function_output(how, 'local.general', 'domain.mesh', namespace, space, mesh, domain_mesh, unit_one, ierr, &
751 pos=ions%pos, atoms=ions%atom)
752
753 safe_deallocate_a(dble_domain_map)
754 safe_deallocate_a(domain_mesh)
755 end if
756
757 pop_sub(local_inside_domain)
758 end subroutine local_inside_domain
759
760 ! ---------------------------------------------------------
761 subroutine create_basins(namespace, space, mesh, ions, ff, basins)
762 type(namespace_t), intent(in) :: namespace
763 class(space_t), intent(in) :: space
764 class(mesh_t), intent(in) :: mesh
765 type(ions_t), intent(in) :: ions
766 real(real64), intent(in) :: ff(:)
767 type(basins_t), intent(inout) :: basins
768
769 logical :: use_atomic_radii
770 integer :: ip, ierr, ia, is, rankmin
771 integer(int64) :: how
772 real(real64) :: bader_threshold, dmin, dd
773 integer, allocatable :: ion_map(:)
774 real(real64), allocatable :: ff2(:,:), ffs(:)
775
776 push_sub(create_basins)
777
778 !%Variable LDBaderThreshold
779 !%Type float
780 !%Default 0.01
781 !%Section Utilities::oct-local_multipoles
782 !%Description
783 !% This variable sets the threshold for the basins calculations. Recommended values:
784 !% 0.01 -> intramolecular volumes; 0.2 -> intermolecular volumes.
785 !%End
786 call parse_variable(namespace, 'LDBaderThreshold', 0.01_real64, bader_threshold)
787
788 ! Copy function used to construct the basins, as we need to modify it
789 safe_allocate(ff2(1:mesh%np,1))
790 ff2(1:mesh%np,1) = ff(1:mesh%np)
791
792 ! Add long range density from atoms
793 safe_allocate(ffs(1:mesh%np))
794 do ia = 1, ions%natoms
795 call species_get_long_range_density(ions%atom(ia)%species, ions%namespace, ions%space, ions%latt, &
796 ions%pos(:, ia), mesh, ffs)
797 do is = 1, ubound(ff2, dim=2)
798 ff2(1:mesh%np, is) = ff2(1:mesh%np, is) - ffs(1:mesh%np)
799 end do
800 end do
801 safe_deallocate_a(ffs)
802
803 call basins_init(basins, namespace, mesh)
804 call basins_analyze(basins, namespace, mesh, ff2(:,1), ff2, bader_threshold)
805
806 if (debug%info) then
807 call parse_variable(namespace, 'LDOutputFormat', 0, how)
808 if (.not. varinfo_valid_option('OutputFormat', how, is_flag=.true.)) then
809 call messages_input_error(namespace, 'LDOutputFormat')
810 end if
811
812 call dio_function_output(how, 'local.general', 'basinsmap', namespace, space, mesh, real(basins%map(1:mesh%np), real64) , &
813 unit_one, ierr, pos=ions%pos, atoms=ions%atom)
814 call dio_function_output(how, 'local.general', 'dens_ff2', namespace, space, mesh, ff2(:,1), &
815 unit_one, ierr, pos=ions%pos, atoms=ions%atom)
816 end if
817 safe_deallocate_a(ff2)
818
819 !%Variable LDUseAtomicRadii
820 !%Type logical
821 !%Default false
822 !%Section Utilities::oct-local_multipoles
823 !%Description
824 !% If set, atomic radii will be used to assign lone pairs to ion.
825 !%End
826 call parse_variable(global_namespace, 'LDUseAtomicRadii', .false., use_atomic_radii)
827
828 if (use_atomic_radii) then
829 safe_allocate(ion_map(1:ions%natoms))
830 ion_map = 0
831 do ia = 1, ions%natoms
832 ion_map(ia) = basins%map(mesh_nearest_point(mesh, ions%pos(:, ia), dmin, rankmin))
833 end do
834
835 do ip = 1, mesh%np
836 ! Check if lonely pair has no ions assigned to it
837 if (all(ion_map(:) /= basins%map(ip))) then
838 ! Assign lonely pair to ion in a atomic radii distance
839 do ia = 1, ions%natoms
840 dd = norm2(mesh%x(ip, :) - ions%pos(:, ia))
841 if (dd <= ions%atom(ia)%species%get_vdw_radius()) basins%map(ip) = ion_map(ia)
842 end do
843 end if
844 end do
845
846 safe_deallocate_a(ion_map)
847 end if
848
849 pop_sub(create_basins)
850 end subroutine create_basins
851
852 ! ---------------------------------------------------------
853 subroutine box_domain_create_mask(domain, mesh)
854 class(local_domain_t), intent(inout) :: domain
855 class(mesh_t), intent(in) :: mesh
856
857 push_sub(box_domain_create_mask)
858
859 domain%mesh_mask = domain%box%contains_points(mesh%np, mesh%x)
860
862 end subroutine box_domain_create_mask
863
864 ! ---------------------------------------------------------
865 subroutine bader_domain_create_mask(domain, basins, mesh, ions)
866 class(local_domain_t), intent(inout) :: domain
867 type(basins_t), intent(in) :: basins
868 class(mesh_t), intent(in) :: mesh
869 type(ions_t), intent(in) :: ions
870
871 integer :: ia, ib, ip, ix, n_basins, rankmin
872 real(real64) :: dmin
873 integer, allocatable :: domain_map(:)
874
876
877 ! Count then number of basins to consider. That is the number of ions specified in the input by the user
878 n_basins = 0
879 do ia = 1, ions%natoms
880 if (loct_isinstringlist(ia, domain%clist)) n_basins = n_basins + 1
881 end do
882
883 safe_allocate(domain_map(1:n_basins))
884
885 ! Assign basins to ions
886 domain_map = 0
887 ib = 0
888 do ia = 1, ions%natoms
889 if (.not. loct_isinstringlist(ia, domain%clist)) cycle
890 ib = ib + 1
891 ix = mesh_nearest_point(mesh, ions%pos(:, ia), dmin, rankmin)
892 domain_map(ib) = basins%map(ix)
893 end do
894
895 ! Now construct the mask: a point belongs to this Bader domain if it is part
896 ! of at least one basin that is part of this domain
897 do ip = 1, mesh%np
898 domain%mesh_mask(ip) = any(domain_map(1:n_basins) == basins%map(ip))
899 end do
900
901 safe_deallocate_a(domain_map)
902
904 end subroutine bader_domain_create_mask
905
906 ! ---------------------------------------------------------
907 !Check for the ions inside each local domain.
908 subroutine local_ions_mask(mesh_mask, ions, mesh, ions_mask)
909 logical, intent(in) :: mesh_mask(:)
910 type(ions_t), intent(in) :: ions
911 class(mesh_t), intent(in) :: mesh
912 logical, intent(out) :: ions_mask(:)
913
914 integer :: ia, ix
915 integer, allocatable :: mask_tmp(:)
916 integer :: rankmin
917 real(real64) :: dmin
918
919 push_sub(local_ions_mask)
920
921 safe_allocate(mask_tmp(1:ions%natoms))
922 mask_tmp = 0
923
924 do ia = 1, ions%natoms
925 ix = mesh_nearest_point(mesh, ions%pos(:, ia), dmin, rankmin)
926 if (rankmin /= mesh%mpi_grp%rank) cycle
927 if (mesh_mask(ix)) mask_tmp(ia) = 1
928 end do
929
930 if (mesh%parallel_in_domains) then
931 call mesh%allreduce(mask_tmp, ions%natoms)
932 end if
933 ions_mask = mask_tmp == 1
934
935 safe_deallocate_a(mask_tmp)
936
937 pop_sub(local_ions_mask)
938 end subroutine local_ions_mask
939
940end program oct_local_multipoles
941
942!! Local Variables:
943!! mode: f90
944!! coding: utf-8
945!! 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:151
subroutine, public basins_analyze(this, namespace, mesh, f, rho, threshold)
Definition: basins.F90:188
subroutine, public basins_end(this)
Definition: basins.F90:170
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:154
real(real64), parameter, public m_two
Definition: global.F90:189
subroutine, public global_end()
Finalise parser varinfo file, and MPI.
Definition: global.F90:381
real(real64), parameter, public m_zero
Definition: global.F90:187
subroutine, public global_init(communicator)
Initialise Octopus.
Definition: global.F90:324
integer, parameter, public length
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:114
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:166
subroutine, public io_close(iunit, grp)
Definition: io.F90:468
subroutine, public io_end()
Definition: io.F90:265
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:395
subroutine, public kick_end(kick)
Definition: kick.F90:795
subroutine, public kick_init(kick, namespace, space, kpoints, nspin)
Definition: kick.F90:223
subroutine, public local_write_init(writ, namespace, lab, iter, dt)
subroutine, public local_write_end(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)
logical function, public loct_isinstringlist(a, s)
Definition: loct.F90:316
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
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:380
subroutine, public messages_end()
Definition: messages.F90:277
subroutine, public messages_print_with_emphasis(msg, iunit, namespace)
Definition: messages.F90:930
character(len=512), private msg
Definition: messages.F90:165
subroutine, public messages_init(output_dir)
Definition: messages.F90:224
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:543
subroutine, public print_date(str)
Definition: messages.F90:1017
subroutine, public messages_info(no_lines, iunit, verbose_limit, stress, all_nodes, namespace)
Definition: messages.F90:624
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:420
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:723
subroutine, public messages_experimental(name, namespace)
Definition: messages.F90:1097
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:266
This module handles the communicators for the various parallelization strategies.
Definition: multicomm.F90:145
type(namespace_t), public global_namespace
Definition: namespace.F90:132
subroutine, public parser_init()
Initialise the Octopus parser.
Definition: parser.F90:450
subroutine, public parser_end()
End the Octopus parser.
Definition: parser.F90:481
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
subroutine, public profiling_end(namespace)
Definition: profiling.F90:413
subroutine, public profiling_init(namespace)
Create profiling subdirectory.
Definition: profiling.F90:255
subroutine, public read_from_default_file(iunit, read_data, spec)
subroutine, public restart_module_init(namespace)
Definition: restart.F90:306
subroutine, public restart_read(restart, iunit, lines, nlines, ierr)
Definition: restart.F90:933
integer, parameter, public restart_undefined
Definition: restart.F90:229
subroutine, public restart_close(restart, iunit)
Close a file previously opened with restart_open.
Definition: restart.F90:950
subroutine, public restart_init(restart, namespace, data_type, type, mc, ierr, mesh, dir, exact)
Initializes a restart object.
Definition: restart.F90:514
integer, parameter, public restart_type_dump
Definition: restart.F90:225
subroutine, public restart_write(restart, iunit, lines, nlines, ierr)
Definition: restart.F90:906
integer function, public restart_open(restart, filename, status, position, silent)
Open file "filename" found inside the current restart directory. Depending on the type of restart,...
Definition: restart.F90:859
integer, parameter, public restart_type_load
Definition: restart.F90:225
subroutine, public restart_end(restart)
Definition: restart.F90:720
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:132
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:118
subroutine, public print_header()
This subroutine prints the logo followed by information about the compilation and the system....
Definition: utils.F90:301
Definition: xc.F90:114
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:141
Class implementing a parallelepiped box. Currently this is restricted to a rectangular cuboid (all th...
Class implementing a spherical box.
Definition: box_sphere.F90:132
Class implementing a box that is an union other boxes.
Definition: box_union.F90:133
Class describing the electron system.
Definition: electrons.F90:214
Describes mesh distribution to nodes.
Definition: mesh.F90:186
Stores all communicators and groups.
Definition: multicomm.F90:206
int true(void)