CONOPT
Loading...
Searching...
No Matches
mp_leastsq12.f90
Go to the documentation of this file.
1!> @file mp_leastsq12.f90
2!! @ingroup FORTOPENMP_EXAMPLES
3!!
4!! The model is similar to leastsq, but we have added non-binding bounds
5!! on the res-variables so they cannot be eliminated and moved to the
6!! post-triangle. The efficiency bottlenecks have therefore been moved
7!! to other parts of CONOPT.
8!!
9!! We solve the following nonlinear least squares model:
10!!
11!! min sum(i, res(i)**2 )
12!!
13!! sum(j, a(i,j)*x(j) + b(i,j)*x(j)**2 ) + res(i) = obs(i)
14!!
15!! where a, b, and obs are known data, and res and x are the
16!! variables of the model.
17!!
18!! We also use an approximate Hessian.
19!!
20!!
21!! For more information about the individual callbacks, please have a look at the source code.
22
23#if defined(_WIN32) && !defined(_WIN64)
24#define dec_directives_win32
25#endif
26
27module lsq12
28 Integer, parameter :: Nobs = 700
29 Integer, parameter :: DimX = 500
30 real*8, dimension(Nobs,DimX) :: a, b
31 real*8, dimension(Nobs) :: obs
32end module lsq12
34
35REAL FUNCTION rndx( )
36!
37! Defines a pseudo random number between 0 and 1
38!
39 IMPLICIT NONE
40
41 Integer, save :: seed = 12359
42
43 seed = mod(seed*1027+25,1048576)
44 rndx = float(seed)/float(1048576)
45
46END FUNCTION rndx
47
48subroutine defdata
49!
50! Define values for A, B, and Obs
51!
52 Use lsq12
53 IMPLICIT NONE
54
55 Integer :: i, j
56 real*8, Parameter :: xtarg = -1.0
57 real*8, Parameter :: noise = 1.0
58 Real, External :: Rndx
59 real*8 :: o
60
61 do i = 1, nobs
62 o = 0.d0
63 do j = 1, dimx
64 a(i,j) = rndx()
65 b(i,j) = rndx()
66 o = o + a(i,j) * xtarg + b(i,j) * xtarg**2
67 enddo
68 obs(i) = o + noise * rndx()
69 enddo
70
71end subroutine defdata
72!
73! Main program.
74!
75!> Main program. A simple setup and call of CONOPT
76!!
77Program leastsquare
78
80 Use conopt
81 Use lsq12
82 Use omp_lib
83 implicit None
84!
85! Declare the user callback routines as Integer, External:
86!
87 Integer, External :: lsq_readmatrix ! Mandatory Matrix definition routine defined below
88 Integer, External :: lsq_fdeval ! Function and Derivative evaluation routine
89 ! needed a nonlinear model.
90 Integer, External :: std_status ! Standard callback for displaying solution status
91 Integer, External :: std_solution ! Standard callback for displaying solution values
92 Integer, External :: std_message ! Standard callback for managing messages
93 Integer, External :: std_errmsg ! Standard callback for managing error messages
94 Integer, External :: lsq_2dlagrstr ! Callback for second derivatives structure
95 Integer, External :: lsq_2dlagrval ! Callback for second derivatives values
96#ifdef dec_directives_win32
97!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
98!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
99!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Status
100!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Solution
101!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
102!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
103!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DLagrStr
104!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DLagrVal
105#endif
106!
107! Control vector
108!
109 INTEGER, Dimension(:), Pointer :: cntvect
110 INTEGER :: coi_error
111!
112 real*8 time0, time1, time2
113
114 call startup
115!
116! Define data
117!
118 Call defdata
119!
120! Create and initialize a Control Vector
121!
122 coi_error = coi_create( cntvect )
123!
124! Tell CONOPT about the size of the model by populating the Control Vector:
125!
126 coi_error = max( coi_error, coidef_numvar( cntvect, nobs + dimx ) )
127 coi_error = max( coi_error, coidef_numcon( cntvect, nobs + 1 ) )
128 coi_error = max( coi_error, coidef_numnz( cntvect, nobs * (dimx+2) ) )
129 coi_error = max( coi_error, coidef_numnlnz( cntvect, nobs * (dimx+1) ) )
130 coi_error = max( coi_error, coidef_numhess( cntvect, nobs + dimx ) )
131 coi_error = max( coi_error, coidef_optdir( cntvect, -1 ) ) ! Minimize
132 coi_error = max( coi_error, coidef_objcon( cntvect, nobs + 1 ) ) ! Objective is last constraint
133 coi_error = max( coi_error, coidef_optfile( cntvect, 'mp_leastsq12.opt' ) )
134!
135! Tell CONOPT about the callback routines:
136!
137 coi_error = max( coi_error, coidef_readmatrix( cntvect, lsq_readmatrix ) )
138 coi_error = max( coi_error, coidef_fdeval( cntvect, lsq_fdeval ) )
139 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
140 coi_error = max( coi_error, coidef_solution( cntvect, std_solution ) )
141 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
142 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
143 coi_error = max( coi_error, coidef_2dlagrstr( cntvect, lsq_2dlagrstr ) )
144 coi_error = max( coi_error, coidef_2dlagrval( cntvect, lsq_2dlagrval ) )
145
146#if defined(CONOPT_LICENSE_INT_1) && defined(CONOPT_LICENSE_INT_2) && defined(CONOPT_LICENSE_INT_3) && defined(CONOPT_LICENSE_TEXT)
147 coi_error = max( coi_error, coidef_license( cntvect, conopt_license_int_1, conopt_license_int_2, conopt_license_int_3, conopt_license_text) )
148#endif
149
150 If ( coi_error .ne. 0 ) THEN
151 write(*,*)
152 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
153 write(*,*)
154 call flog( "Skipping Solve due to setup errors", 1 )
155 ENDIF
156!
157! Start CONOPT -- First time using a single thread:
158!
159 time0 = omp_get_wtime()
160 coi_error = coi_solve( cntvect )
161 time1 = omp_get_wtime() - time0
162
163 write(*,*)
164 write(*,*) 'End of Least Square example 12 with 1 thread. Return code=',coi_error
165
166 If ( coi_error /= 0 ) then
167 call flog( "One Thread: Errors encountered during solution", 1 )
168 elseif ( stacalls == 0 .or. solcalls == 0 ) then
169 call flog( "One Thread: Status or Solution routine was not called", 1 )
170 elseif ( .not. ( sstat == 1 .and. mstat == 2 ) ) then
171 call flog( "One Thread: Solver or Model status was not as expected (1,2)", 1 )
172 endif
173!
174! Start CONOPT again -- this time using multiple threads:
175!
176 coi_error = max( coi_error, coidef_threads( cntvect, 0 ) ) ! 0 means use as many threads as you can
177 time0 = omp_get_wtime()
178 coi_error = coi_solve( cntvect )
179 time2 = omp_get_wtime() - time0
180
181 write(*,*)
182 write(*,*) 'Multi thread: End of Least Square example 12. Return code=',coi_error
183
184 If ( coi_error /= 0 ) then
185 call flog( "Multi thread: Errors encountered during solution", 1 )
186 elseif ( stacalls == 0 .or. solcalls == 0 ) then
187 call flog( "Multi thread: Status or Solution routine was not called", 1 )
188 elseif ( .not. ( sstat == 1 .and. mstat == 2 ) ) then
189 call flog( "Multi thread: Solver or Model status was not as expected (1,2)", 1 )
190 endif
191
192 write(*,*)
193 write(*,*) 'End of Least Square example 12. Return code=',coi_error
194
195 if ( coi_free( cntvect ) /= 0 ) call flog( "Error while freeing control vector", 1 )
196
197 write(*,*)
198 write(*,"('Time for single thread',f10.3)") time1
199 write(*,"('Time for multi thread',f10.3)") time2
200 write(*,"('Speedup ',f10.3)") time1/time2
201 write(*,"('Efficiency ',f10.3)") time1/time2/omp_get_max_threads()
202
203 call flog( "Successful Solve", 0 )
204
205End Program leastsquare
206!
207! ============================================================================
208! Define information about the model:
209!
210
211!> Define information about the model
212!!
213!! @include{doc} readMatrix_params.dox
214Integer Function lsq_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
215 colsta, rowno, value, nlflag, n, m, nz, &
216 usrmem )
217#ifdef dec_directives_win32
218!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
219#endif
220 Use lsq12
221 implicit none
222 integer, intent (in) :: n ! number of variables
223 integer, intent (in) :: m ! number of constraints
224 integer, intent (in) :: nz ! number of nonzeros
225 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
226 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
227 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
228 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
229 ! (not defined here)
230 integer, intent (out), dimension(m) :: type ! vector of equation types
231 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
232 ! (not defined here)
233 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
234 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
235 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
236 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
237 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
238 real*8 usrmem(*) ! optional user memory
239
240 Integer :: i, j, k
241!
242! Information about Variables:
243! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
244! Default: the status information in Vsta is not used.
245!
246 do i = 1, dimx
247 curr(i) = -0.8d0
248 enddo
249 do i = 1, nobs ! prevent post triangle from being identified.
250 lower(dimx+i) = -1000.0d0
251 upper(dimx+i) = +1000.0d0
252 enddo
253!
254! Information about Constraints:
255! Default: Rhs = 0
256! Default: the status information in Esta and the function
257! value in FV are not used.
258! Default: Type: There is no default.
259! 0 = Equality,
260! 1 = Greater than or equal,
261! 2 = Less than or equal,
262! 3 = Non binding.
263!
264! Constraints 1 to Nobs:
265! Rhs = Obs(i) and type Equality
266!
267 do i = 1, nobs
268 rhs(i) = obs(i)
269 type(i) = 0
270 enddo
271!
272! Constraint Nobs + 1 (Objective)
273! Rhs = 0 and type Non binding
274!
275 type(nobs+1) = 3
276!
277! Information about the Jacobian. CONOPT expects a columnwise
278! representation in Rowno, Value, Nlflag and Colsta.
279!
280! Colsta = Start of column indices (No Defaults):
281! Rowno = Row indices
282! Value = Value of derivative (by default only linear
283! derivatives are used)
284! Nlflag = 0 for linear and 1 for nonlinear derivative
285! (not needed for completely linear models)
286!
287!
288! Indices
289! x(i) res(j)
290! j: NL L=1
291! obj: NL
292!
293 k = 1
294 do j = 1, dimx
295 colsta(j) = k
296 do i = 1, nobs
297 rowno(k) = i
298 nlflag(k) = 1
299 k = k + 1
300 enddo
301 enddo
302 do i = 1, nobs
303 colsta(dimx+i) = k
304 rowno(k) = i
305 nlflag(k) = 0
306 value(k) = 1.d0
307 k = k + 1
308 rowno(k) = nobs+1
309 nlflag(k) = 1
310 k = k + 1
311 enddo
312 colsta(dimx+nobs+1) = k
313
314 lsq_readmatrix = 0 ! Return value means OK
315
316end Function lsq_readmatrix
317!
318!==========================================================================
319! Compute nonlinear terms and non-constant Jacobian elements
320!
321
322!> Compute nonlinear terms and non-constant Jacobian elements
323!!
324!! @include{doc} fdeval_params.dox
325Integer Function lsq_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
326 n, nz, thread, usrmem )
327#ifdef dec_directives_win32
328!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
329#endif
330 use lsq12
331 implicit none
332 integer, intent (in) :: n ! number of variables
333 integer, intent (in) :: rowno ! number of the row to be evaluated
334 integer, intent (in) :: nz ! number of nonzeros in this row
335 real*8, intent (in), dimension(n) :: x ! vector of current solution values
336 real*8, intent (in out) :: g ! constraint value
337 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
338 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
339 ! in this row. Ffor information only.
340 integer, intent (in) :: mode ! evaluation mode: 1 = function value
341 ! 2 = derivatives, 3 = both
342 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
343 ! as errcnt is incremented
344 integer, intent (in out) :: errcnt ! error counter to be incremented in case
345 ! of function evaluation errors.
346 integer, intent (in) :: thread
347 real*8 usrmem(*) ! optional user memory
348
349 integer :: i,j
350 real*8 :: s
351!
352! Row Nobs+1: the objective function is nonlinear
353!
354 if ( rowno .eq. nobs+1 ) then
355!
356! Mode = 1 or 3. Function value: G = sum(i, res(i)**2)
357!
358 if ( mode .eq. 1 .or. mode .eq. 3 ) then
359 s = 0.d0
360 do i = 1, nobs
361 s = s + x(dimx+i)**2
362 enddo
363 g = s
364 endif
365!
366! Mode = 2 or 3: Derivative values:
367!
368 if ( mode .eq. 2 .or. mode .eq. 3 ) then
369 do i = 1, nobs
370 jac(dimx+i) = 2*x(dimx+i)
371 enddo
372 endif
373!
374! Row 1 to Nobs: The observation definitions:
375!
376 else
377!
378! Mode = 1 or 3: Function value - x-part only
379!
380 if ( mode .eq. 1 .or. mode .eq. 3 ) then
381 s = 0.d0
382 do j = 1, dimx
383 s = s + a(rowno,j)*x(j) + b(rowno,j)*x(j)**2
384 enddo
385 g = s
386 endif
388! Mode = 2 or 3: Derivatives
389!
390 if ( mode .eq. 2 .or. mode .eq. 3 ) then
391 do j = 1, dimx
392 jac(j) = a(rowno,j) + 2.d0*b(rowno,j)*x(j)
393 enddo
394 endif
395 endif
396 lsq_fdeval = 0
397
398end Function lsq_fdeval
399
400
401!> Specify the structure of the Lagrangian of the Hessian
402!!
403!! @include{doc} 2DLagrStr_params.dox
404Integer Function lsq_2dlagrstr( HSRW, HSCL, NODRV, N, M, NHESS, UsrMem )
405#ifdef dec_directives_win32
406!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DLagrStr
407#endif
408 use lsq12
409 Implicit None
410
411 Integer, Intent (IN) :: n, m, nhess
412 Integer, Intent (IN OUT) :: nodrv
413 Integer, Dimension(Nhess), Intent (Out) :: hsrw, hscl
414 real*8, Intent(IN OUT) :: usrmem(*)
415
416 integer :: i
417!
418! We assume that this is a Large Residual model, where the most
419! important 2nd derivatives are those that appear directly in the
420! objective.
421! We will therefore only define the 2nd derivatives corresponding to
422! the objective constraint where the 2nd derivatives is 2*I.
423! CONOPT will by default complain that some nonlinear variables do
424! not appear in the Hessian. We will therefore define one dummy element
425! on the diagonal of the Hessian for each of the nonlinear X-variables
426! and give it the numerical value 0.
427!
428!
429! Define structure of Hessian
430!
431 do i = 1, dimx+nobs
432 hsrw(i) = i
433 hscl(i) = i
434 enddo
435
436 lsq_2dlagrstr = 0
437
438End Function lsq_2dlagrstr
439
440
441!> Compute the Lagrangian of the Hessian
442!!
443!! @include{doc} 2DLagrVal_params.dox
444Integer Function lsq_2dlagrval( X, U, HSRW, HSCL, HSVL, NODRV, N, M, NHESS, UsrMem )
445#ifdef dec_directives_win32
446!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DLagrVal
447#endif
448 use lsq12
449 Implicit None
450
451 Integer, Intent (IN) :: n, m, nhess
452 Integer, Intent (IN OUT) :: nodrv
453 real*8, Dimension(N), Intent (IN) :: x
454 real*8, Dimension(M), Intent (IN) :: u
455 Integer, Dimension(Nhess), Intent (In) :: hsrw, hscl
456 real*8, Dimension(NHess), Intent (Out) :: hsvl
457 real*8, Intent(IN OUT) :: usrmem(*)
458
459 integer :: i
460!
461! We assume that this is a Large Residual model, where the most
462! important 2nd derivatives are those that appear directly in the
463! objective.
464! We will therefore only define the 2nd derivatives corresponding to
465! the objective constraint where the 2nd derivatives is 2*I.
466! CONOPT will by default complain that some nonlinear variables do
467! not appear in the Hessian. We will therefore define one dummy element
468! on the diagonal of the Hessian for each of the nonlinear X-variables
469! and give it the numerical value 0.
470!
471! Normal Evaluation mode
472!
473 do i = 1, dimx
474 hsvl(i) = 0.d0
475 enddo
476 do i = 1, nobs
477 hsvl(dimx+i) = 2.d0*u(nobs+1)
478 enddo
479
480 lsq_2dlagrval = 0
481
482End Function lsq_2dlagrval
483
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(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_2dlagrstr(cntvect, coi_2dlagrstr)
define callback routine for providing the structure of the second derivatives of the Lagrangian.
Definition conopt.f90:1547
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_2dlagrval(cntvect, coi_2dlagrval)
define callback routine for computing the values of the second derivatives of the Lagrangian.
Definition conopt.f90:1573
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_numhess(cntvect, numhess)
defines the Number of Hessian Nonzeros.
Definition conopt.f90:188
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
integer function lsq_readmatrix(lower, curr, upper, vsta, type, rhs, esta, colsta, rowno, value, nlflag, n, m, nz, usrmem)
Define information about the model.
Definition leastl1.f90:166
integer function lsq_fdeval(x, g, jac, rowno, jcnm, mode, ignerr, errcnt, n, nz, thread, usrmem)
Compute nonlinear terms and non-constant Jacobian elements.
Definition leastl1.f90:287
integer function lsq_2dlagrstr(hsrw, hscl, nodrv, n, m, nhess, usrmem)
Specify the structure of the Lagrangian of the Hessian.
Definition leastsq2.f90:359
integer function lsq_2dlagrval(x, u, hsrw, hscl, hsvl, nodrv, n, m, nhess, usrmem)
Compute the Lagrangian of the Hessian.
Definition leastsq2.f90:396
#define dimx
Definition leastsq.c:16
#define nobs
Definition leastsq.c:15
void defdata()
Defines the data for the problem.
Definition leastsq.c:35
float rndx()
Defines a pseudo random number between 0 and 1.
Definition leastsq.c:22
program leastsquare
Main program. A simple setup and call of CONOPT.
Definition leastsq.f90:68
real *8, dimension(nobs) obs
real *8, dimension(nobs, dimx) b
real *8, dimension(nobs, dimx) a
subroutine flog(msg, code)
Definition comdeclp.f90:48
subroutine startup
Definition comdeclp.f90:31