xref: /petsc/src/dm/impls/composite/packm.c (revision 9566063d113dddea24716c546802770db7481bc0)
1a12302e2SJed Brown 
2ccd284c7SBarry Smith #include <../src/dm/impls/composite/packimpl.h>       /*I  "petscdmcomposite.h"  I*/
3a12302e2SJed Brown 
4b412c318SBarry Smith static PetscErrorCode DMCreateMatrix_Composite_Nest(DM dm,Mat *J)
5a12302e2SJed Brown {
6b989ae6dSJed Brown   const DM_Composite           *com = (DM_Composite*)dm->data;
7b989ae6dSJed Brown   const struct DMCompositeLink *rlink,*clink;
8b989ae6dSJed Brown   IS                           *isg;
9b989ae6dSJed Brown   Mat                          *submats;
10b989ae6dSJed Brown   PetscInt                     i,j,n;
11a12302e2SJed Brown 
12a12302e2SJed Brown   PetscFunctionBegin;
139ae5db72SJed Brown   n = com->nDM;                 /* Total number of entries */
14b989ae6dSJed Brown 
15b989ae6dSJed Brown   /* Explicit index sets are not required for MatCreateNest, but getting them here allows MatNest to do consistency
16b989ae6dSJed Brown    * checking and allows ISEqual to compare by identity instead of by contents. */
17*9566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetGlobalISs(dm,&isg));
18b989ae6dSJed Brown 
19b989ae6dSJed Brown   /* Get submatrices */
20*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n*n,&submats));
21b989ae6dSJed Brown   for (i=0,rlink=com->next; rlink; i++,rlink=rlink->next) {
22b989ae6dSJed Brown     for (j=0,clink=com->next; clink; j++,clink=clink->next) {
230298fd71SBarry Smith       Mat sub = NULL;
24b989ae6dSJed Brown       if (i == j) {
25*9566063dSJacob Faibussowitsch         PetscCall(DMCreateMatrix(rlink->dm,&sub));
2628b400f6SJacob Faibussowitsch       } else PetscCheck(!com->FormCoupleLocations,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Cannot manage off-diagonal parts yet");
27b989ae6dSJed Brown       submats[i*n+j] = sub;
28b989ae6dSJed Brown     }
29b989ae6dSJed Brown   }
30b989ae6dSJed Brown 
31*9566063dSJacob Faibussowitsch   PetscCall(MatCreateNest(PetscObjectComm((PetscObject)dm),n,isg,n,isg,submats,J));
32b989ae6dSJed Brown 
33b989ae6dSJed Brown   /* Disown references */
34*9566063dSJacob Faibussowitsch   for (i=0; i<n; i++) PetscCall(ISDestroy(&isg[i]));
35*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(isg));
36b989ae6dSJed Brown 
37b989ae6dSJed Brown   for (i=0; i<n*n; i++) {
38*9566063dSJacob Faibussowitsch     if (submats[i]) PetscCall(MatDestroy(&submats[i]));
39b989ae6dSJed Brown   }
40*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(submats));
41a12302e2SJed Brown   PetscFunctionReturn(0);
42a12302e2SJed Brown }
43a12302e2SJed Brown 
44b412c318SBarry Smith static PetscErrorCode DMCreateMatrix_Composite_AIJ(DM dm,Mat *J)
45a12302e2SJed Brown {
46a12302e2SJed Brown   PetscErrorCode         ierr;
47a12302e2SJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
483bf036e2SBarry Smith   struct DMCompositeLink *next;
49a12302e2SJed Brown   PetscInt               m,*dnz,*onz,i,j,mA;
50a12302e2SJed Brown   Mat                    Atmp;
51a12302e2SJed Brown   PetscMPIInt            rank;
52fcfd50ebSBarry Smith   PetscBool              dense = PETSC_FALSE;
53a12302e2SJed Brown 
54a12302e2SJed Brown   PetscFunctionBegin;
55a12302e2SJed Brown   /* use global vector to determine layout needed for matrix */
56a12302e2SJed Brown   m = com->n;
57a12302e2SJed Brown 
58*9566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)dm),J));
59*9566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*J,m,m,PETSC_DETERMINE,PETSC_DETERMINE));
60*9566063dSJacob Faibussowitsch   PetscCall(MatSetType(*J,dm->mattype));
61a12302e2SJed Brown 
62a12302e2SJed Brown   /*
63a12302e2SJed Brown      Extremely inefficient but will compute entire Jacobian for testing
64a12302e2SJed Brown   */
65*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(((PetscObject)dm)->options,((PetscObject)dm)->prefix,"-dmcomposite_dense_jacobian",&dense,NULL));
66a12302e2SJed Brown   if (dense) {
67a12302e2SJed Brown     PetscInt    rstart,rend,*indices;
68a12302e2SJed Brown     PetscScalar *values;
69a12302e2SJed Brown 
70a12302e2SJed Brown     mA   = com->N;
71*9566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJSetPreallocation(*J,mA,NULL,mA-m,NULL));
72*9566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(*J,mA,NULL));
73a12302e2SJed Brown 
74*9566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(*J,&rstart,&rend));
75*9566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(mA,&values,mA,&indices));
76*9566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(values,mA));
77a12302e2SJed Brown     for (i=0; i<mA; i++) indices[i] = i;
78a12302e2SJed Brown     for (i=rstart; i<rend; i++) {
79*9566063dSJacob Faibussowitsch       PetscCall(MatSetValues(*J,1,&i,mA,indices,values,INSERT_VALUES));
80a12302e2SJed Brown     }
81*9566063dSJacob Faibussowitsch     PetscCall(PetscFree2(values,indices));
82*9566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY));
83*9566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY));
84a12302e2SJed Brown     PetscFunctionReturn(0);
85a12302e2SJed Brown   }
86a12302e2SJed Brown 
87*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank));
88*9566063dSJacob Faibussowitsch   ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)dm),m,m,dnz,onz);PetscCall(ierr);
89a12302e2SJed Brown   /* loop over packed objects, handling one at at time */
90a12302e2SJed Brown   next = com->next;
91a12302e2SJed Brown   while (next) {
92a12302e2SJed Brown     PetscInt       nc,rstart,*ccols,maxnc;
93a12302e2SJed Brown     const PetscInt *cols,*rstarts;
94a12302e2SJed Brown     PetscMPIInt    proc;
95a12302e2SJed Brown 
96*9566063dSJacob Faibussowitsch     PetscCall(DMCreateMatrix(next->dm,&Atmp));
97*9566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(Atmp,&rstart,NULL));
98*9566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRanges(Atmp,&rstarts));
99*9566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(Atmp,&mA,NULL));
100a12302e2SJed Brown 
101a12302e2SJed Brown     maxnc = 0;
102a12302e2SJed Brown     for (i=0; i<mA; i++) {
103*9566063dSJacob Faibussowitsch       PetscCall(MatGetRow(Atmp,rstart+i,&nc,NULL,NULL));
104a12302e2SJed Brown       maxnc = PetscMax(nc,maxnc);
105*9566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(Atmp,rstart+i,&nc,NULL,NULL));
106a12302e2SJed Brown     }
107*9566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(maxnc,&ccols));
108a12302e2SJed Brown     for (i=0; i<mA; i++) {
109*9566063dSJacob Faibussowitsch       PetscCall(MatGetRow(Atmp,rstart+i,&nc,&cols,NULL));
110a12302e2SJed Brown       /* remap the columns taking into how much they are shifted on each process */
111a12302e2SJed Brown       for (j=0; j<nc; j++) {
112a12302e2SJed Brown         proc = 0;
113a12302e2SJed Brown         while (cols[j] >= rstarts[proc+1]) proc++;
114a12302e2SJed Brown         ccols[j] = cols[j] + next->grstarts[proc] - rstarts[proc];
115a12302e2SJed Brown       }
116*9566063dSJacob Faibussowitsch       PetscCall(MatPreallocateSet(com->rstart+next->rstart+i,nc,ccols,dnz,onz));
117*9566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(Atmp,rstart+i,&nc,&cols,NULL));
118a12302e2SJed Brown     }
119*9566063dSJacob Faibussowitsch     PetscCall(PetscFree(ccols));
120*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Atmp));
121a12302e2SJed Brown     next = next->next;
122a12302e2SJed Brown   }
123a12302e2SJed Brown   if (com->FormCoupleLocations) {
124*9566063dSJacob Faibussowitsch     PetscCall((*com->FormCoupleLocations)(dm,NULL,dnz,onz,__rstart,__nrows,__start,__end));
125a12302e2SJed Brown   }
126*9566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(*J,0,dnz,0,onz));
127*9566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(*J,0,dnz));
128*9566063dSJacob Faibussowitsch   ierr = MatPreallocateFinalize(dnz,onz);PetscCall(ierr);
129a12302e2SJed Brown 
130fcfd50ebSBarry Smith   if (dm->prealloc_only) PetscFunctionReturn(0);
131ff6157d0SJed Brown 
132a12302e2SJed Brown   next = com->next;
133a12302e2SJed Brown   while (next) {
134a12302e2SJed Brown     PetscInt          nc,rstart,row,maxnc,*ccols;
135a12302e2SJed Brown     const PetscInt    *cols,*rstarts;
136a12302e2SJed Brown     const PetscScalar *values;
137a12302e2SJed Brown     PetscMPIInt       proc;
138a12302e2SJed Brown 
139*9566063dSJacob Faibussowitsch     PetscCall(DMCreateMatrix(next->dm,&Atmp));
140*9566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(Atmp,&rstart,NULL));
141*9566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRanges(Atmp,&rstarts));
142*9566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(Atmp,&mA,NULL));
143a12302e2SJed Brown     maxnc = 0;
144a12302e2SJed Brown     for (i=0; i<mA; i++) {
145*9566063dSJacob Faibussowitsch       PetscCall(MatGetRow(Atmp,rstart+i,&nc,NULL,NULL));
146a12302e2SJed Brown       maxnc = PetscMax(nc,maxnc);
147*9566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(Atmp,rstart+i,&nc,NULL,NULL));
148a12302e2SJed Brown     }
149*9566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(maxnc,&ccols));
150a12302e2SJed Brown     for (i=0; i<mA; i++) {
151*9566063dSJacob Faibussowitsch       PetscCall(MatGetRow(Atmp,rstart+i,&nc,(const PetscInt**)&cols,&values));
152a12302e2SJed Brown       for (j=0; j<nc; j++) {
153a12302e2SJed Brown         proc = 0;
154a12302e2SJed Brown         while (cols[j] >= rstarts[proc+1]) proc++;
155a12302e2SJed Brown         ccols[j] = cols[j] + next->grstarts[proc] - rstarts[proc];
156a12302e2SJed Brown       }
157a12302e2SJed Brown       row  = com->rstart+next->rstart+i;
158*9566063dSJacob Faibussowitsch       PetscCall(MatSetValues(*J,1,&row,nc,ccols,values,INSERT_VALUES));
159*9566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(Atmp,rstart+i,&nc,(const PetscInt**)&cols,&values));
160a12302e2SJed Brown     }
161*9566063dSJacob Faibussowitsch     PetscCall(PetscFree(ccols));
162*9566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Atmp));
163a12302e2SJed Brown     next = next->next;
164a12302e2SJed Brown   }
165a12302e2SJed Brown   if (com->FormCoupleLocations) {
166a12302e2SJed Brown     PetscInt __rstart;
167*9566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(*J,&__rstart,NULL));
168*9566063dSJacob Faibussowitsch     PetscCall((*com->FormCoupleLocations)(dm,*J,NULL,NULL,__rstart,0,0,0));
169a12302e2SJed Brown   }
170*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY));
171*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY));
172a12302e2SJed Brown   PetscFunctionReturn(0);
173a12302e2SJed Brown }
174a12302e2SJed Brown 
175b412c318SBarry Smith PetscErrorCode DMCreateMatrix_Composite(DM dm,Mat *J)
176a12302e2SJed Brown {
177a12302e2SJed Brown   PetscBool              usenest;
17845b6f7e9SBarry Smith   ISLocalToGlobalMapping ltogmap;
179a12302e2SJed Brown 
180a12302e2SJed Brown   PetscFunctionBegin;
181*9566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
182*9566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dm));
183*9566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(dm->mattype,MATNEST,&usenest));
184a12302e2SJed Brown   if (usenest) {
185*9566063dSJacob Faibussowitsch     PetscCall(DMCreateMatrix_Composite_Nest(dm,J));
186a12302e2SJed Brown   } else {
187*9566063dSJacob Faibussowitsch     PetscCall(DMCreateMatrix_Composite_AIJ(dm,J));
188a12302e2SJed Brown   }
189a12302e2SJed Brown 
190*9566063dSJacob Faibussowitsch   PetscCall(DMGetLocalToGlobalMapping(dm,&ltogmap));
191*9566063dSJacob Faibussowitsch   PetscCall(MatSetLocalToGlobalMapping(*J,ltogmap,ltogmap));
192*9566063dSJacob Faibussowitsch   PetscCall(MatSetDM(*J,dm));
193a12302e2SJed Brown   PetscFunctionReturn(0);
194a12302e2SJed Brown }
195