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*327415f7SBarry Smith PetscFunctionBeginUser; 169566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 189566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-sigma1",&sigma1,NULL)); 199566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 20c4762a1bSJed Brown dim = n*n; 21c4762a1bSJed Brown 229566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 239566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,dim,dim)); 249566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATAIJ)); 259566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 269566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 27c4762a1bSJed Brown 28c4762a1bSJed Brown sigma2 = 10.0*PETSC_i; 29c4762a1bSJed Brown h2 = 1.0/((n+1)*(n+1)); 30c4762a1bSJed Brown 319566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&rstart,&rend)); 32c4762a1bSJed Brown for (Ii=rstart; Ii<rend; Ii++) { 33c4762a1bSJed Brown v = -1.0; i = Ii/n; j = Ii - i*n; 34c4762a1bSJed Brown if (i>0) { 359566063dSJacob Faibussowitsch J = Ii-n; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES)); 36c4762a1bSJed Brown } 37c4762a1bSJed Brown if (i<n-1) { 389566063dSJacob Faibussowitsch J = Ii+n; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES)); 39c4762a1bSJed Brown } 40c4762a1bSJed Brown if (j>0) { 419566063dSJacob Faibussowitsch J = Ii-1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES)); 42c4762a1bSJed Brown } 43c4762a1bSJed Brown if (j<n-1) { 449566063dSJacob Faibussowitsch J = Ii+1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES)); 45c4762a1bSJed Brown } 46c4762a1bSJed Brown v = 4.0 - sigma1*h2; 479566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&v,ADD_VALUES)); 48c4762a1bSJed Brown } 499566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 509566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 51c4762a1bSJed Brown 52c4762a1bSJed Brown /* Check whether A is symmetric */ 539566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL, "-check_symmetric", &flg)); 54c4762a1bSJed Brown if (flg) { 559566063dSJacob Faibussowitsch PetscCall(MatIsSymmetric(A,0.0,&flg)); 5628b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_USER,"A is not symmetric"); 57c4762a1bSJed Brown } 589566063dSJacob Faibussowitsch PetscCall(MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE)); 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* make A complex Hermitian */ 61c4762a1bSJed Brown Ii = 0; J = dim-1; 62c4762a1bSJed Brown if (Ii >= rstart && Ii < rend) { 63c4762a1bSJed Brown v = sigma2*h2; /* RealPart(v) = 0.0 */ 649566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES)); 65c4762a1bSJed Brown v = -sigma2*h2; 669566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&J,1,&Ii,&v,ADD_VALUES)); 67c4762a1bSJed Brown } 68c4762a1bSJed Brown 69c4762a1bSJed Brown Ii = dim-2; J = dim-1; 70c4762a1bSJed Brown if (Ii >= rstart && Ii < rend) { 71c4762a1bSJed Brown v = sigma2*h2; /* RealPart(v) = 0.0 */ 729566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,ADD_VALUES)); 73c4762a1bSJed Brown v = -sigma2*h2; 749566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&J,1,&Ii,&v,ADD_VALUES)); 75c4762a1bSJed Brown } 76c4762a1bSJed Brown 779566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 789566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 799566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A,NULL,"-disp_mat")); 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* Check whether A is Hermitian, then set A->hermitian flag */ 829566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL, "-check_Hermitian", &flg)); 83c4762a1bSJed Brown if (flg && size == 1) { 849566063dSJacob Faibussowitsch PetscCall(MatIsHermitian(A,0.0,&flg)); 8528b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_USER,"A is not Hermitian"); 86c4762a1bSJed Brown } 879566063dSJacob Faibussowitsch PetscCall(MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE)); 88c4762a1bSJed Brown 89c4762a1bSJed Brown #if defined(PETSC_HAVE_SUPERLU_DIST) 90c4762a1bSJed Brown /* Test Cholesky factorization */ 919566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL, "-test_choleskyfactor", &flg)); 92c4762a1bSJed Brown if (flg) { 93c4762a1bSJed Brown Mat F; 94c4762a1bSJed Brown IS perm,iperm; 95c4762a1bSJed Brown MatFactorInfo info; 96c4762a1bSJed Brown PetscInt nneg,nzero,npos; 97c4762a1bSJed Brown 989566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A,MATSOLVERSUPERLU_DIST,MAT_FACTOR_CHOLESKY,&F)); 999566063dSJacob Faibussowitsch PetscCall(MatGetOrdering(A,MATORDERINGND,&perm,&iperm)); 1009566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(F,A,perm,&info)); 1019566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(F,A,&info)); 102c4762a1bSJed Brown 1039566063dSJacob Faibussowitsch PetscCall(MatGetInertia(F,&nneg,&nzero,&npos)); 1049566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD," MatInertia: nneg: %" PetscInt_FMT ", nzero: %" PetscInt_FMT ", npos: %" PetscInt_FMT "\n",nneg,nzero,npos)); 1059566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 1069566063dSJacob Faibussowitsch PetscCall(ISDestroy(&perm)); 1079566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iperm)); 108c4762a1bSJed Brown } 109c4762a1bSJed Brown #endif 110c4762a1bSJed Brown 111c4762a1bSJed Brown /* Create a Hermitian matrix As in sbaij format */ 1129566063dSJacob Faibussowitsch PetscCall(MatConvert(A,MATSBAIJ,MAT_INITIAL_MATRIX,&As)); 1139566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(As,NULL,"-disp_mat")); 114c4762a1bSJed Brown 115c4762a1bSJed Brown /* Test MatMult */ 1169566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A,As,10,&flg)); 11728b400f6SJacob Faibussowitsch PetscCheck(flg,PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"MatMult not equal"); 1189566063dSJacob Faibussowitsch PetscCall(MatMultAddEqual(A,As,10,&flg)); 11928b400f6SJacob Faibussowitsch PetscCheck(flg,PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"MatMultAdd not equal"); 120c4762a1bSJed Brown 121c4762a1bSJed Brown /* Free spaces */ 1229566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1239566063dSJacob Faibussowitsch PetscCall(MatDestroy(&As)); 1249566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 125b122ec5aSJacob Faibussowitsch return 0; 126c4762a1bSJed Brown } 127c4762a1bSJed Brown 128c4762a1bSJed Brown /*TEST 129c4762a1bSJed Brown 130c4762a1bSJed Brown build: 131c4762a1bSJed Brown requires: complex 132c4762a1bSJed Brown 133c4762a1bSJed Brown test: 134c4762a1bSJed Brown args: -n 1000 135c4762a1bSJed Brown output_file: output/ex127.out 136c4762a1bSJed Brown 137c4762a1bSJed Brown test: 138c4762a1bSJed Brown suffix: 2 139c4762a1bSJed Brown nsize: 3 140c4762a1bSJed Brown args: -n 1000 141c4762a1bSJed Brown output_file: output/ex127.out 142c4762a1bSJed Brown 143c4762a1bSJed Brown test: 144c4762a1bSJed Brown suffix: superlu_dist 145c4762a1bSJed Brown nsize: 3 146c4762a1bSJed Brown requires: superlu_dist 147c4762a1bSJed Brown args: -test_choleskyfactor -mat_superlu_dist_rowperm NOROWPERM 148c4762a1bSJed Brown output_file: output/ex127_superlu_dist.out 149c4762a1bSJed Brown TEST*/ 150