xref: /petsc/src/mat/tests/ex240.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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   PetscInt               N, mx = 5,my = 4,i,j,nc,nrow,n,ncols,rstart,*colors,*map;
10c4762a1bSJed Brown   const PetscInt         *cols;
11c4762a1bSJed Brown   const PetscScalar      *vals;
12c4762a1bSJed Brown   Mat                    A,B;
13c4762a1bSJed Brown   PetscReal              norm;
14c4762a1bSJed Brown   MatFDColoring          fdcoloring;
15c4762a1bSJed Brown   ISColoring             iscoloring;
16c4762a1bSJed Brown   PetscScalar            *cm;
17c4762a1bSJed Brown   const ISColoringValue  *icolors;
18c4762a1bSJed Brown   PetscMPIInt            rank;
19c4762a1bSJed Brown   ISLocalToGlobalMapping ltog;
20c4762a1bSJed Brown   PetscBool              single,two;
21c4762a1bSJed Brown 
22*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,NULL,help));
235f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
245f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,mx,my,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da));
255f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
265f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateMatrix(da,&A));
27c4762a1bSJed Brown 
28c4762a1bSJed 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 */
29c4762a1bSJed Brown   /*    first put the coloring in the global ordering */
305f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateColoring(da,IS_COLORING_LOCAL,&iscoloring));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(ISColoringGetColors(iscoloring,&n,&nc,&icolors));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalToGlobalMapping(da,&ltog));
335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(n,&map));
34c4762a1bSJed Brown   for (i=0; i<n; i++) map[i] = i;
355f80ce2aSJacob Faibussowitsch   CHKERRQ(ISLocalToGlobalMappingApply(ltog,n,map,map));
365f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetSize(A,&N,NULL));
375f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(N,&colors));
38c4762a1bSJed Brown   for (i=0; i<N; i++) colors[i] = -1;
39c4762a1bSJed Brown   for (i=0; i<n; i++) colors[map[i]]= icolors[i];
405f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(map));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSynchronizedPrintf(MPI_COMM_WORLD,"[%d]Global colors \n",rank));
425f80ce2aSJacob Faibussowitsch   for (i=0; i<N; i++) CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD,"%" PetscInt_FMT " %" PetscInt_FMT "\n",i,colors[i]));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD,stdout));
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   /*   second, compress the A matrix */
465f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetRandom(A,NULL));
475f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(A,NULL));
485f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetLocalSize(A,&nrow,NULL));
495f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipRange(A,&rstart,NULL));
505f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc1(nrow*nc,&cm));
51c4762a1bSJed Brown   for (i=0; i<nrow; i++) {
525f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetRow(A,rstart+i,&ncols,&cols,&vals));
53c4762a1bSJed Brown     for (j=0; j<ncols; j++) {
542c71b3e2SJacob Faibussowitsch       PetscCheckFalse(colors[cols[j]] < 0,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Global column %" PetscInt_FMT " had no color",cols[j]);
55c4762a1bSJed Brown       cm[i + nrow*colors[cols[j]]] = vals[j];
56c4762a1bSJed Brown     }
575f80ce2aSJacob Faibussowitsch     CHKERRQ(MatRestoreRow(A,rstart+i,&ncols,&cols,&vals));
58c4762a1bSJed Brown   }
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   /* print compressed matrix */
615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSynchronizedPrintf(MPI_COMM_WORLD,"[%d] Compressed matrix \n",rank));
62c4762a1bSJed Brown   for (i=0; i<nrow; i++) {
63c4762a1bSJed Brown     for (j=0; j<nc; j++) {
645f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSynchronizedPrintf(MPI_COMM_WORLD,"%12.4e  ",(double)PetscAbsScalar(cm[i+nrow*j])));
65c4762a1bSJed Brown     }
665f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSynchronizedPrintf(MPI_COMM_WORLD,"\n"));
67c4762a1bSJed Brown   }
685f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD,stdout));
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   /* put the compressed matrix into the standard matrix */
715f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDuplicate(A,MAT_COPY_VALUES,&B));
725f80ce2aSJacob Faibussowitsch   CHKERRQ(MatZeroEntries(A));
735f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(B,0));
745f80ce2aSJacob Faibussowitsch   CHKERRQ(MatFDColoringCreate(A,iscoloring,&fdcoloring));
755f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-single_block",&single));
76c4762a1bSJed Brown   if (single) {
775f80ce2aSJacob Faibussowitsch     CHKERRQ(MatFDColoringSetBlockSize(fdcoloring,PETSC_DEFAULT,nc));
78c4762a1bSJed Brown   }
795f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-two_block",&two));
80c4762a1bSJed Brown   if (two) {
815f80ce2aSJacob Faibussowitsch     CHKERRQ(MatFDColoringSetBlockSize(fdcoloring,PETSC_DEFAULT,2));
82c4762a1bSJed Brown   }
835f80ce2aSJacob Faibussowitsch   CHKERRQ(MatFDColoringSetFromOptions(fdcoloring));
845f80ce2aSJacob Faibussowitsch   CHKERRQ(MatFDColoringSetUp(A,iscoloring,fdcoloring));
85c4762a1bSJed Brown 
865f80ce2aSJacob Faibussowitsch   CHKERRQ(MatFDColoringSetValues(A,fdcoloring,cm));
875f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(A,NULL));
88c4762a1bSJed Brown 
89c4762a1bSJed Brown   /* check the values were put in the correct locations */
905f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAXPY(A,-1.0,B,SAME_NONZERO_PATTERN));
915f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(A,NULL));
925f80ce2aSJacob Faibussowitsch   CHKERRQ(MatNorm(A,NORM_FROBENIUS,&norm));
93c4762a1bSJed Brown   if (norm > PETSC_MACHINE_EPSILON) {
945f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Matrix is not identical, problem with MatFDColoringSetValues()\n"));
95c4762a1bSJed Brown   }
965f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(colors));
975f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(cm));
985f80ce2aSJacob Faibussowitsch   CHKERRQ(ISColoringDestroy(&iscoloring));
995f80ce2aSJacob Faibussowitsch   CHKERRQ(MatFDColoringDestroy(&fdcoloring));
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
1015f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
1025f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
103*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
104*b122ec5aSJacob Faibussowitsch   return 0;
105c4762a1bSJed Brown }
106c4762a1bSJed Brown 
107c4762a1bSJed Brown /*TEST
108c4762a1bSJed Brown 
109c4762a1bSJed Brown    test:
110c4762a1bSJed Brown       nsize: 2
111c4762a1bSJed Brown       requires: !complex
112c4762a1bSJed Brown 
113c4762a1bSJed Brown    test:
114c4762a1bSJed Brown       suffix: single
115c4762a1bSJed Brown       requires: !complex
116c4762a1bSJed Brown       nsize: 2
117c4762a1bSJed Brown       args: -single_block
118c4762a1bSJed Brown       output_file: output/ex240_1.out
119c4762a1bSJed Brown 
120c4762a1bSJed Brown    test:
121c4762a1bSJed Brown       suffix: two
122c4762a1bSJed Brown       requires: !complex
123c4762a1bSJed Brown       nsize: 2
124c4762a1bSJed Brown       args: -two_block
125c4762a1bSJed Brown       output_file: output/ex240_1.out
126c4762a1bSJed Brown 
127c4762a1bSJed Brown    test:
128c4762a1bSJed Brown       suffix: 2
129c4762a1bSJed Brown       requires: !complex
130c4762a1bSJed Brown       nsize: 5
131c4762a1bSJed Brown 
132c4762a1bSJed Brown TEST*/
133