1c4762a1bSJed Brown 24905a7bcSToby Isaac static char help[] = "Tests LU, Cholesky, and QR factorization and MatMatSolve() for a sequential dense matrix. \n\ 3c4762a1bSJed Brown For MATSEQDENSE matrix, the factorization is just a thin wrapper to LAPACK. \n\ 4c4762a1bSJed Brown For MATSEQDENSECUDA, it uses cusolverDn routines \n\n"; 5c4762a1bSJed Brown 6c4762a1bSJed Brown #include <petscmat.h> 7c4762a1bSJed Brown 84905a7bcSToby Isaac static PetscErrorCode createMatsAndVecs(PetscInt m, PetscInt n, PetscInt nrhs, PetscBool full, Mat *_mat, Mat *_RHS, Mat *_SOLU, Vec *_x, Vec *_y, Vec *_b) 9c4762a1bSJed Brown { 10c4762a1bSJed Brown PetscRandom rand; 114905a7bcSToby Isaac Mat mat,RHS,SOLU; 124905a7bcSToby Isaac PetscInt rstart, rend; 134905a7bcSToby Isaac PetscInt cstart, cend; 144905a7bcSToby Isaac PetscScalar value = 1.0; 154905a7bcSToby Isaac Vec x, y, b; 16c4762a1bSJed Brown 174905a7bcSToby Isaac PetscFunctionBegin; 18c4762a1bSJed Brown /* create multiple vectors RHS and SOLU */ 19*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&RHS)); 20*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(RHS,PETSC_DECIDE,PETSC_DECIDE,m,nrhs)); 21*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(RHS,MATDENSE)); 22*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOptionsPrefix(RHS,"rhs_")); 23*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(RHS)); 24*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqDenseSetPreallocation(RHS,NULL)); 25c4762a1bSJed Brown 26*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rand)); 27*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rand)); 28*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(RHS,rand)); 29c4762a1bSJed Brown 304905a7bcSToby Isaac if (m == n) { 31*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(RHS,MAT_DO_NOT_COPY_VALUES,&SOLU)); 324905a7bcSToby Isaac } else { 33*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&SOLU)); 34*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(SOLU,PETSC_DECIDE,PETSC_DECIDE,n,nrhs)); 35*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(SOLU,MATDENSE)); 36*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqDenseSetPreallocation(SOLU,NULL)); 374905a7bcSToby Isaac } 38*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(SOLU,rand)); 39c4762a1bSJed Brown 40c4762a1bSJed Brown /* create matrix */ 41*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&mat)); 42*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(mat,PETSC_DECIDE,PETSC_DECIDE,m,n)); 43*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(mat,MATDENSE)); 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(mat)); 45*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(mat)); 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(mat,&rstart,&rend)); 47*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRangeColumn(mat,&cstart,&cend)); 48c4762a1bSJed Brown if (!full) { 494905a7bcSToby Isaac for (PetscInt i=rstart; i<rend; i++) { 504905a7bcSToby Isaac if (m == n) { 51c4762a1bSJed Brown value = (PetscReal)i+1; 52*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(mat,1,&i,1,&i,&value,INSERT_VALUES)); 534905a7bcSToby Isaac } else { 544905a7bcSToby Isaac for (PetscInt j = cstart; j < cend; j++) { 554905a7bcSToby Isaac value = ((PetscScalar)i+1.)/(PetscSqr(i - j) + 1.); 56*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(mat,1,&i,1,&j,&value,INSERT_VALUES)); 574905a7bcSToby Isaac } 584905a7bcSToby Isaac } 59c4762a1bSJed Brown } 60*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY)); 61*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY)); 62c4762a1bSJed Brown } else { 63*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetRandom(mat,rand)); 644905a7bcSToby Isaac if (m == n) { 65c4762a1bSJed Brown Mat T; 66c4762a1bSJed Brown 67*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatTransposeMult(mat,mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T)); 68*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&mat)); 69c4762a1bSJed Brown mat = T; 70c4762a1bSJed Brown } 714905a7bcSToby Isaac } 72c4762a1bSJed Brown 73c4762a1bSJed Brown /* create single vectors */ 74*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(mat,&x,&b)); 75*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&y)); 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(x,value)); 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rand)); 784905a7bcSToby Isaac *_mat = mat; 794905a7bcSToby Isaac *_RHS = RHS; 804905a7bcSToby Isaac *_SOLU = SOLU; 814905a7bcSToby Isaac *_x = x; 824905a7bcSToby Isaac *_y = y; 834905a7bcSToby Isaac *_b = b; 844905a7bcSToby Isaac PetscFunctionReturn(0); 854905a7bcSToby Isaac } 864905a7bcSToby Isaac 874905a7bcSToby Isaac int main(int argc,char **argv) 884905a7bcSToby Isaac { 894905a7bcSToby Isaac Mat mat,F,RHS,SOLU; 904905a7bcSToby Isaac MatInfo info; 914905a7bcSToby Isaac PetscErrorCode ierr; 924905a7bcSToby Isaac PetscInt m = 15, n = 10,i,j,nrhs=2; 934905a7bcSToby Isaac Vec x,y,b,ytmp; 944905a7bcSToby Isaac IS perm; 954905a7bcSToby Isaac PetscReal norm,tol=PETSC_SMALL; 964905a7bcSToby Isaac PetscMPIInt size; 974905a7bcSToby Isaac char solver[64]; 987275615fSToby Isaac PetscBool inplace,full = PETSC_FALSE, ldl = PETSC_TRUE, qr = PETSC_TRUE; 994905a7bcSToby Isaac 1004905a7bcSToby Isaac ierr = PetscInitialize(&argc,&argv,(char*) 0,help);if (ierr) return ierr; 101*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 1022c71b3e2SJacob Faibussowitsch PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!"); 103*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrcpy(solver,"petsc")); 104*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 105*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 106*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nrhs",&nrhs,NULL)); 107*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-ldl",&ldl,NULL)); 108*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-qr",&qr,NULL)); 109*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-full",&full,NULL)); 110*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL,NULL,"-solver_type",solver,sizeof(solver),NULL)); 1114905a7bcSToby Isaac 112*5f80ce2aSJacob Faibussowitsch CHKERRQ(createMatsAndVecs(n, n, nrhs, full, &mat, &RHS, &SOLU, &x, &y, &b)); 113*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(y,&ytmp)); 114c4762a1bSJed Brown 115c4762a1bSJed Brown /* Only SeqDense* support in-place factorizations and NULL permutations */ 116*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectBaseTypeCompare((PetscObject)mat,MATSEQDENSE,&inplace)); 117*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(mat,&i,NULL)); 118*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(mat,&j,NULL)); 119*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(PETSC_COMM_WORLD,i,j,1,&perm)); 120c4762a1bSJed Brown 121*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetInfo(mat,MAT_LOCAL,&info)); 1228fffc762SJacob Faibussowitsch ierr = PetscPrintf(PETSC_COMM_WORLD,"matrix nonzeros = %" PetscInt_FMT ", allocated nonzeros = %" PetscInt_FMT "\n", 123c4762a1bSJed Brown (PetscInt)info.nz_used,(PetscInt)info.nz_allocated);CHKERRQ(ierr); 124*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(mat,x,b)); 125c4762a1bSJed Brown 126c4762a1bSJed Brown /* Cholesky factorization - perm and factinfo are ignored by LAPACK */ 127c4762a1bSJed Brown /* in-place Cholesky */ 128c4762a1bSJed Brown if (inplace) { 129c4762a1bSJed Brown Mat RHS2; 130c4762a1bSJed Brown 131*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(mat,MAT_COPY_VALUES,&F)); 132*5f80ce2aSJacob Faibussowitsch if (!ldl) CHKERRQ(MatSetOption(F,MAT_SPD,PETSC_TRUE)); 133*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCholeskyFactor(F,perm,0)); 134*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolve(F,b,y)); 135*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(y,-1.0,x)); 136*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(y,NORM_2,&norm)); 137c4762a1bSJed Brown if (norm > tol) { 138*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Warning: Norm of error for in-place Cholesky %g\n",(double)norm)); 139c4762a1bSJed Brown } 140c4762a1bSJed Brown 141*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatSolve(F,RHS,SOLU)); 142*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(mat,SOLU,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&RHS2)); 143*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(RHS,-1.0,RHS2,SAME_NONZERO_PATTERN)); 144*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(RHS,NORM_FROBENIUS,&norm)); 145c4762a1bSJed Brown if (norm > tol) { 146*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error: Norm of residual for in-place Cholesky (MatMatSolve) %g\n",(double)norm)); 147c4762a1bSJed Brown } 148*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&F)); 149*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&RHS2)); 150c4762a1bSJed Brown } 151c4762a1bSJed Brown 152c4762a1bSJed Brown /* out-of-place Cholesky */ 153*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetFactor(mat,solver,MAT_FACTOR_CHOLESKY,&F)); 154*5f80ce2aSJacob Faibussowitsch if (!ldl) CHKERRQ(MatSetOption(F,MAT_SPD,PETSC_TRUE)); 155*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCholeskyFactorSymbolic(F,mat,perm,0)); 156*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCholeskyFactorNumeric(F,mat,0)); 157*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolve(F,b,y)); 158*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(y,-1.0,x)); 159*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(y,NORM_2,&norm)); 160c4762a1bSJed Brown if (norm > tol) { 161*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Warning: Norm of error for out-of-place Cholesky %g\n",(double)norm)); 162c4762a1bSJed Brown } 163*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&F)); 164c4762a1bSJed Brown 165c4762a1bSJed Brown /* LU factorization - perms and factinfo are ignored by LAPACK */ 166c4762a1bSJed Brown i = n-1; 167*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroRows(mat,1,&i,-1.0,NULL,NULL)); 168*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(mat,x,b)); 169c4762a1bSJed Brown 170c4762a1bSJed Brown /* in-place LU */ 171c4762a1bSJed Brown if (inplace) { 172c4762a1bSJed Brown Mat RHS2; 173c4762a1bSJed Brown 174*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(mat,MAT_COPY_VALUES,&F)); 175*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatLUFactor(F,perm,perm,0)); 176*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolve(F,b,y)); 177*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(y,-1.0,x)); 178*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(y,NORM_2,&norm)); 179c4762a1bSJed Brown if (norm > tol) { 180*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Warning: Norm of error for in-place LU %g\n",(double)norm)); 181c4762a1bSJed Brown } 182*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatSolve(F,RHS,SOLU)); 183*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(mat,SOLU,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&RHS2)); 184*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(RHS,-1.0,RHS2,SAME_NONZERO_PATTERN)); 185*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(RHS,NORM_FROBENIUS,&norm)); 186c4762a1bSJed Brown if (norm > tol) { 187*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error: Norm of residual for in-place LU (MatMatSolve) %g\n",(double)norm)); 188c4762a1bSJed Brown } 189*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&F)); 190*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&RHS2)); 191c4762a1bSJed Brown } 192c4762a1bSJed Brown 193c4762a1bSJed Brown /* out-of-place LU */ 194*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetFactor(mat,solver,MAT_FACTOR_LU,&F)); 195*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatLUFactorSymbolic(F,mat,perm,perm,0)); 196*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatLUFactorNumeric(F,mat,0)); 197*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolve(F,b,y)); 198*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(y,-1.0,x)); 199*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(y,NORM_2,&norm)); 200c4762a1bSJed Brown if (norm > tol) { 201*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Warning: Norm of error for out-of-place LU %g\n",(double)norm)); 202c4762a1bSJed Brown } 203c4762a1bSJed Brown 204c4762a1bSJed Brown /* free space */ 205*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&perm)); 206*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&F)); 207*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&mat)); 208*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&RHS)); 209*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&SOLU)); 210*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 211*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&b)); 212*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&y)); 213*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&ytmp)); 2144905a7bcSToby Isaac 2157275615fSToby Isaac if (qr) { 2164905a7bcSToby Isaac /* setup rectanglar */ 217*5f80ce2aSJacob Faibussowitsch CHKERRQ(createMatsAndVecs(m, n, nrhs, full, &mat, &RHS, &SOLU, &x, &y, &b)); 218*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(y,&ytmp)); 2194905a7bcSToby Isaac 2204905a7bcSToby Isaac /* QR factorization - perms and factinfo are ignored by LAPACK */ 221*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(mat,x,b)); 2224905a7bcSToby Isaac 2234905a7bcSToby Isaac /* in-place QR */ 2244905a7bcSToby Isaac if (inplace) { 2254905a7bcSToby Isaac Mat SOLU2; 2264905a7bcSToby Isaac 227*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(mat,MAT_COPY_VALUES,&F)); 228*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatQRFactor(F,NULL,0)); 229*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolve(F,b,y)); 230*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(y,-1.0,x)); 231*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(y,NORM_2,&norm)); 2324905a7bcSToby Isaac if (norm > tol) { 233*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Warning: Norm of error for in-place QR %g\n",(double)norm)); 2344905a7bcSToby Isaac } 235*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(mat,SOLU,MAT_REUSE_MATRIX,PETSC_DEFAULT,&RHS)); 236*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(SOLU, MAT_DO_NOT_COPY_VALUES, &SOLU2)); 237*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatSolve(F,RHS,SOLU2)); 238*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(SOLU2,-1.0,SOLU,SAME_NONZERO_PATTERN)); 239*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(SOLU2,NORM_FROBENIUS,&norm)); 2404905a7bcSToby Isaac if (norm > tol) { 241*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Error: Norm of error for in-place QR (MatMatSolve) %g\n",(double)norm)); 2424905a7bcSToby Isaac } 243*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&F)); 244*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&SOLU2)); 2454905a7bcSToby Isaac } 2464905a7bcSToby Isaac 2474905a7bcSToby Isaac /* out-of-place QR */ 248*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetFactor(mat,solver,MAT_FACTOR_QR,&F)); 249*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatQRFactorSymbolic(F,mat,NULL,NULL)); 250*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatQRFactorNumeric(F,mat,NULL)); 251*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolve(F,b,y)); 252*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(y,-1.0,x)); 253*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(y,NORM_2,&norm)); 2544905a7bcSToby Isaac if (norm > tol) { 255*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Warning: Norm of error for out-of-place QR %g\n",(double)norm)); 2564905a7bcSToby Isaac } 2574905a7bcSToby Isaac 2584905a7bcSToby Isaac if (m == n) { 2594905a7bcSToby Isaac /* out-of-place MatSolveTranspose */ 260*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(mat,x,b)); 261*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolveTranspose(F,b,y)); 262*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(y,-1.0,x)); 263*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(y,NORM_2,&norm)); 2644905a7bcSToby Isaac if (norm > tol) { 265*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Warning: Norm of error for out-of-place QR %g\n",(double)norm)); 2664905a7bcSToby Isaac } 2674905a7bcSToby Isaac } 2684905a7bcSToby Isaac 2694905a7bcSToby Isaac /* free space */ 270*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&F)); 271*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&mat)); 272*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&RHS)); 273*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&SOLU)); 274*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 275*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&b)); 276*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&y)); 277*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&ytmp)); 2787275615fSToby Isaac } 279c4762a1bSJed Brown ierr = PetscFinalize(); 280c4762a1bSJed Brown return ierr; 281c4762a1bSJed Brown } 282c4762a1bSJed Brown 283c4762a1bSJed Brown /*TEST 284c4762a1bSJed Brown 285c4762a1bSJed Brown test: 286c4762a1bSJed Brown 287c4762a1bSJed Brown test: 288c4762a1bSJed Brown requires: cuda 289c4762a1bSJed Brown suffix: seqdensecuda 290db0483a4SToby Isaac args: -mat_type seqdensecuda -rhs_mat_type seqdensecuda -ldl 0 -solver_type {{petsc cuda}} 2917275615fSToby Isaac output_file: output/ex1_1.out 2927275615fSToby Isaac 2937275615fSToby Isaac test: 2947275615fSToby Isaac requires: cuda 2957275615fSToby Isaac suffix: seqdensecuda_2 296db0483a4SToby Isaac args: -ldl 0 -solver_type cuda 297c4762a1bSJed Brown output_file: output/ex1_1.out 298c4762a1bSJed Brown 299c4762a1bSJed Brown test: 300c4762a1bSJed Brown requires: cuda 301c4762a1bSJed Brown suffix: seqdensecuda_seqaijcusparse 3027275615fSToby Isaac args: -mat_type seqaijcusparse -rhs_mat_type seqdensecuda -qr 0 303c4762a1bSJed Brown output_file: output/ex1_2.out 304c4762a1bSJed Brown 305c4762a1bSJed Brown test: 306c4762a1bSJed Brown requires: cuda viennacl 307c4762a1bSJed Brown suffix: seqdensecuda_seqaijviennacl 3087275615fSToby Isaac args: -mat_type seqaijviennacl -rhs_mat_type seqdensecuda -qr 0 309c4762a1bSJed Brown output_file: output/ex1_2.out 310c4762a1bSJed Brown 3114905a7bcSToby Isaac test: 3124905a7bcSToby Isaac suffix: 4 3134905a7bcSToby Isaac args: -m 10 -n 10 3144905a7bcSToby Isaac output_file: output/ex1_1.out 3154905a7bcSToby Isaac 316c4762a1bSJed Brown TEST*/ 317