forked from m3g/packmol
-
Notifications
You must be signed in to change notification settings - Fork 0
/
checkpoint.f90
116 lines (98 loc) · 3.51 KB
/
checkpoint.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
!
! Written by Leandro Martínez, 2009-2011.
! Copyright (c) 2009-2018, Leandro Martínez, Jose Mario Martinez,
! Ernesto G. Birgin.
!
!
! Subroutine that writes the last point obtained when
! a solution was not found
!
subroutine checkpoint(n,x)
use sizes
use compute_data
use input
use usegencan
use ahestetic
implicit none
integer :: i
integer :: n
double precision :: x(n)
double precision :: fx
logical :: movebadprint
character(len=strl) :: xyzout_forced
! All molecules are important
do i = 1, ntfix
comptype(i) = .true.
end do
! Call the subroutine that computes de function value
call computef(n,x,fx)
write(*,dash3_line)
write(*,"(&
&' Packmol was not able to find a solution to your',/,&
&' packing problem with the desired distance tolerance.',/,/,&
&' First of all, be sure if the molecules fit in the',/,&
&' regions specified and if the constraints were set',/,&
&' correctly. ',/,/,&
&' Secondly, try simply running it again with a different ',/,&
&' seed for the random number generator of the initial ',/,&
&' point. This is done by adding the keyword seed to the',/,&
&' input file, as in: ',/,/,&
&' seed 192911 ',/,/,&
&' The best configuration found has a function value of',/,&
&' f = ', e14.7,/,/,&
&' IMPORTANT: ',/,&
&' If the number of molecules and the restraints are',/,&
&' correct, it is still very likely that the current point',/,&
&' fits your needs if your purpose is to run a MD',/,&
&' simulation.',/,&
&' Therefore, we recommend to minimize the energy of the',/,&
&' solution found, equilibrate it and run with it as well.',/&
&)") fx
write(*,dash3_line)
call output(n,x)
write(*,*) ' The solution with the best function value was '
write(*,*) ' written to the output file: ', trim(adjustl(xyzout))
if ( crd ) write(*,*) ' ... and to CRD file: ', trim(adjustl(crdfile))
write(*,dash1_line)
write(*,*) ' Forcing the solution to fit the constraints...'
! CALL GENCAN
init1 = .true.
do i = 1, nloop
iprint1 = 0
iprint2 = 0
call pgencan(n,x,fx)
movebadprint = .false.
call movebad(n,x,fx,movebadprint)
end do
init1 = .false.
write(*,*)
write(*,dash1_line)
xyzout_forced = trim(adjustl(xyzout))//'_FORCED'
call output(n,x)
write(*,*) ' The forced point was writen to the '
write(*,*) ' output file: ', trim(adjustl(xyzout_forced))
if ( crd ) write(*,*) ' ... and to CRD file: ', trim(adjustl(crdfile))
write(*,*)
write(*,*) ' If you want that the packing procedure continues'
write(*,*) ' for a longer time, add the following keyword '
write(*,*) ' to the input file: '
write(*,*)
write(*,*) ' nloop [integer] (ex: nloop 200) '
write(*,*)
write(*,*) ' The default nloop value is 50 for each molecule.'
write(*,*)
write(*,hash1_line)
write(*,*) ' ENDED WITHOUT PERFECT PACKING: '
write(*,*) ' The output file:'
write(*,*)
write(*,*) ' ', trim(adjustl(xyzout))
if ( crd ) write(*,*) ' (... and to CRD file: ', trim(adjustl(crdfile)), ')'
write(*,*)
write(*,*) ' contains the best solution found. '
write(*,*)
write(*,*) ' Very likely, if the input data was correct, '
write(*,*) ' it is a reasonable starting configuration.'
write(*,*) ' Check commentaries above for more details. '
write(*,hash1_line)
return
end subroutine checkpoint