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