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