1 static char help[] = "Tests MatCreateConstantDiagonal().\n" 2 "\n"; 3 4 #include <petscmat.h> 5 6 /*T 7 Concepts: Mat 8 T*/ 9 10 int main(int argc, char **args) 11 { 12 PetscErrorCode ierr; 13 Vec X, Y; 14 Mat A,B,Af; 15 PetscBool flg; 16 PetscReal xnorm,ynorm; 17 18 ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 19 20 ierr = MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,20,20,3.0,&A);CHKERRQ(ierr); 21 ierr = MatCreateVecs(A,&X,&Y);CHKERRQ(ierr); 22 ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 23 24 ierr = VecSetRandom(X,NULL);CHKERRQ(ierr); 25 ierr = VecNorm(X,NORM_2,&xnorm);CHKERRQ(ierr); 26 ierr = MatMult(A,X,Y);CHKERRQ(ierr); 27 ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr); 28 if (PetscAbsReal(ynorm - 3*xnorm) > PETSC_SMALL) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Expected norm %g actual norm %g\n",(double)3*xnorm,(double)ynorm);CHKERRQ(ierr); 29 ierr = MatShift(A,5.0);CHKERRQ(ierr); 30 ierr = MatScale(A,.5);CHKERRQ(ierr); 31 ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 32 33 /* Convert to AIJ (exercises MatGetRow/MatRestoreRow) */ 34 ierr = MatConvert(A,MATAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr); 35 ierr = MatMultEqual(A,B,10,&flg);CHKERRQ(ierr); 36 if (!flg) { ierr = PetscPrintf(PETSC_COMM_WORLD,"Error MatMult\n");CHKERRQ(ierr); } 37 ierr = MatMultAddEqual(A,B,10,&flg);CHKERRQ(ierr); 38 if (!flg) { ierr = PetscPrintf(PETSC_COMM_WORLD,"Error MatMultAdd\n");CHKERRQ(ierr); } 39 ierr = MatMultTransposeEqual(A,B,10,&flg);CHKERRQ(ierr); 40 if (!flg) { ierr = PetscPrintf(PETSC_COMM_WORLD,"Error MatMultTranspose\n");CHKERRQ(ierr); } 41 ierr = MatMultTransposeAddEqual(A,B,10,&flg);CHKERRQ(ierr); 42 if (!flg) { ierr = PetscPrintf(PETSC_COMM_WORLD,"Error MatMultTransposeAdd\n");CHKERRQ(ierr); } 43 44 ierr = MatGetDiagonal(A,Y);CHKERRQ(ierr); 45 ierr = MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_LU,&Af);CHKERRQ(ierr); 46 ierr = MatLUFactorSymbolic(Af,A,NULL,NULL,NULL);CHKERRQ(ierr); 47 ierr = MatLUFactorNumeric(Af,A,NULL);CHKERRQ(ierr); 48 ierr = MatSolve(Af,X,Y);CHKERRQ(ierr); 49 ierr = VecNorm(Y,NORM_2,&ynorm);CHKERRQ(ierr); 50 if (PetscAbsReal(ynorm - xnorm/4) > PETSC_SMALL) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Expected norm %g actual norm %g\n",(double).25*xnorm,(double)ynorm);CHKERRQ(ierr); 51 52 ierr = MatDestroy(&A);CHKERRQ(ierr); 53 ierr = MatDestroy(&B);CHKERRQ(ierr); 54 ierr = MatDestroy(&Af);CHKERRQ(ierr); 55 ierr = VecDestroy(&X);CHKERRQ(ierr); 56 ierr = VecDestroy(&Y);CHKERRQ(ierr); 57 58 ierr = PetscFinalize(); 59 return ierr; 60 } 61 62 /*TEST 63 64 test: 65 nsize: 2 66 67 TEST*/ 68