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