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; 164905a7bcSToby Isaac PetscErrorCode ierr; 17c4762a1bSJed Brown 184905a7bcSToby Isaac PetscFunctionBegin; 19c4762a1bSJed Brown /* create multiple vectors RHS and SOLU */ 20c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&RHS);CHKERRQ(ierr); 214905a7bcSToby Isaac ierr = MatSetSizes(RHS,PETSC_DECIDE,PETSC_DECIDE,m,nrhs);CHKERRQ(ierr); 22c4762a1bSJed Brown ierr = MatSetType(RHS,MATDENSE);CHKERRQ(ierr); 23c4762a1bSJed Brown ierr = MatSetOptionsPrefix(RHS,"rhs_");CHKERRQ(ierr); 24c4762a1bSJed Brown ierr = MatSetFromOptions(RHS);CHKERRQ(ierr); 25c4762a1bSJed Brown ierr = MatSeqDenseSetPreallocation(RHS,NULL);CHKERRQ(ierr); 26c4762a1bSJed Brown 27c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rand);CHKERRQ(ierr); 28c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(rand);CHKERRQ(ierr); 29c4762a1bSJed Brown ierr = MatSetRandom(RHS,rand);CHKERRQ(ierr); 30c4762a1bSJed Brown 314905a7bcSToby Isaac if (m == n) { 32c4762a1bSJed Brown ierr = MatDuplicate(RHS,MAT_DO_NOT_COPY_VALUES,&SOLU);CHKERRQ(ierr); 334905a7bcSToby Isaac } else { 344905a7bcSToby Isaac ierr = MatCreate(PETSC_COMM_WORLD,&SOLU);CHKERRQ(ierr); 354905a7bcSToby Isaac ierr = MatSetSizes(SOLU,PETSC_DECIDE,PETSC_DECIDE,n,nrhs);CHKERRQ(ierr); 364905a7bcSToby Isaac ierr = MatSetType(SOLU,MATDENSE);CHKERRQ(ierr); 374905a7bcSToby Isaac ierr = MatSeqDenseSetPreallocation(SOLU,NULL);CHKERRQ(ierr); 384905a7bcSToby Isaac } 394905a7bcSToby Isaac ierr = MatSetRandom(SOLU,rand);CHKERRQ(ierr); 40c4762a1bSJed Brown 41c4762a1bSJed Brown /* create matrix */ 42c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&mat);CHKERRQ(ierr); 434905a7bcSToby Isaac ierr = MatSetSizes(mat,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr); 44c4762a1bSJed Brown ierr = MatSetType(mat,MATDENSE);CHKERRQ(ierr); 45c4762a1bSJed Brown ierr = MatSetFromOptions(mat);CHKERRQ(ierr); 46c4762a1bSJed Brown ierr = MatSetUp(mat);CHKERRQ(ierr); 47c4762a1bSJed Brown ierr = MatGetOwnershipRange(mat,&rstart,&rend);CHKERRQ(ierr); 484905a7bcSToby Isaac ierr = MatGetOwnershipRangeColumn(mat,&cstart,&cend);CHKERRQ(ierr); 49c4762a1bSJed Brown if (!full) { 504905a7bcSToby Isaac for (PetscInt i=rstart; i<rend; i++) { 514905a7bcSToby Isaac if (m == n) { 52c4762a1bSJed Brown value = (PetscReal)i+1; 53c4762a1bSJed Brown ierr = MatSetValues(mat,1,&i,1,&i,&value,INSERT_VALUES);CHKERRQ(ierr); 544905a7bcSToby Isaac } else { 554905a7bcSToby Isaac for (PetscInt j = cstart; j < cend; j++) { 564905a7bcSToby Isaac value = ((PetscScalar)i+1.)/(PetscSqr(i - j) + 1.); 574905a7bcSToby Isaac ierr = MatSetValues(mat,1,&i,1,&j,&value,INSERT_VALUES);CHKERRQ(ierr); 584905a7bcSToby Isaac } 594905a7bcSToby Isaac } 60c4762a1bSJed Brown } 61c4762a1bSJed Brown ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 62c4762a1bSJed Brown ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 63c4762a1bSJed Brown } else { 644905a7bcSToby Isaac ierr = MatSetRandom(mat,rand);CHKERRQ(ierr); 654905a7bcSToby Isaac if (m == n) { 66c4762a1bSJed Brown Mat T; 67c4762a1bSJed Brown 68c4762a1bSJed Brown ierr = MatMatTransposeMult(mat,mat,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&T);CHKERRQ(ierr); 69c4762a1bSJed Brown ierr = MatDestroy(&mat);CHKERRQ(ierr); 70c4762a1bSJed Brown mat = T; 71c4762a1bSJed Brown } 724905a7bcSToby Isaac } 73c4762a1bSJed Brown 74c4762a1bSJed Brown /* create single vectors */ 75c4762a1bSJed Brown ierr = MatCreateVecs(mat,&x,&b);CHKERRQ(ierr); 76c4762a1bSJed Brown ierr = VecDuplicate(x,&y);CHKERRQ(ierr); 77c4762a1bSJed Brown ierr = VecSet(x,value);CHKERRQ(ierr); 784905a7bcSToby Isaac ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr); 794905a7bcSToby Isaac *_mat = mat; 804905a7bcSToby Isaac *_RHS = RHS; 814905a7bcSToby Isaac *_SOLU = SOLU; 824905a7bcSToby Isaac *_x = x; 834905a7bcSToby Isaac *_y = y; 844905a7bcSToby Isaac *_b = b; 854905a7bcSToby Isaac PetscFunctionReturn(0); 864905a7bcSToby Isaac } 874905a7bcSToby Isaac 884905a7bcSToby Isaac int main(int argc,char **argv) 894905a7bcSToby Isaac { 904905a7bcSToby Isaac Mat mat,F,RHS,SOLU; 914905a7bcSToby Isaac MatInfo info; 924905a7bcSToby Isaac PetscErrorCode ierr; 934905a7bcSToby Isaac PetscInt m = 15, n = 10,i,j,nrhs=2; 944905a7bcSToby Isaac Vec x,y,b,ytmp; 954905a7bcSToby Isaac IS perm; 964905a7bcSToby Isaac PetscReal norm,tol=PETSC_SMALL; 974905a7bcSToby Isaac PetscMPIInt size; 984905a7bcSToby Isaac char solver[64]; 997275615fSToby Isaac PetscBool inplace,full = PETSC_FALSE, ldl = PETSC_TRUE, qr = PETSC_TRUE; 1004905a7bcSToby Isaac 1014905a7bcSToby Isaac ierr = PetscInitialize(&argc,&argv,(char*) 0,help);if (ierr) return ierr; 1024905a7bcSToby Isaac ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 1034905a7bcSToby Isaac if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This is a uniprocessor example only!"); 1044905a7bcSToby Isaac ierr = PetscStrcpy(solver,"petsc");CHKERRQ(ierr); 1054905a7bcSToby Isaac ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr); 1064905a7bcSToby Isaac ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr); 1074905a7bcSToby Isaac ierr = PetscOptionsGetInt(NULL,NULL,"-nrhs",&nrhs,NULL);CHKERRQ(ierr); 1084905a7bcSToby Isaac ierr = PetscOptionsGetBool(NULL,NULL,"-ldl",&ldl,NULL);CHKERRQ(ierr); 1097275615fSToby Isaac ierr = PetscOptionsGetBool(NULL,NULL,"-qr",&qr,NULL);CHKERRQ(ierr); 1104905a7bcSToby Isaac ierr = PetscOptionsGetBool(NULL,NULL,"-full",&full,NULL);CHKERRQ(ierr); 1114905a7bcSToby Isaac ierr = PetscOptionsGetString(NULL,NULL,"-solver_type",solver,sizeof(solver),NULL);CHKERRQ(ierr); 1124905a7bcSToby Isaac 1134905a7bcSToby Isaac ierr = createMatsAndVecs(n, n, nrhs, full, &mat, &RHS, &SOLU, &x, &y, &b);CHKERRQ(ierr); 1144905a7bcSToby Isaac ierr = VecDuplicate(y,&ytmp);CHKERRQ(ierr); 115c4762a1bSJed Brown 116c4762a1bSJed Brown /* Only SeqDense* support in-place factorizations and NULL permutations */ 117c4762a1bSJed Brown ierr = PetscObjectBaseTypeCompare((PetscObject)mat,MATSEQDENSE,&inplace);CHKERRQ(ierr); 118c4762a1bSJed Brown ierr = MatGetLocalSize(mat,&i,NULL);CHKERRQ(ierr); 119c4762a1bSJed Brown ierr = MatGetOwnershipRange(mat,&j,NULL);CHKERRQ(ierr); 120c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_WORLD,i,j,1,&perm);CHKERRQ(ierr); 121c4762a1bSJed Brown 122c4762a1bSJed Brown ierr = MatGetInfo(mat,MAT_LOCAL,&info);CHKERRQ(ierr); 123*8fffc762SJacob Faibussowitsch ierr = PetscPrintf(PETSC_COMM_WORLD,"matrix nonzeros = %" PetscInt_FMT ", allocated nonzeros = %" PetscInt_FMT "\n", 124c4762a1bSJed Brown (PetscInt)info.nz_used,(PetscInt)info.nz_allocated);CHKERRQ(ierr); 125c4762a1bSJed Brown ierr = MatMult(mat,x,b);CHKERRQ(ierr); 126c4762a1bSJed Brown 127c4762a1bSJed Brown /* Cholesky factorization - perm and factinfo are ignored by LAPACK */ 128c4762a1bSJed Brown /* in-place Cholesky */ 129c4762a1bSJed Brown if (inplace) { 130c4762a1bSJed Brown Mat RHS2; 131c4762a1bSJed Brown 132c4762a1bSJed Brown ierr = MatDuplicate(mat,MAT_COPY_VALUES,&F);CHKERRQ(ierr); 133c4762a1bSJed Brown if (!ldl) { ierr = MatSetOption(F,MAT_SPD,PETSC_TRUE);CHKERRQ(ierr); } 134c4762a1bSJed Brown ierr = MatCholeskyFactor(F,perm,0);CHKERRQ(ierr); 135c4762a1bSJed Brown ierr = MatSolve(F,b,y);CHKERRQ(ierr); 136c4762a1bSJed Brown ierr = VecAXPY(y,-1.0,x);CHKERRQ(ierr); 137c4762a1bSJed Brown ierr = VecNorm(y,NORM_2,&norm);CHKERRQ(ierr); 138c4762a1bSJed Brown if (norm > tol) { 139c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Warning: Norm of error for in-place Cholesky %g\n",(double)norm);CHKERRQ(ierr); 140c4762a1bSJed Brown } 141c4762a1bSJed Brown 142c4762a1bSJed Brown ierr = MatMatSolve(F,RHS,SOLU);CHKERRQ(ierr); 143c4762a1bSJed Brown ierr = MatMatMult(mat,SOLU,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&RHS2);CHKERRQ(ierr); 144c4762a1bSJed Brown ierr = MatAXPY(RHS,-1.0,RHS2,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 145c4762a1bSJed Brown ierr = MatNorm(RHS,NORM_FROBENIUS,&norm);CHKERRQ(ierr); 146c4762a1bSJed Brown if (norm > tol) { 147c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Error: Norm of residual for in-place Cholesky (MatMatSolve) %g\n",(double)norm);CHKERRQ(ierr); 148c4762a1bSJed Brown } 149c4762a1bSJed Brown ierr = MatDestroy(&F);CHKERRQ(ierr); 150c4762a1bSJed Brown ierr = MatDestroy(&RHS2);CHKERRQ(ierr); 151c4762a1bSJed Brown } 152c4762a1bSJed Brown 153c4762a1bSJed Brown /* out-of-place Cholesky */ 154c4762a1bSJed Brown ierr = MatGetFactor(mat,solver,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr); 155c4762a1bSJed Brown if (!ldl) { ierr = MatSetOption(F,MAT_SPD,PETSC_TRUE);CHKERRQ(ierr); } 156c4762a1bSJed Brown ierr = MatCholeskyFactorSymbolic(F,mat,perm,0);CHKERRQ(ierr); 157c4762a1bSJed Brown ierr = MatCholeskyFactorNumeric(F,mat,0);CHKERRQ(ierr); 158c4762a1bSJed Brown ierr = MatSolve(F,b,y);CHKERRQ(ierr); 159c4762a1bSJed Brown ierr = VecAXPY(y,-1.0,x);CHKERRQ(ierr); 160c4762a1bSJed Brown ierr = VecNorm(y,NORM_2,&norm);CHKERRQ(ierr); 161c4762a1bSJed Brown if (norm > tol) { 162c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Warning: Norm of error for out-of-place Cholesky %g\n",(double)norm);CHKERRQ(ierr); 163c4762a1bSJed Brown } 164c4762a1bSJed Brown ierr = MatDestroy(&F);CHKERRQ(ierr); 165c4762a1bSJed Brown 166c4762a1bSJed Brown /* LU factorization - perms and factinfo are ignored by LAPACK */ 167c4762a1bSJed Brown i = n-1; 168c4762a1bSJed Brown ierr = MatZeroRows(mat,1,&i,-1.0,NULL,NULL);CHKERRQ(ierr); 169c4762a1bSJed Brown ierr = MatMult(mat,x,b);CHKERRQ(ierr); 170c4762a1bSJed Brown 171c4762a1bSJed Brown /* in-place LU */ 172c4762a1bSJed Brown if (inplace) { 173c4762a1bSJed Brown Mat RHS2; 174c4762a1bSJed Brown 175c4762a1bSJed Brown ierr = MatDuplicate(mat,MAT_COPY_VALUES,&F);CHKERRQ(ierr); 176c4762a1bSJed Brown ierr = MatLUFactor(F,perm,perm,0);CHKERRQ(ierr); 177c4762a1bSJed Brown ierr = MatSolve(F,b,y);CHKERRQ(ierr); 178c4762a1bSJed Brown ierr = VecAXPY(y,-1.0,x);CHKERRQ(ierr); 179c4762a1bSJed Brown ierr = VecNorm(y,NORM_2,&norm);CHKERRQ(ierr); 180c4762a1bSJed Brown if (norm > tol) { 181c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Warning: Norm of error for in-place LU %g\n",(double)norm);CHKERRQ(ierr); 182c4762a1bSJed Brown } 183c4762a1bSJed Brown ierr = MatMatSolve(F,RHS,SOLU);CHKERRQ(ierr); 184c4762a1bSJed Brown ierr = MatMatMult(mat,SOLU,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&RHS2);CHKERRQ(ierr); 185c4762a1bSJed Brown ierr = MatAXPY(RHS,-1.0,RHS2,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 186c4762a1bSJed Brown ierr = MatNorm(RHS,NORM_FROBENIUS,&norm);CHKERRQ(ierr); 187c4762a1bSJed Brown if (norm > tol) { 188c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Error: Norm of residual for in-place LU (MatMatSolve) %g\n",(double)norm);CHKERRQ(ierr); 189c4762a1bSJed Brown } 190c4762a1bSJed Brown ierr = MatDestroy(&F);CHKERRQ(ierr); 191c4762a1bSJed Brown ierr = MatDestroy(&RHS2);CHKERRQ(ierr); 192c4762a1bSJed Brown } 193c4762a1bSJed Brown 194c4762a1bSJed Brown /* out-of-place LU */ 195c4762a1bSJed Brown ierr = MatGetFactor(mat,solver,MAT_FACTOR_LU,&F);CHKERRQ(ierr); 196c4762a1bSJed Brown ierr = MatLUFactorSymbolic(F,mat,perm,perm,0);CHKERRQ(ierr); 197c4762a1bSJed Brown ierr = MatLUFactorNumeric(F,mat,0);CHKERRQ(ierr); 198c4762a1bSJed Brown ierr = MatSolve(F,b,y);CHKERRQ(ierr); 199c4762a1bSJed Brown ierr = VecAXPY(y,-1.0,x);CHKERRQ(ierr); 200c4762a1bSJed Brown ierr = VecNorm(y,NORM_2,&norm);CHKERRQ(ierr); 201c4762a1bSJed Brown if (norm > tol) { 202c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Warning: Norm of error for out-of-place LU %g\n",(double)norm);CHKERRQ(ierr); 203c4762a1bSJed Brown } 204c4762a1bSJed Brown 205c4762a1bSJed Brown /* free space */ 206c4762a1bSJed Brown ierr = ISDestroy(&perm);CHKERRQ(ierr); 207c4762a1bSJed Brown ierr = MatDestroy(&F);CHKERRQ(ierr); 208c4762a1bSJed Brown ierr = MatDestroy(&mat);CHKERRQ(ierr); 209c4762a1bSJed Brown ierr = MatDestroy(&RHS);CHKERRQ(ierr); 210c4762a1bSJed Brown ierr = MatDestroy(&SOLU);CHKERRQ(ierr); 2114905a7bcSToby Isaac ierr = VecDestroy(&x);CHKERRQ(ierr); 2124905a7bcSToby Isaac ierr = VecDestroy(&b);CHKERRQ(ierr); 2134905a7bcSToby Isaac ierr = VecDestroy(&y);CHKERRQ(ierr); 2144905a7bcSToby Isaac ierr = VecDestroy(&ytmp);CHKERRQ(ierr); 2154905a7bcSToby Isaac 2167275615fSToby Isaac if (qr) { 2174905a7bcSToby Isaac /* setup rectanglar */ 2184905a7bcSToby Isaac ierr = createMatsAndVecs(m, n, nrhs, full, &mat, &RHS, &SOLU, &x, &y, &b);CHKERRQ(ierr); 2194905a7bcSToby Isaac ierr = VecDuplicate(y,&ytmp);CHKERRQ(ierr); 2204905a7bcSToby Isaac 2214905a7bcSToby Isaac /* QR factorization - perms and factinfo are ignored by LAPACK */ 2224905a7bcSToby Isaac ierr = MatMult(mat,x,b);CHKERRQ(ierr); 2234905a7bcSToby Isaac 2244905a7bcSToby Isaac /* in-place QR */ 2254905a7bcSToby Isaac if (inplace) { 2264905a7bcSToby Isaac Mat SOLU2; 2274905a7bcSToby Isaac 2284905a7bcSToby Isaac ierr = MatDuplicate(mat,MAT_COPY_VALUES,&F);CHKERRQ(ierr); 2294905a7bcSToby Isaac ierr = MatQRFactor(F,NULL,0);CHKERRQ(ierr); 2304905a7bcSToby Isaac ierr = MatSolve(F,b,y);CHKERRQ(ierr); 2314905a7bcSToby Isaac ierr = VecAXPY(y,-1.0,x);CHKERRQ(ierr); 2324905a7bcSToby Isaac ierr = VecNorm(y,NORM_2,&norm);CHKERRQ(ierr); 2334905a7bcSToby Isaac if (norm > tol) { 2344905a7bcSToby Isaac ierr = PetscPrintf(PETSC_COMM_WORLD,"Warning: Norm of error for in-place QR %g\n",(double)norm);CHKERRQ(ierr); 2354905a7bcSToby Isaac } 2364905a7bcSToby Isaac ierr = MatMatMult(mat,SOLU,MAT_REUSE_MATRIX,PETSC_DEFAULT,&RHS);CHKERRQ(ierr); 2374905a7bcSToby Isaac ierr = MatDuplicate(SOLU, MAT_DO_NOT_COPY_VALUES, &SOLU2);CHKERRQ(ierr); 2384905a7bcSToby Isaac ierr = MatMatSolve(F,RHS,SOLU2);CHKERRQ(ierr); 2394905a7bcSToby Isaac ierr = MatAXPY(SOLU2,-1.0,SOLU,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 2404905a7bcSToby Isaac ierr = MatNorm(SOLU2,NORM_FROBENIUS,&norm);CHKERRQ(ierr); 2414905a7bcSToby Isaac if (norm > tol) { 2424905a7bcSToby Isaac ierr = PetscPrintf(PETSC_COMM_WORLD,"Error: Norm of error for in-place QR (MatMatSolve) %g\n",(double)norm);CHKERRQ(ierr); 2434905a7bcSToby Isaac } 2444905a7bcSToby Isaac ierr = MatDestroy(&F);CHKERRQ(ierr); 2454905a7bcSToby Isaac ierr = MatDestroy(&SOLU2);CHKERRQ(ierr); 2464905a7bcSToby Isaac } 2474905a7bcSToby Isaac 2484905a7bcSToby Isaac /* out-of-place QR */ 2494905a7bcSToby Isaac ierr = MatGetFactor(mat,solver,MAT_FACTOR_QR,&F);CHKERRQ(ierr); 2504905a7bcSToby Isaac ierr = MatQRFactorSymbolic(F,mat,NULL,NULL);CHKERRQ(ierr); 2514905a7bcSToby Isaac ierr = MatQRFactorNumeric(F,mat,NULL);CHKERRQ(ierr); 2524905a7bcSToby Isaac ierr = MatSolve(F,b,y);CHKERRQ(ierr); 2534905a7bcSToby Isaac ierr = VecAXPY(y,-1.0,x);CHKERRQ(ierr); 2544905a7bcSToby Isaac ierr = VecNorm(y,NORM_2,&norm);CHKERRQ(ierr); 2554905a7bcSToby Isaac if (norm > tol) { 2564905a7bcSToby Isaac ierr = PetscPrintf(PETSC_COMM_WORLD,"Warning: Norm of error for out-of-place QR %g\n",(double)norm);CHKERRQ(ierr); 2574905a7bcSToby Isaac } 2584905a7bcSToby Isaac 2594905a7bcSToby Isaac if (m == n) { 2604905a7bcSToby Isaac /* out-of-place MatSolveTranspose */ 2614905a7bcSToby Isaac ierr = MatMultTranspose(mat,x,b);CHKERRQ(ierr); 2624905a7bcSToby Isaac ierr = MatSolveTranspose(F,b,y);CHKERRQ(ierr); 2634905a7bcSToby Isaac ierr = VecAXPY(y,-1.0,x);CHKERRQ(ierr); 2644905a7bcSToby Isaac ierr = VecNorm(y,NORM_2,&norm);CHKERRQ(ierr); 2654905a7bcSToby Isaac if (norm > tol) { 2664905a7bcSToby Isaac ierr = PetscPrintf(PETSC_COMM_WORLD,"Warning: Norm of error for out-of-place QR %g\n",(double)norm);CHKERRQ(ierr); 2674905a7bcSToby Isaac } 2684905a7bcSToby Isaac } 2694905a7bcSToby Isaac 2704905a7bcSToby Isaac /* free space */ 2714905a7bcSToby Isaac ierr = MatDestroy(&F);CHKERRQ(ierr); 2724905a7bcSToby Isaac ierr = MatDestroy(&mat);CHKERRQ(ierr); 2734905a7bcSToby Isaac ierr = MatDestroy(&RHS);CHKERRQ(ierr); 2744905a7bcSToby Isaac ierr = MatDestroy(&SOLU);CHKERRQ(ierr); 275c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 276c4762a1bSJed Brown ierr = VecDestroy(&b);CHKERRQ(ierr); 277c4762a1bSJed Brown ierr = VecDestroy(&y);CHKERRQ(ierr); 278c4762a1bSJed Brown ierr = VecDestroy(&ytmp);CHKERRQ(ierr); 2797275615fSToby Isaac } 280c4762a1bSJed Brown ierr = PetscFinalize(); 281c4762a1bSJed Brown return ierr; 282c4762a1bSJed Brown } 283c4762a1bSJed Brown 284c4762a1bSJed Brown /*TEST 285c4762a1bSJed Brown 286c4762a1bSJed Brown test: 287c4762a1bSJed Brown 288c4762a1bSJed Brown test: 289c4762a1bSJed Brown requires: cuda 290c4762a1bSJed Brown suffix: seqdensecuda 291db0483a4SToby Isaac args: -mat_type seqdensecuda -rhs_mat_type seqdensecuda -ldl 0 -solver_type {{petsc cuda}} 2927275615fSToby Isaac output_file: output/ex1_1.out 2937275615fSToby Isaac 2947275615fSToby Isaac test: 2957275615fSToby Isaac requires: cuda 2967275615fSToby Isaac suffix: seqdensecuda_2 297db0483a4SToby Isaac args: -ldl 0 -solver_type cuda 298c4762a1bSJed Brown output_file: output/ex1_1.out 299c4762a1bSJed Brown 300c4762a1bSJed Brown test: 301c4762a1bSJed Brown requires: cuda 302c4762a1bSJed Brown suffix: seqdensecuda_seqaijcusparse 3037275615fSToby Isaac args: -mat_type seqaijcusparse -rhs_mat_type seqdensecuda -qr 0 304c4762a1bSJed Brown output_file: output/ex1_2.out 305c4762a1bSJed Brown 306c4762a1bSJed Brown test: 307c4762a1bSJed Brown requires: cuda viennacl 308c4762a1bSJed Brown suffix: seqdensecuda_seqaijviennacl 3097275615fSToby Isaac args: -mat_type seqaijviennacl -rhs_mat_type seqdensecuda -qr 0 310c4762a1bSJed Brown output_file: output/ex1_2.out 311c4762a1bSJed Brown 3124905a7bcSToby Isaac test: 3134905a7bcSToby Isaac suffix: 4 3144905a7bcSToby Isaac args: -m 10 -n 10 3154905a7bcSToby Isaac output_file: output/ex1_1.out 3164905a7bcSToby Isaac 317c4762a1bSJed Brown TEST*/ 318