1c4762a1bSJed Brown static char help[] = "Tests MatSolve(), MatSolveTranspose() and MatMatSolve() with SEQDENSE\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown 5c4762a1bSJed Brown int main(int argc,char **args) 6c4762a1bSJed Brown { 7c4762a1bSJed Brown Mat A,RHS,C,F,X; 8c4762a1bSJed Brown Vec u,x,b; 9c4762a1bSJed Brown PetscMPIInt size; 10c4762a1bSJed Brown PetscInt m,n,nsolve,nrhs; 11c4762a1bSJed Brown PetscReal norm,tol=PETSC_SQRT_MACHINE_EPSILON; 12c4762a1bSJed Brown PetscRandom rand; 13c4762a1bSJed Brown PetscBool data_provided,herm,symm,hpd; 14c4762a1bSJed Brown MatFactorType ftyp; 15c4762a1bSJed Brown PetscViewer fd; 16c4762a1bSJed Brown char file[PETSC_MAX_PATH_LEN]; 17c4762a1bSJed Brown 189566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 20be096a46SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor test"); 21c4762a1bSJed Brown /* Determine which type of solver we want to test for */ 22c4762a1bSJed Brown herm = PETSC_FALSE; 23c4762a1bSJed Brown symm = PETSC_FALSE; 24c4762a1bSJed Brown hpd = PETSC_FALSE; 259566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-symmetric_solve",&symm,NULL)); 269566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-hermitian_solve",&herm,NULL)); 279566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-hpd_solve",&hpd,NULL)); 28c4762a1bSJed Brown 29c4762a1bSJed Brown /* Determine file from which we read the matrix A */ 30c4762a1bSJed Brown ftyp = MAT_FACTOR_LU; 319566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&data_provided)); 32c4762a1bSJed Brown if (!data_provided) { /* get matrices from PETSc distribution */ 339566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(file,"${PETSC_DIR}/share/petsc/datafiles/matrices/")); 34c4762a1bSJed Brown if (hpd) { 35c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) 369566063dSJacob Faibussowitsch PetscCall(PetscStrcat(file,"hpd-complex-")); 37c4762a1bSJed Brown #else 389566063dSJacob Faibussowitsch PetscCall(PetscStrcat(file,"spd-real-")); 39c4762a1bSJed Brown #endif 40c4762a1bSJed Brown ftyp = MAT_FACTOR_CHOLESKY; 41c4762a1bSJed Brown } else { 42c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX) 439566063dSJacob Faibussowitsch PetscCall(PetscStrcat(file,"nh-complex-")); 44c4762a1bSJed Brown #else 459566063dSJacob Faibussowitsch PetscCall(PetscStrcat(file,"ns-real-")); 46c4762a1bSJed Brown #endif 47c4762a1bSJed Brown } 48c4762a1bSJed Brown #if defined(PETSC_USE_64BIT_INDICES) 499566063dSJacob Faibussowitsch PetscCall(PetscStrcat(file,"int64-")); 50c4762a1bSJed Brown #else 519566063dSJacob Faibussowitsch PetscCall(PetscStrcat(file,"int32-")); 52c4762a1bSJed Brown #endif 53c4762a1bSJed Brown #if defined(PETSC_USE_REAL_SINGLE) 549566063dSJacob Faibussowitsch PetscCall(PetscStrcat(file,"float32")); 55c4762a1bSJed Brown #else 569566063dSJacob Faibussowitsch PetscCall(PetscStrcat(file,"float64")); 57c4762a1bSJed Brown #endif 58c4762a1bSJed Brown } 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* Load matrix A */ 61c4762a1bSJed Brown #if defined(PETSC_USE_REAL___FLOAT128) 629566063dSJacob Faibussowitsch PetscCall(PetscOptionsInsertString(NULL,"-binary_read_double")); 63c4762a1bSJed Brown #endif 649566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd)); 659566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 669566063dSJacob Faibussowitsch PetscCall(MatLoad(A,fd)); 679566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd)); 689566063dSJacob Faibussowitsch PetscCall(MatConvert(A,MATSEQDENSE,MAT_INPLACE_MATRIX,&A)); 699566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,&m,&n)); 70*08401ef6SPierre Jolivet PetscCheck(m == n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ, "This example is not intended for rectangular matrices (%" PetscInt_FMT ", %" PetscInt_FMT ")", m, n); 71c4762a1bSJed Brown 72a5b23f4aSJose E. Roman /* Create dense matrix C and X; C holds true solution with identical columns */ 73c4762a1bSJed Brown nrhs = 2; 749566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-nrhs",&nrhs,NULL)); 759566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&C)); 769566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C,m,PETSC_DECIDE,PETSC_DECIDE,nrhs)); 779566063dSJacob Faibussowitsch PetscCall(MatSetType(C,MATDENSE)); 789566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 799566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 80c4762a1bSJed Brown 819566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rand)); 829566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rand)); 839566063dSJacob Faibussowitsch PetscCall(MatSetRandom(C,rand)); 849566063dSJacob Faibussowitsch PetscCall(MatDuplicate(C,MAT_DO_NOT_COPY_VALUES,&X)); 859566063dSJacob Faibussowitsch PetscCall(MatDuplicate(C,MAT_DO_NOT_COPY_VALUES,&RHS)); 86c4762a1bSJed Brown 87c4762a1bSJed Brown /* Create vectors */ 889566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&x)); 899566063dSJacob Faibussowitsch PetscCall(VecSetSizes(x,n,PETSC_DECIDE)); 909566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(x)); 919566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&b)); 929566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&u)); /* save the true solution */ 93c4762a1bSJed Brown 94c4762a1bSJed Brown /* make a symmetric matrix */ 95c4762a1bSJed Brown if (symm) { 96c4762a1bSJed Brown Mat AT; 97c4762a1bSJed Brown 989566063dSJacob Faibussowitsch PetscCall(MatTranspose(A,MAT_INITIAL_MATRIX,&AT)); 999566063dSJacob Faibussowitsch PetscCall(MatAXPY(A,1.0,AT,SAME_NONZERO_PATTERN)); 1009566063dSJacob Faibussowitsch PetscCall(MatDestroy(&AT)); 101c4762a1bSJed Brown ftyp = MAT_FACTOR_CHOLESKY; 102c4762a1bSJed Brown } 103c4762a1bSJed Brown /* make an hermitian matrix */ 104c4762a1bSJed Brown if (herm) { 105c4762a1bSJed Brown Mat AH; 106c4762a1bSJed Brown 1079566063dSJacob Faibussowitsch PetscCall(MatHermitianTranspose(A,MAT_INITIAL_MATRIX,&AH)); 1089566063dSJacob Faibussowitsch PetscCall(MatAXPY(A,1.0,AH,SAME_NONZERO_PATTERN)); 1099566063dSJacob Faibussowitsch PetscCall(MatDestroy(&AH)); 110c4762a1bSJed Brown ftyp = MAT_FACTOR_CHOLESKY; 111c4762a1bSJed Brown } 1129566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)A,"A")); 1139566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A,NULL,"-amat_view")); 114c4762a1bSJed Brown 1159566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&F)); 1169566063dSJacob Faibussowitsch PetscCall(MatSetOption(F,MAT_SYMMETRIC,symm)); 117c4762a1bSJed Brown /* it seems that the SPD concept in PETSc extends naturally to Hermitian Positive definitess */ 1189566063dSJacob Faibussowitsch PetscCall(MatSetOption(F,MAT_HERMITIAN,(PetscBool)(hpd || herm))); 1199566063dSJacob Faibussowitsch PetscCall(MatSetOption(F,MAT_SPD,hpd)); 1204396437dSToby Isaac { 1214396437dSToby Isaac PetscInt iftyp = ftyp; 1229566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetEList(NULL,NULL,"-ftype",MatFactorTypes,MAT_FACTOR_NUM_TYPES,&iftyp,NULL)); 1234396437dSToby Isaac ftyp = (MatFactorType) iftyp; 1244396437dSToby Isaac } 125c4762a1bSJed Brown if (ftyp == MAT_FACTOR_LU) { 1269566063dSJacob Faibussowitsch PetscCall(MatLUFactor(F,NULL,NULL,NULL)); 1274396437dSToby Isaac } else if (ftyp == MAT_FACTOR_CHOLESKY) { 1289566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactor(F,NULL,NULL)); 1294396437dSToby Isaac } else if (ftyp == MAT_FACTOR_QR) { 1309566063dSJacob Faibussowitsch PetscCall(MatQRFactor(F,NULL,NULL)); 13198921bdaSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "Factorization %s not supported in this example", MatFactorTypes[ftyp]); 132c4762a1bSJed Brown 133c4762a1bSJed Brown for (nsolve = 0; nsolve < 2; nsolve++) { 1349566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x,rand)); 1359566063dSJacob Faibussowitsch PetscCall(VecCopy(x,u)); 136c4762a1bSJed Brown if (nsolve) { 1379566063dSJacob Faibussowitsch PetscCall(MatMult(A,x,b)); 1389566063dSJacob Faibussowitsch PetscCall(MatSolve(F,b,x)); 139c4762a1bSJed Brown } else { 1409566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A,x,b)); 1419566063dSJacob Faibussowitsch PetscCall(MatSolveTranspose(F,b,x)); 142c4762a1bSJed Brown } 143c4762a1bSJed Brown /* Check the error */ 1449566063dSJacob Faibussowitsch PetscCall(VecAXPY(u,-1.0,x)); /* u <- (-1.0)x + u */ 1459566063dSJacob Faibussowitsch PetscCall(VecNorm(u,NORM_2,&norm)); 146c4762a1bSJed Brown if (norm > tol) { 147c4762a1bSJed Brown PetscReal resi; 148c4762a1bSJed Brown if (nsolve) { 1499566063dSJacob Faibussowitsch PetscCall(MatMult(A,x,u)); /* u = A*x */ 150c4762a1bSJed Brown } else { 1519566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(A,x,u)); /* u = A*x */ 152c4762a1bSJed Brown } 1539566063dSJacob Faibussowitsch PetscCall(VecAXPY(u,-1.0,b)); /* u <- (-1.0)b + u */ 1549566063dSJacob Faibussowitsch PetscCall(VecNorm(u,NORM_2,&resi)); 155c4762a1bSJed Brown if (nsolve) { 1569566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"MatSolve error: Norm of error %g, residual %g\n",(double)norm,(double)resi)); 157c4762a1bSJed Brown } else { 1589566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"MatSolveTranspose error: Norm of error %g, residual %g\n",(double)norm,(double)resi)); 159c4762a1bSJed Brown } 160c4762a1bSJed Brown } 161c4762a1bSJed Brown } 1629566063dSJacob Faibussowitsch PetscCall(MatMatMult(A,C,MAT_REUSE_MATRIX,2.0,&RHS)); 1639566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F,RHS,X)); 164c4762a1bSJed Brown 165c4762a1bSJed Brown /* Check the error */ 1669566063dSJacob Faibussowitsch PetscCall(MatAXPY(X,-1.0,C,SAME_NONZERO_PATTERN)); 1679566063dSJacob Faibussowitsch PetscCall(MatNorm(X,NORM_FROBENIUS,&norm)); 168c4762a1bSJed Brown if (norm > tol) { 1699566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"MatMatSolve: Norm of error %g\n",(double)norm)); 170c4762a1bSJed Brown } 171c4762a1bSJed Brown 172c4762a1bSJed Brown /* Free data structures */ 1739566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1749566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 1759566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 1769566063dSJacob Faibussowitsch PetscCall(MatDestroy(&X)); 1779566063dSJacob Faibussowitsch PetscCall(MatDestroy(&RHS)); 1789566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 1799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 1819566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 1829566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 183b122ec5aSJacob Faibussowitsch return 0; 184c4762a1bSJed Brown } 185c4762a1bSJed Brown 186c4762a1bSJed Brown /*TEST 187c4762a1bSJed Brown 188c4762a1bSJed Brown testset: 189c4762a1bSJed Brown output_file: output/ex215.out 190c4762a1bSJed Brown test: 191c4762a1bSJed Brown suffix: ns 192c4762a1bSJed Brown test: 193c4762a1bSJed Brown suffix: sym 194c4762a1bSJed Brown args: -symmetric_solve 195c4762a1bSJed Brown test: 196c4762a1bSJed Brown suffix: herm 197c4762a1bSJed Brown args: -hermitian_solve 198c4762a1bSJed Brown test: 199c4762a1bSJed Brown suffix: hpd 200c4762a1bSJed Brown args: -hpd_solve 2014396437dSToby Isaac test: 2024396437dSToby Isaac suffix: qr 2034396437dSToby Isaac args: -ftype qr 204c4762a1bSJed Brown 205c4762a1bSJed Brown TEST*/ 206