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