Octopus
utils.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!! Copyright (C) 2021 S. Ohlmann
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
24
25module utils_oct_m
26 use debug_oct_m
27 use global_oct_m
28 use io_oct_m
29 use iso_c_binding
30 use loct_oct_m
32 use mpi_oct_m
36 use unit_oct_m
38 use string_oct_m
39
40 implicit none
41
42 private
43 public :: &
45 index2axis, &
51 lead_dim, &
55
57 module procedure dleading_dimension_is_known, &
65 end interface leading_dimension_is_known
66
67 interface lead_dim
68 module procedure dlead_dim, zlead_dim
69 module procedure dlead_dim2, zlead_dim2
70 end interface lead_dim
71
72contains
73
74 ! ---------------------------------------------------------
75 subroutine get_divisors(nn, n_divisors, divisors)
76 integer, intent(in) :: nn
77 integer, intent(inout) :: n_divisors
78 integer, intent(out) :: divisors(:)
79
80 integer :: ii, max_d
81
82 push_sub(get_divisors)
83
84 assert(n_divisors > 1)
85 max_d = n_divisors
86
87 n_divisors = 1
88 divisors(n_divisors) = 1
89 do ii = 2, nn / 2
90 if (mod(nn, ii) == 0) then
91 n_divisors = n_divisors + 1
92
93 if (n_divisors > max_d - 1) then
94 message(1) = "Internal error in get_divisors. Please increase n_divisors"
95 call messages_fatal(1)
96 end if
97
98 divisors(n_divisors) = ii
99 end if
100 end do
101 n_divisors = n_divisors + 1
102 divisors(n_divisors) = nn
103
104 pop_sub(get_divisors)
105 end subroutine get_divisors
106
107
108 ! ---------------------------------------------------------
109 character pure function index2axis(idir) result(ch)
110 integer, intent(in) :: idir
111
112 select case (idir)
113 case (1)
114 ch = 'x'
115 case (2)
116 ch = 'y'
117 case (3)
118 ch = 'z'
119 case (4)
120 ch = 'w'
121 case default
122 write(ch,'(i1)') idir
123 end select
124
125 end function index2axis
126
127 ! ---------------------------------------------------------
128 pure function index2axisbz(idir) result(ch)
129 integer, intent(in) :: idir
130 character(len=2) :: ch
131
132 select case (idir)
133 case (1)
134 ch = "kx"
135 case (2)
136 ch = "ky"
137 case (3)
138 ch = "kz"
139 case (4)
140 ch = "kw"
141 case default
142 write(ch,'(i2)') idir
143 end select
144
145 end function index2axisbz
146
147
148 ! ---------------------------------------------------------
149 subroutine output_tensor(tensor, ndim, unit, write_average, iunit, namespace)
150 real(real64), intent(in) :: tensor(:,:)
151 integer, intent(in) :: ndim
152 type(unit_t), intent(in) :: unit
153 logical, optional, intent(in) :: write_average
154 integer, optional, intent(in) :: iunit
155 type(namespace_t), optional, intent(in) :: namespace
156
157 real(real64) :: trace
158 integer :: jj, kk
159 logical :: write_average_
160
161 push_sub(output_tensor)
163 write_average_ = optional_default(write_average, .true.)
164
165 trace = m_zero
166 message(1) = ""
167 do jj = 1, ndim
168 do kk = 1, ndim
169 write(message(1), '(a,f20.6)') trim(message(1)), units_from_atomic(unit, tensor(jj, kk))
170 end do
171 trace = trace + tensor(jj, jj)
172 call messages_info(1, iunit=iunit, namespace=namespace)
173 end do
174
175 if (write_average_) then
176 write(message(1), '(a, f20.6)') 'Isotropic average', units_from_atomic(unit, trace/real(ndim, real64) )
177 call messages_info(1, iunit=iunit, namespace=namespace)
178 end if
179
180 pop_sub(output_tensor)
181 end subroutine output_tensor
182
183
184 ! ---------------------------------------------------------
185 subroutine output_dipole(dipole, ndim, iunit, namespace)
186 real(real64), intent(in) :: dipole(:)
187 integer, intent(in) :: ndim
188 integer, optional, intent(in) :: iunit
189 type(namespace_t), optional, intent(in) :: namespace
190
191 integer :: idir
192
193 push_sub(output_dipole)
194
195 write(message(1), '(a,a20,a17)') 'Dipole:', '[' // trim(units_abbrev(units_out%length)) // ']', &
196 '[' // trim(units_abbrev(unit_debye)) // ']'
197 do idir = 1, ndim
198 write(message(1+idir), '(6x,3a,es14.5,3x,2es14.5)') '<', index2axis(idir), '> = ', &
199 units_from_atomic(units_out%length, dipole(idir)), units_from_atomic(unit_debye, dipole(idir))
200 end do
201 call messages_info(1+ndim, iunit=iunit, namespace=namespace)
202
203 pop_sub(output_dipole)
204 end subroutine output_dipole
205
208 subroutine print_header()
209 use iso_fortran_env
210
211 character(len=256) :: sys_name
212
213 ! Let us print our logo
214 call io_dump_file(stdout, trim(trim(conf%share) // '/logo'))
215
216 ! Let us print the version
217 message(1) = ""
218 message(2) = str_center("Running octopus", 70)
219 message(3) = ""
220 call messages_info(3)
221
222 message(1) = &
223 "Version : " // trim(conf%version)
224 message(2) = &
225 "Commit : "// trim(conf%git_commit)
226 message(3) = &
227 "Configuration time : "// trim(conf%config_time)
228 call messages_info(3)
229
230 message(1) = 'Configuration options : ' // trim(get_config_opts())
231 message(2) = 'Optional libraries :' // trim(get_optional_libraries())
232
233 message(3) = 'Architecture : ' + tostring(oct_arch)
234 call messages_info(3)
235
237
238 message(1) = &
239 "C compiler : "//trim(conf%cc)
240 message(2) = &
241 "C compiler flags : "//trim(conf%cflags)
242 message(3) = &
243 "C++ compiler : "//trim(conf%cxx)
244 message(4) = &
245 "C++ compiler flags : "//trim(conf%cxxflags)
246#ifdef HAVE_FC_COMPILER_VERSION
247 message(5) = "Fortran compiler : "//trim(conf%fc) //" ("//compiler_version()//")"
248#else
249 message(5) = "Fortran compiler : "//trim(conf%fc)
250#endif
251 message(6) = &
252 "Fortran compiler flags : "//trim(conf%fcflags)
253 call messages_info(6)
254
255 message(1) = ""
256 call messages_info(1)
257
258 ! Let us print where we are running
259 call loct_sysname(sys_name)
260 write(message(1), '(a)') str_center("The octopus is swimming in " // trim(sys_name), 70)
261 message(2) = ""
262 call messages_info(2)
263
264 call mpi_world%barrier()
265
266 call print_date("Calculation started on ")
267 end subroutine print_header
268
269 character(len=256) function get_config_opts()
270
271 get_config_opts = ''
272#ifdef HAVE_OPENMP
273 get_config_opts = trim(get_config_opts)//' openmp'
274#endif
275#ifdef HAVE_MPI
276 get_config_opts = trim(get_config_opts)//' mpi'
277#endif
278#ifdef HAVE_CUDA
279 get_config_opts = trim(get_config_opts)//' cuda'
280#endif
281#ifdef HAVE_LIBXC_FXC
282 get_config_opts = trim(get_config_opts)//' libxc_fxc'
283#endif
284#ifdef HAVE_LIBXC_KXC
285 get_config_opts = trim(get_config_opts)//' libxc_kxc'
286#endif
287
288 end function get_config_opts
289
290 character(len=256) function get_optional_libraries()
291
292 ! keep in alphabetical order, for ease in seeing if something is listed
294#ifdef HAVE_ADIOS2
296#endif
297#ifdef HAVE_BERKELEYGW
299#endif
300#ifdef HAVE_CGAL
302#endif
303#ifdef HAVE_DFTBPLUS
305#endif
306#ifdef HAVE_ELPA
308#endif
309#ifdef HAVE_ETSF_IO
311#endif
312#ifdef HAVE_GDLIB
314#endif
315#ifdef HAVE_PSOLVER
317#endif
318#ifdef HAVE_LIBVDWXC
320#endif
321#ifdef HAVE_METIS
323#endif
324#ifdef HAVE_NETCDF
326#endif
327#ifdef HAVE_NFFT
329#endif
330#ifdef HAVE_PARMETIS
332#endif
333#ifdef HAVE_PFFT
335#endif
336#ifdef HAVE_PNFFT
338#endif
339#ifdef HAVE_PSPIO
341#endif
342#ifdef HAVE_SCALAPACK
344#endif
345#ifdef HAVE_SPARSKIT
347#endif
348#ifdef HAVE_NLOPT
350#endif
351
352 end function get_optional_libraries
353
354 ! ---------------------------------------------------------
355
356 logical function dleading_dimension_is_known(array) result(known)
357 real(real64), target, intent(in) :: array(:, :)
358
359 integer(c_intptr_t) :: addr1, addr2
360
361 known = .true.
362
363 if (ubound(array, dim = 2) > 1) then
364 addr1 = transfer(c_loc(array(1,1)), 0_c_intptr_t)
365 addr2 = transfer(c_loc(array(1,2)), 0_c_intptr_t)
366 known = ubound(array, dim = 1) == (addr2-addr1)/c_sizeof(array(1, 1))
367 end if
368
369 end function dleading_dimension_is_known
371
372 ! ---------------------------------------------------------
373
374 logical function zleading_dimension_is_known(array) result(known)
375 complex(real64), target, intent(in) :: array(:, :)
376
377 integer(c_intptr_t) :: addr1, addr2
378
379 known = .true.
380
381 if (ubound(array, dim = 2) > 1) then
382 addr1 = transfer(c_loc(array(1,1)), 0_c_intptr_t)
383 addr2 = transfer(c_loc(array(1,2)), 0_c_intptr_t)
384 known = ubound(array, dim = 1) == (addr2-addr1)/c_sizeof(array(1, 1))
385 end if
386
387 end function zleading_dimension_is_known
388
389 ! ---------------------------------------------------------
390
391 logical function ileading_dimension_is_known(array) result(known)
392 integer, target, intent(in) :: array(:, :)
393
394 integer(c_intptr_t) :: addr1, addr2
395
396 known = .true.
398 if (ubound(array, dim = 2) > 1) then
399 addr1 = transfer(c_loc(array(1,1)), 0_c_intptr_t)
400 addr2 = transfer(c_loc(array(1,2)), 0_c_intptr_t)
401
402 known = ubound(array, dim = 1) == (addr2-addr1)/c_sizeof(array(1, 1))
403 end if
404
405 end function ileading_dimension_is_known
406
407
408 logical function lleading_dimension_is_known(array) result(known)
409 integer(int64), target, intent(in) :: array(:, :)
410
411 integer(c_intptr_t) :: addr1, addr2
412
413 known = .true.
415 if (ubound(array, dim = 2) > 1) then
416 addr1 = transfer(c_loc(array(1,1)), 0_c_intptr_t)
417 addr2 = transfer(c_loc(array(1,2)), 0_c_intptr_t)
418
419 known = ubound(array, dim = 1) == (addr2-addr1)/c_sizeof(array(1, 1))
420 end if
421
422 end function lleading_dimension_is_known
423
424
425 ! ---------------------------------------------------------
426
427 logical function dleading_dimension_is_known2(array) result(known)
428 real(real64), target, intent(in) :: array(:, :, :)
429
430 integer(c_intptr_t) :: addr1, addr2
432 known = .true.
433
434 if (ubound(array, dim = 2) > 1) then
435 addr1 = transfer(c_loc(array(1,1,1)), 0_c_intptr_t)
436 addr2 = transfer(c_loc(array(1,2,1)), 0_c_intptr_t)
437 known = ubound(array, dim = 1) == (addr2 - addr1)/c_sizeof(array(1, 1, 1))
438 end if
439
441
442
443 ! ---------------------------------------------------------
444
445 logical function zleading_dimension_is_known2(array) result(known)
446 complex(real64), target, intent(in) :: array(:, :, :)
447
448 integer(c_intptr_t) :: addr1, addr2
449
450 known = .true.
451
452 if (ubound(array, dim = 2) > 1) then
453 addr1 = transfer(c_loc(array(1,1,1)), 0_c_intptr_t)
454 addr2 = transfer(c_loc(array(1,2,1)), 0_c_intptr_t)
455 known = ubound(array, dim = 1) == (addr2 - addr1)/c_sizeof(array(1, 1, 1))
456 end if
457
459
460 ! ---------------------------------------------------------
461
462 logical function ileading_dimension_is_known2(array) result(known)
463 integer, target, intent(in) :: array(:, :, :)
464
465 integer(c_intptr_t) :: addr1, addr2
466
467 known = .true.
469 if (ubound(array, dim = 2) > 1) then
470 addr1 = transfer(c_loc(array(1,1,1)), 0_c_intptr_t)
471 addr2 = transfer(c_loc(array(1,2,1)), 0_c_intptr_t)
472 known = ubound(array, dim = 1) == (addr2 - addr1)/c_sizeof(array(1, 1, 1))
473 end if
474
476
477
478 logical function lleading_dimension_is_known2(array) result(known)
479 integer(int64), target, intent(in) :: array(:, :, :)
480
481 integer(c_intptr_t) :: addr1, addr2
482
483 known = .true.
484
485 if (ubound(array, dim = 2) > 1) then
486 addr1 = transfer(c_loc(array(1,1,1)), 0_c_intptr_t)
487 addr2 = transfer(c_loc(array(1,2,1)), 0_c_intptr_t)
488 known = ubound(array, dim = 1) == (addr2 - addr1)/c_sizeof(array(1, 1, 1))
489 end if
490
492
493 ! ---------------------------------------------------------
494
495 integer function dlead_dim(array) result(lead_dim)
496 real(real64), intent(in) :: array(:, :)
497
498 assert(leading_dimension_is_known(array))
499
500 lead_dim = ubound(array, dim = 1)
501 end function dlead_dim
502
503 ! ---------------------------------------------------------
504
505 integer function zlead_dim(array) result(lead_dim)
506 complex(real64), intent(in) :: array(:, :)
507
508 assert(leading_dimension_is_known(array))
509
510 lead_dim = ubound(array, dim = 1)
511 end function zlead_dim
512
513 ! ---------------------------------------------------------
514
515 integer function dlead_dim2(array) result(lead_dim)
516 real(real64), intent(in) :: array(:, :, :)
517
519
520 lead_dim = ubound(array, dim = 1) * ubound(array, dim = 2)
521 end function dlead_dim2
522
523 ! ---------------------------------------------------------
524
525 integer function zlead_dim2(array) result(lead_dim)
526 complex(real64), intent(in) :: array(:, :, :)
527
529
530 lead_dim = ubound(array, dim = 1) * ubound(array, dim = 2)
531 end function zlead_dim2
532
533 subroutine make_array_larger(array, new_size)
534 integer(int64), allocatable, intent(inout) :: array(:)
535 integer, intent(in) :: new_size
536
537 integer(int64), allocatable :: tmp(:)
538 integer :: copy_size
539
540 push_sub(make_array_larger)
541
542 ! Use allocate here as move_alloc deallocate internally the from array
543 allocate(tmp(1:new_size))
544 copy_size = min(new_size, size(array))
545 tmp(1:copy_size) = array(1:copy_size)
546 safe_deallocate_a(array)
547 call move_alloc(tmp, array)
549 pop_sub(make_array_larger)
550 end subroutine make_array_larger
551
553 subroutine write_vectorization_level()
554 character(len=32) :: vec
555 character(kind=c_char) :: c_str(33)
557 c_str = c_null_char
558 call get_vectorization_level(c_str)
559 call string_c_to_f(c_str, vec)
560 message(1) = 'Vectorization level : ' // trim(vec)
561 call messages_info(1)
562
563 end subroutine write_vectorization_level
564
565end module utils_oct_m
566
567!! Local Variables:
568!! mode: f90
569!! coding: utf-8
570!! End:
Definition: io.F90:116
System information (time, memory, sysname)
Definition: loct.F90:117
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:410
This module contains interfaces for routines in operate.c.
Definition: operate_f.F90:119
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.
This module is intended to contain simple general-purpose utility functions and procedures.
Definition: utils.F90:120
subroutine, public output_tensor(tensor, ndim, unit, write_average, iunit, namespace)
Definition: utils.F90:245
character(len=256) function, public get_config_opts()
Definition: utils.F90:365
integer function zlead_dim2(array)
Definition: utils.F90:549
subroutine write_vectorization_level()
Prints the level of vectorization used for the vectorized finite differences.
Definition: utils.F90:577
subroutine, public output_dipole(dipole, ndim, iunit, namespace)
Definition: utils.F90:281
logical function zleading_dimension_is_known2(array)
Definition: utils.F90:469
logical function lleading_dimension_is_known(array)
Definition: utils.F90:432
logical function zleading_dimension_is_known(array)
Definition: utils.F90:398
subroutine, public get_divisors(nn, n_divisors, divisors)
Definition: utils.F90:171
integer function dlead_dim(array)
Definition: utils.F90:519
subroutine, public make_array_larger(array, new_size)
Definition: utils.F90:557
logical function dleading_dimension_is_known2(array)
Definition: utils.F90:451
character(len=256) function, public get_optional_libraries()
Definition: utils.F90:371
logical function ileading_dimension_is_known(array)
Definition: utils.F90:415
integer function dlead_dim2(array)
Definition: utils.F90:539
character pure function, public index2axis(idir)
Definition: utils.F90:205
logical function lleading_dimension_is_known2(array)
Definition: utils.F90:502
logical function ileading_dimension_is_known2(array)
Definition: utils.F90:486
pure character(len=2) function, public index2axisbz(idir)
Definition: utils.F90:224
subroutine, public print_header()
This subroutine prints the logo followed by information about the compilation and the system....
Definition: utils.F90:304
integer function zlead_dim(array)
Definition: utils.F90:529
logical function dleading_dimension_is_known(array)
Definition: utils.F90:380
void get_vectorization_level(char *level)
Definition: operate.c:733
int true(void)