1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatTranspose(), MatNorm(), and MatAXPY().\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown int main(int argc,char **argv) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown Mat mat,tmat = 0; 9c4762a1bSJed Brown PetscInt m = 4,n,i,j; 10c4762a1bSJed Brown PetscMPIInt size,rank; 11c4762a1bSJed Brown PetscInt rstart,rend,rect = 0; 12c4762a1bSJed Brown PetscBool flg; 13c4762a1bSJed Brown PetscScalar v; 14c4762a1bSJed Brown PetscReal normf,normi,norm1; 15c4762a1bSJed Brown MatInfo info; 16c4762a1bSJed Brown 17*327415f7SBarry Smith PetscFunctionBeginUser; 189566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 199566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 22c4762a1bSJed Brown n = m; 239566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-rect1",&flg)); 24c4762a1bSJed Brown if (flg) {n += 2; rect = 1;} 259566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-rect2",&flg)); 26c4762a1bSJed Brown if (flg) {n -= 2; rect = 1;} 27c4762a1bSJed Brown 28c4762a1bSJed Brown /* Create and assemble matrix */ 299566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&mat)); 309566063dSJacob Faibussowitsch PetscCall(MatSetSizes(mat,PETSC_DECIDE,PETSC_DECIDE,m,n)); 319566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(mat)); 329566063dSJacob Faibussowitsch PetscCall(MatSetUp(mat)); 339566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(mat,&rstart,&rend)); 34c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 35c4762a1bSJed Brown for (j=0; j<n; j++) { 36c4762a1bSJed Brown v = 10*i+j; 379566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat,1,&i,1,&j,&v,INSERT_VALUES)); 38c4762a1bSJed Brown } 39c4762a1bSJed Brown } 409566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY)); 419566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY)); 42c4762a1bSJed Brown 43c4762a1bSJed Brown /* Print info about original matrix */ 449566063dSJacob Faibussowitsch PetscCall(MatGetInfo(mat,MAT_GLOBAL_SUM,&info)); 45d0609cedSBarry Smith PetscCall(PetscPrintf(PETSC_COMM_WORLD,"original matrix nonzeros = %" PetscInt_FMT ", allocated nonzeros = %" PetscInt_FMT "\n",(PetscInt)info.nz_used,(PetscInt)info.nz_allocated)); 469566063dSJacob Faibussowitsch PetscCall(MatNorm(mat,NORM_FROBENIUS,&normf)); 479566063dSJacob Faibussowitsch PetscCall(MatNorm(mat,NORM_1,&norm1)); 489566063dSJacob Faibussowitsch PetscCall(MatNorm(mat,NORM_INFINITY,&normi)); 499566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"original: Frobenious norm = %g, one norm = %g, infinity norm = %g\n",(double)normf,(double)norm1,(double)normi)); 509566063dSJacob Faibussowitsch PetscCall(MatView(mat,PETSC_VIEWER_STDOUT_WORLD)); 51c4762a1bSJed Brown 52c4762a1bSJed Brown /* Form matrix transpose */ 539566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-in_place",&flg)); 54c4762a1bSJed Brown if (flg) { 559566063dSJacob Faibussowitsch PetscCall(MatTranspose(mat,MAT_INPLACE_MATRIX,&mat)); /* in-place transpose */ 56c4762a1bSJed Brown tmat = mat; mat = 0; 57c4762a1bSJed Brown } else { /* out-of-place transpose */ 589566063dSJacob Faibussowitsch PetscCall(MatTranspose(mat,MAT_INITIAL_MATRIX,&tmat)); 59c4762a1bSJed Brown } 60c4762a1bSJed Brown 61c4762a1bSJed Brown /* Print info about transpose matrix */ 629566063dSJacob Faibussowitsch PetscCall(MatGetInfo(tmat,MAT_GLOBAL_SUM,&info)); 63d0609cedSBarry Smith PetscCall(PetscPrintf(PETSC_COMM_WORLD,"transpose matrix nonzeros = %" PetscInt_FMT ", allocated nonzeros = %" PetscInt_FMT "\n",(PetscInt)info.nz_used,(PetscInt)info.nz_allocated)); 649566063dSJacob Faibussowitsch PetscCall(MatNorm(tmat,NORM_FROBENIUS,&normf)); 659566063dSJacob Faibussowitsch PetscCall(MatNorm(tmat,NORM_1,&norm1)); 669566063dSJacob Faibussowitsch PetscCall(MatNorm(tmat,NORM_INFINITY,&normi)); 679566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"transpose: Frobenious norm = %g, one norm = %g, infinity norm = %g\n",(double)normf,(double)norm1,(double)normi)); 689566063dSJacob Faibussowitsch PetscCall(MatView(tmat,PETSC_VIEWER_STDOUT_WORLD)); 69c4762a1bSJed Brown 70c4762a1bSJed Brown /* Test MatAXPY */ 71c4762a1bSJed Brown if (mat && !rect) { 72c4762a1bSJed Brown PetscScalar alpha = 1.0; 739566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL,NULL,"-alpha",&alpha,NULL)); 749566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"matrix addition: B = B + alpha * A\n")); 759566063dSJacob Faibussowitsch PetscCall(MatAXPY(tmat,alpha,mat,DIFFERENT_NONZERO_PATTERN)); 769566063dSJacob Faibussowitsch PetscCall(MatView(tmat,PETSC_VIEWER_STDOUT_WORLD)); 77c4762a1bSJed Brown } 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* Free data structures */ 809566063dSJacob Faibussowitsch PetscCall(MatDestroy(&tmat)); 819566063dSJacob Faibussowitsch if (mat) PetscCall(MatDestroy(&mat)); 82c4762a1bSJed Brown 839566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 84b122ec5aSJacob Faibussowitsch return 0; 85c4762a1bSJed Brown } 86c4762a1bSJed Brown 87c4762a1bSJed Brown /*TEST 88c4762a1bSJed Brown 89c4762a1bSJed Brown test: 90c4762a1bSJed Brown 91c4762a1bSJed Brown testset: 92c4762a1bSJed Brown args: -rect1 93c4762a1bSJed Brown test: 94c4762a1bSJed Brown suffix: r1 95c4762a1bSJed Brown output_file: output/ex49_r1.out 96c4762a1bSJed Brown test: 97c4762a1bSJed Brown suffix: r1_inplace 98c4762a1bSJed Brown args: -in_place 99c4762a1bSJed Brown output_file: output/ex49_r1.out 100c4762a1bSJed Brown test: 101c4762a1bSJed Brown suffix: r1_par 102c4762a1bSJed Brown nsize: 2 103c4762a1bSJed Brown output_file: output/ex49_r1_par.out 104c4762a1bSJed Brown test: 105c4762a1bSJed Brown suffix: r1_par_inplace 106c4762a1bSJed Brown args: -in_place 107c4762a1bSJed Brown nsize: 2 108c4762a1bSJed Brown output_file: output/ex49_r1_par.out 109c4762a1bSJed Brown 110c4762a1bSJed Brown TEST*/ 111