CONOPT
Loading...
Searching...
No Matches
mp_qpbandb.f90
Go to the documentation of this file.
1!> @file mp_qpbandb.f90
2!! @ingroup FORTOPENMP_EXAMPLES
3!!
4!!
5!! QPBandB: Multi Threaded implementation of QP4 MINLP model
6!! The implementation shows how to implement a parallel Branch & Bound
7!! procedure that still is deterministic independent of the order in
8!! which the sub-models are solved (within the framework defined here
9!! with Work-pools).
10!!
11!! When testing that the sequential and parallel algorithms uses the same
12!! number of nodes we have a test that the solves are identical even
13!! when run with different threads. This is some kind of test that
14!! CONOPT Solves are reproducable and that there are no un-initialized
15!! values used in the algorithm.
16!!
17!!
18!! For more information about the individual callbacks, please have a look at the source code.
19
20Module qp_data
21 Implicit None
22 Integer, Parameter :: nn = 32 ! Dimension of Q
23 Integer, Parameter :: nvar = 2*nn ! Number of variables in model
24 Integer, Parameter :: nbin = nn ! Number of binary variables in model
25
26 Integer, Parameter :: nobs = 60 ! Observations used to estimate Q
27 Real(8), dimension(NN,NN) :: q
28 Real(8), dimension(nn) :: mean
29 Real(8) :: target
30 Integer, Parameter :: maxbin = 15
31 Integer, Dimension(Nbin) :: binmap ! Mapping from binary to index in all variables
32End Module qp_data
33
34Module bb_data
35!
36! Data used to manage the Branch and Bound Process
37!
38 Use qp_data
39 Use omp_lib
40 Implicit None
41!
42! Branch and bound node data
43!
44 Type bbnode
45 Integer :: number ! Sequense number of the node.
46 real(8), dimension(NVar) :: values ! Solution value or initial values
47 Integer, dimension(Nbin) :: fixed ! Binary: 0=fixed zero, 1=fixed 1, -1 not fixed.
48 Type (bbnode), Pointer :: prvnode, nxtnode ! Pointer forward and backward in
49 ! double-linked list of active nodes
50 real(8) :: object ! Objective value (if solved) or lower bound if
51 ! not yet solved
52 Integer :: node_stat ! Node Status:
53 ! 0 - created but not solved
54 ! 1 - solved
55 Integer :: solstat ! Solution status. Initialized to -1 and
56 ! redefined when the model has been solved
57 Integer :: modstat ! Model status. Initialized to -1 and
58 ! redefined when the model has been solved
59 Integer :: depth ! Depth of the node in the B&B tree
60 Logical :: isinteger ! Is true iff the solution is integer
61 Integer :: sthread ! The thread being used when solving
62 End type bbnode
63
65 Type (bbnode), Pointer :: bnode
66 End Type bbnodeptr
67
68 Real(8), Parameter :: bintol = 1.d-7
69 Integer :: totnodes ! Total number of nodes in use
70 Integer :: freenodes ! Number of active nodes
71 Integer :: actnodes ! Number of free nodes in a free list
72 Integer :: numsolve ! Counter for finished Solves.
73!
74! Data for best solution
75!
76 Logical :: intfound ! Global flag indicating whether we already have an integer
77 ! feasible solution.
78 Integer :: numintfound ! Number of integer solutions found
79 Real(8) :: bestint ! Objective value for the best integer solution found so far.
80 Real(8), dimension(NVar) :: solution ! The solution corresponding to BestInt
81!
82! Data for managing the Branch and Bound Tree
83!
84 Type (bbnode), Pointer :: actnode ! Pointer to first node in a list of active untested noded
85 Type (bbnode), Pointer :: freenode ! Pointer to first node in a freelist
86!
87! Data for managing a WorkPool with the nodes selected for the next
88! set of parallel Solves
89!
90 Integer, Parameter :: maxmaxpool = 64 ! Allocation size for WorkPool
91 Integer :: maxpool ! Largest size of WorkPool selected by user
92 Logical :: maxdyn ! If .true. then MaxPool above is dynamic
93 Integer :: workpool ! Size of the actual workpool
94 Type (bbnodeptr) :: work_node(maxmaxpool) ! Pointers to the selected nodes in the WorkPool
95!
96 integer :: maxthread ! Maximum number of threads
97 Logical :: parallel ! If true then run the COI_Solve loop over the nodes
98 ! in the WorkPool in parallel
100 Integer, Dimension(:), Pointer :: cntvect
101 End Type controlv
102 Type (controlv), Dimension(:), Pointer :: cntvect1 ! CONOPT Control vectors, one for each thread
103 Real(8) :: stime ! Accumulated time in Solve Loop
104
105 Contains
106!
107! Routines used to manage the Branch and Bound Tree
108!
109Subroutine remove_activenode(Node)
110!
111! Remove a node from the list of active nodes.
112!
113 Implicit None
114 Type (BBNode), Pointer :: Node
115
116 if ( associated(node%prvnode) ) then
117 node%prvnode%nxtnode => node%nxtnode
118 else
119 actnode => node%nxtnode
120 endif
121 if ( associated(node%nxtnode) ) then
122 node%nxtnode%prvnode => node%prvnode
123 endif
124 actnodes = actnodes - 1
125
126End Subroutine remove_activenode
127
128Subroutine return_freenode(Node)
129!
130! Link the node into the free list. The free list does not use prvnode
131!
132 Implicit None
133 Type (BBNode), Pointer :: Node
134
135 node%nxtnode => freenode
136 freenode => node
137 freenodes = freenodes + 1
138
139End Subroutine return_freenode
140
141Subroutine get_freenode(Node)
142!
143! Get a free node to be used in the list of active nodes and initialize
144! the status.
145!
146 Implicit None
147 Type (BBNode), Pointer :: Node
148
149 if ( associated(freenode) ) then
150 node => freenode
151 freenode => freenode%Nxtnode
152 freenodes = freenodes - 1
153 else
154 Allocate(node)
155 Nullify(node%NxtNode);
156 Nullify(node%PrvNode);
157 totnodes = totnodes + 1
158 node%Number = totnodes
159 endif
160 node%Node_Stat = 0
161 node%Solstat =-1
162 node%Modstat =-1
163 node%Sthread =-1
164
165End Subroutine get_freenode
166
168!
169! Deallocates the nodes in the free list between solves so the next
170! solve starts without any nodes.
171!
172 Implicit None
173 Type (BBNode), Pointer :: Node, Next
174
175 node => freenode
176 do while ( Associated(node) )
177 next => node%Nxtnode
178 Deallocate(node);
179 node => next
180 enddo
181 Nullify(freenode)
182
183End Subroutine deallocate_nodes
184
185Subroutine add_activenode(Node)
186!
187! Add a node to the front of the list of Active nodes
188!
189 Implicit None
190 Type (BBNode), pointer :: Node
191
192 actnodes = actnodes + 1
193 if ( associated(actnode) ) then
194 node%NxtNode => actnode
195 actnode%PrvNode => node;
196 else
197 Nullify(node%NxtNode);
198 endif
199 Nullify(node%PrvNode);
200 actnode => node
201
202End Subroutine add_activenode
203
204Subroutine add2workpool( Node )
205!
206! Add a node to the Workpool list
207!
208 Implicit None
209 Type (BBNode), pointer :: Node
210
211 workpool = workpool + 1
212! write(*,*) 'Work_node(',Workpool,') has number',node%number
213 work_node(workpool)%Bnode => node
214
215end Subroutine add2workpool
216
217Subroutine list_node(node)
218
219 Implicit None
220 Type (BBNode), pointer :: Node
221
222 If ( Associated(node) ) Then
223 write(*,*) 'List_Node: Number=',node%Number,' Node_Stat=',node%node_stat
224 Else
225 write(*,*) 'List_Node: Node not associated.'
226 Endif
227End Subroutine list_node
228
229Subroutine list_nodes
230
231 Implicit None
232 Type (BBNode), Pointer :: Node
233 Integer :: N
234
235 n = 0
236 node => actnode
237 do while ( Associated(node) )
238 write(*,"('List_nodes: Node number=',i4,' Node_Stat=',I2,' Object=',1pe15.8)") node%number, node%node_stat, node%object
239 node => node%nxtnode
240 n = n + 1
241 enddo
242 write(*,*) 'List_nodes: Total=',n
243
244End Subroutine list_nodes
245
246Subroutine select_var(Node, IB)
247!
248! Select variable to branch on, counted in the space of binary variables.
249! Currently there is only one rule, namely to select based on max
250! infeasibility. Other Variable selection rules could be implemented
251! in this routine.
252!
253 Implicit None
254 Type (BBNode), Pointer :: Node
255 integer :: IB, i, j
256 Real(8) :: Dist, Maxdist
257
258 maxdist = 0.d0
259 ib = 0
260 do i = 1, nn
261 j = binmap(i)
262 dist = min( node%values(j), 1.d0-node%values(j) )
263 if ( dist > maxdist ) then
264 maxdist = dist
265 ib = i
266 endif
267 enddo
268
269End Subroutine select_var
270
272!
273! Create a workpool of up to MaxPool nodes to be processed in
274! parallel.
275! The nodes are removed from the active list.
276!
277 Implicit None
278 Type (BBNode), Pointer :: Node, Curr, Next
279 Real(8) :: Cutoff
280 integer :: IB ! Branching variable
281
282 Nullify(node)
283 workpool = 0
284 write(*,*) 'Create_WorkPool'
285!
286! If the max size of the workpool is dynamic then use the formula
287! MaxThread*min(8,1+ActNodes/(MaxThread*5)). With few nodes
288! (<MaxThread*5) we will select MaxThread nodes. With more nodes we
289! will select gradually more up to a maximum of MaxThread*8.
290!
291 if ( maxdyn ) maxpool = maxthread*min(8,1+actnodes/(maxthread*5))
292 do while ( workpool < maxpool )
293 cutoff = 1.d30
294 curr => actnode
295!
296! Run through the linked list of active nodes and find the best.
297! We use the objective but other criteria could be used
298!
299 Nullify(node)
300 do while (associated(curr))
301 next => curr%NxtNode
302 if ( curr%Object < cutoff ) then
303 node => curr
304 cutoff = curr%Object
305 endif
306 curr => next
307 enddo
308 if ( .not. associated(node) ) then
309 exit
310 endif
311 if ( node%Node_Stat == 0 ) then
312!
313! The node has been created but not solved.
314! Take it immediately: remove it from the active list and add it
315! to the work_node vector
316!
317 write(*,"('Unsolved node',i5,' selected. Obj=',e24.17)") node%Number,node%Object
318 Call remove_activenode(node)
319 call add2workpool(node)
320 else
321!
322! The node has been solved. Split it in two and take one or both
323! depending on the space left in the workpool. The total number
324! of nodes will grow by one so a dynamic MaxPool must be adjusted
325! before the space test below
326!
327 if ( maxdyn ) maxpool = maxthread*min(8,1+(actnodes+1)/(maxthread*5))
328!
329! Get a new node and, select a branching variable, and copy the
330! old node into the new node.
331!
332 call get_freenode(next)
333 call select_var( node, ib )
334 next%Values = node%Values
335 next%Fixed = node%Fixed
336 next%Object = node%Object
337!
338! Make the two node unsolved and fix the branching variable down and up
339! respectively.
340!
341 node%Node_Stat = 0
342 node%Solstat = -1
343 node%Modstat = -1
344 node%Depth = node%Depth + 1
345 node%Fixed(ib) = 0 ! Branch down
346
347 next%Node_Stat = 0
348 next%Solstat = -1
349 next%Modstat = -1
350 next%Depth = node%Depth
351 next%Fixed(ib) = 1 ! Branch up
352!
353! Add the first node to the workpool
354!
355 Call remove_activenode(node)
356 call add2workpool(node)
357 write(*,"('New-Down node',i5,' selected. Obj=',e24.17)") node%Number,node%Object
358!
359! Add the second node to the workpool or to the pool of active
360! nodes, depending on space in the workpool
361!
362 if ( workpool < maxpool ) then
363 call add2workpool(next)
364 write(*,"('New-Up Next',i5,' selected. Obj=',e24.17)") next%Number,next%Object
365 else
366 call add_activenode(next)
367 endif
368 endif
369 enddo
370
371end subroutine create_workpool
372
374!
375! Solve the models in the WorkPool and integrate the solutions into
376! the active list of nodes. Update the best solution if one is found.
377!
378 Use proginfop
379 Use coidef
380 implicit none
381
382 integer nw ! index for active node in the WorkPool
383 integer Thread
384 Integer, Dimension(:), pointer :: CntVect ! CONOPT Control vector for current node
385 integer COI_Error
386 Type (BBNode), Pointer :: Node, Curr, Next
387 Logical :: NewInt
388 Real(8) Time0
389!
390! COIDEF_UsrMem and COIDEF_WorkMem are not in the Coidef module and must be
391! explicitly declared if used.
392!
393 INTEGER, External :: COIDEF_UsrMem
394#if defined(itl)
395!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: COIDEF_UsrMem
396!DEC$ ATTRIBUTES DECORATE, ALIAS: 'COIDEF_USRMEM'::COIDEF_UsrMem
397#endif
398
399! write(*,*) 'Calling Solve_WorkPool with',WorkPool,' nodes.'
400 time0 = omp_get_wtime()
401 if ( workpool .le. 1 .OR. .NOT. parallel ) Then
402! This is a serial version
403 thread = 1
404 do nw = 1, workpool
405 node => work_node(nw)%Bnode
406 if ( debug ) then
407 write(16,*) 'Starting to solve node',node%number,' using thread',thread
408 endif
409 cntvect => cntvect1(thread)%CntVect
410 coi_error = coidef_usrmem( cntvect, node ) ! Pass the current node through Usrmem
411 If ( coi_error .ne. 0 ) THEN
412 write(*,*)
413 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
414 write(*,*)
415 stop
416 ENDIF
417 node%Sthread = thread
418 coi_error = coi_solve( cntvect )
419 If ( coi_error .ne. 0 ) THEN
420 write(*,*)
421 write(*,*) '**** Fatal Error while solving.'
422 write(*,*)
423 stop
424 ENDIF
425 enddo
426 Else
427! And this is the corresponding Parallel version of the SOLVE loop.
428!$OMP PARALLEL DEFAULT(SHARED) Private(COI_Error,Thread, Node)
429!$OMP DO PRIVATE(Thread, COI_Error, Node, CntVect) SCHEDULE(Guided)
430 do nw = 1, workpool
431 node => work_node(nw)%Bnode
432 thread = 1+omp_get_thread_num()
433 if ( debug ) then
434 write(16,*) 'Starting to solve node',node%number,' using thread',thread
435 endif
436 cntvect => cntvect1(thread)%CntVect
437 coi_error = coidef_usrmem( cntvect, node ) ! Pass the current node through Usrmem
438 If ( coi_error .ne. 0 ) THEN
439 write(*,*)
440 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
441 write(*,*)
442 call flog( "Errors encountered during setup for Solve.", 1 )
443 ENDIF
444 node%Sthread = thread
445 coi_error = coi_solve( cntvect )
446 If ( coi_error .ne. 0 ) THEN
447 write(*,*)
448 write(*,*) '**** Fatal Error while solving.'
449 write(*,*)
450 call flog( "Errors encountered during Solve.", 1 )
451 ENDIF
452 enddo
453!$OMP END PARALLEL
454 Endif
455 stime = stime + omp_get_wtime() - time0
456 numsolve = numsolve + workpool
457!
458! Process the solutions. The loop is not parallelized so the order here
459! is fixed, independent of the order in which the nodes finished
460! solving.
461!
462 newint = .false. ! No new integer solution in this pass
463 do nw = 1, workpool
464 node => work_node(nw)%Bnode
465 write(*,'("Node=",I5," Solstat=",I3," Modstat=",I3," Object=",F24.17)') &
466 node%number, node%Solstat, node%modstat, node%Object
467 if ( node%Solstat .ne. 1 .or. node%Modstat .gt. 2 ) then
468!
469! Node was not solved properly
470!
471 call flog( "Solve did not return an normal solver status.", 1 )
472 Else if ( node%Modstat .gt. 2 ) then
473!
474! Node was not solved to local optimality. Fathom.
475!
476 write(*,*) 'Not optimal -- fathom'
477 call return_freenode(node)
478 else if ( intfound ) then
479 if ( node%object >= bestint ) then
480!
481! A better solution has been found. Fathom.
482!
483 write(*,*) 'Worse than optimal -- fathom'
484 call return_freenode(node)
485 else
486 if ( node%IsInteger ) then
487 write(*,*) 'Integer Solution Stored'
488 numintfound = numintfound + 1
489 newint = .true.
490 bestint = node%object
491 solution = node%Values
492 call return_freenode(node)
493 Else
494! write(*,*) 'Non-Integer. Kept'
495 call add_activenode(node)
496 Endif
497 endif
498 else
499!
500! An Integer solution has not yet been found. We cannot fathom it.
501!
502 if ( node%IsInteger ) then
503 write(*,*) 'Integer Solution Stored'
504 numintfound = numintfound + 1
505 intfound = .true.
506 newint = .true.
507 bestint = node%object
508 solution = node%Values
509 call return_freenode(node)
510 Else
511! write(*,*) 'Non-Integer. Kept'
512 call add_activenode(node)
513 Endif
514 endif
515 enddo
516!
517! If we did find a new integer solution then scan the list of active
518! nodes to check if some of them should be fathomed.
519!
520 if ( newint ) then
521 curr => actnode
522 do while (associated(curr))
523 next => curr%NxtNode
524 if ( curr%Object > bestint ) then
525 write(*,"('Node',I5,' with Object=',F18.12,' fathomed by new objective')") curr%Number, curr%Object
526 call remove_activenode(curr)
527 call return_freenode(curr)
528 endif
529 curr => next
530 enddo
531 endif
532
533end subroutine solve_workpool
534
535End Module bb_data
536
537REAL FUNCTION rndx( )
538!
539! Defines a pseudo random number between 0 and 1
540!
541 IMPLICIT NONE
542
543 Integer, save :: seed = 12359
544
545 seed = mod(seed*1027+25,1048576)
546 rndx = float(seed)/float(1048576)
547
548END FUNCTION rndx
549
550subroutine defdata
551!
552! Define values for A, B, and Obs
553!
554 Use qp_data
555 IMPLICIT NONE
556
557 Integer :: i, j, k
558 Real, External :: Rndx
559 Real(8) :: Sum
560 Real(8), dimension(NN,Nobs) :: Xobs
561
562 Target = 10.0d0
563 do i = 1, nn
564 mean(i) = Target + (rndx()*4.d0-2.d0)
565 enddo
566!
567! Create a positive definite matrix computed as the sum of some
568! outer products
569!
570 do i = 1, nn
571 do j = 1, nn
572 q(i,j) = 0.d0
573 enddo
574 enddo
575 do i = 1, nn
576 sum = 0.d0
577 do k = 1, nobs
578 xobs(i,k) = 100.d0*rndx()
579 sum = sum + xobs(i,k)
580 enddo
581 sum = sum / nobs
582 do k = 1, nobs
583 xobs(i,k) = xobs(i,k) - sum
584 enddo
585 enddo
586 do i = 1, nn
587 do j = 1, i
588 sum = 0.d0
589 do k = 1, nobs
590 sum = sum + xobs(i,k) * xobs(j,k)
591 enddo
592 q(i,j) = sum / (nobs-1)
593 q(j,i) = q(i,j)
594 enddo
595 enddo
596!
597! Create BinMap
598!
599 do i = 1, nn
600 binmap(i) = nn+i
601 enddo
602
603end subroutine defdata
604!
605! Main program. A simple setup and call of CONOPT for a QP model
606!
607!> Main program. A simple setup and call of CONOPT
608!!
610
611 Use proginfop
612 Use coidef
613 use bb_data
614 use qp_data
615 Use omp_lib
616 Implicit None
617!
618! Data needed for Branch and bound
619!
620 Integer :: thread
621!
622! Declare the user callback routines as Integer, External:
623!
624 Integer, External :: qp_readmatrix ! Mandatory Matrix definition routine defined below
625 Integer, External :: qp_fdeval ! Function and Derivative evaluation routine
626 ! needed a nonlinear model.
627 Integer, External :: qp_2dlagrstr ! 2nd derivative routine for structure
628 Integer, External :: qp_2dlagrval ! 2nd derivative routine for values
629 Integer, External :: qp_2ddir ! Directional 2nd derivative routine
630 Integer, External :: bb_status ! callback for displaying solution status
631 Integer, External :: bb_solution ! callback for displaying solution values
632 Integer, External :: bb_message ! callback for managing messages
633 Integer, External :: bb_errmsg ! callback for managing error messages
634 Integer, External :: bb_option ! Define options using a callback. Files cannot be used
635 ! for multithreading since they can be shared and not
636 ! read in the proper order.
637#if defined(itl)
638!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_ReadMatrix
639!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_FDEval
640!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_2DLagrStr
641!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_2DLagrVal
642!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_2DDir
643!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: BB_Status
644!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: BB_Solution
645!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: BB_Message
646!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: BB_ErrMsg
647!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: BB_Option
648#endif
649!
650! Control vector
651!
652 INTEGER :: coi_error
653 Integer, Dimension(:), Pointer :: cntvect ! CONOPT Control vector for current node
654
655 Real(8) time0, time1, stime1
656 Real(8) time8s, time16s, time32s, timeds
657 Real(8) time8p, time16p, time32p, timedp
658 Real(8) stime8s, stime16s, stime32s, stimeds
659 Real(8) stime8p, stime16p, stime32p, stimedp
660 Integer numsolve1
661 Integer numsolve8s, numsolve16s, numsolve32s, numsolveds
662 Integer numsolve8p, numsolve16p, numsolve32p, numsolvedp
663
664 call startup
665 maxthread = omp_get_max_threads()
666 write(11,*) 'Initial MaxThread=',maxthread
667 write(*,*) 'Initial MaxThread=',maxthread
668 If ( maxthread .lt. 2 ) then
669 maxthread = 2
670 call omp_set_num_threads(maxthread)
671 write(11,*) 'Revised MaxThread=',maxthread
672 write(*,*) 'Revised MaxThread=',maxthread
673 Endif
674! Debug = .true. ! Turn on if debugging output is needed
675 If ( debug ) Then
676 Open(12,file='threadA1.txt',status='unknown')
677 Open(16,file='nodesA.txt',status='unknown')
678 Endif
679!
680 call defdata
681!
682! Create a vector of Control Vectors in CntVect1 and allocate and
683! initialize the individual vectors for each thread.
684!
685 Allocate( cntvect1( maxthread ) )
686 coi_error = 0
687 do thread = 1, maxthread
688 coi_error = coi_createfort( cntvect )
689 cntvect1(thread)%CntVect => cntvect
690!
691! Tell CONOPT about the size of the model by populating the Control Vector:
692! We will create a QP model with one constraint, sum(i, x(i) ) = 1
693! and NN variables. NN and the Q matrix are defined in Module QP_Data.
694!
695 coi_error = max( coi_error, coidef_numvar( cntvect, 2*nn ) ) ! NN variables
696 coi_error = max( coi_error, coidef_numcon( cntvect, 4+nn ) ) ! 3+NN constraint + 1 objective
697 coi_error = max( coi_error, coidef_numnz( cntvect, 6*nn ) ) ! 6*NN nonzeros in the Jacobian
698 coi_error = max( coi_error, coidef_numnlnz( cntvect, nn ) ) ! NN of which are nonlinear
699 coi_error = max( coi_error, coidef_numhess( cntvect, nn*(nn+1)/2 ) ) ! Hessian is dense in the NN nonlinear variables
700 coi_error = max( coi_error, coidef_optdir( cntvect, -1 ) ) ! Minimize
701 coi_error = max( coi_error, coidef_objcon( cntvect, 1 ) ) ! Objective is constraint 1
702!
703! Tell CONOPT about the callback routines:
704!
705 coi_error = max( coi_error, coidef_readmatrix( cntvect, qp_readmatrix ) )
706 coi_error = max( coi_error, coidef_fdeval( cntvect, qp_fdeval ) )
707 coi_error = max( coi_error, coidef_2dlagrstr( cntvect, qp_2dlagrstr ) )
708 coi_error = max( coi_error, coidef_2dlagrval( cntvect, qp_2dlagrval ) )
709 coi_error = max( coi_error, coidef_2ddir( cntvect, qp_2ddir ) )
710 coi_error = max( coi_error, coidef_status( cntvect, bb_status ) )
711 coi_error = max( coi_error, coidef_solution( cntvect, bb_solution ) )
712 coi_error = max( coi_error, coidef_message( cntvect, bb_message ) )
713 coi_error = max( coi_error, coidef_errmsg( cntvect, bb_errmsg ) )
714 If ( debug ) &
715 coi_error = max( coi_error, coidef_option( cntvect, bb_option ) )
716
717#if defined(LICENSE_INT_1) && defined(LICENSE_INT_2) && defined(LICENSE_INT_3) && defined(LICENSE_TEXT)
718 coi_error = max( coi_error, coidef_license( cntvect, license_int_1, license_int_2, license_int_3, license_text) )
719#endif
720 enddo
721 If ( coi_error .ne. 0 ) THEN
722 write(*,*)
723 write(*,*) '**** Fatal Error while loading CONOPT Callback routines.'
724 write(*,*)
725 call flog( "Skipping First Solve due to setup errors", 1 )
726 ENDIF
727 write(11,*) 'Finished Setup of model'
728 write(*,*) 'Finished Setup of model'
729!
730! Solve the B and B process with MaxPool = 1
731!
732 parallel = .false.
733 maxdyn = .false.
734 maxpool = 1
735 time0 = omp_get_wtime()
736 Call bandb_solve
737 time1 = omp_get_wtime() - time0
738 stime1 = stime
739 numsolve1 = numsolve
740 write(11,*) 'Finished sequential with Maxpool=',maxpool
741 write(*,*) 'Finished sequential with Maxpool=',maxpool
742 If ( debug ) Then
743 close(12); close(16)
744 Open(12,file='threadB1.txt',status='unknown')
745 Open(16,file='nodesB.txt',status='unknown')
746 Endif
747!
748! Redo the B and B process with MaxPool = 8
749!
750 maxpool = 8
751 time0 = omp_get_wtime()
752 Call bandb_solve
753 time8s = omp_get_wtime() - time0
754 stime8s = stime
755 numsolve8s = numsolve
756 write(11,*) 'Finished sequential with Maxpool=',maxpool
757 write(*,*) 'Finished sequential with Maxpool=',maxpool
758!
759! Redo the B and B process with MaxPool = 16
760!
761 maxpool = 16
762 time0 = omp_get_wtime()
763 Call bandb_solve
764 time16s = omp_get_wtime() - time0
765 stime16s = stime
766 numsolve16s = numsolve
767 write(11,*) 'Finished sequential with Maxpool=',maxpool
768 write(*,*) 'Finished sequential with Maxpool=',maxpool
769!
770! Redo the B and B process with MaxPool = 32
771!
772 maxpool = 32
773 time0 = omp_get_wtime()
774 Call bandb_solve
775 time32s = omp_get_wtime() - time0
776 stime32s = stime
777 numsolve32s = numsolve
778 write(11,*) 'Finished sequential with Maxpool=',maxpool
779 write(*,*) 'Finished sequential with Maxpool=',maxpool
780!
781! Redo the B and B process with Dynamic MaxPool
782!
783 maxdyn = .true.
784 time0 = omp_get_wtime()
785 Call bandb_solve
786 timeds = omp_get_wtime() - time0
787 stimeds = stime
788 numsolveds = numsolve
789 write(11,*) 'Finished sequential with Maxpool=',maxpool
790 write(*,*) 'Finished sequential with Maxpool=',maxpool
791 If ( debug ) Then
792 close(12); close(16)
793 Open(12,file='threadC1.txt',status='unknown')
794 Open(13,file='threadC2.txt',status='unknown')
795 Open(14,file='threadC3.txt',status='unknown')
796 Open(15,file='threadC4.txt',status='unknown')
797 Open(16,file='nodesC.txt',status='unknown')
798 Endif
799!
800! Redo again with same three MaxPool sizes but with a parallel solve loop.
801!
802 parallel = .true.
803 maxdyn = .false.
804 maxpool = 8
805 time0 = omp_get_wtime()
806 Call bandb_solve
807 time8p = omp_get_wtime() - time0
808 stime8p = stime
809 numsolve8p = numsolve
810 write(11,*) 'Finished parallel with Maxpool=',maxpool
811 write(*,*) 'Finished parallel with Maxpool=',maxpool
812 maxpool = 16
813 time0 = omp_get_wtime()
814 Call bandb_solve
815 time16p = omp_get_wtime() - time0
816 stime16p = stime
817 numsolve16p = numsolve
818 write(11,*) 'Finished parallel with Maxpool=',maxpool
819 write(*,*) 'Finished parallel with Maxpool=',maxpool
820 maxpool = 32
821 time0 = omp_get_wtime()
822 Call bandb_solve
823 time32p = omp_get_wtime() - time0
824 stime32p = stime
825 numsolve32p = numsolve
826 maxdyn = .true.
827 time0 = omp_get_wtime()
828 Call bandb_solve
829 timedp = omp_get_wtime() - time0
830 stimedp = stime
831 numsolvedp = numsolve
832 write(11,*) 'Finished parallel with Maxpool=',maxpool
833 write(*,*) 'Finished parallel with Maxpool=',maxpool
834!
835! Clean up
836!
837 DO thread = 1, maxthread
838 coi_error = max( coi_error, coi_free( cntvect1(thread)%CntVect ) )
839 Enddo
840 Deallocate( cntvect1 )
841
842 write(*,*) 'Number of Binary variables in model=',nn
843 write(*,*) 'Number of threads in parallel experiments=',maxthread
844 write(*,*)
845 write(*,"('Loop type MaxPool Total time % Solve Time % Number of Solves Parallel Speedup')")
846 write(*,*)
847 write(*,1000) 'Sequential',' 1',time1, 100.d0, stime1, 100.d0, numsolve1
848 write(*,1000) 'Sequential',' 8',time8s ,time8s /time1*100.d0, stime8s ,stime8s /stime1*100.d0, numsolve8s
849 write(*,1000) 'Sequential','16',time16s,time16s/time1*100.d0, stime16s,stime16s/stime1*100.d0, numsolve16s
850 write(*,1000) 'Sequential','32',time32s,time32s/time1*100.d0, stime32s,stime32s/stime1*100.d0, numsolve32s
851 write(*,1000) 'Sequential',' D',timeds ,timeds /time1*100.d0, stimeds ,stimeds /stime1*100.d0, numsolveds
852 write(*,1000) 'Parallel ',' 8',time8p ,time8p /time1*100.d0, stime8p ,stime8p /stime1*100.d0, numsolve8p , time8s/time8p
853 write(*,1000) 'Parallel ','16',time16p,time16p/time1*100.d0, stime16p,stime16p/stime1*100.d0, numsolve16p, time16s/time16p
854 write(*,1000) 'Parallel ','32',time32p,time32p/time1*100.d0, stime32p,stime32p/stime1*100.d0, numsolve32p, time32s/time32p
855 write(*,1000) 'Parallel ',' D',timedp ,timedp /time1*100.d0, stimedp ,stimedp /stime1*100.d0, numsolvedp , timeds/timedp
856
8571000 Format(a,5x,a2,f12.3,f9.1,f12.3,f9.1,i12,f19.2)
858
859 write(*,*)
860 write(*,*) 'End of QP Branch & Bound example.'
861
862#if defined (nocomp)
863 call flog( "Successful Solve", 0 )
864#else
865 if ( numsolve8s .ne. numsolve8p ) then
866 call flog( "Difference between number of nodes in sequential and parallel solve for MaxPool=8", 1 )
867 else if ( numsolve16s .ne. numsolve16p ) then
868 call flog( "Difference between number of nodes in sequential and parallel solve for MaxPool=16", 1 )
869 else if ( numsolve32s .ne. numsolve32p ) then
870 call flog( "Difference between number of nodes in sequential and parallel solve for MaxPool=32", 1 )
871 else if ( numsolveds .ne. numsolvedp ) then
872 call flog( "Difference between number of nodes in sequential and parallel solve for Dynamic MaxPool", 1 )
873 else
874 call flog( "Successful Solve", 0 )
875 endif
876#endif
877
878End Program qp_threads
879
880
881!> Sets runtime options
882!!
883!! @include{doc} option_params.dox
884Integer Function bb_option( ncall, rval, ival, lval, usrmem, name )
885#if defined(itl)
886!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: BB_Option
887#endif
888 integer ncall, ival, lval
889 character(Len=*) :: name
890 real*8 rval
891 real*8 usrmem(*) ! optional user memory
892
893 Select Case (ncall)
894 case (1)
895 name = 'locon3'
896 ival = 0
897 case (2)
898 name = 'lotria'
899 ival = 0
900 case default
901 name = ' '
902 end Select
903 bb_option = 0
904
905end Function bb_option
906
907Subroutine bandb_solve
908
909 Use bb_data
910 Use qp_data
911 implicit none
912 Integer i
913 Type (BBNode), Pointer :: Node
914!
915! Initialize the first node for B&B tree
916!
917 numsolve = 0
918 numintfound = 0
919 Nullify(freenode)
920 Nullify(actnode)
921 intfound = .false.
922 totnodes = 0 ! Total number of nodes in use
923 actnodes = 0 ! Number of active nodes
924 freenodes = 0 ! Number of free nodes in a free list
925 stime = 0 ! Reset Solve loop timer
926!
927 Call get_freenode(node) ! Create a root node with no binaries fixed
928 node%Depth = 0;
929 do i = 1, nvar
930 node%Values(i) = 0.d0
931 enddo
932 do i = 1, nbin
933 node%Fixed(i) = -1
934 enddo
935 node%Object = 1.e29 ! We need a value for the first selection before
936 ! the node has been solved
937 Call add_activenode(node)
938 do while ( actnodes > 0 )
939 write(*,*) 'Creating WorkPool. Active nodes=',actnodes
940 Call create_workpool
941 write(*,*) 'WorkPool with ',workpool,' nodes selected. Remaining Active nodes=',actnodes
942 if ( workpool > 0 ) then
943 Call solve_workpool
944 endif
945 Enddo
946 Call deallocate_nodes ! Clean up again in the node list
947
948 if ( intfound ) then
949 write(*,*)
950 write(*,*) 'An integer solution with objective',bestint,' has been found.'
951 write(*,*)
952 write(*,*) 'Solution values:'
953 do i = 1, 2*nn
954 write(*,"(I5,1pe20.12)") i, solution(i)
955 enddo
956 else
957 write(*,*) 'No integer solution has been found'
958 endif
959 write(*,*) 'Number of Solves: ', numsolve
960 write(*,*) 'Number of Integer Solutions:', numintfound
961 write(*,*)
962
963End Subroutine bandb_solve
964!
965! ============================================================================
966! Define information about the model:
967!
968
969!> Define information about the model
970!!
971!! @include{doc} readMatrix_params.dox
972Integer Function qp_readmatrix( lower, curr, upper, vsta, type, rhs, esta, &
973 colsta, rowno, value, nlflag, n, m, nz, &
974 usrmem )
975#if defined(itl)
976!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_ReadMatrix
977#endif
978 Use bb_data
979 Use qp_data
980 implicit none
981 integer, intent (in) :: n ! number of variables
982 integer, intent (in) :: m ! number of constraints
983 integer, intent (in) :: nz ! number of nonzeros
984 real*8, intent (in out), dimension(n) :: lower ! vector of lower bounds
985 real*8, intent (in out), dimension(n) :: curr ! vector of initial values
986 real*8, intent (in out), dimension(n) :: upper ! vector of upper bounds
987 integer, intent (in out), dimension(n) :: vsta ! vector of initial variable status
988 ! (not defined here)
989 integer, intent (out), dimension(m) :: type ! vector of equation types
990 integer, intent (in out), dimension(m) :: esta ! vector of initial equation status
991 ! (not defined here)
992 real*8, intent (in out), dimension(m) :: rhs ! vector of right hand sides
993 integer, intent (in out), dimension(n+1) :: colsta ! vector with start of column indices
994 integer, intent (out), dimension(nz) :: rowno ! vector of row numbers
995 integer, intent (in out), dimension(nz) :: nlflag ! vector of nonlinearity flags
996 real*8, intent (in out), dimension(nz) :: value ! vector of matrix values
997 Type (bbnode) :: usrmem
998 integer :: i, j
999!
1000! Information about Variables:
1001! Default: Lower = -Inf, Curr = 0, and Upper = +inf.
1002! Default: the status information in Vsta is not used.
1003!
1004! Nondefault: Lower bound = 0 for all variables
1005!
1006 do i = 1, nvar
1007 lower(i) = 0.0d0
1008 enddo
1009 do i = 1, nn
1010 upper(nn+i) = 1.d0
1011 enddo
1012 do i = 1, nvar
1013 curr(i) = usrmem%Values(i)
1014 enddo
1015!
1016! Fix the binary variables if indicated by Fixed
1017!
1018 do i = 1, nbin
1019 if ( usrmem%fixed(i) == 0 ) then
1020 j = binmap(i)
1021 lower(j) = 0.d0
1022 upper(j) = 0.d0
1023 elseif ( usrmem%fixed(i) == 1 ) then
1024 j = binmap(i)
1025 lower(j) = 1.d0
1026 upper(j) = 1.d0
1027 endif
1028 enddo
1029!
1030! Information about Constraints:
1031! Default: Rhs = 0
1032! Default: the status information in Esta and the function
1033! value in FV are not used.
1034! Default: Type: There is no default.
1035! 0 = Equality,
1036! 1 = Greater than or equal,
1037! 2 = Less than or equal,
1038! 3 = Non binding.
1039!
1040! Constraint 1 (Objective)
1041! Rhs = 0.0 (default) and type Non binding
1042!
1043 type(1) = 3
1044!
1045! Constraint 2 (Sum of all X-vars = 1)
1046! Rhs = 1 and type Equality
1047!
1048 rhs(2) = 1.d0
1049 type(2) = 0
1050!
1051! Constraint 3 (average of all mean = target)
1052! Rhs = Target and type Eguality
1053!
1054 rhs(3) = target
1055 type(3) = 0
1056!
1057! Constraint 4 (sum of all Bin-vars =L= Maxbin)
1058! Rhs = Maxbin and type L
1059!
1060 rhs(4) = maxbin
1061 type(4) = 2
1062!
1063! Next N constraints ( x(i) - bin(i) =L= 0 )
1064! Rhs = 0 (default) and type L
1065!
1066 do i = 1,nn
1067 type(4+i) = 2
1068 enddo
1069!
1070! Information about the Jacobian.
1071!
1072! Colsta = Start of column indices (No Defaults):
1073! Rowno = Row indices
1074! Value = Value of derivative (by default only linear
1075! derivatives are used)
1076! Nlflag = 0 for linear and 1 for nonlinear derivative
1077! (not needed for completely linear models)
1078!
1079 j = 1 ! counter for current nonzero
1080 do i = 1, nn
1081 colsta(i) = j
1082 rowno(j) = 1 ! Objective - nonlinear
1083 nlflag(j) = 1
1084 j = j + 1
1085 rowno(j) = 2 ! sum(x) = 1
1086 value(j) = 1.d0
1087 nlflag(j) = 0
1088 j = j + 1
1089 rowno(j) = 3 ! sum(mean*x) = target
1090 value(j) = mean(i)
1091 nlflag(j) = 0
1092 j = j + 1
1093 rowno(j) = 4+i ! x(i) - bin(i) =l= 0
1094 value(j) = 1.d0
1095 nlflag(j) = 0
1096 j = j + 1
1097 enddo
1098 do i = 1, nn
1099 colsta(nn+i) = j
1100 rowno(j) = 4 ! sum(bin) =l= Maxbin
1101 value(j) = 1.d0
1102 nlflag(j) = 0
1103 j = j + 1
1104 rowno(j) = 4+i ! x(i) - bin(i) =l= 0
1105 value(j) = -1.d0
1106 nlflag(j) = 0
1107 j = j + 1
1108 enddo
1109 colsta(2*nn+1) = j
1110
1111 qp_readmatrix = 0 ! Return value means OK
1112
1113end Function qp_readmatrix
1114!
1115!==========================================================================
1116! Compute nonlinear terms and non-constant Jacobian elements
1117!
1118
1119!> Compute nonlinear terms and non-constant Jacobian elements
1120!!
1121!! @include{doc} fdeval_params.dox
1122Integer Function qp_fdeval( x, g, jac, rowno, jcnm, mode, ignerr, errcnt, &
1123 n, nz, thread, usrmem )
1124#if defined(itl)
1125!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_FDEval
1126#endif
1127 Use bb_data
1128 Use qp_data
1129 implicit none
1130 integer, intent (in) :: n ! number of variables
1131 integer, intent (in) :: rowno ! number of the row to be evaluated
1132 integer, intent (in) :: nz ! number of nonzeros in this row
1133 real*8, intent (in), dimension(n) :: x ! vector of current solution values
1134 real*8, intent (in out) :: g ! constraint value
1135 real*8, intent (in out), dimension(n) :: jac ! vector of derivatives for current constraint
1136 integer, intent (in), dimension(nz) :: jcnm ! list of variables that appear nonlinearly
1137 ! in this row. Ffor information only.
1138 integer, intent (in) :: mode ! evaluation mode: 1 = function value
1139 ! 2 = derivatives, 3 = both
1140 integer, intent (in) :: ignerr ! if 1 then errors can be ignored as long
1141 ! as errcnt is incremented
1142 integer, intent (in out) :: errcnt ! error counter to be incremented in case
1143 ! of function evaluation errors.
1144 integer, intent (in) :: thread !
1145 Type (bbnode) :: usrmem
1146 Integer :: i, j
1147 real*8 :: sum
1148!
1149! Row 1: The objective function is nonlinear
1150!
1151 if ( rowno .eq. 1 ) then
1152!
1153! Mode = 1 or 3. Function value: G = x * Q * x / 2
1154!
1155 if ( mode .eq. 1 .or. mode .eq. 3 ) then
1156 sum = 0.d0
1157 do i = 1, nn
1158 do j = 1, nn
1159 sum = sum + x(i)*q(i,j)*x(j)
1160 enddo
1161 enddo
1162 g = sum / 2.d0
1163 endif
1164!
1165! Mode = 2 or 3: Derivative values: Q*x
1166!
1167 if ( mode .eq. 2 .or. mode .eq. 3 ) then
1168 do i = 1, nn
1169 sum = 0
1170 do j = 1, nn
1171 sum = sum + q(i,j) * x(j)
1172 enddo
1173 jac(i) = sum
1174 enddo
1175 endif
1176 qp_fdeval = 0
1177 Else
1178!
1179! Other rows are linear and will not be called.
1180!
1181 qp_fdeval = 1
1182 endif
1183 qp_fdeval = 0
1184
1185end Function qp_fdeval
1186
1187
1188!> Computes the second derivative of a constraint in a direction
1189!!
1190!! @include{doc} 2DDir_params.dox
1191Integer Function qp_2ddir( X, DX, D2G, ROWNO, JCNM, NODRV, N, NJ, THREAD, USRMEM )
1192#if defined(itl)
1193!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_2DDir
1194#endif
1195 Use qp_data
1196 Use bb_data
1197 implicit none
1198 INTEGER, Intent(IN) :: rowno, n, nj, thread
1199 INTEGER, Intent(IN OUT) :: nodrv
1200 INTEGER, Dimension(NJ), Intent(IN) :: jcnm
1201 real*8, Dimension(N), Intent(IN) :: x
1202 real*8, Dimension(N), Intent(IN) :: dx
1203 real*8, Dimension(N), Intent(OUT) :: d2g
1204 Type (bbnode) :: usrmem
1205 integer :: i,j
1206 real(8) :: sum
1207!
1208! Is only called for Irow = 1, the nonlinear objective
1209! Return H*Dx = Q*Dx in D2G
1210!
1211 if ( rowno .eq. 1 ) then
1212 do i = 1, nn
1213 sum = 0.d0
1214 do j = 1, nn
1215 sum = sum + q(i,j) * dx(j)
1216 enddo
1217 d2g(i) = sum
1218 enddo
1219 qp_2ddir = 0
1220 else
1221 qp_2ddir = 1
1222 endif
1223
1224End Function qp_2ddir
1225
1226
1227!> Specify the structure of the Lagrangian of the Hessian
1228!!
1229!! @include{doc} 2DLagrStr_params.dox
1230Integer Function qp_2dlagrstr( ROWNO, COLNO, NODRV, N, M, NHESS, USRMEM )
1231#if defined(itl)
1232!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_2DLagrStr
1233#endif
1234 Use qp_data
1235 Use bb_data
1236 implicit none
1237 INTEGER n, m, nhess, nodrv
1238 INTEGER rowno(nhess), colno(nhess)
1239 Type (bbnode) :: usrmem
1240
1241 Integer :: i,j,k
1242
1243 k = 1
1244 Do i = 1, nn
1245 do j = i,nn
1246 rowno(k) = j
1247 colno(k) = i
1248 k = k + 1
1249 Enddo
1250 Enddo
1251 qp_2dlagrstr = 0
1252
1253End function qp_2dlagrstr
1254
1255
1256!> Compute the Lagrangian of the Hessian
1257!!
1258!! @include{doc} 2DLagrVal_params.dox
1259Integer Function qp_2dlagrval( X, U, ROWNO, COLNO, VALUE, NODRV, N, M, NHESS, USRMEM )
1260#if defined(itl)
1261!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: QP_2DLagrVal
1262#endif
1263 Use qp_data
1264 Use bb_data
1265 implicit none
1266 INTEGER n, m, nhess, nodrv
1267 real*8 x(n), u(m), value(nhess)
1268 INTEGER rowno(nhess), colno(nhess)
1269 Type (bbnode) :: usrmem
1270
1271 Integer :: i,j,k
1272
1273 k = 1
1274 Do i = 1, nn
1275 do j = i,nn
1276 value(k) = q(i,j)*u(1)
1277 k = k + 1
1278 Enddo
1279 Enddo
1280 qp_2dlagrval = 0
1281
1282End function qp_2dlagrval
1283!
1284! Toutines for Message, Status, Solution, and ErrMsg called
1285! BB_Message, BB_Status, BB_Solution, and BB_ErrMsg.
1286! They are very simple and we avoid written output since overlapping
1287! threads can give unreadable output.
1288!
1289Integer Function bb_status( MODSTA, SOLSTA, ITER, OBJVAL, USRMEM )
1290#if defined(itl)
1291!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: BB_Status
1292#endif
1293!
1294! Simple implementation in which we just store the status and
1295! objective value in the current node.
1296!
1297 Use bb_data
1298 IMPLICIT NONE
1299 INTEGER, Intent(IN) :: modsta, solsta, iter
1300 real*8, Intent(IN) :: objval
1301 Type (bbnode) :: usrmem
1302
1303 usrmem%Solstat = solsta
1304 usrmem%Modstat = modsta
1305 usrmem%Object = objval
1306 usrmem%Node_Stat = 1
1307
1308 bb_status = 0
1309 return
1310
1311END Function bb_status
1312
1313Integer Function bb_solution( XVAL, XMAR, XBAS, XSTA, YVAL, YMAR, YBAS, YSTA, N, M, USRMEM )
1314#if defined(itl)
1315!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: BB_Solution
1316#endif
1317!
1318! Store the solution in the Values field of the solved node.
1319! Also check for integrality. If done here it is inside the
1320! parallel loop.
1321!
1322 Use bb_data
1323 IMPLICIT NONE
1324 INTEGER, Intent(IN) :: n, m
1325 INTEGER, Intent(IN), Dimension(N) :: xbas, xsta
1326 INTEGER, Intent(IN), Dimension(M) :: ybas, ysta
1327 real*8, Intent(IN), Dimension(N) :: xval, xmar
1328 real*8, Intent(IN), Dimension(M) :: yval, ymar
1329 Type (bbnode) :: usrmem
1330
1331 Integer :: numintvar, i, j
1332 Real(8) :: dist
1333
1334 DO i = 1, n
1335 usrmem%Values(i) = xval(i)
1336 ENDDO
1337 bb_solution = 0
1338
1339 numintvar = 0
1340 do i = 1, nn
1341 j = binmap(i)
1342 dist = min( usrmem%values(j), 1.d0-usrmem%values(j) )
1343 if ( dist <= bintol ) numintvar = numintvar + 1
1344 enddo
1345 usrmem%IsInteger = (numintvar == nn)
1346
1347END Function bb_solution
1348
1349Integer Function bb_message( SMSG, DMSG, NMSG, LLEN, USRMEM, MSGV )
1350#if defined(itl)
1351!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: BB_Message
1352#endif
1353!
1354! Simple implementation in which does nothing. Messages from different
1355! threads can give a confusing output.
1356!
1357 Use bb_data
1358 Use proginfop
1359 IMPLICIT NONE
1360
1361 INTEGER, Intent(IN) :: smsg, dmsg, nmsg
1362 INTEGER, Intent(IN) , Dimension(*) :: llen
1363 CHARACTER(Len=133), Intent(IN), Dimension(*) :: msgv
1364 Type (bbnode) :: usrmem
1365
1366 Integer :: i, thread, num
1367
1368 if ( debug ) then
1369!
1370! write to document file
1371!
1372 thread = usrmem%Sthread
1373 num = usrmem%Number
1374 do i = 1, dmsg
1375 write(11+thread,"(i4,': ',(A))") num,msgv(i)(1:llen(i))
1376 enddo
1377 endif
1378 bb_message = 0
1379 return
1380
1381END Function bb_message
1382
1383Integer Function bb_errmsg( ROWNO, COLNO, POSNO, MSGLEN, USRMEM, MSG )
1384#if defined(itl)
1385!DEC$ ATTRIBUTES STDCALL, REFERENCE, NOMIXED_STR_LEN_ARG :: BB_ErrMsg
1386#endif
1387!
1388! Simple implementation in which does nothing. Messages from different
1389! threads can give a confusing output.
1390!
1391 Use bb_data
1392 IMPLICIT NONE
1393
1394 INTEGER, Intent(IN) :: rowno, colno, posno, msglen
1395 CHARACTER(Len=*), Intent(IN) :: msg
1396 Type (bbnode) :: usrmem
1397 bb_errmsg = 0
1398 return
1399
1400END Function bb_errmsg
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_2dlagrstr(cntvect, coi_2dlagrstr)
define callback routine for providing the structure of the second derivatives of the Lagrangian.
integer function coidef_2ddir(cntvect, coi_2ddir)
define callback routine for computing the second derivative for a constraint in a direction.
integer function coidef_option(cntvect, coi_option)
define callback routine for defining runtime options.
integer function coidef_2dlagrval(cntvect, coi_2dlagrval)
define callback routine for computing the values of the second derivatives of the Lagrangian.
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_numhess(cntvect, numhess)
defines the Number of Hessian Nonzeros.
Definition coistart.f90:513
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
#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
int select_var(struct BBNode *node)
Definition mp_qpbandb.c:240
void add_activenode(struct BBNode *node)
Definition mp_qpbandb.c:219
void return_freenode(struct BBNode *node)
Definition mp_qpbandb.c:181
void solve_workpool()
Definition mp_qpbandb.c:360
void bandb_solve()
Definition mp_qpbandb.c:502
void add2workpool(struct BBNode *node)
Definition mp_qpbandb.c:232
void create_workpool()
Definition mp_qpbandb.c:266
#define maxmaxpool
Definition mp_qpbandb.c:131
struct BBNode * get_freenode()
Definition mp_qpbandb.c:191
void remove_activenode(struct BBNode *node)
Definition mp_qpbandb.c:168
void list_nodes()
Definition mp_qpbandb.c:153
void list_node(struct BBNode *node)
Definition mp_qpbandb.c:144
integer function bb_message(smsg, dmsg, nmsg, llen, usrmem, msgv)
integer function bb_solution(xval, xmar, xbas, xsta, yval, ymar, ybas, ysta, n, m, usrmem)
program qp_threads
Main program. A simple setup and call of CONOPT.
integer function bb_status(modsta, solsta, iter, objval, usrmem)
integer function bb_option(ncall, rval, ival, lval, usrmem, name)
Sets runtime options.
integer function bb_errmsg(rowno, colno, posno, msglen, usrmem, msg)
#define nj
Definition mp_trans.c:46
type(controlv), dimension(:), pointer cntvect1
logical maxdyn
integer maxthread
subroutine deallocate_nodes
subroutine flog(msg, code)
Definition comdeclp.f90:42
logical debug
Definition comdeclp.f90:17
subroutine startup
Definition comdeclp.f90:25
integer, parameter maxbin
real(8) target
integer, parameter nn
real(8), dimension(nn) mean
integer, parameter nbin
integer, parameter nvar
real(8), dimension(nn, nn) q
integer function qp_readmatrix(lower, curr, upper, vsta, type, rhs, esta, colsta, rowno, value, nlflag, n, m, nz, usrmem)
Define information about the model.
Definition qp1.f90:209
integer function qp_fdeval(x, g, jac, rowno, jcnm, mode, ignerr, errcnt, n, nz, thread, usrmem)
Compute nonlinear terms and non-constant Jacobian elements.
Definition qp1.f90:296
integer function qp_2ddir(x, dx, d2g, rowno, jcnm, nodrv, n, nj, thread, usrmem)
Computes the second derivative of a constraint in a direction.
Definition qp2.f90:300
integer function qp_2dlagrval(x, u, rowno, colno, value, nodrv, n, m, nhess, usrmem)
Compute the Lagrangian of the Hessian.
Definition qp3.f90:317
integer function qp_2dlagrstr(rowno, colno, nodrv, n, m, nhess, usrmem)
Specify the structure of the Lagrangian of the Hessian.
Definition qp3.f90:294