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