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