1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatMult(), MatMultAdd(), MatMultTranspose().\n\ 3c4762a1bSJed Brown Also MatMultTransposeAdd(), MatScale(), MatGetDiagonal(), and MatDiagonalScale().\n\n"; 4c4762a1bSJed Brown 5c4762a1bSJed Brown #include <petscmat.h> 6c4762a1bSJed Brown 7c4762a1bSJed Brown int main(int argc,char **args) 8c4762a1bSJed Brown { 9c4762a1bSJed Brown Mat C; 10c4762a1bSJed Brown Vec s,u,w,x,y,z; 11c4762a1bSJed Brown PetscInt i,j,m = 8,n,rstart,rend,vstart,vend; 12c4762a1bSJed Brown PetscScalar one = 1.0,negone = -1.0,v,alpha=0.1; 13c4762a1bSJed Brown PetscReal norm, tol = PETSC_SQRT_MACHINE_EPSILON; 14c4762a1bSJed Brown PetscBool flg; 15c4762a1bSJed Brown 16*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 175f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_COMMON)); 185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 19c4762a1bSJed Brown n = m; 205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-rectA",&flg)); 21c4762a1bSJed Brown if (flg) n += 2; 225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-rectB",&flg)); 23c4762a1bSJed Brown if (flg) n -= 2; 24c4762a1bSJed Brown 25c4762a1bSJed Brown /* ---------- Assemble matrix and vectors ----------- */ 26c4762a1bSJed Brown 275f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m,n)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(C)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(C)); 315f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(C,&rstart,&rend)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&x)); 335f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(x,PETSC_DECIDE,m)); 345f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(x)); 355f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&z)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&w)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&y)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(y,PETSC_DECIDE,n)); 395f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(y)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(y,&u)); 415f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(y,&s)); 425f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetOwnershipRange(y,&vstart,&vend)); 43c4762a1bSJed Brown 44c4762a1bSJed Brown /* Assembly */ 45c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 46c4762a1bSJed Brown v = 100*(i+1); 475f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(z,1,&i,&v,INSERT_VALUES)); 48c4762a1bSJed Brown for (j=0; j<n; j++) { 49c4762a1bSJed Brown v = 10*(i+1)+j+1; 505f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES)); 51c4762a1bSJed Brown } 52c4762a1bSJed Brown } 53c4762a1bSJed Brown 54c4762a1bSJed Brown /* Flush off proc Vec values and do more assembly */ 555f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(z)); 56c4762a1bSJed Brown for (i=vstart; i<vend; i++) { 57c4762a1bSJed Brown v = one*((PetscReal)i); 585f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(y,1,&i,&v,INSERT_VALUES)); 59c4762a1bSJed Brown v = 100.0*i; 605f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(u,1,&i,&v,INSERT_VALUES)); 61c4762a1bSJed Brown } 62c4762a1bSJed Brown 63c4762a1bSJed Brown /* Flush off proc Mat values and do more assembly */ 645f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(C,MAT_FLUSH_ASSEMBLY)); 65c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 66c4762a1bSJed Brown for (j=0; j<n; j++) { 67c4762a1bSJed Brown v = 10*(i+1)+j+1; 685f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES)); 69c4762a1bSJed Brown } 70c4762a1bSJed Brown } 71c4762a1bSJed Brown /* Try overlap Coomunication with the next stage XXXSetValues */ 725f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(z)); 73c4762a1bSJed Brown 745f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(C,MAT_FLUSH_ASSEMBLY)); 75c4762a1bSJed Brown CHKMEMQ; 76c4762a1bSJed Brown /* The Assembly for the second Stage */ 775f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 785f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 795f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(y)); 805f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(y)); 815f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(C,alpha)); 825f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(u)); 835f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(u)); 845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"testing MatMult()\n")); 85c4762a1bSJed Brown CHKMEMQ; 865f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(C,y,x)); 87c4762a1bSJed Brown CHKMEMQ; 885f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 895f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"testing MatMultAdd()\n")); 905f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAdd(C,y,z,w)); 915f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(x,one,z)); 925f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(x,negone,w)); 935f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(x,NORM_2,&norm)); 94c4762a1bSJed Brown if (norm > tol) { 955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Norm of error difference = %g\n",(double)norm)); 96c4762a1bSJed Brown } 97c4762a1bSJed Brown 98c4762a1bSJed Brown /* ------- Test MatMultTranspose(), MatMultTransposeAdd() ------- */ 99c4762a1bSJed Brown 100c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 101c4762a1bSJed Brown v = one*((PetscReal)i); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(x,1,&i,&v,INSERT_VALUES)); 103c4762a1bSJed Brown } 1045f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(x)); 1055f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(x)); 1065f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"testing MatMultTranspose()\n")); 1075f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(C,x,y)); 1085f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(y,PETSC_VIEWER_STDOUT_WORLD)); 109c4762a1bSJed Brown 1105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"testing MatMultTransposeAdd()\n")); 1115f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTransposeAdd(C,x,u,s)); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(y,one,u)); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(y,negone,s)); 1145f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(y,NORM_2,&norm)); 115c4762a1bSJed Brown if (norm > tol) { 1165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Norm of error difference = %g\n",(double)norm)); 117c4762a1bSJed Brown } 118c4762a1bSJed Brown 119c4762a1bSJed Brown /* -------------------- Test MatGetDiagonal() ------------------ */ 120c4762a1bSJed Brown 1215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"testing MatGetDiagonal(), MatDiagonalScale()\n")); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 1235f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(x,one)); 1245f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetDiagonal(C,x)); 1255f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 126c4762a1bSJed Brown for (i=vstart; i<vend; i++) { 127c4762a1bSJed Brown v = one*((PetscReal)(i+1)); 1285f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetValues(y,1,&i,&v,INSERT_VALUES)); 129c4762a1bSJed Brown } 130c4762a1bSJed Brown 131c4762a1bSJed Brown /* -------------------- Test () MatDiagonalScale ------------------ */ 1325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-test_diagonalscale",&flg)); 133c4762a1bSJed Brown if (flg) { 1345f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(C,x,y)); 1355f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 136c4762a1bSJed Brown } 137c4762a1bSJed Brown /* Free data structures */ 1385f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&u)); CHKERRQ(VecDestroy(&s)); 1395f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&w)); CHKERRQ(VecDestroy(&x)); 1405f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&y)); CHKERRQ(VecDestroy(&z)); 1415f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 142c4762a1bSJed Brown 143*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 144*b122ec5aSJacob Faibussowitsch return 0; 145c4762a1bSJed Brown } 146c4762a1bSJed Brown 147c4762a1bSJed Brown /*TEST 148c4762a1bSJed Brown 149c4762a1bSJed Brown test: 150c4762a1bSJed Brown suffix: 11_A 151c4762a1bSJed Brown args: -mat_type seqaij -rectA 152c4762a1bSJed Brown filter: grep -v type 153c4762a1bSJed Brown 154c4762a1bSJed Brown test: 155c4762a1bSJed Brown args: -mat_type seqdense -rectA 156c4762a1bSJed Brown suffix: 12_A 157c4762a1bSJed Brown 158c4762a1bSJed Brown test: 159c4762a1bSJed Brown args: -mat_type seqaij -rectB 160c4762a1bSJed Brown suffix: 11_B 161c4762a1bSJed Brown filter: grep -v type 162c4762a1bSJed Brown 163c4762a1bSJed Brown test: 164c4762a1bSJed Brown args: -mat_type seqdense -rectB 165c4762a1bSJed Brown suffix: 12_B 166c4762a1bSJed Brown 167c4762a1bSJed Brown test: 168c4762a1bSJed Brown suffix: 21 169c4762a1bSJed Brown args: -mat_type mpiaij 170c4762a1bSJed Brown filter: grep -v type 171c4762a1bSJed Brown 172c4762a1bSJed Brown test: 173c4762a1bSJed Brown suffix: 22 174c4762a1bSJed Brown args: -mat_type mpidense 175c4762a1bSJed Brown 176c4762a1bSJed Brown test: 177c4762a1bSJed Brown suffix: 23 178c4762a1bSJed Brown nsize: 3 179c4762a1bSJed Brown args: -mat_type mpiaij 180c4762a1bSJed Brown filter: grep -v type 181c4762a1bSJed Brown 182c4762a1bSJed Brown test: 183c4762a1bSJed Brown suffix: 24 184c4762a1bSJed Brown nsize: 3 185c4762a1bSJed Brown args: -mat_type mpidense 186c4762a1bSJed Brown 187c4762a1bSJed Brown test: 188c4762a1bSJed Brown suffix: 2_aijcusparse_1 189c4762a1bSJed Brown args: -mat_type mpiaijcusparse -vec_type cuda 190c4762a1bSJed Brown filter: grep -v type 191c4762a1bSJed Brown output_file: output/ex5_21.out 192c4762a1bSJed Brown requires: cuda 193c4762a1bSJed Brown 194c4762a1bSJed Brown test: 195c4762a1bSJed Brown nsize: 3 196c4762a1bSJed Brown suffix: 2_aijcusparse_2 197c4762a1bSJed Brown filter: grep -v type 198c4762a1bSJed Brown args: -mat_type mpiaijcusparse -vec_type cuda 199bd46da1dSJunchao Zhang args: -sf_type {{basic neighbor}} 200c4762a1bSJed Brown output_file: output/ex5_23.out 201c4762a1bSJed Brown requires: cuda 202c4762a1bSJed Brown 203c4762a1bSJed Brown test: 204c4762a1bSJed Brown nsize: 3 205c4762a1bSJed Brown suffix: 2_aijcusparse_3 206c4762a1bSJed Brown filter: grep -v type 207c4762a1bSJed Brown args: -mat_type mpiaijcusparse -vec_type cuda 208c20d7725SJed Brown args: -sf_type {{basic neighbor}} 209c4762a1bSJed Brown output_file: output/ex5_23.out 210dfd57a17SPierre Jolivet requires: cuda defined(PETSC_HAVE_MPI_GPU_AWARE) 211c4762a1bSJed Brown 212c4762a1bSJed Brown test: 213c4762a1bSJed Brown suffix: 31 214c4762a1bSJed Brown args: -mat_type mpiaij -test_diagonalscale 215c4762a1bSJed Brown filter: grep -v type 216c4762a1bSJed Brown 217c4762a1bSJed Brown test: 218c4762a1bSJed Brown suffix: 32 219c4762a1bSJed Brown args: -mat_type mpibaij -test_diagonalscale 220c4762a1bSJed Brown filter: grep -v Mat_ 221c4762a1bSJed Brown 222c4762a1bSJed Brown test: 223c4762a1bSJed Brown suffix: 33 224c4762a1bSJed Brown nsize: 3 225c4762a1bSJed Brown args: -mat_type mpiaij -test_diagonalscale 226c4762a1bSJed Brown filter: grep -v type 227c4762a1bSJed Brown 228c4762a1bSJed Brown test: 229c4762a1bSJed Brown suffix: 34 230c4762a1bSJed Brown nsize: 3 231c4762a1bSJed Brown args: -mat_type mpibaij -test_diagonalscale 232c4762a1bSJed Brown filter: grep -v Mat_ 233c4762a1bSJed Brown 234c4762a1bSJed Brown test: 235c4762a1bSJed Brown suffix: 3_aijcusparse_1 236c4762a1bSJed Brown args: -mat_type mpiaijcusparse -vec_type cuda -test_diagonalscale 237c4762a1bSJed Brown filter: grep -v type 238c4762a1bSJed Brown output_file: output/ex5_31.out 239c4762a1bSJed Brown requires: cuda 240c4762a1bSJed Brown 241c4762a1bSJed Brown test: 242c4762a1bSJed Brown suffix: 3_aijcusparse_2 243c4762a1bSJed Brown nsize: 3 244c4762a1bSJed Brown args: -mat_type mpiaijcusparse -vec_type cuda -test_diagonalscale 245c4762a1bSJed Brown filter: grep -v type 246c4762a1bSJed Brown output_file: output/ex5_33.out 247c4762a1bSJed Brown requires: cuda 248c4762a1bSJed Brown 249c4762a1bSJed Brown test: 25035990778SJunchao Zhang suffix: 3_kokkos 25135990778SJunchao Zhang nsize: 3 25235990778SJunchao Zhang args: -mat_type mpiaijkokkos -vec_type kokkos -test_diagonalscale 25335990778SJunchao Zhang filter: grep -v type 25435990778SJunchao Zhang output_file: output/ex5_33.out 2553078479eSJunchao Zhang requires: !sycl kokkos_kernels 25635990778SJunchao Zhang 25735990778SJunchao Zhang test: 258c4762a1bSJed Brown suffix: aijcusparse_1 259c4762a1bSJed Brown args: -mat_type seqaijcusparse -vec_type cuda -rectA 260c4762a1bSJed Brown filter: grep -v type 261c4762a1bSJed Brown output_file: output/ex5_11_A.out 262c4762a1bSJed Brown requires: cuda 263c4762a1bSJed Brown 264c4762a1bSJed Brown test: 265c4762a1bSJed Brown suffix: aijcusparse_2 266c4762a1bSJed Brown args: -mat_type seqaijcusparse -vec_type cuda -rectB 267c4762a1bSJed Brown filter: grep -v type 268c4762a1bSJed Brown output_file: output/ex5_11_B.out 269c4762a1bSJed Brown requires: cuda 270c4762a1bSJed Brown 271c4762a1bSJed Brown test: 272c4762a1bSJed Brown suffix: sell_1 273c4762a1bSJed Brown args: -mat_type sell 274c4762a1bSJed Brown output_file: output/ex5_41.out 275c4762a1bSJed Brown 276c4762a1bSJed Brown test: 277c4762a1bSJed Brown suffix: sell_2 278c4762a1bSJed Brown nsize: 3 279c4762a1bSJed Brown args: -mat_type sell 280c4762a1bSJed Brown output_file: output/ex5_43.out 281c4762a1bSJed Brown 282c4762a1bSJed Brown test: 283c4762a1bSJed Brown suffix: sell_3 284c4762a1bSJed Brown args: -mat_type sell -test_diagonalscale 285c4762a1bSJed Brown output_file: output/ex5_51.out 286c4762a1bSJed Brown 287c4762a1bSJed Brown test: 288c4762a1bSJed Brown suffix: sell_4 289c4762a1bSJed Brown nsize: 3 290c4762a1bSJed Brown args: -mat_type sell -test_diagonalscale 291c4762a1bSJed Brown output_file: output/ex5_53.out 292c4762a1bSJed Brown 293c4762a1bSJed Brown TEST*/ 294