xref: /petsc/src/dm/impls/composite/packm.c (revision a12302e241325cab813a3e2ea9582901af8a2164)
1*a12302e2SJed Brown #define PETSCDM_DLL
2*a12302e2SJed Brown 
3*a12302e2SJed Brown #include "packimpl.h"
4*a12302e2SJed Brown 
5*a12302e2SJed Brown #undef __FUNCT__
6*a12302e2SJed Brown #define __FUNCT__ "DMGetMatrix_Composite_Nest"
7*a12302e2SJed Brown static PetscErrorCode DMGetMatrix_Composite_Nest(DM dm,const MatType mtype,Mat *J)
8*a12302e2SJed Brown {
9*a12302e2SJed Brown 
10*a12302e2SJed Brown   PetscFunctionBegin;
11*a12302e2SJed Brown   SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"not implemented");
12*a12302e2SJed Brown   PetscFunctionReturn(0);
13*a12302e2SJed Brown }
14*a12302e2SJed Brown 
15*a12302e2SJed Brown #undef __FUNCT__
16*a12302e2SJed Brown #define __FUNCT__ "DMGetMatrix_Composite_AIJ"
17*a12302e2SJed Brown static PetscErrorCode DMGetMatrix_Composite_AIJ(DM dm,const MatType mtype,Mat *J)
18*a12302e2SJed Brown {
19*a12302e2SJed Brown   PetscErrorCode         ierr;
20*a12302e2SJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
21*a12302e2SJed Brown   struct DMCompositeLink *next = com->next;
22*a12302e2SJed Brown   PetscInt               m,*dnz,*onz,i,j,mA;
23*a12302e2SJed Brown   Mat                    Atmp;
24*a12302e2SJed Brown   PetscMPIInt            rank;
25*a12302e2SJed Brown   PetscBool              dense = PETSC_FALSE;
26*a12302e2SJed Brown 
27*a12302e2SJed Brown   PetscFunctionBegin;
28*a12302e2SJed Brown   /* use global vector to determine layout needed for matrix */
29*a12302e2SJed Brown   m = com->n;
30*a12302e2SJed Brown 
31*a12302e2SJed Brown   ierr = MatCreate(((PetscObject)dm)->comm,J);CHKERRQ(ierr);
32*a12302e2SJed Brown   ierr = MatSetSizes(*J,m,m,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
33*a12302e2SJed Brown   ierr = MatSetType(*J,mtype);CHKERRQ(ierr);
34*a12302e2SJed Brown 
35*a12302e2SJed Brown   /*
36*a12302e2SJed Brown      Extremely inefficient but will compute entire Jacobian for testing
37*a12302e2SJed Brown   */
38*a12302e2SJed Brown   ierr = PetscOptionsGetBool(PETSC_NULL,"-dmcomposite_dense_jacobian",&dense,PETSC_NULL);CHKERRQ(ierr);
39*a12302e2SJed Brown   if (dense) {
40*a12302e2SJed Brown     PetscInt    rstart,rend,*indices;
41*a12302e2SJed Brown     PetscScalar *values;
42*a12302e2SJed Brown 
43*a12302e2SJed Brown     mA    = com->N;
44*a12302e2SJed Brown     ierr = MatMPIAIJSetPreallocation(*J,mA,PETSC_NULL,mA-m,PETSC_NULL);CHKERRQ(ierr);
45*a12302e2SJed Brown     ierr = MatSeqAIJSetPreallocation(*J,mA,PETSC_NULL);CHKERRQ(ierr);
46*a12302e2SJed Brown 
47*a12302e2SJed Brown     ierr = MatGetOwnershipRange(*J,&rstart,&rend);CHKERRQ(ierr);
48*a12302e2SJed Brown     ierr = PetscMalloc2(mA,PetscScalar,&values,mA,PetscInt,&indices);CHKERRQ(ierr);
49*a12302e2SJed Brown     ierr = PetscMemzero(values,mA*sizeof(PetscScalar));CHKERRQ(ierr);
50*a12302e2SJed Brown     for (i=0; i<mA; i++) indices[i] = i;
51*a12302e2SJed Brown     for (i=rstart; i<rend; i++) {
52*a12302e2SJed Brown       ierr = MatSetValues(*J,1,&i,mA,indices,values,INSERT_VALUES);CHKERRQ(ierr);
53*a12302e2SJed Brown     }
54*a12302e2SJed Brown     ierr = PetscFree2(values,indices);CHKERRQ(ierr);
55*a12302e2SJed Brown     ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
56*a12302e2SJed Brown     ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
57*a12302e2SJed Brown     PetscFunctionReturn(0);
58*a12302e2SJed Brown   }
59*a12302e2SJed Brown 
60*a12302e2SJed Brown   ierr = MPI_Comm_rank(((PetscObject)dm)->comm,&rank);CHKERRQ(ierr);
61*a12302e2SJed Brown   ierr = MatPreallocateInitialize(((PetscObject)dm)->comm,m,m,dnz,onz);CHKERRQ(ierr);
62*a12302e2SJed Brown   /* loop over packed objects, handling one at at time */
63*a12302e2SJed Brown   next = com->next;
64*a12302e2SJed Brown   while (next) {
65*a12302e2SJed Brown     if (next->type == DMCOMPOSITE_ARRAY) {
66*a12302e2SJed Brown       if (rank == next->rank) {  /* zero the "little" block */
67*a12302e2SJed Brown         for (j=com->rstart+next->rstart; j<com->rstart+next->rstart+next->n; j++) {
68*a12302e2SJed Brown           for (i=com->rstart+next->rstart; i<com->rstart+next->rstart+next->n; i++) {
69*a12302e2SJed Brown             ierr = MatPreallocateSet(j,1,&i,dnz,onz);CHKERRQ(ierr);
70*a12302e2SJed Brown           }
71*a12302e2SJed Brown         }
72*a12302e2SJed Brown       }
73*a12302e2SJed Brown     } else if (next->type == DMCOMPOSITE_DM) {
74*a12302e2SJed Brown       PetscInt       nc,rstart,*ccols,maxnc;
75*a12302e2SJed Brown       const PetscInt *cols,*rstarts;
76*a12302e2SJed Brown       PetscMPIInt    proc;
77*a12302e2SJed Brown 
78*a12302e2SJed Brown       ierr = DMGetMatrix(next->dm,mtype,&Atmp);CHKERRQ(ierr);
79*a12302e2SJed Brown       ierr = MatGetOwnershipRange(Atmp,&rstart,PETSC_NULL);CHKERRQ(ierr);
80*a12302e2SJed Brown       ierr = MatGetOwnershipRanges(Atmp,&rstarts);CHKERRQ(ierr);
81*a12302e2SJed Brown       ierr = MatGetLocalSize(Atmp,&mA,PETSC_NULL);CHKERRQ(ierr);
82*a12302e2SJed Brown 
83*a12302e2SJed Brown       maxnc = 0;
84*a12302e2SJed Brown       for (i=0; i<mA; i++) {
85*a12302e2SJed Brown         ierr  = MatGetRow(Atmp,rstart+i,&nc,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
86*a12302e2SJed Brown         ierr  = MatRestoreRow(Atmp,rstart+i,&nc,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
87*a12302e2SJed Brown         maxnc = PetscMax(nc,maxnc);
88*a12302e2SJed Brown       }
89*a12302e2SJed Brown       ierr = PetscMalloc(maxnc*sizeof(PetscInt),&ccols);CHKERRQ(ierr);
90*a12302e2SJed Brown       for (i=0; i<mA; i++) {
91*a12302e2SJed Brown         ierr = MatGetRow(Atmp,rstart+i,&nc,&cols,PETSC_NULL);CHKERRQ(ierr);
92*a12302e2SJed Brown         /* remap the columns taking into how much they are shifted on each process */
93*a12302e2SJed Brown         for (j=0; j<nc; j++) {
94*a12302e2SJed Brown           proc = 0;
95*a12302e2SJed Brown           while (cols[j] >= rstarts[proc+1]) proc++;
96*a12302e2SJed Brown           ccols[j] = cols[j] + next->grstarts[proc] - rstarts[proc];
97*a12302e2SJed Brown         }
98*a12302e2SJed Brown         ierr = MatPreallocateSet(com->rstart+next->rstart+i,nc,ccols,dnz,onz);CHKERRQ(ierr);
99*a12302e2SJed Brown         ierr = MatRestoreRow(Atmp,rstart+i,&nc,&cols,PETSC_NULL);CHKERRQ(ierr);
100*a12302e2SJed Brown       }
101*a12302e2SJed Brown       ierr = PetscFree(ccols);CHKERRQ(ierr);
102*a12302e2SJed Brown       ierr = MatDestroy(Atmp);CHKERRQ(ierr);
103*a12302e2SJed Brown     } else {
104*a12302e2SJed Brown       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
105*a12302e2SJed Brown     }
106*a12302e2SJed Brown     next = next->next;
107*a12302e2SJed Brown   }
108*a12302e2SJed Brown   if (com->FormCoupleLocations) {
109*a12302e2SJed Brown     ierr = (*com->FormCoupleLocations)(dm,PETSC_NULL,dnz,onz,__rstart,__nrows,__start,__end);CHKERRQ(ierr);
110*a12302e2SJed Brown   }
111*a12302e2SJed Brown   ierr = MatMPIAIJSetPreallocation(*J,0,dnz,0,onz);CHKERRQ(ierr);
112*a12302e2SJed Brown   ierr = MatSeqAIJSetPreallocation(*J,0,dnz);CHKERRQ(ierr);
113*a12302e2SJed Brown   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
114*a12302e2SJed Brown 
115*a12302e2SJed Brown   next = com->next;
116*a12302e2SJed Brown   while (next) {
117*a12302e2SJed Brown     if (next->type == DMCOMPOSITE_ARRAY) {
118*a12302e2SJed Brown       if (rank == next->rank) {
119*a12302e2SJed Brown         for (j=com->rstart+next->rstart; j<com->rstart+next->rstart+next->n; j++) {
120*a12302e2SJed Brown           for (i=com->rstart+next->rstart; i<com->rstart+next->rstart+next->n; i++) {
121*a12302e2SJed Brown             ierr = MatSetValue(*J,j,i,0.0,INSERT_VALUES);CHKERRQ(ierr);
122*a12302e2SJed Brown           }
123*a12302e2SJed Brown         }
124*a12302e2SJed Brown       }
125*a12302e2SJed Brown     } else if (next->type == DMCOMPOSITE_DM) {
126*a12302e2SJed Brown       PetscInt          nc,rstart,row,maxnc,*ccols;
127*a12302e2SJed Brown       const PetscInt    *cols,*rstarts;
128*a12302e2SJed Brown       const PetscScalar *values;
129*a12302e2SJed Brown       PetscMPIInt       proc;
130*a12302e2SJed Brown 
131*a12302e2SJed Brown       ierr = DMGetMatrix(next->dm,mtype,&Atmp);CHKERRQ(ierr);
132*a12302e2SJed Brown       ierr = MatGetOwnershipRange(Atmp,&rstart,PETSC_NULL);CHKERRQ(ierr);
133*a12302e2SJed Brown       ierr = MatGetOwnershipRanges(Atmp,&rstarts);CHKERRQ(ierr);
134*a12302e2SJed Brown       ierr = MatGetLocalSize(Atmp,&mA,PETSC_NULL);CHKERRQ(ierr);
135*a12302e2SJed Brown       maxnc = 0;
136*a12302e2SJed Brown       for (i=0; i<mA; i++) {
137*a12302e2SJed Brown         ierr  = MatGetRow(Atmp,rstart+i,&nc,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
138*a12302e2SJed Brown         ierr  = MatRestoreRow(Atmp,rstart+i,&nc,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
139*a12302e2SJed Brown         maxnc = PetscMax(nc,maxnc);
140*a12302e2SJed Brown       }
141*a12302e2SJed Brown       ierr = PetscMalloc(maxnc*sizeof(PetscInt),&ccols);CHKERRQ(ierr);
142*a12302e2SJed Brown       for (i=0; i<mA; i++) {
143*a12302e2SJed Brown         ierr = MatGetRow(Atmp,rstart+i,&nc,(const PetscInt **)&cols,&values);CHKERRQ(ierr);
144*a12302e2SJed Brown         for (j=0; j<nc; j++) {
145*a12302e2SJed Brown           proc = 0;
146*a12302e2SJed Brown           while (cols[j] >= rstarts[proc+1]) proc++;
147*a12302e2SJed Brown           ccols[j] = cols[j] + next->grstarts[proc] - rstarts[proc];
148*a12302e2SJed Brown         }
149*a12302e2SJed Brown         row  = com->rstart+next->rstart+i;
150*a12302e2SJed Brown         ierr = MatSetValues(*J,1,&row,nc,ccols,values,INSERT_VALUES);CHKERRQ(ierr);
151*a12302e2SJed Brown         ierr = MatRestoreRow(Atmp,rstart+i,&nc,(const PetscInt **)&cols,&values);CHKERRQ(ierr);
152*a12302e2SJed Brown       }
153*a12302e2SJed Brown       ierr = PetscFree(ccols);CHKERRQ(ierr);
154*a12302e2SJed Brown       ierr = MatDestroy(Atmp);CHKERRQ(ierr);
155*a12302e2SJed Brown     } else {
156*a12302e2SJed Brown       SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_SUP,"Cannot handle that object type yet");
157*a12302e2SJed Brown     }
158*a12302e2SJed Brown     next = next->next;
159*a12302e2SJed Brown   }
160*a12302e2SJed Brown   if (com->FormCoupleLocations) {
161*a12302e2SJed Brown     PetscInt __rstart;
162*a12302e2SJed Brown     ierr = MatGetOwnershipRange(*J,&__rstart,PETSC_NULL);CHKERRQ(ierr);
163*a12302e2SJed Brown     ierr = (*com->FormCoupleLocations)(dm,*J,PETSC_NULL,PETSC_NULL,__rstart,0,0,0);CHKERRQ(ierr);
164*a12302e2SJed Brown   }
165*a12302e2SJed Brown   ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
166*a12302e2SJed Brown   ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
167*a12302e2SJed Brown   PetscFunctionReturn(0);
168*a12302e2SJed Brown }
169*a12302e2SJed Brown 
170*a12302e2SJed Brown #undef __FUNCT__
171*a12302e2SJed Brown #define __FUNCT__ "DMGetMatrix_Composite"
172*a12302e2SJed Brown PetscErrorCode PETSCDM_DLLEXPORT DMGetMatrix_Composite(DM dm,const MatType mtype,Mat *J)
173*a12302e2SJed Brown {
174*a12302e2SJed Brown   PetscErrorCode         ierr;
175*a12302e2SJed Brown   PetscBool              usenest;
176*a12302e2SJed Brown   ISLocalToGlobalMapping ltogmap,ltogmapb;
177*a12302e2SJed Brown 
178*a12302e2SJed Brown   PetscFunctionBegin;
179*a12302e2SJed Brown   ierr = PetscStrcmp(mtype,MATNEST,&usenest);CHKERRQ(ierr);
180*a12302e2SJed Brown   if (usenest) {
181*a12302e2SJed Brown     ierr = DMGetMatrix_Composite_Nest(dm,mtype,J);CHKERRQ(ierr);
182*a12302e2SJed Brown   } else {
183*a12302e2SJed Brown     ierr = DMGetMatrix_Composite_AIJ(dm,mtype?mtype:MATAIJ,J);CHKERRQ(ierr);
184*a12302e2SJed Brown   }
185*a12302e2SJed Brown 
186*a12302e2SJed Brown   ierr = DMGetLocalToGlobalMapping(dm,&ltogmap);CHKERRQ(ierr);
187*a12302e2SJed Brown   ierr = DMGetLocalToGlobalMappingBlock(dm,&ltogmapb);CHKERRQ(ierr);
188*a12302e2SJed Brown   ierr = MatSetLocalToGlobalMapping(*J,ltogmap,ltogmap);CHKERRQ(ierr);
189*a12302e2SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(*J,ltogmapb,ltogmapb);CHKERRQ(ierr);
190*a12302e2SJed Brown   PetscFunctionReturn(0);
191*a12302e2SJed Brown }
192