ogeval Subroutine

private subroutine ogeval(me, Valvar, Valcon, Varder, Dercon)

COMPUTES SCALED CONTRAINTS+MERIT AND DERIVATIVES FROM SCALED VARIABLES

2008/01/16 | J. SCHOENMAEKERS | NEW

Type Bound

optgra

Arguments

Type IntentOptional Attributes Name
class(optgra), intent(inout) :: me
real(kind=wp), intent(in) :: Valvar(me%Numvar)

SCALED VARIABLES: valvar(var) = x_var / varsca(var)

real(kind=wp), intent(out) :: Valcon(me%Numcon+1)

SCALED CONTRAINTS+MERIT AND DERIVATIVES:

 CONSTRAINTS VALUE (1:NUMCON)
 MERIT       VALUE (NUMCON+1)
  • UNSCALED INSIDE CALVAL/CALDER (physical units)
  • SCALED ON EXIT FROM ogeval: valcon(con) = raw_value / consca(con)
integer(kind=ip), intent(in) :: Varder

DERIVATIVES COMPUTATION MODE

  • 0: VALUES ONLY (DERCON UNCHANGED)
  • 1: USER DEFINED (CALL CALDER) *-1: USER DEFINED, SPECIAL MODE
  • 2: NUMERIC WITH DOUBLE DIFFERENCING
  • 3: NUMERIC WITH SINGLE DIFFERENCING
real(kind=wp), intent(out) :: Dercon(me%Numcon+1,me%Numvar)

derivatives of scaled valcon w.r.t. scaled valvar: dercon(con,var) = d(valcon(con))/d(valvar(var))

  • set by calder if dervar=1 or -1
  • set by finite diff. if dervar=2 or 3
  • not modified if dervar=0

Calls

proc~~ogeval~~CallsGraph proc~ogeval optgra%ogeval proc~ogpwri optgra%ogpwri proc~ogeval->proc~ogpwri proc~ogwrit optgra%ogwrit proc~ogeval->proc~ogwrit proc~ogpwri_end optgra%ogpwri_end proc~ogpwri->proc~ogpwri_end proc~ogpwri_start optgra%ogpwri_start proc~ogpwri->proc~ogpwri_start

Called by

proc~~ogeval~~CalledByGraph proc~ogeval optgra%ogeval proc~ogcorr optgra%ogcorr proc~ogcorr->proc~ogeval proc~ogexec optgra%ogexec proc~ogexec->proc~ogeval proc~ogexec->proc~ogcorr proc~ogopti optgra%ogopti proc~ogexec->proc~ogopti proc~ogopti->proc~ogeval

