69 integer,
intent(in) :: nn
70 integer,
intent(inout) :: n_divisors
71 integer,
intent(out) :: divisors(:)
77 assert(n_divisors > 1)
81 divisors(n_divisors) = 1
83 if (mod(nn, ii) == 0)
then
84 n_divisors = n_divisors + 1
86 if (n_divisors > max_d - 1)
then
87 message(1) =
"Internal error in get_divisors. Please increase n_divisors"
91 divisors(n_divisors) = ii
94 n_divisors = n_divisors + 1
95 divisors(n_divisors) = nn
102 character pure function index2axis(idir) result(ch)
103 integer,
intent(in) :: idir
115 write(ch,
'(i1)') idir
122 integer,
intent(in) :: idir
123 character(len=2) :: ch
135 write(ch,
'(i2)') idir
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
150 real(real64) :: trace
152 logical :: write_average_
156 write_average_ = optional_default(write_average, .
true.)
162 write(message(1),
'(a,f20.6)') trim(message(1)), units_from_atomic(unit, tensor(jj, kk))
164 trace = trace + tensor(jj, jj)
165 call messages_info(1, iunit=iunit, namespace=namespace)
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)
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
188 write(message(1),
'(a,a20,a17)')
'Dipole:',
'[' // trim(units_abbrev(units_out%length)) //
']', &
189 '[' // trim(units_abbrev(unit_debye)) //
']'
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))
194 call messages_info(1+ndim, iunit=iunit, namespace=namespace)
204 character(len=256) :: sys_name
207 call io_dump_file(stdout, trim(trim(conf%share) //
'/logo'))
211 message(2) = str_center(
"Running octopus", 70)
213 call messages_info(3)
216 "Version : " // trim(conf%version)
218 "Commit : "// trim(conf%git_commit)
220 "Configuration time : "// trim(conf%config_time)
221 call messages_info(3)
226 message(3) =
'Architecture : ' + tostring(oct_arch)
227 call messages_info(3)
230 "C compiler : "//trim(conf%cc)
232 "C compiler flags : "//trim(conf%cflags)
234 "C++ compiler : "//trim(conf%cxx)
236 "C++ compiler flags : "//trim(conf%cxxflags)
237#ifdef HAVE_FC_COMPILER_VERSION
238 message(5) =
"Fortran compiler : "//trim(conf%fc) //
" ("//compiler_version()//
")"
240 message(5) =
"Fortran compiler : "//trim(conf%fc)
243 "Fortran compiler flags : "//trim(conf%fcflags)
244 call messages_info(6)
247 call messages_info(1)
250 call loct_sysname(sys_name)
251 write(message(1),
'(a)') str_center(
"The octopus is swimming in " // trim(sys_name), 70)
253 call messages_info(2)
255 call mpi_world%barrier()
257 call print_date(
"Calculation started on ")
284#ifdef HAVE_BLUE_GENE_Q
300#ifdef HAVE_BERKELEYGW
309#if (defined(HAVE_CLBLAS)) || (defined(HAVE_CLBLAST))
369 real(real64),
intent(in) :: array(:, :)
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))
385 complex(real64),
intent(in) :: array(:, :)
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))
400 integer,
intent(in) :: array(:, :)
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))
416 real(real64),
intent(in) :: array(:, :, :)
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))
432 complex(real64),
intent(in) :: array(:, :, :)
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))
447 integer,
intent(in) :: array(:, :, :)
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))
462 integer function dlead_dim(array)
result(lead_dim)
463 real(real64),
intent(in) :: array(:, :)
472 integer function zlead_dim(array)
result(lead_dim)
473 complex(real64),
intent(in) :: array(:, :)
482 integer function dlead_dim2(array)
result(lead_dim)
483 real(real64),
intent(in) :: array(:, :, :)
487 lead_dim = ubound(array, dim = 1) * ubound(array, dim = 2)
492 integer function zlead_dim2(array)
result(lead_dim)
493 complex(real64),
intent(in) :: array(:, :, :)
497 lead_dim = ubound(array, dim = 1) * ubound(array, dim = 2)
501 integer(int64),
allocatable,
intent(inout) :: array(:)
502 integer,
intent(in) :: new_size
504 integer(int64),
allocatable :: tmp(:)
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)
!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
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
This module defines the unit system, used for input and output.
This module is intended to contain simple general-purpose utility functions and procedures.
subroutine, public output_tensor(tensor, ndim, unit, write_average, iunit, namespace)
character(len=256) function, public get_config_opts()
integer function zlead_dim2(array)
subroutine, public output_dipole(dipole, ndim, iunit, namespace)
logical function zleading_dimension_is_known2(array)
logical function zleading_dimension_is_known(array)
subroutine, public get_divisors(nn, n_divisors, divisors)
integer function dlead_dim(array)
subroutine, public make_array_larger(array, new_size)
logical function dleading_dimension_is_known2(array)
character(len=256) function, public get_optional_libraries()
logical function ileading_dimension_is_known(array)
integer function dlead_dim2(array)
character pure function, public index2axis(idir)
logical function ileading_dimension_is_known2(array)
pure character(len=2) function, public index2axisbz(idir)
subroutine, public print_header()
This subroutine prints the logo followed by information about the compilation and the system....
integer function zlead_dim(array)
logical function dleading_dimension_is_known(array)