1 subroutine minimize(etot,x,iretcode,nfun)
2 implicit real*8 (a-h,o-z)
4 parameter (liv=60,lv=(77+maxvar*(maxvar+17)/2))
5 *********************************************************************
6 * OPTIMIZE sets up SUMSL or DFP and provides a simple interface for *
7 * the calling subprogram. *
8 * when d(i)=1.0, then v(35) is the length of the initial step, *
9 * calculated in the usual pythagorean way. *
10 * absolute convergence occurs when the function is within v(31) of *
11 * zero. unless you know the minimum value in advance, abs convg *
12 * is probably not useful. *
13 * relative convergence is when the model predicts that the function *
14 * will decrease by less than v(32)*abs(fun). *
15 *********************************************************************
16 include 'COMMON.IOUNITS'
19 include 'COMMON.MINIM'
22 double precision minval,x(maxvar),d(maxvar),v(1:lv),xx(maxvar)
23 double precision energia(0:n_ene)
24 external func,gradient,fdum
25 external func_restr,grad_restr
26 logical not_done,change,reduce
27 c common /przechowalnia/ v
30 C Following 3 line for diagnostics; comment out if not needed
31 write(iout,*) "Enter MINIMIZE liv",liv," lv",lv
32 write (iout,*) "Coordinates before minimization"
34 call deflt(2,iv,liv,lv,v)
35 * 12 means fresh start, dont call deflt
37 * max num of fun calls
38 if (maxfun.eq.0) maxfun=500
40 * max num of iterations
41 if (maxmin.eq.0) maxmin=1000
47 if (print_min_ini+print_min_stat+print_min_res.gt.0) iv(21)=iout
48 * 1 means to print out result
50 * 1 means to print out summary stats
52 * 1 means to print initial x and d
54 * min val for v(radfac) default is 0.1
56 * max val for v(radfac) default is 4.0
59 * check false conv if (act fnctn decrease) .lt. v(26)*(exp decrease)
60 * the sumsl default is 0.1
62 * false conv if (act fnctn decrease) .lt. v(34)
63 * the sumsl default is 100*machep
65 * absolute convergence
66 if (tolf.eq.0.0D0) tolf=1.0D-4
68 * relative convergence
69 if (rtolf.eq.0.0D0) rtolf=1.0D-4
71 * controls initial step size
73 * large vals of d correspond to small components of step
80 cd print *,'Calling SUMSL'
82 call x2xx(x,xx,nvar_restr)
83 call sumsl(nvar_restr,d,xx,func_restr,grad_restr,
84 & iv,liv,lv,v,idum,rdum,fdum)
87 call sumsl(nvar,d,x,func,gradient,iv,liv,lv,v,idum,rdum,fdum)
91 cd print *,'Exit SUMSL; return code:',iretcode,' energy:',etot
92 cd write (iout,'(/a,i4,a,f15.5/)') 'SUMSL return code:',iv(1),
94 call var_to_geom(nvar,x)
96 c call etotal(energia(0))
98 c call enerprint(energia(0))
104 c----------------------------------------------------------------------------
105 subroutine ergastulum
106 implicit real*8 (a-h,o-z)
111 include 'COMMON.SETUP'
112 include 'COMMON.DERIV'
114 include 'COMMON.IOUNITS'
115 include 'COMMON.FFIELD'
116 include 'COMMON.INTERACT'
118 include 'COMMON.TIME1'
119 double precision z(maxres6),d_a_tmp(maxres6)
120 double precision edum(0:n_ene),time_order(0:10)
121 double precision Gcopy(maxres2,maxres2)
122 common /przechowalnia/ Gcopy
124 C Workers wait for variables and NF, and NFL from the boss
126 do while (iorder.ge.0)
127 c write (*,*) 'Processor',fg_rank,' CG group',kolor,
128 c & ' receives order from Master'
130 call MPI_Bcast(iorder,1,MPI_INTEGER,king,FG_COMM,IERR)
131 time_Bcast=time_Bcast+MPI_Wtime()-time00
132 if (icall.gt.4 .and. iorder.ge.0)
133 & time_order(iorder)=time_order(iorder)+MPI_Wtime()-time00
136 c & 'Processor',fg_rank,' completed receive MPI_BCAST order',iorder
137 if (iorder.eq.0) then
140 c write (2,*) "After etotal"
141 c write (2,*) "dimen",dimen," dimen3",dimen3
143 else if (iorder.eq.2) then
145 cmd call etotal_short(edum)
146 c write (2,*) "After etotal_short"
147 c write (2,*) "dimen",dimen," dimen3",dimen3
149 else if (iorder.eq.3) then
151 cmd call etotal_long(edum)
152 c write (2,*) "After etotal_long"
153 c write (2,*) "dimen",dimen," dimen3",dimen3
155 else if (iorder.eq.1) then
157 c write (2,*) "After sum_gradient"
158 c write (2,*) "dimen",dimen," dimen3",dimen3
160 else if (iorder.eq.4) then
161 cmd call ginv_mult(z,d_a_tmp)
162 else if (iorder.eq.5) then
163 c Setup MD things for a slave
164 dimen=(nct-nnt+1)+nside
165 dimen1=(nct-nnt)+(nct-nnt+1)
167 c write (2,*) "dimen",dimen," dimen3",dimen3
169 call int_bounds(dimen,igmult_start,igmult_end)
170 igmult_start=igmult_start-1
171 call MPI_Allgather(3*igmult_start,1,MPI_INTEGER,
172 & ng_start(0),1,MPI_INTEGER,FG_COMM,IERROR)
173 my_ng_count=igmult_end-igmult_start
174 call MPI_Allgather(3*my_ng_count,1,MPI_INTEGER,ng_counts(0),1,
175 & MPI_INTEGER,FG_COMM,IERROR)
176 c write (2,*) "ng_start",(ng_start(i),i=0,nfgtasks-1)
177 c write (2,*) "ng_counts",(ng_counts(i),i=0,nfgtasks-1)
178 myginv_ng_count=maxres2*my_ng_count
179 c write (2,*) "igmult_start",igmult_start," igmult_end",
180 c & igmult_end," my_ng_count",my_ng_count
182 call MPI_Allgather(maxres2*igmult_start,1,MPI_INTEGER,
183 & nginv_start(0),1,MPI_INTEGER,FG_COMM,IERROR)
184 call MPI_Allgather(myginv_ng_count,1,MPI_INTEGER,
185 & nginv_counts(0),1,MPI_INTEGER,FG_COMM,IERROR)
186 c write (2,*) "nginv_start",(nginv_start(i),i=0,nfgtasks-1)
187 c write (2,*) "nginv_counts",(nginv_counts(i),i=0,nfgtasks-1)
189 c call MPI_Barrier(FG_COMM,IERROR)
191 call MPI_Scatterv(ginv(1,1),nginv_counts(0),
192 & nginv_start(0),MPI_DOUBLE_PRECISION,gcopy(1,1),
193 & myginv_ng_count,MPI_DOUBLE_PRECISION,king,FG_COMM,IERR)
195 time_scatter_ginv=time_scatter_ginv+MPI_Wtime()-time00
202 c write (2,*) "dimen",dimen," dimen3",dimen3
203 c write (2,*) "End MD setup"
205 c write (iout,*) "My chunk of ginv_block"
206 c call MATOUT2(my_ng_count,dimen3,maxres2,maxers2,ginv_block)
207 else if (iorder.eq.6) then
208 call int_from_cart1(.false.)
209 else if (iorder.eq.7) then
211 else if (iorder.eq.8) then
213 else if (iorder.eq.9) then
214 cmd call fricmat_mult(z,d_a_tmp)
215 else if (iorder.eq.10) then
216 cmd call setup_fricmat
219 write (*,*) 'Processor',fg_rank,' CG group',kolor,
220 & ' absolute rank',myrank,' leves ERGASTULUM.'
221 write(*,*)'Processor',fg_rank,' wait times for respective orders',
222 & (' order[',i,']',time_order(i),i=0,10)
226 ************************************************************************
227 subroutine func(n,x,nf,f,uiparm,urparm,ufparm)
228 implicit real*8 (a-h,o-z)
230 include 'COMMON.DERIV'
231 include 'COMMON.IOUNITS'
234 double precision energia(0:n_ene)
236 double precision ufparm
242 c write (iout,'(10f8.3)') (rad2deg*x(i),i=1,n)
246 cd print *,'func',nf,nfl,icg
247 call var_to_geom(n,x)
250 cd write (iout,*) 'ETOTAL called from FUNC'
251 call etotal(energia(0))
255 c write (iout,'(10f8.3)') (rad2deg*x(i),i=1,n)
256 c write (iout,*) 'f=',etot
261 ************************************************************************
262 subroutine func_restr(n,x,nf,f,uiparm,urparm,ufparm)
263 implicit real*8 (a-h,o-z)
265 include 'COMMON.DERIV'
266 include 'COMMON.IOUNITS'
269 double precision energia(0:n_ene)
271 double precision ufparm
277 c write (iout,'(10f8.3)') (rad2deg*x(i),i=1,n)
281 call var_to_geom_restr(n,x)
284 cd write (iout,*) 'ETOTAL called from FUNC'
285 call etotal(energia(0))
289 c write (iout,'(10f8.3)') (rad2deg*x(i),i=1,n)
290 c write (iout,*) 'f=',etot
295 c-------------------------------------------------------
296 subroutine x2xx(x,xx,n)
297 implicit real*8 (a-h,o-z)
300 include 'COMMON.CHAIN'
301 include 'COMMON.INTERACT'
302 double precision xx(maxvar),x(maxvar)
312 if (mask_phi(i).eq.1) then
320 if (mask_theta(i).eq.1) then
328 if (itype(i).ne.10) then
330 if (mask_side(i).eq.1) then
343 subroutine xx2x(x,xx)
344 implicit real*8 (a-h,o-z)
347 include 'COMMON.CHAIN'
348 include 'COMMON.INTERACT'
349 double precision xx(maxvar),x(maxvar)
359 if (mask_phi(i).eq.1) then
367 if (mask_theta(i).eq.1) then
375 if (itype(i).ne.10) then
377 if (mask_side(i).eq.1) then
388 c----------------------------------------------------------
389 subroutine minim_dc(etot,iretcode,nfun)
390 implicit real*8 (a-h,o-z)
392 parameter (liv=60,lv=(77+maxvar*(maxvar+17)/2))
396 include 'COMMON.SETUP'
397 include 'COMMON.IOUNITS'
400 include 'COMMON.MINIM'
401 include 'COMMON.CHAIN'
403 double precision minval,x(maxvar),d(maxvar),v(1:lv),xx(maxvar)
404 c common /przechowalnia/ v
406 double precision energia(0:n_ene)
407 external func_dc,grad_dc,fdum
408 logical not_done,change,reduce
409 double precision g(maxvar),f1
411 call deflt(2,iv,liv,lv,v)
412 * 12 means fresh start, dont call deflt
414 * max num of fun calls
415 if (maxfun.eq.0) maxfun=500
417 * max num of iterations
418 if (maxmin.eq.0) maxmin=1000
422 * selects output unit
423 if (print_min_ini+print_min_stat+print_min_res.gt.0) iv(21)=iout
424 * 1 means to print out result
426 * 1 means to print out summary stats
427 iv(23)=print_min_stat
428 * 1 means to print initial x and d
430 * min val for v(radfac) default is 0.1
432 * max val for v(radfac) default is 4.0
435 * check false conv if (act fnctn decrease) .lt. v(26)*(exp decrease)
436 * the sumsl default is 0.1
438 * false conv if (act fnctn decrease) .lt. v(34)
439 * the sumsl default is 100*machep
441 * absolute convergence
442 if (tolf.eq.0.0D0) tolf=1.0D-4
444 * relative convergence
445 if (rtolf.eq.0.0D0) rtolf=1.0D-4
447 * controls initial step size
449 * large vals of d correspond to small components of step
462 if (ialph(i,1).gt.0) then
470 call sumsl(k,d,x,func_dc,grad_dc,iv,liv,lv,v,idum,rdum,fdum)
480 if (ialph(i,1).gt.0) then
491 cd call func_dc(k,x,nf,f,idum,rdum,fdum)
492 cd call grad_dc(k,x,nf,g,idum,rdum,fdum)
496 cd call func_dc(k,x,nf,f1,idum,rdum,fdum)
498 cd print '(i5,2f15.5)',i,g(i),(f1-f)/1.0D-5
507 subroutine func_dc(n,x,nf,f,uiparm,urparm,ufparm)
508 implicit real*8 (a-h,o-z)
513 include 'COMMON.SETUP'
514 include 'COMMON.DERIV'
515 include 'COMMON.IOUNITS'
517 include 'COMMON.CHAIN'
519 double precision energia(0:n_ene)
520 double precision ufparm
537 if (ialph(i,1).gt.0) then
547 call etotal(energia(0))
550 cd print *,'func_dc ',nf,nfl,f
555 subroutine grad_dc(n,x,nf,g,uiparm,urparm,ufparm)
556 implicit real*8 (a-h,o-z)
561 include 'COMMON.SETUP'
562 include 'COMMON.CHAIN'
563 include 'COMMON.DERIV'
565 include 'COMMON.INTERACT'
566 include 'COMMON.FFIELD'
568 include 'COMMON.IOUNITS'
571 double precision urparm(1)
572 dimension x(maxvar),g(maxvar)
578 cd print *,'grad_dc ',nf,nfl,nf-nfl+1,icg
579 if (nf-nfl+1) 20,30,40
580 20 call func_dc(n,x,nf,f,uiparm,urparm,ufparm)
594 if (ialph(i,1).gt.0) then
604 C Evaluate the derivatives of virtual bond lengths and SC vectors in variables.
616 if (ialph(i,1).gt.0) then