Source Code

   subroutine ogeval(me,Valvar,Valcon,Varder,Dercon)

      !! COMPUTES SCALED CONTRAINTS+MERIT AND DERIVATIVES
      !! FROM     SCALED VARIABLES
      !!
      !! 2008/01/16 | J. SCHOENMAEKERS | NEW

      class(optgra),intent(inout) :: me
      real(wp),intent(in) :: Valvar(me%Numvar) !! SCALED VARIABLES:
                                               !! `valvar(var) = x_var / varsca(var)`
      real(wp),intent(out) :: Valcon(me%Numcon+1) !! SCALED CONTRAINTS+MERIT AND DERIVATIVES:
                                                  !!```
                                                  !! CONSTRAINTS VALUE (1:NUMCON)
                                                  !! MERIT       VALUE (NUMCON+1)
                                                  !!```
                                                  !!
                                                  !!  * UNSCALED INSIDE CALVAL/CALDER (physical units)
                                                  !!  * SCALED ON EXIT FROM [[ogeval]]: `valcon(con) = raw_value / consca(con)`
      integer(ip),intent(in) :: Varder !! DERIVATIVES COMPUTATION MODE
                                       !!
                                       !!  * 0: VALUES ONLY (DERCON UNCHANGED)
                                       !!  * 1: USER DEFINED      (CALL CALDER)
                                       !!  *-1: USER DEFINED, SPECIAL MODE
                                       !!  * 2: NUMERIC WITH DOUBLE DIFFERENCING
                                       !!  * 3: NUMERIC WITH SINGLE DIFFERENCING
      real(wp),intent(out) :: Dercon(me%Numcon+1,me%Numvar)
                                       !! derivatives of scaled valcon w.r.t. scaled valvar:
                                       !! `dercon(con,var) = d(valcon(con))/d(valvar(var))`
                                       !!
                                       !!  * set by calder if dervar=1 or -1
                                       !!  * set by finite diff. if dervar=2 or 3
                                       !!  * not modified if dervar=0

      integer(ip) :: var , con , cod , len , ind , numvio
      real(wp) :: val , sca , fac , per , sav , der , err , conerr , convio
      character(len=3) :: typ
      character(len=3) :: sta
      character(len=name_len) :: nam
      character(len=str_len) :: str

      !real(wp) :: ggg(4,4) , bbb(4) , vvv(4)  ! JW : not sure what this was for
      real(wp) :: objval
      real(wp) , dimension(:) , allocatable :: varvec
      real(wp) , dimension(:) , allocatable :: convec

      allocate (varvec(me%Numvar))
      allocate (convec(me%Numcon+1))

      ! ======================================================================
      ! GENERAL
      ! ----------------------------------------------------------------------
      write (str,'()')
      call me%ogwrit(3,str)
      select case (Varder)
       case ( 0 );    write (str,'("COMPUTE RESULTS")')
       case ( 1, -1); write (str,'("COMPUTE RESULTS",'//'   " AND DERIVATIVES USER DEFINED")')
       case ( 2 );    write (str,'("COMPUTE RESULTS",'//'   " AND DERIVATIVES BY DOUBLE DIFFERENCING")')
       case ( 3 );    write (str,'("COMPUTE RESULTS",'//'   " AND DERIVATIVES BY SINGLE DIFFERENCING")')
      end select
      call me%ogwrit(3,str)

      ! ======================================================================
      ! WRITE VARIABLES
      ! ----------------------------------------------------------------------
      write (str,'()')
      call me%ogwrit(3,str)
      write (str,'("VARIABLES NOT SCALED:")')
      call me%ogwrit(3,str)
      write (str,'()')
      call me%ogwrit(3,str)
      do var = 1 , me%Numvar
         val = Valvar(var)
         sca = me%Varsca(var)
         cod = me%Vartyp(var)
         select case (cod)
          case ( 0 ); typ = "FRE"
          case ( 1 ); typ = "PAR"
         end select
         nam = me%Varstr(var)
         len = me%Varlen(var)
         write (str,'("VAR/VAL/SCA/TYP/NAM=",'//'  I5,D14.6,D9.1,1X,A3,1X,A)') &
            var , val*sca , sca , typ , nam(1:len)
         call me%ogwrit(3,str)
      end do

      ! ======================================================================
      ! DE-SCALE VARIABLES
      ! ----------------------------------------------------------------------
      do var = 1 , me%Numvar
         varvec(var) = Valvar(var)*me%Varsca(var)
      end do

      ! ======================================================================
      ! GET RESULTS
      ! GET DERIVATIVES IF USER DEFINED
      ! ----------------------------------------------------------------------
      select case (Varder)
       case ( 0 );     call me%calval(varvec,Valcon,0)
       case ( 1, -1 ); call me%calval(varvec,Valcon,1)
       case ( 2 );     call me%calval(varvec,Valcon,1)
       case ( 3 );     call me%calval(varvec,Valcon,1)
      end select

! ======================================================================
!    IF ( 1==2 ) THEN
!       ggg(1,1) = +1D+01
!       ggg(2,1) = +1D+00
!       ggg(3,1) = +2D+00
!       ggg(4,1) = +3D+00
!       ggg(1,2) = ggg(2,1)
!       ggg(2,2) = +1D+01
!       ggg(3,2) = +4D+00
!       ggg(4,2) = +5D+00
!       ggg(1,3) = ggg(3,1)
!       ggg(2,3) = ggg(3,2)
!       ggg(3,3) = +1D+01
!       ggg(4,3) = +6D+00
!       ggg(1,4) = ggg(4,1)
!       ggg(2,4) = ggg(4,2)
!       ggg(3,4) = ggg(4,3)
!       ggg(4,4) = +1D+01
!       ! ----------------------------------------------------------------------
!       bbb(1) = +1D+01
!       bbb(2) = +1D+01
!       bbb(3) = +1D+01
!       bbb(4) = +1D+01
!       ! ----------------------------------------------------------------------
!       CALL mul2m(ggg,4,1,1,4,Valvar,4,1,1,4,vvv,4,1,1,1)
!       CALL mul2m(Valvar,1,1,1,1,vvv,4,1,1,4,Valcon,1,1,1,1)
!       CALL mul2m(bbb,1,1,1,1,Valvar,4,1,1,4,vvv,1,1,1,1)
!       CALL mulvs(Valcon,0.5_wp,Valcon,1)
!       CALL sum2v(Valcon,vvv,Valcon,1)
!       CALL mulvs(Valcon,-1.0_wp,Valcon,1)
!       WRITE (str,*) "VALCON=" , (Valcon(ind),ind=1,1)
!       CALL me%ogwrit(3,str)
!    end if

      ! ----------------------------------------------------------------------
      ! SCALE RESULTS
      ! ----------------------------------------------------------------------
      do con = 1 , me%Numcon + 1
         convec(con) = Valcon(con)
         sca = me%Consca(con)
         cod = me%Contyp(con)
         if ( cod==-1 ) sca = -sca
         Valcon(con) = Valcon(con)/sca
      end do

      ! ======================================================================
      ! WRITE RESULTS
      ! ----------------------------------------------------------------------
      write (str,'()')
      call me%ogwrit(3,str)
      write (str,'("RESULTS NOT SCALED:")')
      call me%ogwrit(3,str)
      write (str,'()')
      call me%ogwrit(3,str)
      conerr = 0.0_wp  ! total constraint error (scaled to constr. threshod)
      convio = 0.0_wp  ! total constaint error norm (unscaled)
      ind = 0          ! index of largest constraint violation
      fac = 0.0_wp     ! value of largest constraint violation
      numvio = 0       ! number of violated constraints
      do con = 1 , me%Numcon + 1
         val = Valcon(con)
         sca = me%Consca(con)
         cod = me%Contyp(con)
         sta = "   "
         err = 0.0_wp
         if ( cod==-1 ) sca = -sca
         if ( cod==-2 ) typ = "DER"
         if ( cod==-1 .and. con<=me%Numcon ) typ = "LTE"
         if ( cod==0 .and. con<=me%Numcon ) typ = "EQU"
         if ( cod==1 .and. con<=me%Numcon ) typ = "GTE"
         if ( cod==-1 .and. con>me%Numcon ) typ = "MIN"
         if ( cod==1 .and. con>me%Numcon ) typ = "MAX"
         if ( cod==0 .and. con<=me%Numcon .and. abs(val)>1.0_wp ) then
            sta = "VIO"
            err = abs(val)
            numvio = numvio + 1
         end if
         if ( cod/=0 .and. con<=me%Numcon .and. cod/=-2 .and. -val>1.0_wp ) then
            sta = "VIO"
            err = abs(val)
            numvio = numvio + 1
         end if
         conerr = conerr + err
         convio = convio + (err*sca)**2
         if ( err>fac ) ind = con
         if ( err>fac ) fac = err
         nam = me%Constr(con)
         len = me%Conlen(con)
         write (str,'("CON/VAL/SCA/TYP/STA/NAM=",'//'  I5,D14.6,D9.1,1X,A3,1X,A3,1X,A)') &
            con , val*sca , sca , typ , sta , nam(1:len)
         call me%ogwrit(3,str)
      end do
      write (str,'()')
      call me%ogwrit(3,str)
      write (str,'("CONSTRAINT ERROR.:",2(1X,D13.6),I6)') conerr , fac , ind
      call me%ogwrit(3,str)
      write (str,'()')
      call me%ogwrit(3,str)
      ! write pygmo-style log output
      objval = -Valcon(me%Numcon+1)
      convio = sqrt(convio)
      call me%ogpwri(objval,numvio,convio)

      ! ======================================================================
      ! NO DERIVATIVES
      ! ----------------------------------------------------------------------
      if ( Varder==0 ) then
         return
      elseif ( Varder==1 .or. Varder==-1 ) then
         call me%calder(varvec,convec,Dercon)
      end if

! ----------------------------------------------------------------------
!    IF ( 1==2 ) THEN
!       CALL mul2m(Valvar,1,1,1,1,ggg,4,1,1,4,Dercon,1,1,1,4)
!       CALL sum2v(Dercon,bbb,Dercon,4)
!       CALL mulvs(Dercon,-1.0_wp,Dercon,4)
!       WRITE (str,*) "DERCON=" , (Dercon(1,ind),ind=1,4)
!       CALL me%ogwrit(3,str)
!    end if

      ! ======================================================================
      ! WRITE DERIVATIVES
      ! ----------------------------------------------------------------------
      write (str,'()')
      call me%ogwrit(3,str)
      write (str,'("DERIVATIVES SCALED:")')
      call me%ogwrit(3,str)
      write (str,'()')
      call me%ogwrit(3,str)

      do var = 1 , me%Numvar

         ! ----------------------------------------------------------------------
         ! WRITE VARIABLE
         ! ----------------------------------------------------------------------
         val = Valvar(var)
         sca = me%Varsca(var)
         cod = me%Vartyp(var)
         if ( cod==0 ) typ = "FRE"
         if ( cod==1 ) typ = "PAR"
         nam = me%Varstr(var)
         len = me%Varlen(var)
         write (str,'("VAR/VAL/SCA/TYP/NAM=",'//'  I5,D14.6,D9.1,1X,A3,1X,A)') &
                  var , val*sca , sca , typ , nam(1:len)
         call me%ogwrit(4,str)
         write (str,'()')
         call me%ogwrit(4,str)

         select case (Varder)
          case ( 2 )  ! DERIVATIVES BY DOUBLE DIFFERENCING
            per = me%Varper(var)
            sav = varvec(var)
            varvec(var) = sav + per
            call me%calval(varvec,Dercon(1:,var),0)  ! JW added : here
            varvec(var) = sav - per
            call me%calval(varvec,convec,0)
            fac = 0.5_wp/per
            do con = 1 , me%Numcon + 1
               Dercon(con,var) = (Dercon(con,var)-convec(con))*fac
            end do
            varvec(var) = sav
          case ( 3 )  ! DERIVATIVES BY SINGLE DIFFERENCING
            per = me%Varper(var)
            sav = varvec(var)
            varvec(var) = sav + per
            call me%calval(varvec,Dercon(1:,var),0)  ! JW added : here
            fac = 1.0_wp/per
            do con = 1 , me%Numcon + 1
               Dercon(con,var) = (Dercon(con,var)-convec(con))*fac
            end do
            varvec(var) = sav
         end select

         ! ----------------------------------------------------------------------
         ! SCALE DERIVATIVES
         ! ----------------------------------------------------------------------
         do con = 1 , me%Numcon + 1
            fac = me%Varsca(var)/me%Consca(con)
            if ( me%Contyp(con)==-1 ) fac = -fac
            Dercon(con,var) = Dercon(con,var)*fac
         end do

         ! ----------------------------------------------------------------------
         ! WRITE DERIVATIVES
         ! ======================================================================
         do con = 1 , me%Numcon + 1
            der = Dercon(con,var)
            if ( der/=0.0_wp ) then
               sca = me%Consca(con)
               cod = me%Contyp(con)
               if ( cod==-1 ) sca = -sca
               if ( cod==-2 ) typ = "DER"
               if ( cod==-1 .and. con<=me%Numcon ) typ = "LTE"
               if ( cod==0 .and. con<=me%Numcon ) typ = "EQU"
               if ( cod==1 .and. con<=me%Numcon ) typ = "GTE"
               if ( cod==-1 .and. con>me%Numcon ) typ = "MIN"
               if ( cod==1 .and. con>me%Numcon ) typ = "MAX"
               nam = me%Constr(con)
               len = me%Conlen(con)
               write (str,'("CON/DER/SCA/TYP/NAM=",'//'  I5,D14.6,D9.1,1X,A3,1X,A)') &
                  con , der*sca/me%Varsca(var) , sca , typ , nam(1:len)
               call me%ogwrit(4,str)
            end if
         end do
         write (str,'()')
         call me%ogwrit(4,str)

      end do

      deallocate (varvec)
      deallocate (convec)

   end subroutine ogeval