1*c4762a1bSJed Brown static char help[] = "Test MatMult() for Hermitian matrix.\n\n"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown #include <petscmat.h> 4*c4762a1bSJed Brown 5*c4762a1bSJed Brown int main(int argc,char **args) 6*c4762a1bSJed Brown { 7*c4762a1bSJed Brown Mat A,As; 8*c4762a1bSJed Brown PetscBool flg; 9*c4762a1bSJed Brown PetscErrorCode ierr; 10*c4762a1bSJed Brown PetscMPIInt size; 11*c4762a1bSJed Brown PetscInt i,j; 12*c4762a1bSJed Brown PetscScalar v,sigma2; 13*c4762a1bSJed Brown PetscReal h2,sigma1=100.0; 14*c4762a1bSJed Brown PetscInt dim,Ii,J,n = 3,rstart,rend; 15*c4762a1bSJed Brown 16*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 17*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 18*c4762a1bSJed Brown ierr = PetscOptionsGetReal(NULL,NULL,"-sigma1",&sigma1,NULL);CHKERRQ(ierr); 19*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr); 20*c4762a1bSJed Brown dim = n*n; 21*c4762a1bSJed Brown 22*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 23*c4762a1bSJed Brown ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,dim,dim);CHKERRQ(ierr); 24*c4762a1bSJed Brown ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr); 25*c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 26*c4762a1bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 27*c4762a1bSJed Brown 28*c4762a1bSJed Brown sigma2 = 10.0*PETSC_i; 29*c4762a1bSJed Brown h2 = 1.0/((n+1)*(n+1)); 30*c4762a1bSJed Brown 31*c4762a1bSJed Brown ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 32*c4762a1bSJed Brown for (Ii=rstart; Ii<rend; Ii++) { 33*c4762a1bSJed Brown v = -1.0; i = Ii/n; j = Ii - i*n; 34*c4762a1bSJed Brown if (i>0) { 35*c4762a1bSJed Brown J = Ii-n; ierr = MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr); 36*c4762a1bSJed Brown } 37*c4762a1bSJed Brown if (i<n-1) { 38*c4762a1bSJed Brown J = Ii+n; ierr = MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr); 39*c4762a1bSJed Brown } 40*c4762a1bSJed Brown if (j>0) { 41*c4762a1bSJed Brown J = Ii-1; ierr = MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr); 42*c4762a1bSJed Brown } 43*c4762a1bSJed Brown if (j<n-1) { 44*c4762a1bSJed Brown J = Ii+1; ierr = MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr); 45*c4762a1bSJed Brown } 46*c4762a1bSJed Brown v = 4.0 - sigma1*h2; 47*c4762a1bSJed Brown ierr = MatSetValues(A,1,&Ii,1,&Ii,&v,ADD_VALUES);CHKERRQ(ierr); 48*c4762a1bSJed Brown } 49*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 50*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 51*c4762a1bSJed Brown 52*c4762a1bSJed Brown /* Check whether A is symmetric */ 53*c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL, "-check_symmetric", &flg);CHKERRQ(ierr); 54*c4762a1bSJed Brown if (flg) { 55*c4762a1bSJed Brown ierr = MatIsSymmetric(A,0.0,&flg);CHKERRQ(ierr); 56*c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"A is not symmetric"); 57*c4762a1bSJed Brown } 58*c4762a1bSJed Brown ierr = MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr); 59*c4762a1bSJed Brown 60*c4762a1bSJed Brown /* make A complex Hermitian */ 61*c4762a1bSJed Brown Ii = 0; J = dim-1; 62*c4762a1bSJed Brown if (Ii >= rstart && Ii < rend) { 63*c4762a1bSJed Brown v = sigma2*h2; /* RealPart(v) = 0.0 */ 64*c4762a1bSJed Brown ierr = MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr); 65*c4762a1bSJed Brown v = -sigma2*h2; 66*c4762a1bSJed Brown ierr = MatSetValues(A,1,&J,1,&Ii,&v,ADD_VALUES);CHKERRQ(ierr); 67*c4762a1bSJed Brown } 68*c4762a1bSJed Brown 69*c4762a1bSJed Brown Ii = dim-2; J = dim-1; 70*c4762a1bSJed Brown if (Ii >= rstart && Ii < rend) { 71*c4762a1bSJed Brown v = sigma2*h2; /* RealPart(v) = 0.0 */ 72*c4762a1bSJed Brown ierr = MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr); 73*c4762a1bSJed Brown v = -sigma2*h2; 74*c4762a1bSJed Brown ierr = MatSetValues(A,1,&J,1,&Ii,&v,ADD_VALUES);CHKERRQ(ierr); 75*c4762a1bSJed Brown } 76*c4762a1bSJed Brown 77*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 78*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 79*c4762a1bSJed Brown ierr = MatViewFromOptions(A,NULL,"-disp_mat");CHKERRQ(ierr); 80*c4762a1bSJed Brown 81*c4762a1bSJed Brown /* Check whether A is Hermitian, then set A->hermitian flag */ 82*c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL, "-check_Hermitian", &flg);CHKERRQ(ierr); 83*c4762a1bSJed Brown if (flg && size == 1) { 84*c4762a1bSJed Brown ierr = MatIsHermitian(A,0.0,&flg);CHKERRQ(ierr); 85*c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"A is not Hermitian"); 86*c4762a1bSJed Brown } 87*c4762a1bSJed Brown ierr = MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);CHKERRQ(ierr); 88*c4762a1bSJed Brown 89*c4762a1bSJed Brown #if defined(PETSC_HAVE_SUPERLU_DIST) 90*c4762a1bSJed Brown /* Test Cholesky factorization */ 91*c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL, "-test_choleskyfactor", &flg);CHKERRQ(ierr); 92*c4762a1bSJed Brown if (flg) { 93*c4762a1bSJed Brown Mat F; 94*c4762a1bSJed Brown IS perm,iperm; 95*c4762a1bSJed Brown MatFactorInfo info; 96*c4762a1bSJed Brown PetscInt nneg,nzero,npos; 97*c4762a1bSJed Brown 98*c4762a1bSJed Brown ierr = MatGetFactor(A,MATSOLVERSUPERLU_DIST,MAT_FACTOR_CHOLESKY,&F);CHKERRQ(ierr); 99*c4762a1bSJed Brown ierr = MatGetOrdering(A,MATORDERINGND,&perm,&iperm);CHKERRQ(ierr); 100*c4762a1bSJed Brown ierr = MatCholeskyFactorSymbolic(F,A,perm,&info);CHKERRQ(ierr); 101*c4762a1bSJed Brown ierr = MatCholeskyFactorNumeric(F,A,&info);CHKERRQ(ierr); 102*c4762a1bSJed Brown 103*c4762a1bSJed Brown ierr = MatGetInertia(F,&nneg,&nzero,&npos);CHKERRQ(ierr); 104*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," MatInertia: nneg: %D, nzero: %D, npos: %D\n",nneg,nzero,npos);CHKERRQ(ierr); 105*c4762a1bSJed Brown ierr = MatDestroy(&F);CHKERRQ(ierr); 106*c4762a1bSJed Brown ierr = ISDestroy(&perm);CHKERRQ(ierr); 107*c4762a1bSJed Brown ierr = ISDestroy(&iperm);CHKERRQ(ierr); 108*c4762a1bSJed Brown } 109*c4762a1bSJed Brown #endif 110*c4762a1bSJed Brown 111*c4762a1bSJed Brown /* Create a Hermitian matrix As in sbaij format */ 112*c4762a1bSJed Brown ierr = MatConvert(A,MATSBAIJ,MAT_INITIAL_MATRIX,&As);CHKERRQ(ierr); 113*c4762a1bSJed Brown ierr = MatViewFromOptions(As,NULL,"-disp_mat");CHKERRQ(ierr); 114*c4762a1bSJed Brown 115*c4762a1bSJed Brown /* Test MatMult */ 116*c4762a1bSJed Brown ierr = MatMultEqual(A,As,10,&flg);CHKERRQ(ierr); 117*c4762a1bSJed Brown if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"MatMult not equal"); 118*c4762a1bSJed Brown ierr = MatMultAddEqual(A,As,10,&flg);CHKERRQ(ierr); 119*c4762a1bSJed Brown if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"MatMultAdd not equal"); 120*c4762a1bSJed Brown 121*c4762a1bSJed Brown /* Free spaces */ 122*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 123*c4762a1bSJed Brown ierr = MatDestroy(&As);CHKERRQ(ierr); 124*c4762a1bSJed Brown ierr = PetscFinalize(); 125*c4762a1bSJed Brown return ierr; 126*c4762a1bSJed Brown } 127*c4762a1bSJed Brown 128*c4762a1bSJed Brown 129*c4762a1bSJed Brown /*TEST 130*c4762a1bSJed Brown 131*c4762a1bSJed Brown build: 132*c4762a1bSJed Brown requires: complex 133*c4762a1bSJed Brown 134*c4762a1bSJed Brown test: 135*c4762a1bSJed Brown args: -n 1000 136*c4762a1bSJed Brown output_file: output/ex127.out 137*c4762a1bSJed Brown 138*c4762a1bSJed Brown test: 139*c4762a1bSJed Brown suffix: 2 140*c4762a1bSJed Brown nsize: 3 141*c4762a1bSJed Brown args: -n 1000 142*c4762a1bSJed Brown output_file: output/ex127.out 143*c4762a1bSJed Brown 144*c4762a1bSJed Brown 145*c4762a1bSJed Brown test: 146*c4762a1bSJed Brown suffix: superlu_dist 147*c4762a1bSJed Brown nsize: 3 148*c4762a1bSJed Brown requires: superlu_dist 149*c4762a1bSJed Brown args: -test_choleskyfactor -mat_superlu_dist_rowperm NOROWPERM 150*c4762a1bSJed Brown output_file: output/ex127_superlu_dist.out 151*c4762a1bSJed Brown TEST*/ 152