xref: /petsc/src/mat/tests/ex103.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
1c4762a1bSJed Brown static char help[] = "Test MatSetValues() by converting MATDENSE to MATELEMENTAL. \n\
2c4762a1bSJed Brown Modified from the code contributed by Yaning Liu @lbl.gov \n\n";
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown  Example:
5c4762a1bSJed Brown    mpiexec -n <np> ./ex103
6c4762a1bSJed Brown    mpiexec -n <np> ./ex103 -mat_type elemental -mat_view
7c4762a1bSJed Brown    mpiexec -n <np> ./ex103 -mat_type aij
8c4762a1bSJed Brown */
9c4762a1bSJed Brown 
10c4762a1bSJed Brown #include <petscmat.h>
11c4762a1bSJed Brown 
12c4762a1bSJed Brown int main(int argc, char** argv)
13c4762a1bSJed Brown {
14c4762a1bSJed Brown   Mat            A,A_elemental;
15c4762a1bSJed Brown   PetscInt       i,j,M=10,N=5,nrows,ncols;
16c4762a1bSJed Brown   PetscErrorCode ierr;
17c4762a1bSJed Brown   PetscMPIInt    rank,size;
18c4762a1bSJed Brown   IS             isrows,iscols;
19c4762a1bSJed Brown   const PetscInt *rows,*cols;
20c4762a1bSJed Brown   PetscScalar    *v;
21c4762a1bSJed Brown   MatType        type;
22c4762a1bSJed Brown   PetscBool      isDense,isAIJ,flg;
23c4762a1bSJed Brown 
24c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, (char*)0, help);if (ierr) return ierr;
255f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
265f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   /* Creat a matrix */
295f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD, &A));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N));
335f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATDENSE));
345f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
355f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(A));
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   /* Set local matrix entries */
385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipIS(A,&isrows,&iscols));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetLocalSize(isrows,&nrows));
405f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetIndices(isrows,&rows));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetLocalSize(iscols,&ncols));
425f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetIndices(iscols,&cols));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(nrows*ncols,&v));
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   for (i=0; i<nrows; i++) {
46c4762a1bSJed Brown     for (j=0; j<ncols; j++) {
47c4762a1bSJed Brown       if (size == 1) {
48c4762a1bSJed Brown         v[i*ncols+j] = (PetscScalar)(i+j);
49c4762a1bSJed Brown       } else {
50c4762a1bSJed Brown         v[i*ncols+j] = (PetscScalar)rank+j*0.1;
51c4762a1bSJed Brown       }
52c4762a1bSJed Brown     }
53c4762a1bSJed Brown   }
545f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,nrows,rows,ncols,cols,v,INSERT_VALUES));
555f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
575f80ce2aSJacob Faibussowitsch   //CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%" PetscInt_FMT "] local nrows %" PetscInt_FMT ", ncols %" PetscInt_FMT "\n",rank,nrows,ncols));
585f80ce2aSJacob Faibussowitsch   //CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   /* Test MatSetValues() by converting A to A_elemental */
615f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetType(A,&type));
62c4762a1bSJed Brown   if (size == 1) {
635f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATSEQDENSE,&isDense));
645f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isAIJ));
65c4762a1bSJed Brown   } else {
665f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATMPIDENSE,&isDense));
675f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isAIJ));
68c4762a1bSJed Brown   }
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   if (isDense || isAIJ) {
71c4762a1bSJed Brown     Mat Aexplicit;
725f80ce2aSJacob Faibussowitsch     CHKERRQ(MatConvert(A, MATELEMENTAL, MAT_INITIAL_MATRIX, &A_elemental));
735f80ce2aSJacob Faibussowitsch     CHKERRQ(MatComputeOperator(A_elemental,isAIJ ? MATAIJ : MATDENSE,&Aexplicit));
745f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultEqual(Aexplicit,A_elemental,5,&flg));
75*28b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Aexplicit != A_elemental.");
765f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&Aexplicit));
77c4762a1bSJed Brown 
78c4762a1bSJed Brown     /* Test MAT_REUSE_MATRIX which is only supported for inplace conversion */
795f80ce2aSJacob Faibussowitsch     CHKERRQ(MatConvert(A, MATELEMENTAL, MAT_INPLACE_MATRIX, &A));
805f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultEqual(A_elemental,A,5,&flg));
81*28b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"A_elemental != A.");
825f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&A_elemental));
83c4762a1bSJed Brown   }
84c4762a1bSJed Brown 
855f80ce2aSJacob Faibussowitsch   CHKERRQ(ISRestoreIndices(isrows,&rows));
865f80ce2aSJacob Faibussowitsch   CHKERRQ(ISRestoreIndices(iscols,&cols));
875f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isrows));
885f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&iscols));
895f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(v));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
91c4762a1bSJed Brown   ierr = PetscFinalize();
92c4762a1bSJed Brown   return ierr;
93c4762a1bSJed Brown }
94c4762a1bSJed Brown 
95c4762a1bSJed Brown /*TEST
96c4762a1bSJed Brown 
97c4762a1bSJed Brown    build:
98c4762a1bSJed Brown       requires: elemental
99c4762a1bSJed Brown 
100c4762a1bSJed Brown    test:
101c4762a1bSJed Brown       nsize: 6
102c4762a1bSJed Brown 
103c4762a1bSJed Brown    test:
104c4762a1bSJed Brown       suffix: 2
105c4762a1bSJed Brown       nsize: 6
106c4762a1bSJed Brown       args: -mat_type aij
107c4762a1bSJed Brown       output_file: output/ex103_1.out
108c4762a1bSJed Brown 
109c4762a1bSJed Brown    test:
110c4762a1bSJed Brown       suffix: 3
111c4762a1bSJed Brown       nsize: 6
112c4762a1bSJed Brown       args: -mat_type elemental
113c4762a1bSJed Brown       output_file: output/ex103_1.out
114c4762a1bSJed Brown 
115c4762a1bSJed Brown TEST*/
116