1*c4762a1bSJed Brown static char help[] = "Test MatSetValues() by converting MATDENSE to MATELEMENTAL. \n\ 2*c4762a1bSJed Brown Modified from the code contributed by Yaning Liu @lbl.gov \n\n"; 3*c4762a1bSJed Brown /* 4*c4762a1bSJed Brown Example: 5*c4762a1bSJed Brown mpiexec -n <np> ./ex103 6*c4762a1bSJed Brown mpiexec -n <np> ./ex103 -mat_type elemental -mat_view 7*c4762a1bSJed Brown mpiexec -n <np> ./ex103 -mat_type aij 8*c4762a1bSJed Brown */ 9*c4762a1bSJed Brown 10*c4762a1bSJed Brown #include <petscmat.h> 11*c4762a1bSJed Brown 12*c4762a1bSJed Brown int main(int argc, char** argv) 13*c4762a1bSJed Brown { 14*c4762a1bSJed Brown Mat A,A_elemental; 15*c4762a1bSJed Brown PetscInt i,j,M=10,N=5,nrows,ncols; 16*c4762a1bSJed Brown PetscErrorCode ierr; 17*c4762a1bSJed Brown PetscMPIInt rank,size; 18*c4762a1bSJed Brown IS isrows,iscols; 19*c4762a1bSJed Brown const PetscInt *rows,*cols; 20*c4762a1bSJed Brown PetscScalar *v; 21*c4762a1bSJed Brown MatType type; 22*c4762a1bSJed Brown PetscBool isDense,isAIJ,flg; 23*c4762a1bSJed Brown 24*c4762a1bSJed Brown ierr = PetscInitialize(&argc, &argv, (char*)0, help);if (ierr) return ierr; 25*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 26*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 27*c4762a1bSJed Brown 28*c4762a1bSJed Brown /* Creat a matrix */ 29*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);CHKERRQ(ierr); 30*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL);CHKERRQ(ierr); 31*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD, &A);CHKERRQ(ierr); 32*c4762a1bSJed Brown ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr); 33*c4762a1bSJed Brown ierr = MatSetType(A,MATDENSE);CHKERRQ(ierr); 34*c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 35*c4762a1bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 36*c4762a1bSJed Brown 37*c4762a1bSJed Brown /* Set local matrix entries */ 38*c4762a1bSJed Brown ierr = MatGetOwnershipIS(A,&isrows,&iscols);CHKERRQ(ierr); 39*c4762a1bSJed Brown ierr = ISGetLocalSize(isrows,&nrows);CHKERRQ(ierr); 40*c4762a1bSJed Brown ierr = ISGetIndices(isrows,&rows);CHKERRQ(ierr); 41*c4762a1bSJed Brown ierr = ISGetLocalSize(iscols,&ncols);CHKERRQ(ierr); 42*c4762a1bSJed Brown ierr = ISGetIndices(iscols,&cols);CHKERRQ(ierr); 43*c4762a1bSJed Brown ierr = PetscMalloc1(nrows*ncols,&v);CHKERRQ(ierr); 44*c4762a1bSJed Brown 45*c4762a1bSJed Brown for (i=0; i<nrows; i++) { 46*c4762a1bSJed Brown for (j=0; j<ncols; j++) { 47*c4762a1bSJed Brown if (size == 1) { 48*c4762a1bSJed Brown v[i*ncols+j] = (PetscScalar)(i+j); 49*c4762a1bSJed Brown } else { 50*c4762a1bSJed Brown v[i*ncols+j] = (PetscScalar)rank+j*0.1; 51*c4762a1bSJed Brown } 52*c4762a1bSJed Brown } 53*c4762a1bSJed Brown } 54*c4762a1bSJed Brown ierr = MatSetValues(A,nrows,rows,ncols,cols,v,INSERT_VALUES);CHKERRQ(ierr); 55*c4762a1bSJed Brown ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 56*c4762a1bSJed Brown ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 57*c4762a1bSJed Brown //ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%D] local nrows %D, ncols %D\n",rank,nrows,ncols);CHKERRQ(ierr); 58*c4762a1bSJed Brown //ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr); 59*c4762a1bSJed Brown 60*c4762a1bSJed Brown /* Test MatSetValues() by converting A to A_elemental */ 61*c4762a1bSJed Brown ierr = MatGetType(A,&type);CHKERRQ(ierr); 62*c4762a1bSJed Brown if (size == 1) { 63*c4762a1bSJed Brown ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQDENSE,&isDense);CHKERRQ(ierr); 64*c4762a1bSJed Brown ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isAIJ);CHKERRQ(ierr); 65*c4762a1bSJed Brown } else { 66*c4762a1bSJed Brown ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIDENSE,&isDense);CHKERRQ(ierr); 67*c4762a1bSJed Brown ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isAIJ);CHKERRQ(ierr); 68*c4762a1bSJed Brown } 69*c4762a1bSJed Brown 70*c4762a1bSJed Brown if (isDense || isAIJ) { 71*c4762a1bSJed Brown Mat Aexplicit; 72*c4762a1bSJed Brown ierr = MatConvert(A, MATELEMENTAL, MAT_INITIAL_MATRIX, &A_elemental);CHKERRQ(ierr); 73*c4762a1bSJed Brown ierr = MatComputeOperator(A_elemental,isAIJ ? MATAIJ : MATDENSE,&Aexplicit);CHKERRQ(ierr); 74*c4762a1bSJed Brown ierr = MatMultEqual(Aexplicit,A_elemental,5,&flg);CHKERRQ(ierr); 75*c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Aexplicit != A_elemental."); 76*c4762a1bSJed Brown ierr = MatDestroy(&Aexplicit);CHKERRQ(ierr); 77*c4762a1bSJed Brown 78*c4762a1bSJed Brown /* Test MAT_REUSE_MATRIX which is only supported for inplace conversion */ 79*c4762a1bSJed Brown ierr = MatConvert(A, MATELEMENTAL, MAT_INPLACE_MATRIX, &A);CHKERRQ(ierr); 80*c4762a1bSJed Brown ierr = MatMultEqual(A_elemental,A,5,&flg);CHKERRQ(ierr); 81*c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"A_elemental != A."); 82*c4762a1bSJed Brown ierr = MatDestroy(&A_elemental);CHKERRQ(ierr); 83*c4762a1bSJed Brown } 84*c4762a1bSJed Brown 85*c4762a1bSJed Brown ierr = ISRestoreIndices(isrows,&rows);CHKERRQ(ierr); 86*c4762a1bSJed Brown ierr = ISRestoreIndices(iscols,&cols);CHKERRQ(ierr); 87*c4762a1bSJed Brown ierr = ISDestroy(&isrows);CHKERRQ(ierr); 88*c4762a1bSJed Brown ierr = ISDestroy(&iscols);CHKERRQ(ierr); 89*c4762a1bSJed Brown ierr = PetscFree(v);CHKERRQ(ierr); 90*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 91*c4762a1bSJed Brown ierr = PetscFinalize(); 92*c4762a1bSJed Brown return ierr; 93*c4762a1bSJed Brown } 94*c4762a1bSJed Brown 95*c4762a1bSJed Brown 96*c4762a1bSJed Brown /*TEST 97*c4762a1bSJed Brown 98*c4762a1bSJed Brown build: 99*c4762a1bSJed Brown requires: elemental 100*c4762a1bSJed Brown 101*c4762a1bSJed Brown test: 102*c4762a1bSJed Brown nsize: 6 103*c4762a1bSJed Brown 104*c4762a1bSJed Brown test: 105*c4762a1bSJed Brown suffix: 2 106*c4762a1bSJed Brown nsize: 6 107*c4762a1bSJed Brown args: -mat_type aij 108*c4762a1bSJed Brown output_file: output/ex103_1.out 109*c4762a1bSJed Brown 110*c4762a1bSJed Brown test: 111*c4762a1bSJed Brown suffix: 3 112*c4762a1bSJed Brown nsize: 6 113*c4762a1bSJed Brown args: -mat_type elemental 114*c4762a1bSJed Brown output_file: output/ex103_1.out 115*c4762a1bSJed Brown 116*c4762a1bSJed Brown TEST*/ 117