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