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
35 call deflt(2,iv,liv,lv,v)
36 * 12 means fresh start, dont call deflt
38 * max num of fun calls
39 if (maxfun.eq.0) maxfun=500
41 * max num of iterations
42 if (maxmin.eq.0) maxmin=1000
49 * 1 means to print out result
51 * 1 means to print out summary stats
53 * 1 means to print initial x and d
55 * min val for v(radfac) default is 0.1
57 * max val for v(radfac) default is 4.0
60 * check false conv if (act fnctn decrease) .lt. v(26)*(exp decrease)
61 * the sumsl default is 0.1
63 * false conv if (act fnctn decrease) .lt. v(34)
64 * the sumsl default is 100*machep
66 * absolute convergence
67 if (tolf.eq.0.0D0) tolf=1.0D-4
69 * relative convergence
70 if (rtolf.eq.0.0D0) rtolf=1.0D-4
72 * controls initial step size
74 * large vals of d correspond to small components of step
81 cd print *,'Calling SUMSL'
82 c call var_to_geom(nvar,x)
84 c call etotal(energia(0))
87 call x2xx(x,xx,nvar_restr)
88 call sumsl(nvar_restr,d,xx,func_restr,grad_restr,
89 & iv,liv,lv,v,idum,rdum,fdum)
92 call sumsl(nvar,d,x,func,gradient,iv,liv,lv,v,idum,rdum,fdum)
96 cd print *,'Exit SUMSL; return code:',iretcode,' energy:',etot
97 cd write (iout,'(/a,i4/)') 'SUMSL return code:',iv(1)
100 call var_to_geom(nvar,x)
102 c write (iout,'(a)') 'Reduction worked, minimizing again...'
107 c call etotal(energia(0))
109 c call enerprint(energia(0))
112 c write (*,*) 'Processor',MyID,' leaves MINIMIZE.'
119 c----------------------------------------------------------------------------
120 subroutine ergastulum
121 implicit real*8 (a-h,o-z)
126 include 'COMMON.SETUP'
127 include 'COMMON.DERIV'
129 include 'COMMON.IOUNITS'
130 include 'COMMON.FFIELD'
131 include 'COMMON.INTERACT'
133 include 'COMMON.TIME1'
134 double precision z(maxres6),d_a_tmp(maxres6)
135 double precision edum(0:n_ene),time_order(0:10)
136 double precision Gcopy(maxres2,maxres2)
137 common /przechowalnia/ Gcopy
139 C Workers wait for variables and NF, and NFL from the boss
141 do while (iorder.ge.0)
142 c write (*,*) 'Processor',fg_rank,' CG group',kolor,
143 c & ' receives order from Master'
145 call MPI_Bcast(iorder,1,MPI_INTEGER,king,FG_COMM,IERR)
146 time_Bcast=time_Bcast+MPI_Wtime()-time00
147 if (icall.gt.4 .and. iorder.ge.0)
148 & time_order(iorder)=time_order(iorder)+MPI_Wtime()-time00
151 c & 'Processor',fg_rank,' completed receive MPI_BCAST order',iorder
152 if (iorder.eq.0) then
155 c write (2,*) "After etotal"
156 c write (2,*) "dimen",dimen," dimen3",dimen3
158 else if (iorder.eq.2) then
160 cmd call etotal_short(edum)
161 c write (2,*) "After etotal_short"
162 c write (2,*) "dimen",dimen," dimen3",dimen3
164 else if (iorder.eq.3) then
166 cmd call etotal_long(edum)
167 c write (2,*) "After etotal_long"
168 c write (2,*) "dimen",dimen," dimen3",dimen3
170 else if (iorder.eq.1) then
172 c write (2,*) "After sum_gradient"
173 c write (2,*) "dimen",dimen," dimen3",dimen3
175 else if (iorder.eq.4) then
176 cmd call ginv_mult(z,d_a_tmp)
177 else if (iorder.eq.5) then
178 c Setup MD things for a slave
179 dimen=(nct-nnt+1)+nside
180 dimen1=(nct-nnt)+(nct-nnt+1)
182 c write (2,*) "dimen",dimen," dimen3",dimen3
184 call int_bounds(dimen,igmult_start,igmult_end)
185 igmult_start=igmult_start-1
186 call MPI_Allgather(3*igmult_start,1,MPI_INTEGER,
187 & ng_start(0),1,MPI_INTEGER,FG_COMM,IERROR)
188 my_ng_count=igmult_end-igmult_start
189 call MPI_Allgather(3*my_ng_count,1,MPI_INTEGER,ng_counts(0),1,
190 & MPI_INTEGER,FG_COMM,IERROR)
191 c write (2,*) "ng_start",(ng_start(i),i=0,nfgtasks-1)
192 c write (2,*) "ng_counts",(ng_counts(i),i=0,nfgtasks-1)
193 myginv_ng_count=maxres2*my_ng_count
194 c write (2,*) "igmult_start",igmult_start," igmult_end",
195 c & igmult_end," my_ng_count",my_ng_count
197 call MPI_Allgather(maxres2*igmult_start,1,MPI_INTEGER,
198 & nginv_start(0),1,MPI_INTEGER,FG_COMM,IERROR)
199 call MPI_Allgather(myginv_ng_count,1,MPI_INTEGER,
200 & nginv_counts(0),1,MPI_INTEGER,FG_COMM,IERROR)
201 c write (2,*) "nginv_start",(nginv_start(i),i=0,nfgtasks-1)
202 c write (2,*) "nginv_counts",(nginv_counts(i),i=0,nfgtasks-1)
204 c call MPI_Barrier(FG_COMM,IERROR)
206 call MPI_Scatterv(ginv(1,1),nginv_counts(0),
207 & nginv_start(0),MPI_DOUBLE_PRECISION,gcopy(1,1),
208 & myginv_ng_count,MPI_DOUBLE_PRECISION,king,FG_COMM,IERR)
210 time_scatter_ginv=time_scatter_ginv+MPI_Wtime()-time00
217 c write (2,*) "dimen",dimen," dimen3",dimen3
218 c write (2,*) "End MD setup"
220 c write (iout,*) "My chunk of ginv_block"
221 c call MATOUT2(my_ng_count,dimen3,maxres2,maxers2,ginv_block)
222 else if (iorder.eq.6) then
223 call int_from_cart1(.false.)
224 else if (iorder.eq.7) then
226 else if (iorder.eq.8) then
228 else if (iorder.eq.9) then
229 cmd call fricmat_mult(z,d_a_tmp)
230 else if (iorder.eq.10) then
231 cmd call setup_fricmat
234 write (*,*) 'Processor',fg_rank,' CG group',kolor,
235 & ' absolute rank',myrank,' leves ERGASTULUM.'
236 write(*,*)'Processor',fg_rank,' wait times for respective orders',
237 & (' order[',i,']',time_order(i),i=0,10)
241 ************************************************************************
242 subroutine func(n,x,nf,f,uiparm,urparm,ufparm)
243 implicit real*8 (a-h,o-z)
245 include 'COMMON.DERIV'
246 include 'COMMON.IOUNITS'
249 double precision energia(0:n_ene)
251 double precision ufparm
257 c write (iout,'(10f8.3)') (rad2deg*x(i),i=1,n)
261 cd print *,'func',nf,nfl,icg
262 call var_to_geom(n,x)
265 cd write (iout,*) 'ETOTAL called from FUNC'
266 call etotal(energia(0))
270 c write (iout,'(10f8.3)') (rad2deg*x(i),i=1,n)
271 c write (iout,*) 'f=',etot
276 ************************************************************************
277 subroutine func_restr(n,x,nf,f,uiparm,urparm,ufparm)
278 implicit real*8 (a-h,o-z)
280 include 'COMMON.DERIV'
281 include 'COMMON.IOUNITS'
284 double precision energia(0:n_ene)
286 double precision ufparm
292 c write (iout,'(10f8.3)') (rad2deg*x(i),i=1,n)
296 call var_to_geom_restr(n,x)
299 cd write (iout,*) 'ETOTAL called from FUNC'
300 call etotal(energia(0))
304 c write (iout,'(10f8.3)') (rad2deg*x(i),i=1,n)
305 c write (iout,*) 'f=',etot
310 c-------------------------------------------------------
311 subroutine x2xx(x,xx,n)
312 implicit real*8 (a-h,o-z)
315 include 'COMMON.CHAIN'
316 include 'COMMON.INTERACT'
317 double precision xx(maxvar),x(maxvar)
327 if (mask_phi(i).eq.1) then
335 if (mask_theta(i).eq.1) then
343 if (itype(i).ne.10) then
345 if (mask_side(i).eq.1) then
358 subroutine xx2x(x,xx)
359 implicit real*8 (a-h,o-z)
362 include 'COMMON.CHAIN'
363 include 'COMMON.INTERACT'
364 double precision xx(maxvar),x(maxvar)
374 if (mask_phi(i).eq.1) then
382 if (mask_theta(i).eq.1) then
390 if (itype(i).ne.10) then
392 if (mask_side(i).eq.1) then
403 c----------------------------------------------------------
404 subroutine minim_dc(etot,iretcode,nfun)
405 implicit real*8 (a-h,o-z)
407 parameter (liv=60,lv=(77+maxvar*(maxvar+17)/2))
411 include 'COMMON.SETUP'
412 include 'COMMON.IOUNITS'
415 include 'COMMON.MINIM'
416 include 'COMMON.CHAIN'
418 double precision minval,x(maxvar),d(maxvar),v(1:lv),xx(maxvar)
419 c common /przechowalnia/ v
421 double precision energia(0:n_ene)
422 external func_dc,grad_dc,fdum
423 logical not_done,change,reduce
424 double precision g(maxvar),f1
426 call deflt(2,iv,liv,lv,v)
427 * 12 means fresh start, dont call deflt
429 * max num of fun calls
430 if (maxfun.eq.0) maxfun=500
432 * max num of iterations
433 if (maxmin.eq.0) maxmin=1000
437 * selects output unit
440 * 1 means to print out result
442 * 1 means to print out summary stats
444 * 1 means to print initial x and d
446 * min val for v(radfac) default is 0.1
448 * max val for v(radfac) default is 4.0
451 * check false conv if (act fnctn decrease) .lt. v(26)*(exp decrease)
452 * the sumsl default is 0.1
454 * false conv if (act fnctn decrease) .lt. v(34)
455 * the sumsl default is 100*machep
457 * absolute convergence
458 if (tolf.eq.0.0D0) tolf=1.0D-4
460 * relative convergence
461 if (rtolf.eq.0.0D0) rtolf=1.0D-4
463 * controls initial step size
465 * large vals of d correspond to small components of step
478 if (ialph(i,1).gt.0) then
486 call sumsl(k,d,x,func_dc,grad_dc,iv,liv,lv,v,idum,rdum,fdum)
496 if (ialph(i,1).gt.0) then
507 cd call func_dc(k,x,nf,f,idum,rdum,fdum)
508 cd call grad_dc(k,x,nf,g,idum,rdum,fdum)
512 cd call func_dc(k,x,nf,f1,idum,rdum,fdum)
514 cd print '(i5,2f15.5)',i,g(i),(f1-f)/1.0D-5
523 subroutine func_dc(n,x,nf,f,uiparm,urparm,ufparm)
524 implicit real*8 (a-h,o-z)
529 include 'COMMON.SETUP'
530 include 'COMMON.DERIV'
531 include 'COMMON.IOUNITS'
533 include 'COMMON.CHAIN'
535 double precision energia(0:n_ene)
536 double precision ufparm
553 if (ialph(i,1).gt.0) then
563 call etotal(energia(0))
566 cd print *,'func_dc ',nf,nfl,f
571 subroutine grad_dc(n,x,nf,g,uiparm,urparm,ufparm)
572 implicit real*8 (a-h,o-z)
577 include 'COMMON.SETUP'
578 include 'COMMON.CHAIN'
579 include 'COMMON.DERIV'
581 include 'COMMON.INTERACT'
582 include 'COMMON.FFIELD'
584 include 'COMMON.IOUNITS'
587 double precision urparm(1)
588 dimension x(maxvar),g(maxvar)
594 cd print *,'grad_dc ',nf,nfl,nf-nfl+1,icg
595 if (nf-nfl+1) 20,30,40
596 20 call func_dc(n,x,nf,f,uiparm,urparm,ufparm)
610 if (ialph(i,1).gt.0) then
620 C Evaluate the derivatives of virtual bond lengths and SC vectors in variables.
632 if (ialph(i,1).gt.0) then