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