(git:e5b1968)
Loading...
Searching...
No Matches
environment.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2025 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! **************************************************************************************************
18 USE bibliography, ONLY: frigo2005,&
19 marek2014,&
20 solca2024,&
21 cite_reference
22 USE cp2k_info, ONLY: &
26 USE cp_files, ONLY: close_file,&
38 diag_init,&
44 USE cp_log_handling, ONLY: &
56 USE fft_tools, ONLY: fwfft,&
57 fft3d,&
63 USE grid_api, ONLY: grid_backend_auto,&
69 USE header, ONLY: cp2k_footer,&
71 USE input_constants, ONLY: &
80 USE input_section_types, ONLY: &
84 USE kinds, ONLY: default_path_length,&
86 dp,&
87 int_8,&
90 USE machine, ONLY: &
95 USE mp_perf_env, ONLY: add_mp_perf_env,&
102 USE parallel_rng_types, ONLY: gaussian,&
103 check_rng,&
106 USE physcon, ONLY: write_physcon
112 USE timings, ONLY: add_timer_env,&
122
123!$ USE OMP_LIB, ONLY: omp_get_max_threads, omp_get_thread_num, omp_get_num_threads
124#include "./base/base_uses.f90"
125
126 IMPLICIT NONE
127
128 PRIVATE
129
130 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'environment'
131
132 ! Public subroutines
133
135
136CONTAINS
137
138! **************************************************************************************************
139!> \brief Initializes a CP2K run (setting of the global environment variables)
140!> \param para_env ...
141!> \param output_unit ...
142!> \param globenv ...
143!> \param input_file_name ...
144!> \param wdir ...
145!> \par History
146!> JGH (28.11.2001) : default for pp_library_path
147!> - print keys added (17.01.2002, MK)
148!> - merged with cp2k_input (30.10.2002,MK)
149!> \author JGH,MK
150! **************************************************************************************************
151 SUBROUTINE cp2k_init(para_env, output_unit, globenv, input_file_name, wdir)
152
153 TYPE(mp_para_env_type), POINTER :: para_env
154 INTEGER :: output_unit
155 TYPE(global_environment_type), POINTER :: globenv
156 CHARACTER(LEN=*) :: input_file_name
157 CHARACTER(LEN=*), OPTIONAL :: wdir
158
159 CHARACTER(LEN=10*default_string_length) :: cp_flags
160 INTEGER :: i, ilen, my_output_unit
161 TYPE(cp_logger_type), POINTER :: logger
162
163 ! create a timer_env
164
165 CALL add_timer_env()
166
167 ! Message passing performance
168 CALL add_mp_perf_env()
169
170 ! Init the default logger
171 IF (para_env%is_source()) THEN
172 my_output_unit = output_unit
173 ELSE
174 my_output_unit = -1
175 END IF
176 NULLIFY (logger)
177 CALL cp_logger_create(logger, para_env=para_env, &
178 default_global_unit_nr=output_unit, &
179 close_global_unit_on_dealloc=.false.)
180 CALL cp_add_default_logger(logger)
181 CALL cp_logger_release(logger)
182
183 ! Initialize timing
184 CALL timeset(root_cp2k_name, globenv%handle)
185
186 ! Print header
187 CALL cp2k_header(my_output_unit, wdir)
188
189 IF (my_output_unit > 0) THEN
190 WRITE (unit=my_output_unit, fmt="(/,T2,A,T31,A50)") &
191 "CP2K| version string: ", adjustr(trim(cp2k_version))
192 WRITE (unit=my_output_unit, fmt="(T2,A,T41,A40)") &
193 "CP2K| source code revision number:", &
194 adjustr(compile_revision)
195 cp_flags = cp2k_flags()
196 ilen = len_trim(cp_flags)
197 WRITE (unit=my_output_unit, fmt="(T2,A)") &
198 "CP2K| "//cp_flags(1:73)
199 IF (ilen > 73) THEN
200 DO i = 0, (ilen - 75)/61
201 WRITE (unit=my_output_unit, fmt="(T2,A)") &
202 "CP2K| "//trim(cp_flags(74 + i*61:min(74 + (i + 1)*61, ilen)))
203 END DO
204 END IF
205 WRITE (unit=my_output_unit, fmt="(T2,A,T41,A40)") &
206 "CP2K| is freely available from ", &
207 adjustr(trim(cp2k_home))
208 WRITE (unit=my_output_unit, fmt="(T2,A,T31,A50)") &
209 "CP2K| Program compiled at", &
210 adjustr(compile_date(1:min(50, len(compile_date))))
211 WRITE (unit=my_output_unit, fmt="(T2,A,T31,A50)") &
212 "CP2K| Program compiled on", &
213 adjustr(compile_host(1:min(50, len(compile_host))))
214 WRITE (unit=my_output_unit, fmt="(T2,A,T31,A50)") &
215 "CP2K| Program compiled for", &
216 adjustr(compile_arch(1:min(50, len(compile_arch))))
217 WRITE (unit=my_output_unit, fmt="(T2,A,T31,A50)") &
218 "CP2K| Data directory path", &
219 adjustr(trim(get_data_dir()))
220 WRITE (unit=my_output_unit, fmt="(T2,A,T31,A50)") &
221 "CP2K| Input file name", &
222 adjustr(trim(input_file_name))
223 FLUSH (my_output_unit) ! ignore &GLOBAL / FLUSH_SHOULD_FLUSH
224 END IF
225
226#if defined(__FAST_MATH__)
227 CALL cp_warn(__location__, &
228 "During compilation one of the following flags was active:"// &
229 " `-ffast-math` (GCC)"// &
230 " `-hfpN` (Cray, N > 0, default N=2)"// &
231 " This can lead to wrong results and numerical instabilities"// &
232 " and is therefore no longer supported.")
233
234#if !defined(__FORCE_USE_FAST_MATH)
235#error "-ffast-math (GCC) or -hfpN (N>0, Cray) can lead to wrong results and numerical instabilities and are therefore no longer supported"
236#endif
237#endif
238
239#if defined(NDEBUG)
240#error "Please do not build CP2K with NDEBUG. There is no performance advantage and asserts will save your neck."
241#endif
242
243 END SUBROUTINE cp2k_init
244
245! **************************************************************************************************
246!> \brief echoes the list of host names and pids
247!> \param para_env ...
248!> \param output_unit ...
249! **************************************************************************************************
250 SUBROUTINE echo_all_hosts(para_env, output_unit)
251 TYPE(mp_para_env_type), POINTER :: para_env
252 INTEGER :: output_unit
253
254 CHARACTER(LEN=default_string_length) :: string
255 INTEGER :: ipe
256 INTEGER, ALLOCATABLE, DIMENSION(:) :: all_pid
257 INTEGER, ALLOCATABLE, DIMENSION(:, :) :: all_host
258
259 ! Print a list of all started processes
260
261 ALLOCATE (all_pid(para_env%num_pe))
262 all_pid(:) = 0
263 all_pid(para_env%mepos + 1) = r_pid
264
265 CALL para_env%sum(all_pid)
266 ALLOCATE (all_host(30, para_env%num_pe))
267 all_host(:, :) = 0
268 CALL string_to_ascii(r_host_name, all_host(:, para_env%mepos + 1))
269 CALL para_env%sum(all_host)
270 IF (output_unit > 0) THEN
271
272 WRITE (unit=output_unit, fmt="(T2,A)") ""
273 DO ipe = 1, para_env%num_pe
274 CALL ascii_to_string(all_host(:, ipe), string)
275 WRITE (unit=output_unit, fmt="(T2,A,T63,I8,T71,I10)") &
276 trim(r_user_name)//"@"//trim(string)// &
277 " has created rank and process ", ipe - 1, all_pid(ipe)
278 END DO
279 WRITE (unit=output_unit, fmt="(T2,A)") ""
280 END IF
281 DEALLOCATE (all_pid)
282 DEALLOCATE (all_host)
283
284 END SUBROUTINE echo_all_hosts
285
286! **************************************************************************************************
287!> \brief echoes the list the number of process per host
288!> \param para_env ...
289!> \param output_unit ...
290! **************************************************************************************************
291 SUBROUTINE echo_all_process_host(para_env, output_unit)
292 TYPE(mp_para_env_type), POINTER :: para_env
293 INTEGER :: output_unit
294
295 CHARACTER(LEN=default_string_length) :: string, string_sec
296 INTEGER :: ipe, jpe, nr_occu
297 INTEGER, ALLOCATABLE, DIMENSION(:) :: all_pid
298 INTEGER, ALLOCATABLE, DIMENSION(:, :) :: all_host
299
300 ALLOCATE (all_host(30, para_env%num_pe))
301 all_host(:, :) = 0
302
303 IF (m_procrun(r_pid) .EQ. 1) THEN
304 CALL string_to_ascii(r_host_name, all_host(:, para_env%mepos + 1))
305 CALL para_env%sum(all_host)
306 END IF
307
308 IF (output_unit > 0) THEN
309 ALLOCATE (all_pid(para_env%num_pe))
310 all_pid(:) = 0
311
312 WRITE (unit=output_unit, fmt="(T2,A)") ""
313 DO ipe = 1, para_env%num_pe
314 nr_occu = 0
315 IF (all_pid(ipe) .NE. -1) THEN
316 CALL ascii_to_string(all_host(:, ipe), string)
317 DO jpe = 1, para_env%num_pe
318 CALL ascii_to_string(all_host(:, jpe), string_sec)
319 IF (string .EQ. string_sec) THEN
320 nr_occu = nr_occu + 1
321 all_pid(jpe) = -1
322 END IF
323 END DO
324 WRITE (unit=output_unit, fmt="(T2,A,T63,I8,A)") &
325 trim(r_user_name)//"@"//trim(string)// &
326 " is running ", nr_occu, " processes"
327 WRITE (unit=output_unit, fmt="(T2,A)") ""
328 END IF
329 END DO
330 DEALLOCATE (all_pid)
331
332 END IF
333
334 DEALLOCATE (all_host)
335
336 END SUBROUTINE echo_all_process_host
337
338! **************************************************************************************************
339!> \brief read part of cp2k_init
340!> \param root_section ...
341!> \param para_env ...
342!> \param globenv the globenv
343!> \author fawzi
344!> \note
345!> The following routines need to be synchronized wrt. adding/removing
346!> of the default environments (logging, performance,error):
347!> environment:cp2k_init, environment:cp2k_finalize,
348!> f77_interface:f_env_add_defaults, f77_interface:f_env_rm_defaults,
349!> f77_interface:create_force_env, f77_interface:destroy_force_env
350! **************************************************************************************************
351 SUBROUTINE cp2k_read(root_section, para_env, globenv)
352
353 TYPE(section_vals_type), POINTER :: root_section
354 TYPE(mp_para_env_type), POINTER :: para_env
355 TYPE(global_environment_type), POINTER :: globenv
356
357 CHARACTER(LEN=3*default_string_length) :: message
358 CHARACTER(LEN=default_string_length) :: c_val
359 INTEGER :: i, iw
360 TYPE(cp_logger_type), POINTER :: logger
361
362 ! Read the input/output section
363
364 logger => cp_get_default_logger()
365
366 ! try to use better names for the local log if it is not too late
367 CALL section_vals_val_get(root_section, "GLOBAL%OUTPUT_FILE_NAME", &
368 c_val=c_val)
369 IF (c_val /= "") THEN
370 CALL cp_logger_set(logger, &
371 local_filename=trim(c_val)//"_localLog")
372 END IF
373
374 ! Process project name
375 CALL section_vals_val_get(root_section, "GLOBAL%PROJECT", c_val=c_val)
376 IF (index(c_val(:len_trim(c_val)), " ") > 0) THEN
377 message = "Project name <"//trim(c_val)// &
378 "> contains spaces which will be replaced with underscores"
379 cpwarn(trim(message))
380 DO i = 1, len_trim(c_val)
381 ! Replace space with underscore
382 IF (c_val(i:i) == " ") c_val(i:i) = "_"
383 END DO
384 CALL section_vals_val_set(root_section, "GLOBAL%PROJECT", c_val=trim(c_val))
385 END IF
386 IF (c_val /= "") THEN
387 CALL cp_logger_set(logger, local_filename=trim(c_val)//"_localLog")
388 END IF
389 logger%iter_info%project_name = c_val
390
391 CALL section_vals_val_get(root_section, "GLOBAL%PRINT_LEVEL", i_val=logger%iter_info%print_level)
392
393 ! Read the CP2K section
394 CALL read_cp2k_section(root_section, para_env, globenv)
395
396 iw = cp_print_key_unit_nr(logger, root_section, "GLOBAL%PRINT/BASIC_DATA_TYPES", &
397 extension=".Log")
398 IF (iw > 0) CALL print_kind_info(iw)
399 CALL cp_print_key_finished_output(iw, logger, root_section, &
400 "GLOBAL%PRINT/BASIC_DATA_TYPES")
401
402 iw = cp_print_key_unit_nr(logger, root_section, "GLOBAL%PRINT/PHYSCON", &
403 extension=".Log")
404 IF (iw > 0) CALL write_physcon(iw)
405 CALL cp_print_key_finished_output(iw, logger, root_section, &
406 "GLOBAL%PRINT/PHYSCON")
407
408 END SUBROUTINE cp2k_read
409
410! **************************************************************************************************
411!> \brief globenv initializations that need the input and error
412!> \param root_section ...
413!> \param para_env ...
414!> \param globenv the global environment to initialize
415!> \author fawzi
416!> \note
417!> if possible do the initializations here as the environment
418!> (error,...) is setup, instead of cp2k_init
419! **************************************************************************************************
420 SUBROUTINE cp2k_setup(root_section, para_env, globenv)
421
422 TYPE(section_vals_type), POINTER :: root_section
423 TYPE(mp_para_env_type), POINTER :: para_env
424 TYPE(global_environment_type), POINTER :: globenv
425
426 INTEGER :: iw, maxl
427 INTEGER, DIMENSION(:), POINTER :: seed_vals
428 REAL(kind=dp), DIMENSION(3, 2) :: initial_seed
429 TYPE(cp_logger_type), POINTER :: logger
430
431 NULLIFY (logger)
432 logger => cp_get_default_logger()
433
434 ! Initialize the parallel random number generator
435
436 iw = cp_print_key_unit_nr(logger, root_section, "GLOBAL%PRINT/RNG_MATRICES", &
437 extension=".Log")
438 IF (iw > 0) THEN
439 CALL write_rng_matrices(iw)
440 END IF
441
442 CALL cp_print_key_finished_output(iw, logger, root_section, &
443 "GLOBAL%PRINT/RNG_MATRICES")
444
445 ! Initialize a global normally Gaussian distributed (pseudo)random number stream
446
447 CALL section_vals_val_get(root_section, "GLOBAL%SEED", i_vals=seed_vals)
448 IF (SIZE(seed_vals) == 1) THEN
449 initial_seed(:, :) = real(seed_vals(1), kind=dp)
450 ELSE IF (SIZE(seed_vals) == 6) THEN
451 initial_seed(1:3, 1:2) = reshape(real(seed_vals(:), kind=dp), (/3, 2/))
452 ELSE
453 cpabort("Supply exactly 1 or 6 arguments for SEED in &GLOBAL only!")
454 END IF
455
456 globenv%gaussian_rng_stream = rng_stream_type( &
457 name="Global Gaussian random numbers", &
458 distribution_type=gaussian, &
459 seed=initial_seed, &
460 extended_precision=.true.)
461
462 iw = cp_print_key_unit_nr(logger, root_section, "GLOBAL%PRINT/RNG_CHECK", &
463 extension=".Log")
464 IF (iw > 0) THEN
465 CALL check_rng(iw, para_env%is_source())
466 END IF
467
468 CALL cp_print_key_finished_output(iw, logger, root_section, &
469 "GLOBAL%PRINT/RNG_CHECK")
470
471 iw = cp_print_key_unit_nr(logger, root_section, "GLOBAL%PRINT/GLOBAL_GAUSSIAN_RNG", &
472 extension=".Log")
473 IF (iw > 0) &
474 CALL globenv%gaussian_rng_stream%write(iw, write_all=.true.)
475
476 CALL cp_print_key_finished_output(iw, logger, root_section, &
477 "GLOBAL%PRINT/GLOBAL_GAUSSIAN_RNG")
478
479 CALL section_vals_val_get(root_section, "GLOBAL%PRINT%SPHERICAL_HARMONICS", i_val=maxl)
480 IF (maxl >= 0) THEN
481 iw = cp_print_key_unit_nr(logger, root_section, "GLOBAL%PRINT", &
482 extension=".Log")
483 CALL init_orbital_pointers(maxl)
484 CALL init_spherical_harmonics(maxl, iw)
487 CALL cp_print_key_finished_output(iw, logger, root_section, &
488 "GLOBAL%PRINT")
489 END IF
490
491 END SUBROUTINE cp2k_setup
492
493! **************************************************************************************************
494!> \brief read the global section of new input
495!> \param root_section ...
496!> \param para_env ...
497!> \param globenv ...
498!> \par History
499!> 06-2005 [created]
500!> \author MI
501!> \note
502!> Should not be required anymore once everything is converted
503!> to get information directly from the input structure
504! **************************************************************************************************
505 SUBROUTINE read_global_section(root_section, para_env, globenv)
506
507 TYPE(section_vals_type), POINTER :: root_section
508 TYPE(mp_para_env_type), POINTER :: para_env
509 TYPE(global_environment_type), POINTER :: globenv
510
511 CHARACTER(LEN=6), PARAMETER :: start_section_label = "GLOBAL"
512
513 CHARACTER(LEN=13) :: omp_stacksize, tracing_string
514 CHARACTER(LEN=6) :: print_level_string
515 CHARACTER(LEN=default_path_length) :: basis_set_file_name, coord_file_name, &
516 mm_potential_file_name, &
517 potential_file_name
518 CHARACTER(LEN=default_string_length) :: env_num, model_name, project_name
519 CHARACTER(LEN=default_string_length), &
520 DIMENSION(:), POINTER :: trace_routines
521 INTEGER :: cpuid, cpuid_static, i_cholesky, i_dgemm, i_diag, i_fft, i_grid_backend, &
522 iforce_eval, method_name_id, n_rep_val, nforce_eval, num_threads, output_unit, &
523 print_level, trace_max, unit_nr
524 INTEGER(kind=int_8) :: buffers, buffers_avr, buffers_max, buffers_min, cached, cached_avr, &
525 cached_max, cached_min, memfree, memfree_avr, memfree_max, memfree_min, memlikelyfree, &
526 memlikelyfree_avr, memlikelyfree_max, memlikelyfree_min, memtotal, memtotal_avr, &
527 memtotal_max, memtotal_min, slab, slab_avr, slab_max, slab_min, sreclaimable, &
528 sreclaimable_avr, sreclaimable_max, sreclaimable_min
529 INTEGER, DIMENSION(:), POINTER :: i_force_eval
530 LOGICAL :: ata, do_echo_all_hosts, efl, explicit, &
531 flag, report_maxloc, trace, &
532 trace_master
533 TYPE(cp_logger_type), POINTER :: logger
534 TYPE(enumeration_type), POINTER :: enum1, enum2
535 TYPE(keyword_type), POINTER :: keyword
536 TYPE(section_type), POINTER :: section
537 TYPE(section_vals_type), POINTER :: dft_section, force_env_sections, &
538 global_section, qmmm_section, &
539 subsys_section
540
541 NULLIFY (dft_section, global_section, i_force_eval)
542
543 logger => cp_get_default_logger()
544 global_section => section_vals_get_subs_vals(root_section, "GLOBAL")
545 CALL section_vals_val_get(global_section, "BLACS_GRID", i_val=globenv%blacs_grid_layout)
546 CALL section_vals_val_get(global_section, "BLACS_REPEATABLE", l_val=globenv%blacs_repeatable)
547 CALL section_vals_val_get(global_section, "PREFERRED_DIAG_LIBRARY", i_val=i_diag)
548 CALL section_vals_val_get(global_section, "PREFERRED_CHOLESKY_LIBRARY", i_val=i_cholesky)
549 CALL section_vals_val_get(global_section, "PREFERRED_DGEMM_LIBRARY", i_val=i_dgemm)
550 CALL section_vals_val_get(global_section, "EPS_CHECK_DIAG", r_val=globenv%eps_check_diag)
551 CALL section_vals_val_get(global_section, "ENABLE_MPI_IO", l_val=flag)
552 CALL cp_mpi_io_set(flag)
553 CALL section_vals_val_get(global_section, "ELPA_KERNEL", i_val=globenv%k_elpa)
554 CALL section_vals_val_get(global_section, "ELPA_NEIGVEC_MIN", i_val=globenv%elpa_neigvec_min)
555 CALL section_vals_val_get(global_section, "ELPA_QR", l_val=globenv%elpa_qr)
556 CALL section_vals_val_get(global_section, "ELPA_QR_UNSAFE", l_val=globenv%elpa_qr_unsafe)
557 unit_nr = cp_print_key_unit_nr(logger, global_section, "PRINT_ELPA", extension=".Log")
558 IF (unit_nr > 0) globenv%elpa_print = .true.
559 CALL cp_print_key_finished_output(unit_nr, logger, global_section, "PRINT_ELPA")
560 CALL section_vals_val_get(global_section, "DLAF_NEIGVEC_MIN", i_val=globenv%dlaf_neigvec_min)
561 CALL section_vals_val_get(global_section, "DLAF_CHOLESKY_N_MIN", i_val=globenv%dlaf_cholesky_n_min)
562 CALL section_vals_val_get(global_section, "PREFERRED_FFT_LIBRARY", i_val=i_fft)
563 CALL section_vals_val_get(global_section, "PRINT_LEVEL", i_val=print_level)
564 CALL section_vals_val_get(global_section, "PROGRAM_NAME", i_val=globenv%prog_name_id)
565 CALL section_vals_val_get(global_section, "FFT_POOL_SCRATCH_LIMIT", i_val=globenv%fft_pool_scratch_limit)
566 CALL section_vals_val_get(global_section, "FFTW_PLAN_TYPE", i_val=globenv%fftw_plan_type)
567 CALL section_vals_val_get(global_section, "PROJECT_NAME", c_val=project_name)
568 CALL section_vals_val_get(global_section, "FFTW_WISDOM_FILE_NAME", c_val=globenv%fftw_wisdom_file_name)
569 CALL section_vals_val_get(global_section, "RUN_TYPE", i_val=globenv%run_type_id)
570 CALL cp2k_get_walltime(section=global_section, keyword_name="WALLTIME", &
571 walltime=globenv%cp2k_target_time)
572 CALL section_vals_val_get(global_section, "TRACE", l_val=trace)
573 CALL section_vals_val_get(global_section, "TRACE_MASTER", l_val=trace_master)
574 CALL section_vals_val_get(global_section, "TRACE_MAX", i_val=trace_max)
575 CALL section_vals_val_get(global_section, "TRACE_ROUTINES", explicit=explicit)
576 IF (explicit) THEN
577 CALL section_vals_val_get(global_section, "TRACE_ROUTINES", c_vals=trace_routines)
578 ELSE
579 NULLIFY (trace_routines)
580 END IF
581 CALL section_vals_val_get(global_section, "FLUSH_SHOULD_FLUSH", l_val=flush_should_flush)
582 CALL section_vals_val_get(global_section, "ECHO_ALL_HOSTS", l_val=do_echo_all_hosts)
583 report_maxloc = section_get_lval(global_section, "TIMINGS%REPORT_MAXLOC")
584 global_timings_level = section_get_ival(global_section, "TIMINGS%TIMINGS_LEVEL")
585 do_echo_all_hosts = do_echo_all_hosts .OR. report_maxloc
586 force_env_sections => section_vals_get_subs_vals(root_section, "FORCE_EVAL")
587 CALL section_vals_get(force_env_sections, n_repetition=nforce_eval)
588 output_unit = cp_print_key_unit_nr(logger, global_section, "PROGRAM_RUN_INFO", &
589 extension=".log")
590
591 CALL fm_setup(global_section)
592 CALL fm_diag_rules_setup(global_section)
593 CALL dgemm_setup(global_section)
594
595 IF (trace .AND. (.NOT. trace_master .OR. para_env%mepos == 0)) THEN
596 unit_nr = -1
597 IF (logger%para_env%is_source() .OR. .NOT. trace_master) &
598 unit_nr = cp_logger_get_default_unit_nr(logger, local=.true.)
599 WRITE (tracing_string, "(I6.6,A1,I6.6)") para_env%mepos, ":", para_env%num_pe
600 IF (ASSOCIATED(trace_routines)) THEN
601 CALL timings_setup_tracing(trace_max, unit_nr, tracing_string, trace_routines)
602 ELSE
603 CALL timings_setup_tracing(trace_max, unit_nr, tracing_string)
604 END IF
605 END IF
606
607 CALL section_vals_val_get(global_section, "TIMINGS%TIME_MPI", l_val=mp_collect_timings)
608
609 SELECT CASE (i_diag)
611 globenv%diag_library = "ScaLAPACK"
612 CASE (fm_diag_type_elpa)
613 globenv%diag_library = "ELPA"
616 globenv%diag_library = "cuSOLVER"
617 CASE (fm_diag_type_dlaf)
618 globenv%diag_library = "DLAF"
620 CASE DEFAULT
621 cpabort("Unknown diagonalization library specified")
622 END SELECT
623
624 SELECT CASE (i_cholesky)
626 globenv%cholesky_library = "ScaLAPACK"
629 globenv%cholesky_library = "DLAF"
631 dlaf_cholesky_n_min = globenv%dlaf_cholesky_n_min
633 CASE DEFAULT
634 cpabort("Unknown Cholesky decomposition library specified")
635 END SELECT
636
637 SELECT CASE (i_fft)
638 CASE (do_fft_sg)
639 globenv%default_fft_library = "FFTSG"
640 CASE (do_fft_fftw3)
641 globenv%default_fft_library = "FFTW3"
643 CASE DEFAULT
644 cpabort("Unknown FFT library specified")
645 END SELECT
646
647 SELECT CASE (i_dgemm)
648 CASE (do_dgemm_spla)
649 globenv%default_dgemm_library = "SPLA"
650 CASE (do_dgemm_blas)
651 globenv%default_dgemm_library = "BLAS"
652 CASE DEFAULT
653 cpabort("Unknown DGEMM library specified")
654 END SELECT
655
656 IF (globenv%run_type_id == 0) THEN
657 SELECT CASE (globenv%prog_name_id)
658 CASE (do_farming, do_test)
659 globenv%run_type_id = none_run
660 CASE (do_cp2k)
661 IF (nforce_eval /= 1) THEN
662 ! multiple force_eval corresponds at the moment to RESPA calculations only
663 ! default MD
664 globenv%run_type_id = mol_dyn_run
665 ELSE
666 CALL section_vals_val_get(force_env_sections, "METHOD", i_val=method_name_id)
667 SELECT CASE (method_name_id)
668 CASE (do_fist)
669 globenv%run_type_id = mol_dyn_run
670 CASE (do_eip)
671 globenv%run_type_id = mol_dyn_run
672 CASE (do_qs)
673 globenv%run_type_id = energy_run
674 CASE (do_sirius)
675 globenv%run_type_id = energy_run
676 END SELECT
677 END IF
678 END SELECT
679 END IF
680
681 IF (globenv%prog_name_id == do_farming .AND. globenv%run_type_id /= none_run) THEN
682 cpabort("FARMING program supports only NONE as run type")
683 END IF
684
685 IF (globenv%prog_name_id == do_test .AND. globenv%run_type_id /= none_run) &
686 cpabort("TEST program supports only NONE as run type")
687
688 CALL m_memory_details(memtotal, memfree, buffers, cached, slab, sreclaimable, memlikelyfree)
689 memtotal_avr = memtotal
690 memfree_avr = memfree
691 buffers_avr = buffers
692 cached_avr = cached
693 slab_avr = slab
694 sreclaimable_avr = sreclaimable
695 memlikelyfree_avr = memlikelyfree
696 CALL para_env%sum(memtotal_avr); memtotal_avr = memtotal_avr/para_env%num_pe/1024
697 CALL para_env%sum(memfree_avr); memfree_avr = memfree_avr/para_env%num_pe/1024
698 CALL para_env%sum(buffers_avr); buffers_avr = buffers_avr/para_env%num_pe/1024
699 CALL para_env%sum(cached_avr); cached_avr = cached_avr/para_env%num_pe/1024
700 CALL para_env%sum(slab_avr); slab_avr = slab_avr/para_env%num_pe/1024
701 CALL para_env%sum(sreclaimable_avr); sreclaimable_avr = sreclaimable_avr/para_env%num_pe/1024
702 CALL para_env%sum(memlikelyfree_avr); memlikelyfree_avr = memlikelyfree_avr/para_env%num_pe/1024
703
704 memtotal_min = -memtotal
705 memfree_min = -memfree
706 buffers_min = -buffers
707 cached_min = -cached
708 slab_min = -slab
709 sreclaimable_min = -sreclaimable
710 memlikelyfree_min = -memlikelyfree
711 CALL para_env%max(memtotal_min); memtotal_min = -memtotal_min/1024
712 CALL para_env%max(memfree_min); memfree_min = -memfree_min/1024
713 CALL para_env%max(buffers_min); buffers_min = -buffers_min/1024
714 CALL para_env%max(cached_min); cached_min = -cached_min/1024
715 CALL para_env%max(slab_min); slab_min = -slab_min/1024
716 CALL para_env%max(sreclaimable_min); sreclaimable_min = -sreclaimable_min/1024
717 CALL para_env%max(memlikelyfree_min); memlikelyfree_min = -memlikelyfree_min/1024
718
719 memtotal_max = memtotal
720 memfree_max = memfree
721 buffers_max = buffers
722 cached_max = cached
723 slab_max = slab
724 sreclaimable_max = sreclaimable
725 memlikelyfree_max = memlikelyfree
726 CALL para_env%max(memtotal_max); memtotal_max = memtotal_max/1024
727 CALL para_env%max(memfree_max); memfree_max = memfree_max/1024
728 CALL para_env%max(buffers_max); buffers_max = buffers_max/1024
729 CALL para_env%max(cached_max); cached_max = cached_max/1024
730 CALL para_env%max(slab_max); slab_max = slab_max/1024
731 CALL para_env%max(sreclaimable_max); sreclaimable_max = sreclaimable_max/1024
732 CALL para_env%max(memlikelyfree_max); memlikelyfree_max = memlikelyfree_max/1024
733
734 memtotal = memtotal/1024
735 memfree = memfree/1024
736 buffers = buffers/1024
737 cached = cached/1024
738 slab = slab/1024
739 sreclaimable = sreclaimable/1024
740 memlikelyfree = memlikelyfree/1024
741
742 ! Print a list of all started processes
743 IF (do_echo_all_hosts) THEN
744 CALL echo_all_hosts(para_env, output_unit)
745
746 ! Print the number of processes per host
747 CALL echo_all_process_host(para_env, output_unit)
748 END IF
749
750 num_threads = 1
751!$ num_threads = omp_get_max_threads()
752 IF (output_unit > 0) THEN
753 WRITE (unit=output_unit, fmt=*)
754 CALL multiple_fe_list(force_env_sections, root_section, i_force_eval, nforce_eval)
755 DO iforce_eval = 1, nforce_eval
756 dft_section => section_vals_get_subs_vals3(force_env_sections, "DFT", &
757 i_rep_section=i_force_eval(iforce_eval))
758 qmmm_section => section_vals_get_subs_vals3(force_env_sections, "QMMM", &
759 i_rep_section=i_force_eval(iforce_eval))
760 CALL section_vals_val_get(dft_section, "BASIS_SET_FILE_NAME", &
761 c_val=basis_set_file_name)
762 CALL section_vals_val_get(dft_section, "POTENTIAL_FILE_NAME", &
763 c_val=potential_file_name)
764
765 CALL section_vals_val_get(qmmm_section, "MM_POTENTIAL_FILE_NAME", &
766 c_val=mm_potential_file_name)
767 ! SUBSYS - If any
768 subsys_section => section_vals_get_subs_vals3(force_env_sections, "SUBSYS", &
769 i_rep_section=i_force_eval(iforce_eval))
770 CALL section_vals_get(subsys_section, explicit=explicit)
771 coord_file_name = "__STD_INPUT__"
772 IF (explicit) THEN
773 CALL section_vals_val_get(subsys_section, "TOPOLOGY%COORD_FILE_NAME", &
774 n_rep_val=n_rep_val)
775 IF (n_rep_val == 1) THEN
776 CALL section_vals_val_get(subsys_section, "TOPOLOGY%COORD_FILE_NAME", &
777 c_val=coord_file_name)
778 END IF
779 END IF
780 CALL integer_to_string(i_force_eval(iforce_eval), env_num)
781
782 WRITE (unit=output_unit, fmt="(T2,A,T41,A)") &
783 start_section_label//"| Force Environment number", &
784 adjustr(env_num(:40)), &
785 start_section_label//"| Basis set file name", &
786 adjustr(basis_set_file_name(:40)), &
787 start_section_label//"| Potential file name", &
788 adjustr(potential_file_name(:40)), &
789 start_section_label//"| MM Potential file name", &
790 adjustr(mm_potential_file_name(:40)), &
791 start_section_label//"| Coordinate file name", &
792 adjustr(coord_file_name(:40))
793 END DO
794 DEALLOCATE (i_force_eval)
795
796 NULLIFY (enum1, enum2, keyword, section)
797 CALL create_global_section(section)
798 keyword => section_get_keyword(section, "PROGRAM_NAME")
799 CALL keyword_get(keyword, enum=enum1)
800 keyword => section_get_keyword(section, "RUN_TYPE")
801 CALL keyword_get(keyword, enum=enum2)
802
803 WRITE (unit=output_unit, fmt="(T2,A,T41,A40)") &
804 start_section_label//"| Method name", &
805 adjustr(trim(enum_i2c(enum1, globenv%prog_name_id))), &
806 start_section_label//"| Project name", &
807 adjustr(project_name(:40)), &
808 start_section_label//"| Run type", &
809 adjustr(trim(enum_i2c(enum2, globenv%run_type_id))), &
810 start_section_label//"| FFT library", &
811 adjustr(globenv%default_fft_library(:40)), &
812 start_section_label//"| Diagonalization library", &
813 adjustr(globenv%diag_library(:40)), &
814 start_section_label//"| Cholesky decomposition library", &
815 adjustr(globenv%cholesky_library(:40)), &
816 start_section_label//"| DGEMM library", &
817 adjustr(globenv%default_dgemm_library(:40))
818
819 IF (globenv%diag_library == "ELPA") THEN
820 WRITE (unit=output_unit, fmt="(T2,A,T71,I10)") &
821 start_section_label//"| Minimum number of eigenvectors for ELPA usage", &
822 globenv%elpa_neigvec_min
823 END IF
824
825 IF (globenv%diag_library == "DLAF") THEN
826 WRITE (unit=output_unit, fmt="(T2,A,T71,I10)") &
827 start_section_label//"| Minimum number of eigenvectors for DLAF usage", &
828 globenv%dlaf_neigvec_min
829 END IF
830
831 IF (globenv%cholesky_library == "DLAF") THEN
832 WRITE (unit=output_unit, fmt="(T2,A,T71,I10)") &
833 start_section_label//"| Minimum matrix size for Cholesky decomposition with DLAF", &
834 globenv%dlaf_cholesky_n_min
835 END IF
836
837#if defined(__CHECK_DIAG)
838 ! Perform default check if no threshold value has been specified explicitly
839 IF (globenv%eps_check_diag < 0.0_dp) THEN
840 WRITE (unit=output_unit, fmt="(T2,A,T71,ES10.3)") &
841 start_section_label//"| Orthonormality check for eigenvectors enabled", &
843 ELSE
844 WRITE (unit=output_unit, fmt="(T2,A,T71,ES10.3)") &
845 start_section_label//"| Orthonormality check for eigenvectors enabled", &
846 globenv%eps_check_diag
847 END IF
848#else
849 IF (globenv%eps_check_diag < 0.0_dp) THEN
850 WRITE (unit=output_unit, fmt="(T2,A,T73,A)") &
851 start_section_label//"| Orthonormality check for eigenvectors", &
852 "DISABLED"
853 ELSE
854 WRITE (unit=output_unit, fmt="(T2,A,T71,ES10.3)") &
855 start_section_label//"| Orthonormality check for eigenvectors enabled", &
856 globenv%eps_check_diag
857 END IF
858#endif
859 CALL section_release(section)
860
861 SELECT CASE (cp_fm_get_mm_type())
862 CASE (do_scalapack)
863 WRITE (unit=output_unit, fmt="(T2,A,T72,A)") &
864 start_section_label//"| Matrix multiplication library", "ScaLAPACK"
865 CASE (do_cosma)
866 WRITE (unit=output_unit, fmt="(T2,A,T76,A)") &
867 start_section_label//"| Matrix multiplication library", "COSMA"
868 END SELECT
869
870 CALL section_vals_val_get(global_section, "ALLTOALL_SGL", l_val=ata)
871 WRITE (unit=output_unit, fmt="(T2,A,T80,L1)") &
872 start_section_label//"| All-to-all communication in single precision", ata
873 CALL section_vals_val_get(global_section, "EXTENDED_FFT_LENGTHS", l_val=efl)
874 WRITE (unit=output_unit, fmt="(T2,A,T80,L1)") &
875 start_section_label//"| FFTs using library dependent lengths", efl
876
877 SELECT CASE (print_level)
878 CASE (silent_print_level)
879 print_level_string = "SILENT"
880 CASE (low_print_level)
881 print_level_string = " LOW"
882 CASE (medium_print_level)
883 print_level_string = "MEDIUM"
884 CASE (high_print_level)
885 print_level_string = " HIGH"
886 CASE (debug_print_level)
887 print_level_string = " DEBUG"
888 CASE DEFAULT
889 cpabort("Unknown print_level")
890 END SELECT
891
892 CALL section_vals_val_get(global_section, "GRID%BACKEND", i_val=i_grid_backend)
893 SELECT CASE (i_grid_backend)
894 CASE (grid_backend_auto)
895 WRITE (unit=output_unit, fmt="(T2,A,T75,A6)") &
896 start_section_label//"| Grid backend", "AUTO"
897 CASE (grid_backend_cpu)
898 WRITE (unit=output_unit, fmt="(T2,A,T75,A6)") &
899 start_section_label//"| Grid backend", "CPU"
900 CASE (grid_backend_dgemm)
901 WRITE (unit=output_unit, fmt="(T2,A,T75,A6)") &
902 start_section_label//"| Grid backend", "DGEMM"
903 CASE (grid_backend_gpu)
904 WRITE (unit=output_unit, fmt="(T2,A,T75,A6)") &
905 start_section_label//"| Grid backend", "GPU"
906 CASE (grid_backend_hip)
907 WRITE (unit=output_unit, fmt="(T2,A,T75,A6)") &
908 start_section_label//"| Grid backend", "HIP"
909 CASE (grid_backend_ref)
910 WRITE (unit=output_unit, fmt="(T2,A,T75,A6)") &
911 start_section_label//"| Grid backend", "REF"
912 END SELECT
913
914 WRITE (unit=output_unit, fmt="(T2,A,T75,A6)") &
915 start_section_label//"| Global print level", print_level_string
916 WRITE (unit=output_unit, fmt="(T2,A,T75,L6)") &
917 start_section_label//"| MPI I/O enabled", flag
918 WRITE (unit=output_unit, fmt="(T2,A,T75,I6)") &
919 start_section_label//"| Total number of message passing processes", &
920 para_env%num_pe, &
921 start_section_label//"| Number of threads for this process", &
922 num_threads, &
923 start_section_label//"| This output is from process", para_env%mepos
924
925 CALL m_omp_get_stacksize(omp_stacksize)
926 WRITE (unit=output_unit, fmt="(T2,A,T68,A13)") &
927 start_section_label//"| OpenMP stack size per thread (OMP_STACKSIZE)", &
928 adjustr(omp_stacksize)
929
930 IF (0 <= m_omp_trace_issues()) THEN ! only show in header if enabled
931 WRITE (unit=output_unit, fmt="(T2,A,T68,A13)") &
932 start_section_label//"| OpenMP issue trace (CP2K_OMP_TRACE)", &
933 "enabled"
934 END IF
935
936 CALL m_cpuinfo(model_name)
937 WRITE (unit=output_unit, fmt="(T2,A,T30,A51)") &
938 start_section_label//"| CPU model name", adjustr(trim(model_name))
939
940 cpuid = m_cpuid()
941 cpuid_static = m_cpuid_static()
942
943 IF ((cpuid > 0) .OR. (cpuid_static > 0)) THEN
944 WRITE (unit=output_unit, fmt="(T2,A,T75,I6)") &
945 start_section_label//"| CPUID", cpuid
946 IF (cpuid /= cpuid_static) THEN
947 WRITE (unit=output_unit, fmt="(T2,A,T75,I6)") &
948 start_section_label//"| Compiled for CPUID", cpuid_static
949 END IF
950 END IF
951
952 ! filter cpuids by vlen to show more relevant information
953 IF (m_cpuid_vlen(cpuid_static) < m_cpuid_vlen(cpuid)) THEN
954 ! base/machine_cpuid.c relies on the (same) target flags as the Fortran code
955 CALL cp_hint(__location__, "The compiler target flags ("// &
956 trim(m_cpuid_name(cpuid_static))//") used to build this binary cannot exploit "// &
957 "all extensions of this CPU model ("//trim(m_cpuid_name(cpuid))//"). "// &
958 "Consider compiler target flags as part of FCFLAGS and CFLAGS (ARCH file).")
959 END IF
960
961 WRITE (unit=output_unit, fmt="()")
962 WRITE (unit=output_unit, fmt="(T2,A)") "MEMORY| system memory details [Kb]"
963 WRITE (unit=output_unit, fmt="(T2,A23,4A14)") "MEMORY| ", "rank 0", "min", "max", "average"
964 WRITE (unit=output_unit, fmt="(T2,A23,4I14)") "MEMORY| MemTotal ", memtotal, memtotal_min, memtotal_max, memtotal_avr
965 WRITE (unit=output_unit, fmt="(T2,A23,4I14)") "MEMORY| MemFree ", memfree, memfree_min, memfree_max, memfree_avr
966 WRITE (unit=output_unit, fmt="(T2,A23,4I14)") "MEMORY| Buffers ", buffers, buffers_min, buffers_max, buffers_avr
967 WRITE (unit=output_unit, fmt="(T2,A23,4I14)") "MEMORY| Cached ", cached, cached_min, cached_max, cached_avr
968 WRITE (unit=output_unit, fmt="(T2,A23,4I14)") "MEMORY| Slab ", slab, slab_min, slab_max, slab_avr
969 WRITE (unit=output_unit, fmt="(T2,A23,4I14)") &
970 "MEMORY| SReclaimable ", sreclaimable, sreclaimable_min, sreclaimable_max, &
971 sreclaimable_avr
972 WRITE (unit=output_unit, fmt="(T2,A23,4I14)") &
973 "MEMORY| MemLikelyFree ", memlikelyfree, memlikelyfree_min, memlikelyfree_max, &
974 memlikelyfree_avr
975 WRITE (unit=output_unit, fmt='()')
976
977 END IF
978
979 CALL cp_print_key_finished_output(output_unit, logger, global_section, &
980 "PROGRAM_RUN_INFO")
981
982 END SUBROUTINE read_global_section
983
984! **************************************************************************************************
985!> \brief ...
986!> \param root_section ...
987!> \param para_env ...
988!> \param globenv ...
989!> \par History
990!> 2-Dec-2000 (JGH) added default fft library
991!> \author JGH,MK
992! **************************************************************************************************
993 SUBROUTINE read_cp2k_section(root_section, para_env, globenv)
994
995 TYPE(section_vals_type), POINTER :: root_section
996 TYPE(mp_para_env_type), POINTER :: para_env
997 TYPE(global_environment_type), POINTER :: globenv
998
999 INTEGER :: output_unit
1000 TYPE(cp_logger_type), POINTER :: logger
1001 TYPE(section_vals_type), POINTER :: global_section
1002
1003 global_section => section_vals_get_subs_vals(root_section, "GLOBAL")
1004 CALL read_global_section(root_section, para_env, globenv)
1005 logger => cp_get_default_logger()
1006 output_unit = cp_print_key_unit_nr(logger, global_section, "PROGRAM_RUN_INFO", &
1007 extension=".log")
1008
1009 CALL fft_setup_library(globenv, global_section)
1010 CALL diag_setup_library(globenv)
1011
1012 CALL cp_print_key_finished_output(output_unit, logger, global_section, &
1013 "PROGRAM_RUN_INFO")
1014
1015 END SUBROUTINE read_cp2k_section
1016
1017! **************************************************************************************************
1018!> \brief check FFT preferred library availability, if not switch
1019!> \param globenv ...
1020!> \param global_section ...
1021!> \par History
1022!> 2-Dec-2000 (JGH) added default fft library
1023!> Nov-2013 (MI) refactoring
1024!> \author JGH,MK
1025! **************************************************************************************************
1026 SUBROUTINE fft_setup_library(globenv, global_section)
1027
1028 TYPE(global_environment_type), POINTER :: globenv
1029 TYPE(section_vals_type), POINTER :: global_section
1030
1031 CHARACTER(LEN=3*default_string_length) :: message
1032 COMPLEX(KIND=dp), DIMENSION(4, 4, 4) :: zz
1033 INTEGER :: stat
1034 INTEGER, DIMENSION(3) :: n
1035 LOGICAL :: try_fftw
1036
1037 n(:) = 4
1038 zz(:, :, :) = 0.0_dp
1039
1040 ! Setup the FFT library
1041 ! If the user has specified PREFERRED_FFT_LIBRARY try that first (default FFTW3)
1042 ! If that one is not available, try FFTW3 (unless it has been tried already)
1043 ! If FFTW3 is not available use FFTSG
1044
1045 IF (globenv%default_fft_library == "FFTW3") THEN
1046 try_fftw = .false.
1047 ELSE
1048 try_fftw = .true.
1049 END IF
1050
1051 ! Initialize FFT library with the user's preferred FFT library
1052 CALL init_fft(fftlib=trim(globenv%default_fft_library), &
1053 alltoall=section_get_lval(global_section, "ALLTOALL_SGL"), &
1054 fftsg_sizes=.NOT. section_get_lval(global_section, "EXTENDED_FFT_LENGTHS"), &
1055 pool_limit=globenv%fft_pool_scratch_limit, &
1056 wisdom_file=globenv%fftw_wisdom_file_name, &
1057 plan_style=globenv%fftw_plan_type)
1058
1059 ! Check for FFT library
1060 CALL fft3d(fwfft, n, zz, status=stat)
1061 IF (stat /= 0) THEN
1062 IF (try_fftw) THEN
1063 message = "FFT library "//trim(globenv%default_fft_library)// &
1064 " is not available. Trying FFT library FFTW3."
1065 cpwarn(trim(message))
1066 globenv%default_fft_library = "FFTW3"
1067 CALL init_fft(fftlib=trim(globenv%default_fft_library), &
1068 alltoall=section_get_lval(global_section, "ALLTOALL_SGL"), &
1069 fftsg_sizes=.NOT. section_get_lval(global_section, "EXTENDED_FFT_LENGTHS"), &
1070 pool_limit=globenv%fft_pool_scratch_limit, &
1071 wisdom_file=globenv%fftw_wisdom_file_name, &
1072 plan_style=globenv%fftw_plan_type)
1073
1074 CALL fft3d(fwfft, n, zz, status=stat)
1075 END IF
1076 IF (stat /= 0) THEN
1077 message = "FFT library "//trim(globenv%default_fft_library)// &
1078 " is not available. Trying FFT library FFTSG."
1079 cpwarn(trim(message))
1080 globenv%default_fft_library = "FFTSG"
1081 CALL init_fft(fftlib=trim(globenv%default_fft_library), &
1082 alltoall=section_get_lval(global_section, "ALLTOALL_SGL"), &
1083 fftsg_sizes=.NOT. section_get_lval(global_section, "EXTENDED_FFT_LENGTHS"), &
1084 pool_limit=globenv%fft_pool_scratch_limit, &
1085 wisdom_file=globenv%fftw_wisdom_file_name, &
1086 plan_style=globenv%fftw_plan_type)
1087
1088 CALL fft3d(fwfft, n, zz, status=stat)
1089 IF (stat /= 0) THEN
1090 cpabort("FFT library FFTSG does not work. No FFT library available.")
1091 END IF
1092 END IF
1093 END IF
1094
1095 END SUBROUTINE fft_setup_library
1096
1097! **************************************************************************************************
1098!> \brief availability diagonalizatioon library
1099!>
1100!> \param globenv ...
1101!> \author MI
1102! **************************************************************************************************
1103 SUBROUTINE diag_setup_library(globenv)
1104 TYPE(global_environment_type), POINTER :: globenv
1105
1106 CHARACTER(LEN=3*default_string_length) :: message
1107 LOGICAL :: fallback_applied
1108
1109 CALL diag_init(diag_lib=trim(globenv%diag_library), &
1110 fallback_applied=fallback_applied, &
1111 elpa_kernel=globenv%k_elpa, &
1112 elpa_neigvec_min_input=globenv%elpa_neigvec_min, &
1113 elpa_qr=globenv%elpa_qr, &
1114 elpa_print=globenv%elpa_print, &
1115 elpa_qr_unsafe=globenv%elpa_qr_unsafe, &
1116 dlaf_neigvec_min_input=globenv%dlaf_neigvec_min, &
1117 eps_check_diag_input=globenv%eps_check_diag)
1118
1119 IF (fallback_applied) THEN
1120 message = "Diagonalization library "//trim(globenv%diag_library)// &
1121 " is not available. The ScaLAPACK library is used as fallback."
1122 cpwarn(trim(message))
1123 END IF
1124
1125 END SUBROUTINE diag_setup_library
1126
1127! **************************************************************************************************
1128!> \brief ...
1129!> \param glob_section ...
1130! **************************************************************************************************
1131 SUBROUTINE fm_setup(glob_section)
1132 TYPE(section_vals_type), POINTER :: glob_section
1133
1134 INTEGER :: mm_type, ncb, nrb
1135 LOGICAL :: force_me
1136 TYPE(section_vals_type), POINTER :: fm_section
1137
1138 fm_section => section_vals_get_subs_vals(glob_section, "FM")
1139
1140 CALL section_vals_val_get(fm_section, "NROW_BLOCKS", i_val=nrb)
1141 CALL section_vals_val_get(fm_section, "NCOL_BLOCKS", i_val=ncb)
1142 CALL section_vals_val_get(fm_section, "FORCE_BLOCK_SIZE", l_val=force_me)
1143 CALL cp_fm_struct_config(nrow_block=nrb, ncol_block=ncb, force_block=force_me)
1144
1145 CALL section_vals_val_get(fm_section, "TYPE_OF_MATRIX_MULTIPLICATION", i_val=mm_type)
1146 CALL cp_fm_setup(mm_type)
1147
1148 END SUBROUTINE fm_setup
1149
1150! **************************************************************************************************
1151!> \brief ...
1152!> \param glob_section ...
1153! **************************************************************************************************
1154 SUBROUTINE dgemm_setup(glob_section)
1155 TYPE(section_vals_type), POINTER :: glob_section
1156
1157 INTEGER :: dgemm_type
1158
1159 CALL section_vals_val_get(glob_section, "PREFERRED_DGEMM_LIBRARY", i_val=dgemm_type)
1160
1161 CALL local_gemm_set_library(dgemm_type)
1162
1163 END SUBROUTINE dgemm_setup
1164
1165! **************************************************************************************************
1166!> \brief Parses the input section used to define the heuristic rules which determine if
1167!> a FM matrix should be redistributed before diagonalizing it.
1168!> \param glob_section the global input section
1169!> \author Nico Holmberg [01.2018]
1170! **************************************************************************************************
1171 SUBROUTINE fm_diag_rules_setup(glob_section)
1172 TYPE(section_vals_type), POINTER :: glob_section
1173
1174 INTEGER :: a, x
1175 LOGICAL :: elpa_force_redistribute, should_print
1176 TYPE(section_vals_type), POINTER :: section
1177
1178 section => section_vals_get_subs_vals(glob_section, "FM_DIAG_SETTINGS")
1179
1180 CALL section_vals_val_get(section, "PARAMETER_A", i_val=a)
1181 CALL section_vals_val_get(section, "PARAMETER_X", i_val=x)
1182 CALL section_vals_val_get(section, "PRINT_FM_REDISTRIBUTE", l_val=should_print)
1183 CALL section_vals_val_get(section, "ELPA_FORCE_REDISTRIBUTE", l_val=elpa_force_redistribute)
1184
1185 CALL cp_fm_redistribute_init(a, x, should_print, elpa_force_redistribute)
1186
1187 END SUBROUTINE fm_diag_rules_setup
1188! **************************************************************************************************
1189!> \brief reads the Walltime also in format HH:MM:SS
1190!> \param section ...
1191!> \param keyword_name ...
1192!> \param walltime ...
1193!> \par History
1194!> none
1195!> \author Mandes
1196! **************************************************************************************************
1197 SUBROUTINE cp2k_get_walltime(section, keyword_name, walltime)
1198 TYPE(section_vals_type), POINTER :: section
1199 CHARACTER(LEN=*), INTENT(in) :: keyword_name
1200 REAL(kind=dp), INTENT(out) :: walltime
1201
1202 CHARACTER(LEN=1) :: c1, c2
1203 CHARACTER(LEN=100) :: txt
1204 INTEGER :: hours, ierr, minutes, n, seconds
1205
1206 CALL section_vals_val_get(section, keyword_name, c_val=txt)
1207 n = len_trim(txt)
1208
1209 IF (n == 0) THEN
1210 walltime = -1.0_dp
1211 ELSE IF (index(txt, ":") == 0) THEN
1212 READ (txt(1:n), fmt=*, iostat=ierr) walltime
1213 IF (ierr /= 0) cpabort('Could not parse WALLTIME: "'//txt(1:n)//'"')
1214 ELSE
1215 READ (txt(1:n), fmt="(I2,A1,I2,A1,I2)", iostat=ierr) hours, c1, minutes, c2, seconds
1216 IF (n /= 8 .OR. ierr /= 0 .OR. c1 .NE. ":" .OR. c2 .NE. ":") &
1217 cpabort('Could not parse WALLTIME: "'//txt(1:n)//'"')
1218 walltime = 3600.0_dp*real(hours, dp) + 60.0_dp*real(minutes, dp) + real(seconds, dp)
1219 END IF
1220 END SUBROUTINE cp2k_get_walltime
1221
1222! **************************************************************************************************
1223!> \brief Writes final timings and banner for CP2K
1224!> \param root_section ...
1225!> \param para_env ...
1226!> \param globenv ...
1227!> \param wdir ...
1228!> \param q_finalize ...
1229!> \par History
1230!> none
1231!> \author JGH,MK
1232!> \note
1233!> The following routines need to be synchronized wrt. adding/removing
1234!> of the default environments (logging, performance,error):
1235!> environment:cp2k_init, environment:cp2k_finalize,
1236!> f77_interface:f_env_add_defaults, f77_interface:f_env_rm_defaults,
1237!> f77_interface:create_force_env, f77_interface:destroy_force_env
1238! **************************************************************************************************
1239 SUBROUTINE cp2k_finalize(root_section, para_env, globenv, wdir, q_finalize)
1240
1241 TYPE(section_vals_type), POINTER :: root_section
1242 TYPE(mp_para_env_type), POINTER :: para_env
1243 TYPE(global_environment_type), POINTER :: globenv
1244 CHARACTER(LEN=*), OPTIONAL :: wdir
1245 LOGICAL, INTENT(IN), OPTIONAL :: q_finalize
1246
1247 CHARACTER(LEN=default_path_length) :: cg_filename
1248 INTEGER :: cg_mode, iw, unit_exit
1249 LOGICAL :: delete_it, do_finalize, report_maxloc, &
1250 sort_by_self_time
1251 REAL(kind=dp) :: r_timings
1252 TYPE(cp_logger_type), POINTER :: logger
1253
1254 ! Look if we inherited a failure, more care is needed if so
1255 ! i.e. the input is most likely not available
1256 ! Set flag if this is a development version
1257
1258 do_finalize = .true.
1259 IF (PRESENT(q_finalize)) do_finalize = q_finalize
1260 ! Clean up
1261 NULLIFY (logger)
1262 logger => cp_get_default_logger()
1263 IF (do_finalize) THEN
1267 CALL diag_finalize()
1268 ! finalize the fft (i.e. writes the wisdom if FFTW3 )
1269 CALL finalize_fft(para_env, globenv%fftw_wisdom_file_name)
1270 CALL finalize_libvori()
1271 END IF
1272
1273 ! Write message passing performance info
1274
1275 iw = cp_print_key_unit_nr(logger, root_section, "GLOBAL%PROGRAM_RUN_INFO", &
1276 extension=".log")
1277 CALL describe_mp_perf_env(iw)
1278 CALL cp_print_key_finished_output(iw, logger, root_section, &
1279 "GLOBAL%PROGRAM_RUN_INFO")
1280
1281 CALL collect_citations_from_ranks(para_env)
1282 iw = cp_print_key_unit_nr(logger, root_section, "GLOBAL%REFERENCES", &
1283 extension=".Log")
1284 IF (iw > 0) THEN
1285 WRITE (unit=iw, fmt="(/,T2,A)") repeat("-", 79)
1286 WRITE (unit=iw, fmt="(T2,A,T80,A)") "-", "-"
1287 WRITE (unit=iw, fmt="(T2,A,T30,A,T80,A)") "-", "R E F E R E N C E S", "-"
1288 WRITE (unit=iw, fmt="(T2,A,T80,A)") "-", "-"
1289 WRITE (unit=iw, fmt="(T2,A)") repeat("-", 79)
1290 WRITE (unit=iw, fmt="(T2,A)") ""
1291 WRITE (unit=iw, fmt="(T2,A)") trim(cp2k_version)//", the CP2K developers group ("//trim(cp2k_year)//")."
1292 WRITE (unit=iw, fmt="(T2,A)") "CP2K is freely available from "//trim(cp2k_home)//" ."
1293 WRITE (unit=iw, fmt="(T2,A)") ""
1294 CALL print_cited_references(unit=iw)
1295 END IF
1296 CALL cp_print_key_finished_output(iw, logger, root_section, &
1297 "GLOBAL%REFERENCES")
1298
1299 CALL timestop(globenv%handle) ! corresponding the "CP2K" in cp2k_init
1300
1301 iw = cp_print_key_unit_nr(logger, root_section, "GLOBAL%TIMINGS", &
1302 extension=".Log")
1303 r_timings = section_get_rval(root_section, "GLOBAL%TIMINGS%THRESHOLD")
1304 sort_by_self_time = section_get_lval(root_section, "GLOBAL%TIMINGS%SORT_BY_SELF_TIME")
1305 report_maxloc = section_get_lval(root_section, "GLOBAL%TIMINGS%REPORT_MAXLOC")
1306 IF (m_energy() .NE. 0.0_dp) THEN
1307 CALL timings_report_print(iw, r_timings, sort_by_self_time, cost_type_energy, report_maxloc, para_env)
1308 END IF
1309 CALL timings_report_print(iw, r_timings, sort_by_self_time, cost_type_time, report_maxloc, para_env)
1310
1311 ! Write the callgraph, if desired by user
1312 CALL section_vals_val_get(root_section, "GLOBAL%CALLGRAPH", i_val=cg_mode)
1313 IF (cg_mode /= callgraph_none) THEN
1314 CALL section_vals_val_get(root_section, "GLOBAL%CALLGRAPH_FILE_NAME", c_val=cg_filename)
1315 IF (len_trim(cg_filename) == 0) cg_filename = trim(logger%iter_info%project_name)
1316 IF (cg_mode == callgraph_all) & !incorporate mpi-rank into filename
1317 cg_filename = trim(cg_filename)//"_"//trim(adjustl(cp_to_string(para_env%mepos)))
1318 IF (iw > 0) THEN
1319 WRITE (unit=iw, fmt="(T2,3X,A)") "Writing callgraph to: "//trim(cg_filename)//".callgraph"
1320 WRITE (unit=iw, fmt="()")
1321 WRITE (unit=iw, fmt="(T2,A)") "-------------------------------------------------------------------------------"
1322 END IF
1323 IF (cg_mode == callgraph_all .OR. para_env%is_source()) &
1324 CALL timings_report_callgraph(trim(cg_filename)//".callgraph")
1325 END IF
1326
1327 CALL cp_print_key_finished_output(iw, logger, root_section, &
1328 "GLOBAL%TIMINGS")
1329
1330 CALL rm_mp_perf_env()
1331 CALL rm_timer_env()
1332
1333 IF (para_env%is_source()) THEN
1334 iw = cp_print_key_unit_nr(logger, root_section, "GLOBAL%PROGRAM_RUN_INFO", &
1335 extension=".log")
1336
1337 ! Deleting (if existing) the external EXIT files
1338 delete_it = .false.
1339 INQUIRE (file="EXIT", exist=delete_it)
1340 IF (delete_it) THEN
1341 CALL open_file(file_name="EXIT", unit_number=unit_exit)
1342 CALL close_file(unit_number=unit_exit, file_status="DELETE")
1343 END IF
1344
1345 delete_it = .false.
1346 INQUIRE (file=trim(logger%iter_info%project_name)//".EXIT", exist=delete_it)
1347 IF (delete_it) THEN
1348 CALL open_file(file_name=trim(logger%iter_info%project_name)//".EXIT", unit_number=unit_exit)
1349 CALL close_file(unit_number=unit_exit, file_status="DELETE")
1350 END IF
1351
1352 ! Print OpenMP issue counter and number of warnings for this workload
1353 IF (iw > 0) THEN
1354 IF (0 <= m_omp_trace_issues()) THEN
1355 WRITE (iw, "(T2,A,I0)") "The number of traced issues for OpenMP : ", m_omp_trace_issues()
1356 END IF
1357 WRITE (iw, "(T2,A,I0)") "The number of warnings for this run is : ", warning_counter
1358 WRITE (iw, *) ""
1359 WRITE (unit=iw, fmt="(T2,A)") repeat("-", 79)
1360 END IF
1361
1362 ! Update the runtime environment variables
1363 CALL get_runtime_info()
1364
1365 ! Just a choice, do not print the CP2K footer if there is a failure
1366 CALL cp2k_footer(iw, wdir)
1367 IF (iw > 0) FLUSH (iw) ! ignore &GLOBAL / FLUSH_SHOULD_FLUSH
1368
1369 CALL cp_print_key_finished_output(iw, logger, root_section, &
1370 "GLOBAL%PROGRAM_RUN_INFO")
1371 END IF
1372
1373 ! Release message passing environment
1375
1376 END SUBROUTINE cp2k_finalize
1377
1378END MODULE environment
Target architecture or instruction set extension according to compiler target flags.
Definition machine.F:96
Trace OpenMP constructs if ennvironment variable CP2K_OMP_TRACE=1.
Definition machine.F:107
collects all references to literature in CP2K as new algorithms / method are included from literature...
integer, save, public marek2014
integer, save, public solca2024
integer, save, public frigo2005
some minimal info about CP2K, including its version and license
Definition cp2k_info.F:16
character(len=default_string_length), public r_host_name
Definition cp2k_info.F:68
character(len= *), parameter, public cp2k_home
Definition cp2k_info.F:44
character(len= *), parameter, public compile_host
Definition cp2k_info.F:62
character(len= *), parameter, public compile_arch
Definition cp2k_info.F:50
integer, public r_pid
Definition cp2k_info.F:69
character(len= *), parameter, public compile_revision
Definition cp2k_info.F:38
character(len= *), parameter, public compile_date
Definition cp2k_info.F:56
character(len= *), parameter, public cp2k_year
Definition cp2k_info.F:43
character(len=10 *default_string_length) function, public cp2k_flags()
list all compile time options that influence the capabilities of cp2k. All new flags should be added ...
Definition cp2k_info.F:80
character(len= *), parameter, public cp2k_version
Definition cp2k_info.F:42
subroutine, public get_runtime_info()
...
Definition cp2k_info.F:336
character(len=default_string_length), public r_user_name
Definition cp2k_info.F:68
Module that contains the routines for error handling.
integer, save, public warning_counter
Utility routines to open and close files. Tracking of preconnections.
Definition cp_files.F:16
subroutine, public open_file(file_name, file_status, file_form, file_action, file_position, file_pad, unit_number, debug, skip_get_unit_number, file_access)
Opens the requested file using a free unit number.
Definition cp_files.F:308
subroutine, public close_file(unit_number, file_status, keep_preconnection)
Close an open file given by its logical unit number. Optionally, keep the file and unit preconnected.
Definition cp_files.F:119
character(len=default_path_length) function, public get_data_dir()
Returns path of data directory if set, otherwise an empty string.
Definition cp_files.F:555
various cholesky decomposition related routines
integer, parameter, public fm_cholesky_type_dlaf
integer, parameter, public fm_cholesky_type_scalapack
integer, save, public dlaf_cholesky_n_min
integer, save, public cholesky_type
Auxiliary tools to redistribute cp_fm_type matrices before and after diagonalization....
subroutine, public cp_fm_redistribute_init(a, x, should_print, elpa_force_redistribute)
Initializes the parameters that determine how to calculate the optimal number of CPUs for diagonalizi...
used for collecting some of the diagonalization schemes available for cp_fm_type. cp_fm_power also mo...
Definition cp_fm_diag.F:17
real(kind=dp), parameter, public eps_check_diag_default
Definition cp_fm_diag.F:87
integer, parameter, public fm_diag_type_cusolver
Definition cp_fm_diag.F:105
integer, parameter, public fm_diag_type_dlaf
Definition cp_fm_diag.F:105
integer, parameter, public fm_diag_type_scalapack
Definition cp_fm_diag.F:105
subroutine, public diag_finalize()
Finalize the diagonalization library.
Definition cp_fm_diag.F:214
subroutine, public diag_init(diag_lib, fallback_applied, elpa_kernel, elpa_neigvec_min_input, elpa_qr, elpa_print, elpa_qr_unsafe, dlaf_neigvec_min_input, eps_check_diag_input)
Setup the diagonalization library to be used.
Definition cp_fm_diag.F:153
integer, parameter, public fm_diag_type_elpa
Definition cp_fm_diag.F:105
represent the structure of a full matrix
subroutine, public cp_fm_struct_config(nrow_block, ncol_block, force_block)
allows to modify the default settings for matrix creation
represent a full matrix distributed on many processors
Definition cp_fm_types.F:15
integer function, public cp_fm_get_mm_type()
...
subroutine, public cp_fm_setup(mm_type)
...
various routines to log and control the output. The idea is that decisions about where to log should ...
recursive integer function, public cp_logger_get_default_unit_nr(logger, local, skip_not_ionode)
asks the default unit number of the given logger. try to use cp_logger_get_unit_nr
subroutine, public cp_logger_set(logger, local_filename, global_filename)
sets various attributes of the given logger
subroutine, public cp_rm_default_logger()
the cousin of cp_add_default_logger, decrements the stack, so that the default logger is what it has ...
subroutine, public cp_logger_release(logger)
releases this logger
subroutine, public cp_logger_create(logger, para_env, print_level, default_global_unit_nr, default_local_unit_nr, global_filename, local_filename, close_global_unit_on_dealloc, iter_info, close_local_unit_on_dealloc, suffix, template_logger)
initializes a logger
subroutine, public cp_add_default_logger(logger)
adds a default logger. MUST be called before logging occours
type(cp_logger_type) function, pointer, public cp_get_default_logger()
returns the default logger
routines to handle the output, The idea is to remove the decision of wheter to output and what to out...
subroutine, public cp_mpi_io_set(flag)
Sets flag which determines whether or not to use MPI I/O for I/O routines that have been parallized w...
integer, parameter, public debug_print_level
integer function, public cp_print_key_unit_nr(logger, basis_section, print_key_path, extension, middle_name, local, log_filename, ignore_should_output, file_form, file_position, file_action, file_status, do_backup, on_file, is_new_file, mpi_io, fout)
...
integer, parameter, public low_print_level
integer, parameter, public medium_print_level
subroutine, public cp_print_key_finished_output(unit_nr, logger, basis_section, print_key_path, local, ignore_should_output, on_file, mpi_io)
should be called after you finish working with a unit obtained with cp_print_key_unit_nr,...
integer, parameter, public high_print_level
integer, parameter, public silent_print_level
Sets up and terminates the global environment variables.
Definition environment.F:17
subroutine, public cp2k_finalize(root_section, para_env, globenv, wdir, q_finalize)
Writes final timings and banner for CP2K.
subroutine, public cp2k_read(root_section, para_env, globenv)
read part of cp2k_init
subroutine, public cp2k_get_walltime(section, keyword_name, walltime)
reads the Walltime also in format HH:MM:SS
subroutine, public cp2k_init(para_env, output_unit, globenv, input_file_name, wdir)
Initializes a CP2K run (setting of the global environment variables)
subroutine, public cp2k_setup(root_section, para_env, globenv)
globenv initializations that need the input and error
subroutine, public init_fft(fftlib, alltoall, fftsg_sizes, pool_limit, wisdom_file, plan_style)
...
Definition fft_tools.F:186
integer, parameter, public fwfft
Definition fft_tools.F:146
subroutine, public finalize_fft(para_env, wisdom_file)
does whatever is needed to finalize the current fft setup
Definition fft_tools.F:216
Interface for the force calculations.
subroutine, public multiple_fe_list(force_env_sections, root_section, i_force_eval, nforce_eval)
returns the order of the multiple force_env
Calculation of the incomplete Gamma function F_n(t) for multi-center integrals over Cartesian Gaussia...
Definition gamma.F:15
subroutine, public deallocate_md_ftable()
Deallocate the table of F_n(t) values.
Definition gamma.F:121
Define type storing the global information of a run. Keep the amount of stored data small....
Fortran API for the grid package, which is written in C.
Definition grid_api.F:12
integer, parameter, public grid_backend_auto
Definition grid_api.F:66
integer, parameter, public grid_backend_gpu
Definition grid_api.F:70
integer, parameter, public grid_backend_hip
Definition grid_api.F:71
integer, parameter, public grid_backend_dgemm
Definition grid_api.F:69
integer, parameter, public grid_backend_cpu
Definition grid_api.F:68
integer, parameter, public grid_backend_ref
Definition grid_api.F:67
subroutine, public cp2k_header(iw, wdir)
...
Definition header.F:40
subroutine, public cp2k_footer(iw, wdir)
...
Definition header.F:69
collects all constants needed in input so that they can be used without circular dependencies
integer, parameter, public energy_run
integer, parameter, public do_fft_sg
integer, parameter, public callgraph_all
integer, parameter, public do_farming
integer, parameter, public do_cosma
integer, parameter, public do_cp2k
integer, parameter, public do_scalapack
integer, parameter, public do_eip
integer, parameter, public do_test
integer, parameter, public do_fist
integer, parameter, public do_sirius
integer, parameter, public do_dgemm_blas
integer, parameter, public callgraph_none
integer, parameter, public mol_dyn_run
integer, parameter, public do_fft_fftw3
integer, parameter, public none_run
integer, parameter, public do_qs
integer, parameter, public do_dgemm_spla
builds the global input section for cp2k
subroutine, public create_global_section(section)
section to hold global settings for the whole program
represents an enumeration, i.e. a mapping between integers and strings
character(len=default_string_length) function, public enum_i2c(enum, i)
maps an integer to a string
represents keywords in an input
subroutine, public keyword_get(keyword, names, usage, description, type_of_var, n_var, default_value, lone_keyword_value, repeats, enum, citations)
...
objects that represent the structure of input sections and the data contained in an input section
real(kind=dp) function, public section_get_rval(section_vals, keyword_name)
...
subroutine, public section_vals_val_set(section_vals, keyword_name, i_rep_section, i_rep_val, val, l_val, i_val, r_val, c_val, l_vals_ptr, i_vals_ptr, r_vals_ptr, c_vals_ptr)
sets the requested value
integer function, public section_get_ival(section_vals, keyword_name)
...
recursive type(section_vals_type) function, pointer, public section_vals_get_subs_vals(section_vals, subsection_name, i_rep_section, can_return_null)
returns the values of the requested subsection
recursive subroutine, public section_release(section)
releases the given keyword list (see doc/ReferenceCounting.html)
recursive type(keyword_type) function, pointer, public section_get_keyword(section, keyword_name)
returns the requested keyword
subroutine, public section_vals_get(section_vals, ref_count, n_repetition, n_subs_vals_rep, section, explicit)
returns various attributes about the section_vals
type(section_vals_type) function, pointer, public section_vals_get_subs_vals3(section_vals, subsection_name, i_rep_section)
returns the values of the n-th non default subsection (null if no such section exists (not so many no...
subroutine, public section_vals_val_get(section_vals, keyword_name, i_rep_section, i_rep_val, n_rep_val, val, l_val, i_val, r_val, c_val, l_vals, i_vals, r_vals, c_vals, explicit)
returns the requested value
logical function, public section_get_lval(section_vals, keyword_name)
...
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public int_8
Definition kinds.F:54
integer, parameter, public dp
Definition kinds.F:34
integer, parameter, public default_string_length
Definition kinds.F:57
integer, parameter, public default_path_length
Definition kinds.F:58
subroutine, public print_kind_info(iw)
Print informations about the used data types.
Definition kinds.F:72
subroutine, public local_gemm_set_library(dgemm_library)
...
Machine interface based on Fortran 2003 and POSIX.
Definition machine.F:17
logical, save, public flush_should_flush
Definition machine.F:115
integer function, public m_procrun(pid)
Returns if a process is running on the local machine 1 if yes and 0 if not.
Definition machine.F:407
subroutine, public m_memory_details(memtotal, memfree, buffers, cached, slab, sreclaimable, memlikelyfree)
get more detailed memory info, all units are bytes. the only 'useful' option is MemLikelyFree which i...
Definition machine.F:510
pure integer function, public m_cpuid()
Target architecture or instruction set extension according to CPU-check at runtime.
Definition machine.F:194
pure integer function, public m_cpuid_vlen(cpuid, typesize)
Determine vector-length for a given CPUID.
Definition machine.F:276
real(kind=dp) function, public m_energy()
returns the energy used since some time in the past. The precise meaning depends on the infrastructur...
Definition machine.F:327
subroutine, public m_cpuinfo(model_name)
reads /proc/cpuinfo if it exists (i.e. Linux) to return relevant info
Definition machine.F:161
pure character(len=default_string_length) function, public m_cpuid_name(cpuid)
Determine name of target architecture for a given CPUID.
Definition machine.F:232
subroutine, public m_omp_get_stacksize(omp_stacksize)
Retrieve environment variable OMP_STACKSIZE.
Definition machine.F:762
Interface to the message passing library MPI.
logical, save, public mp_collect_timings
Defines all routines to deal with the performance of MPI routines.
Definition mp_perf_env.F:11
subroutine, public rm_mp_perf_env()
...
subroutine, public describe_mp_perf_env(scr)
...
subroutine, public add_mp_perf_env(perf_env)
start and stop the performance indicators for every call to start there has to be (exactly) one call ...
Definition mp_perf_env.F:76
Provides Cartesian and spherical orbital pointers and indices.
subroutine, public init_orbital_pointers(maxl)
Initialize or update the orbital pointers.
subroutine, public deallocate_orbital_pointers()
Deallocate the orbital pointers.
Calculation of the spherical harmonics and the corresponding orbital transformation matrices.
subroutine, public init_spherical_harmonics(maxl, output_unit)
Initialize or update the orbital transformation matrices.
subroutine, public deallocate_spherical_harmonics()
Deallocate the orbital transformation matrices.
Parallel (pseudo)random number generator (RNG) for multiple streams and substreams of random numbers.
subroutine, public check_rng(output_unit, ionode)
...
subroutine, public write_rng_matrices(output_unit)
Write the transformation matrices of the two MRG components (raised to the specified output)
integer, parameter, public gaussian
Definition of physical constants:
Definition physcon.F:68
real(kind=dp), parameter, public seconds
Definition physcon.F:150
subroutine, public write_physcon(output_unit)
Write all basic physical constants used by CP2K to a logical output unit.
Definition physcon.F:217
provides a uniform framework to add references to CP2K cite and output these
subroutine, public collect_citations_from_ranks(para_env)
Checks for each reference if any mpi-rank has marked it for citation.
subroutine, public cite_reference(key)
marks a given reference as cited.
subroutine, public print_cited_references(unit)
printout of all cited references in the journal format sorted by publication year
Utilities for string manipulations.
subroutine, public integer_to_string(inumber, string)
Converts an integer number to a string. The WRITE statement will return an error message,...
subroutine, public string_to_ascii(string, nascii)
Convert a string to sequence of integer numbers.
subroutine, public ascii_to_string(nascii, string)
Convert a sequence of integer numbers (ASCII code) to a string. Blanks are inserted for invalid ASCII...
Timing routines for accounting.
integer, parameter, public cost_type_energy
subroutine, public timings_report_callgraph(filename)
Write accumulated callgraph information as cachegrind-file. http://kcachegrind.sourceforge....
integer, parameter, public cost_type_time
subroutine, public timings_report_print(iw, r_timings, sort_by_self_time, cost_type, report_maxloc, para_env)
Print accumulated information on timers.
Timing routines for accounting.
Definition timings.F:17
integer, save, public global_timings_level
Definition timings.F:68
subroutine, public timings_setup_tracing(trace_max, unit_nr, trace_str, routine_names)
Set routine tracer.
Definition timings.F:398
subroutine, public add_timer_env(timer_env)
adds the given timer_env to the top of the stack
Definition timings.F:93
subroutine, public rm_timer_env()
removes the current timer env from the stack
Definition timings.F:134
character(len=default_string_length), parameter, public root_cp2k_name
Definition timings.F:70
Interface for Voronoi Integration and output of BQB files.
subroutine, public finalize_libvori()
Call libvori's finalize if support is compiled in.
type of a logger, at the moment it contains just a print level starting at which level it should be l...
contains the initially parsed file and the initial parallel environment
represent a keyword in the input
represent a section of the input file
stores all the informations relevant to an mpi environment