CONOPT
Loading...
Searching...
No Matches
mp_leastsq5.f90
Go to the documentation of this file.
1!> @file mp_leastsq5.f90
2!! @ingroup FORTOPENMP_EXAMPLES
3!!
4!! @copydoc leastsq5.f90
5
6module lsq5
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 lsq5
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 lsq5
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 lsq5
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 :: lsq_2ddir ! Directional 2nd derivative routine
70 Integer, External :: std_status ! Standard callback for displaying solution status
71 Integer, External :: std_solution ! Standard callback for displaying solution values
72 Integer, External :: std_message ! Standard callback for managing messages
73 Integer, External :: std_errmsg ! Standard callback for managing error messages
74#if defined(itl)
75!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
76!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
77!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DDir
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#endif
83!
84! Control vector
85!
86 INTEGER, Dimension(:), Pointer :: cntvect
87 INTEGER :: coi_error
88
89 real*8 time0, time1, time2, time3
90
91 call startup
92!
93! Define data
94!
95 Call defdata
96!
97! Create and initialize a Control Vector
98!
99 coi_error = coi_createfort( cntvect )
100!
101! Tell CONOPT about the size of the model by populating the Control Vector:
102!
103 coi_error = max( coi_error, coidef_numvar( cntvect, nobs + dimx ) )
104 coi_error = max( coi_error, coidef_numcon( cntvect, nobs + 1 ) )
105 coi_error = max( coi_error, coidef_numnz( cntvect, nobs * (dimx+2) ) )
106 coi_error = max( coi_error, coidef_numnlnz( cntvect, nobs * (dimx+1) ) )
107 coi_error = max( coi_error, coidef_optdir( cntvect, -1 ) ) ! Minimize
108 coi_error = max( coi_error, coidef_objcon( cntvect, nobs + 1 ) ) ! Objective is last constraint
109 coi_error = max( coi_error, coidef_optfile( cntvect, 'mp_leastsq5.opt' ) )
110!
111! Tell CONOPT about the callback routines:
112!
113 coi_error = max( coi_error, coidef_readmatrix( cntvect, lsq_readmatrix ) )
114 coi_error = max( coi_error, coidef_fdeval( cntvect, lsq_fdeval ) )
115 coi_error = max( coi_error, coidef_2ddir( cntvect, lsq_2ddir ) )
116 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
117 coi_error = max( coi_error, coidef_solution( cntvect, std_solution ) )
118 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
119 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
120 coi_error = max( coi_error, coidef_debugfv( cntvect, 0 ) ) ! Debug Fdeval on or off
121
122#if defined(LICENSE_INT_1) && defined(LICENSE_INT_2) && defined(LICENSE_INT_3) && defined(LICENSE_TEXT)
123 coi_error = max( coi_error, coidef_license( cntvect, license_int_1, license_int_2, license_int_3, license_text) )
124#endif
125
126 If ( coi_error .ne. 0 ) THEN
127 write(*,*)
128 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
129 write(*,*)
130 call flog( "Skipping Solve due to setup errors", 1 )
131 ENDIF
132!
133! Start CONOPT with a single thread:
134!
135 time0 = omp_get_wtime()
136 coi_error = coi_solve( cntvect )
137 time1 = omp_get_wtime() - time0
138
139 write(*,*)
140 write(*,*) 'Single Thread: End of Least Square example 5. Return code=',coi_error
141
142 If ( coi_error /= 0 ) then
143 call flog( "Single Thread: Errors encountered during solution", 1 )
144 elseif ( stacalls == 0 .or. solcalls == 0 ) then
145 call flog( "Single Thread: Status or Solution routine was not called", 1 )
146 elseif ( sstat /= 1 .or. mstat /= 2 ) then
147 call flog( "Single Thread: Solver and Model Status was not as expected (1,2)", 1 )
148 elseif ( abs( obj - 19.44434311d0 ) > 1.d-7 ) then
149 call flog( "Single Thread: Incorrect objective returned", 1 )
150 endif
151!
152! Start CONOPT again this time with multiple threads:
153!
154 coi_error = max( coi_error, coidef_threads( cntvect, 0 ) ) ! 0 means use as many threads as you can
155 time0 = omp_get_wtime()
156 coi_error = coi_solve( cntvect )
157 time2 = omp_get_wtime() - time0
158
159 write(*,*)
160 write(*,*) 'Multiple Threads: End of Least Square example 5. Return code=',coi_error
161
162 If ( coi_error /= 0 ) then
163 call flog( "Multiple Threads: Errors encountered during solution", 1 )
164 elseif ( stacalls == 0 .or. solcalls == 0 ) then
165 call flog( "Multiple Threads: Status or Solution routine was not called", 1 )
166 elseif ( sstat /= 1 .or. mstat /= 2 ) then
167 call flog( "Multiple Threads: Solver and Model Status was not as expected (1,2)", 1 )
168 elseif ( abs( obj - 19.44434311d0 ) > 1.d-7 ) then
169 call flog( "Multiple Threads: Incorrect objective returned", 1 )
170 endif
171!
172! Start CONOPT again this time with multiple threads except for 2DDir:
173!
174 coi_error = max( coi_error, coidef_thread2d( cntvect, 1 ) ) ! 1 means use only 1 thread for 2DDir
175 time0 = omp_get_wtime()
176 coi_error = coi_solve( cntvect )
177 time3 = omp_get_wtime() - time0
178
179 write(*,*)
180 write(*,*) 'Multiple Threads -2DDir: End of Least Square example 5. Return code=',coi_error
181
182 If ( coi_error /= 0 ) then
183 call flog( "Multiple Threads -2DDir: Errors encountered during solution", 1 )
184 elseif ( stacalls == 0 .or. solcalls == 0 ) then
185 call flog( "Multiple Threads -2DDir: Status or Solution routine was not called", 1 )
186 elseif ( sstat /= 1 .or. mstat /= 2 ) then
187 call flog( "Multiple Threads -2DDir: Solver and Model Status was not as expected (1,2)", 1 )
188 elseif ( abs( obj - 19.44434311d0 ) > 1.d-7 ) then
189 call flog( "Multiple Threads -2DDir: Incorrect objective returned", 1 )
190 endif
191
192 if ( coi_free( cntvect ) /= 0 ) call flog( "Error while freeing control vector", 1 )
193
194 write(*,*)
195 write(*,"('Time for single thread ',f10.3)") time1
196 write(*,"('Time for multi thread incl 2DDir',f10.3)") time2
197 write(*,"('Time for multi thread excl 2DDir',f10.3)") time3
198 write(*,"('Speedup incl 2DDir',f10.3)") time1/time2
199 write(*,"('Efficiency incl 2DDir',f10.3)") time1/time2/omp_get_max_threads()
200 write(*,"('Speedup excl 2DDir',f10.3)") time1/time3
201 write(*,"('Efficiency excl 2DDir',f10.3)") time1/time3/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#if defined(itl)
218!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
219#endif
220 Use lsq5
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!
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! 3:
289!
290 k = 1
291 do j = 1, dimx
292 colsta(j) = k
293 do i = 1, nobs
294 rowno(k) = i
295 nlflag(k) = 1
296 k = k + 1
297 enddo
298 enddo
299 do i = 1, nobs
300 colsta(dimx+i) = k
301 rowno(k) = i
302 nlflag(k) = 0
303 value(k) = 1.d0
304 k = k + 1
305 rowno(k) = nobs+1
306 nlflag(k) = 1
307 k = k + 1
308 enddo
309 colsta(dimx+nobs+1) = k
310
311 lsq_readmatrix = 0 ! Return value means OK
312
313end Function lsq_readmatrix
314!
315!==========================================================================
316! Compute nonlinear terms and non-constant Jacobian elements
317!
318
319!> Compute nonlinear terms and non-constant Jacobian elements
320!!
321!! @include{doc} fdeval_params.dox
322Integer Function lsq_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
323 n, nz, thread, usrmem )
324#if defined(itl)
325!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
326#endif
327 use lsq5
328 implicit none
329 integer, intent (in) :: n ! number of variables
330 integer, intent (in) :: rowno ! number of the row to be evaluated
331 integer, intent (in) :: nz ! number of nonzeros in this row
332 real*8, intent (in), dimension(n) :: x ! vector of current solution values
333 real*8, intent (in out) :: g ! constraint value
334 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
335 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
336 ! in this row. Ffor information only.
337 integer, intent (in) :: mode ! evaluation mode: 1 = function value
338 ! 2 = derivatives, 3 = both
339 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
340 ! as errcnt is incremented
341 integer, intent (in out) :: errcnt ! error counter to be incremented in case
342 ! of function evaluation errors.
343 integer, intent (in) :: thread
344 real*8 usrmem(*) ! optional user memory
345
346 integer :: i,j
347 real*8 :: s
348!
349! Row 1: the objective function is nonlinear
350!
351 if ( rowno .eq. nobs+1 ) then
352!
353! Mode = 1 or 3. Function value: G = sum(i, res(i)**2)
354!
355 if ( mode .eq. 1 .or. mode .eq. 3 ) then
356 s = 0.d0
357 do i = 1, nobs
358 s = s + x(dimx+i)**2
359 enddo
360 g = s
361 endif
362!
363! Mode = 2 or 3: Derivative values:
364!
365 if ( mode .eq. 2 .or. mode .eq. 3 ) then
366 do i = 1, nobs
367 jac(dimx+i) = 2*x(dimx+i)
368 enddo
369 endif
370!
371! Row 1 to Nobs: The observation definitions:
372!
373 else
374!
375! Mode = 1 or 3: Function value - x-part only
376!
377 if ( mode .eq. 1 .or. mode .eq. 3 ) then
378 s = 0.d0
379 do j = 1, dimx
380 s = s + a(rowno,j)*x(j) + b(rowno,j)*x(j)**2
381 enddo
382 g = s
383 endif
384!
385! Mode = 2 or 3: Derivatives
386!
387 if ( mode .eq. 2 .or. mode .eq. 3 ) then
388 do j = 1, dimx
389 jac(j) = a(rowno,j) + 2.d0*b(rowno,j)*x(j)
390 enddo
391 endif
392 endif
393 lsq_fdeval = 0
394
395end Function lsq_fdeval
396
397
398!> Computes the second derivative of a constraint in a direction
399!!
400!! @include{doc} 2DDir_params.dox
401Integer Function lsq_2ddir( X, DX, D2G, ROWNO, JCNM, NODRV, N, NJ, THREAD, USRMEM )
402#if defined(itl)
403!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_2DDir
404#endif
405 Use lsq5
406 implicit none
407 INTEGER, Intent(IN) :: rowno, n, nj, thread
408 INTEGER, Intent(IN OUT) :: nodrv
409 INTEGER, Dimension(NJ), Intent(IN) :: jcnm
410 real*8, Dimension(N), Intent(IN) :: x
411 real*8, Dimension(N), Intent(IN) :: dx
412 real*8, Dimension(N), Intent(OUT) :: d2g
413 real*8, Intent(IN OUT) :: usrmem(*)
414 integer :: i, j
415
416 if ( rowno .eq. nobs+1 ) then
417!
418! Objective row: sum(i, res(i)**2 )
419!
420 do j = 1, dimx
421 d2g(j) = 0.d0
422 enddo
423 do i = 1, nobs
424 d2g(dimx+i) = 2.d0*dx(dimx+i)
425 enddo
426 else
427!
428! Constraint row with b(i,j)*x(j)**2 as only nonlinear term
429!
430 do j = 1, dimx
431 d2g(j) = 2.d0*b(rowno,j)*dx(j)
432 enddo
433 do i = 1, nobs
434 d2g(dimx+i) = 0.d0
435 enddo
436 endif
437
438 lsq_2ddir = 0
439
440End Function lsq_2ddir
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_2ddir(cntvect, coi_2ddir)
define callback routine for computing the second derivative for a constraint in a direction.
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_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_thread2d(cntvect, thread2d)
number of threads allowed for simultaneous 2DDir calls.
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
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_2ddir(x, dx, d2g, rowno, jcnm, nodrv, n, nj, thread, usrmem)
Computes the second derivative of a constraint in a direction.
#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
#define nj
Definition mp_trans.c:46
real *8, dimension(nobs) obs
real *8, dimension(nobs, dimx) a
real *8, dimension(nobs, dimx) b
subroutine flog(msg, code)
Definition comdeclp.f90:42
real *8 obj
Definition comdeclp.f90:13
subroutine startup
Definition comdeclp.f90:25