xref: /petsc/src/mat/tests/ex240.c (revision 1e1ea65d8de51fde77ce8a787efbef25e407badc)
1c4762a1bSJed Brown static char help[] ="Tests MatFDColoringSetValues()\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscdm.h>
4c4762a1bSJed Brown #include <petscdmda.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown int main(int argc,char **argv)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   DM                     da;
9c4762a1bSJed Brown   PetscErrorCode         ierr;
10c4762a1bSJed Brown   PetscInt               N, mx = 5,my = 4,i,j,nc,nrow,n,ncols,rstart,*colors,*map;
11c4762a1bSJed Brown   const PetscInt         *cols;
12c4762a1bSJed Brown   const PetscScalar      *vals;
13c4762a1bSJed Brown   Mat                    A,B;
14c4762a1bSJed Brown   PetscReal              norm;
15c4762a1bSJed Brown   MatFDColoring          fdcoloring;
16c4762a1bSJed Brown   ISColoring             iscoloring;
17c4762a1bSJed Brown   PetscScalar            *cm;
18c4762a1bSJed Brown   const ISColoringValue  *icolors;
19c4762a1bSJed Brown   PetscMPIInt            rank;
20c4762a1bSJed Brown   ISLocalToGlobalMapping ltog;
21c4762a1bSJed Brown   PetscBool              single,two;
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
24ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
25c4762a1bSJed Brown   ierr = DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,mx,my,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da);CHKERRQ(ierr);
26c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
27c4762a1bSJed Brown   ierr = DMCreateMatrix(da,&A);CHKERRQ(ierr);
28c4762a1bSJed Brown 
29c4762a1bSJed Brown   /* as a test copy the matrices from the standard format to the compressed format; this is not scalable but is ok because just for testing */
30c4762a1bSJed Brown   /*    first put the coloring in the global ordering */
31c4762a1bSJed Brown   ierr = DMCreateColoring(da,IS_COLORING_LOCAL,&iscoloring);CHKERRQ(ierr);
32c4762a1bSJed Brown   ierr = ISColoringGetColors(iscoloring,&n,&nc,&icolors);CHKERRQ(ierr);
33c4762a1bSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
34c4762a1bSJed Brown   ierr = PetscMalloc1(n,&map);CHKERRQ(ierr);
35c4762a1bSJed Brown   for (i=0; i<n; i++) map[i] = i;
36c4762a1bSJed Brown   ierr = ISLocalToGlobalMappingApply(ltog,n,map,map);CHKERRQ(ierr);
37c4762a1bSJed Brown   ierr = MatGetSize(A,&N,NULL);CHKERRQ(ierr);
38c4762a1bSJed Brown   ierr = PetscMalloc1(N,&colors);CHKERRQ(ierr);
39c4762a1bSJed Brown   for (i=0; i<N; i++) colors[i] = -1;
40c4762a1bSJed Brown   for (i=0; i<n; i++) colors[map[i]]= icolors[i];
41c4762a1bSJed Brown   ierr = PetscFree(map);CHKERRQ(ierr);
42*1e1ea65dSPierre Jolivet   ierr = PetscSynchronizedPrintf(MPI_COMM_WORLD,"[%d]Global colors \n",rank);CHKERRQ(ierr);
43c4762a1bSJed Brown   for (i=0; i<N; i++) {ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%D %D\n",i,colors[i]);CHKERRQ(ierr);}
44c4762a1bSJed Brown   ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,stdout);CHKERRQ(ierr);
45c4762a1bSJed Brown 
46c4762a1bSJed Brown   /*   second, compress the A matrix */
47c4762a1bSJed Brown   ierr = MatSetRandom(A,NULL);CHKERRQ(ierr);
48c4762a1bSJed Brown   ierr = MatView(A,NULL);CHKERRQ(ierr);
49c4762a1bSJed Brown   ierr = MatGetLocalSize(A,&nrow,NULL);CHKERRQ(ierr);
50c4762a1bSJed Brown   ierr = MatGetOwnershipRange(A,&rstart,NULL);CHKERRQ(ierr);
51c4762a1bSJed Brown   ierr = PetscCalloc1(nrow*nc,&cm);CHKERRQ(ierr);
52c4762a1bSJed Brown   for (i=0; i<nrow; i++) {
53c4762a1bSJed Brown     ierr = MatGetRow(A,rstart+i,&ncols,&cols,&vals);CHKERRQ(ierr);
54c4762a1bSJed Brown     for (j=0; j<ncols; j++) {
55c4762a1bSJed Brown       if (colors[cols[j]] < 0) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Global column %D had no color",cols[j]);
56c4762a1bSJed Brown       cm[i + nrow*colors[cols[j]]] = vals[j];
57c4762a1bSJed Brown     }
58c4762a1bSJed Brown     ierr = MatRestoreRow(A,rstart+i,&ncols,&cols,&vals);CHKERRQ(ierr);
59c4762a1bSJed Brown   }
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   /* print compressed matrix */
62*1e1ea65dSPierre Jolivet   ierr = PetscSynchronizedPrintf(MPI_COMM_WORLD,"[%d] Compressed matrix \n",rank);CHKERRQ(ierr);
63c4762a1bSJed Brown   for (i=0; i<nrow; i++) {
64c4762a1bSJed Brown     for (j=0; j<nc; j++) {
65c4762a1bSJed Brown       ierr = PetscSynchronizedPrintf(MPI_COMM_WORLD,"%12.4e  ",cm[i+nrow*j]);CHKERRQ(ierr);
66c4762a1bSJed Brown     }
67*1e1ea65dSPierre Jolivet     ierr = PetscSynchronizedPrintf(MPI_COMM_WORLD,"\n");CHKERRQ(ierr);
68c4762a1bSJed Brown   }
69c4762a1bSJed Brown   ierr = PetscSynchronizedFlush(PETSC_COMM_WORLD,stdout);CHKERRQ(ierr);
70c4762a1bSJed Brown 
71c4762a1bSJed Brown   /* put the compressed matrix into the standard matrix */
72c4762a1bSJed Brown   ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
73c4762a1bSJed Brown   ierr = MatZeroEntries(A);CHKERRQ(ierr);
74c4762a1bSJed Brown   ierr = MatView(B,0);CHKERRQ(ierr);
75c4762a1bSJed Brown   ierr = MatFDColoringCreate(A,iscoloring,&fdcoloring);CHKERRQ(ierr);
76c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-single_block",&single);CHKERRQ(ierr);
77c4762a1bSJed Brown   if (single) {
78c4762a1bSJed Brown     ierr = MatFDColoringSetBlockSize(fdcoloring,PETSC_DEFAULT,nc);CHKERRQ(ierr);
79c4762a1bSJed Brown   }
80c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-two_block",&two);CHKERRQ(ierr);
81c4762a1bSJed Brown   if (two) {
82c4762a1bSJed Brown     ierr = MatFDColoringSetBlockSize(fdcoloring,PETSC_DEFAULT,2);CHKERRQ(ierr);
83c4762a1bSJed Brown   }
84c4762a1bSJed Brown   ierr = MatFDColoringSetFromOptions(fdcoloring);CHKERRQ(ierr);
85c4762a1bSJed Brown   ierr = MatFDColoringSetUp(A,iscoloring,fdcoloring);CHKERRQ(ierr);
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   ierr = MatFDColoringSetValues(A,fdcoloring,cm);CHKERRQ(ierr);
88c4762a1bSJed Brown   ierr = MatView(A,NULL);CHKERRQ(ierr);
89c4762a1bSJed Brown 
90c4762a1bSJed Brown   /* check the values were put in the correct locations */
91c4762a1bSJed Brown   ierr = MatAXPY(A,-1.0,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
92c4762a1bSJed Brown   ierr = MatView(A,NULL);CHKERRQ(ierr);
93c4762a1bSJed Brown   ierr = MatNorm(A,NORM_FROBENIUS,&norm);CHKERRQ(ierr);
94c4762a1bSJed Brown   if (norm > PETSC_MACHINE_EPSILON) {
95*1e1ea65dSPierre Jolivet     ierr = PetscPrintf(PETSC_COMM_WORLD,"Matrix is not identical, problem with MatFDColoringSetValues()\n");CHKERRQ(ierr);
96c4762a1bSJed Brown   }
97c4762a1bSJed Brown   ierr = PetscFree(colors);CHKERRQ(ierr);
98c4762a1bSJed Brown   ierr = PetscFree(cm);CHKERRQ(ierr);
99c4762a1bSJed Brown   ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr);
100*1e1ea65dSPierre Jolivet   ierr = MatFDColoringDestroy(&fdcoloring);CHKERRQ(ierr);
101c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
102c4762a1bSJed Brown   ierr = MatDestroy(&B);CHKERRQ(ierr);
103c4762a1bSJed Brown   ierr = DMDestroy(&da);CHKERRQ(ierr);
104c4762a1bSJed Brown   ierr = PetscFinalize();
105c4762a1bSJed Brown   return ierr;
106c4762a1bSJed Brown }
107c4762a1bSJed Brown 
108c4762a1bSJed Brown /*TEST
109c4762a1bSJed Brown 
110c4762a1bSJed Brown    test:
111c4762a1bSJed Brown       nsize: 2
112c4762a1bSJed Brown       requires: !complex
113c4762a1bSJed Brown 
114c4762a1bSJed Brown    test:
115c4762a1bSJed Brown       suffix: single
116c4762a1bSJed Brown       requires: !complex
117c4762a1bSJed Brown       nsize: 2
118c4762a1bSJed Brown       args: -single_block
119c4762a1bSJed Brown       output_file: output/ex240_1.out
120c4762a1bSJed Brown 
121c4762a1bSJed Brown    test:
122c4762a1bSJed Brown       suffix: two
123c4762a1bSJed Brown       requires: !complex
124c4762a1bSJed Brown       nsize: 2
125c4762a1bSJed Brown       args: -two_block
126c4762a1bSJed Brown       output_file: output/ex240_1.out
127c4762a1bSJed Brown 
128c4762a1bSJed Brown    test:
129c4762a1bSJed Brown       suffix: 2
130c4762a1bSJed Brown       requires: !complex
131c4762a1bSJed Brown       nsize: 5
132c4762a1bSJed Brown 
133c4762a1bSJed Brown TEST*/
134