1c4762a1bSJed Brown static char help[] = "Test MatMult() for Hermitian matrix.\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown 5c4762a1bSJed Brown int main(int argc,char **args) 6c4762a1bSJed Brown { 7c4762a1bSJed Brown Mat A,As; 8c4762a1bSJed Brown PetscBool flg; 9c4762a1bSJed Brown PetscMPIInt size; 10c4762a1bSJed Brown PetscInt i,j; 11c4762a1bSJed Brown PetscScalar v,sigma2; 12c4762a1bSJed Brown PetscReal h2,sigma1=100.0; 13c4762a1bSJed Brown PetscInt dim,Ii,J,n = 3,rstart,rend; 14c4762a1bSJed Brown 15*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 165f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 175f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-sigma1",&sigma1,NULL)); 185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 19c4762a1bSJed Brown dim = n*n; 20c4762a1bSJed Brown 215f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 225f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,dim,dim)); 235f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATAIJ)); 245f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 255f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 26c4762a1bSJed Brown 27c4762a1bSJed Brown sigma2 = 10.0*PETSC_i; 28c4762a1bSJed Brown h2 = 1.0/((n+1)*(n+1)); 29c4762a1bSJed Brown 305f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(A,&rstart,&rend)); 31c4762a1bSJed Brown for (Ii=rstart; Ii<rend; Ii++) { 32c4762a1bSJed Brown v = -1.0; i = Ii/n; j = Ii - i*n; 33c4762a1bSJed Brown if (i>0) { 345f80ce2aSJacob Faibussowitsch J = Ii-n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES)); 35c4762a1bSJed Brown } 36c4762a1bSJed Brown if (i<n-1) { 375f80ce2aSJacob Faibussowitsch J = Ii+n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES)); 38c4762a1bSJed Brown } 39c4762a1bSJed Brown if (j>0) { 405f80ce2aSJacob Faibussowitsch J = Ii-1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES)); 41c4762a1bSJed Brown } 42c4762a1bSJed Brown if (j<n-1) { 435f80ce2aSJacob Faibussowitsch J = Ii+1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES)); 44c4762a1bSJed Brown } 45c4762a1bSJed Brown v = 4.0 - sigma1*h2; 465f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&Ii,1,&Ii,&v,ADD_VALUES)); 47c4762a1bSJed Brown } 485f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 495f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 50c4762a1bSJed Brown 51c4762a1bSJed Brown /* Check whether A is symmetric */ 525f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL, "-check_symmetric", &flg)); 53c4762a1bSJed Brown if (flg) { 545f80ce2aSJacob Faibussowitsch CHKERRQ(MatIsSymmetric(A,0.0,&flg)); 5528b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_USER,"A is not symmetric"); 56c4762a1bSJed Brown } 575f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE)); 58c4762a1bSJed Brown 59c4762a1bSJed Brown /* make A complex Hermitian */ 60c4762a1bSJed Brown Ii = 0; J = dim-1; 61c4762a1bSJed Brown if (Ii >= rstart && Ii < rend) { 62c4762a1bSJed Brown v = sigma2*h2; /* RealPart(v) = 0.0 */ 635f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES)); 64c4762a1bSJed Brown v = -sigma2*h2; 655f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&J,1,&Ii,&v,ADD_VALUES)); 66c4762a1bSJed Brown } 67c4762a1bSJed Brown 68c4762a1bSJed Brown Ii = dim-2; J = dim-1; 69c4762a1bSJed Brown if (Ii >= rstart && Ii < rend) { 70c4762a1bSJed Brown v = sigma2*h2; /* RealPart(v) = 0.0 */ 715f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES)); 72c4762a1bSJed Brown v = -sigma2*h2; 735f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,1,&J,1,&Ii,&v,ADD_VALUES)); 74c4762a1bSJed Brown } 75c4762a1bSJed Brown 765f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 775f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 785f80ce2aSJacob Faibussowitsch CHKERRQ(MatViewFromOptions(A,NULL,"-disp_mat")); 79c4762a1bSJed Brown 80c4762a1bSJed Brown /* Check whether A is Hermitian, then set A->hermitian flag */ 815f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL, "-check_Hermitian", &flg)); 82c4762a1bSJed Brown if (flg && size == 1) { 835f80ce2aSJacob Faibussowitsch CHKERRQ(MatIsHermitian(A,0.0,&flg)); 8428b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_USER,"A is not Hermitian"); 85c4762a1bSJed Brown } 865f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE)); 87c4762a1bSJed Brown 88c4762a1bSJed Brown #if defined(PETSC_HAVE_SUPERLU_DIST) 89c4762a1bSJed Brown /* Test Cholesky factorization */ 905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL, "-test_choleskyfactor", &flg)); 91c4762a1bSJed Brown if (flg) { 92c4762a1bSJed Brown Mat F; 93c4762a1bSJed Brown IS perm,iperm; 94c4762a1bSJed Brown MatFactorInfo info; 95c4762a1bSJed Brown PetscInt nneg,nzero,npos; 96c4762a1bSJed Brown 975f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetFactor(A,MATSOLVERSUPERLU_DIST,MAT_FACTOR_CHOLESKY,&F)); 985f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOrdering(A,MATORDERINGND,&perm,&iperm)); 995f80ce2aSJacob Faibussowitsch CHKERRQ(MatCholeskyFactorSymbolic(F,A,perm,&info)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(MatCholeskyFactorNumeric(F,A,&info)); 101c4762a1bSJed Brown 1025f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetInertia(F,&nneg,&nzero,&npos)); 1035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," MatInertia: nneg: %" PetscInt_FMT ", nzero: %" PetscInt_FMT ", npos: %" PetscInt_FMT "\n",nneg,nzero,npos)); 1045f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&F)); 1055f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&perm)); 1065f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&iperm)); 107c4762a1bSJed Brown } 108c4762a1bSJed Brown #endif 109c4762a1bSJed Brown 110c4762a1bSJed Brown /* Create a Hermitian matrix As in sbaij format */ 1115f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(A,MATSBAIJ,MAT_INITIAL_MATRIX,&As)); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(MatViewFromOptions(As,NULL,"-disp_mat")); 113c4762a1bSJed Brown 114c4762a1bSJed Brown /* Test MatMult */ 1155f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(A,As,10,&flg)); 11628b400f6SJacob Faibussowitsch PetscCheck(flg,PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"MatMult not equal"); 1175f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAddEqual(A,As,10,&flg)); 11828b400f6SJacob Faibussowitsch PetscCheck(flg,PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"MatMultAdd not equal"); 119c4762a1bSJed Brown 120c4762a1bSJed Brown /* Free spaces */ 1215f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&As)); 123*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 124*b122ec5aSJacob Faibussowitsch return 0; 125c4762a1bSJed Brown } 126c4762a1bSJed Brown 127c4762a1bSJed Brown /*TEST 128c4762a1bSJed Brown 129c4762a1bSJed Brown build: 130c4762a1bSJed Brown requires: complex 131c4762a1bSJed Brown 132c4762a1bSJed Brown test: 133c4762a1bSJed Brown args: -n 1000 134c4762a1bSJed Brown output_file: output/ex127.out 135c4762a1bSJed Brown 136c4762a1bSJed Brown test: 137c4762a1bSJed Brown suffix: 2 138c4762a1bSJed Brown nsize: 3 139c4762a1bSJed Brown args: -n 1000 140c4762a1bSJed Brown output_file: output/ex127.out 141c4762a1bSJed Brown 142c4762a1bSJed Brown test: 143c4762a1bSJed Brown suffix: superlu_dist 144c4762a1bSJed Brown nsize: 3 145c4762a1bSJed Brown requires: superlu_dist 146c4762a1bSJed Brown args: -test_choleskyfactor -mat_superlu_dist_rowperm NOROWPERM 147c4762a1bSJed Brown output_file: output/ex127_superlu_dist.out 148c4762a1bSJed Brown TEST*/ 149