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