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