1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "Tests MatTranspose(), MatNorm(), and MatAXPY().\n\n"; 3*c4762a1bSJed Brown 4*c4762a1bSJed Brown #include <petscmat.h> 5*c4762a1bSJed Brown 6*c4762a1bSJed Brown int main(int argc,char **argv) 7*c4762a1bSJed Brown { 8*c4762a1bSJed Brown Mat mat,tmat = 0; 9*c4762a1bSJed Brown PetscInt m = 4,n,i,j; 10*c4762a1bSJed Brown PetscErrorCode ierr; 11*c4762a1bSJed Brown PetscMPIInt size,rank; 12*c4762a1bSJed Brown PetscInt rstart,rend,rect = 0; 13*c4762a1bSJed Brown PetscBool flg; 14*c4762a1bSJed Brown PetscScalar v; 15*c4762a1bSJed Brown PetscReal normf,normi,norm1; 16*c4762a1bSJed Brown MatInfo info; 17*c4762a1bSJed Brown 18*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 19*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr); 20*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 21*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 22*c4762a1bSJed Brown n = m; 23*c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL,"-rect1",&flg);CHKERRQ(ierr); 24*c4762a1bSJed Brown if (flg) {n += 2; rect = 1;} 25*c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL,"-rect2",&flg);CHKERRQ(ierr); 26*c4762a1bSJed Brown if (flg) {n -= 2; rect = 1;} 27*c4762a1bSJed Brown 28*c4762a1bSJed Brown /* Create and assemble matrix */ 29*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&mat);CHKERRQ(ierr); 30*c4762a1bSJed Brown ierr = MatSetSizes(mat,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr); 31*c4762a1bSJed Brown ierr = MatSetFromOptions(mat);CHKERRQ(ierr); 32*c4762a1bSJed Brown ierr = MatSetUp(mat);CHKERRQ(ierr); 33*c4762a1bSJed Brown ierr = MatGetOwnershipRange(mat,&rstart,&rend);CHKERRQ(ierr); 34*c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 35*c4762a1bSJed Brown for (j=0; j<n; j++) { 36*c4762a1bSJed Brown v = 10*i+j; 37*c4762a1bSJed Brown ierr = MatSetValues(mat,1,&i,1,&j,&v,INSERT_VALUES);CHKERRQ(ierr); 38*c4762a1bSJed Brown } 39*c4762a1bSJed Brown } 40*c4762a1bSJed Brown ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 41*c4762a1bSJed Brown ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 42*c4762a1bSJed Brown 43*c4762a1bSJed Brown /* Print info about original matrix */ 44*c4762a1bSJed Brown ierr = MatGetInfo(mat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); 45*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"original matrix nonzeros = %D, allocated nonzeros = %D\n", 46*c4762a1bSJed Brown (PetscInt)info.nz_used,(PetscInt)info.nz_allocated);CHKERRQ(ierr); 47*c4762a1bSJed Brown ierr = MatNorm(mat,NORM_FROBENIUS,&normf);CHKERRQ(ierr); 48*c4762a1bSJed Brown ierr = MatNorm(mat,NORM_1,&norm1);CHKERRQ(ierr); 49*c4762a1bSJed Brown ierr = MatNorm(mat,NORM_INFINITY,&normi);CHKERRQ(ierr); 50*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"original: Frobenious norm = %g, one norm = %g, infinity norm = %g\n",(double)normf,(double)norm1,(double)normi);CHKERRQ(ierr); 51*c4762a1bSJed Brown ierr = MatView(mat,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 52*c4762a1bSJed Brown 53*c4762a1bSJed Brown /* Form matrix transpose */ 54*c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL,"-in_place",&flg);CHKERRQ(ierr); 55*c4762a1bSJed Brown if (flg) { 56*c4762a1bSJed Brown ierr = MatTranspose(mat,MAT_INPLACE_MATRIX,&mat);CHKERRQ(ierr); /* in-place transpose */ 57*c4762a1bSJed Brown tmat = mat; mat = 0; 58*c4762a1bSJed Brown } else { /* out-of-place transpose */ 59*c4762a1bSJed Brown ierr = MatTranspose(mat,MAT_INITIAL_MATRIX,&tmat);CHKERRQ(ierr); 60*c4762a1bSJed Brown } 61*c4762a1bSJed Brown 62*c4762a1bSJed Brown /* Print info about transpose matrix */ 63*c4762a1bSJed Brown ierr = MatGetInfo(tmat,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); 64*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"transpose matrix nonzeros = %D, allocated nonzeros = %D\n", 65*c4762a1bSJed Brown (PetscInt)info.nz_used,(PetscInt)info.nz_allocated);CHKERRQ(ierr); 66*c4762a1bSJed Brown ierr = MatNorm(tmat,NORM_FROBENIUS,&normf);CHKERRQ(ierr); 67*c4762a1bSJed Brown ierr = MatNorm(tmat,NORM_1,&norm1);CHKERRQ(ierr); 68*c4762a1bSJed Brown ierr = MatNorm(tmat,NORM_INFINITY,&normi);CHKERRQ(ierr); 69*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"transpose: Frobenious norm = %g, one norm = %g, infinity norm = %g\n",(double)normf,(double)norm1,(double)normi);CHKERRQ(ierr); 70*c4762a1bSJed Brown ierr = MatView(tmat,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 71*c4762a1bSJed Brown 72*c4762a1bSJed Brown 73*c4762a1bSJed Brown /* Test MatAXPY */ 74*c4762a1bSJed Brown if (mat && !rect) { 75*c4762a1bSJed Brown PetscScalar alpha = 1.0; 76*c4762a1bSJed Brown ierr = PetscOptionsGetScalar(NULL,NULL,"-alpha",&alpha,NULL);CHKERRQ(ierr); 77*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"matrix addition: B = B + alpha * A\n");CHKERRQ(ierr); 78*c4762a1bSJed Brown ierr = MatAXPY(tmat,alpha,mat,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 79*c4762a1bSJed Brown ierr = MatView(tmat,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 80*c4762a1bSJed Brown } 81*c4762a1bSJed Brown 82*c4762a1bSJed Brown /* Free data structures */ 83*c4762a1bSJed Brown ierr = MatDestroy(&tmat);CHKERRQ(ierr); 84*c4762a1bSJed Brown if (mat) {ierr = MatDestroy(&mat);CHKERRQ(ierr);} 85*c4762a1bSJed Brown 86*c4762a1bSJed Brown ierr = PetscFinalize(); 87*c4762a1bSJed Brown return ierr; 88*c4762a1bSJed Brown } 89*c4762a1bSJed Brown 90*c4762a1bSJed Brown /*TEST 91*c4762a1bSJed Brown 92*c4762a1bSJed Brown test: 93*c4762a1bSJed Brown 94*c4762a1bSJed Brown testset: 95*c4762a1bSJed Brown args: -rect1 96*c4762a1bSJed Brown test: 97*c4762a1bSJed Brown suffix: r1 98*c4762a1bSJed Brown output_file: output/ex49_r1.out 99*c4762a1bSJed Brown test: 100*c4762a1bSJed Brown suffix: r1_inplace 101*c4762a1bSJed Brown args: -in_place 102*c4762a1bSJed Brown output_file: output/ex49_r1.out 103*c4762a1bSJed Brown test: 104*c4762a1bSJed Brown suffix: r1_par 105*c4762a1bSJed Brown nsize: 2 106*c4762a1bSJed Brown output_file: output/ex49_r1_par.out 107*c4762a1bSJed Brown test: 108*c4762a1bSJed Brown suffix: r1_par_inplace 109*c4762a1bSJed Brown args: -in_place 110*c4762a1bSJed Brown nsize: 2 111*c4762a1bSJed Brown output_file: output/ex49_r1_par.out 112*c4762a1bSJed Brown 113*c4762a1bSJed Brown TEST*/ 114