xref: /petsc/src/mat/tests/ex103.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
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;
25ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
26ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   /* Creat a matrix */
29c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);CHKERRQ(ierr);
30c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL);CHKERRQ(ierr);
31c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD, &A);CHKERRQ(ierr);
32c4762a1bSJed Brown   ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N);CHKERRQ(ierr);
33c4762a1bSJed Brown   ierr = MatSetType(A,MATDENSE);CHKERRQ(ierr);
34c4762a1bSJed Brown   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
35c4762a1bSJed Brown   ierr = MatSetUp(A);CHKERRQ(ierr);
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   /* Set local matrix entries */
38c4762a1bSJed Brown   ierr = MatGetOwnershipIS(A,&isrows,&iscols);CHKERRQ(ierr);
39c4762a1bSJed Brown   ierr = ISGetLocalSize(isrows,&nrows);CHKERRQ(ierr);
40c4762a1bSJed Brown   ierr = ISGetIndices(isrows,&rows);CHKERRQ(ierr);
41c4762a1bSJed Brown   ierr = ISGetLocalSize(iscols,&ncols);CHKERRQ(ierr);
42c4762a1bSJed Brown   ierr = ISGetIndices(iscols,&cols);CHKERRQ(ierr);
43c4762a1bSJed Brown   ierr = PetscMalloc1(nrows*ncols,&v);CHKERRQ(ierr);
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   }
54c4762a1bSJed Brown   ierr = MatSetValues(A,nrows,rows,ncols,cols,v,INSERT_VALUES);CHKERRQ(ierr);
55c4762a1bSJed Brown   ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
56c4762a1bSJed Brown   ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
578fffc762SJacob Faibussowitsch   //ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%" PetscInt_FMT "] local nrows %" PetscInt_FMT ", ncols %" PetscInt_FMT "\n",rank,nrows,ncols);CHKERRQ(ierr);
58c4762a1bSJed Brown   //ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);CHKERRQ(ierr);
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   /* Test MatSetValues() by converting A to A_elemental */
61c4762a1bSJed Brown   ierr = MatGetType(A,&type);CHKERRQ(ierr);
62c4762a1bSJed Brown   if (size == 1) {
63c4762a1bSJed Brown     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQDENSE,&isDense);CHKERRQ(ierr);
64c4762a1bSJed Brown     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isAIJ);CHKERRQ(ierr);
65c4762a1bSJed Brown   } else {
66c4762a1bSJed Brown     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIDENSE,&isDense);CHKERRQ(ierr);
67c4762a1bSJed Brown     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&isAIJ);CHKERRQ(ierr);
68c4762a1bSJed Brown   }
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   if (isDense || isAIJ) {
71c4762a1bSJed Brown     Mat Aexplicit;
72c4762a1bSJed Brown     ierr = MatConvert(A, MATELEMENTAL, MAT_INITIAL_MATRIX, &A_elemental);CHKERRQ(ierr);
73c4762a1bSJed Brown     ierr = MatComputeOperator(A_elemental,isAIJ ? MATAIJ : MATDENSE,&Aexplicit);CHKERRQ(ierr);
74c4762a1bSJed Brown     ierr = MatMultEqual(Aexplicit,A_elemental,5,&flg);CHKERRQ(ierr);
75*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Aexplicit != A_elemental.");
76c4762a1bSJed Brown     ierr = MatDestroy(&Aexplicit);CHKERRQ(ierr);
77c4762a1bSJed Brown 
78c4762a1bSJed Brown     /* Test MAT_REUSE_MATRIX which is only supported for inplace conversion */
79c4762a1bSJed Brown     ierr = MatConvert(A, MATELEMENTAL, MAT_INPLACE_MATRIX, &A);CHKERRQ(ierr);
80c4762a1bSJed Brown     ierr = MatMultEqual(A_elemental,A,5,&flg);CHKERRQ(ierr);
81*2c71b3e2SJacob Faibussowitsch     PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"A_elemental != A.");
82c4762a1bSJed Brown     ierr = MatDestroy(&A_elemental);CHKERRQ(ierr);
83c4762a1bSJed Brown   }
84c4762a1bSJed Brown 
85c4762a1bSJed Brown   ierr = ISRestoreIndices(isrows,&rows);CHKERRQ(ierr);
86c4762a1bSJed Brown   ierr = ISRestoreIndices(iscols,&cols);CHKERRQ(ierr);
87c4762a1bSJed Brown   ierr = ISDestroy(&isrows);CHKERRQ(ierr);
88c4762a1bSJed Brown   ierr = ISDestroy(&iscols);CHKERRQ(ierr);
89c4762a1bSJed Brown   ierr = PetscFree(v);CHKERRQ(ierr);
90c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
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