subroutine maxlik_init(nvarr,xrange,comm) c Optimize the UNRES energy function by minimization of a quartic target c function or by the VMC method. implicit none #ifndef ISNAN external proc_proc #endif #ifdef WINPGI cMS$ATTRIBUTES C :: proc_proc #endif include "DIMENSIONS" include "DIMENSIONS.ZSCOPT" #ifdef MPI include "mpif.h" integer IERROR,ERRCODE,kolor,key,comm include "COMMON.MPI" #endif include "COMMON.IOUNITS" include "COMMON.OPTIM" include "COMMON.XBOUND" integer nvarr,iparm double precision rr,x(max_paropt) integer idumm integer i #ifdef PYTHON double precision xrange(maxvar,2) integer number_of_variables common /patch/ number_of_variables Cf2py intent(in) comm Cf2py intent(out) nvarr Cf2Py intent(out) xrange #endif c print *,"Starting..." #ifdef MPI #ifndef PYTHON c print *,"Initializing MPI..." call MPI_Init( IERROR ) ALL_COMM = MPI_COMM_WORLD #else ALL_COMM = comm #endif call MPI_Comm_rank( ALL_COMM, me, IERROR ) call MPI_Comm_size( ALL_COMM, nprocs, IERROR ) c print *,"Finished initializing MPI..." Master = 0 Master1 = 1 c print *,"Me",me," Master",master," Ierror",ierror #ifndef PYTHON if (ierror.gt.0) then write(iout,*) "SEVERE ERROR - Can't initialize MPI." call mpi_finalize(ierror) stop endif #endif #endif #ifndef ISNAN c NaNQ initialization i=-1 rr=dacos(100.0d0) #ifdef WINPGI idumm=proc_proc(rr,i) #else call proc_proc(rr,i) #endif #endif call initialize c print *,"calling openunits" call openunits c print *,"openunits called" call read_general_data(*10) write (iout,'(80(1h-)/10x, & "Maximum likelihood optimization of UNRES energy function", & " v. 05/10/16"/80(1h-))') call flush(iout) call cinfo c call promienie(*10) write (iout,*) "Finished READ_GENERAL_DATA" call flush(iout) do iparm=1,nparmset call parmread(iparm,*10) enddo write (iout,*) "Finished parmread" call flush(iout) call read_optim_parm(*10) call print_general_data(*10) call read_protein_data(*10) write (iout,*) "Finished READ_PROTEIN_DATA" call flush(iout) call read_database(*10) write (iout,*) "Finished READ_DATABASE" call flush(iout) #ifdef MPI c write (iout,*) Me,' calling PROC_GROUPS' call proc_groups c write (iout,*) Me,' calling WORK_PARTITION_MAP' c call work_partition_map(nvarr) #endif call proc_data(nvarr,x,*10) #ifdef PYTHON number_of_variables=nvarr do i=1,nvarr xrange(i,1)=x_low(i) xrange(i,2)=x_up(i) enddo write (iout,*) "xrange from MAXLIK_INIT" do i=1,nvarr write (iout,*) i,xrange(i,1),xrange(i,2) enddo #endif call read_thermal return 10 write (iout,*) "Error termination of the program" call MPI_Finalize( IERROR ) return end c------------------------------------------------------------------------ subroutine maxlik_optim(x,fmin) c Optimize the UNRES energy function by minimization of a quartic target c function or by the VMC method. implicit none #ifndef ISNAN external proc_proc #endif #ifdef WINPGI cMS$ATTRIBUTES C :: proc_proc #endif include "DIMENSIONS" include "DIMENSIONS.ZSCOPT" #ifdef MPI include "mpif.h" integer IERROR,ERRCODE,kolor,key include "COMMON.MPI" #endif include "COMMON.IOUNITS" include "COMMON.OPTIM" integer nvarr,iparm,i double precision rr,x(max_paropt),fmin #ifdef PYTHON integer number_of_variables common /patch/ number_of_variables Cf2py intent(inout) x Cf2py intent(out) fmin nvarr=number_of_variables write (*,*) "MAXLIK_OPTIM: Variables from PYTHON:",nvarr write (iout,*) "MAXLIK_OPTIM: Variables from PYTHON:",nvarr do i=1,nvarr write (iout,*) i,x(i) enddo #endif #ifdef MPI if (me.eq.Master) then #endif call maxlikopt(nvarr,x,fmin) write (iout,*) "fmin from MAXLIK_OPTIM:",fmin #ifdef MPI call jebadelko(nvarr) else call jebadelko(nvarr) endif call bilans call MPI_Finalize( IERROR ) #else call bilans #endif return end