CONOPT
Loading...
Searching...
No Matches
mp_leastsq13.f90
Go to the documentation of this file.
1!> @file mp_leastsq13.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!! This implementation has an approximate 2DLagr routine.
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 lsq13
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 lsq13
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 lsq13
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 lsq13
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#ifdef dec_directives_win32
95!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
96!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
97!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Status
98!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Solution
99!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_Message
100!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Std_ErrMsg
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_create( 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_optdir( cntvect, -1 ) ) ! Minimize
127 coi_error = max( coi_error, coidef_objcon( cntvect, nobs + 1 ) ) ! Objective is last constraint
128 coi_error = max( coi_error, coidef_optfile( cntvect, 'mp_leastsq13.opt' ) )
129!
130! Tell CONOPT about the callback routines:
131!
132 coi_error = max( coi_error, coidef_readmatrix( cntvect, lsq_readmatrix ) )
133 coi_error = max( coi_error, coidef_fdeval( cntvect, lsq_fdeval ) )
134 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
135 coi_error = max( coi_error, coidef_solution( cntvect, std_solution ) )
136 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
137 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
138
139#if defined(CONOPT_LICENSE_INT_1) && defined(CONOPT_LICENSE_INT_2) && defined(CONOPT_LICENSE_INT_3) && defined(CONOPT_LICENSE_TEXT)
140 coi_error = max( coi_error, coidef_license( cntvect, conopt_license_int_1, conopt_license_int_2, conopt_license_int_3, conopt_license_text) )
141#endif
142
143 If ( coi_error .ne. 0 ) THEN
144 write(*,*)
145 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
146 write(*,*)
147 call flog( "Skipping Solve due to setup errors", 1 )
148 ENDIF
149!
150! Start CONOPT -- First time using a single thread:
151!
152 time0 = omp_get_wtime()
153 coi_error = coi_solve( cntvect )
154 time1 = omp_get_wtime() - time0
155
156 write(*,*)
157 write(*,*) 'End of Least Square example 13 with 1 thread. Return code=',coi_error
158
159 If ( coi_error /= 0 ) then
160 call flog( "One Thread: Errors encountered during solution", 1 )
161 elseif ( stacalls == 0 .or. solcalls == 0 ) then
162 call flog( "One Thread: Status or Solution routine was not called", 1 )
163 elseif ( .not. ( sstat == 1 .and. mstat == 2 ) ) then
164 call flog( "One Thread: Solver or Model status was not as expected (1,2)", 1 )
165 endif
166!
167! Start CONOPT again -- this time using multiple threads:
168!
169 coi_error = max( coi_error, coidef_threads( cntvect, 0 ) ) ! 0 means use as many threads as you can
170 time0 = omp_get_wtime()
171 coi_error = coi_solve( cntvect )
172 time2 = omp_get_wtime() - time0
173
174 write(*,*)
175 write(*,*) 'Multi thread: End of Least Square example 13. Return code=',coi_error
176
177 If ( coi_error /= 0 ) then
178 call flog( "Multi thread: Errors encountered during solution", 1 )
179 elseif ( stacalls == 0 .or. solcalls == 0 ) then
180 call flog( "Multi thread: Status or Solution routine was not called", 1 )
181 elseif ( .not. ( sstat == 1 .and. mstat == 2 ) ) then
182 call flog( "Multi thread: Solver or Model status was not as expected (1,2)", 1 )
183 endif
184
185 write(*,*)
186 write(*,*) 'End of Least Square example 13. Return code=',coi_error
187
188 if ( coi_free( cntvect ) /= 0 ) call flog( "Error while freeing control vector", 1 )
189
190 write(*,*)
191 write(*,"('Time for single thread',f10.3)") time1
192 write(*,"('Time for multi thread',f10.3)") time2
193 write(*,"('Speedup ',f10.3)") time1/time2
194 write(*,"('Efficiency ',f10.3)") time1/time2/omp_get_max_threads()
195
196 call flog( "Successful Solve", 0 )
197
198End Program leastsquare
199!
200! ============================================================================
201! Define information about the model:
202!
203
204!> Define information about the model
205!!
206!! @include{doc} readMatrix_params.dox
207Integer Function lsq_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
208 colsta, rowno, value, nlflag, n, m, nz, &
209 usrmem )
210#ifdef dec_directives_win32
211!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
212#endif
213 Use lsq13
214 implicit none
215 integer, intent (in) :: n ! number of variables
216 integer, intent (in) :: m ! number of constraints
217 integer, intent (in) :: nz ! number of nonzeros
218 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
219 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
220 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
221 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
222 ! (not defined here)
223 integer, intent (out), dimension(m) :: type ! vector of equation types
224 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
225 ! (not defined here)
226 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
227 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
228 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
229 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
230 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
231 real*8 usrmem(*) ! optional user memory
232
233 Integer :: i, j, k
234!
235! Information about Variables:
236! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
237! Default: the status information in Vsta is not used.
238!
239 do i = 1, dimx
240 curr(i) = -0.8d0
241 enddo
242!
243! Information about Constraints:
244! Default: Rhs = 0
245! Default: the status information in Esta and the function
246! value in FV are not used.
247! Default: Type: There is no default.
248! 0 = Equality,
249! 1 = Greater than or equal,
250! 2 = Less than or equal,
251! 3 = Non binding.
252!
253! Constraints 1 to Nobs:
254! Rhs = Obs(i) and type Equality
255!
256 do i = 1, nobs
257 rhs(i) = obs(i)
258 type(i) = 0
259 enddo
260!
261! Constraint Nobs + 1 (Objective)
262! Rhs = 0 and type Non binding
263!
264 type(nobs+1) = 3
265!
266! Information about the Jacobian. CONOPT expects a columnwise
267! representation in Rowno, Value, Nlflag and Colsta.
268!
269! Colsta = Start of column indices (No Defaults):
270! Rowno = Row indices
271! Value = Value of derivative (by default only linear
272! derivatives are used)
273! Nlflag = 0 for linear and 1 for nonlinear derivative
274! (not needed for completely linear models)
275!
276!
277! Indices
278! x(i) res(j)
279! j: NL L=1
280! obj: NL
281!
282 k = 1
283 do j = 1, dimx
284 colsta(j) = k
285 do i = 1, nobs
286 rowno(k) = i
287 nlflag(k) = 1
288 k = k + 1
289 enddo
290 enddo
291 do i = 1, nobs
292 colsta(dimx+i) = k
293 rowno(k) = i
294 nlflag(k) = 0
295 value(k) = 1.d0
296 k = k + 1
297 rowno(k) = nobs+1
298 nlflag(k) = 1
299 k = k + 1
300 enddo
301 colsta(dimx+nobs+1) = k
303 lsq_readmatrix = 0 ! Return value means OK
304
305end Function lsq_readmatrix
306!
307!==========================================================================
308! Compute nonlinear terms and non-constant Jacobian elements
309!
310
311!> Compute nonlinear terms and non-constant Jacobian elements
312!!
313!! @include{doc} fdeval_params.dox
314Integer Function lsq_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
315 n, nz, thread, usrmem )
316#ifdef dec_directives_win32
317!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
318#endif
319 use lsq13
320 implicit none
321 integer, intent (in) :: n ! number of variables
322 integer, intent (in) :: rowno ! number of the row to be evaluated
323 integer, intent (in) :: nz ! number of nonzeros in this row
324 real*8, intent (in), dimension(n) :: x ! vector of current solution values
325 real*8, intent (in out) :: g ! constraint value
326 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
327 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
328 ! in this row. Ffor information only.
329 integer, intent (in) :: mode ! evaluation mode: 1 = function value
330 ! 2 = derivatives, 3 = both
331 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
332 ! as errcnt is incremented
333 integer, intent (in out) :: errcnt ! error counter to be incremented in case
334 ! of function evaluation errors.
335 integer, intent (in) :: thread
336 real*8 usrmem(*) ! optional user memory
337
338 integer :: i,j
339 real*8 :: s
340!
341! Row Nobs+1: the objective function is nonlinear
342!
343 if ( rowno .eq. nobs+1 ) then
344!
345! Mode = 1 or 3. Function value: G = sum(i, res(i)**2)
346!
347 if ( mode .eq. 1 .or. mode .eq. 3 ) then
348 s = 0.d0
349 do i = 1, nobs
350 s = s + x(dimx+i)**2
351 enddo
352 g = s
353 endif
354!
355! Mode = 2 or 3: Derivative values:
356!
357 if ( mode .eq. 2 .or. mode .eq. 3 ) then
358 do i = 1, nobs
359 jac(dimx+i) = 2*x(dimx+i)
360 enddo
361 endif
362!
363! Row 1 to Nobs: The observation definitions:
364!
365 else
366!
367! Mode = 1 or 3: Function value - x-part only
368!
369 if ( mode .eq. 1 .or. mode .eq. 3 ) then
370 s = 0.d0
371 do j = 1, dimx
372 s = s + a(rowno,j)*x(j) + b(rowno,j)*x(j)**2
373 enddo
374 g = s
375 endif
376!
377! Mode = 2 or 3: Derivatives
378!
379 if ( mode .eq. 2 .or. mode .eq. 3 ) then
380 do j = 1, dimx
381 jac(j) = a(rowno,j) + 2.d0*b(rowno,j)*x(j)
382 enddo
383 endif
384 endif
385 lsq_fdeval = 0
386
387end Function lsq_fdeval
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_optfile(cntvect, optfile)
define callback routine for defining an options file.
Definition conopt.f90:928
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_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
#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