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
35 use unit_oct_m
37 use string_oct_m
38
39 implicit none
40
41 private
42 public :: &
44 index2axis, &
50 lead_dim, &
54
58 end interface leading_dimension_is_known
59
60 interface lead_dim
61 module procedure dlead_dim, zlead_dim
62 module procedure dlead_dim2, zlead_dim2
63 end interface lead_dim
64
65contains
66
67 ! ---------------------------------------------------------
68 subroutine get_divisors(nn, n_divisors, divisors)
69 integer, intent(in) :: nn
70 integer, intent(inout) :: n_divisors
71 integer, intent(out) :: divisors(:)
72
73 integer :: ii, max_d
74
75 push_sub(get_divisors)
76
77 assert(n_divisors > 1)
78 max_d = n_divisors
79
80 n_divisors = 1
81 divisors(n_divisors) = 1
82 do ii = 2, nn / 2
83 if (mod(nn, ii) == 0) then
84 n_divisors = n_divisors + 1
85
86 if (n_divisors > max_d - 1) then
87 message(1) = "Internal error in get_divisors. Please increase n_divisors"
88 call messages_fatal(1)
89 end if
90
91 divisors(n_divisors) = ii
92 end if
93 end do
94 n_divisors = n_divisors + 1
95 divisors(n_divisors) = nn
96
97 pop_sub(get_divisors)
98 end subroutine get_divisors
99
100
101 ! ---------------------------------------------------------
102 character pure function index2axis(idir) result(ch)
103 integer, intent(in) :: idir
104
105 select case (idir)
106 case (1)
107 ch = 'x'
108 case (2)
109 ch = 'y'
110 case (3)
111 ch = 'z'
112 case (4)
113 ch = 'w'
114 case default
115 write(ch,'(i1)') idir
116 end select
117
118 end function index2axis
119
120 ! ---------------------------------------------------------
121 pure function index2axisbz(idir) result(ch)
122 integer, intent(in) :: idir
123 character(len=2) :: ch
124
125 select case (idir)
126 case (1)
127 ch = "kx"
128 case (2)
129 ch = "ky"
130 case (3)
131 ch = "kz"
132 case (4)
133 ch = "kw"
134 case default
135 write(ch,'(i2)') idir
136 end select
137
138 end function index2axisbz
139
140
141 ! ---------------------------------------------------------
142 subroutine output_tensor(tensor, ndim, unit, write_average, iunit, namespace)
143 real(real64), intent(in) :: tensor(:,:)
144 integer, intent(in) :: ndim
145 type(unit_t), intent(in) :: unit
146 logical, optional, intent(in) :: write_average
147 integer, optional, intent(in) :: iunit
148 type(namespace_t), optional, intent(in) :: namespace
149
150 real(real64) :: trace
151 integer :: jj, kk
152 logical :: write_average_
154 push_sub(output_tensor)
155
156 write_average_ = optional_default(write_average, .true.)
157
158 trace = m_zero
159 message(1) = ""
160 do jj = 1, ndim
161 do kk = 1, ndim
162 write(message(1), '(a,f20.6)') trim(message(1)), units_from_atomic(unit, tensor(jj, kk))
163 end do
164 trace = trace + tensor(jj, jj)
165 call messages_info(1, iunit=iunit, namespace=namespace)
166 end do
167
168 if (write_average_) then
169 write(message(1), '(a, f20.6)') 'Isotropic average', units_from_atomic(unit, trace/real(ndim, real64) )
170 call messages_info(1, iunit=iunit, namespace=namespace)
171 end if
172
173 pop_sub(output_tensor)
174 end subroutine output_tensor
175
176
177 ! ---------------------------------------------------------
178 subroutine output_dipole(dipole, ndim, iunit, namespace)
179 real(real64), intent(in) :: dipole(:)
180 integer, intent(in) :: ndim
181 integer, optional, intent(in) :: iunit
182 type(namespace_t), optional, intent(in) :: namespace
183
184 integer :: idir
185
186 push_sub(output_dipole)
187
188 write(message(1), '(a,a20,a17)') 'Dipole:', '[' // trim(units_abbrev(units_out%length)) // ']', &
189 '[' // trim(units_abbrev(unit_debye)) // ']'
190 do idir = 1, ndim
191 write(message(1+idir), '(6x,3a,es14.5,3x,2es14.5)') '<', index2axis(idir), '> = ', &
192 units_from_atomic(units_out%length, dipole(idir)), units_from_atomic(unit_debye, dipole(idir))
193 end do
194 call messages_info(1+ndim, iunit=iunit, namespace=namespace)
196 pop_sub(output_dipole)
197 end subroutine output_dipole
198
201 subroutine print_header()
202 use iso_fortran_env
203
204 character(len=256) :: sys_name
205
206 ! Let us print our logo
207 call io_dump_file(stdout, trim(trim(conf%share) // '/logo'))
208
209 ! Let us print the version
210 message(1) = ""
211 message(2) = str_center("Running octopus", 70)
212 message(3) = ""
213 call messages_info(3)
215 message(1) = &
216 "Version : " // trim(conf%version)
217 message(2) = &
218 "Commit : "// trim(conf%git_commit)
219 message(3) = &
220 "Configuration time : "// trim(conf%config_time)
221 call messages_info(3)
222
223 message(1) = 'Configuration options : ' // trim(get_config_opts())
224 message(2) = 'Optional libraries :' // trim(get_optional_libraries())
225
226 message(3) = 'Architecture : ' + tostring(oct_arch)
227 call messages_info(3)
228
229 message(1) = &
230 "C compiler : "//trim(conf%cc)
231 message(2) = &
232 "C compiler flags : "//trim(conf%cflags)
233 message(3) = &
234 "C++ compiler : "//trim(conf%cxx)
235 message(4) = &
236 "C++ compiler flags : "//trim(conf%cxxflags)
237#ifdef HAVE_FC_COMPILER_VERSION
238 message(5) = "Fortran compiler : "//trim(conf%fc) //" ("//compiler_version()//")"
239#else
240 message(5) = "Fortran compiler : "//trim(conf%fc)
241#endif
242 message(6) = &
243 "Fortran compiler flags : "//trim(conf%fcflags)
244 call messages_info(6)
245
246 message(1) = ""
247 call messages_info(1)
248
249 ! Let us print where we are running
250 call loct_sysname(sys_name)
251 write(message(1), '(a)') str_center("The octopus is swimming in " // trim(sys_name), 70)
252 message(2) = ""
253 call messages_info(2)
254
255 call mpi_world%barrier()
256
257 call print_date("Calculation started on ")
258 end subroutine print_header
259
260 character(len=256) function get_config_opts()
261
262 get_config_opts = ''
263#ifdef HAVE_OPENMP
264 get_config_opts = trim(get_config_opts)//' openmp'
265#endif
266#ifdef HAVE_MPI
267 get_config_opts = trim(get_config_opts)//' mpi'
268#endif
269#ifdef HAVE_OPENCL
270 get_config_opts = trim(get_config_opts)//' opencl'
271#endif
272#ifdef HAVE_CUDA
273 get_config_opts = trim(get_config_opts)//' cuda'
274#endif
275#ifdef HAVE_M128D
276 get_config_opts = trim(get_config_opts)//' sse2'
277#endif
278#ifdef HAVE_M256D
279 get_config_opts = trim(get_config_opts)//' avx'
280#endif
281#ifdef HAVE_BLUE_GENE
282 get_config_opts = trim(get_config_opts)//' bluegene/p'
283#endif
284#ifdef HAVE_BLUE_GENE_Q
285 get_config_opts = trim(get_config_opts)//' bluegene/q'
286#endif
287#ifdef HAVE_LIBXC_FXC
288 get_config_opts = trim(get_config_opts)//' libxc_fxc'
289#endif
290#ifdef HAVE_LIBXC_KXC
291 get_config_opts = trim(get_config_opts)//' libxc_kxc'
292#endif
293
294 end function get_config_opts
295
296 character(len=256) function get_optional_libraries()
297
298 ! keep in alphabetical order, for ease in seeing if something is listed
300#ifdef HAVE_BERKELEYGW
302#endif
303#ifdef HAVE_CGAL
305#endif
306#ifdef HAVE_CLFFT
308#endif
309#if (defined(HAVE_CLBLAS)) || (defined(HAVE_CLBLAST))
311#endif
312#ifdef HAVE_DFTBPLUS
314#endif
315#ifdef HAVE_ELPA
317#endif
318#ifdef HAVE_ETSF_IO
320#endif
321#ifdef HAVE_GDLIB
323#endif
324#ifdef HAVE_LIBFM
326#endif
327#ifdef HAVE_PSOLVER
329#endif
330#ifdef HAVE_LIBVDWXC
332#endif
333#ifdef HAVE_METIS
335#endif
336#ifdef HAVE_NETCDF
338#endif
339#ifdef HAVE_NFFT
341#endif
342#ifdef HAVE_PARMETIS
344#endif
345#ifdef HAVE_PFFT
347#endif
348#ifdef HAVE_PNFFT
350#endif
351#ifdef HAVE_PSPIO
353#endif
354#ifdef HAVE_SCALAPACK
356#endif
357#ifdef HAVE_SPARSKIT
359#endif
360#ifdef HAVE_NLOPT
362#endif
363
364 end function get_optional_libraries
365
366 ! ---------------------------------------------------------
367
368 logical function dleading_dimension_is_known(array) result(known)
369 real(real64), intent(in) :: array(:, :)
370
371 known = .true.
372
373#if defined(HAVE_FORTRAN_LOC)
374 if (ubound(array, dim = 2) > 1) then
375 known = ubound(array, dim = 1) == (loc(array(1, 2)) - loc(array(1, 1)))/c_sizeof(array(1, 1))
376 end if
377#endif
378
379 end function dleading_dimension_is_known
380
381
382 ! ---------------------------------------------------------
383
384 logical function zleading_dimension_is_known(array) result(known)
385 complex(real64), intent(in) :: array(:, :)
386
387 known = .true.
388
389#if defined(HAVE_FORTRAN_LOC)
390 if (ubound(array, dim = 2) > 1) then
391 known = ubound(array, dim = 1) == (loc(array(1, 2)) - loc(array(1, 1)))/c_sizeof(array(1, 1))
392 end if
393#endif
394
395 end function zleading_dimension_is_known
396
397 ! ---------------------------------------------------------
398
399 logical function ileading_dimension_is_known(array) result(known)
400 integer, intent(in) :: array(:, :)
401
402 known = .true.
403
404#if defined(HAVE_FORTRAN_LOC) && defined(HAVE_FC_SIZEOF)
405 if (ubound(array, dim = 2) > 1) then
406 known = ubound(array, dim = 1) == (loc(array(1, 2)) - loc(array(1, 1)))/sizeof(array(1, 1))
407 end if
408#endif
409
410 end function ileading_dimension_is_known
411
412
413 ! ---------------------------------------------------------
414
415 logical function dleading_dimension_is_known2(array) result(known)
416 real(real64), intent(in) :: array(:, :, :)
417
418 known = .true.
419
420#if defined(HAVE_FORTRAN_LOC)
421 if (ubound(array, dim = 2) > 1) then
422 known = ubound(array, dim = 1) == (loc(array(1, 2, 1)) - loc(array(1, 1, 1)))/c_sizeof(array(1, 1, 1))
423 end if
424#endif
425
427
428
429 ! ---------------------------------------------------------
430
431 logical function zleading_dimension_is_known2(array) result(known)
432 complex(real64), intent(in) :: array(:, :, :)
433
434 known = .true.
435
436#if defined(HAVE_FORTRAN_LOC)
437 if (ubound(array, dim = 2) > 1) then
438 known = ubound(array, dim = 1) == (loc(array(1, 2, 1)) - loc(array(1, 1, 1)))/c_sizeof(array(1, 1, 1))
439 end if
440#endif
441
443
444 ! ---------------------------------------------------------
445
446 logical function ileading_dimension_is_known2(array) result(known)
447 integer, intent(in) :: array(:, :, :)
449 known = .true.
450
451#if defined(HAVE_FORTRAN_LOC) && defined(HAVE_FC_SIZEOF)
452 if (ubound(array, dim = 2) > 1) then
453 known = ubound(array, dim = 1) == (loc(array(1, 2, 1)) - loc(array(1, 1, 1)))/sizeof(array(1, 1, 1))
454 end if
455#endif
456
458
459
460 ! ---------------------------------------------------------
461
462 integer function dlead_dim(array) result(lead_dim)
463 real(real64), intent(in) :: array(:, :)
465 assert(leading_dimension_is_known(array))
466
467 lead_dim = ubound(array, dim = 1)
468 end function dlead_dim
469
470 ! ---------------------------------------------------------
471
472 integer function zlead_dim(array) result(lead_dim)
473 complex(real64), intent(in) :: array(:, :)
474
475 assert(leading_dimension_is_known(array))
476
477 lead_dim = ubound(array, dim = 1)
478 end function zlead_dim
480 ! ---------------------------------------------------------
481
482 integer function dlead_dim2(array) result(lead_dim)
483 real(real64), intent(in) :: array(:, :, :)
484
485 assert(leading_dimension_is_known(array))
486
487 lead_dim = ubound(array, dim = 1) * ubound(array, dim = 2)
488 end function dlead_dim2
489
490 ! ---------------------------------------------------------
491
492 integer function zlead_dim2(array) result(lead_dim)
493 complex(real64), intent(in) :: array(:, :, :)
494
496
497 lead_dim = ubound(array, dim = 1) * ubound(array, dim = 2)
498 end function zlead_dim2
499
500 subroutine make_array_larger(array, new_size)
501 integer(int64), allocatable, intent(inout) :: array(:)
502 integer, intent(in) :: new_size
503
504 integer(int64), allocatable :: tmp(:)
505 integer :: copy_size
506
507 push_sub(make_array_larger)
508
509 safe_allocate(tmp(1:new_size))
510 copy_size = min(new_size, size(array))
511 tmp(1:copy_size) = array(1:copy_size)
512 safe_deallocate_a(array)
513 call move_alloc(tmp, array)
514
515 pop_sub(make_array_larger)
516 end subroutine make_array_larger
517
518end module utils_oct_m
519
520!! Local Variables:
521!! mode: f90
522!! coding: utf-8
523!! End:
!in the one should use SAFE_DEALLOCATE_P for pointers and SAFE_DEALLOCATE_A !for arrays !A special version of the SAFE_ALLOCATE macro named SAFE_ALLOCATE_TYPE is also !provided to allocate a polymorphic variable This is necessary because of the !special Fortran syntax type::var !Some versions of GCC have a bug in the sizeof() function such that the compiler crashes with a ICE ! when passing a polymorphic variable to the function and explicit array bounds are given. ! The workaround is not to pass the bounds to sizeof. Otherwise we could just use SAFE_ALLOCATE_TYPE. ! the TOSTRING macro converts a macro into a string ! do not use the STRINGIFY macro ! Whenever a procedure is not called too many times
Definition: io.F90:114
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
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.
This module is intended to contain simple general-purpose utility functions and procedures.
Definition: utils.F90:118
subroutine, public output_tensor(tensor, ndim, unit, write_average, iunit, namespace)
Definition: utils.F90:236
character(len=256) function, public get_config_opts()
Definition: utils.F90:354
integer function zlead_dim2(array)
Definition: utils.F90:573
subroutine, public output_dipole(dipole, ndim, iunit, namespace)
Definition: utils.F90:272
logical function zleading_dimension_is_known2(array)
Definition: utils.F90:512
logical function zleading_dimension_is_known(array)
Definition: utils.F90:465
subroutine, public get_divisors(nn, n_divisors, divisors)
Definition: utils.F90:162
integer function dlead_dim(array)
Definition: utils.F90:543
subroutine, public make_array_larger(array, new_size)
Definition: utils.F90:581
logical function dleading_dimension_is_known2(array)
Definition: utils.F90:496
character(len=256) function, public get_optional_libraries()
Definition: utils.F90:377
logical function ileading_dimension_is_known(array)
Definition: utils.F90:480
integer function dlead_dim2(array)
Definition: utils.F90:563
character pure function, public index2axis(idir)
Definition: utils.F90:196
logical function ileading_dimension_is_known2(array)
Definition: utils.F90:527
pure character(len=2) function, public index2axisbz(idir)
Definition: utils.F90:215
subroutine, public print_header()
This subroutine prints the logo followed by information about the compilation and the system....
Definition: utils.F90:295
integer function zlead_dim(array)
Definition: utils.F90:553
logical function dleading_dimension_is_known(array)
Definition: utils.F90:449
int true(void)