(git:374b731)
Loading...
Searching...
No Matches
dm_ls_chebyshev.F
Go to the documentation of this file.
1!--------------------------------------------------------------------------------------------------!
2! CP2K: A general program to perform molecular dynamics simulations !
3! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
4! !
5! SPDX-License-Identifier: GPL-2.0-or-later !
6!--------------------------------------------------------------------------------------------------!
7
8! **************************************************************************************************
9!> \brief Routines using linear scaling chebyshev methods
10!> \par History
11!> 2012.10 created [Jinwoong Cha]
12!> \author Jinwoong Cha
13! **************************************************************************************************
19 USE cp_output_handling, ONLY: cp_p_file,&
23 USE dbcsr_api, ONLY: &
24 dbcsr_add, dbcsr_add_on_diag, dbcsr_copy, dbcsr_create, dbcsr_frobenius_norm, &
25 dbcsr_get_info, dbcsr_get_occupation, dbcsr_multiply, dbcsr_release, dbcsr_scale, &
26 dbcsr_set, dbcsr_trace, dbcsr_type, dbcsr_type_no_symmetry
31 USE kinds, ONLY: default_string_length,&
32 dp
33 USE machine, ONLY: m_flush,&
35 USE mathconstants, ONLY: pi
37#include "./base/base_uses.f90"
38
39 IMPLICIT NONE
40
41 PRIVATE
42
43 CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'dm_ls_chebyshev'
44
45 PUBLIC :: compute_chebyshev
46
47CONTAINS
48
49! **************************************************************************************************
50!> \brief compute chebyshev polynomials up to order n for a given value of x
51!> \param value ...
52!> \param x ...
53!> \param n ...
54!> \par History
55!> 2012.11 created [Jinwoong Cha]
56!> \author Jinwoong Cha
57! **************************************************************************************************
58 SUBROUTINE chebyshev_poly(value, x, n)
59 REAL(KIND=dp), INTENT(OUT) :: value
60 REAL(KIND=dp), INTENT(IN) :: x
61 INTEGER, INTENT(IN) :: n
62
63!polynomial values
64!number of chev polynomials
65
66 value = cos((n - 1)*acos(x))
67
68 END SUBROUTINE chebyshev_poly
69
70! **************************************************************************************************
71!> \brief kernel for chebyshev polynomials expansion (Jackson kernel)
72!> \param value ...
73!> \param n ...
74!> \param nc ...
75!> \par History
76!> 2012.11 created [Jinwoong Cha]
77!> \author Jinwoong Cha
78! **************************************************************************************************
79 SUBROUTINE kernel(value, n, nc)
80 REAL(KIND=dp), INTENT(OUT) :: value
81 INTEGER, INTENT(IN) :: n, nc
82
83!kernel at n
84!n-1 order of chebyshev polynomials
85!number of total chebyshev polynomials
86!Kernel define
87
88 value = 1.0_dp/(nc + 1.0_dp)*((nc - (n - 1) + 1.0_dp)* &
89 cos(pi*(n - 1)/(nc + 1.0_dp)) + sin(pi*(n - 1)/(nc + 1.0_dp))*1.0_dp/tan(pi/(nc + 1.0_dp)))
90
91 END SUBROUTINE kernel
92
93! **************************************************************************************************
94!> \brief compute properties based on chebyshev expansion
95!> \param qs_env ...
96!> \param ls_scf_env ...
97!> \par History
98!> 2012.10 created [Jinwoong Cha]
99!> \author Jinwoong Cha
100! **************************************************************************************************
101 SUBROUTINE compute_chebyshev(qs_env, ls_scf_env)
102 TYPE(qs_environment_type), POINTER :: qs_env
103 TYPE(ls_scf_env_type) :: ls_scf_env
104
105 CHARACTER(len=*), PARAMETER :: routinen = 'compute_chebyshev'
106 REAL(kind=dp), PARAMETER :: scale_evals = 1.01_dp
107
108 CHARACTER(LEN=30) :: middle_name
109 CHARACTER(LEN=default_string_length) :: title
110 INTEGER :: handle, icheb, igrid, iinte, ispin, &
111 iwindow, n_gridpoint_dos, ncheb, &
112 ninte, nrows, nwindow, unit_cube, &
113 unit_dos, unit_nr
114 LOGICAL :: converged, write_cubes
115 REAL(kind=dp) :: chev_t, chev_t_dos, dummy1, final, frob_matrix, initial, interval_a, &
116 interval_b, max_ev, min_ev, occ, orbital_occ, summa, t1, t2
117 REAL(kind=dp), ALLOCATABLE, DIMENSION(:) :: chev_e, chev_es_dos, dos, dummy2, ev1, &
118 ev2, kernel_g, mu, sev1, sev2, trace_dm
119 REAL(kind=dp), ALLOCATABLE, DIMENSION(:, :) :: aitchev_t, e_inte, gdensity, sqrt_vec
120 REAL(kind=dp), DIMENSION(:), POINTER :: tmp_r
121 TYPE(cp_logger_type), POINTER :: logger
122 TYPE(dbcsr_type) :: matrix_dummy1, matrix_f, matrix_tmp1, &
123 matrix_tmp2, matrix_tmp3
124 TYPE(dbcsr_type), DIMENSION(:), POINTER :: matrix_dummy2
125
126 IF (.NOT. ls_scf_env%chebyshev%compute_chebyshev) RETURN
127
128 CALL timeset(routinen, handle)
129
130 ! get a useful output_unit
131 logger => cp_get_default_logger()
132 IF (logger%para_env%is_source()) THEN
133 unit_nr = cp_logger_get_default_unit_nr(logger, local=.true.)
134 ELSE
135 unit_nr = -1
136 END IF
137
138 ncheb = ls_scf_env%chebyshev%n_chebyshev
139 ninte = 2*ncheb
140 n_gridpoint_dos = ls_scf_env%chebyshev%n_gridpoint_dos
141
142 write_cubes = btest(cp_print_key_should_output(logger%iter_info, ls_scf_env%chebyshev%print_key_cube), cp_p_file)
143 IF (write_cubes) THEN
144 IF (ASSOCIATED(ls_scf_env%chebyshev%min_energy)) DEALLOCATE (ls_scf_env%chebyshev%min_energy)
145 CALL section_vals_val_get(ls_scf_env%chebyshev%print_key_cube, "MIN_ENERGY", r_vals=tmp_r)
146 ALLOCATE (ls_scf_env%chebyshev%min_energy(SIZE(tmp_r)))
147 ls_scf_env%chebyshev%min_energy = tmp_r
148
149 IF (ASSOCIATED(ls_scf_env%chebyshev%max_energy)) DEALLOCATE (ls_scf_env%chebyshev%max_energy)
150 CALL section_vals_val_get(ls_scf_env%chebyshev%print_key_cube, "MAX_ENERGY", r_vals=tmp_r)
151 ALLOCATE (ls_scf_env%chebyshev%max_energy(SIZE(tmp_r)))
152 ls_scf_env%chebyshev%max_energy = tmp_r
153
154 nwindow = SIZE(ls_scf_env%chebyshev%min_energy)
155 ELSE
156 nwindow = 0
157 END IF
158
159 ALLOCATE (ev1(1:nwindow))
160 ALLOCATE (ev2(1:nwindow))
161 ALLOCATE (sev1(1:nwindow))
162 ALLOCATE (sev2(1:nwindow))
163 ALLOCATE (trace_dm(1:nwindow))
164 ALLOCATE (matrix_dummy2(1:nwindow))
165
166 DO iwindow = 1, nwindow
167 ev1(iwindow) = ls_scf_env%chebyshev%min_energy(iwindow)
168 ev2(iwindow) = ls_scf_env%chebyshev%max_energy(iwindow)
169 END DO
170
171 IF (unit_nr > 0) THEN
172 WRITE (unit_nr, '()')
173 WRITE (unit_nr, '(T2,A)') "STARTING CHEBYSHEV CALCULATION"
174 END IF
175
176 ! create 3 temporary matrices
177 CALL dbcsr_create(matrix_tmp1, template=ls_scf_env%matrix_s, matrix_type=dbcsr_type_no_symmetry)
178 CALL dbcsr_create(matrix_tmp2, template=ls_scf_env%matrix_s, matrix_type=dbcsr_type_no_symmetry)
179 CALL dbcsr_create(matrix_tmp3, template=ls_scf_env%matrix_s, matrix_type=dbcsr_type_no_symmetry)
180 CALL dbcsr_create(matrix_f, template=ls_scf_env%matrix_s, matrix_type=dbcsr_type_no_symmetry)
181 CALL dbcsr_create(matrix_dummy1, template=ls_scf_env%matrix_s, matrix_type=dbcsr_type_no_symmetry)
182
183 DO iwindow = 1, nwindow
184 CALL dbcsr_create(matrix_dummy2(iwindow), template=ls_scf_env%matrix_s, &
185 matrix_type=dbcsr_type_no_symmetry)
186 END DO
187
188 DO ispin = 1, SIZE(ls_scf_env%matrix_ks)
189 ! create matrix_F=inv(sqrt(S))*H*inv(sqrt(S))
190 CALL dbcsr_multiply("N", "N", 1.0_dp, ls_scf_env%matrix_s_sqrt_inv, ls_scf_env%matrix_ks(ispin), &
191 0.0_dp, matrix_tmp1, filter_eps=ls_scf_env%eps_filter)
192 CALL dbcsr_multiply("N", "N", 1.0_dp, matrix_tmp1, ls_scf_env%matrix_s_sqrt_inv, &
193 0.0_dp, matrix_f, filter_eps=ls_scf_env%eps_filter)
194
195 ! find largest and smallest eigenvalues
196 CALL arnoldi_extremal(matrix_f, max_ev, min_ev, converged=converged, max_iter=ls_scf_env%max_iter_lanczos, &
197 threshold=ls_scf_env%eps_lanczos) !Lanczos algorithm to calculate eigenvalue
198 IF (unit_nr > 0) WRITE (unit_nr, '(T2,A,2F16.8,A,L2)') &
199 "smallest largest eigenvalue", min_ev, max_ev, " converged ", converged
200 IF (nwindow > 0) THEN
201 IF (unit_nr > 0) WRITE (unit_nr, '(T2,A,1000F16.8)') "requested interval-min_energy", ev1(:)
202 IF (unit_nr > 0) WRITE (unit_nr, '(T2,A,1000F16.8)') "requested interval-max_energy", ev2(:)
203 END IF
204 interval_a = (max_ev - min_ev)*scale_evals/2
205 interval_b = (max_ev + min_ev)/2
206
207 sev1(:) = (ev1(:) - interval_b)/interval_a !scaled ev1 vector
208 sev2(:) = (ev2(:) - interval_b)/interval_a !scaled ev2 vector
209
210 !chebyshev domain,pi*sqrt(1-x^2) vector construction and chebyshev polynomials for integration (for g(E))
211 ALLOCATE (e_inte(1:ninte + 1, 1:nwindow))
212 ALLOCATE (sqrt_vec(1:ninte + 1, 1:nwindow))
213
214 DO iwindow = 1, nwindow
215 DO iinte = 1, ninte + 1
216 e_inte(iinte, iwindow) = sev1(iwindow) + ((sev2(iwindow) - sev1(iwindow))/ninte)*(iinte - 1)
217 sqrt_vec(iinte, iwindow) = pi*sqrt(1.0_dp - e_inte(iinte, iwindow)*e_inte(iinte, iwindow))
218 END DO
219 END DO
220
221 !integral.. (identical to the coefficient for g(E))
222
223 ALLOCATE (aitchev_t(1:ncheb, 1:nwindow)) !after intergral. =>ainte
224
225 DO iwindow = 1, nwindow
226 DO icheb = 1, ncheb
227 CALL chebyshev_poly(initial, e_inte(1, iwindow), icheb)
228 CALL chebyshev_poly(final, e_inte(1, iwindow), icheb)
229 summa = (sev2(iwindow) - sev1(iwindow))/(2.0_dp*ninte)*(initial/sqrt_vec(1, iwindow) + final/sqrt_vec(ninte + 1, iwindow))
230 DO iinte = 2, ninte
231 CALL chebyshev_poly(chev_t, e_inte(iinte, iwindow), icheb)
232 summa = summa + ((sev2(iwindow) - sev1(iwindow))/ninte)*(chev_t/sqrt_vec(iinte, iwindow))
233 END DO
234 aitchev_t(icheb, iwindow) = summa
235 summa = 0
236 END DO
237 END DO
238
239 ! scale the matrix to get evals in the interval -1,1
240 CALL dbcsr_add_on_diag(matrix_f, -interval_b)
241 CALL dbcsr_scale(matrix_f, 1/interval_a)
242
243 ! compute chebyshev matrix recursion
244 CALL dbcsr_get_info(matrix=matrix_f, nfullrows_total=nrows) !get information about a matrix
245 CALL dbcsr_set(matrix_dummy1, 0.0_dp) !empty matrix creation(for density matrix)
246
247 DO iwindow = 1, nwindow
248 CALL dbcsr_set(matrix_dummy2(iwindow), 0.0_dp) !empty matrix creation(for density matrix)
249 END DO
250
251 ALLOCATE (mu(1:ncheb))
252 ALLOCATE (kernel_g(1:ncheb))
253 CALL kernel(kernel_g(1), 1, ncheb)
254 CALL kernel(kernel_g(2), 2, ncheb)
255
256 CALL dbcsr_set(matrix_tmp1, 0.0_dp) !matrix creation
257 CALL dbcsr_add_on_diag(matrix_tmp1, 1.0_dp) !add a only number to diagonal elements
258 CALL dbcsr_trace(matrix_tmp1, trace=mu(1))
259 CALL dbcsr_copy(matrix_tmp2, matrix_f) !make matrix_tmp2 = matrix_F
260 CALL dbcsr_trace(matrix_tmp2, trace=mu(2))
261
262 DO iwindow = 1, nwindow
263 CALL dbcsr_copy(matrix_dummy1, matrix_tmp1)
264 CALL dbcsr_copy(matrix_dummy2(iwindow), matrix_tmp2) !matrix_dummy2=
265 CALL dbcsr_scale(matrix_dummy1, kernel_g(1)*aitchev_t(1, iwindow)) !first term of chebyshev poly(matrix)
266 CALL dbcsr_scale(matrix_dummy2(iwindow), 2.0_dp*kernel_g(2)*aitchev_t(2, iwindow)) !second term of chebyshev poly(matrix)
267
268 CALL dbcsr_add(matrix_dummy2(iwindow), matrix_dummy1, 1.0_dp, 1.0_dp)
269 END DO
270
271 DO icheb = 2, ncheb - 1
272 t1 = m_walltime()
273 CALL dbcsr_multiply("N", "N", 2.0_dp, matrix_f, matrix_tmp2, &
274 -1.0_dp, matrix_tmp1, filter_eps=ls_scf_env%eps_filter) !matrix multiplication(Recursion)
275 CALL dbcsr_copy(matrix_tmp3, matrix_tmp1)
276 CALL dbcsr_copy(matrix_tmp1, matrix_tmp2)
277 CALL dbcsr_copy(matrix_tmp2, matrix_tmp3)
278 CALL dbcsr_trace(matrix_tmp2, trace=mu(icheb + 1)) !icheb+1 th coefficient
279 CALL kernel(kernel_g(icheb + 1), icheb + 1, ncheb)
280
281 DO iwindow = 1, nwindow
282
283 CALL dbcsr_copy(matrix_dummy1, matrix_tmp2)
284 CALL dbcsr_scale(matrix_dummy1, 2.0_dp*kernel_g(icheb + 1)*aitchev_t(icheb + 1, iwindow)) !second term of chebyshev poly(matrix)
285 CALL dbcsr_add(matrix_dummy2(iwindow), matrix_dummy1, 1.0_dp, 1.0_dp)
286 CALL dbcsr_trace(matrix_dummy2(iwindow), trace=trace_dm(iwindow)) !icheb+1 th coefficient
287
288 END DO
289
290 occ = dbcsr_get_occupation(matrix_tmp1)
291 t2 = m_walltime()
292 IF (unit_nr > 0 .AND. mod(icheb, 20) == 0) THEN
293 CALL m_flush(unit_nr)
294 IF (nwindow > 0) THEN
295 WRITE (unit_nr, '(T2,A,I5,1X,A,1X,F8.3,1X,A,1X,F8.6,1X,A,1X,1000F16.8)') &
296 "Iter.", icheb, "time=", t2 - t1, "occ=", occ, "traces=", trace_dm(:)
297 ELSE
298 WRITE (unit_nr, '(T2,A,I5,1X,A,1X,F8.3,1X,A,1X,F8.6)') &
299 "Iter.", icheb, "time=", t2 - t1, "occ=", occ
300 END IF
301 END IF
302 END DO
303
304 DO iwindow = 1, nwindow
305 IF (SIZE(ls_scf_env%matrix_ks) == 1) THEN
306 orbital_occ = 2.0_dp
307 ELSE
308 orbital_occ = 1.0_dp
309 END IF
310 CALL dbcsr_multiply("N", "N", 1.0_dp, ls_scf_env%matrix_s_sqrt_inv, matrix_dummy2(iwindow), &
311 0.0_dp, matrix_tmp1, filter_eps=ls_scf_env%eps_filter)
312 CALL dbcsr_multiply("N", "N", orbital_occ, matrix_tmp1, ls_scf_env%matrix_s_sqrt_inv, &
313 0.0_dp, matrix_tmp2, filter_eps=ls_scf_env%eps_filter)
314 CALL dbcsr_copy(matrix_dummy2(iwindow), matrix_tmp2)
315
316 ! look at the difference with the density matrix from the ls routines
317 IF (.false.) THEN
318 CALL dbcsr_copy(matrix_tmp1, matrix_tmp2)
319 CALL dbcsr_add(matrix_tmp1, ls_scf_env%matrix_p(ispin), 1.0_dp, -1.0_dp) !comparison
320 frob_matrix = dbcsr_frobenius_norm(matrix_tmp1)
321 IF (unit_nr > 0) WRITE (unit_nr, *) "Difference between Chebyshev DM and LS DM", frob_matrix
322 END IF
323 END DO
324
325 write_cubes = btest(cp_print_key_should_output(logger%iter_info, &
326 ls_scf_env%chebyshev%print_key_cube), cp_p_file)
327 IF (write_cubes) THEN
328 DO iwindow = 1, nwindow
329 WRITE (middle_name, "(A,I0)") "E_DENSITY_WINDOW_", iwindow
330 WRITE (title, "(A,1X,F16.8,1X,A,1X,F16.8)") "Energy range : ", ev1(iwindow), "to", ev2(iwindow)
331 unit_cube = cp_print_key_unit_nr(logger, ls_scf_env%chebyshev%print_key_cube, &
332 "", extension=".cube", & !added 01/22/2012
333 middle_name=trim(middle_name), log_filename=.false.)
334 CALL write_matrix_to_cube(qs_env, ls_scf_env, matrix_dummy2(iwindow), unit_cube, title, &
335 section_get_ivals(ls_scf_env%chebyshev%print_key_cube, "STRIDE"))
336 CALL cp_print_key_finished_output(unit_cube, logger, ls_scf_env%chebyshev%print_key_cube, "")
337 END DO
338 END IF
339
340 END DO
341
342 ! Chebyshev expansion with calculated coefficient
343 ! grid construction and rescaling (by J)
344 unit_dos = cp_print_key_unit_nr(logger, ls_scf_env%chebyshev%print_key_dos, "", extension=".xy", &
345 middle_name="DOS", log_filename=.false.)
346
347 IF (unit_dos > 0) THEN
348 ALLOCATE (dos(1:n_gridpoint_dos))
349 ALLOCATE (gdensity(1:n_gridpoint_dos, 1:nwindow))
350 ALLOCATE (chev_e(1:n_gridpoint_dos))
351 ALLOCATE (chev_es_dos(1:n_gridpoint_dos))
352 ALLOCATE (dummy2(1:nwindow))
353 DO igrid = 1, n_gridpoint_dos
354 chev_e(igrid) = min_ev + (igrid - 1)*(max_ev - min_ev)/(n_gridpoint_dos - 1)
355 chev_es_dos(igrid) = (chev_e(igrid) - interval_b)/interval_a
356 END DO
357 DO igrid = 1, n_gridpoint_dos
358 dummy1 = 0.0_dp !summation of polynomials
359 dummy2(:) = 0.0_dp !summation of polynomials
360 DO icheb = 2, ncheb
361 CALL chebyshev_poly(chev_t_dos, chev_es_dos(igrid), icheb)
362 dummy1 = dummy1 + kernel_g(icheb)*mu(icheb)*chev_t_dos
363 DO iwindow = 1, nwindow
364 dummy2(iwindow) = dummy2(iwindow) + kernel_g(icheb)*aitchev_t(icheb, iwindow)*chev_t_dos
365 END DO
366 END DO
367 dos(igrid) = 1.0_dp/(interval_a*nrows* &
368 (pi*sqrt(1.0_dp - chev_es_dos(igrid)*chev_es_dos(igrid))))*(kernel_g(1)*mu(1) + 2.0_dp*dummy1)
369 DO iwindow = 1, nwindow
370 gdensity(igrid, iwindow) = kernel_g(1)*aitchev_t(1, iwindow) + 2.0_dp*dummy2(iwindow)
371 END DO
372 WRITE (unit_dos, '(1000F16.8)') chev_e(igrid), dos(igrid), gdensity(igrid, :)
373 END DO
374 DEALLOCATE (chev_es_dos, chev_e, dos, gdensity)
375 END IF
376 CALL cp_print_key_finished_output(unit_dos, logger, ls_scf_env%chebyshev%print_key_dos, "")
377
378 ! free the matrices
379 CALL dbcsr_release(matrix_tmp1)
380 CALL dbcsr_release(matrix_tmp2)
381 CALL dbcsr_release(matrix_tmp3)
382 CALL dbcsr_release(matrix_f)
383 CALL dbcsr_release(matrix_dummy1)
384
385 DO iwindow = 1, nwindow
386 CALL dbcsr_release(matrix_dummy2(iwindow))
387 END DO
388
389 DEALLOCATE (ev1, ev2, sev1, sev2, matrix_dummy2)
390
391 !Need deallocation
392 DEALLOCATE (mu, kernel_g, aitchev_t, e_inte, sqrt_vec)
393
394 IF (unit_nr > 0) WRITE (unit_nr, '(T2,A)') "ENDING CHEBYSHEV CALCULATION"
395
396 CALL timestop(handle)
397
398 END SUBROUTINE compute_chebyshev
399
400END MODULE dm_ls_chebyshev
arnoldi iteration using dbcsr
Definition arnoldi_api.F:15
subroutine, public arnoldi_extremal(matrix_a, max_ev, min_ev, converged, threshold, max_iter)
simple wrapper to estimate extremal eigenvalues with arnoldi, using the old lanczos interface this hi...
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
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...
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)
...
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 cp_p_file
integer function, public cp_print_key_should_output(iteration_info, basis_section, print_key_path, used_print_key, first_time)
returns what should be done with the given property if btest(res,cp_p_store) then the property should...
Routines using linear scaling chebyshev methods.
subroutine, public compute_chebyshev(qs_env, ls_scf_env)
compute properties based on chebyshev expansion
Routines for a linear scaling quickstep SCF run based on the density matrix, with a focus on the inte...
subroutine, public write_matrix_to_cube(qs_env, ls_scf_env, matrix_p_ls, unit_nr, title, stride)
...
Types needed for a linear scaling quickstep SCF run based on the density matrix.
objects that represent the structure of input sections and the data contained in an input section
integer function, dimension(:), pointer, public section_get_ivals(section_vals, keyword_name)
...
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
Defines the basic variable types.
Definition kinds.F:23
integer, parameter, public dp
Definition kinds.F:34
integer, parameter, public default_string_length
Definition kinds.F:57
Machine interface based on Fortran 2003 and POSIX.
Definition machine.F:17
subroutine, public m_flush(lunit)
flushes units if the &GLOBAL flag is set accordingly
Definition machine.F:106
real(kind=dp) function, public m_walltime()
returns time from a real-time clock, protected against rolling early/easily
Definition machine.F:123
Definition of mathematical constants and functions.
real(kind=dp), parameter, public pi
type of a logger, at the moment it contains just a print level starting at which level it should be l...