34#include "../base/base_uses.f90"
38 CHARACTER(LEN=80) :: file, out1, out2, out3, wq_char, &
39 path_file, out3_stride
40 CHARACTER(LEN=480) :: a, b
41 CHARACTER(LEN=default_string_length) :: active_label, per_label
42 INTEGER :: istat, coor, i, id, ip, &
43 it, iw, ix, j, ncount, ndim, &
44 ndw, nf, nfes, ngauss, nh, &
45 nprd, nt, nt_p, nwr, p, q, &
46 stat, unit_nr, unit_nr2
47 INTEGER,
POINTER :: i_map(:), idw(:), ind(:), &
48 inds(:), iperd(:), iprd(:), &
49 ngrid(:), nn(:, :), nn_max(:), &
51 LOGICAL :: fix, l_cp2k, l_cpmd, &
52 l_dp, l_fes_int, l_fmin, &
53 l_grid, l_math, l_orac, &
54 l_pmin, lstride, l_popt, l_int, &
56 REAL(kind=
dp) :: delta_s_glob, diff, dp2, dum, &
57 eps_cut, sc, ss, x0w(3), &
59 REAL(kind=
dp),
POINTER :: delta_s(:, :), dp_cut(:), dp_grid(:), fes(:), &
60 gauss(:, :), ss0(:, :), tmpr(:), ww(:), x0(:), xf(:)
88 IF (command_argument_count() == 0)
THEN
89 WRITE (iw, *)
'USAGE:'
91 WRITE (iw, *)
'[-ngrid 50 .. ..] (Mesh dimension. Default :: 100)'
92 WRITE (iw, *)
'[-dp 0.05 .. ..] (Alternative to -ngrid, allows the specification of the mesh dx)'
93 WRITE (iw, *)
'[-ndim 3 ] (Number of collective variables NCV)'
94 WRITE (iw, *)
'[-ndw 1 3 .. ] (CVs for the free energy surface)'
95 WRITE (iw, *)
'[-periodic 2 3 ..] (CVs with periodic boundary conditions (-pi,pi] )'
96 WRITE (iw, *)
'[-stride 10 ] (How often the FES is written)'
97 WRITE (iw, *)
'[-fix 1.1 .. ..] (Define the region for the FES)'
98 WRITE (iw, *)
' (If omitted this is automatically calculated)'
99 WRITE (iw, *)
'[-cutoff 2. ] (The hills are cutoffed at 2)'
100 WRITE (iw, *)
'[-file filename]'
101 WRITE (iw, *)
'[-out filename]'
102 WRITE (iw, *)
'[-integrated_fes] (When projecting the FES print the integrated value, '
103 WRITE (iw, *)
' rather then the minimum value (minimum value is default))'
104 WRITE (iw, *)
'[-orac] (If energies are written in orac intern units)'
105 WRITE (iw, *)
'[-cp2k] (Specify if a CP2K restart file is provided)'
106 WRITE (iw, *)
'[-cpmd] (Specify if CPMD colvar_mtd and parvar_mtd are provided)'
107 WRITE (iw, *)
' (With CPMD you do not need to specify -file, parvar_mtd and'
108 WRITE (iw, *)
' colvar_mtd are expected to be present in the working directory)'
109 WRITE (iw, *)
'[-mathlab] (File storing FES in Mathlab format. Default format Gnuplot)'
110 WRITE (iw, *)
'[-cube] (File storing FES in GAUSSIAN CUBE format. Default format Gnuplot)'
111 WRITE (iw, *)
'[-find-minima] (Tries to finds all minima in the computed FES)'
112 WRITE (iw, *)
'[-find-path] (Finds MEP between all minima (found) in the computed FES)'
113 WRITE (iw, *)
'[-point-a] (Specifies point (a) when using -find-path option)'
114 WRITE (iw, *)
'[-point-b] (Specifies point (b) when using -find-path option)'
115 WRITE (iw, *)
'[-opt-path filename] (Optimize initial MEP of mep-nreplica points in the same format as mep.data)'
116 WRITE (iw, *)
'[-mep-kb] (Specifies the value of the force constant for the MEP: default 0.1_dp)'
117 WRITE (iw, *)
'[-mep-nreplica] (Specifies the number of replica points used in the MEP: default 8)'
118 WRITE (iw, *)
'[-mep-iter] (Specifies the maximum number of iterations used in the MEP: default 10000)'
120 WRITE (iw, *)
'DEFAULT OUTPUT: fes.dat'
122 cpabort(
"Please provide arguments to run FES!")
125 DO i = 1, command_argument_count()
126 CALL get_command_argument(i, wq_char, status=istat)
129 IF (index(wq_char,
'-file') .NE. 0)
THEN
130 CALL get_command_argument(i + 1, wq_char, status=istat)
132 READ (wq_char, *) file
135 IF (index(wq_char,
'-out') .NE. 0)
THEN
136 CALL get_command_argument(i + 1, wq_char, status=istat)
138 READ (wq_char, *) out1
144 IF (index(wq_char,
'-ndim') .NE. 0)
THEN
145 CALL get_command_argument(i + 1, wq_char, status=istat)
147 READ (wq_char, *) ndim
150 IF (index(wq_char,
'-stride') .NE. 0)
THEN
151 CALL get_command_argument(i + 1, wq_char, status=istat)
153 READ (wq_char, *) nt_p
157 IF (index(wq_char,
'-cutoff') .NE. 0)
THEN
158 CALL get_command_argument(i + 1, wq_char, status=istat)
160 READ (wq_char, *) eps_cut
163 IF (index(wq_char,
'-integrated_fes') .NE. 0)
THEN
167 IF (index(wq_char,
'-orac') .NE. 0)
THEN
171 IF (index(wq_char,
'-cp2k') .NE. 0)
THEN
175 IF (index(wq_char,
'-cpmd') .NE. 0)
THEN
179 IF (index(wq_char,
'-find-minima') .NE. 0)
THEN
183 IF (index(wq_char,
'-find-path') .NE. 0)
THEN
187 IF (index(wq_char,
'-mathlab') .NE. 0)
THEN
191 IF (index(wq_char,
'-cube') .NE. 0)
THEN
195 IF (index(wq_char,
'-opt-path') .NE. 0)
THEN
197 CALL get_command_argument(i + 1, wq_char, status=istat)
199 READ (wq_char, *) path_file
203 IF (count((/l_orac, l_cp2k, l_cpmd/)) /= 1) &
204 cpabort(
"Error! You've to specify either ORAC, CP2K or CPMD!")
213 CALL random_number(dum)
217 ALLOCATE (ngrid(ndim), stat=stat)
218 IF (stat /= 0) cpabort(
"Allocation Error")
219 ALLOCATE (dp_grid(ndim), stat=stat)
220 IF (stat /= 0) cpabort(
"Allocation Error")
221 ALLOCATE (idw(ndw), stat=stat)
222 IF (stat /= 0) cpabort(
"Allocation Error")
223 ALLOCATE (iperd(ndim), stat=stat)
224 IF (stat /= 0) cpabort(
"Allocation Error")
225 ALLOCATE (iprd(nprd), stat=stat)
226 IF (stat /= 0) cpabort(
"Allocation Error")
232 DO i = 1, command_argument_count()
233 CALL get_command_argument(i, wq_char, status=istat)
236 IF (index(wq_char,
'-ndw') .NE. 0)
THEN
240 ndw_loop:
DO ix = i + 1, command_argument_count()
241 CALL get_command_argument(ix, wq_char, status=istat)
243 IF (index(wq_char,
'-') .EQ. 0)
THEN
250 ALLOCATE (idw(ndw), stat=stat)
251 IF (stat /= 0) cpabort(
"Allocation Error")
254 CALL get_command_argument(i + id, wq_char, status=istat)
256 READ (wq_char, *) idw(id)
260 IF (index(wq_char,
'-periodic') .NE. 0)
THEN
262 nprd_loop:
DO ix = i + 1, command_argument_count()
263 CALL get_command_argument(ix, wq_char, status=istat)
265 IF (index(wq_char,
'-') .EQ. 0)
THEN
273 ALLOCATE (iprd(nprd), stat=stat)
274 IF (stat /= 0) cpabort(
"Allocation Error")
277 CALL get_command_argument(i + id, wq_char, status=istat)
279 READ (wq_char, *) iprd(id)
283 IF (index(wq_char,
'-ngrid') .NE. 0)
THEN
285 CALL get_command_argument(i + ix, wq_char, status=istat)
287 READ (wq_char, *) ngrid(ix)
292 IF (index(wq_char,
'-dp') .NE. 0)
THEN
296 CALL get_command_argument(i + ix, wq_char, status=istat)
298 READ (wq_char, *) dp_grid(ix)
302 IF (index(wq_char,
'-fix') .NE. 0)
THEN
305 CALL get_command_argument(i + 2*(id - 1) + 1, wq_char, status=istat)
307 READ (wq_char, *) x0w(id)
308 CALL get_command_argument(i + 2*(id - 1) + 2, wq_char, status=istat)
310 READ (wq_char, *) xfw(id)
316 ALLOCATE (mep_input_data%minima(ndw, 2))
317 mep_input_data%minima = huge(0.0_dp)
318 mep_input_data%max_iter = 10000
319 mep_input_data%kb = 0.1_dp
320 mep_input_data%nreplica = 8
322 DO i = 1, command_argument_count()
323 CALL get_command_argument(i, wq_char, status=istat)
326 IF (index(wq_char,
'-point-a') .NE. 0)
THEN
328 CALL get_command_argument(i + id, wq_char, status=istat)
330 READ (wq_char, *) mep_input_data%minima(id, 1)
334 IF (index(wq_char,
'-point-b') .NE. 0)
THEN
336 CALL get_command_argument(i + id, wq_char, status=istat)
338 READ (wq_char, *) mep_input_data%minima(id, 2)
342 IF (index(wq_char,
'-mep-iter') .NE. 0)
THEN
343 CALL get_command_argument(i + 1, wq_char, status=istat)
345 READ (wq_char, *) mep_input_data%max_iter
348 IF (index(wq_char,
'-mep-kb') .NE. 0)
THEN
349 CALL get_command_argument(i + 1, wq_char, status=istat)
351 READ (wq_char, *) mep_input_data%kb
354 IF (index(wq_char,
'-mep-nreplica') .NE. 0)
THEN
355 CALL get_command_argument(i + 1, wq_char, status=istat)
357 READ (wq_char, *) mep_input_data%nreplica
361 IF (any(mep_input_data%minima == huge(0.0_dp))) &
362 cpabort(
"-find-path requires the specification of -point-a and -point-b !")
364 ALLOCATE (mep_input_data%minima(0, 0))
369 mep_input_data%nreplica = 0
370 mep_input_data%max_iter = 10000
371 mep_input_data%kb = 0.1_dp
373 DO i = 1, command_argument_count()
374 CALL get_command_argument(i, wq_char, status=istat)
377 IF (index(wq_char,
'-mep-kb') .NE. 0)
THEN
378 CALL get_command_argument(i + 1, wq_char, status=istat)
380 READ (wq_char, *) mep_input_data%kb
383 IF (index(wq_char,
'-mep-iter') .NE. 0)
THEN
384 CALL get_command_argument(i + 1, wq_char, status=istat)
386 READ (wq_char, *) mep_input_data%max_iter
389 IF (index(wq_char,
'-mep-nreplica') .NE. 0)
THEN
390 CALL get_command_argument(i + 1, wq_char, status=istat)
392 READ (wq_char, *) mep_input_data%nreplica
396 ALLOCATE (mep_input_data%minima(ndw, mep_input_data%nreplica))
398 CALL open_file(unit_number=unit_nr, file_name=path_file, file_status=
"OLD")
399 DO id = 1, mep_input_data%nreplica
400 READ (unit_nr, *) j, mep_input_data%minima(:, id)
404 DO id = 1, mep_input_data%nreplica
405 WRITE (*, *) mep_input_data%minima(:, id)
410 ALLOCATE (i_map(ndim), stat=stat)
411 IF (stat /= 0) cpabort(
"Allocation Error")
419 IF (i_map(id) == 0)
THEN
427 i_map = ndim - i_map + 1
431 iperd(i_map(iprd(id))) = 1
436 ALLOCATE (tmp(ndim), stat=stat)
437 IF (stat /= 0) cpabort(
"Allocation Error")
440 ngrid(i_map(i)) = tmp(i)
447 WRITE (iw,
'(/,70("*"))')
448 WRITE (iw,
'("FES|",T7,A,/)')
"Parsing file: <"//trim(file)//
">"
450 CALL open_file(unit_number=unit_nr, file_name=file, file_status=
"OLD")
452 CALL get_val_res(unit=unit_nr, section=
"&METADYN", keyword=
"NHILLS_START_VAL", i_val=nt)
455 CALL get_val_res(unit=unit_nr, section=
"&METADYN", keyword=
"HILL_TAIL_CUTOFF", r_val=sc)
456 CALL get_val_res(unit=unit_nr, section=
"&METADYN", keyword=
"P_EXPONENT", i_val=p)
457 CALL get_val_res(unit=unit_nr, section=
"&METADYN", keyword=
"Q_EXPONENT", i_val=q)
458 ELSE IF (l_orac .OR. l_cpmd)
THEN
461 READ (unit_nr, *,
END=100, ERR=100) dum
467 ALLOCATE (x0(ndim), stat=stat)
468 IF (stat /= 0) cpabort(
"Allocation Error")
469 ALLOCATE (xf(ndim), stat=stat)
470 IF (stat /= 0) cpabort(
"Allocation Error")
471 ALLOCATE (ss0(ndim, nt), stat=stat)
472 IF (stat /= 0) cpabort(
"Allocation Error")
473 ALLOCATE (delta_s(ndim, nt), stat=stat)
474 IF (stat /= 0) cpabort(
"Allocation Error")
475 ALLOCATE (ww(nt), stat=stat)
476 IF (stat /= 0) cpabort(
"Allocation Error")
477 ALLOCATE (ind(ndim), stat=stat)
478 IF (stat /= 0) cpabort(
"Allocation Error")
479 ALLOCATE (inds(ndim), stat=stat)
480 IF (stat /= 0) cpabort(
"Allocation Error")
481 ALLOCATE (nn(ndim, nt), stat=stat)
482 IF (stat /= 0) cpabort(
"Allocation Error")
483 ALLOCATE (nn_max(ndim), stat=stat)
484 IF (stat /= 0) cpabort(
"Allocation Error")
485 ALLOCATE (dp_cut(ndim), stat=stat)
486 IF (stat /= 0) cpabort(
"Allocation Error")
489 CALL get_val_res(unit=unit_nr, section=
"&METADYN", subsection=
"&SPAWNED_HILLS_POS")
491 READ (unit_nr,
'(A120)') a
492 DO WHILE (a(len_trim(a):len_trim(a)) ==
'\\')
493 READ (unit_nr,
'(A120)') b
494 a = a(1:len_trim(a) - 1)//b(1:len_trim(b))
496 READ (a, *) (ss0(i_map(id), i), id=1, ndim)
498 CALL get_val_res(unit=unit_nr, section=
"&METADYN", subsection=
"&SPAWNED_HILLS_SCALE")
500 READ (unit_nr,
'(A120)') a
501 DO WHILE (a(len_trim(a):len_trim(a)) ==
'\\')
502 READ (unit_nr,
'(A120)') b
503 a = a(1:len_trim(a) - 1)//b(1:len_trim(b))
505 READ (a, *) (delta_s(i_map(id), i), id=1, ndim)
507 CALL get_val_res(unit=unit_nr, section=
"&METADYN", subsection=
"&SPAWNED_HILLS_HEIGHT")
509 READ (unit_nr, *) ww(i)
511 ELSE IF (l_orac)
THEN
513 READ (unit_nr, *) dum, (ss0(i_map(id), i), id=1, ndim), (delta_s(i_map(id), i), id=1, ndim), ww(i)
515 ELSE IF (l_cpmd)
THEN
516 CALL open_file(unit_number=unit_nr2, file_name=
"parvar_mtd", file_status=
"OLD")
518 READ (unit_nr, *) dum, (ss0(i_map(id), i), id=1, ndim), (delta_s(id, i), id=1, ndim)
519 READ (unit_nr2, *) dum, dum, delta_s_glob, ww(i)
520 delta_s(1:ndim, i) = delta_s_glob*delta_s(1:ndim, i)
522 CALL close_file(unit_nr2)
524 CALL close_file(unit_nr)
527 IF (l_orac) ww = ww*10000._dp/4.187_dp
531 x0(id) = huge(1.0_dp)
532 xf(id) = -huge(1.0_dp)
536 DO id = 1, ndim - ndw
537 x0(id) = min(x0(id), ss0(id, it) - 3.*delta_s(id, it))
538 xf(id) = max(xf(id), ss0(id, it) + 3.*delta_s(id, it))
542 DO id = ndim, ndim - ndw + 1, -1
550 IF (iperd(id) == 1)
THEN
554 x0(id) = min(x0(id), ss0(id, it) - 3.*delta_s(id, it))
555 xf(id) = max(xf(id), ss0(id, it) + 3.*delta_s(id, it))
562 ALLOCATE (tmpr(ndim))
565 dp_grid(i_map(i)) = tmpr(i)
568 ngrid = int((xf - x0)/dp_grid) + 1
570 dp_grid = (xf - x0)/real(ngrid - 1, kind=dp)
573 WRITE (iw,
'(70("*"))')
574 WRITE (iw,
'("FES|",T7,A,/)')
"Parameters for FES:"
575 WRITE (iw,
'("FES|",T7,A15,5x,i7)')
"NDIM ::", ndim
576 WRITE (iw,
'("FES|",T7,A15,5x,i7)')
"NWD ::", ndw
577 WRITE (iw,
'("FES|",T7,A15,5x,i7)')
"HILLS ::", nt
582 active_label =
"(NO MAPPED)"
583 IF (iperd(i) /= 0) per_label =
"(PERIODIC)"
584 IF (it <= ndw) active_label =
"( MAPPED)"
585 j = minloc((i_map - i)**2, 1)
586 WRITE (iw,
'("FES|",T7,"COLVAR # ",i3," ::",5x,"(",f7.3," ,",f7.3,")",T48,A,T60,A)') &
587 j, x0(i), xf(i), trim(per_label), trim(active_label)
589 WRITE (iw,
'("FES|",T7,a15,5x,7i7)')
"NGRID ::", (ngrid(id), id=ndim, ndim - ndw + 1, -1)
590 WRITE (iw,
'("FES|",T7,a15,5x,5f7.3)')
"DX ::", (dp_grid(id), id=ndim, ndim - ndw + 1, -1)
591 WRITE (iw,
'("FES|",T7,a15,5x,g10.5)')
"CUTOFF ::", eps_cut
592 WRITE (iw,
'(70("*"),/)')
596 dp_cut = sqrt(log(abs(ww(i))/eps_cut))*2.0_dp*delta_s(:, i)
597 nn(:, i) = int(dp_cut/dp_grid)
598 ww(i) = ww(i)**(1.0_dp/real(ndim, kind=dp))
601 nn_max = maxval(nn, dim=2)
602 ngauss = maxval(nn_max)*2 + 1
603 nfes = product(ngrid)
605 ALLOCATE (gauss(-maxval(nn_max):maxval(nn_max), ndim))
610 nf = min(nh + nt_p - 1, nt)
619 WRITE (iw,
'(/,"FES|",T7,A)')
"Computing Free Energy Surface"
621 stride:
DO WHILE (nh <= nt)
622 hills:
DO it = nh, nf
623 ind = int((ss0(:, it) - x0)/dp_grid) + 1
627 coor = ind(i) - nn(i, it) - 1
628 ss = x0(i) + coor*dp_grid(i) - dp_grid(i)
629 DO ip = -nn(i, it), nn(i, it)
632 IF (iperd(i) == 0)
THEN
633 IF (coor .GT. ngrid(i)) cycle
634 IF (coor .LT. 1) cycle
636 diff = ss - ss0(i, it)
637 dp2 = (diff/delta_s(i, it))**2
638 gauss(ip, i) = ww(it)*exp(-0.5_dp*dp2)
639 IF (sc > 0.0_dp .AND. p > 0.0_dp .AND. q > 0.0_dp .AND. q > p)
THEN
640 gauss(ip, i) = gauss(ip, i)*(1 - (diff/sc*delta_s(i, it))**p)/(1 - (diff/sc*delta_s(i, it))**q)
645 CALL fes_compute_low(ndim, nn(:, it), fes, gauss, ind, inds, nfes, ndim, ngauss, ngrid, iperd)
647 IF (.NOT. lstride .AND. mod(it, nwr) == 0)
THEN
648 WRITE (iw,
'("FES|",T7,a,i4,a2)')
"Mapping Gaussians ::", int(10*anint(10.*it/nt)),
" %"
649 ELSEIF (.NOT. lstride .AND. it == nt)
THEN
650 WRITE (iw,
'("FES|",T7,a,i4,a2)')
"Mapping Gaussians ::", int(10*anint(10.*it/nt)),
" %"
656 WRITE (iw,
'("FES|",T7,a13,i5," |-| Gaussians from ",i6," to",i6)')
"Done frame ::", ncount, nh, nf
663 IF (ncount < 10)
THEN
664 WRITE (out3_stride,
'(A,i1)') trim(out3), ncount
665 ELSEIF (ncount < 100)
THEN
666 WRITE (out3_stride,
'(A,i2)') trim(out3), ncount
668 WRITE (out3_stride,
'(A,i3)') trim(out3), ncount
670 CALL open_file(unit_number=unit_nr, file_name=out3_stride, file_action=
"WRITE", file_status=
"UNKNOWN", file_form=
"FORMATTED")
672 CALL fes_only_write(ndim, fes, ind, ndim, ngrid, dp_grid, ndw, l_fes_int, unit_nr)
673 CALL close_file(unit_nr)
677 nf = min(nh + nt_p - 1, nt)
687 WRITE (iw,
'("FES|",T7,A)')
"Dumping FES structure in file: < "//trim(out3)//
" >"
688 CALL open_file(unit_number=unit_nr, file_name=out3, file_action=
"WRITE", file_status=
"UNKNOWN", file_form=
"FORMATTED")
691 CALL fes_cube_write(ndim, fes, ind, ndim, ngrid, dp_grid, x0, ndw, l_fes_int, file)
694 IF (l_math)
WRITE (unit_nr,
'(10g12.5)') (ngrid(id), id=ndim, ndim - ndw + 1, -1), ix
696 CALL fes_write(unit_nr, ndim, fes, ind, ndim, ngrid, dp_grid, x0, ndw, l_fes_int)
698 CALL close_file(unit_nr)
701 IF (l_fmin)
CALL fes_min(fes, ndim, iperd, ngrid, dp_grid, x0, ndw)
704 IF ((l_pmin) .AND. (l_popt)) cpabort(
"USE EITHER -find-path OR -opt-path")
705 IF (l_pmin) l_int = .true.
706 IF (l_popt) l_int = .false.
708 IF ((l_pmin) .OR. (l_popt))
CALL fes_path(fes, ndim, ngrid, dp_grid, iperd, x0, ndw, mep_input_data, l_int)
728 DEALLOCATE (mep_input_data%minima)
731 WRITE (iw,
'(/,A,/)')
"FES| NORMAL FES TERMINATION."
program graph
Program to Map on grid the hills spawned during a metadynamics run.
Utility routines to open and close files. Tracking of preconnections.
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.
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.
Program to Map on grid the hills spawned during a metadynamics run.
subroutine, public fes_min(fes, ndim, iperd, ngrid, dp_grid, x0, ndw)
Finds minima of the FES.
recursive subroutine, public fes_compute_low(idim, nn, fes, gauss, ind, ind0, nfes, ndim, ngauss, ngrid, iperd)
Efficiently map the gaussians on the grid.
subroutine, public fes_cube_write(idim, fes, pos, ndim, ngrid, dp_grid, x0, ndw, l_fes_int, file)
Dump FES with a GAUSSIAN cube format - Useful for multidimensional FES.
recursive subroutine, public fes_write(unit_nr, idim, fes, pos, ndim, ngrid, dp_grid, x0, ndw, l_fes_int, array)
Writes the FES on the file.
recursive subroutine, public fes_only_write(idim, fes, pos, ndim, ngrid, dp_grid, ndw, l_fes_int, unit_nr)
Writes the FES on the file when stride is requested.
subroutine, public fes_path(fes, ndim, ngrid, dp_grid, iperd, x0, ndw, mep_input_data, l_int)
Finds path between two points (a) and (b)
Module containing utils for mapping FESs.
subroutine, public get_val_res(unit, section, keyword, subsection, i_val, r_val)
Parser informations from the cp2k input/restart.
Defines the basic variable types.
integer, parameter, public dp
integer, parameter, public default_string_length
Definition of mathematical constants and functions.
real(kind=dp), parameter, public pi