/usr/include/linbox/algorithms/hybrid-det.h is in liblinbox-dev 1.3.2-1.1.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 | /* linbox/algorithms/hybrid-det.h
* Copyright (C) 2005 Anna Urbanska
*
* Written by Anna Urbanska <aniau@astronet.pl>
* Modified by JGD <Jean-Guillaume.Dumas@imag.fr>
*
* ========LICENCE========
* This file is part of the library LinBox.
*
* LinBox is free software: you can redistribute it and/or modify
* it under the terms of the GNU Lesser General Public
* License as published by the Free Software Foundation; either
* version 2.1 of the License, or (at your option) any later version.
*
* This library is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
* Lesser General Public License for more details.
*
* You should have received a copy of the GNU Lesser General Public
* License along with this library; if not, write to the Free Software
* Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
* ========LICENCE========
*/
#ifndef __LINBOX_hybrid_det_H
#define __LINBOX_hybrid_det_H
#include "linbox/linbox-config.h"
#include "linbox/util/debug.h"
#include "linbox/blackbox/sparse.h"
#include "linbox/field/PID-integer.h"
#include "linbox/randiter/random-prime.h"
#include "linbox/algorithms/rational-solver.h"
#include "linbox/algorithms/last-invariant-factor.h"
#include "linbox/field/modular.h"
#include "linbox/algorithms/cra-domain.h"
#include "linbox/randiter/random-prime.h"
#include "linbox/algorithms/matrix-hom.h"
#include "linbox/solutions/det.h"
// #define _LB_H_DET_TIMING
namespace LinBox
{
template <class Blackbox, class MyMethod>
struct IntegerModularDetReduced {
private:
const Blackbox &A;
const MyMethod &M;
/* contains the factor by which to divide det */
integer beta;
size_t factor;
size_t iter_count;
size_t iter_count2;
typename Vector<PID_integer>::Dense moduli;
public:
typename Vector<PID_integer>::Dense primes;
size_t iterations()
{
return iter_count;
}
size_t iterations2()
{
return iter_count2;
}
#if 0
int iter_count;
int iter_count2;
Vector <PID_integer>:: Dense moduli;
Vector <PID_integer>:: Dense primes;
#endif
IntegerModularDetReduced(const Blackbox& b, const MyMethod& n, const integer& divisor, const size_t& f) :
A(b), M(n), beta(divisor), factor(f)
{
moduli.resize(factor);
primes.resize(factor);
iter_count = 0;
iter_count2 = 0;
}
template<typename Field>
typename Field::Element& operator()(typename Field::Element& d, const Field& F)
{
if (beta > 1) {
if (iter_count2 < factor) {
Field D(primes[iter_count2]);
typename Field::Element pbeta;
typename Field::Element kbeta;
typename Field::Element current_moduli;
D.init(pbeta, beta);
D.init(current_moduli, moduli[iter_count2]);
D.div(kbeta, current_moduli,pbeta);
++this->iter_count2;
return d=kbeta;
}
}
typedef typename Blackbox::template rebind<Field>::other FBlackbox;
FBlackbox Ap(A,F);
detin( d, Ap, M);
if (beta > 1) {
typename Field::Element y;
F.init(y,beta);
F.div(d,d,y);
}
if (iter_count < factor) {
moduli[iter_count] = d;
}
++this->iter_count;
return d;
}
void Beta(Integer& b) { beta = b; iter_count2=0;}
};
/** \brief Compute the determinant of A over the integers
*
* The determinant of a linear operator A, represented as a
* black box, is computed over the integers.
*
* This variant is a hybrid between Last Invariant Factor and Chinese Remaindering
* It performs several determinants mod primes
* Then switches to the LIF method, producing a factor of det.
* It then comes back to the CRA if necessary to compute
* the remaining (usually small) factor of the determinant.
*
* @param d Field element into which to store the result
* @param A Black box of which to compute the determinant
* @param tag explicit over the integers
* @param M may be a Method::BlasElimination (default) or a Method::Wiedemann.
\ingroup solutions
*/
template <class Blackbox, class MyMethod>
typename Blackbox::Field::Element & lif_cra_det (typename Blackbox::Field::Element &d,
const Blackbox &A,
const RingCategories::IntegerTag &tag,
const MyMethod &M)
{
//commentator().setReportStream(std::cout);
typedef Modular<double> myModular;
typedef typename Blackbox::Field Integers;
typedef typename Integers::Element Integer_t;
commentator().start ("Integer Determinant - hybrid version ", "det");
size_t myfactor=5;
size_t early_counter=0;
//double a = 0.0;//0.28;//0.0013//if time(lif)/time(10lu) < a * log(lif) then calculate bonus
Integer_t lif = 1;
Integer_t bonus = 1;
Integer_t beta = 1;
d=1;
double p_size = 26-(int)ceil(log((double)A.rowdim())*0.7213475205);
RandomPrimeIterator genprime( (Integer_t)p_size );
//cout << "prime size: " << p_size << "\n";
EarlySingleCRA<myModular> cra(4UL);
IntegerModularDetReduced<Blackbox,MyMethod> iteration(A, M, beta,myfactor);
#if 0
if (A.rowdim() < 200 ) {
cra(d,iteration,genprime);
commentator().stop ( "first step", NULL, "det");
commentator().report (Commentator::LEVEL_NORMAL, INTERNAL_DESCRIPTION)
<< "Iterations done " << iteration.iterations() << "\n";
}
else {}
#endif
Integer_t res;
#ifdef _LB_H_DET_TIMING
Timer BT;
double time1, time2;
BT.start();
#endif
{
++genprime;
myModular D(*genprime);
iteration.primes[early_counter] = *genprime;
myModular::Element r;
D.init(r,0);
cra.initialize( D, iteration(r, D));
++early_counter;
}
while ( early_counter < myfactor && !cra.terminated() ) {
++genprime;
while (cra.noncoprime(*genprime)) ++genprime;
myModular D(*genprime);
iteration.primes[early_counter] = *genprime;
// prime(p, early_counter);
myModular::Element r;
D.init(r,0);
cra.progress( D, iteration(r, D));
++early_counter;
}
#ifdef _LB_H_DET_TIMING
BT.stop();
time1 = BT.usertime()/myfactor;
if (time1 < 0) time1 = 0;
#endif
cra.result(res);
if (early_counter < myfactor) {
/* determinant found */
commentator().report (Commentator::LEVEL_NORMAL, INTERNAL_DESCRIPTION) << myfactor << "\n";
commentator().stop ( "first step", NULL, "det");
commentator().report (Commentator::LEVEL_NORMAL, INTERNAL_DESCRIPTION)
<< "Iterations done " << iteration.iterations()<< "\n";
return d=res;
}
commentator().report (Commentator::LEVEL_NORMAL, INTERNAL_DESCRIPTION) << "no very early termination \n";
/* turn to LU when matrix size small - not to be used at the moment */
#if 0
if (A.rowdim() < 50 ) {
while ( !cra.terminated() ) {
genprime.randomPrime(p);
myModular D(p);
iteration.primes[early_counter] = p;
// prime(p, early_counter);
myModular::Element r;
D.init(r,0);
cra.progress( D, iteration(r, D));
++early_counter;
}
commentator().stop ( "zero step", NULL, "det");
commentator().report (Commentator::LEVEL_NORMAL, INTERNAL_DESCRIPTION)
<< "Iterations done " << iteration.iterations() << "\n";
return d;
}
RandomPrime genprime1( 26-(int)ceil(log((double)A.rowdim())*0.7213475205));
Integers ZZ;
RationalSolver < Integers , myModular, RandomPrime, DixonTraits > RSolver(A. field(), genprime);
#endif
RationalSolver < Integers , myModular, RandomPrimeIterator, DixonTraits > RSolver;
typename Vector<Integers>:: Dense r_num1 (A. coldim());
LastInvariantFactor < Integers ,RationalSolver < Integers, myModular, RandomPrimeIterator, DixonTraits > > LIF(RSolver);
#ifdef _LB_H_DET_TIMING
BT.start();
#endif
if (LIF.lastInvariantFactor1(lif, r_num1, A)==0) {
//if (lif==0)
d = 0;
commentator().stop ("is 0", NULL, "det");
return d;
}
#ifdef _LB_H_DET_TIMING
BT.stop();
time2 = BT.usertime();
if (time2 < 0) time2 =0;
#endif
//commentator().report (Commentator::LEVEL_NORMAL, INTERNAL_DESCRIPTION)
// << "5 LU time: " << time1 << " LIF time: " << time2 << ")\n";
commentator().report (Commentator::LEVEL_NORMAL, INTERNAL_DESCRIPTION) << "lif calculated\n";
#if 0
LIF.lastInvariantFactor_Bonus(lif, bonus, A);
if (lif==0) {
d = 0;
commentator().stop ("done", NULL, "det");
return d;
}
if (bonus == 1) {
d = lif;
commentator().stop ("done", NULL, "det");
return d;
}
#endif
beta = lif*bonus;
iteration.Beta(beta);
//RandomPrime genprime2( 26-(int)ceil(log((double)A.rowdim())*0.7213475205));
EarlySingleCRA< Modular<double> > cra2(4UL);
Integer_t k = 1;
early_counter = 0;
while ( early_counter < myfactor && !cra2.terminated() ) {
myModular D(iteration.primes[early_counter]);
myModular::Element r;
D.init(r,0);
cra2.progress( D, iteration(r, D) );
++early_counter;
}
if (early_counter < myfactor) {
/* determinant found */
k = cra2.result(res);
commentator().stop ("second step ", NULL, "det");
commentator().report (Commentator::LEVEL_NORMAL, INTERNAL_DESCRIPTION)
<< "Iterations done " << iteration.iterations()<< "(" << iteration.iterations2() << ")\n";
}
else if (0/* time2 < a*log2(lif)*time1/p_size*/) {
typename Vector<Integers>:: Dense r_num2 (A. coldim());
Integer_t lif2=1;
LIF.lastInvariantFactor1(lif2,r_num2,A);
LIF.bonus(bonus,lif, lif2, r_num1,r_num2);
//lif2 = lif;
if ((bonus > 1) || (lif2 != lif)) {
//cout << "lif: " <<lif <<",\n "<< lif2<< "przed\n";
lif = lcm(lif, lif2);
//cout << "lif: " <<lif << "po";
beta=lif*bonus;
iteration.Beta(beta);
//iteration.Restart2() // included in Beta();
//iteration.Moduli(moduli);
//iteration.Primes(primes);
k=1;
EarlySingleCRA< Modular<double> > cra3(4UL);
early_counter = 0;
while ( (early_counter < myfactor) && (!cra3.terminated() )) {
myModular D(iteration.primes[early_counter]);
myModular::Element r;
D.init(r,0);
cra3.progress( D, iteration(r, D));
++early_counter;
//iteration.Inc();
}
if (early_counter < myfactor) {
/* determinant found based on the initial LU */
cra3.result(k);
commentator().stop ("third step - recalc", NULL, "det");
commentator().report (Commentator::LEVEL_NORMAL, INTERNAL_DESCRIPTION)
<< "Iterations done " << iteration.iterations()<<"(" << iteration.iterations2() << ")\n";
//<< "bonus size " << log2(bonus) << "\n";
}
else {
/* enter the cra loop */
//cra3(k,iteration, genprime);
while (!cra3.terminated()) {
++genprime;
while (cra3.noncoprime(*genprime)) ++genprime;
myModular D(*genprime);
myModular::Element r;
D.init(r,0);
cra3.progress( D, iteration(r, D));
}
cra3.result(k);
commentator().stop ("third step, bonus > 1", NULL, "det");
commentator().report (Commentator::LEVEL_NORMAL, INTERNAL_DESCRIPTION)
<< "Iterations done " << iteration.iterations()<< "(" << iteration.iterations2() << ")\n";
//<< "bonus size " << log2(bonus) << "\n";
}
}
else {
//cra2(k,iteration, genprime);
while (!cra2.terminated()) {
++genprime;
while (cra2.noncoprime(*genprime)) ++genprime;
myModular D(*genprime);
myModular::Element r;
D.init(r,0);
cra2.progress( D, iteration(r, D));
}
cra2.result(k);
commentator().stop ("third step, bonus = 1", NULL, "det");
commentator().report (Commentator::LEVEL_NORMAL, INTERNAL_DESCRIPTION)
<< "Iterations done " << iteration.iterations()<< "(" << iteration.iterations2() << ")\n";
}
}
else {
while (!cra2.terminated()) {
++genprime;
while (cra2.noncoprime(*genprime)) ++genprime;
myModular D(*genprime);
myModular::Element r;
D.init(r,0);
cra2.progress( D, iteration(r, D));
}
cra2.result(k);
commentator().stop ("second step+", NULL, "det");
commentator().report (Commentator::LEVEL_NORMAL, INTERNAL_DESCRIPTION)
<< "Iterations done " << iteration.iterations()<< "(" << iteration.iterations2() << ")\n";
/* enter the cra loop */
//cra2(k,iteration, genprime);
}
//commentator().stop ("second step", NULL, "det");
//commentator().report (Commentator::LEVEL_NORMAL, INTERNAL_DESCRIPTION)
// << "Iterations done " << iteration.iterations() << "("
// << iteration.iterations2() << " )\n";
d = k*beta;
Integer_t tmp;
#ifdef _LB_H_DET_TIMING
commentator().report (Commentator::LEVEL_NORMAL, INTERNAL_DESCRIPTION)
<< "1 LU time: " << time1 << " LIF time: " << time2 << ")\n";
#endif
commentator().report(Commentator::LEVEL_NORMAL, INTERNAL_DESCRIPTION)
<< "det/lif " << k<< "\n";
//commentator().report(Commentator::LEVEL_NORMAL, INTERNAL_DESCRIPTION)
// << "determinant size " << Integers::log2(tmp,abs(d))<< " lif size "<< Integers::log2(tmp,beta) << "\n";
return d ;
}
#if 0
template <class Integers, class MyMethod>
typename Integers::Element & lif_cra_det (typename Integers::Element &d,
const SparseMatrix<Integers> &A,
const RingCategories::IntegerTag &tag,
const MyMethod &M)
{
//commentator().setReportStream(std::cout);
typedef Modular<double> myModular;
//typedef PID_integer Integers;
typedef typename Integers::Element Integer;
commentator().start ("Integer Determinant - hybrid version for sparse matrices", "det");
size_t myfactor=5;
size_t early_counter=0;
//double a = 0.0;//0.28;//0.0013i //if time(lif)/time(10lu) < a * log(lif) then calculate bonus
Integer lif = 1;
Integer bonus = 1;
Integer beta = 1;
d=1;
double p_size = 26-(int)ceil(log((double)A.rowdim())*0.7213475205);
RandomPrime genprime( (Integer)p_size );
commentator().report (Commentator::LEVEL_NORMAL, INTERNAL_DESCRIPTION) << "prime size: " << p_size << "\n";
ChineseRemainder< myModular > cra(3UL);
IntegerModularDetReduced<SparseMatrix<Integers >,MyMethod> iteration(A, M, beta,myfactor);
#if 0
if (A.rowdim() < 200 ) {
cra(d,iteration,genprime);
commentator().stop ( "first step", NULL, "det");
commentator().report (Commentator::LEVEL_NORMAL, INTERNAL_DESCRIPTION)
<< "Iterations done " << iteration.iterations() << "\n";
}
else
#endif
Integer p;
Integer res;
Timer BT;
double time1, time2;
BT.start();
genprime.randomPrime(p);
myModular D(p);
iteration.primes[early_counter] = p;
myModular::Element r;
D.init(r,0);
cra.initialize( D, iteration(r, D));
++early_counter;
while ( early_counter < myfactor && !cra.terminated() ) {
genprime.randomPrime(p);
while (cra.noncoprime(p)) genprime.randomPrime(p);
myModular D(p);
iteration.primes[early_counter] = p;
// prime(p, early_counter);
myModular::Element r;
D.init(r,0);
cra.progress( D, iteration(r, D));
++early_counter;
}
BT.stop();
time1 = BT.usertime()/myfactor;
if (time1 < 0) time1 = 0;
cra.result(res);
if (early_counter < myfactor) {
//determinant found
commentator().stop ( "first step", NULL, "det");
commentator().report (Commentator::LEVEL_NORMAL, INTERNAL_DESCRIPTION)
<< "Iterations done " << iteration.iterations()<< "\n";
return d=res;
}
#if 0
cout << "no very early termination \n";
turn to LU when matrix size small - not to be used at the moment
if (A.rowdim() < 50 ) {
while ( !cra.terminated() ) {
genprime.randomPrime(p);
myModular D(p);
iteration.primes[early_counter] = p;
// prime(p, early_counter);
myModular::Element r;
D.init(r,0);
cra.progress( D, iteration(r, D));
++early_counter;
}
commentator().stop ( "zero step", NULL, "det");
commentator().report (Commentator::LEVEL_NORMAL, INTERNAL_DESCRIPTION)
<< "Iterations done " << iteration.iterations() << "\n";
return d;
}
RandomPrime genprime1( 26-(int)ceil(log((double)A.rowdim())*0.7213475205));
Integers ZZ;
#endif
typedef RationalSolver < Integers , myModular, RandomPrime, BlockHankelTraits > Solver;
Solver RSolver(A. field(), genprime);
typename Vector<Integers>:: Dense r_num1 (A. coldim());
LastInvariantFactor < Integers ,Solver > LIF(RSolver);
BT.start();
if (LIF.lastInvariantFactor1(lif, r_num1, A)==0) {
//if (lif==0)
d = 0;
commentator().stop ("is 0", NULL, "det");
return d;
}
BT.stop();
time2 = BT.usertime();
if (time2 < 0) time2 =0;
#if 0
commentator().report (Commentator::LEVEL_NORMAL, INTERNAL_DESCRIPTION)
<< "5 LU time: " << time1 << " LIF time: " << time2 << ")\n";
cout << "lif calculated\n";
LIF.lastInvariantFactor_Bonus(lif, bonus, A);
if (lif==0) {
d = 0;
commentator().stop ("done", NULL, "det");
return d;
}
if (bonus == 1) {
d = lif;
commentator().stop ("done", NULL, "det");
return d;
}
#endif
beta = lif*bonus;
iteration.Beta(beta);
//RandomPrime genprime2( 26-(int)ceil(log((double)A.rowdim())*0.7213475205));
ChineseRemainder< Modular<double> > cra2(3UL);
Integer k = 1;
early_counter = 0;
while ( early_counter < myfactor && !cra2.terminated() ) {
myModular D(iteration.primes[early_counter]);
myModular::Element r;
D.init(r,0);
cra2.progress( D, iteration(r, D));
++early_counter;
}
if (early_counter < myfactor) {
// determinant found
k = cra2.result(res);
commentator().stop ("second step ", NULL, "det");
commentator().report (Commentator::LEVEL_NORMAL, INTERNAL_DESCRIPTION)
<< "Iterations done " << iteration.iterations()<< "(" << iteration.iterations2() << ")\n";
}
else if (0) { //time2 < a*log2(lif)*time1/p_size)
typename Vector<Integers>:: Dense r_num2 (A. coldim());
Integer lif2=1;
LIF.lastInvariantFactor1(lif2,r_num2,A);
LIF.bonus(bonus,lif, lif2, r_num1,r_num2);
//lif2 = lif;
if ((bonus > 1) || (lif2 != lif)) {
//cout << "lif: " <<lif <<",\n "<< lif2<< "przed\n";
lif = lcm(lif, lif2);
//cout << "lif: " <<lif << "po";
beta=lif*bonus;
iteration.Beta(beta);
//iteration.Restart2() // included in Beta();
//iteration.Moduli(moduli);
//iteration.Primes(primes);
k=1;
ChineseRemainder< Modular<double> > cra3(3UL);
early_counter = 0;
while ( (early_counter < myfactor) && (!cra3.terminated() )) {
myModular D(iteration.primes[early_counter]);
myModular::Element r;
D.init(r,0);
cra3.progress( D, iteration(r, D));
++early_counter;
//iteration.Inc();
}
if (early_counter < myfactor) {
// determinant found based on the initial LU
cra3.result(k);
commentator().stop ("third step - recalc", NULL, "det");
commentator().report (Commentator::LEVEL_NORMAL, INTERNAL_DESCRIPTION)
<< "Iterations done " << iteration.iterations()<<"(" << iteration.iterations2() << ")\n";
//<< "bonus size " << log2(bonus) << "\n";
}
else {
// enter the cra loop
//cra3(k,iteration, genprime);
while (!cra3.terminated()) {
genprime.randomPrime(p);
while (cra3.noncoprime(p)) genprime.randomPrime(p);
myModular D(p);
myModular::Element r;
D.init(r,0);
cra3.progress( D, iteration(r, D));
}
cra3.result(k);
commentator().stop ("third step, bonus > 1", NULL, "det");
commentator().report (Commentator::LEVEL_NORMAL, INTERNAL_DESCRIPTION)
<< "Iterations done " << iteration.iterations()<< "(" << iteration.iterations2() << ")\n";
//<< "bonus size " << log2(bonus) << "\n";
}
}
else {
//cra2(k,iteration, genprime);
while (!cra2.terminated()) {
genprime.randomPrime(p);
while (cra2.noncoprime(p)) genprime.randomPrime(p);
myModular D(p);
myModular::Element r;
D.init(r,0);
cra2.progress( D, iteration(r, D));
}
cra2.result(k);
commentator().stop ("third step, bonus = 1", NULL, "det");
commentator().report (Commentator::LEVEL_NORMAL, INTERNAL_DESCRIPTION)
<< "Iterations done " << iteration.iterations()<< "(" << iteration.iterations2() << ")\n";
}
}
else {
while (!cra2.terminated()) {
genprime.randomPrime(p);
while (cra2.noncoprime(p)) genprime.randomPrime(p);
myModular D(p);
myModular::Element r;
D.init(r,0);
cra2.progress( D, iteration(r, D));
}
cra2.result(k);
commentator().stop ("second step+", NULL, "det");
commentator().report (Commentator::LEVEL_NORMAL, INTERNAL_DESCRIPTION)
<< "Iterations done " << iteration.iterations()<< "(" << iteration.iterations2() << ")\n";
// enter the cra loop
//cra2(k,iteration, genprime);
}
//commentator().stop ("second step", NULL, "det");
//commentator().report (Commentator::LEVEL_NORMAL, INTERNAL_DESCRIPTION)
// << "Iterations done " << iteration.iterations() << "("
// << iteration.iterations2() << " )\n";
d = k*beta;
Integer tmp;
commentator().report (Commentator::LEVEL_NORMAL, INTERNAL_DESCRIPTION)
<< "1 LU time: " << time1 << " LIF time: " << time2 << ")\n";
commentator().report(Commentator::LEVEL_NORMAL, INTERNAL_DESCRIPTION)
<< "det/lif " << k<< "\n";
//commentator().report(Commentator::LEVEL_NORMAL, INTERNAL_DESCRIPTION)
// << "determinant size " << Integers::log2(tmp,abs(d))<< " lif size "<< Integers::log2(tmp,beta) << "\n";
return d ;
}
#endif
} // end of LinBox namespace
#undef _LB_H_DET_TIMING
#endif // __LINBOX_hybrid_det_H
// vim:sts=8:sw=8:ts=8:noet:sr:cino=>s,f0,{0,g0,(0,:0,t0,+0,=s
// Local Variables:
// mode: C++
// tab-width: 8
// indent-tabs-mode: nil
// c-basic-offset: 8
// End:
|