CONOPT
Loading...
Searching...
No Matches
mp_elec.f90
Go to the documentation of this file.
1!> @file mp_elec.f90
2!! @ingroup FORTOPENMP_EXAMPLES
3!!
4!! @copydoc elec.f90
5
6#if defined(_WIN32) && !defined(_WIN64)
7#define dec_directives_win32
8#endif
9
10Module random
11 Integer :: seed
12End Module random
13
14!> Main program. A simple setup and call of CONOPT
15!!
16Program electron
17
19 Use conopt
20 Use omp_lib
21 Use random
22 implicit None
23!
24! Declare the user callback routines as Integer, External:
25!
26 Integer, External :: elec_readmatrix ! Mandatory Matrix definition routine defined below
27 Integer, External :: elec_fdeval ! Function and Derivative evaluation routine
28 ! needed a nonlinear model.
29 Integer, External :: std_status ! Standard callback for displaying solution status
30 Integer, External :: std_solution ! Standard callback for displaying solution values
31 Integer, External :: std_message ! Standard callback for managing messages
32 Integer, External :: std_errmsg ! Standard callback for managing error messages
33#ifdef dec_directives_win32
34!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Elec_ReadMatrix
35!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Elec_FDEval
36!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Status
37!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Solution
38!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
39!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
40#endif
41!
42! Control vector
43!
44 INTEGER, Dimension(:), Pointer :: cntvect
45 INTEGER :: coi_error
46!
47! The model shown here in GAMS format
48!
49!
50!Set i electrons /i1 * i%np%/
51! ut(i,i) upper triangular part;
52!
53!Alias (i,j);
54!ut(i,j)$(ord(j) > ord(i)) = yes;
55!
56!Variables x(i) x-coordinate of the electron
57! y(i) y-coordinate of the electron
58! z(i) z-coordinate of the electron
59! potential Coulomb potential;
60!
61!Equations obj objective
62! ball(i) points on unit ball;
63!
64!obj.. potential =e=
65! sum{ut(i,j), 1.0/sqrt(sqr(x[i]-x[j]) + sqr(y[i]-y[j]) + sqr(z[i]-z[j]))};
66!
67!ball(i).. sqr(x(i)) + sqr(y(i)) + sqr(z(i)) =e= 1;
68!
69!
70!* Set the starting point to a quasi-uniform distribution
71!* of electrons on a unit sphere
72!
73!scalar pi a famous constant;
74!pi = 2*arctan(inf);
75!
76!parameter theta(i), phi(i);
77!theta(i) = 2*pi*uniform(0,1);
78!phi(i) = pi*uniform(0,1);
79!
80!x.l(i) = cos(theta(i))*sin(phi(i));
81!y.l(i) = sin(theta(i))*sin(phi(i));
82!z.l(i) = cos(phi(i));
83!
84 Integer :: ne ! Number of Electrons
85
86 real*8 time0, time1, time2
87!
88! Create and initialize a Control Vector
89!
90 call startup
91
92 coi_error = coi_create( cntvect )
93!
94! Tell CONOPT about the size of the model by populating the Control Vector:
95!
96 ne = 200 ! Number of electrons
97 coi_error = max( coi_error, coidef_numvar( cntvect, 3*ne ) ) ! # variables
98 coi_error = max( coi_error, coidef_numcon( cntvect, ne+1 ) ) ! # constraints
99 coi_error = max( coi_error, coidef_numnz( cntvect, 6*ne ) ) ! # nonzeros in the Jacobian
100 coi_error = max( coi_error, coidef_numnlnz( cntvect, 6*ne ) ) ! # of which are nonlinear
101 coi_error = max( coi_error, coidef_optdir( cntvect, -1 ) ) ! Minimize
102 coi_error = max( coi_error, coidef_objcon( cntvect, ne+1 ) ) ! Objective is constraint NE+1
103 coi_error = max( coi_error, coidef_optfile( cntvect, 'mp_elec.opt' ) )
104!
105! Tell CONOPT about the callback routines:
106!
107 coi_error = max( coi_error, coidef_readmatrix( cntvect, elec_readmatrix ) )
108 coi_error = max( coi_error, coidef_fdeval( cntvect, elec_fdeval ) )
109 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
110 coi_error = max( coi_error, coidef_solution( cntvect, std_solution ) )
111 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
112 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
113
114#if defined(CONOPT_LICENSE_INT_1) && defined(CONOPT_LICENSE_INT_2) && defined(CONOPT_LICENSE_INT_3) && defined(CONOPT_LICENSE_TEXT)
115 coi_error = max( coi_error, coidef_license( cntvect, conopt_license_int_1, conopt_license_int_2, conopt_license_int_3, conopt_license_text) )
116#endif
117
118 If ( coi_error .ne. 0 ) THEN
119 write(*,*)
120 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
121 write(*,*)
122 call flog( "Skipping Solve due to setup errors", 1 )
123 ENDIF
124!
125! Start CONOPT -- First time using a single thread:
126!
127 seed = 12345
128 time0 = omp_get_wtime()
129 coi_error = coi_solve( cntvect )
130 time1 = omp_get_wtime() - time0
131
132 write(*,*)
133 write(*,*) 'End of Electron example single thread. Return code=',coi_error
134
135 If ( coi_error /= 0 ) then
136 call flog( "One Thread: Errors encountered during solution", 1 )
137 elseif ( stacalls == 0 .or. solcalls == 0 ) then
138 call flog( "One Thread: Status or Solution routine was not called", 1 )
139 elseif ( sstat /= 1 .or. mstat /= 2 ) then
140 call flog( "One Thread: Solver and Model Status was not as expected (1,2)", 1 )
141 endif
142!
143! Start CONOPT again -- this time using multiple threads:
144!
145 seed = 12345 ! Make sure we get the same model the second time
146 coi_error = max( coi_error, coidef_threads( cntvect, 0 ) ) ! 0 means use as many threads as you can
147 time0 = omp_get_wtime()
148 coi_error = coi_solve( cntvect )
149 time2 = omp_get_wtime() - time0
150
151 write(*,*)
152 write(*,*) 'End of Electron example multi thread. Return code=',coi_error
153
154 If ( coi_error /= 0 ) then
155 call flog( "Multi Thread: Errors encountered during solution", 1 )
156 elseif ( stacalls == 0 .or. solcalls == 0 ) then
157 call flog( "Multi Thread: Status or Solution routine was not called", 1 )
158 elseif ( sstat /= 1 .or. mstat /= 2 ) then
159 call flog( "Multi Thread: Solver and Model Status was not as expected (1,2)", 1 )
160 endif
161
162 if ( coi_free( cntvect ) /= 0 ) call flog( "Error while freeing control vector", 1 )
163
164 write(*,*)
165 write(*,"('Time for single thread',f10.3)") time1
166 write(*,"('Time for multi thread',f10.3)") time2
167 write(*,"('Speedup ',f10.3)") time1/time2
168 write(*,"('Efficiency ',f10.3)") time1/time2/omp_get_max_threads()
169
170 call flog( "Successful Solve", 0 )
171
172End Program electron
173
174REAL FUNCTION rndx( )
175!
176! Defines a pseudo random number between 0 and 1
177!
178 Use random
179 IMPLICIT NONE
180
181 seed = mod(seed*1027+25,1048576)
182 rndx = float(seed)/float(1048576)
183
184END FUNCTION rndx
185!
186! ============================================================================
187! Define information about the model:
188!
189
190!> Define information about the model
191!!
192!! @include{doc} readMatrix_params.dox
193Integer Function elec_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
194 colsta, rowno, value, nlflag, n, m, nz, &
195 usrmem )
196#ifdef dec_directives_win32
197!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Elec_ReadMatrix
198#endif
199 implicit none
200 integer, intent (in) :: n ! number of variables
201 integer, intent (in) :: m ! number of constraints
202 integer, intent (in) :: nz ! number of nonzeros
203 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
204 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
205 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
206 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
207 ! (not defined here)
208 integer, intent (out), dimension(m) :: type ! vector of equation types
209 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
210 ! (not defined here)
211 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
212 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
213 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
214 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
215 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
216 real*8 usrmem(*) ! optional user memory
217
218 Integer :: ne
219 Integer :: i, k
220 real*8, parameter :: pi = 3.141592
221 real*8 :: theta, phi
222 Real, External :: rndx
223
224 ne = n / 3
225!
226! Information about Variables:
227! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
228! Default: the status information in Vsta is not used.
229!
230 DO i = 1, ne
231 theta = rndx()
232 phi = rndx()
233 curr(i ) = cos(theta)*sin(phi)
234 curr(i+ ne) = sin(theta)*sin(phi)
235 curr(i+2*ne) = cos(phi)
236 enddo
237!
238! Information about Constraints:
239! Default: Rhs = 0
240! Default: the status information in Esta and the function
241! value in FV are not used.
242! Default: Type: There is no default.
243! 0 = Equality,
244! 1 = Greater than or equal,
245! 2 = Less than or equal,
246! 3 = Non binding.
247!
248! Constraint 1 (Objective)
249! Rhs = 0.0 and type Non binding
250!
251 DO i = 1, ne
252 Type(i) = 0
253 rhs(i) = 1.d0
254 Enddo
255 type(ne+1) = 3
256!
257! Information about the Jacobian. We have to define Rowno, Value,
258! Nlflag and Colsta.
259!
260! Colsta = Start of column indices (No Defaults):
261! Rowno = Row indices
262! Value = Value of derivative (by default only linear
263! derivatives are used)
264! Nlflag = 0 for linear and 1 for nonlinear derivative
265! (not needed for completely linear models)
266!
267! Indices
268! x(1) x(2) .. y(1) y(2) .. z(1) z(2) ..
269! 1: 1 2*Ne+1 3*Ne+1
270! 2: 3 2*Ne+3 3*Ne+3
271! ..
272! Obj: 2 4 2*Ne+2 2*Ne+4 3*Ne+2 3*Ne+4
273!
274! Nonlinearity Structure: L = 0 are linear and NL = 1 are nonlinear
275! All nonzeros are nonlinear
276!
277 DO i = 1, n+1
278 colsta(i) = 2*i-1
279 enddo
280 k = 0
281 DO i = 1, ne ! x variablex
282 k = k + 1
283 rowno(k) = i ! Ball constraint
284 nlflag(k) = 1
285 k = k + 1
286 rowno(k) = ne+1 ! Objective
287 nlflag(k) = 1
288 enddo
289 DO i = 1, ne ! y variablex
290 k = k + 1
291 rowno(k) = i ! Ball constraint
292 nlflag(k) = 1
293 k = k + 1
294 rowno(k) = ne+1 ! Objective
295 nlflag(k) = 1
296 enddo
297 DO i = 1, ne ! z variablex
298 k = k + 1
299 rowno(k) = i ! Ball constraint
300 nlflag(k) = 1
301 k = k + 1
302 rowno(k) = ne+1 ! Objective
303 nlflag(k) = 1
304 enddo
306 elec_readmatrix = 0 ! Return value means OK
307
308end Function elec_readmatrix
309!
310!==========================================================================
311! Compute nonlinear terms and non-constant Jacobian elements
312!
313
314!> Compute nonlinear terms and non-constant Jacobian elements
315!!
316!! @include{doc} fdeval_params.dox
317Integer Function elec_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
318 n, nz, thread, usrmem )
319#ifdef dec_directives_win32
320!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Elec_FDEval
321#endif
322 implicit none
323 integer, intent (in) :: n ! number of variables
324 integer, intent (in) :: rowno ! number of the row to be evaluated
325 integer, intent (in) :: nz ! number of nonzeros in this row
326 real*8, intent (in), dimension(n) :: x ! vector of current solution values
327 real*8, intent (in out) :: g ! constraint value
328 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
329 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
330 ! in this row. Ffor information only.
331 integer, intent (in) :: mode ! evaluation mode: 1 = function value
332 ! 2 = derivatives, 3 = both
333 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
334 ! as errcnt is incremented
335 integer, intent (in out) :: errcnt ! error counter to be incremented in case
336 ! of function evaluation errors.
337 integer, intent (in) :: thread
338 real*8 usrmem(*) ! optional user memory
339
340 Integer :: ne
341 Integer :: i, j
342 real*8 :: dist, dist32, tx, ty, tz
343
344 ne = n / 3
345!
346! Row NE+1: the objective function is nonlinear
347!
348 if ( rowno == ne+1 ) then
349!
350! Mode = 1 or 3. Function value:
351!
352 if ( mode .eq. 1 .or. mode .eq. 3 ) then
353 g = 0.d0
354 do i = 2, ne
355 do j = 1, i-1
356 dist = (x(i)-x(j))**2 + (x(i+ne)-x(j+ne))**2 + (x(i+2*ne)-x(j+2*ne))**2
357 g = g + 1.d0/sqrt(dist)
358 enddo
359 enddo
360 endif
361!
362! Mode = 2 or 3: Derivative values: w.r.t. x: 2*(x(i)-x(j))*(-1/2)*dist(-3/2)
363!
364 if ( mode .eq. 2 .or. mode .eq. 3 ) then
365 do i = 1, 3*ne
366 jac(i) = 0.d0
367 enddo
368 do i = 2, ne
369 do j = 1, i-1
370 dist = (x(i)-x(j))**2 + (x(i+ne)-x(j+ne))**2 + (x(i+2*ne)-x(j+2*ne))**2
371 dist32 = (1.d0/sqrt(dist))**3
372 tx = -(x(i)-x(j))*dist32
373 ty = -(x(i+ne)-x(j+ne))*dist32
374 tz = -(x(i+2*ne)-x(j+2*ne))*dist32
375 jac(i) = jac(i) + tx
376 jac(j) = jac(j) - tx
377 jac(i+ne) = jac(i+ne) + ty
378 jac(j+ne) = jac(j+ne) - ty
379 jac(i+2*ne) = jac(i+2*ne) + tz
380 jac(j+2*ne) = jac(j+2*ne) - tz
381 enddo
382 enddo
383 endif
384!
385 Else ! this is ball constraint rowno
386!
387! Mode = 1 or 3. Function value:
388!
389 i = rowno
390 if ( mode .eq. 1 .or. mode .eq. 3 ) then
391 g = x(i)**2 + x(i+ne)**2 + x(i+2*ne)**2
392 endif
393!
394! Mode = 2 or 3: Derivative values:
395!
396 if ( mode .eq. 2 .or. mode .eq. 3 ) then
397 jac(i) = 2.d0*x(i)
398 jac(i+ne) = 2.d0*x(i+ne)
399 jac(i+2*ne) = 2.d0*x(i+2*ne)
400 endif
401
402 endif
403 elec_fdeval = 0
404
405end Function elec_fdeval
integer function std_solution(xval, xmar, xbas, xsta, yval, ymar, ybas, ysta, n, m, usrmem)
Definition comdecl.f90:132
integer function std_status(modsta, solsta, iter, objval, usrmem)
Definition comdecl.f90:88
integer function std_message(smsg, dmsg, nmsg, llen, usrmem, msgv)
Definition comdecl.f90:205
integer function std_errmsg(rowno, colno, posno, msglen, usrmem, msg)
Definition comdecl.f90:248
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:277
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:157
integer(c_int) function coidef_message(cntvect, coi_message)
define callback routine for handling messages returned during the solution process.
Definition conopt.f90:1265
integer(c_int) function coidef_solution(cntvect, coi_solution)
define callback routine for returning the final solution values.
Definition conopt.f90:1238
integer(c_int) function coidef_status(cntvect, coi_status)
define callback routine for returning the completion status.
Definition conopt.f90:1212
integer(c_int) function coidef_readmatrix(cntvect, coi_readmatrix)
define callback routine for providing the matrix data to CONOPT.
Definition conopt.f90:1111
integer(c_int) function coidef_errmsg(cntvect, coi_errmsg)
define callback routine for returning error messages for row, column or Jacobian elements.
Definition conopt.f90:1291
integer(c_int) function coidef_fdeval(cntvect, coi_fdeval)
define callback routine for performing function and derivative evaluations.
Definition conopt.f90:1135
integer(c_int) function coidef_optfile(cntvect, optfile)
define callback routine for defining an options file.
Definition conopt.f90:928
integer(c_int) function coidef_license(cntvect, licint1, licint2, licint3, licstring)
define the License Information.
Definition conopt.f90:293
integer(c_int) function coidef_threads(cntvect, threads)
number of threads allowed internally in CONOPT.
Definition conopt.f90:627
integer(c_int) function coidef_numvar(cntvect, numvar)
defines the number of variables in the model.
Definition conopt.f90:97
integer(c_int) function coidef_numcon(cntvect, numcon)
defines the number of constraints in the model.
Definition conopt.f90:121
integer(c_int) function coidef_numnlnz(cntvect, numnlnz)
defines the Number of Nonlinear Nonzeros.
Definition conopt.f90:167
integer(c_int) function coidef_optdir(cntvect, optdir)
defines the Optimization Direction.
Definition conopt.f90:213
integer(c_int) function coidef_numnz(cntvect, numnz)
defines the number of nonzero elements in the Jacobian.
Definition conopt.f90:144
integer(c_int) function coidef_objcon(cntvect, objcon)
defines the Objective Constraint.
Definition conopt.f90:239
integer(c_int) function coi_create(cntvect)
initializes CONOPT and creates the control vector.
Definition conopt.f90:1726
integer(c_int) function coi_free(cntvect)
frees the control vector.
Definition conopt.f90:1749
integer(c_int) function coi_solve(cntvect)
method for starting the solving process of CONOPT.
Definition conopt.f90:1625
float rndx()
Defines a pseudo random number between 0 and 1.
Definition leastsq.c:22
program electron
Main program. A simple setup and call of CONOPT.
Definition mp_elec.f90:18
subroutine flog(msg, code)
Definition comdeclp.f90:48
subroutine startup
Definition comdeclp.f90:31