CONOPT
Loading...
Searching...
No Matches
mp_leastsq11.f90
Go to the documentation of this file.
1!> @file mp_leastsq11.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!!
19!! For more information about the individual callbacks, please have a look at the source code.
20
21#if defined(_WIN32) && !defined(_WIN64)
22#define dec_directives_win32
23#endif
24
25module lsq11
26 Integer, parameter :: Nobs = 700
27 Integer, parameter :: DimX = 500
28 real*8, dimension(Nobs,DimX) :: a, b
29 real*8, dimension(Nobs) :: obs
30end module lsq11
32
33REAL FUNCTION rndx( )
34!
35! Defines a pseudo random number between 0 and 1
36!
37 IMPLICIT NONE
38
39 Integer, save :: seed = 12359
40
41 seed = mod(seed*1027+25,1048576)
42 rndx = float(seed)/float(1048576)
43
44END FUNCTION rndx
45
46subroutine defdata
47!
48! Define values for A, B, and Obs
49!
50 Use lsq11
51 IMPLICIT NONE
52
53 Integer :: i, j
54 real*8, Parameter :: xtarg = -1.0
55 real*8, Parameter :: noise = 1.0
56 Real, External :: Rndx
57 real*8 :: o
58
59 do i = 1, nobs
60 o = 0.d0
61 do j = 1, dimx
62 a(i,j) = rndx()
63 b(i,j) = rndx()
64 o = o + a(i,j) * xtarg + b(i,j) * xtarg**2
65 enddo
66 obs(i) = o + noise * rndx()
67 enddo
68
69end subroutine defdata
70!
71! Main program.
72!
73!> Main program. A simple setup and call of CONOPT
74!!
75Program leastsquare
76
78 Use conopt
79 Use lsq11
80 Use omp_lib
81 implicit None
82!
83! Declare the user callback routines as Integer, External:
84!
85 Integer, External :: lsq_readmatrix ! Mandatory Matrix definition routine defined below
86 Integer, External :: lsq_fdeval ! Function and Derivative evaluation routine
87 ! needed a nonlinear model.
88 Integer, External :: std_status ! Standard callback for displaying solution status
89 Integer, External :: std_solution ! Standard callback for displaying solution values
90 Integer, External :: std_message ! Standard callback for managing messages
91 Integer, External :: std_errmsg ! Standard callback for managing error messages
92#ifdef dec_directives_win32
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#endif
100!
101! Control vector
102!
103 INTEGER, Dimension(:), Pointer :: cntvect
104 INTEGER :: coi_error
105!
106 real*8 time0, time1, time2
107
108 call startup
109!
110! Define data
111!
112 Call defdata
113!
114! Create and initialize a Control Vector
115!
116 coi_error = coi_create( cntvect )
117!
118! Tell CONOPT about the size of the model by populating the Control Vector:
119!
120 coi_error = max( coi_error, coidef_numvar( cntvect, nobs + dimx ) )
121 coi_error = max( coi_error, coidef_numcon( cntvect, nobs + 1 ) )
122 coi_error = max( coi_error, coidef_numnz( cntvect, nobs * (dimx+2) ) )
123 coi_error = max( coi_error, coidef_numnlnz( cntvect, nobs * (dimx+1) ) )
124 coi_error = max( coi_error, coidef_optdir( cntvect, -1 ) ) ! Minimize
125 coi_error = max( coi_error, coidef_objcon( cntvect, nobs + 1 ) ) ! Objective is last constraint
126 coi_error = max( coi_error, coidef_optfile( cntvect, 'mp_leastsq11.opt' ) )
127!
128! Tell CONOPT about the callback routines:
129!
130 coi_error = max( coi_error, coidef_readmatrix( cntvect, lsq_readmatrix ) )
131 coi_error = max( coi_error, coidef_fdeval( cntvect, lsq_fdeval ) )
132 coi_error = max( coi_error, coidef_status( cntvect, std_status ) )
133 coi_error = max( coi_error, coidef_solution( cntvect, std_solution ) )
134 coi_error = max( coi_error, coidef_message( cntvect, std_message ) )
135 coi_error = max( coi_error, coidef_errmsg( cntvect, std_errmsg ) )
136
137#if defined(CONOPT_LICENSE_INT_1) && defined(CONOPT_LICENSE_INT_2) && defined(CONOPT_LICENSE_INT_3) && defined(CONOPT_LICENSE_TEXT)
138 coi_error = max( coi_error, coidef_license( cntvect, conopt_license_int_1, conopt_license_int_2, conopt_license_int_3, conopt_license_text) )
139#endif
140
141 If ( coi_error .ne. 0 ) THEN
142 write(*,*)
143 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
144 write(*,*)
145 call flog( "Skipping Solve due to setup errors", 1 )
146 ENDIF
147!
148! Start CONOPT -- First time using a single thread:
149!
150 time0 = omp_get_wtime()
151 coi_error = coi_solve( cntvect )
152 time1 = omp_get_wtime() - time0
153
154 write(*,*)
155 write(*,*) 'End of Least Square example 11 with 1 thread. Return code=',coi_error
156
157 If ( coi_error /= 0 ) then
158 call flog( "One Thread: Errors encountered during solution", 1 )
159 elseif ( stacalls == 0 .or. solcalls == 0 ) then
160 call flog( "One Thread: Status or Solution routine was not called", 1 )
161 elseif ( .not. ( sstat == 1 .and. mstat == 2 ) ) then
162 call flog( "One Thread: Solver or Model status was not as expected (1,2)", 1 )
163 endif
164!
165! Start CONOPT again -- this time using multiple threads:
166!
167 coi_error = max( coi_error, coidef_threads( cntvect, 0 ) ) ! 0 means use as many threads as you can
168 time0 = omp_get_wtime()
169 coi_error = coi_solve( cntvect )
170 time2 = omp_get_wtime() - time0
171
172 write(*,*)
173 write(*,*) 'Multi thread: End of Least Square example 11. Return code=',coi_error
174
175 If ( coi_error /= 0 ) then
176 call flog( "Multi thread: Errors encountered during solution", 1 )
177 elseif ( stacalls == 0 .or. solcalls == 0 ) then
178 call flog( "Multi thread: Status or Solution routine was not called", 1 )
179 elseif ( .not. ( sstat == 1 .and. mstat == 2 ) ) then
180 call flog( "Multi thread: Solver or Model status was not as expected (1,2)", 1 )
181 endif
182
183 write(*,*)
184 write(*,*) 'End of Least Square example 11. Return code=',coi_error
185
186 if ( coi_free( cntvect ) /= 0 ) call flog( "Error while freeing control vector", 1 )
187
188 write(*,*)
189 write(*,"('Time for single thread',f10.3)") time1
190 write(*,"('Time for multi thread',f10.3)") time2
191 write(*,"('Speedup ',f10.3)") time1/time2
192 write(*,"('Efficiency ',f10.3)") time1/time2/omp_get_max_threads()
193
194 call flog( "Successful Solve", 0 )
195
196End Program leastsquare
197!
198! ============================================================================
199! Define information about the model:
200!
201
202!> Define information about the model
203!!
204!! @include{doc} readMatrix_params.dox
205Integer Function lsq_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
206 colsta, rowno, value, nlflag, n, m, nz, &
207 usrmem )
208#ifdef dec_directives_win32
209!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_ReadMatrix
210#endif
211 Use lsq11
212 implicit none
213 integer, intent (in) :: n ! number of variables
214 integer, intent (in) :: m ! number of constraints
215 integer, intent (in) :: nz ! number of nonzeros
216 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
217 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
218 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
219 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
220 ! (not defined here)
221 integer, intent (out), dimension(m) :: type ! vector of equation types
222 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
223 ! (not defined here)
224 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
225 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
226 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
227 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
228 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
229 real*8 usrmem(*) ! optional user memory
230
231 Integer :: i, j, k
232!
233! Information about Variables:
234! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
235! Default: the status information in Vsta is not used.
236!
237 do i = 1, dimx
238 curr(i) = -0.8d0
239 enddo
240 do i = 1, nobs ! prevent post triangle from being identified.
241 lower(dimx+i) = -1000.0d0
242 upper(dimx+i) = +1000.0d0
243 enddo
244!
245! Information about Constraints:
246! Default: Rhs = 0
247! Default: the status information in Esta and the function
248! value in FV are not used.
249! Default: Type: There is no default.
250! 0 = Equality,
251! 1 = Greater than or equal,
252! 2 = Less than or equal,
253! 3 = Non binding.
254!
255! Constraints 1 to Nobs:
256! Rhs = Obs(i) and type Equality
257!
258 do i = 1, nobs
259 rhs(i) = obs(i)
260 type(i) = 0
261 enddo
262!
263! Constraint Nobs + 1 (Objective)
264! Rhs = 0 and type Non binding
265!
266 type(nobs+1) = 3
267!
268! Information about the Jacobian. CONOPT expects a columnwise
269! representation in Rowno, Value, Nlflag and Colsta.
270!
271! Colsta = Start of column indices (No Defaults):
272! Rowno = Row indices
273! Value = Value of derivative (by default only linear
274! derivatives are used)
275! Nlflag = 0 for linear and 1 for nonlinear derivative
276! (not needed for completely linear models)
277!
278!
279! Indices
280! x(i) res(j)
281! j: NL L=1
282! obj: NL
283!
284 k = 1
285 do j = 1, dimx
286 colsta(j) = k
287 do i = 1, nobs
288 rowno(k) = i
289 nlflag(k) = 1
290 k = k + 1
291 enddo
292 enddo
293 do i = 1, nobs
294 colsta(dimx+i) = k
295 rowno(k) = i
296 nlflag(k) = 0
297 value(k) = 1.d0
298 k = k + 1
299 rowno(k) = nobs+1
300 nlflag(k) = 1
301 k = k + 1
302 enddo
303 colsta(dimx+nobs+1) = k
305 lsq_readmatrix = 0 ! Return value means OK
306
307end Function lsq_readmatrix
308!
309!==========================================================================
310! Compute nonlinear terms and non-constant Jacobian elements
311!
312
313!> Compute nonlinear terms and non-constant Jacobian elements
314!!
315!! @include{doc} fdeval_params.dox
316Integer Function lsq_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
317 n, nz, thread, usrmem )
318#ifdef dec_directives_win32
319!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: Lsq_FDEval
320#endif
321 use lsq11
322 implicit none
323 integer, intent (in) :: n ! number of variables
324 integer, intent (in) :: rowno ! number of the row to be evaluated
325 integer, intent (in) :: nz ! number of nonzeros in this row
326 real*8, intent (in), dimension(n) :: x ! vector of current solution values
327 real*8, intent (in out) :: g ! constraint value
328 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
329 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
330 ! in this row. Ffor information only.
331 integer, intent (in) :: mode ! evaluation mode: 1 = function value
332 ! 2 = derivatives, 3 = both
333 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
334 ! as errcnt is incremented
335 integer, intent (in out) :: errcnt ! error counter to be incremented in case
336 ! of function evaluation errors.
337 integer, intent (in) :: thread
338 real*8 usrmem(*) ! optional user memory
339
340 integer :: i,j
341 real*8 :: s
342!
343! Row Nobs+1: the objective function is nonlinear
344!
345 if ( rowno .eq. nobs+1 ) then
346!
347! Mode = 1 or 3. Function value: G = sum(i, res(i)**2)
348!
349 if ( mode .eq. 1 .or. mode .eq. 3 ) then
350 s = 0.d0
351 do i = 1, nobs
352 s = s + x(dimx+i)**2
353 enddo
354 g = s
355 endif
356!
357! Mode = 2 or 3: Derivative values:
358!
359 if ( mode .eq. 2 .or. mode .eq. 3 ) then
360 do i = 1, nobs
361 jac(dimx+i) = 2*x(dimx+i)
362 enddo
363 endif
364!
365! Row 1 to Nobs: The observation definitions:
366!
367 else
368!
369! Mode = 1 or 3: Function value - x-part only
370!
371 if ( mode .eq. 1 .or. mode .eq. 3 ) then
372 s = 0.d0
373 do j = 1, dimx
374 s = s + a(rowno,j)*x(j) + b(rowno,j)*x(j)**2
375 enddo
376 g = s
377 endif
378!
379! Mode = 2 or 3: Derivatives
380!
381 if ( mode .eq. 2 .or. mode .eq. 3 ) then
382 do j = 1, dimx
383 jac(j) = a(rowno,j) + 2.d0*b(rowno,j)*x(j)
384 enddo
385 endif
386 endif
387 lsq_fdeval = 0
388
389end 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, dimx) b
real *8, dimension(nobs) obs
real *8, dimension(nobs, dimx) a
subroutine flog(msg, code)
Definition comdeclp.f90:48
subroutine startup
Definition comdeclp.f90:31