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