CONOPT
Loading...
Searching...
No Matches
mp_elec2.f90
Go to the documentation of this file.
1!> @file mp_elec2.f90
2!! @ingroup FORTOPENMP_EXAMPLES
3!!
4!! @copydoc elec2.f90
5
6Module hesscalls
7 Integer :: Hcalls
8End Module hesscalls
9
10Module random
11 Integer :: seed
12End Module random
13
14!> Main program. A simple setup and call of CONOPT
15!!
16Program electron
17
18 Use proginfop
19 Use coidef
20 Use hesscalls
21 Use omp_lib
22 Use random
23 implicit None
24!
25! Declare the user callback routines as Integer, External:
26!
27 Integer, External :: elec_readmatrix ! Mandatory Matrix definition routine defined below
28 Integer, External :: elec_fdeval ! Function and Derivative evaluation routine
29 ! needed a nonlinear model.
30 Integer, External :: std_status ! Standard callback for displaying solution status
31 Integer, External :: std_solution ! Standard callback for displaying solution values
32 Integer, External :: std_message ! Standard callback for managing messages
33 Integer, External :: std_errmsg ! Standard callback for managing error messages
34 Integer, External :: elec_2dlagrsize ! Callback for second derivatives size
35 Integer, External :: elec_2dlagrstr ! Callback for second derivatives structure
36 Integer, External :: elec_2dlagrval ! Callback for second derivatives values
37#if defined(itl)
38!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Elec_ReadMatrix
39!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Elec_FDEval
40!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Status
41!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Solution
42!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
43!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
44!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Elec_2DLagrSize
45!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Elec_2DLagrStr
46!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Elec_2DLagrVal
47#endif
48!
49! Control vector
50!
51 INTEGER, Dimension(:), Pointer :: cntvect
52 INTEGER :: coi_error
53!
54! The model shown here in GAMS format
55!
56!
57!Set i electrons /i1 * i%np%/
58! ut(i,i) upper triangular part;
59!
60!Alias (i,j);
61!ut(i,j)$(ord(j) > ord(i)) = yes;
62!
63!Variables x(i) x-coordinate of the electron
64! y(i) y-coordinate of the electron
65! z(i) z-coordinate of the electron
66! potential Coulomb potential;
67!
68!Equations obj objective
69! ball(i) points on unit ball;
70!
71!obj.. potential =e=
72! sum{ut(i,j), 1.0/sqrt(sqr(x[i]-x[j]) + sqr(y[i]-y[j]) + sqr(z[i]-z[j]))};
73!
74!ball(i).. sqr(x(i)) + sqr(y(i)) + sqr(z(i)) =e= 1;
75!
76!
77!* Set the starting point to a quasi-uniform distribution
78!* of electrons on a unit sphere
79!
80!scalar pi a famous constant;
81!pi = 2*arctan(inf);
82!
83!parameter theta(i), phi(i);
84!theta(i) = 2*pi*uniform(0,1);
85!phi(i) = pi*uniform(0,1);
86!
87!x.l(i) = cos(theta(i))*sin(phi(i));
88!y.l(i) = sin(theta(i))*sin(phi(i));
89!z.l(i) = cos(phi(i));
90!
91 Integer :: ne ! Number of Electrons
92
93 real*8 time0, time1, time2
94!
95! Create and initialize a Control Vector
96!
97 call startup
98
99 coi_error = coi_createfort( cntvect )
100!
101! Tell CONOPT about the size of the model by populating the Control Vector:
102!
103 ne = 150 ! Number of electrons
104 coi_error = max( coi_error, coidef_numvar( cntvect, 3*ne ) ) ! # variables
105 coi_error = max( coi_error, coidef_numcon( cntvect, ne+1 ) ) ! # constraints
106 coi_error = max( coi_error, coidef_numnz( cntvect, 6*ne ) ) ! # nonzeros in the Jacobian
107 coi_error = max( coi_error, coidef_numnlnz( cntvect, 6*ne ) ) ! # of which are nonlinear
108 coi_error = max( coi_error, coidef_optdir( cntvect, -1 ) ) ! Minimize
109 coi_error = max( coi_error, coidef_objcon( cntvect, ne+1 ) ) ! Objective is constraint NE+1
110 coi_error = max( coi_error, coidef_debug2d( cntvect, 0 ) ) ! turn 2nd derivative debugging on or off
111 coi_error = max( coi_error, coidef_optfile( cntvect, 'mp_elec2.opt' ) )
112!
113! Tell CONOPT about the callback routines:
114!
115 coi_error = max( coi_error, coidef_readmatrix( cntvect, elec_readmatrix ) )
116 coi_error = max( coi_error, coidef_fdeval( cntvect, elec_fdeval ) )
117 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
118 coi_error = max( coi_error, coidef_solution( cntvect, std_solution ) )
119 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
120 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
121 coi_error = max( coi_error, coidef_2dlagrsize( cntvect, elec_2dlagrsize ) )
122 coi_error = max( coi_error, coidef_2dlagrstr( cntvect, elec_2dlagrstr ) )
123 coi_error = max( coi_error, coidef_2dlagrval( cntvect, elec_2dlagrval ) )
124
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) )
127#endif
128
129 If ( coi_error .ne. 0 ) THEN
130 write(*,*)
131 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
132 write(*,*)
133 call flog( "Skipping Solve due to setup errors", 1 )
134 ENDIF
135!
136! Start CONOPT with enough space for the Hessian of the Lagrangian -- first time single thread
137!
138 coi_error = max( coi_error, coidef_hessfac( cntvect, 3.01d0*(ne+1)/4.d0 ) ) ! enough to allow Hessian to be use
139 hcalls = 0
140 seed = 12345
141 time0 = omp_get_wtime()
142 coi_error = coi_solve( cntvect )
143 time1 = omp_get_wtime() - time0
144
145 write(*,*)
146 write(*,*) 'End of Electron example one thread. Return code=',coi_error
147
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 )
156 endif
157!
158! Start CONOPT again -- this time using multiple threads:
159!
160 seed = 12345 ! Make sure we get the same model the second time
161 coi_error = max( coi_error, coidef_threads( cntvect, 0 ) ) ! 0 means use as many threads as you can
162 time0 = omp_get_wtime()
163 coi_error = coi_solve( cntvect )
164 time2 = omp_get_wtime() - time0
165
166 write(*,*)
167 write(*,*) 'End of Electron example multi threads. Return code=',coi_error
168
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 )
177 endif
178
179 if ( coi_free( cntvect ) /= 0 ) call flog( "Error while freeing control vector", 1 )
180
181 write(*,*)
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()
186
187 call flog( "Successful Solve", 0 )
188
189End Program electron
190
191REAL FUNCTION rndx( )
192!
193! Defines a pseudo random number between 0 and 1
194!
195 Use random
196 IMPLICIT NONE
197
198 seed = mod(seed*1027+25,1048576)
199 rndx = float(seed)/float(1048576)
200
201END FUNCTION rndx
202!
203! ============================================================================
204! Define information about the model:
205!
206
207!> Define information about the model
208!!
209!! @include{doc} readMatrix_params.dox
210Integer Function elec_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
211 colsta, rowno, value, nlflag, n, m, nz, &
212 usrmem )
213#if defined(itl)
214!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Elec_ReadMatrix
215#endif
216 implicit none
217 integer, intent (in) :: n ! number of variables
218 integer, intent (in) :: m ! number of constraints
219 integer, intent (in) :: nz ! number of nonzeros
220 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
221 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
222 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
223 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
224 ! (not defined here)
225 integer, intent (out), dimension(m) :: type ! vector of equation types
226 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
227 ! (not defined here)
228 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
229 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
230 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
231 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
232 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
233 real*8 usrmem(*) ! optional user memory
234
235 Integer :: ne
236 Integer :: i, k
237 real*8, parameter :: pi = 3.141592
238 real*8 :: theta, phi
239 Real, External :: rndx
240
241 ne = n / 3
242!
243! Information about Variables:
244! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
245! Default: the status information in Vsta is not used.
246!
247 DO i = 1, ne
248 theta = rndx()
249 phi = rndx()
250 curr(i ) = cos(theta)*sin(phi)
251 curr(i+ ne) = sin(theta)*sin(phi)
252 curr(i+2*ne) = cos(phi)
253 enddo
254!
255! Information about Constraints:
256! Default: Rhs = 0
257! Default: the status information in Esta and the function
258! value in FV are not used.
259! Default: Type: There is no default.
260! 0 = Equality,
261! 1 = Greater than or equal,
262! 2 = Less than or equal,
263! 3 = Non binding.
264!
265! Constraint 1 (Objective)
266! Rhs = 0.0 and type Non binding
267!
268 DO i = 1, ne
269 Type(i) = 0
270 rhs(i) = 1.d0
271 Enddo
272 type(ne+1) = 3
273!
274! Information about the Jacobian. We have to define Rowno, Value,
275! Nlflag and Colsta.
276!
277! Colsta = Start of column indices (No Defaults):
278! Rowno = Row indices
279! Value = Value of derivative (by default only linear
280! derivatives are used)
281! Nlflag = 0 for linear and 1 for nonlinear derivative
282! (not needed for completely linear models)
283!
284! Indices
285! x(1) x(2) .. y(1) y(2) .. z(1) z(2) ..
286! 1: 1 2*Ne+1 3*Ne+1
287! 2: 3 2*Ne+3 3*Ne+3
288! ..
289! Obj: 2 4 2*Ne+2 2*Ne+4 3*Ne+2 3*Ne+4
290!
291! Nonlinearity Structure: L = 0 are linear and NL = 1 are nonlinear
292! All nonzeros are nonlinear
293!
294 DO i = 1, n+1
295 colsta(i) = 2*i-1
296 enddo
297 k = 0
298 DO i = 1, ne ! x variablex
299 k = k + 1
300 rowno(k) = i ! Ball constraint
301 nlflag(k) = 1
302 k = k + 1
303 rowno(k) = ne+1 ! Objective
304 nlflag(k) = 1
305 enddo
306 DO i = 1, ne ! y variablex
307 k = k + 1
308 rowno(k) = i ! Ball constraint
309 nlflag(k) = 1
310 k = k + 1
311 rowno(k) = ne+1 ! Objective
312 nlflag(k) = 1
313 enddo
314 DO i = 1, ne ! z variablex
315 k = k + 1
316 rowno(k) = i ! Ball constraint
317 nlflag(k) = 1
318 k = k + 1
319 rowno(k) = ne+1 ! Objective
320 nlflag(k) = 1
321 enddo
322
323 elec_readmatrix = 0 ! Return value means OK
324
325end Function elec_readmatrix
326!
327!==========================================================================
328! Compute nonlinear terms and non-constant Jacobian elements
329!
330
331!> Compute nonlinear terms and non-constant Jacobian elements
332!!
333!! @include{doc} fdeval_params.dox
334Integer Function elec_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
335 n, nz, thread, usrmem )
336#if defined(itl)
337!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Elec_FDEval
338#endif
339 implicit none
340 integer, intent (in) :: n ! number of variables
341 integer, intent (in) :: rowno ! number of the row to be evaluated
342 integer, intent (in) :: nz ! number of nonzeros in this row
343 real*8, intent (in), dimension(n) :: x ! vector of current solution values
344 real*8, intent (in out) :: g ! constraint value
345 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
346 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
347 ! in this row. Ffor information only.
348 integer, intent (in) :: mode ! evaluation mode: 1 = function value
349 ! 2 = derivatives, 3 = both
350 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
351 ! as errcnt is incremented
352 integer, intent (in out) :: errcnt ! error counter to be incremented in case
353 ! of function evaluation errors.
354 integer, intent (in) :: thread
355 real*8 usrmem(*) ! optional user memory
356
357 Integer :: ne
358 Integer :: i, j
359 real*8 :: dist, dist32, tx, ty, tz
360
361 ne = n / 3
362!
363! Row NE+1: the objective function is nonlinear
364!
365 if ( rowno == ne+1 ) then
366!
367! Mode = 1 or 3. Function value:
368!
369 if ( mode .eq. 1 .or. mode .eq. 3 ) then
370 g = 0.d0
371 do i = 1, ne-1
372 do j = i+1,ne
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)
375 enddo
376 enddo
377 endif
378!
379! Mode = 2 or 3: Derivative values: w.r.t. x: 2*(x(i)-x(j))*(-1/2)*dist(-3/2)
380!
381 if ( mode .eq. 2 .or. mode .eq. 3 ) then
382 do i = 1, 3*ne
383 jac(i) = 0.d0
384 enddo
385 do i = 1, ne-1
386 do j = i+1,ne
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
392 jac(i) = jac(i) + tx
393 jac(j) = jac(j) - tx
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
398 enddo
399 enddo
400 endif
401!
402 Else ! this is ball constraint rowno
403!
404! Mode = 1 or 3. Function value:
405!
406 i = rowno
407 if ( mode .eq. 1 .or. mode .eq. 3 ) then
408 g = x(i)**2 + x(i+ne)**2 + x(i+2*ne)**2
409 endif
410!
411! Mode = 2 or 3: Derivative values:
412!
413 if ( mode .eq. 2 .or. mode .eq. 3 ) then
414 jac(i) = 2.d0*x(i)
415 jac(i+ne) = 2.d0*x(i+ne)
416 jac(i+2*ne) = 2.d0*x(i+2*ne)
417 endif
418
419 endif
420 elec_fdeval = 0
421
422end Function elec_fdeval
423
424Integer Function elec_2dlagrsize( NODRV, NumVar, NumCon, NumHess, MaxHess, UsrMem )
425#if defined(itl)
426!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Elec_2DLagrSize
427#endif
428 Use hesscalls
429 Implicit None
430
431 Integer, Intent (IN) :: numvar, numcon, maxhess
432 Integer, Intent (IN OUT) :: nodrv, numhess
433 real*8, Intent(IN OUT) :: usrmem(*)
434!
435! Sizing call -- the hessian is completely dense
436!
437 numhess = numvar * (numvar+1)/2
439
440End Function elec_2dlagrsize
441
442
443!> Specify the structure of the Lagrangian of the Hessian
444!!
445!! @include{doc} 2DLagrStr_params.dox
446Integer Function elec_2dlagrstr( HSRW, HSCL, NODRV, NumVar, NumCon, NumHess, UsrMem )
447#if defined(itl)
448!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Elec_2DLagrStr
449#endif
450 Use hesscalls
451 Implicit None
452
453 Integer, Intent (IN) :: numvar, numcon, numhess
454 Integer, Intent (IN) :: nodrv
455 Integer, Dimension(*) :: hsrw, hscl
456 real*8, Intent(IN OUT) :: usrmem(*)
457
458 integer :: i, j, k
459!
460! Define structure of Hessian: A dense lower-triangular matrix
461!
462 k = 0
463 do i = 1, numvar
464 do j = i, numvar
465 k = k + 1
466 hsrw(k) = j
467 hscl(k) = i
468 enddo
469 enddo
471
472End Function elec_2dlagrstr
473
474
475!> Compute the Lagrangian of the Hessian
476!!
477!! @include{doc} 2DLagrVal_params.dox
478Integer Function elec_2dlagrval( X, U, HSRW, HSCL, HSVL, NODRV, NumVar, NumCon, NumHess, UsrMem )
479#if defined(itl)
480!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Elec_2DLagrVal
481#endif
482 Use hesscalls
483 Implicit None
484
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(*)
492
493 Integer :: ne ! Number of electrons
494 integer :: k, l, ix, iy, iz, jx, jy, jz
495 real*8 :: tx, ty, tz, dist, distsq, dist3, dist5
496!
497! Evaluation of 2nd derivatives. Hsvl is zero on entry:
498!
499 ne = numvar / 3
500 hcalls = hcalls + 1
501 k = 0
502! do i = 1, N
503! do j = i, N
504! k = k + 1
505! Hsvl(k) = 0.d0
506! enddo
507! enddo
508 do l = 1, ne ! constraint no l with 3 elements on the diagonal, each with value 2.0
509 k = pos(l,l) ! 1+(l-1)*(2*N+2-l)/2 ! variable x(l)
510 hsvl(k) = hsvl(k) + 2.d0*u(l)
511 k = pos(l+ne,l+ne) ! 1+(l+NE-1)*(2*N+2-(l+NE))/2 ! variable y(l)
512 hsvl(k) = hsvl(k) + 2.d0*u(l)
513 k = pos(l+2*ne,l+2*ne) ! 1+(l+2*NE-1)*(2*N+2-(l+2*NE))/2 ! variable z(l)
514 hsvl(k) = hsvl(k) + 2.d0*u(l)
515 enddo
516 if ( u(ne+1) /= 0.0d0 ) Then ! Nonzero dual on the objective
517 do ix = 1, ne-1
518 iy = ix + ne
519 iz = iy + ne
520 do jx = ix+1, ne
521 jy = jx + ne
522 jz = jy + ne
523! Objective computation:
524! dist = (x(i)-x(j))**2 + (x(i+NE)-x(j+NE))**2 + (x(i+2*NE)-x(j+2*NE))**2
525! g = g + 1.d0/sqrt(dist)
526 tx = x(ix)-x(jx) ! x1 - x2
527 ty = x(iy)-x(jy) ! y1 - y2
528 tz = x(iz)-x(jz) ! z1 - z2
529! The first derivatives of tx, ty, and tz w.r.t. the x variables are +1 or -1
530! The second derivatives of tx, ty, and tz are zero
531 dist = tx**2 + ty**2 + tz**2
532! The first derivatives of dist are
533! d(dist)/dx = d(dist)/d(tx)*d(tx)/dx = 2*tx*(+/-) (+1 for x1 and -1 for x2)
534! and similarly for d(dist)/dy = 2*ty*(+/-1) and d(dist)/dz = 2*tz*(+/-1)
535! The second derivatives of dist are
536! d2(dist)/dxidxj = 2 if xi=xj and 0 otherwise
537! d2(dist)/dyidyj = 2 if yi=yj and 0 otherwise
538! d2(dist)/dzidzj = 2 if zi=zj and 0 otherwise
539! All cross 2nd derivatives dxdy, dxdz, and dydz are zero
540!
541 distsq = 1/sqrt(dist) ! = dist**(-0.5)
542 dist3 = distsq**3
543 dist5 = distsq**5
544! The first derivatives of distsq are
545! d(distsq)/dx = -0.5*dist**(-1.5)*d(dist)/dx = -dist**(-1.5)*tx*(+/-1) = -distsq**3*tx*(+/-1)
546! d(distsq)/dy = -0.5*dist**(-1.5)*d(dist)/dy = -dist**(-1.5)*ty*(+/-1) = -distsq**3*ty*(+/-1)
547! d(distsq)/dz = -0.5*dist**(-1.5)*d(dist)/dz = -dist**(-1.5)*tz*(+/-1) = -distsq**3*tz*(+/-1)
548! The second derivatives of distsq are
549! d2(distsq)/dxidxj = -1.5*dist**(-2.5)*d(dist)/dxj*tx*(+/-1) - dist**(-1.5)*d(tx)/dxj*(+/-1)
550! = -1.5*distsq**5*2*tx*tx*(+/-1) - distsq**3*(+/-1)
551! d2(distsq)/dxidyj = -1.5*dist**(-2.5)*d(dist)/dyj*tx*(+/-1) - dist**(-1.5)*d(tx)/dyj*(+/-1)
552! = -1.5*distsq**5*2*ty*tx*(+/-1)
553 k = pos(ix,ix)
554 hsvl(k) = hsvl(k) + 3.0d0*dist5*tx**2 - dist3
555 k = pos(jx,jx)
556 hsvl(k) = hsvl(k) + 3.0d0*dist5*tx**2 - dist3
557 k = pos(ix,jx)
558 hsvl(k) = hsvl(k) - 3.0d0*dist5*tx**2 + dist3
559 k = pos(iy,iy)
560 hsvl(k) = hsvl(k) + 3.0d0*dist5*ty**2 - dist3
561 k = pos(jy,jy)
562 hsvl(k) = hsvl(k) + 3.0d0*dist5*ty**2 - dist3
563 k = pos(iy,jy)
564 hsvl(k) = hsvl(k) - 3.0d0*dist5*ty**2 + dist3
565 k = pos(iz,iz)
566 hsvl(k) = hsvl(k) + 3.0d0*dist5*tz**2 - dist3
567 k = pos(jz,jz)
568 hsvl(k) = hsvl(k) + 3.0d0*dist5*tz**2 - dist3
569 k = pos(iz,jz)
570 hsvl(k) = hsvl(k) - 3.0d0*dist5*tz**2 + dist3
571!
572! Off-diagonal blocks
573!
574 k = pos(ix,iy)
575 hsvl(k) = hsvl(k) + 3.0d0*dist5*tx*ty
576 k = pos(ix,jy)
577 hsvl(k) = hsvl(k) - 3.0d0*dist5*tx*ty
578 k = pos(jx,iy)
579 hsvl(k) = hsvl(k) - 3.0d0*dist5*tx*ty
580 k = pos(jx,jy)
581 hsvl(k) = hsvl(k) + 3.0d0*dist5*tx*ty
582 k = pos(ix,iz)
583 hsvl(k) = hsvl(k) + 3.0d0*dist5*tx*tz
584 k = pos(ix,jz)
585 hsvl(k) = hsvl(k) - 3.0d0*dist5*tx*tz
586 k = pos(jx,iz)
587 hsvl(k) = hsvl(k) - 3.0d0*dist5*tx*tz
588 k = pos(jx,jz)
589 hsvl(k) = hsvl(k) + 3.0d0*dist5*tx*tz
590 k = pos(iy,iz)
591 hsvl(k) = hsvl(k) + 3.0d0*dist5*ty*tz
592 k = pos(iy,jz)
593 hsvl(k) = hsvl(k) - 3.0d0*dist5*ty*tz
594 k = pos(jy,iz)
595 hsvl(k) = hsvl(k) - 3.0d0*dist5*ty*tz
596 k = pos(jy,jz)
597 hsvl(k) = hsvl(k) + 3.0d0*dist5*ty*tz
598 enddo
599 enddo
600 endif
601
603
604 Contains
605
606 Integer function pos(i,j)
607 Integer, Intent(in) :: i, j
608 pos = 1+(i-1)*(2*numvar+2-i)/2 + (j-i) ! Hessian position for variable (j,i), j >= i
609 End function pos
610
611End Function elec_2dlagrval
integer function coidef_2dlagrsize(cntvect, coi_2dlagrsize)
integer function std_solution(xval, xmar, xbas, xsta, yval, ymar, ybas, ysta, n, m, usrmem)
Definition comdecl.f90:128
integer function std_status(modsta, solsta, iter, objval, usrmem)
Definition comdecl.f90:82
integer function std_message(smsg, dmsg, nmsg, llen, usrmem, msgv)
Definition comdecl.f90:203
integer function std_errmsg(rowno, colno, posno, msglen, usrmem, msg)
Definition comdecl.f90:248
integer function elec_2dlagrval(x, u, hsrw, hscl, hsvl, nodrv, numvar, numcon, numhess, usrmem)
Compute the Lagrangian of the Hessian.
Definition elec2.f90:483
integer function elec_2dlagrsize(nodrv, numvar, numcon, numhess, maxhess, usrmem)
Definition elec2.f90:429
integer function pos(i, j)
Definition elec2.f90:611
integer function elec_2dlagrstr(hsrw, hscl, nodrv, numvar, numcon, numhess, usrmem)
Specify the structure of the Lagrangian of the Hessian.
Definition elec2.f90:451
integer function elec_fdeval(x, g, jac, rowno, jcnm, mode, ignerr, errcnt, n, nz, thread, usrmem)
Compute nonlinear terms and non-constant Jacobian elements.
Definition elec.f90:288
integer function elec_readmatrix(lower, curr, upper, vsta, type, rhs, esta, colsta, rowno, value, nlflag, n, m, nz, usrmem)
Define information about the model.
Definition elec.f90:165
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.
Definition coistart.f90:680
integer function coidef_numvar(cntvect, numvar)
defines the number of variables in the model.
Definition coistart.f90:358
integer function coidef_objcon(cntvect, objcon)
defines the Objective Constraint.
Definition coistart.f90:629
integer function coidef_numnz(cntvect, numnz)
defines the number of nonzero elements in the Jacobian.
Definition coistart.f90:437
integer function coidef_optdir(cntvect, optdir)
defines the Optimization Direction.
Definition coistart.f90:552
integer function coidef_numnlnz(cntvect, numnlnz)
defines the Number of Nonlinear Nonzeros.
Definition coistart.f90:476
integer function coidef_numcon(cntvect, numcon)
defines the number of constraints in the model.
Definition coistart.f90:398
integer function coi_solve(cntvect)
method for starting the solving process of CONOPT.
Definition coistart.f90:14
float rndx()
Defines a pseudo random number between 0 and 1.
Definition leastsq.c:23
program electron
Main program. A simple setup and call of CONOPT.
Definition mp_elec.f90:12
integer hcalls
Definition elec2.f90:62
subroutine flog(msg, code)
Definition comdeclp.f90:42
subroutine startup
Definition comdeclp.f90:25