xref: /petsc/src/mat/tests/ex103.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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   PetscMPIInt    rank,size;
17c4762a1bSJed Brown   IS             isrows,iscols;
18c4762a1bSJed Brown   const PetscInt *rows,*cols;
19c4762a1bSJed Brown   PetscScalar    *v;
20c4762a1bSJed Brown   MatType        type;
21c4762a1bSJed Brown   PetscBool      isDense,isAIJ,flg;
22c4762a1bSJed Brown 
23*327415f7SBarry Smith   PetscFunctionBeginUser;
249566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char*)0, help));
259566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
269566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   /* Creat a matrix */
299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL));
309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL));
319566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
329566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N));
339566063dSJacob Faibussowitsch   PetscCall(MatSetType(A,MATDENSE));
349566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
359566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   /* Set local matrix entries */
389566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipIS(A,&isrows,&iscols));
399566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(isrows,&nrows));
409566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(isrows,&rows));
419566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(iscols,&ncols));
429566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(iscols,&cols));
439566063dSJacob Faibussowitsch   PetscCall(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   }
549566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A,nrows,rows,ncols,cols,v,INSERT_VALUES));
559566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
569566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
579566063dSJacob Faibussowitsch   //PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%" PetscInt_FMT "] local nrows %" PetscInt_FMT ", ncols %" PetscInt_FMT "\n",rank,nrows,ncols));
589566063dSJacob Faibussowitsch   //PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT));
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   /* Test MatSetValues() by converting A to A_elemental */
619566063dSJacob Faibussowitsch   PetscCall(MatGetType(A,&type));
62c4762a1bSJed Brown   if (size == 1) {
639566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)A,MATSEQDENSE,&isDense));
649566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isAIJ));
65c4762a1bSJed Brown   } else {
669566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)A,MATMPIDENSE,&isDense));
679566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isAIJ));
68c4762a1bSJed Brown   }
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   if (isDense || isAIJ) {
71c4762a1bSJed Brown     Mat Aexplicit;
729566063dSJacob Faibussowitsch     PetscCall(MatConvert(A, MATELEMENTAL, MAT_INITIAL_MATRIX, &A_elemental));
739566063dSJacob Faibussowitsch     PetscCall(MatComputeOperator(A_elemental,isAIJ ? MATAIJ : MATDENSE,&Aexplicit));
749566063dSJacob Faibussowitsch     PetscCall(MatMultEqual(Aexplicit,A_elemental,5,&flg));
7528b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Aexplicit != A_elemental.");
769566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Aexplicit));
77c4762a1bSJed Brown 
78c4762a1bSJed Brown     /* Test MAT_REUSE_MATRIX which is only supported for inplace conversion */
799566063dSJacob Faibussowitsch     PetscCall(MatConvert(A, MATELEMENTAL, MAT_INPLACE_MATRIX, &A));
809566063dSJacob Faibussowitsch     PetscCall(MatMultEqual(A_elemental,A,5,&flg));
8128b400f6SJacob Faibussowitsch     PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"A_elemental != A.");
829566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&A_elemental));
83c4762a1bSJed Brown   }
84c4762a1bSJed Brown 
859566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(isrows,&rows));
869566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(iscols,&cols));
879566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isrows));
889566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iscols));
899566063dSJacob Faibussowitsch   PetscCall(PetscFree(v));
909566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
919566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
92b122ec5aSJacob Faibussowitsch   return 0;
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