Beispiel #1
0
void testTiming(Blackbox & A)
{
	typedef typename Blackbox::MatrixDomain Dom;
	typedef typename Dom::Block Block;

	Dom MD = A.domain();
	size_t m = A.rowdim(), n = A.coldim();
	size_t k = (m + n)/2;

	LinBox::UserTimer timer;

	Block B(n,k), C(m,k), D(k,m), E(k,n), F(k,k);
	MD.random(B); MD.random(D);

	vector<typename Dom::Element> v1, v2(m);
	typename Dom::RandIter r(MD);
	typename Dom::Element x;
	for(size_t i = 0; i != n; ++i){
		r.random(x);
		v1.push_back(x);
	}


	//Tests:
	cout << "Timing tests:" << endl << endl;

	timer.clear(); timer.start();
	for(size_t j = 0; j != m; ++j) A.apply(v2,v1);
	timer.stop();
	cout << "apply using vectors time: " << timer << endl;

	timer.clear(); timer.start();
	A.applyTranspose(C,B);
	timer.stop();
	cout << "apply using row addin time: " << timer << endl;

	timer.clear(); timer.start();
	A.unpackingApplyTranspose(C,B);
	timer.stop();
	cout << "apply using block axpy time: " << timer << endl;

	timer.clear(); timer.start();
	MD.mul(F, D, C);
	timer.stop();
	cout << "Matrix Domain mul time: " << timer << endl;

	cout << "End of timing tests" << endl << endl;

} // testTiming
Beispiel #2
0
void largeTest (Blackbox & A)
{
	//Use for large blackboxes
	typedef typename Blackbox::MatrixDomain Dom;
	typedef typename Dom::Block Block;

	Dom MD = A.domain();
	size_t m = A.coldim();
	size_t n = 2000;

	LinBox::UserTimer timer;

	Block B(m,n), C(m,n);
	MD.random(B);

	cout << "Test: " << A.rowdim() << "x" << m << "blackbox multiplied by " << m << "x" << n << "block\nblock size: 2048\n\n";

	timer.clear(); timer.start();
	A.unpackingApply(C,B,2048);
	timer.stop();
	cout << "unpacking apply time: " << timer << endl;
}  //end largeTest
Beispiel #3
0
bool testAssociativity(Blackbox& A)
{
	typedef typename Blackbox::MatrixDomain Dom;
	Dom MD = A.domain();
	size_t m = A.rowdim(), n = A.coldim() - 100;
	size_t k = (m + n)/2;
	typename Dom::Block B(A.field(),k,m), C(A.field(),m,n);
	MD.random(B); MD.random(C);

	typename Dom::Block D(A.field(),m,n), E(A.field(),k,n);

	A.apply(D, C); // D = AC
	MD.mul(E,B,D); // E = B(AC)

	typename Dom::Block F(A.field(),k,m), G(A.field(),k,n);

	A.unpackingApplyTranspose(F,B); // F = BA
	MD.mul(G,F,C); // G = (BA)C
	return MD.areEqual(E,G);


} // testAssociativity
Beispiel #4
0
static bool testTransposeBlackbox(Blackbox & A)
{
	typedef typename Blackbox::Field Field;
	commentator().start ("Testing Transpose", "testTranspose", 1);

	Transpose<Blackbox> B(A);

	bool ret = true, ret1;

	size_t m = A.rowdim(), n = A.coldim();
	const Field & F = A.field();
	VectorDomain<Field> VD (F);
	BlasVector<Field> x(F,n), y(F,m), z(F,n), w(F,m);

	VD.random(x);
	A.apply(y, x);
	B.applyTranspose(w, x);
	ret1 = VD.areEqual(y, w);
	if (not ret1) commentator().report() << "A and B^T disagree, FAIL" << std::endl;
	ret = ret and ret1;

	VD.random(y);
	A.applyTranspose(x, y);
	B.apply(z, y);
	ret1 = VD.areEqual(x, z);
	if (not ret1) commentator().report() << "A^T and B disagree, FAIL" << std::endl;
	ret = ret and ret1;

	ret1 = testBlackboxNoRW(B);
	if (not ret1) commentator().report() << "testBlackbox A^T FAIL" << std::endl;
	ret = ret and ret1;

	commentator().stop (MSG_STATUS (ret), (const char *) 0, "testTranspose");

	return ret;
}
Beispiel #5
0
int main (int argc, char **argv)
{
	//     commentator().setMaxDetailLevel (-1);
	//     commentator().setMaxDepth (-1);
	//     commentator().setReportStream (std::cerr);


	if (argc < 2 || argc > 4) {
		std::cerr << "Usage: omp_smithvalence <matrix-file-in-supported-format> [-ata|-aat|valence] [coprime]" << std::endl;
        std::cerr << "       Optional parameters valence and coprime are integers." << std::endl;
        std::cerr << "       Prime factors of valence will be used for local computation." << std::endl;
        std::cerr << "       coprime will be used for overall rank computation." << std::endl;
		return -1;
	}

	std::ifstream input (argv[1]);
	if (!input) { std::cerr << "Error opening matrix file " << argv[1] << std::endl; return -1; }

	Givaro::ZRing<Integer> ZZ;
	MatrixStream< Givaro::ZRing<Integer> > ms( ZZ, input );
	typedef SparseMatrix<Givaro::ZRing<Integer>>  Blackbox;
	Blackbox A (ms);
	input.close();

	std::cout << "A is " << A.rowdim() << " by " << A.coldim() << std::endl;

	Givaro::ZRing<Integer>::Element val_A;

	LinBox::Timer chrono; chrono.start();
	if (argc >= 3) {
		Transpose<Blackbox> T(&A);
		if (strcmp(argv[2],"-ata") == 0) {
			Compose< Transpose<Blackbox>, Blackbox > C (&T, &A);
			std::cout << "A^T A is " << C.rowdim() << " by " << C.coldim() << std::endl;
			valence(val_A, C);
		}
		else if (strcmp(argv[2],"-aat") == 0) {
			Compose< Blackbox, Transpose<Blackbox> > C (&A, &T);
			std::cout << "A A^T is " << C.rowdim() << " by " << C.coldim() << std::endl;
			valence(val_A, C);
		}
		else {
			std::cout << "Suppose primes are contained in " << argv[2] << std::endl;
			val_A = LinBox::Integer(argv[2]);
		}
	}
	else {
		if (A.rowdim() != A.coldim()) {
			std::cerr << "Valence works only on square matrices, try either to change the dimension in the matrix file, or to compute the valence of A A^T or A^T A, via the -aat or -ata options."  << std::endl;
			exit(0);
		}
		else
			valence (val_A, A);
	}

	std::cout << "Valence is " << val_A << std::endl;

	std::vector<Givaro::Integer> Moduli;
	std::vector<size_t> exponents;
	Givaro::IntFactorDom<> FTD;

	typedef std::pair<Givaro::Integer,unsigned long> PairIntRk;
	std::vector< PairIntRk > smith;


Givaro::Integer coprimeV=2;
	if (argc >= 4) {
		coprimeV =Givaro::Integer(argv[3]);
	}
	while ( gcd(val_A,coprimeV) > 1 ) {
		FTD.nextprimein(coprimeV);
	}

	if (argc >= 4) {
		std::cout << "Suppose " << argv[3] << " is coprime with Smith form" << std::endl;
	}

	std::cout << "Integer rank: " << std::endl;

	unsigned long coprimeR; LRank(coprimeR, argv[1], coprimeV);
	smith.push_back(PairIntRk(coprimeV, coprimeR));
	//         std::cerr << "Rank mod " << coprimeV << " is " << coprimeR << std::endl;

	std::cout << "Some factors (50000 factoring loop bound): ";
	FTD.set(Moduli, exponents, val_A, 50000);
	std::vector<size_t>::const_iterator eit=exponents.begin();
	for(std::vector<Givaro::Integer>::const_iterator mit=Moduli.begin();
	    mit != Moduli.end(); ++mit,++eit)
		std::cout << *mit << '^' << *eit << ' ';
	std::cout << std::endl;

	std::vector<Givaro::Integer> SmithDiagonal(coprimeR,Givaro::Integer(1));

	std::cout << "num procs: " << omp_get_num_procs() << std::endl;
	std::cout << "max threads: " << omp_get_max_threads() << std::endl;
#pragma omp parallel for shared(SmithDiagonal, Moduli, coprimeR)
	for(size_t j=0; j<Moduli.size(); ++j) {
		unsigned long r; LRank(r, argv[1], Moduli[j]);
		std::cerr << "Rank mod " << Moduli[j] << " is " << r << " on thread: " << omp_get_thread_num() << std::endl;
		smith.push_back(PairIntRk( Moduli[j], r));
		for(size_t i=r; i < coprimeR; ++i)
			SmithDiagonal[i] *= Moduli[j];
	}


	/*
	   for(std::vector<Givaro::Integer>::const_iterator mit=Moduli.begin();
	   mit != Moduli.end(); ++mit) {
	   unsigned long r; LRank(r, argv[1], *mit);
	   std::cerr << "Rank mod " << *mit << " is " << r << std::endl;
	   smith.push_back(PairIntRk(*mit, r));
	   for(size_t i=r; i < coprimeR; ++i)
	   SmithDiagonal[i] *= *mit;
	   }
	   */

	eit=exponents.begin();
	std::vector<PairIntRk>::const_iterator sit=smith.begin();
	for( ++sit; sit != smith.end(); ++sit, ++eit) {
            if (sit->second != coprimeR) {
                std::vector<size_t> ranks;
                ranks.push_back(sit->second);
                size_t effexp;
                if (*eit > 1) {
                    if (sit->first == 2)
                        PRankPowerOfTwo(ranks, effexp, argv[1], *eit, coprimeR);
                    else
                        PRank(ranks, effexp, argv[1], sit->first, *eit, coprimeR);
                }
                else {
			// if (sit->first == 2)
			PRank(ranks, effexp, argv[1], sit->first, 2, coprimeR);
			// else
			// PRank(ranks, effexp, argv[1], sit->first, 2, coprimeR);
		}
                if (ranks.size() == 1) ranks.push_back(coprimeR);

                if (effexp < *eit) {
                    for(size_t expo = effexp<<1; ranks.back() < coprimeR; expo<<=1) {
                        if (sit->first == 2)
                            PRankIntegerPowerOfTwo(ranks, argv[1], expo, coprimeR);
                        else
                            PRankInteger(ranks, argv[1], sit->first, expo, coprimeR);
                    }
                } else {

                    for(size_t expo = (*eit)<<1; ranks.back() < coprimeR; expo<<=1) {
                        if (sit->first == 2)
                            PRankPowerOfTwo(ranks, effexp, argv[1], expo, coprimeR);
                        else
                            PRank(ranks, effexp, argv[1], sit->first, expo, coprimeR);
                        if (ranks.size() < expo) {
                            std::cerr << "It seems we need a larger prime power, it will take longer ..." << std::endl;
                                // break;
                            if (sit->first == 2)
                                PRankIntegerPowerOfTwo(ranks, argv[1], expo, coprimeR);
                            else
                                PRankInteger(ranks, argv[1], sit->first, expo, coprimeR);
                        }
                    }
                }

                std::vector<size_t>::const_iterator rit=ranks.begin();
// 			unsigned long modrank = *rit;
                for(++rit; rit!= ranks.end(); ++rit) {
                    if ((*rit)>= coprimeR) break;
                    for(size_t i=(*rit); i < coprimeR; ++i)
                        SmithDiagonal[i] *= sit->first;
// 				modrank = *rit;
                }
            }
	}

Givaro::Integer si=1;
	size_t num=0;
	for( std::vector<Givaro::Integer>::const_iterator dit=SmithDiagonal.begin();
	     dit != SmithDiagonal.end(); ++dit) {
		if (*dit == si) ++num;
		else {
			std::cerr << '[' << si << ',' << num << "] ";
			num=1;
			si = *dit;
		}
	}
	std::cerr << '[' << si << ',' << num << "] " << std::endl;
	chrono.stop();
	std::cerr << chrono << std::endl;


	return 0;
}
Beispiel #6
0
bool testQLUP(const Field &F, size_t n, unsigned int iterations, int rseed, double sparsity = 0.05)
{
	bool res = true;

	commentator().start ("Testing Sparse elimination qlup", "testQLUP", iterations);

	size_t Ni = n;
	size_t Nj = n;
	integer card; F.cardinality(card);
	typename Field::RandIter generator (F,card,rseed);
	RandStream stream (F, generator, sparsity, n, n);

	for (size_t i = 0; i < iterations; ++i) {
		commentator().startIteration ((unsigned)i);


		stream.reset();

		Blackbox A (F, stream);

		std::ostream & report = commentator().report (Commentator::LEVEL_UNIMPORTANT, INTERNAL_DESCRIPTION);

		F.write( report ) << endl;
		A.write( report,Tag::FileFormat::Maple ) << endl;

		DenseVector<Field> u(F,Nj), v(F,Ni), w1(F,Nj), w2(F,Ni), w3(F,Ni), w(F,Ni);
		for(auto it=u.begin();it!=u.end();++it)
			generator.random (*it);


		A.apply(v,u);


		unsigned long rank;

		Method::SparseElimination SE;
		SE.strategy(Specifier::PIVOT_LINEAR);
		GaussDomain<Field> GD ( F );
		typename Field::Element determinant;
		Blackbox L(F, A.rowdim(), A.coldim());
		Permutation<Field> Q((int)A.rowdim(),F);
		Permutation<Field> P((int)A.coldim(),F);

		GD.QLUPin(rank, determinant,
			  Q, L, A, P,
			  A.rowdim(), A.coldim() );

		Q.apply(w, L.apply(w3, A.apply(w2, P.apply(w1,u) ) ) );

		bool error = false;
		auto itv=v.begin();
		auto itw=w.begin();
		for( ; itw!=w.end();++itw,++itv) {
			if (! F.areEqual(*itw,*itv) ) {
				error = true;
			}
		}

		if (error) {
			res = false;

			report << "ERROR : matrix(" << u.size() << ",1,[";
			for(auto itu=u.begin(); itu!=u.end();++itu)
				report << *itu << ',';
			report << "]);\n[";
			for(auto itv2=v.begin(); itv2!=v.end();++itv2)
				report << *itv2 << ' ';
			report << "]  !=  [";
			for(auto itw2=w.begin(); itw2!=w.end();++itw2)
				report << *itw2 << ' ';
			report << "]" << std::endl;


			report << "w1: [";
			for(auto itw2=w1.begin(); itw2!=w1.end();++itw2)
				report << *itw2 << ' ';
			report << "]" << std::endl;
			report << "w2: [";
			for(auto itw2=w2.begin(); itw2!=w2.end();++itw2)
				report << *itw2 << ' ';
			report << "]" << std::endl;
			report << "w3: [";
			for(auto itw2=w3.begin(); itw2!=w3.end();++itw2)
				report << *itw2 << ' ';
			report << "]" << std::endl;
		}

		commentator().stop ("done");
		commentator().progress ();
	}

	commentator().stop (MSG_STATUS (res), (const char *) 0, "testQLUP");

	return res;
}
Beispiel #7
0
bool testQLUPnullspace(const Field &F, size_t n, unsigned int iterations, int rseed, double sparsity = 0.05)
{
	bool res = true;

	commentator().start ("Testing Sparse elimination qlup nullspacebasis", "testQLUPnullspace", iterations);

	size_t Ni = n;
	size_t Nj = n;
	integer card; F.cardinality(card);
	typename Field::RandIter generator (F,card,rseed);
	RandStream stream (F, generator, sparsity, n, n, rseed);

	for (size_t i = 0; i < iterations; ++i) {
		commentator().startIteration ((unsigned)i);

		stream.reset();
		Blackbox A (F, stream);

		std::ostream & report = commentator().report (Commentator::LEVEL_UNIMPORTANT, INTERNAL_DESCRIPTION);

		F.write( report ) << endl;
		A.write( report, Tag::FileFormat::Maple ) << endl;


		Method::SparseElimination SE;
		SE.strategy(Specifier::PIVOT_LINEAR);
		GaussDomain<Field> GD ( F );

		Blackbox CopyA ( A );
		Blackbox X(F, A.coldim(), A.coldim() );

		GD.nullspacebasisin(X, CopyA );

		size_t nullity = X.coldim();

		DenseVector<Field> u(F,nullity);
		for(auto it=u.begin();it!=u.end();++it)
			generator.random (*it);
		DenseVector<Field> v(F,Nj);
		X.apply(v,u);
		report << "Random combination of the rows of the NullSpace basis" << std::endl;

		DenseVector<Field> w(F,Ni);
		A.apply(w, v);

		VectorDomain<Field> VD(F);

		if (! VD.isZero(w)) {
			res=false;
			A.write( report, Tag::FileFormat::Maple ) << endl;

			report << "ERROR u: matrix(" << u.size() << ",1,[";
			for(auto itu=u.begin(); itu!=u.end();++itu)
				report << *itu << ',';
			report << "]);\n[";
			report << "ERROR v: matrix(" << v.size() << ",1,[";
			for(auto itu=v.begin(); itu!=v.end();++itu)
				report << *itu << ',';
			report << "]);\n[";
			for(auto itv=w.begin(); itv!=w.end();++itv)
				report << *itv << ' ';
			report << "]  !=  0" << std::endl;

		}

		commentator().stop ("done");
		commentator().progress ();
	}

	commentator().stop (MSG_STATUS (res), (const char *) 0, "testQLUPnullspace");

	return res;
}