51 INTEGER,
Dimension(:),
Pointer :: cntvect
93 real*8 time0, time1, time2
99 coi_error = coi_createfort( cntvect )
106 coi_error = max( coi_error,
coidef_numnz( cntvect, 6*ne ) )
111 coi_error = max( coi_error,
coidef_optfile( cntvect,
'mp_elec2.opt' ) )
125#if defined(LICENSE_INT_1) && defined(LICENSE_INT_2) && defined(LICENSE_INT_3) && defined(LICENSE_TEXT)
126 coi_error = max( coi_error,
coidef_license( cntvect, license_int_1, license_int_2, license_int_3, license_text) )
129 If ( coi_error .ne. 0 )
THEN
131 write(*,*)
'**** Fatal Error while loading CONOPT Callback routines.'
133 call flog(
"Skipping Solve due to setup errors", 1 )
138 coi_error = max( coi_error,
coidef_hessfac( cntvect, 3.01d0*(ne+1)/4.d0 ) )
141 time0 = omp_get_wtime()
143 time1 = omp_get_wtime() - time0
146 write(*,*)
'End of Electron example one thread. Return code=',coi_error
148 If ( coi_error /= 0 )
then
149 call flog(
"One Thread: Errors encountered during solution", 1 )
150 elseif ( stacalls == 0 .or. solcalls == 0 )
then
151 call flog(
"One Thread: Status or Solution routine was not called", 1 )
152 elseif ( sstat /= 1 .or. mstat /= 2 )
then
153 call flog(
"One Thread: Solver and Model Status was not as expected (1,2)", 1 )
154 elseif (
hcalls == 0 )
then
155 call flog(
"One Thread: Solve did not use 2nd derivatives but should have enough memory", 1 )
162 time0 = omp_get_wtime()
164 time2 = omp_get_wtime() - time0
167 write(*,*)
'End of Electron example multi threads. Return code=',coi_error
169 If ( coi_error /= 0 )
then
170 call flog(
"Multi thread: Errors encountered during solution", 1 )
171 elseif ( stacalls == 0 .or. solcalls == 0 )
then
172 call flog(
"Multi thread: Status or Solution routine was not called", 1 )
173 elseif ( sstat /= 1 .or. mstat /= 2 )
then
174 call flog(
"Multi thread: Solver and Model Status was not as expected (1,2)", 1 )
175 elseif (
hcalls == 0 )
then
176 call flog(
"Multi thread: Solve did not use 2nd derivatives but should have enough memory", 1 )
179 if ( coi_free( cntvect ) /= 0 )
call flog(
"Error while freeing control vector", 1 )
182 write(*,
"('Time for single thread',f10.3)") time1
183 write(*,
"('Time for multi thread',f10.3)") time2
184 write(*,
"('Speedup ',f10.3)") time1/time2
185 write(*,
"('Efficiency ',f10.3)") time1/time2/omp_get_max_threads()
187 call flog(
"Successful Solve", 0 )
198 seed = mod(seed*1027+25,1048576)
199 rndx = float(seed)/float(1048576)
211 colsta, rowno, value, nlflag, n, m, nz, &
217 integer,
intent (in) :: n
218 integer,
intent (in) :: m
219 integer,
intent (in) :: nz
220 real*8,
intent (in out),
dimension(n) :: lower
221 real*8,
intent (in out),
dimension(n) :: curr
222 real*8,
intent (in out),
dimension(n) :: upper
223 integer,
intent (in out),
dimension(n) :: vsta
225 integer,
intent (out),
dimension(m) ::
type
226 integer,
intent (in out),
dimension(m) :: esta
228 real*8,
intent (in out),
dimension(m) :: rhs
229 integer,
intent (in out),
dimension(n+1) :: colsta
230 integer,
intent (out),
dimension(nz) :: rowno
231 integer,
intent (in out),
dimension(nz) :: nlflag
232 real*8,
intent (in out),
dimension(nz) ::
value
237 real*8,
parameter :: pi = 3.141592
239 Real,
External ::
rndx
250 curr(i ) = cos(theta)*sin(phi)
251 curr(i+ ne) = sin(theta)*sin(phi)
252 curr(i+2*ne) = cos(phi)
334Integer Function elec_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
335 n, nz, thread, usrmem )
340 integer,
intent (in) :: n
341 integer,
intent (in) :: rowno
342 integer,
intent (in) :: nz
343 real*8,
intent (in),
dimension(n) :: x
344 real*8,
intent (in out) :: g
345 real*8,
intent (in out),
dimension(n) :: jac
346 integer,
intent (in),
dimension(nz) :: jcnm
348 integer,
intent (in) :: mode
350 integer,
intent (in) :: ignerr
352 integer,
intent (in out) :: errcnt
354 integer,
intent (in) :: thread
359 real*8 :: dist, dist32, tx, ty, tz
365 if ( rowno == ne+1 )
then
369 if ( mode .eq. 1 .or. mode .eq. 3 )
then
373 dist = (x(i)-x(j))**2 + (x(i+ne)-x(j+ne))**2 + (x(i+2*ne)-x(j+2*ne))**2
374 g = g + 1.d0/sqrt(dist)
381 if ( mode .eq. 2 .or. mode .eq. 3 )
then
387 dist = (x(i)-x(j))**2 + (x(i+ne)-x(j+ne))**2 + (x(i+2*ne)-x(j+2*ne))**2
388 dist32 = (1.d0/sqrt(dist))**3
389 tx = -(x(i)-x(j))*dist32
390 ty = -(x(i+ne)-x(j+ne))*dist32
391 tz = -(x(i+2*ne)-x(j+2*ne))*dist32
394 jac(i+ne) = jac(i+ne) + ty
395 jac(j+ne) = jac(j+ne) - ty
396 jac(i+2*ne) = jac(i+2*ne) + tz
397 jac(j+2*ne) = jac(j+2*ne) - tz
407 if ( mode .eq. 1 .or. mode .eq. 3 )
then
408 g = x(i)**2 + x(i+ne)**2 + x(i+2*ne)**2
413 if ( mode .eq. 2 .or. mode .eq. 3 )
then
415 jac(i+ne) = 2.d0*x(i+ne)
416 jac(i+2*ne) = 2.d0*x(i+2*ne)
431 Integer,
Intent (IN) :: numvar, numcon, maxhess
432 Integer,
Intent (IN OUT) :: nodrv, numhess
433 real*8,
Intent(IN OUT) :: usrmem(*)
437 numhess = numvar * (numvar+1)/2
446Integer Function elec_2dlagrstr( HSRW, HSCL, NODRV, NumVar, NumCon, NumHess, UsrMem )
453 Integer,
Intent (IN) :: numvar, numcon, numhess
454 Integer,
Intent (IN) :: nodrv
455 Integer,
Dimension(*) :: hsrw, hscl
456 real*8,
Intent(IN OUT) :: usrmem(*)
478Integer Function elec_2dlagrval( X, U, HSRW, HSCL, HSVL, NODRV, NumVar, NumCon, NumHess, UsrMem )
485 Integer,
Intent (IN) :: numvar, numcon, numhess
486 Integer,
Intent (IN OUT) :: nodrv
487 real*8,
Dimension(NumVar),
Intent (IN) :: x
488 real*8,
Dimension(NumCon),
Intent (IN) :: u
489 Integer,
Dimension(NumHess),
Intent (IN) :: hsrw, hscl
490 real*8,
Dimension(NumHess),
Intent (Out) :: hsvl
491 real*8,
Intent(IN OUT) :: usrmem(*)
494 integer :: k, l, ix, iy, iz, jx, jy, jz
495 real*8 :: tx, ty, tz, dist, distsq, dist3, dist5
510 hsvl(k) = hsvl(k) + 2.d0*u(l)
512 hsvl(k) = hsvl(k) + 2.d0*u(l)
513 k =
pos(l+2*ne,l+2*ne)
514 hsvl(k) = hsvl(k) + 2.d0*u(l)
516 if ( u(ne+1) /= 0.0d0 )
Then
531 dist = tx**2 + ty**2 + tz**2
541 distsq = 1/sqrt(dist)
554 hsvl(k) = hsvl(k) + 3.0d0*dist5*tx**2 - dist3
556 hsvl(k) = hsvl(k) + 3.0d0*dist5*tx**2 - dist3
558 hsvl(k) = hsvl(k) - 3.0d0*dist5*tx**2 + dist3
560 hsvl(k) = hsvl(k) + 3.0d0*dist5*ty**2 - dist3
562 hsvl(k) = hsvl(k) + 3.0d0*dist5*ty**2 - dist3
564 hsvl(k) = hsvl(k) - 3.0d0*dist5*ty**2 + dist3
566 hsvl(k) = hsvl(k) + 3.0d0*dist5*tz**2 - dist3
568 hsvl(k) = hsvl(k) + 3.0d0*dist5*tz**2 - dist3
570 hsvl(k) = hsvl(k) - 3.0d0*dist5*tz**2 + dist3
575 hsvl(k) = hsvl(k) + 3.0d0*dist5*tx*ty
577 hsvl(k) = hsvl(k) - 3.0d0*dist5*tx*ty
579 hsvl(k) = hsvl(k) - 3.0d0*dist5*tx*ty
581 hsvl(k) = hsvl(k) + 3.0d0*dist5*tx*ty
583 hsvl(k) = hsvl(k) + 3.0d0*dist5*tx*tz
585 hsvl(k) = hsvl(k) - 3.0d0*dist5*tx*tz
587 hsvl(k) = hsvl(k) - 3.0d0*dist5*tx*tz
589 hsvl(k) = hsvl(k) + 3.0d0*dist5*tx*tz
591 hsvl(k) = hsvl(k) + 3.0d0*dist5*ty*tz
593 hsvl(k) = hsvl(k) - 3.0d0*dist5*ty*tz
595 hsvl(k) = hsvl(k) - 3.0d0*dist5*ty*tz
597 hsvl(k) = hsvl(k) + 3.0d0*dist5*ty*tz
607 Integer,
Intent(in) :: i, j
608 pos = 1+(i-1)*(2*numvar+2-i)/2 + (j-i)
integer function coidef_2dlagrsize(cntvect, coi_2dlagrsize)
integer function std_solution(xval, xmar, xbas, xsta, yval, ymar, ybas, ysta, n, m, usrmem)
integer function std_status(modsta, solsta, iter, objval, usrmem)
integer function std_message(smsg, dmsg, nmsg, llen, usrmem, msgv)
integer function std_errmsg(rowno, colno, posno, msglen, usrmem, msg)
integer function elec_2dlagrval(x, u, hsrw, hscl, hsvl, nodrv, numvar, numcon, numhess, usrmem)
Compute the Lagrangian of the Hessian.
integer function elec_2dlagrsize(nodrv, numvar, numcon, numhess, maxhess, usrmem)
integer function pos(i, j)
integer function elec_2dlagrstr(hsrw, hscl, nodrv, numvar, numcon, numhess, usrmem)
Specify the structure of the Lagrangian of the Hessian.
integer function elec_fdeval(x, g, jac, rowno, jcnm, mode, ignerr, errcnt, n, nz, thread, usrmem)
Compute nonlinear terms and non-constant Jacobian elements.
integer function elec_readmatrix(lower, curr, upper, vsta, type, rhs, esta, colsta, rowno, value, nlflag, n, m, nz, usrmem)
Define information about the model.
integer function coidef_fdeval(cntvect, coi_fdeval)
define callback routine for performing function and derivative evaluations.
integer function coidef_errmsg(cntvect, coi_errmsg)
define callback routine for returning error messages for row, column or Jacobian elements.
integer function coidef_message(cntvect, coi_message)
define callback routine for handling messages returned during the solution process.
integer function coidef_readmatrix(cntvect, coi_readmatrix)
define callback routine for providing the matrix data to CONOPT.
integer function coidef_status(cntvect, coi_status)
define callback routine for returning the completion status.
integer function coidef_solution(cntvect, coi_solution)
define callback routine for returning the final solution values.
integer function coidef_2dlagrstr(cntvect, coi_2dlagrstr)
define callback routine for providing the structure of the second derivatives of the Lagrangian.
integer function coidef_optfile(cntvect, optfile)
define callback routine for defining an options file.
integer function coidef_2dlagrval(cntvect, coi_2dlagrval)
define callback routine for computing the values of the second derivatives of the Lagrangian.
integer function coidef_threads(cntvect, threads)
number of threads allowed internally in CONOPT.
integer function coidef_hessfac(cntvect, hessfac)
factor for Hessian density relative to Jacobian density HessFac.
integer function coidef_debug2d(cntvect, debug2d)
turn debugging of 2nd derivatives on and off.
integer function coidef_license(cntvect, licint1, licint2, licint3, licstring)
define the License Information.
integer function coidef_numvar(cntvect, numvar)
defines the number of variables in the model.
integer function coidef_objcon(cntvect, objcon)
defines the Objective Constraint.
integer function coidef_numnz(cntvect, numnz)
defines the number of nonzero elements in the Jacobian.
integer function coidef_optdir(cntvect, optdir)
defines the Optimization Direction.
integer function coidef_numnlnz(cntvect, numnlnz)
defines the Number of Nonlinear Nonzeros.
integer function coidef_numcon(cntvect, numcon)
defines the number of constraints in the model.
integer function coi_solve(cntvect)
method for starting the solving process of CONOPT.
float rndx()
Defines a pseudo random number between 0 and 1.
program electron
Main program. A simple setup and call of CONOPT.
subroutine flog(msg, code)