LCOV - code coverage report
Current view: top level - src - environment.F (source / functions) Coverage Total Hit
Test: CP2K Regtests (git:cccd2f3) Lines: 87.7 % 568 498
Test Date: 2026-05-06 07:07:47 Functions: 100.0 % 14 14

            Line data    Source code
       1              : !--------------------------------------------------------------------------------------------------!
       2              : !   CP2K: A general program to perform molecular dynamics simulations                              !
       3              : !   Copyright 2000-2026 CP2K developers group <https://cp2k.org>                                   !
       4              : !                                                                                                  !
       5              : !   SPDX-License-Identifier: GPL-2.0-or-later                                                      !
       6              : !--------------------------------------------------------------------------------------------------!
       7              : 
       8              : ! **************************************************************************************************
       9              : !> \brief Sets up and terminates the global environment variables
      10              : !> \par History
      11              : !>      - Merged with Quickstep MODULE start_program_run (17.01.2002,MK)
      12              : !>      - Compile information added (16.01.2002,MK)
      13              : !>      - Merged with MODULE cp2k_input, some rearrangements (30.10.2002,MK)
      14              : !>      - Update seed input (24.10.2016,MK)
      15              : !> \author JGH,MK
      16              : ! **************************************************************************************************
      17              : MODULE environment
      18              :    USE bibliography,                    ONLY: Frigo2005,&
      19              :                                               Marek2014,&
      20              :                                               Solca2024,&
      21              :                                               cite_reference
      22              :    USE cp2k_info,                       ONLY: &
      23              :         compile_arch, compile_date, compile_host, compile_revision, cp2k_flags, cp2k_home, &
      24              :         cp2k_version, cp2k_year, get_runtime_info, r_host_name, r_pid, r_user_name
      25              :    USE cp_error_handling,               ONLY: warning_counter
      26              :    USE cp_files,                        ONLY: close_file,&
      27              :                                               get_data_dir,&
      28              :                                               open_file
      29              :    USE cp_fm_cholesky,                  ONLY: FM_CHOLESKY_TYPE_DLAF,&
      30              :                                               FM_CHOLESKY_TYPE_SCALAPACK,&
      31              :                                               cholesky_type,&
      32              :                                               dlaf_cholesky_n_min
      33              :    USE cp_fm_diag,                      ONLY: FM_DIAG_TYPE_CUSOLVER,&
      34              :                                               FM_DIAG_TYPE_DLAF,&
      35              :                                               FM_DIAG_TYPE_ELPA,&
      36              :                                               FM_DIAG_TYPE_SCALAPACK,&
      37              :                                               cusolver_n_min,&
      38              :                                               diag_finalize,&
      39              :                                               diag_init,&
      40              :                                               eps_check_diag_default
      41              :    USE cp_fm_diag_utils,                ONLY: cp_fm_redistribute_init
      42              :    USE cp_fm_struct,                    ONLY: cp_fm_struct_config
      43              :    USE cp_fm_types,                     ONLY: cp_fm_get_mm_type,&
      44              :                                               cp_fm_setup
      45              :    USE cp_log_handling,                 ONLY: &
      46              :         cp_add_default_logger, cp_get_default_logger, cp_logger_create, &
      47              :         cp_logger_get_default_unit_nr, cp_logger_release, cp_logger_set, cp_logger_type, &
      48              :         cp_rm_default_logger, cp_to_string
      49              :    USE cp_output_handling,              ONLY: cp_mpi_io_set,&
      50              :                                               cp_print_key_finished_output,&
      51              :                                               cp_print_key_unit_nr,&
      52              :                                               debug_print_level,&
      53              :                                               high_print_level,&
      54              :                                               low_print_level,&
      55              :                                               medium_print_level,&
      56              :                                               silent_print_level
      57              :    USE fft_tools,                       ONLY: FWFFT,&
      58              :                                               fft3d,&
      59              :                                               finalize_fft,&
      60              :                                               init_fft
      61              :    USE force_env_types,                 ONLY: multiple_fe_list
      62              :    USE gamma,                           ONLY: deallocate_md_ftable
      63              :    USE global_types,                    ONLY: global_environment_type
      64              :    USE grid_api,                        ONLY: GRID_BACKEND_AUTO,&
      65              :                                               GRID_BACKEND_CPU,&
      66              :                                               GRID_BACKEND_DGEMM,&
      67              :                                               GRID_BACKEND_GPU,&
      68              :                                               GRID_BACKEND_HIP,&
      69              :                                               GRID_BACKEND_REF
      70              :    USE header,                          ONLY: cp2k_footer,&
      71              :                                               cp2k_header
      72              :    USE input_constants,                 ONLY: &
      73              :         callgraph_all, callgraph_none, do_cosma, do_cp2k, do_dgemm_blas, do_dgemm_spla, do_eip, &
      74              :         do_farming, do_fft_fftw3, do_fft_sg, do_fist, do_qs, do_scalapack, do_sirius, do_test, &
      75              :         energy_run, mol_dyn_run, none_run
      76              :    USE input_cp2k_global,               ONLY: create_global_section
      77              :    USE input_enumeration_types,         ONLY: enum_i2c,&
      78              :                                               enumeration_type
      79              :    USE input_keyword_types,             ONLY: keyword_get,&
      80              :                                               keyword_type
      81              :    USE input_section_types,             ONLY: &
      82              :         section_get_ival, section_get_keyword, section_get_lval, section_get_rval, &
      83              :         section_release, section_type, section_vals_get, section_vals_get_subs_vals, &
      84              :         section_vals_get_subs_vals3, section_vals_type, section_vals_val_get, section_vals_val_set
      85              :    USE kinds,                           ONLY: default_path_length,&
      86              :                                               default_string_length,&
      87              :                                               dp,&
      88              :                                               int_8,&
      89              :                                               print_kind_info
      90              :    USE local_gemm_api,                  ONLY: local_gemm_set_library
      91              :    USE machine,                         ONLY: &
      92              :         flush_should_flush, m_cpuid, m_cpuid_name, m_cpuid_static, m_cpuid_vlen, m_cpuinfo, &
      93              :         m_energy, m_memory_details, m_omp_get_stacksize, m_omp_trace_issues, m_procrun
      94              :    USE message_passing,                 ONLY: mp_collect_timings,&
      95              :                                               mp_para_env_type
      96              :    USE mp_perf_env,                     ONLY: add_mp_perf_env,&
      97              :                                               describe_mp_perf_env,&
      98              :                                               rm_mp_perf_env
      99              :    USE orbital_pointers,                ONLY: deallocate_orbital_pointers,&
     100              :                                               init_orbital_pointers
     101              :    USE orbital_transformation_matrices, ONLY: deallocate_spherical_harmonics,&
     102              :                                               init_spherical_harmonics
     103              :    USE parallel_rng_types,              ONLY: GAUSSIAN,&
     104              :                                               check_rng,&
     105              :                                               rng_stream_type,&
     106              :                                               write_rng_matrices
     107              :    USE physcon,                         ONLY: write_physcon
     108              :    USE reference_manager,               ONLY: collect_citations_from_ranks,&
     109              :                                               print_cited_references
     110              :    USE string_utilities,                ONLY: ascii_to_string,&
     111              :                                               integer_to_string,&
     112              :                                               string_to_ascii
     113              :    USE timings,                         ONLY: add_timer_env,&
     114              :                                               global_timings_level,&
     115              :                                               rm_timer_env,&
     116              :                                               root_cp2k_name,&
     117              :                                               timings_setup_tracing
     118              :    USE timings_report,                  ONLY: cost_type_energy,&
     119              :                                               cost_type_time,&
     120              :                                               timings_report_callgraph,&
     121              :                                               timings_report_print
     122              :    USE voronoi_interface,               ONLY: finalize_libvori
     123              : 
     124              : !$ USE OMP_LIB, ONLY: omp_get_max_threads, omp_get_thread_num, omp_get_num_threads
     125              : #include "./base/base_uses.f90"
     126              : 
     127              :    IMPLICIT NONE
     128              : 
     129              :    PRIVATE
     130              : 
     131              :    CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'environment'
     132              : 
     133              :    ! Public subroutines
     134              : 
     135              :    PUBLIC :: cp2k_finalize, cp2k_init, cp2k_read, cp2k_setup, cp2k_get_walltime
     136              : 
     137              : CONTAINS
     138              : 
     139              : ! **************************************************************************************************
     140              : !> \brief Initializes a CP2K run (setting of the global environment variables)
     141              : !> \param para_env ...
     142              : !> \param output_unit ...
     143              : !> \param globenv ...
     144              : !> \param input_file_name ...
     145              : !> \param wdir ...
     146              : !> \par History
     147              : !>      JGH (28.11.2001) : default for pp_library_path
     148              : !>      - print keys added (17.01.2002, MK)
     149              : !>      - merged with cp2k_input (30.10.2002,MK)
     150              : !> \author JGH,MK
     151              : ! **************************************************************************************************
     152        10420 :    SUBROUTINE cp2k_init(para_env, output_unit, globenv, input_file_name, wdir)
     153              : 
     154              :       TYPE(mp_para_env_type), POINTER                    :: para_env
     155              :       INTEGER                                            :: output_unit
     156              :       TYPE(global_environment_type), POINTER             :: globenv
     157              :       CHARACTER(LEN=*)                                   :: input_file_name
     158              :       CHARACTER(LEN=*), OPTIONAL                         :: wdir
     159              : 
     160              :       CHARACTER(LEN=10*default_string_length)            :: cp_flags
     161              :       INTEGER                                            :: i, ilen, my_output_unit
     162              :       TYPE(cp_logger_type), POINTER                      :: logger
     163              : 
     164              :       ! create a timer_env
     165              : 
     166        10420 :       CALL add_timer_env()
     167              : 
     168              :       ! Message passing performance
     169        10420 :       CALL add_mp_perf_env()
     170              : 
     171              :       ! Init the default logger
     172        10420 :       IF (para_env%is_source()) THEN
     173         5314 :          my_output_unit = output_unit
     174              :       ELSE
     175         5106 :          my_output_unit = -1
     176              :       END IF
     177        10420 :       NULLIFY (logger)
     178              :       CALL cp_logger_create(logger, para_env=para_env, &
     179              :                             default_global_unit_nr=output_unit, &
     180        10420 :                             close_global_unit_on_dealloc=.FALSE.)
     181        10420 :       CALL cp_add_default_logger(logger)
     182        10420 :       CALL cp_logger_release(logger)
     183              : 
     184              :       ! Initialize timing
     185        10420 :       CALL timeset(root_cp2k_name, globenv%handle)
     186              : 
     187              :       ! Print header
     188        10944 :       CALL cp2k_header(my_output_unit, wdir)
     189              : 
     190        10420 :       IF (my_output_unit > 0) THEN
     191              :          WRITE (UNIT=my_output_unit, FMT="(/,T2,A,T31,A50)") &
     192         5314 :             "CP2K| version string: ", ADJUSTR(TRIM(cp2k_version))
     193              :          WRITE (UNIT=my_output_unit, FMT="(T2,A,T41,A40)") &
     194         5314 :             "CP2K| source code revision number:", &
     195        10628 :             ADJUSTR(compile_revision)
     196         5314 :          cp_flags = cp2k_flags()
     197         5314 :          ilen = LEN_TRIM(cp_flags)
     198              :          WRITE (UNIT=my_output_unit, FMT="(T2,A)") &
     199         5314 :             "CP2K| "//cp_flags(1:73)
     200         5314 :          IF (ilen > 73) THEN
     201        21256 :             DO i = 0, (ilen - 75)/61
     202              :                WRITE (UNIT=my_output_unit, FMT="(T2,A)") &
     203        21256 :                   "CP2K|            "//TRIM(cp_flags(74 + i*61:MIN(74 + (i + 1)*61, ilen)))
     204              :             END DO
     205              :          END IF
     206              :          WRITE (UNIT=my_output_unit, FMT="(T2,A,T41,A40)") &
     207         5314 :             "CP2K| is freely available from ", &
     208        10628 :             ADJUSTR(TRIM(cp2k_home))
     209              :          WRITE (UNIT=my_output_unit, FMT="(T2,A,T31,A50)") &
     210         5314 :             "CP2K| Program compiled at", &
     211        10628 :             ADJUSTR(compile_date(1:MIN(50, LEN(compile_date))))
     212              :          WRITE (UNIT=my_output_unit, FMT="(T2,A,T31,A50)") &
     213         5314 :             "CP2K| Program compiled on", &
     214        10628 :             ADJUSTR(compile_host(1:MIN(50, LEN(compile_host))))
     215              :          WRITE (UNIT=my_output_unit, FMT="(T2,A,T31,A50)") &
     216         5314 :             "CP2K| Program compiled for", &
     217        10628 :             ADJUSTR(compile_arch(1:MIN(50, LEN(compile_arch))))
     218              :          WRITE (UNIT=my_output_unit, FMT="(T2,A,T31,A50)") &
     219         5314 :             "CP2K| Data directory path", &
     220        10628 :             ADJUSTR(TRIM(get_data_dir()))
     221              :          WRITE (UNIT=my_output_unit, FMT="(T2,A,T31,A50)") &
     222         5314 :             "CP2K| Input file name", &
     223        10628 :             ADJUSTR(TRIM(input_file_name))
     224         5314 :          FLUSH (my_output_unit) ! ignore &GLOBAL / FLUSH_SHOULD_FLUSH
     225              :       END IF
     226              : 
     227              : #if defined(__FAST_MATH__)
     228              :       CALL cp_warn(__LOCATION__, &
     229              :                    "During compilation one of the following flags was active:"// &
     230              :                    "   `-ffast-math` (GCC)"// &
     231              :                    "   `-hfpN` (Cray, N > 0, default N=2)"// &
     232              :                    " This can lead to wrong results and numerical instabilities"// &
     233              :                    " and is therefore no longer supported.")
     234              : 
     235              : #if !defined(__FORCE_USE_FAST_MATH)
     236              : #error "-ffast-math (GCC) or -hfpN (N>0, Cray) can lead to wrong results and numerical instabilities and are therefore no longer supported"
     237              : #endif
     238              : #endif
     239              : 
     240              : #if defined(NDEBUG)
     241              : #error "Please do not build CP2K with NDEBUG. There is no performance advantage and asserts will save your neck."
     242              : #endif
     243              : 
     244        10420 :    END SUBROUTINE cp2k_init
     245              : 
     246              : ! **************************************************************************************************
     247              : !> \brief echoes the list of host names and pids
     248              : !> \param para_env ...
     249              : !> \param output_unit ...
     250              : ! **************************************************************************************************
     251            2 :    SUBROUTINE echo_all_hosts(para_env, output_unit)
     252              :       TYPE(mp_para_env_type), POINTER                    :: para_env
     253              :       INTEGER, INTENT(IN)                                :: output_unit
     254              : 
     255              :       CHARACTER(LEN=default_string_length)               :: string
     256              :       INTEGER                                            :: ipe
     257            2 :       INTEGER, ALLOCATABLE, DIMENSION(:)                 :: all_pid
     258              :       INTEGER, ALLOCATABLE, DIMENSION(:, :)              :: all_host
     259              : 
     260              :       ! Print a list of all started processes
     261              : 
     262           10 :       ALLOCATE (all_pid(para_env%num_pe), SOURCE=0)
     263            2 :       all_pid(para_env%mepos + 1) = r_pid
     264              : 
     265            2 :       CALL para_env%sum(all_pid)
     266          130 :       ALLOCATE (all_host(30, para_env%num_pe), SOURCE=0)
     267              : 
     268            2 :       CALL string_to_ascii(r_host_name, all_host(:, para_env%mepos + 1))
     269            2 :       CALL para_env%sum(all_host)
     270            2 :       IF (output_unit > 0) THEN
     271            1 :          WRITE (UNIT=output_unit, FMT="(T2,A)") ""
     272            3 :          DO ipe = 1, para_env%num_pe
     273            2 :             CALL ascii_to_string(all_host(:, ipe), string)
     274              :             WRITE (UNIT=output_unit, FMT="(T2,A,T63,I8,T71,I10)") &
     275              :                TRIM(r_user_name)//"@"//TRIM(string)// &
     276            3 :                " has created rank and process ", ipe - 1, all_pid(ipe)
     277              :          END DO
     278              :       END IF
     279              : 
     280            2 :       DEALLOCATE (all_pid)
     281            2 :       DEALLOCATE (all_host)
     282              : 
     283            2 :    END SUBROUTINE echo_all_hosts
     284              : 
     285              : ! **************************************************************************************************
     286              : !> \brief echoes the list the number of process per host
     287              : !> \param para_env ...
     288              : !> \param output_unit ...
     289              : !> \param node_count Count number of distributed systems (nodes)
     290              : ! **************************************************************************************************
     291        10420 :    SUBROUTINE echo_all_process_host(para_env, output_unit, node_count)
     292              :       TYPE(mp_para_env_type), POINTER                    :: para_env
     293              :       INTEGER, INTENT(IN)                                :: output_unit
     294              :       INTEGER, INTENT(OUT), OPTIONAL                     :: node_count
     295              : 
     296              :       CHARACTER(LEN=default_string_length)               :: string, string_sec
     297              :       INTEGER                                            :: ipe, jpe, nr_dist, nr_occu
     298        10420 :       INTEGER, ALLOCATABLE, DIMENSION(:)                 :: all_pid
     299              :       INTEGER, ALLOCATABLE, DIMENSION(:, :)              :: all_host
     300              : 
     301       670883 :       ALLOCATE (all_host(30, para_env%num_pe), SOURCE=0)
     302        51893 :       ALLOCATE (all_pid(para_env%num_pe), SOURCE=0)
     303              : 
     304        10420 :       IF (m_procrun(r_pid) == 1) THEN
     305        10420 :          CALL string_to_ascii(r_host_name, all_host(:, para_env%mepos + 1))
     306        10420 :          CALL para_env%sum(all_host)
     307              :       END IF
     308              : 
     309        10420 :       nr_dist = 0
     310        10420 :       IF (output_unit > 0) WRITE (UNIT=output_unit, FMT="(T2,A)") ""
     311        31053 :       DO ipe = 1, para_env%num_pe
     312        20633 :          nr_occu = 0
     313        31053 :          IF (all_pid(ipe) /= -1) THEN
     314        10420 :             CALL ascii_to_string(all_host(:, ipe), string)
     315        31053 :             DO jpe = 1, para_env%num_pe
     316        20633 :                CALL ascii_to_string(all_host(:, jpe), string_sec)
     317        31053 :                IF (string == string_sec) THEN
     318        20633 :                   nr_occu = nr_occu + 1
     319        20633 :                   all_pid(jpe) = -1
     320              :                END IF
     321              :             END DO
     322        10420 :             IF (output_unit > 0) THEN
     323              :                WRITE (UNIT=output_unit, FMT="(T2,A,T63,I8,A)") &
     324              :                   TRIM(r_user_name)//"@"//TRIM(string)// &
     325            1 :                   " is running ", nr_occu, " processes"
     326              :             END IF
     327        10420 :             nr_dist = nr_dist + 1
     328              :          END IF
     329              :       END DO
     330              : 
     331        10420 :       DEALLOCATE (all_pid)
     332        10420 :       DEALLOCATE (all_host)
     333              : 
     334        10420 :       CPASSERT(0 < nr_dist)
     335        10420 :       IF (PRESENT(node_count)) node_count = nr_dist
     336              : 
     337        10420 :    END SUBROUTINE echo_all_process_host
     338              : 
     339              : ! **************************************************************************************************
     340              : !> \brief read part of cp2k_init
     341              : !> \param root_section ...
     342              : !> \param para_env ...
     343              : !> \param globenv the globenv
     344              : !> \author fawzi
     345              : !> \note
     346              : !>      The following routines need to be synchronized wrt. adding/removing
     347              : !>      of the default environments (logging, performance,error):
     348              : !>      environment:cp2k_init, environment:cp2k_finalize,
     349              : !>      f77_interface:f_env_add_defaults, f77_interface:f_env_rm_defaults,
     350              : !>      f77_interface:create_force_env, f77_interface:destroy_force_env
     351              : ! **************************************************************************************************
     352        10420 :    SUBROUTINE cp2k_read(root_section, para_env, globenv)
     353              : 
     354              :       TYPE(section_vals_type), POINTER                   :: root_section
     355              :       TYPE(mp_para_env_type), POINTER                    :: para_env
     356              :       TYPE(global_environment_type), POINTER             :: globenv
     357              : 
     358              :       CHARACTER(LEN=3*default_string_length)             :: message
     359              :       CHARACTER(LEN=default_string_length)               :: c_val
     360              :       INTEGER                                            :: i, iw
     361              :       TYPE(cp_logger_type), POINTER                      :: logger
     362              : 
     363              :       ! Read the input/output section
     364              : 
     365        10420 :       logger => cp_get_default_logger()
     366              : 
     367              :       ! try to use better names for the local log if it is not too late
     368              :       CALL section_vals_val_get(root_section, "GLOBAL%OUTPUT_FILE_NAME", &
     369        10420 :                                 c_val=c_val)
     370        10420 :       IF (c_val /= "") THEN
     371              :          CALL cp_logger_set(logger, &
     372          144 :                             local_filename=TRIM(c_val)//"_localLog")
     373              :       END IF
     374              : 
     375              :       ! Process project name
     376        10420 :       CALL section_vals_val_get(root_section, "GLOBAL%PROJECT", c_val=c_val)
     377        10420 :       IF (INDEX(c_val(:LEN_TRIM(c_val)), " ") > 0) THEN
     378              :          message = "Project name <"//TRIM(c_val)// &
     379            2 :                    "> contains spaces which will be replaced with underscores"
     380            2 :          CPWARN(TRIM(message))
     381           24 :          DO i = 1, LEN_TRIM(c_val)
     382              :             ! Replace space with underscore
     383           24 :             IF (c_val(i:i) == " ") c_val(i:i) = "_"
     384              :          END DO
     385            2 :          CALL section_vals_val_set(root_section, "GLOBAL%PROJECT", c_val=TRIM(c_val))
     386              :       END IF
     387        10420 :       IF (c_val /= "") THEN
     388        10420 :          CALL cp_logger_set(logger, local_filename=TRIM(c_val)//"_localLog")
     389              :       END IF
     390        10420 :       logger%iter_info%project_name = c_val
     391              : 
     392        10420 :       CALL section_vals_val_get(root_section, "GLOBAL%PRINT_LEVEL", i_val=logger%iter_info%print_level)
     393              : 
     394              :       ! Read the CP2K section
     395        10420 :       CALL read_cp2k_section(root_section, para_env, globenv)
     396              : 
     397              :       iw = cp_print_key_unit_nr(logger, root_section, "GLOBAL%PRINT/BASIC_DATA_TYPES", &
     398        10420 :                                 extension=".Log")
     399        10420 :       IF (iw > 0) CALL print_kind_info(iw)
     400              :       CALL cp_print_key_finished_output(iw, logger, root_section, &
     401        10420 :                                         "GLOBAL%PRINT/BASIC_DATA_TYPES")
     402              : 
     403              :       iw = cp_print_key_unit_nr(logger, root_section, "GLOBAL%PRINT/PHYSCON", &
     404        10420 :                                 extension=".Log")
     405        10420 :       IF (iw > 0) CALL write_physcon(iw)
     406              :       CALL cp_print_key_finished_output(iw, logger, root_section, &
     407        10420 :                                         "GLOBAL%PRINT/PHYSCON")
     408              : 
     409        10420 :    END SUBROUTINE cp2k_read
     410              : 
     411              : ! **************************************************************************************************
     412              : !> \brief globenv initializations that need the input and error
     413              : !> \param root_section ...
     414              : !> \param para_env ...
     415              : !> \param globenv the global environment to initialize
     416              : !> \author fawzi
     417              : !> \note
     418              : !>      if possible do the initializations here as the environment
     419              : !>      (error,...) is setup, instead of cp2k_init
     420              : ! **************************************************************************************************
     421        20840 :    SUBROUTINE cp2k_setup(root_section, para_env, globenv)
     422              : 
     423              :       TYPE(section_vals_type), POINTER                   :: root_section
     424              :       TYPE(mp_para_env_type), POINTER                    :: para_env
     425              :       TYPE(global_environment_type), POINTER             :: globenv
     426              : 
     427              :       INTEGER                                            :: iw, maxl
     428        10420 :       INTEGER, DIMENSION(:), POINTER                     :: seed_vals
     429              :       REAL(KIND=dp), DIMENSION(3, 2)                     :: initial_seed
     430              :       TYPE(cp_logger_type), POINTER                      :: logger
     431              : 
     432        10420 :       NULLIFY (logger)
     433        20840 :       logger => cp_get_default_logger()
     434              : 
     435              :       ! Initialize the parallel random number generator
     436              : 
     437              :       iw = cp_print_key_unit_nr(logger, root_section, "GLOBAL%PRINT/RNG_MATRICES", &
     438        10420 :                                 extension=".Log")
     439        10420 :       IF (iw > 0) THEN
     440            1 :          CALL write_rng_matrices(iw)
     441              :       END IF
     442              : 
     443              :       CALL cp_print_key_finished_output(iw, logger, root_section, &
     444        10420 :                                         "GLOBAL%PRINT/RNG_MATRICES")
     445              : 
     446              :       ! Initialize a global normally Gaussian distributed (pseudo)random number stream
     447              : 
     448        10420 :       CALL section_vals_val_get(root_section, "GLOBAL%SEED", i_vals=seed_vals)
     449        10420 :       IF (SIZE(seed_vals) == 1) THEN
     450        93780 :          initial_seed(:, :) = REAL(seed_vals(1), KIND=dp)
     451            0 :       ELSE IF (SIZE(seed_vals) == 6) THEN
     452            0 :          initial_seed(1:3, 1:2) = RESHAPE(REAL(seed_vals(:), KIND=dp), [3, 2])
     453              :       ELSE
     454            0 :          CPABORT("Supply exactly 1 or 6 arguments for SEED in &GLOBAL only!")
     455              :       END IF
     456              : 
     457              :       globenv%gaussian_rng_stream = rng_stream_type( &
     458              :                                     name="Global Gaussian random numbers", &
     459              :                                     distribution_type=GAUSSIAN, &
     460              :                                     seed=initial_seed, &
     461        10420 :                                     extended_precision=.TRUE.)
     462              : 
     463              :       iw = cp_print_key_unit_nr(logger, root_section, "GLOBAL%PRINT/RNG_CHECK", &
     464        10420 :                                 extension=".Log")
     465        10420 :       IF (iw > 0) THEN
     466            1 :          CALL check_rng(iw, para_env%is_source())
     467              :       END IF
     468              : 
     469              :       CALL cp_print_key_finished_output(iw, logger, root_section, &
     470        10420 :                                         "GLOBAL%PRINT/RNG_CHECK")
     471              : 
     472              :       iw = cp_print_key_unit_nr(logger, root_section, "GLOBAL%PRINT/GLOBAL_GAUSSIAN_RNG", &
     473        10420 :                                 extension=".Log")
     474        10420 :       IF (iw > 0) &
     475            1 :          CALL globenv%gaussian_rng_stream%write(iw, write_all=.TRUE.)
     476              : 
     477              :       CALL cp_print_key_finished_output(iw, logger, root_section, &
     478        10420 :                                         "GLOBAL%PRINT/GLOBAL_GAUSSIAN_RNG")
     479              : 
     480        10420 :       CALL section_vals_val_get(root_section, "GLOBAL%PRINT%SPHERICAL_HARMONICS", i_val=maxl)
     481        10420 :       IF (maxl >= 0) THEN
     482              :          iw = cp_print_key_unit_nr(logger, root_section, "GLOBAL%PRINT", &
     483            2 :                                    extension=".Log")
     484            2 :          CALL init_orbital_pointers(maxl)
     485            2 :          CALL init_spherical_harmonics(maxl, iw)
     486            2 :          CALL deallocate_spherical_harmonics()
     487            2 :          CALL deallocate_orbital_pointers()
     488              :          CALL cp_print_key_finished_output(iw, logger, root_section, &
     489            2 :                                            "GLOBAL%PRINT")
     490              :       END IF
     491              : 
     492        10420 :    END SUBROUTINE cp2k_setup
     493              : 
     494              : ! **************************************************************************************************
     495              : !> \brief read the global section of new input
     496              : !> \param root_section ...
     497              : !> \param para_env ...
     498              : !> \param globenv ...
     499              : !> \par History
     500              : !>      06-2005 [created]
     501              : !> \author MI
     502              : !> \note
     503              : !>      Should not be required anymore once everything is converted
     504              : !>      to get information directly from the input structure
     505              : ! **************************************************************************************************
     506       156300 :    SUBROUTINE read_global_section(root_section, para_env, globenv)
     507              : 
     508              :       TYPE(section_vals_type), POINTER                   :: root_section
     509              :       TYPE(mp_para_env_type), POINTER                    :: para_env
     510              :       TYPE(global_environment_type), POINTER             :: globenv
     511              : 
     512              :       CHARACTER(LEN=6), PARAMETER                        :: start_section_label = "GLOBAL"
     513              : 
     514              :       CHARACTER(LEN=13)                                  :: omp_stacksize, tracing_string
     515              :       CHARACTER(LEN=6)                                   :: print_level_string
     516              :       CHARACTER(LEN=default_path_length)                 :: basis_set_file_name, coord_file_name, &
     517              :                                                             mm_potential_file_name, &
     518              :                                                             potential_file_name
     519              :       CHARACTER(LEN=default_string_length)               :: env_num, model_name, project_name
     520              :       CHARACTER(LEN=default_string_length), &
     521        10420 :          DIMENSION(:), POINTER                           :: trace_routines
     522              :       INTEGER :: cpuid, cpuid_static, i_cholesky, i_dgemm, i_diag, i_fft, i_grid_backend, &
     523              :          iforce_eval, method_name_id, n_rep_val, nforce_eval, node_count, num_threads, &
     524              :          output_unit, print_level, trace_max, unit_nr
     525              :       INTEGER(kind=int_8) :: Buffers, Buffers_avr, Buffers_max, Buffers_min, Cached, Cached_avr, &
     526              :          Cached_max, Cached_min, MemFree, MemFree_avr, MemFree_max, MemFree_min, MemLikelyFree, &
     527              :          MemLikelyFree_avr, MemLikelyFree_max, MemLikelyFree_min, MemTotal, MemTotal_avr, &
     528              :          MemTotal_max, MemTotal_min, Slab, Slab_avr, Slab_max, Slab_min, SReclaimable, &
     529              :          SReclaimable_avr, SReclaimable_max, SReclaimable_min
     530        10420 :       INTEGER, DIMENSION(:), POINTER                     :: i_force_eval
     531              :       LOGICAL                                            :: ata, do_echo_all_hosts, efl, explicit, &
     532              :                                                             flag, report_maxloc, trace, &
     533              :                                                             trace_master
     534              :       TYPE(cp_logger_type), POINTER                      :: logger
     535              :       TYPE(enumeration_type), POINTER                    :: enum1, enum2
     536              :       TYPE(keyword_type), POINTER                        :: keyword
     537              :       TYPE(section_type), POINTER                        :: section
     538              :       TYPE(section_vals_type), POINTER                   :: dft_section, force_env_sections, &
     539              :                                                             global_section, qmmm_section, &
     540              :                                                             subsys_section
     541              : 
     542        10420 :       NULLIFY (dft_section, global_section, i_force_eval)
     543              : 
     544        20840 :       logger => cp_get_default_logger()
     545        10420 :       global_section => section_vals_get_subs_vals(root_section, "GLOBAL")
     546        10420 :       CALL section_vals_val_get(global_section, "BLACS_GRID", i_val=globenv%blacs_grid_layout)
     547        10420 :       CALL section_vals_val_get(global_section, "BLACS_REPEATABLE", l_val=globenv%blacs_repeatable)
     548        10420 :       CALL section_vals_val_get(global_section, "PREFERRED_DIAG_LIBRARY", i_val=i_diag)
     549              :       CALL section_vals_val_get(global_section, "DIRECT_GENERALIZED_DIAGONALIZATION", &
     550        10420 :                                 l_val=globenv%direct_generalized_diagonalization)
     551        10420 :       CALL section_vals_val_get(global_section, "PREFERRED_CHOLESKY_LIBRARY", i_val=i_cholesky)
     552        10420 :       CALL section_vals_val_get(global_section, "PREFERRED_DGEMM_LIBRARY", i_val=i_dgemm)
     553        10420 :       CALL section_vals_val_get(global_section, "EPS_CHECK_DIAG", r_val=globenv%eps_check_diag)
     554        10420 :       CALL section_vals_val_get(global_section, "ENABLE_MPI_IO", l_val=flag)
     555        10420 :       CALL cp_mpi_io_set(flag)
     556        10420 :       CALL section_vals_val_get(global_section, "ELPA_KERNEL", i_val=globenv%k_elpa)
     557        10420 :       CALL section_vals_val_get(global_section, "ELPA_NEIGVEC_MIN", i_val=globenv%elpa_neigvec_min)
     558        10420 :       CALL section_vals_val_get(global_section, "ELPA_QR", l_val=globenv%elpa_qr)
     559        10420 :       CALL section_vals_val_get(global_section, "ELPA_ONE_STAGE", l_val=globenv%elpa_one_stage)
     560        10420 :       CALL section_vals_val_get(global_section, "ELPA_PRINT", l_val=globenv%elpa_print)
     561        10420 :       CALL section_vals_val_get(global_section, "DLAF_NEIGVEC_MIN", i_val=globenv%dlaf_neigvec_min)
     562        10420 :       CALL section_vals_val_get(global_section, "DLAF_CHOLESKY_N_MIN", i_val=globenv%dlaf_cholesky_n_min)
     563        10420 :       CALL section_vals_val_get(global_section, "PREFERRED_FFT_LIBRARY", i_val=i_fft)
     564        10420 :       CALL section_vals_val_get(global_section, "PRINT_LEVEL", i_val=print_level)
     565        10420 :       CALL section_vals_val_get(global_section, "PROGRAM_NAME", i_val=globenv%prog_name_id)
     566        10420 :       CALL section_vals_val_get(global_section, "FFT_POOL_SCRATCH_LIMIT", i_val=globenv%fft_pool_scratch_limit)
     567        10420 :       CALL section_vals_val_get(global_section, "FFTW_PLAN_TYPE", i_val=globenv%fftw_plan_type)
     568        10420 :       CALL section_vals_val_get(global_section, "PROJECT_NAME", c_val=project_name)
     569        10420 :       CALL section_vals_val_get(global_section, "FFTW_WISDOM_FILE_NAME", c_val=globenv%fftw_wisdom_file_name)
     570        10420 :       CALL section_vals_val_get(global_section, "RUN_TYPE", i_val=globenv%run_type_id)
     571              :       CALL cp2k_get_walltime(section=global_section, keyword_name="WALLTIME", &
     572        10420 :                              walltime=globenv%cp2k_target_time)
     573        10420 :       CALL section_vals_val_get(global_section, "TRACE", l_val=trace)
     574        10420 :       CALL section_vals_val_get(global_section, "TRACE_MASTER", l_val=trace_MASTER)
     575        10420 :       CALL section_vals_val_get(global_section, "TRACE_MAX", i_val=trace_max)
     576        10420 :       CALL section_vals_val_get(global_section, "TRACE_ROUTINES", explicit=explicit)
     577        10420 :       IF (explicit) THEN
     578            0 :          CALL section_vals_val_get(global_section, "TRACE_ROUTINES", c_vals=trace_routines)
     579              :       ELSE
     580        10420 :          NULLIFY (trace_routines)
     581              :       END IF
     582        10420 :       CALL section_vals_val_get(global_section, "FLUSH_SHOULD_FLUSH", l_val=flush_should_flush)
     583        10420 :       CALL section_vals_val_get(global_section, "ECHO_ALL_HOSTS", l_val=do_echo_all_hosts)
     584        10420 :       report_maxloc = section_get_lval(global_section, "TIMINGS%REPORT_MAXLOC")
     585        10420 :       global_timings_level = section_get_ival(global_section, "TIMINGS%TIMINGS_LEVEL")
     586        10420 :       do_echo_all_hosts = do_echo_all_hosts .OR. report_maxloc
     587        10420 :       force_env_sections => section_vals_get_subs_vals(root_section, "FORCE_EVAL")
     588        10420 :       CALL section_vals_get(force_env_sections, n_repetition=nforce_eval)
     589              :       output_unit = cp_print_key_unit_nr(logger, global_section, "PROGRAM_RUN_INFO", &
     590        10420 :                                          extension=".log")
     591              : 
     592        10420 :       CALL fm_setup(global_section)
     593        10420 :       CALL fm_diag_rules_setup(global_section)
     594        10420 :       CALL dgemm_setup(global_section)
     595              : 
     596        10420 :       IF (trace .AND. (.NOT. trace_master .OR. para_env%mepos == 0)) THEN
     597            0 :          unit_nr = -1
     598            0 :          IF (logger%para_env%is_source() .OR. .NOT. trace_master) &
     599            0 :             unit_nr = cp_logger_get_default_unit_nr(logger, local=.TRUE.)
     600            0 :          WRITE (tracing_string, "(I6.6,A1,I6.6)") para_env%mepos, ":", para_env%num_pe
     601            0 :          IF (ASSOCIATED(trace_routines)) THEN
     602            0 :             CALL timings_setup_tracing(trace_max, unit_nr, tracing_string, trace_routines)
     603              :          ELSE
     604            0 :             CALL timings_setup_tracing(trace_max, unit_nr, tracing_string)
     605              :          END IF
     606              :       END IF
     607              : 
     608        10420 :       CALL section_vals_val_get(global_section, "TIMINGS%TIME_MPI", l_val=mp_collect_timings)
     609              : 
     610          198 :       SELECT CASE (i_diag)
     611              :       CASE (FM_DIAG_TYPE_SCALAPACK)
     612          198 :          globenv%diag_library = "ScaLAPACK"
     613              :       CASE (FM_DIAG_TYPE_ELPA)
     614        10222 :          globenv%diag_library = "ELPA"
     615        10222 :          CALL cite_reference(Marek2014)
     616              :       CASE (FM_DIAG_TYPE_CUSOLVER)
     617            0 :          globenv%diag_library = "cuSOLVER"
     618              :       CASE (FM_DIAG_TYPE_DLAF)
     619            0 :          globenv%diag_library = "DLAF"
     620            0 :          CALL cite_reference(Solca2024)
     621              :       CASE DEFAULT
     622        10420 :          CPABORT("Unknown diagonalization library specified")
     623              :       END SELECT
     624              : 
     625        10420 :       SELECT CASE (i_cholesky)
     626              :       CASE (FM_CHOLESKY_TYPE_SCALAPACK)
     627        10420 :          globenv%cholesky_library = "ScaLAPACK"
     628        10420 :          cholesky_type = FM_CHOLESKY_TYPE_SCALAPACK
     629              :       CASE (FM_CHOLESKY_TYPE_DLAF)
     630            0 :          globenv%cholesky_library = "DLAF"
     631            0 :          cholesky_type = FM_CHOLESKY_TYPE_DLAF
     632            0 :          dlaf_cholesky_n_min = globenv%dlaf_cholesky_n_min
     633            0 :          CALL cite_reference(Solca2024)
     634              :       CASE DEFAULT
     635        10420 :          CPABORT("Unknown Cholesky decomposition library specified")
     636              :       END SELECT
     637              : 
     638           10 :       SELECT CASE (i_fft)
     639              :       CASE (do_fft_sg)
     640           10 :          globenv%default_fft_library = "FFTSG"
     641              :       CASE (do_fft_fftw3)
     642        10410 :          globenv%default_fft_library = "FFTW3"
     643        10410 :          CALL cite_reference(Frigo2005)
     644              :       CASE DEFAULT
     645        10420 :          CPABORT("Unknown FFT library specified")
     646              :       END SELECT
     647              : 
     648            0 :       SELECT CASE (i_dgemm)
     649              :       CASE (do_dgemm_spla)
     650            0 :          globenv%default_dgemm_library = "SPLA"
     651              :       CASE (do_dgemm_blas)
     652        10420 :          globenv%default_dgemm_library = "BLAS"
     653              :       CASE DEFAULT
     654        10420 :          CPABORT("Unknown DGEMM library specified")
     655              :       END SELECT
     656              : 
     657        10420 :       IF (globenv%run_type_id == 0) THEN
     658            0 :          SELECT CASE (globenv%prog_name_id)
     659              :          CASE (do_farming, do_test)
     660            0 :             globenv%run_type_id = none_run
     661              :          CASE (do_cp2k)
     662            0 :             IF (nforce_eval /= 1) THEN
     663              :                ! multiple force_eval corresponds at the moment to RESPA calculations only
     664              :                ! default MD
     665            0 :                globenv%run_type_id = mol_dyn_run
     666              :             ELSE
     667            0 :                CALL section_vals_val_get(force_env_sections, "METHOD", i_val=method_name_id)
     668            0 :                SELECT CASE (method_name_id)
     669              :                CASE (do_fist)
     670            0 :                   globenv%run_type_id = mol_dyn_run
     671              :                CASE (do_eip)
     672            0 :                   globenv%run_type_id = mol_dyn_run
     673              :                CASE (do_qs)
     674            0 :                   globenv%run_type_id = energy_run
     675              :                CASE (do_sirius)
     676            0 :                   globenv%run_type_id = energy_run
     677              :                END SELECT
     678              :             END IF
     679              :          END SELECT
     680              :       END IF
     681              : 
     682        10420 :       IF (globenv%prog_name_id == do_farming .AND. globenv%run_type_id /= none_run) THEN
     683            0 :          CPABORT("FARMING program supports only NONE as run type")
     684              :       END IF
     685              : 
     686        10420 :       IF (globenv%prog_name_id == do_test .AND. globenv%run_type_id /= none_run) &
     687            0 :          CPABORT("TEST program supports only NONE as run type")
     688              : 
     689        10420 :       CALL m_memory_details(MemTotal, MemFree, Buffers, Cached, Slab, SReclaimable, MemLikelyFree)
     690        10420 :       MemTotal_avr = MemTotal
     691        10420 :       MemFree_avr = MemFree
     692        10420 :       Buffers_avr = Buffers
     693        10420 :       Cached_avr = Cached
     694        10420 :       Slab_avr = Slab
     695        10420 :       SReclaimable_avr = SReclaimable
     696        10420 :       MemLikelyFree_avr = MemLikelyFree
     697        10420 :       CALL para_env%sum(MemTotal_avr); MemTotal_avr = MemTotal_avr/para_env%num_pe/1024
     698        10420 :       CALL para_env%sum(MemFree_avr); MemFree_avr = MemFree_avr/para_env%num_pe/1024
     699        10420 :       CALL para_env%sum(Buffers_avr); Buffers_avr = Buffers_avr/para_env%num_pe/1024
     700        10420 :       CALL para_env%sum(Cached_avr); Cached_avr = Cached_avr/para_env%num_pe/1024
     701        10420 :       CALL para_env%sum(Slab_avr); Slab_avr = Slab_avr/para_env%num_pe/1024
     702        10420 :       CALL para_env%sum(SReclaimable_avr); SReclaimable_avr = SReclaimable_avr/para_env%num_pe/1024
     703        10420 :       CALL para_env%sum(MemLikelyFree_avr); MemLikelyFree_avr = MemLikelyFree_avr/para_env%num_pe/1024
     704              : 
     705        10420 :       MemTotal_min = -MemTotal
     706        10420 :       MemFree_min = -MemFree
     707        10420 :       Buffers_min = -Buffers
     708        10420 :       Cached_min = -Cached
     709        10420 :       Slab_min = -Slab
     710        10420 :       SReclaimable_min = -SReclaimable
     711        10420 :       MemLikelyFree_min = -MemLikelyFree
     712        10420 :       CALL para_env%max(MemTotal_min); MemTotal_min = -MemTotal_min/1024
     713        10420 :       CALL para_env%max(MemFree_min); MemFree_min = -MemFree_min/1024
     714        10420 :       CALL para_env%max(Buffers_min); Buffers_min = -Buffers_min/1024
     715        10420 :       CALL para_env%max(Cached_min); Cached_min = -Cached_min/1024
     716        10420 :       CALL para_env%max(Slab_min); Slab_min = -Slab_min/1024
     717        10420 :       CALL para_env%max(SReclaimable_min); SReclaimable_min = -SReclaimable_min/1024
     718        10420 :       CALL para_env%max(MemLikelyFree_min); MemLikelyFree_min = -MemLikelyFree_min/1024
     719              : 
     720        10420 :       MemTotal_max = MemTotal
     721        10420 :       MemFree_max = MemFree
     722        10420 :       Buffers_max = Buffers
     723        10420 :       Cached_max = Cached
     724        10420 :       Slab_max = Slab
     725        10420 :       SReclaimable_max = SReclaimable
     726        10420 :       MemLikelyFree_max = MemLikelyFree
     727        10420 :       CALL para_env%max(MemTotal_max); MemTotal_max = MemTotal_max/1024
     728        10420 :       CALL para_env%max(MemFree_max); MemFree_max = MemFree_max/1024
     729        10420 :       CALL para_env%max(Buffers_max); Buffers_max = Buffers_max/1024
     730        10420 :       CALL para_env%max(Cached_max); Cached_max = Cached_max/1024
     731        10420 :       CALL para_env%max(Slab_max); Slab_max = Slab_max/1024
     732        10420 :       CALL para_env%max(SReclaimable_max); SReclaimable_max = SReclaimable_max/1024
     733        10420 :       CALL para_env%max(MemLikelyFree_max); MemLikelyFree_max = MemLikelyFree_max/1024
     734              : 
     735        10420 :       MemTotal = MemTotal/1024
     736        10420 :       MemFree = MemFree/1024
     737        10420 :       Buffers = Buffers/1024
     738        10420 :       Cached = Cached/1024
     739        10420 :       Slab = Slab/1024
     740        10420 :       SReclaimable = SReclaimable/1024
     741        10420 :       MemLikelyFree = MemLikelyFree/1024
     742              : 
     743        10420 :       node_count = 1
     744              :       ! Print a list of all started processes
     745        10420 :       IF (do_echo_all_hosts) THEN
     746            2 :          CALL echo_all_hosts(para_env, output_unit)
     747              : 
     748              :          ! Print the number of processes per host
     749            2 :          CALL echo_all_process_host(para_env, output_unit, node_count)
     750              :       ELSE ! no echo
     751        10418 :          CALL echo_all_process_host(para_env, 0, node_count)
     752              :       END IF
     753              : 
     754        10420 :       num_threads = 1
     755        10420 : !$    num_threads = omp_get_max_threads()
     756        10420 :       IF (output_unit > 0) THEN
     757         5314 :          WRITE (UNIT=output_unit, FMT=*)
     758         5314 :          CALL multiple_fe_list(force_env_sections, root_section, i_force_eval, nforce_eval)
     759        10568 :          DO iforce_eval = 1, nforce_eval
     760              :             dft_section => section_vals_get_subs_vals3(force_env_sections, "DFT", &
     761         5254 :                                                        i_rep_section=i_force_eval(iforce_eval))
     762              :             qmmm_section => section_vals_get_subs_vals3(force_env_sections, "QMMM", &
     763         5254 :                                                         i_rep_section=i_force_eval(iforce_eval))
     764              :             CALL section_vals_val_get(dft_section, "BASIS_SET_FILE_NAME", &
     765         5254 :                                       c_val=basis_set_file_name)
     766              :             CALL section_vals_val_get(dft_section, "POTENTIAL_FILE_NAME", &
     767         5254 :                                       c_val=potential_file_name)
     768              : 
     769              :             CALL section_vals_val_get(qmmm_section, "MM_POTENTIAL_FILE_NAME", &
     770         5254 :                                       c_val=mm_potential_file_name)
     771              :             ! SUBSYS - If any
     772              :             subsys_section => section_vals_get_subs_vals3(force_env_sections, "SUBSYS", &
     773         5254 :                                                           i_rep_section=i_force_eval(iforce_eval))
     774         5254 :             CALL section_vals_get(subsys_section, explicit=explicit)
     775         5254 :             coord_file_name = "__STD_INPUT__"
     776         5254 :             IF (explicit) THEN
     777              :                CALL section_vals_val_get(subsys_section, "TOPOLOGY%COORD_FILE_NAME", &
     778         5147 :                                          n_rep_val=n_rep_val)
     779         5147 :                IF (n_rep_val == 1) THEN
     780              :                   CALL section_vals_val_get(subsys_section, "TOPOLOGY%COORD_FILE_NAME", &
     781          922 :                                             c_val=coord_file_name)
     782              :                END IF
     783              :             END IF
     784         5254 :             CALL integer_to_string(i_force_eval(iforce_eval), env_num)
     785              : 
     786              :             WRITE (UNIT=output_unit, FMT="(T2,A,T41,A)") &
     787         5254 :                start_section_label//"| Force Environment number", &
     788         5254 :                ADJUSTR(env_num(:40)), &
     789         5254 :                start_section_label//"| Basis set file name", &
     790         5254 :                ADJUSTR(basis_set_file_name(:40)), &
     791         5254 :                start_section_label//"| Potential file name", &
     792         5254 :                ADJUSTR(potential_file_name(:40)), &
     793         5254 :                start_section_label//"| MM Potential file name", &
     794         5254 :                ADJUSTR(mm_potential_file_name(:40)), &
     795         5254 :                start_section_label//"| Coordinate file name", &
     796        21076 :                ADJUSTR(coord_file_name(:40))
     797              :          END DO
     798         5314 :          DEALLOCATE (i_force_eval)
     799              : 
     800         5314 :          NULLIFY (enum1, enum2, keyword, section)
     801         5314 :          CALL create_global_section(section)
     802         5314 :          keyword => section_get_keyword(section, "PROGRAM_NAME")
     803         5314 :          CALL keyword_get(keyword, enum=enum1)
     804         5314 :          keyword => section_get_keyword(section, "RUN_TYPE")
     805         5314 :          CALL keyword_get(keyword, enum=enum2)
     806              : 
     807              :          WRITE (UNIT=output_unit, FMT="(T2,A,T41,A40)") &
     808         5314 :             start_section_label//"| Method name", &
     809         5314 :             ADJUSTR(TRIM(enum_i2c(enum1, globenv%prog_name_id))), &
     810         5314 :             start_section_label//"| Project name", &
     811         5314 :             ADJUSTR(project_name(:40)), &
     812         5314 :             start_section_label//"| Run type", &
     813         5314 :             ADJUSTR(TRIM(enum_i2c(enum2, globenv%run_type_id))), &
     814         5314 :             start_section_label//"| FFT library", &
     815         5314 :             ADJUSTR(globenv%default_fft_library(:40)), &
     816         5314 :             start_section_label//"| Diagonalization library", &
     817         5314 :             ADJUSTR(globenv%diag_library(:40)), &
     818         5314 :             start_section_label//"| Cholesky decomposition library", &
     819         5314 :             ADJUSTR(globenv%cholesky_library(:40)), &
     820         5314 :             start_section_label//"| DGEMM library", &
     821        10628 :             ADJUSTR(globenv%default_dgemm_library(:40))
     822              : 
     823         5314 :          IF (globenv%diag_library == "ELPA") THEN
     824              :             WRITE (UNIT=output_unit, FMT="(T2,A,T71,I10)") &
     825         5214 :                start_section_label//"| Minimum number of eigenvectors for ELPA usage", &
     826        10428 :                globenv%elpa_neigvec_min
     827              :          END IF
     828              : 
     829         5314 :          IF (globenv%diag_library == "DLAF") THEN
     830              :             WRITE (UNIT=output_unit, FMT="(T2,A,T71,I10)") &
     831            0 :                start_section_label//"| Minimum number of eigenvectors for DLAF usage", &
     832            0 :                globenv%dlaf_neigvec_min
     833              :          END IF
     834              : 
     835         5314 :          IF (globenv%diag_library == "cuSOLVER" .OR. globenv%diag_library == "ScaLAPACK" .OR. &
     836              :              globenv%diag_library == "DLAF") THEN
     837              :             WRITE (UNIT=output_unit, FMT="(T2,A,T71,L10)") &
     838          100 :                start_section_label//"| Direct generalized diagonalization requested", &
     839          200 :                globenv%direct_generalized_diagonalization
     840              :          END IF
     841              : 
     842         5314 :          IF (globenv%diag_library == "cuSOLVER") THEN
     843              :             WRITE (UNIT=output_unit, FMT="(T2,A,T71,I10)") &
     844            0 :                start_section_label//"| Minimum matrix size for cuSOLVER diagonalization", &
     845            0 :                cusolver_n_min
     846              :          END IF
     847              : 
     848         5314 :          IF (globenv%cholesky_library == "DLAF") THEN
     849              :             WRITE (UNIT=output_unit, FMT="(T2,A,T71,I10)") &
     850            0 :                start_section_label//"| Minimum matrix size for Cholesky decomposition with DLAF", &
     851            0 :                globenv%dlaf_cholesky_n_min
     852              :          END IF
     853              : 
     854              : #if defined(__CHECK_DIAG)
     855              :          ! Perform default check if no threshold value has been specified explicitly
     856              :          IF (globenv%eps_check_diag < 0.0_dp) THEN
     857              :             WRITE (UNIT=output_unit, FMT="(T2,A,T71,ES10.3)") &
     858              :                start_section_label//"| Orthonormality check for eigenvectors enabled", &
     859              :                eps_check_diag_default
     860              :          ELSE
     861              :             WRITE (UNIT=output_unit, FMT="(T2,A,T71,ES10.3)") &
     862              :                start_section_label//"| Orthonormality check for eigenvectors enabled", &
     863              :                globenv%eps_check_diag
     864              :          END IF
     865              : #else
     866         5314 :          IF (globenv%eps_check_diag < 0.0_dp) THEN
     867              :             WRITE (UNIT=output_unit, FMT="(T2,A,T73,A)") &
     868         5306 :                start_section_label//"| Orthonormality check for eigenvectors", &
     869        10612 :                "DISABLED"
     870              :          ELSE
     871              :             WRITE (UNIT=output_unit, FMT="(T2,A,T71,ES10.3)") &
     872            8 :                start_section_label//"| Orthonormality check for eigenvectors enabled", &
     873           16 :                globenv%eps_check_diag
     874              :          END IF
     875              : #endif
     876         5314 :          CALL section_release(section)
     877              : 
     878            0 :          SELECT CASE (cp_fm_get_mm_type())
     879              :          CASE (do_scalapack)
     880              :             WRITE (UNIT=output_unit, FMT="(T2,A,T72,A)") &
     881            0 :                start_section_label//"| Matrix multiplication library", "ScaLAPACK"
     882              :          CASE (do_cosma)
     883              :             WRITE (UNIT=output_unit, FMT="(T2,A,T76,A)") &
     884         5314 :                start_section_label//"| Matrix multiplication library", "COSMA"
     885              :          END SELECT
     886              : 
     887         5314 :          CALL section_vals_val_get(global_section, "ALLTOALL_SGL", l_val=ata)
     888              :          WRITE (UNIT=output_unit, FMT="(T2,A,T80,L1)") &
     889         5314 :             start_section_label//"| All-to-all communication in single precision", ata
     890         5314 :          CALL section_vals_val_get(global_section, "EXTENDED_FFT_LENGTHS", l_val=efl)
     891              :          WRITE (UNIT=output_unit, FMT="(T2,A,T80,L1)") &
     892         5314 :             start_section_label//"| FFTs using library dependent lengths", efl
     893              : 
     894          136 :          SELECT CASE (print_level)
     895              :          CASE (silent_print_level)
     896          136 :             print_level_string = "SILENT"
     897              :          CASE (low_print_level)
     898         2536 :             print_level_string = "   LOW"
     899              :          CASE (medium_print_level)
     900         2577 :             print_level_string = "MEDIUM"
     901              :          CASE (high_print_level)
     902           39 :             print_level_string = "  HIGH"
     903              :          CASE (debug_print_level)
     904           26 :             print_level_string = " DEBUG"
     905              :          CASE DEFAULT
     906         5314 :             CPABORT("Unknown print_level")
     907              :          END SELECT
     908              : 
     909         5314 :          CALL section_vals_val_get(global_section, "GRID%BACKEND", i_val=i_grid_backend)
     910         5304 :          SELECT CASE (i_grid_backend)
     911              :          CASE (GRID_BACKEND_AUTO)
     912              :             WRITE (UNIT=output_unit, FMT="(T2,A,T75,A6)") &
     913         5304 :                start_section_label//"| Grid backend", "AUTO"
     914              :          CASE (GRID_BACKEND_CPU)
     915              :             WRITE (UNIT=output_unit, FMT="(T2,A,T75,A6)") &
     916            4 :                start_section_label//"| Grid backend", "CPU"
     917              :          CASE (GRID_BACKEND_DGEMM)
     918              :             WRITE (UNIT=output_unit, FMT="(T2,A,T75,A6)") &
     919            4 :                start_section_label//"| Grid backend", "DGEMM"
     920              :          CASE (GRID_BACKEND_GPU)
     921              :             WRITE (UNIT=output_unit, FMT="(T2,A,T75,A6)") &
     922            0 :                start_section_label//"| Grid backend", "GPU"
     923              :          CASE (GRID_BACKEND_HIP)
     924              :             WRITE (UNIT=output_unit, FMT="(T2,A,T75,A6)") &
     925            0 :                start_section_label//"| Grid backend", "HIP"
     926              :          CASE (GRID_BACKEND_REF)
     927              :             WRITE (UNIT=output_unit, FMT="(T2,A,T75,A6)") &
     928         5314 :                start_section_label//"| Grid backend", "REF"
     929              :          END SELECT
     930              : 
     931              :          WRITE (UNIT=output_unit, FMT="(T2,A,T75,A6)") &
     932         5314 :             start_section_label//"| Global print level", print_level_string
     933              :          WRITE (UNIT=output_unit, FMT="(T2,A,T75,L6)") &
     934         5314 :             start_section_label//"| MPI I/O enabled", flag
     935              :          WRITE (UNIT=output_unit, FMT="(T2,A,T75,I6)") &
     936         5314 :             start_section_label//"| Total number of message passing processes", &
     937         5314 :             para_env%num_pe, &
     938         5314 :             start_section_label//"| Number of distributed systems (nodes)", &
     939         5314 :             node_count, &
     940         5314 :             start_section_label//"| Number of threads for this process", &
     941         5314 :             num_threads, &
     942        10628 :             start_section_label//"| This output is from process", para_env%mepos
     943              : 
     944         5314 :          CALL m_omp_get_stacksize(omp_stacksize)
     945              :          WRITE (UNIT=output_unit, FMT="(T2,A,T68,A13)") &
     946         5314 :             start_section_label//"| OpenMP stack size per thread (OMP_STACKSIZE)", &
     947        10628 :             ADJUSTR(omp_stacksize)
     948              : 
     949         5314 :          IF (0 <= m_omp_trace_issues()) THEN ! only show in header if enabled
     950              :             WRITE (UNIT=output_unit, FMT="(T2,A,T68,A13)") &
     951            0 :                start_section_label//"| OpenMP issue trace (CP2K_OMP_TRACE)", &
     952            0 :                "enabled"
     953              :          END IF
     954              : 
     955         5314 :          CALL m_cpuinfo(model_name)
     956              :          WRITE (UNIT=output_unit, FMT="(T2,A,T30,A51)") &
     957         5314 :             start_section_label//"| CPU model name", ADJUSTR(TRIM(model_name))
     958              : 
     959         5314 :          cpuid = m_cpuid()
     960         5314 :          cpuid_static = m_cpuid_static()
     961              : 
     962         5314 :          IF ((cpuid > 0) .OR. (cpuid_static > 0)) THEN
     963              :             WRITE (UNIT=output_unit, FMT="(T2,A,T75,I6)") &
     964         5314 :                start_section_label//"| CPUID", cpuid
     965         5314 :             IF (cpuid /= cpuid_static) THEN
     966              :                WRITE (UNIT=output_unit, FMT="(T2,A,T75,I6)") &
     967            0 :                   start_section_label//"| Compiled for CPUID", cpuid_static
     968              :             END IF
     969              :          END IF
     970              : 
     971              :          ! filter cpuids by vlen to show more relevant information
     972         5314 :          IF (m_cpuid_vlen(cpuid_static) < m_cpuid_vlen(cpuid)) THEN
     973              :             ! base/machine_cpuid.c relies on the (same) target flags as the Fortran code
     974              :             CALL cp_hint(__LOCATION__, "The compiler target flags ("// &
     975              :                          TRIM(m_cpuid_name(cpuid_static))//") used to build this binary cannot exploit "// &
     976              :                          "all extensions of this CPU model ("//TRIM(m_cpuid_name(cpuid))//"). "// &
     977            0 :                          "Consider compiler target flags as part of FCFLAGS and CFLAGS (ARCH file).")
     978              :          END IF
     979              : 
     980         5314 :          WRITE (UNIT=output_unit, FMT="()")
     981         5314 :          WRITE (UNIT=output_unit, FMT="(T2,A)") "MEMORY| system memory details [Kb]"
     982         5314 :          WRITE (UNIT=output_unit, FMT="(T2,A23,4A14)") "MEMORY|                ", "rank 0", "min", "max", "average"
     983         5314 :          WRITE (UNIT=output_unit, FMT="(T2,A23,4I14)") "MEMORY| MemTotal       ", memtotal, memtotal_min, memtotal_max, memtotal_avr
     984         5314 :          WRITE (UNIT=output_unit, FMT="(T2,A23,4I14)") "MEMORY| MemFree        ", memFree, memfree_min, memfree_max, memfree_avr
     985         5314 :          WRITE (UNIT=output_unit, FMT="(T2,A23,4I14)") "MEMORY| Buffers        ", Buffers, Buffers_min, Buffers_max, Buffers_avr
     986         5314 :          WRITE (UNIT=output_unit, FMT="(T2,A23,4I14)") "MEMORY| Cached         ", Cached, Cached_min, Cached_max, Cached_avr
     987         5314 :          WRITE (UNIT=output_unit, FMT="(T2,A23,4I14)") "MEMORY| Slab           ", Slab, Slab_min, Slab_max, Slab_avr
     988              :          WRITE (UNIT=output_unit, FMT="(T2,A23,4I14)") &
     989         5314 :             "MEMORY| SReclaimable   ", SReclaimable, SReclaimable_min, SReclaimable_max, &
     990        10628 :             SReclaimable_avr
     991              :          WRITE (UNIT=output_unit, FMT="(T2,A23,4I14)") &
     992         5314 :             "MEMORY| MemLikelyFree  ", MemLikelyFree, MemLikelyFree_min, MemLikelyFree_max, &
     993        10628 :             MemLikelyFree_avr
     994         5314 :          WRITE (UNIT=output_unit, FMT='()')
     995              : 
     996              :       END IF
     997              : 
     998              :       CALL cp_print_key_finished_output(output_unit, logger, global_section, &
     999        10420 :                                         "PROGRAM_RUN_INFO")
    1000              : 
    1001        10420 :    END SUBROUTINE read_global_section
    1002              : 
    1003              : ! **************************************************************************************************
    1004              : !> \brief ...
    1005              : !> \param root_section ...
    1006              : !> \param para_env ...
    1007              : !> \param globenv ...
    1008              : !> \par History
    1009              : !>      2-Dec-2000 (JGH) added default fft library
    1010              : !> \author JGH,MK
    1011              : ! **************************************************************************************************
    1012        10420 :    SUBROUTINE read_cp2k_section(root_section, para_env, globenv)
    1013              : 
    1014              :       TYPE(section_vals_type), POINTER                   :: root_section
    1015              :       TYPE(mp_para_env_type), POINTER                    :: para_env
    1016              :       TYPE(global_environment_type), POINTER             :: globenv
    1017              : 
    1018              :       INTEGER                                            :: output_unit
    1019              :       TYPE(cp_logger_type), POINTER                      :: logger
    1020              :       TYPE(section_vals_type), POINTER                   :: global_section
    1021              : 
    1022        10420 :       global_section => section_vals_get_subs_vals(root_section, "GLOBAL")
    1023        10420 :       CALL read_global_section(root_section, para_env, globenv)
    1024        10420 :       logger => cp_get_default_logger()
    1025              :       output_unit = cp_print_key_unit_nr(logger, global_section, "PROGRAM_RUN_INFO", &
    1026        10420 :                                          extension=".log")
    1027              : 
    1028        10420 :       CALL fft_setup_library(globenv, global_section)
    1029        10420 :       CALL diag_setup_library(globenv)
    1030              : 
    1031              :       CALL cp_print_key_finished_output(output_unit, logger, global_section, &
    1032        10420 :                                         "PROGRAM_RUN_INFO")
    1033              : 
    1034        10420 :    END SUBROUTINE read_cp2k_section
    1035              : 
    1036              : ! **************************************************************************************************
    1037              : !> \brief check FFT preferred library availability, if not switch
    1038              : !> \param globenv ...
    1039              : !> \param global_section ...
    1040              : !> \par History
    1041              : !>      2-Dec-2000 (JGH) added default fft library
    1042              : !>      Nov-2013 (MI) refactoring
    1043              : !> \author JGH,MK
    1044              : ! **************************************************************************************************
    1045        10420 :    SUBROUTINE fft_setup_library(globenv, global_section)
    1046              : 
    1047              :       TYPE(global_environment_type), POINTER             :: globenv
    1048              :       TYPE(section_vals_type), POINTER                   :: global_section
    1049              : 
    1050              :       CHARACTER(LEN=3*default_string_length)             :: message
    1051              :       COMPLEX(KIND=dp), DIMENSION(4, 4, 4)               :: zz
    1052              :       INTEGER                                            :: stat
    1053              :       INTEGER, DIMENSION(3)                              :: n
    1054              :       LOGICAL                                            :: try_fftw
    1055              : 
    1056        41680 :       n(:) = 4
    1057        10420 :       zz(:, :, :) = 0.0_dp
    1058              : 
    1059              :       ! Setup the FFT library
    1060              :       ! If the user has specified PREFERRED_FFT_LIBRARY try that first (default FFTW3)
    1061              :       ! If that one is not available, try FFTW3 (unless it has been tried already)
    1062              :       ! If FFTW3 is not available use FFTSG
    1063              : 
    1064        10420 :       IF (globenv%default_fft_library == "FFTW3") THEN
    1065              :          try_fftw = .FALSE.
    1066              :       ELSE
    1067           10 :          try_fftw = .TRUE.
    1068              :       END IF
    1069              : 
    1070              :       ! Initialize FFT library with the user's preferred FFT library
    1071              :       CALL init_fft(fftlib=TRIM(globenv%default_fft_library), &
    1072              :                     alltoall=section_get_lval(global_section, "ALLTOALL_SGL"), &
    1073              :                     fftsg_sizes=.NOT. section_get_lval(global_section, "EXTENDED_FFT_LENGTHS"), &
    1074              :                     pool_limit=globenv%fft_pool_scratch_limit, &
    1075              :                     wisdom_file=globenv%fftw_wisdom_file_name, &
    1076        10420 :                     plan_style=globenv%fftw_plan_type)
    1077              : 
    1078              :       ! Check for FFT library
    1079        10420 :       CALL fft3d(FWFFT, n, zz, status=stat)
    1080        10420 :       IF (stat /= 0) THEN
    1081            0 :          IF (try_fftw) THEN
    1082              :             message = "FFT library "//TRIM(globenv%default_fft_library)// &
    1083            0 :                       " is not available. Trying FFT library FFTW3."
    1084            0 :             CPWARN(TRIM(message))
    1085            0 :             globenv%default_fft_library = "FFTW3"
    1086              :             CALL init_fft(fftlib=TRIM(globenv%default_fft_library), &
    1087              :                           alltoall=section_get_lval(global_section, "ALLTOALL_SGL"), &
    1088              :                           fftsg_sizes=.NOT. section_get_lval(global_section, "EXTENDED_FFT_LENGTHS"), &
    1089              :                           pool_limit=globenv%fft_pool_scratch_limit, &
    1090              :                           wisdom_file=globenv%fftw_wisdom_file_name, &
    1091            0 :                           plan_style=globenv%fftw_plan_type)
    1092              : 
    1093            0 :             CALL fft3d(FWFFT, n, zz, status=stat)
    1094              :          END IF
    1095            0 :          IF (stat /= 0) THEN
    1096              :             message = "FFT library "//TRIM(globenv%default_fft_library)// &
    1097            0 :                       " is not available. Trying FFT library FFTSG."
    1098            0 :             CPWARN(TRIM(message))
    1099            0 :             globenv%default_fft_library = "FFTSG"
    1100              :             CALL init_fft(fftlib=TRIM(globenv%default_fft_library), &
    1101              :                           alltoall=section_get_lval(global_section, "ALLTOALL_SGL"), &
    1102              :                           fftsg_sizes=.NOT. section_get_lval(global_section, "EXTENDED_FFT_LENGTHS"), &
    1103              :                           pool_limit=globenv%fft_pool_scratch_limit, &
    1104              :                           wisdom_file=globenv%fftw_wisdom_file_name, &
    1105            0 :                           plan_style=globenv%fftw_plan_type)
    1106              : 
    1107            0 :             CALL fft3d(FWFFT, n, zz, status=stat)
    1108            0 :             IF (stat /= 0) THEN
    1109            0 :                CPABORT("FFT library FFTSG does not work. No FFT library available.")
    1110              :             END IF
    1111              :          END IF
    1112              :       END IF
    1113              : 
    1114        10420 :    END SUBROUTINE fft_setup_library
    1115              : 
    1116              : ! **************************************************************************************************
    1117              : !> \brief availability diagonalizatioon library
    1118              : !>
    1119              : !> \param globenv ...
    1120              : !> \author MI
    1121              : ! **************************************************************************************************
    1122        10420 :    SUBROUTINE diag_setup_library(globenv)
    1123              :       TYPE(global_environment_type), POINTER             :: globenv
    1124              : 
    1125              :       CHARACTER(LEN=3*default_string_length)             :: message
    1126              :       LOGICAL                                            :: fallback_applied
    1127              : 
    1128              :       CALL diag_init(diag_lib=TRIM(globenv%diag_library), &
    1129              :                      fallback_applied=fallback_applied, &
    1130              :                      elpa_kernel=globenv%k_elpa, &
    1131              :                      elpa_neigvec_min_input=globenv%elpa_neigvec_min, &
    1132              :                      elpa_qr=globenv%elpa_qr, &
    1133              :                      elpa_print=globenv%elpa_print, &
    1134              :                      elpa_one_stage=globenv%elpa_one_stage, &
    1135              :                      dlaf_neigvec_min_input=globenv%dlaf_neigvec_min, &
    1136              :                      eps_check_diag_input=globenv%eps_check_diag, &
    1137              :                      direct_generalized_diagonalization_input= &
    1138        10420 :                      globenv%direct_generalized_diagonalization)
    1139              : 
    1140        10420 :       IF (fallback_applied) THEN
    1141              :          message = "Diagonalization library "//TRIM(globenv%diag_library)// &
    1142            0 :                    " is not available. The ScaLAPACK library is used as fallback."
    1143            0 :          CPWARN(TRIM(message))
    1144              :       END IF
    1145              : 
    1146        10420 :    END SUBROUTINE diag_setup_library
    1147              : 
    1148              : ! **************************************************************************************************
    1149              : !> \brief ...
    1150              : !> \param glob_section ...
    1151              : ! **************************************************************************************************
    1152        52100 :    SUBROUTINE fm_setup(glob_section)
    1153              :       TYPE(section_vals_type), POINTER                   :: glob_section
    1154              : 
    1155              :       INTEGER                                            :: mm_type, ncb, nrb
    1156              :       LOGICAL                                            :: force_me
    1157              :       TYPE(section_vals_type), POINTER                   :: fm_section
    1158              : 
    1159        10420 :       fm_section => section_vals_get_subs_vals(glob_section, "FM")
    1160              : 
    1161        10420 :       CALL section_vals_val_get(fm_section, "NROW_BLOCKS", i_val=nrb)
    1162        10420 :       CALL section_vals_val_get(fm_section, "NCOL_BLOCKS", i_val=ncb)
    1163        10420 :       CALL section_vals_val_get(fm_section, "FORCE_BLOCK_SIZE", l_val=force_me)
    1164        10420 :       CALL cp_fm_struct_config(nrow_block=nrb, ncol_block=ncb, force_block=force_me)
    1165              : 
    1166        10420 :       CALL section_vals_val_get(fm_section, "TYPE_OF_MATRIX_MULTIPLICATION", i_val=mm_type)
    1167        10420 :       CALL cp_fm_setup(mm_type)
    1168              : 
    1169        10420 :    END SUBROUTINE fm_setup
    1170              : 
    1171              : ! **************************************************************************************************
    1172              : !> \brief ...
    1173              : !> \param glob_section ...
    1174              : ! **************************************************************************************************
    1175        10420 :    SUBROUTINE dgemm_setup(glob_section)
    1176              :       TYPE(section_vals_type), POINTER                   :: glob_section
    1177              : 
    1178              :       INTEGER                                            :: dgemm_type
    1179              : 
    1180        10420 :       CALL section_vals_val_get(glob_section, "PREFERRED_DGEMM_LIBRARY", i_val=dgemm_type)
    1181              : 
    1182        10420 :       CALL local_gemm_set_library(dgemm_type)
    1183              : 
    1184        10420 :    END SUBROUTINE dgemm_setup
    1185              : 
    1186              : ! **************************************************************************************************
    1187              : !> \brief   Parses the input section used to define the heuristic rules which determine if
    1188              : !>          a FM matrix should be redistributed before diagonalizing it.
    1189              : !> \param glob_section the global input section
    1190              : !> \author Nico Holmberg [01.2018]
    1191              : ! **************************************************************************************************
    1192        52100 :    SUBROUTINE fm_diag_rules_setup(glob_section)
    1193              :       TYPE(section_vals_type), POINTER                   :: glob_section
    1194              : 
    1195              :       INTEGER                                            :: a, x
    1196              :       LOGICAL                                            :: elpa_force_redistribute, should_print
    1197              :       TYPE(section_vals_type), POINTER                   :: section
    1198              : 
    1199        10420 :       section => section_vals_get_subs_vals(glob_section, "FM_DIAG_SETTINGS")
    1200              : 
    1201        10420 :       CALL section_vals_val_get(section, "PARAMETER_A", i_val=a)
    1202        10420 :       CALL section_vals_val_get(section, "PARAMETER_X", i_val=x)
    1203        10420 :       CALL section_vals_val_get(section, "PRINT_FM_REDISTRIBUTE", l_val=should_print)
    1204        10420 :       CALL section_vals_val_get(section, "ELPA_FORCE_REDISTRIBUTE", l_val=elpa_force_redistribute)
    1205              : 
    1206        10420 :       CALL cp_fm_redistribute_init(a, x, should_print, elpa_force_redistribute)
    1207              : 
    1208        10420 :    END SUBROUTINE fm_diag_rules_setup
    1209              : ! **************************************************************************************************
    1210              : !> \brief reads the Walltime also in format HH:MM:SS
    1211              : !> \param section ...
    1212              : !> \param keyword_name ...
    1213              : !> \param walltime ...
    1214              : !> \par History
    1215              : !>      none
    1216              : !> \author Mandes
    1217              : ! **************************************************************************************************
    1218        10440 :    SUBROUTINE cp2k_get_walltime(section, keyword_name, walltime)
    1219              :       TYPE(section_vals_type), POINTER                   :: section
    1220              :       CHARACTER(LEN=*), INTENT(in)                       :: keyword_name
    1221              :       REAL(KIND=dp), INTENT(out)                         :: walltime
    1222              : 
    1223              :       CHARACTER(LEN=1)                                   :: c1, c2
    1224              :       CHARACTER(LEN=100)                                 :: txt
    1225              :       INTEGER                                            :: hours, ierr, minutes, n, seconds
    1226              : 
    1227        10440 :       CALL section_vals_val_get(section, keyword_name, c_val=txt)
    1228        10440 :       n = LEN_TRIM(txt)
    1229              : 
    1230        10440 :       IF (n == 0) THEN
    1231        10199 :          walltime = -1.0_dp
    1232          241 :       ELSE IF (INDEX(txt, ":") == 0) THEN
    1233          189 :          READ (txt(1:n), FMT=*, IOSTAT=ierr) walltime
    1234          189 :          IF (ierr /= 0) CPABORT('Could not parse WALLTIME: "'//txt(1:n)//'"')
    1235              :       ELSE
    1236           52 :          READ (txt(1:n), FMT="(I2,A1,I2,A1,I2)", IOSTAT=ierr) hours, c1, minutes, c2, seconds
    1237           52 :          IF (n /= 8 .OR. ierr /= 0 .OR. c1 /= ":" .OR. c2 /= ":") &
    1238            0 :             CPABORT('Could not parse WALLTIME: "'//txt(1:n)//'"')
    1239           52 :          walltime = 3600.0_dp*REAL(hours, dp) + 60.0_dp*REAL(minutes, dp) + REAL(seconds, dp)
    1240              :       END IF
    1241        10440 :    END SUBROUTINE cp2k_get_walltime
    1242              : 
    1243              : ! **************************************************************************************************
    1244              : !> \brief Writes final timings and banner for CP2K
    1245              : !> \param root_section ...
    1246              : !> \param para_env ...
    1247              : !> \param globenv ...
    1248              : !> \param wdir ...
    1249              : !> \param q_finalize ...
    1250              : !> \par History
    1251              : !>      none
    1252              : !> \author JGH,MK
    1253              : !> \note
    1254              : !>      The following routines need to be synchronized wrt. adding/removing
    1255              : !>      of the default environments (logging, performance,error):
    1256              : !>      environment:cp2k_init, environment:cp2k_finalize,
    1257              : !>      f77_interface:f_env_add_defaults, f77_interface:f_env_rm_defaults,
    1258              : !>      f77_interface:create_force_env, f77_interface:destroy_force_env
    1259              : ! **************************************************************************************************
    1260        20838 :    SUBROUTINE cp2k_finalize(root_section, para_env, globenv, wdir, q_finalize)
    1261              : 
    1262              :       TYPE(section_vals_type), POINTER                   :: root_section
    1263              :       TYPE(mp_para_env_type), POINTER                    :: para_env
    1264              :       TYPE(global_environment_type), POINTER             :: globenv
    1265              :       CHARACTER(LEN=*), OPTIONAL                         :: wdir
    1266              :       LOGICAL, INTENT(IN), OPTIONAL                      :: q_finalize
    1267              : 
    1268              :       CHARACTER(LEN=default_path_length)                 :: cg_filename
    1269              :       INTEGER                                            :: cg_mode, iw, unit_exit
    1270              :       LOGICAL                                            :: delete_it, do_finalize, report_maxloc, &
    1271              :                                                             sort_by_self_time
    1272              :       REAL(KIND=dp)                                      :: r_timings
    1273              :       TYPE(cp_logger_type), POINTER                      :: logger
    1274              : 
    1275              :       ! Look if we inherited a failure, more care is needed if so
    1276              :       ! i.e. the input is most likely not available
    1277              :       ! Set flag if this is a development version
    1278              : 
    1279        10419 :       do_finalize = .TRUE.
    1280        10419 :       IF (PRESENT(q_finalize)) do_finalize = q_finalize
    1281              :       ! Clean up
    1282        10419 :       NULLIFY (logger)
    1283        10419 :       logger => cp_get_default_logger()
    1284        10419 :       IF (do_finalize) THEN
    1285        10209 :          CALL deallocate_spherical_harmonics()
    1286        10209 :          CALL deallocate_orbital_pointers()
    1287        10209 :          CALL deallocate_md_ftable()
    1288        10209 :          CALL diag_finalize()
    1289              :          ! finalize the fft (i.e. writes the wisdom if FFTW3 )
    1290        10209 :          CALL finalize_fft(para_env, globenv%fftw_wisdom_file_name)
    1291        10209 :          CALL finalize_libvori()
    1292              :       END IF
    1293              : 
    1294              :       ! Write message passing performance info
    1295              : 
    1296              :       iw = cp_print_key_unit_nr(logger, root_section, "GLOBAL%PROGRAM_RUN_INFO", &
    1297        10419 :                                 extension=".log")
    1298        10419 :       CALL describe_mp_perf_env(iw)
    1299              :       CALL cp_print_key_finished_output(iw, logger, root_section, &
    1300        10419 :                                         "GLOBAL%PROGRAM_RUN_INFO")
    1301              : 
    1302        10419 :       CALL collect_citations_from_ranks(para_env)
    1303              :       iw = cp_print_key_unit_nr(logger, root_section, "GLOBAL%REFERENCES", &
    1304        10419 :                                 extension=".Log")
    1305        10419 :       IF (iw > 0) THEN
    1306         5305 :          WRITE (UNIT=iw, FMT="(/,T2,A)") REPEAT("-", 79)
    1307         5305 :          WRITE (UNIT=iw, FMT="(T2,A,T80,A)") "-", "-"
    1308         5305 :          WRITE (UNIT=iw, FMT="(T2,A,T30,A,T80,A)") "-", "R E F E R E N C E S", "-"
    1309         5305 :          WRITE (UNIT=iw, FMT="(T2,A,T80,A)") "-", "-"
    1310         5305 :          WRITE (UNIT=iw, FMT="(T2,A)") REPEAT("-", 79)
    1311         5305 :          WRITE (UNIT=iw, FMT="(T2,A)") ""
    1312         5305 :          WRITE (UNIT=iw, FMT="(T2,A)") TRIM(cp2k_version)//", the CP2K developers group ("//TRIM(cp2k_year)//")."
    1313         5305 :          WRITE (UNIT=iw, FMT="(T2,A)") "CP2K is freely available from "//TRIM(cp2k_home)//" ."
    1314         5305 :          WRITE (UNIT=iw, FMT="(T2,A)") ""
    1315         5305 :          CALL print_cited_references(unit=iw)
    1316              :       END IF
    1317              :       CALL cp_print_key_finished_output(iw, logger, root_section, &
    1318        10419 :                                         "GLOBAL%REFERENCES")
    1319              : 
    1320        10419 :       CALL timestop(globenv%handle) ! corresponding the "CP2K" in cp2k_init
    1321              : 
    1322              :       iw = cp_print_key_unit_nr(logger, root_section, "GLOBAL%TIMINGS", &
    1323        10419 :                                 extension=".Log")
    1324        10419 :       r_timings = section_get_rval(root_section, "GLOBAL%TIMINGS%THRESHOLD")
    1325        10419 :       sort_by_self_time = section_get_lval(root_section, "GLOBAL%TIMINGS%SORT_BY_SELF_TIME")
    1326        10419 :       report_maxloc = section_get_lval(root_section, "GLOBAL%TIMINGS%REPORT_MAXLOC")
    1327        10419 :       IF (m_energy() /= 0.0_dp) THEN
    1328            0 :          CALL timings_report_print(iw, r_timings, sort_by_self_time, cost_type_energy, report_maxloc, para_env)
    1329              :       END IF
    1330        10419 :       CALL timings_report_print(iw, r_timings, sort_by_self_time, cost_type_time, report_maxloc, para_env)
    1331              : 
    1332              :       ! Write the callgraph, if desired by user
    1333        10419 :       CALL section_vals_val_get(root_section, "GLOBAL%CALLGRAPH", i_val=cg_mode)
    1334        10419 :       IF (cg_mode /= CALLGRAPH_NONE) THEN
    1335            2 :          CALL section_vals_val_get(root_section, "GLOBAL%CALLGRAPH_FILE_NAME", c_val=cg_filename)
    1336            2 :          IF (LEN_TRIM(cg_filename) == 0) cg_filename = TRIM(logger%iter_info%project_name)
    1337            2 :          IF (cg_mode == CALLGRAPH_ALL) & !incorporate mpi-rank into filename
    1338            0 :             cg_filename = TRIM(cg_filename)//"_"//TRIM(ADJUSTL(cp_to_string(para_env%mepos)))
    1339            2 :          IF (iw > 0) THEN
    1340            1 :             WRITE (UNIT=iw, FMT="(T2,3X,A)") "Writing callgraph to: "//TRIM(cg_filename)//".callgraph"
    1341            1 :             WRITE (UNIT=iw, FMT="()")
    1342            1 :             WRITE (UNIT=iw, FMT="(T2,A)") "-------------------------------------------------------------------------------"
    1343              :          END IF
    1344            2 :          IF (cg_mode == CALLGRAPH_ALL .OR. para_env%is_source()) &
    1345            1 :             CALL timings_report_callgraph(TRIM(cg_filename)//".callgraph")
    1346              :       END IF
    1347              : 
    1348              :       CALL cp_print_key_finished_output(iw, logger, root_section, &
    1349        10419 :                                         "GLOBAL%TIMINGS")
    1350              : 
    1351        10419 :       CALL rm_mp_perf_env()
    1352        10419 :       CALL rm_timer_env()
    1353              : 
    1354        10419 :       IF (para_env%is_source()) THEN
    1355              :          iw = cp_print_key_unit_nr(logger, root_section, "GLOBAL%PROGRAM_RUN_INFO", &
    1356         5313 :                                    extension=".log")
    1357              : 
    1358              :          ! Deleting (if existing) the external EXIT files
    1359         5313 :          delete_it = .FALSE.
    1360         5313 :          INQUIRE (FILE="EXIT", EXIST=delete_it)
    1361         5313 :          IF (delete_it) THEN
    1362            0 :             CALL open_file(file_name="EXIT", unit_number=unit_exit)
    1363            0 :             CALL close_file(unit_number=unit_exit, file_status="DELETE")
    1364              :          END IF
    1365              : 
    1366         5313 :          delete_it = .FALSE.
    1367         5313 :          INQUIRE (FILE=TRIM(logger%iter_info%project_name)//".EXIT", EXIST=delete_it)
    1368         5313 :          IF (delete_it) THEN
    1369            0 :             CALL open_file(file_name=TRIM(logger%iter_info%project_name)//".EXIT", unit_number=unit_exit)
    1370            0 :             CALL close_file(unit_number=unit_exit, file_status="DELETE")
    1371              :          END IF
    1372              : 
    1373              :          ! Print OpenMP issue counter and number of warnings for this workload
    1374         5313 :          IF (iw > 0) THEN
    1375         5313 :             IF (0 <= m_omp_trace_issues()) THEN
    1376            0 :                WRITE (iw, "(T2,A,I0)") "The number of traced issues for OpenMP : ", m_omp_trace_issues()
    1377              :             END IF
    1378         5313 :             WRITE (iw, "(T2,A,I0)") "The number of warnings for this run is : ", warning_counter
    1379         5313 :             WRITE (iw, *) ""
    1380         5313 :             WRITE (UNIT=iw, FMT="(T2,A)") REPEAT("-", 79)
    1381              :          END IF
    1382              : 
    1383              :          ! Update the runtime environment variables
    1384         5313 :          CALL get_runtime_info()
    1385              : 
    1386              :          ! Just a choice, do not print the CP2K footer if there is a failure
    1387         5575 :          CALL cp2k_footer(iw, wdir)
    1388         5313 :          IF (iw > 0) FLUSH (iw)  ! ignore &GLOBAL / FLUSH_SHOULD_FLUSH
    1389              : 
    1390              :          CALL cp_print_key_finished_output(iw, logger, root_section, &
    1391         5313 :                                            "GLOBAL%PROGRAM_RUN_INFO")
    1392              :       END IF
    1393              : 
    1394              :       ! Release message passing environment
    1395        10419 :       CALL cp_rm_default_logger()
    1396              : 
    1397        10419 :    END SUBROUTINE cp2k_finalize
    1398              : 
    1399              : END MODULE environment
        

Generated by: LCOV version 2.0-1