1 2 static char help[] = "Tests ScaLAPACK interface.\n\n"; 3 4 #include <petscmat.h> 5 6 int main(int argc,char **args) 7 { 8 Mat Cdense,Caij,B,C,Ct; 9 Vec d; 10 PetscInt i,j,M = 5,N,mb = 2,nb,nrows,ncols,mloc,nloc; 11 const PetscInt *rows,*cols; 12 IS isrows,iscols; 13 PetscErrorCode ierr; 14 PetscScalar *v; 15 PetscMPIInt rank; 16 PetscReal Cnorm; 17 PetscBool flg,mats_view=PETSC_FALSE; 18 19 ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 20 ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 21 ierr = PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);CHKERRQ(ierr); 22 N = M; 23 ierr = PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL);CHKERRQ(ierr); 24 ierr = PetscOptionsGetInt(NULL,NULL,"-mb",&mb,NULL);CHKERRQ(ierr); 25 nb = mb; 26 ierr = PetscOptionsGetInt(NULL,NULL,"-nb",&nb,NULL);CHKERRQ(ierr); 27 ierr = PetscOptionsHasName(NULL,NULL,"-mats_view",&mats_view);CHKERRQ(ierr); 28 29 ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr); 30 ierr = MatSetType(C,MATSCALAPACK);CHKERRQ(ierr); 31 mloc = PETSC_DECIDE; 32 ierr = PetscSplitOwnershipEqual(PETSC_COMM_WORLD,&mloc,&M);CHKERRQ(ierr); 33 nloc = PETSC_DECIDE; 34 ierr = PetscSplitOwnershipEqual(PETSC_COMM_WORLD,&nloc,&N);CHKERRQ(ierr); 35 ierr = MatSetSizes(C,mloc,nloc,M,N);CHKERRQ(ierr); 36 ierr = MatScaLAPACKSetBlockSizes(C,mb,nb);CHKERRQ(ierr); 37 ierr = MatSetFromOptions(C);CHKERRQ(ierr); 38 ierr = MatSetUp(C);CHKERRQ(ierr); 39 /*ierr = MatCreateScaLAPACK(PETSC_COMM_WORLD,mb,nb,M,N,0,0,&C);CHKERRQ(ierr); */ 40 41 ierr = MatGetOwnershipIS(C,&isrows,&iscols);CHKERRQ(ierr); 42 ierr = ISGetLocalSize(isrows,&nrows);CHKERRQ(ierr); 43 ierr = ISGetIndices(isrows,&rows);CHKERRQ(ierr); 44 ierr = ISGetLocalSize(iscols,&ncols);CHKERRQ(ierr); 45 ierr = ISGetIndices(iscols,&cols);CHKERRQ(ierr); 46 ierr = PetscMalloc1(nrows*ncols,&v);CHKERRQ(ierr); 47 for (i=0;i<nrows;i++) { 48 for (j=0;j<ncols;j++) v[i*ncols+j] = (PetscReal)(rows[i]+1+(cols[j]+1)*0.01); 49 } 50 ierr = MatSetValues(C,nrows,rows,ncols,cols,v,INSERT_VALUES);CHKERRQ(ierr); 51 ierr = PetscFree(v);CHKERRQ(ierr); 52 ierr = ISRestoreIndices(isrows,&rows);CHKERRQ(ierr); 53 ierr = ISRestoreIndices(iscols,&cols);CHKERRQ(ierr); 54 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 55 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 56 ierr = ISDestroy(&isrows);CHKERRQ(ierr); 57 ierr = ISDestroy(&iscols);CHKERRQ(ierr); 58 59 /* Test MatView(), MatDuplicate() and out-of-place MatConvert() */ 60 ierr = MatDuplicate(C,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 61 if (mats_view) { 62 ierr = PetscPrintf(PETSC_COMM_WORLD,"Duplicated C:\n");CHKERRQ(ierr); 63 ierr = MatView(B,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 64 } 65 ierr = MatDestroy(&B);CHKERRQ(ierr); 66 ierr = MatConvert(C,MATDENSE,MAT_INITIAL_MATRIX,&Cdense);CHKERRQ(ierr); 67 ierr = MatMultEqual(C,Cdense,5,&flg);CHKERRQ(ierr); 68 if (!flg) SETERRQ(PETSC_COMM_WORLD,1,"Check fails: Cdense != C"); 69 70 /* Test MatNorm() */ 71 ierr = MatNorm(C,NORM_1,&Cnorm);CHKERRQ(ierr); 72 73 /* Test MatTranspose(), MatZeroEntries() and MatGetDiagonal() */ 74 ierr = MatTranspose(C,MAT_INITIAL_MATRIX,&Ct);CHKERRQ(ierr); 75 ierr = MatConjugate(Ct);CHKERRQ(ierr); 76 if (mats_view) { 77 ierr = PetscPrintf(PETSC_COMM_WORLD,"C's Transpose Conjugate:\n");CHKERRQ(ierr); 78 ierr = MatView(Ct,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 79 } 80 ierr = MatZeroEntries(Ct);CHKERRQ(ierr); 81 if (M>N) { ierr = MatCreateVecs(C,&d,NULL);CHKERRQ(ierr); } 82 else { ierr = MatCreateVecs(C,NULL,&d);CHKERRQ(ierr); } 83 ierr = MatGetDiagonal(C,d);CHKERRQ(ierr); 84 if (mats_view) { 85 ierr = PetscPrintf(PETSC_COMM_WORLD,"Diagonal of C:\n");CHKERRQ(ierr); 86 ierr = VecView(d,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 87 } 88 if (M>N) { 89 ierr = MatDiagonalScale(C,NULL,d);CHKERRQ(ierr); 90 } else { 91 ierr = MatDiagonalScale(C,d,NULL);CHKERRQ(ierr); 92 } 93 if (mats_view) { 94 ierr = PetscPrintf(PETSC_COMM_WORLD,"Diagonal Scaled C:\n");CHKERRQ(ierr); 95 ierr = MatView(C,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 96 } 97 98 /* Test MatAXPY(), MatAYPX() and in-place MatConvert() */ 99 ierr = MatCreate(PETSC_COMM_WORLD,&B);CHKERRQ(ierr); 100 ierr = MatSetType(B,MATSCALAPACK);CHKERRQ(ierr); 101 ierr = MatSetSizes(B,mloc,nloc,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 102 ierr = MatScaLAPACKSetBlockSizes(B,mb,nb);CHKERRQ(ierr); 103 ierr = MatSetFromOptions(B);CHKERRQ(ierr); 104 ierr = MatSetUp(B);CHKERRQ(ierr); 105 /* ierr = MatCreateScaLAPACK(PETSC_COMM_WORLD,mb,nb,M,N,0,0,&B);CHKERRQ(ierr); */ 106 ierr = MatGetOwnershipIS(B,&isrows,&iscols);CHKERRQ(ierr); 107 ierr = ISGetLocalSize(isrows,&nrows);CHKERRQ(ierr); 108 ierr = ISGetIndices(isrows,&rows);CHKERRQ(ierr); 109 ierr = ISGetLocalSize(iscols,&ncols);CHKERRQ(ierr); 110 ierr = ISGetIndices(iscols,&cols);CHKERRQ(ierr); 111 ierr = PetscMalloc1(nrows*ncols,&v);CHKERRQ(ierr); 112 for (i=0;i<nrows;i++) { 113 for (j=0;j<ncols;j++) v[i*ncols+j] = (PetscReal)(1000*rows[i]+cols[j]); 114 } 115 ierr = MatSetValues(B,nrows,rows,ncols,cols,v,INSERT_VALUES);CHKERRQ(ierr); 116 ierr = PetscFree(v);CHKERRQ(ierr); 117 ierr = ISRestoreIndices(isrows,&rows);CHKERRQ(ierr); 118 ierr = ISRestoreIndices(iscols,&cols);CHKERRQ(ierr); 119 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 120 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 121 if (mats_view) { 122 ierr = PetscPrintf(PETSC_COMM_WORLD,"B:\n");CHKERRQ(ierr); 123 ierr = MatView(B,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 124 } 125 ierr = MatAXPY(B,2.5,C,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 126 ierr = MatAYPX(B,3.75,C,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 127 ierr = MatConvert(B,MATDENSE,MAT_INPLACE_MATRIX,&B);CHKERRQ(ierr); 128 if (mats_view) { 129 ierr = PetscPrintf(PETSC_COMM_WORLD,"B after MatAXPY and MatAYPX:\n");CHKERRQ(ierr); 130 ierr = MatView(B,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 131 } 132 ierr = ISDestroy(&isrows);CHKERRQ(ierr); 133 ierr = ISDestroy(&iscols);CHKERRQ(ierr); 134 ierr = MatDestroy(&B);CHKERRQ(ierr); 135 136 /* Test MatMatTransposeMult(): B = C*C^T */ 137 ierr = MatMatTransposeMult(C,C,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&B);CHKERRQ(ierr); 138 ierr = MatScale(C,2.0);CHKERRQ(ierr); 139 ierr = MatMatTransposeMult(C,C,MAT_REUSE_MATRIX,PETSC_DEFAULT,&B);CHKERRQ(ierr); 140 ierr = MatMatTransposeMultEqual(C,C,B,10,&flg);CHKERRQ(ierr); 141 if (!flg) SETERRQ(PETSC_COMM_WORLD,1,"Chech fails: B != C*C^T"); 142 143 if (mats_view) { 144 ierr = PetscPrintf(PETSC_COMM_WORLD,"C MatMatTransposeMult C:\n");CHKERRQ(ierr); 145 ierr = MatView(B,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 146 } 147 148 /* Test MatMult() */ 149 ierr = MatComputeOperator(C,MATAIJ,&Caij);CHKERRQ(ierr); 150 ierr = MatMultEqual(C,Caij,5,&flg);CHKERRQ(ierr); 151 if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_NOTSAMETYPE,"C != Caij. MatMultEqual() fails"); 152 ierr = MatMultTransposeEqual(C,Caij,5,&flg);CHKERRQ(ierr); 153 if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_NOTSAMETYPE,"C != Caij. MatMultTransposeEqual() fails"); 154 155 /* Test MatMultAdd() and MatMultTransposeAddEqual() */ 156 ierr = MatMultAddEqual(C,Caij,5,&flg);CHKERRQ(ierr); 157 if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_NOTSAMETYPE,"C != Caij. MatMultAddEqual() fails"); 158 ierr = MatMultTransposeAddEqual(C,Caij,5,&flg);CHKERRQ(ierr); 159 if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_NOTSAMETYPE,"C != Caij. MatMultTransposeAddEqual() fails"); 160 161 /* Test MatMatMult() */ 162 ierr = PetscOptionsHasName(NULL,NULL,"-test_matmatmult",&flg);CHKERRQ(ierr); 163 if (flg) { 164 Mat CC,CCaij; 165 ierr = MatMatMult(C,C,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&CC);CHKERRQ(ierr); 166 ierr = MatMatMult(Caij,Caij,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&CCaij);CHKERRQ(ierr); 167 ierr = MatMultEqual(CC,CCaij,5,&flg);CHKERRQ(ierr); 168 if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_NOTSAMETYPE,"CC != CCaij. MatMatMult() fails"); 169 ierr = MatDestroy(&CCaij);CHKERRQ(ierr); 170 ierr = MatDestroy(&CC);CHKERRQ(ierr); 171 } 172 173 ierr = MatDestroy(&Cdense);CHKERRQ(ierr); 174 ierr = MatDestroy(&Caij);CHKERRQ(ierr); 175 ierr = MatDestroy(&B);CHKERRQ(ierr); 176 ierr = MatDestroy(&C);CHKERRQ(ierr); 177 ierr = MatDestroy(&Ct);CHKERRQ(ierr); 178 ierr = VecDestroy(&d);CHKERRQ(ierr); 179 ierr = PetscFinalize(); 180 return ierr; 181 } 182 183 /*TEST 184 185 build: 186 requires: scalapack 187 188 test: 189 nsize: 2 190 args: -mb 5 -nb 5 -M 12 -N 10 191 requires: scalapack 192 193 test: 194 suffix: 2 195 nsize: 6 196 args: -mb 8 -nb 6 -M 20 -N 50 197 requires: scalapack 198 output_file: output/ex242_1.out 199 200 test: 201 suffix: 3 202 nsize: 3 203 args: -mb 2 -nb 2 -M 20 -N 20 -test_matmatmult 204 requires: scalapack 205 output_file: output/ex242_1.out 206 207 TEST*/ 208