1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "Test MatAXPY(), and illustrate how to reduce number of mallocs used during MatSetValues() calls \n\ 3*c4762a1bSJed Brown Matrix C is copied from ~petsc/src/ksp/ksp/tutorials/ex5.c\n\n"; 4*c4762a1bSJed Brown /* 5*c4762a1bSJed Brown Example: ./ex132 -mat_view ascii::ascii_info 6*c4762a1bSJed Brown */ 7*c4762a1bSJed Brown 8*c4762a1bSJed Brown #include <petscmat.h> 9*c4762a1bSJed Brown 10*c4762a1bSJed Brown int main(int argc,char **args) 11*c4762a1bSJed Brown { 12*c4762a1bSJed Brown Mat C,C1,C2; /* matrix */ 13*c4762a1bSJed Brown PetscScalar v; 14*c4762a1bSJed Brown PetscInt Ii,J,Istart,Iend; 15*c4762a1bSJed Brown PetscErrorCode ierr; 16*c4762a1bSJed Brown PetscInt i,j,m = 3,n = 2; 17*c4762a1bSJed Brown PetscMPIInt size,rank; 18*c4762a1bSJed Brown PetscBool mat_nonsymmetric = PETSC_FALSE; 19*c4762a1bSJed Brown MatInfo info; 20*c4762a1bSJed Brown 21*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 22*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr); 23*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 24*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 25*c4762a1bSJed Brown n = 2*size; 26*c4762a1bSJed Brown 27*c4762a1bSJed Brown /* Set flag if we are doing a nonsymmetric problem; the default is symmetric. */ 28*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-mat_nonsym",&mat_nonsymmetric,NULL);CHKERRQ(ierr); 29*c4762a1bSJed Brown 30*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr); 31*c4762a1bSJed Brown ierr = MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);CHKERRQ(ierr); 32*c4762a1bSJed Brown ierr = MatSetFromOptions(C);CHKERRQ(ierr); 33*c4762a1bSJed Brown ierr = MatSeqAIJSetPreallocation(C,5,NULL);CHKERRQ(ierr); 34*c4762a1bSJed Brown 35*c4762a1bSJed Brown ierr = MatGetOwnershipRange(C,&Istart,&Iend);CHKERRQ(ierr); 36*c4762a1bSJed Brown for (Ii=Istart; Ii<Iend; Ii++) { 37*c4762a1bSJed Brown v = -1.0; i = Ii/n; j = Ii - i*n; 38*c4762a1bSJed Brown if (i>0) {J = Ii - n; ierr = MatSetValues(C,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);} 39*c4762a1bSJed Brown if (i<m-1) {J = Ii + n; ierr = MatSetValues(C,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);} 40*c4762a1bSJed Brown if (j>0) {J = Ii - 1; ierr = MatSetValues(C,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);} 41*c4762a1bSJed Brown if (j<n-1) {J = Ii + 1; ierr = MatSetValues(C,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);} 42*c4762a1bSJed Brown v = 4.0; ierr = MatSetValues(C,1,&Ii,1,&Ii,&v,ADD_VALUES);CHKERRQ(ierr); 43*c4762a1bSJed Brown } 44*c4762a1bSJed Brown 45*c4762a1bSJed Brown /* Make the matrix nonsymmetric if desired */ 46*c4762a1bSJed Brown if (mat_nonsymmetric) { 47*c4762a1bSJed Brown for (Ii=Istart; Ii<Iend; Ii++) { 48*c4762a1bSJed Brown v = -1.5; i = Ii/n; 49*c4762a1bSJed Brown if (i>1) {J = Ii-n-1; ierr = MatSetValues(C,1,&Ii,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);} 50*c4762a1bSJed Brown } 51*c4762a1bSJed Brown } else { 52*c4762a1bSJed Brown ierr = MatSetOption(C,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr); 53*c4762a1bSJed Brown ierr = MatSetOption(C,MAT_SYMMETRY_ETERNAL,PETSC_TRUE);CHKERRQ(ierr); 54*c4762a1bSJed Brown } 55*c4762a1bSJed Brown ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 56*c4762a1bSJed Brown ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 57*c4762a1bSJed Brown 58*c4762a1bSJed Brown /* First, create C1 = 2.0*C1 + C, C1 has less non-zeros than C */ 59*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "\ncreate C1 = 2.0*C1 + C, C1 has less non-zeros than C \n");CHKERRQ(ierr); 60*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&C1);CHKERRQ(ierr); 61*c4762a1bSJed Brown ierr = MatSetSizes(C1,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);CHKERRQ(ierr); 62*c4762a1bSJed Brown ierr = MatSetFromOptions(C1);CHKERRQ(ierr); 63*c4762a1bSJed Brown ierr = MatSeqAIJSetPreallocation(C1,1,NULL);CHKERRQ(ierr); 64*c4762a1bSJed Brown for (Ii=Istart; Ii<Iend; Ii++) { 65*c4762a1bSJed Brown ierr = MatSetValues(C1,1,&Ii,1,&Ii,&v,ADD_VALUES);CHKERRQ(ierr); 66*c4762a1bSJed Brown } 67*c4762a1bSJed Brown ierr = MatAssemblyBegin(C1,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 68*c4762a1bSJed Brown ierr = MatAssemblyEnd(C1,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 69*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," MatAXPY(C1,2.0,C,DIFFERENT_NONZERO_PATTERN)...\n");CHKERRQ(ierr); 70*c4762a1bSJed Brown ierr = MatAXPY(C1,2.0,C,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 71*c4762a1bSJed Brown ierr = MatGetInfo(C1,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); 72*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," C1: nz_allocated = %g; nz_used = %g; nz_unneeded = %g\n",info.nz_allocated,info.nz_used, info.nz_unneeded);CHKERRQ(ierr); 73*c4762a1bSJed Brown 74*c4762a1bSJed Brown /* Secondly, create C2 = 2.0*C2 + C, C2 has non-zero pattern of C2 + C */ 75*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD, "\ncreate C2 = 2.0*C2 + C, C2 has non-zero pattern of C2 + C \n");CHKERRQ(ierr); 76*c4762a1bSJed Brown ierr = MatDuplicate(C,MAT_DO_NOT_COPY_VALUES,&C2);CHKERRQ(ierr); 77*c4762a1bSJed Brown 78*c4762a1bSJed Brown for (Ii=Istart; Ii<Iend; Ii++) { 79*c4762a1bSJed Brown v = 1.0; 80*c4762a1bSJed Brown ierr = MatSetValues(C2,1,&Ii,1,&Ii,&v,ADD_VALUES);CHKERRQ(ierr); 81*c4762a1bSJed Brown } 82*c4762a1bSJed Brown ierr = MatAssemblyBegin(C2,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 83*c4762a1bSJed Brown ierr = MatAssemblyEnd(C2,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 84*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," MatAXPY(C2,2.0,C,SUBSET_NONZERO_PATTERN)...\n");CHKERRQ(ierr); 85*c4762a1bSJed Brown ierr = MatAXPY(C2,2.0,C,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 86*c4762a1bSJed Brown ierr = MatGetInfo(C2,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); 87*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD," C2: nz_allocated = %g; nz_used = %g; nz_unneeded = %g\n",info.nz_allocated,info.nz_used, info.nz_unneeded);CHKERRQ(ierr); 88*c4762a1bSJed Brown 89*c4762a1bSJed Brown ierr = MatDestroy(&C1);CHKERRQ(ierr); 90*c4762a1bSJed Brown ierr = MatDestroy(&C2);CHKERRQ(ierr); 91*c4762a1bSJed Brown ierr = MatDestroy(&C);CHKERRQ(ierr); 92*c4762a1bSJed Brown 93*c4762a1bSJed Brown ierr = PetscFinalize(); 94*c4762a1bSJed Brown return ierr; 95*c4762a1bSJed Brown } 96*c4762a1bSJed Brown 97*c4762a1bSJed Brown /*TEST 98*c4762a1bSJed Brown 99*c4762a1bSJed Brown test: 100*c4762a1bSJed Brown 101*c4762a1bSJed Brown TEST*/ 102