xref: /petsc/src/dm/impls/composite/packm.c (revision dcca6d9d80ebd869fe6029bd05a3aa9faafef49e)
1a12302e2SJed Brown 
2ccd284c7SBarry Smith #include <../src/dm/impls/composite/packimpl.h>       /*I  "petscdmcomposite.h"  I*/
3a12302e2SJed Brown 
4a12302e2SJed Brown #undef __FUNCT__
5950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_Composite_Nest"
6b412c318SBarry Smith static PetscErrorCode DMCreateMatrix_Composite_Nest(DM dm,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) {
260298fd71SBarry Smith       Mat sub = NULL;
27b989ae6dSJed Brown       if (i == j) {
28b412c318SBarry Smith         ierr = DMCreateMatrix(rlink->dm,&sub);CHKERRQ(ierr);
29ce94432eSBarry Smith       } else if (com->FormCoupleLocations) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Cannot manage off-diagonal parts yet");
30b989ae6dSJed Brown       submats[i*n+j] = sub;
31b989ae6dSJed Brown     }
32b989ae6dSJed Brown   }
33b989ae6dSJed Brown 
34ce94432eSBarry Smith   ierr = MatCreateNest(PetscObjectComm((PetscObject)dm),n,isg,n,isg,submats,J);CHKERRQ(ierr);
35b989ae6dSJed Brown 
36b989ae6dSJed Brown   /* Disown references */
37fcfd50ebSBarry Smith   for (i=0; i<n; i++) {ierr = ISDestroy(&isg[i]);CHKERRQ(ierr);}
38b989ae6dSJed Brown   ierr = PetscFree(isg);CHKERRQ(ierr);
39b989ae6dSJed Brown 
40b989ae6dSJed Brown   for (i=0; i<n*n; i++) {
41fcfd50ebSBarry Smith     if (submats[i]) {ierr = MatDestroy(&submats[i]);CHKERRQ(ierr);}
42b989ae6dSJed Brown   }
43b989ae6dSJed Brown   ierr = PetscFree(submats);CHKERRQ(ierr);
44a12302e2SJed Brown   PetscFunctionReturn(0);
45a12302e2SJed Brown }
46a12302e2SJed Brown 
47a12302e2SJed Brown #undef __FUNCT__
48950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_Composite_AIJ"
49b412c318SBarry Smith static PetscErrorCode DMCreateMatrix_Composite_AIJ(DM dm,Mat *J)
50a12302e2SJed Brown {
51a12302e2SJed Brown   PetscErrorCode         ierr;
52a12302e2SJed Brown   DM_Composite           *com = (DM_Composite*)dm->data;
533bf036e2SBarry Smith   struct DMCompositeLink *next;
54a12302e2SJed Brown   PetscInt               m,*dnz,*onz,i,j,mA;
55a12302e2SJed Brown   Mat                    Atmp;
56a12302e2SJed Brown   PetscMPIInt            rank;
57fcfd50ebSBarry Smith   PetscBool              dense = PETSC_FALSE;
58a12302e2SJed Brown 
59a12302e2SJed Brown   PetscFunctionBegin;
60a12302e2SJed Brown   /* use global vector to determine layout needed for matrix */
61a12302e2SJed Brown   m = com->n;
62a12302e2SJed Brown 
63ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr);
64a12302e2SJed Brown   ierr = MatSetSizes(*J,m,m,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
65b412c318SBarry Smith   ierr = MatSetType(*J,dm->mattype);CHKERRQ(ierr);
66a12302e2SJed Brown 
67a12302e2SJed Brown   /*
68a12302e2SJed Brown      Extremely inefficient but will compute entire Jacobian for testing
69a12302e2SJed Brown   */
700298fd71SBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)dm)->prefix,"-dmcomposite_dense_jacobian",&dense,NULL);CHKERRQ(ierr);
71a12302e2SJed Brown   if (dense) {
72a12302e2SJed Brown     PetscInt    rstart,rend,*indices;
73a12302e2SJed Brown     PetscScalar *values;
74a12302e2SJed Brown 
75a12302e2SJed Brown     mA   = com->N;
760298fd71SBarry Smith     ierr = MatMPIAIJSetPreallocation(*J,mA,NULL,mA-m,NULL);CHKERRQ(ierr);
770298fd71SBarry Smith     ierr = MatSeqAIJSetPreallocation(*J,mA,NULL);CHKERRQ(ierr);
78a12302e2SJed Brown 
79a12302e2SJed Brown     ierr = MatGetOwnershipRange(*J,&rstart,&rend);CHKERRQ(ierr);
80*dcca6d9dSJed Brown     ierr = PetscMalloc2(mA,&values,mA,&indices);CHKERRQ(ierr);
81a12302e2SJed Brown     ierr = PetscMemzero(values,mA*sizeof(PetscScalar));CHKERRQ(ierr);
82a12302e2SJed Brown     for (i=0; i<mA; i++) indices[i] = i;
83a12302e2SJed Brown     for (i=rstart; i<rend; i++) {
84a12302e2SJed Brown       ierr = MatSetValues(*J,1,&i,mA,indices,values,INSERT_VALUES);CHKERRQ(ierr);
85a12302e2SJed Brown     }
86a12302e2SJed Brown     ierr = PetscFree2(values,indices);CHKERRQ(ierr);
87a12302e2SJed Brown     ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
88a12302e2SJed Brown     ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
89a12302e2SJed Brown     PetscFunctionReturn(0);
90a12302e2SJed Brown   }
91a12302e2SJed Brown 
92ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
93ce94432eSBarry Smith   ierr = MatPreallocateInitialize(PetscObjectComm((PetscObject)dm),m,m,dnz,onz);CHKERRQ(ierr);
94a12302e2SJed Brown   /* loop over packed objects, handling one at at time */
95a12302e2SJed Brown   next = com->next;
96a12302e2SJed Brown   while (next) {
97a12302e2SJed Brown     PetscInt       nc,rstart,*ccols,maxnc;
98a12302e2SJed Brown     const PetscInt *cols,*rstarts;
99a12302e2SJed Brown     PetscMPIInt    proc;
100a12302e2SJed Brown 
101b412c318SBarry Smith     ierr = DMCreateMatrix(next->dm,&Atmp);CHKERRQ(ierr);
1020298fd71SBarry Smith     ierr = MatGetOwnershipRange(Atmp,&rstart,NULL);CHKERRQ(ierr);
103a12302e2SJed Brown     ierr = MatGetOwnershipRanges(Atmp,&rstarts);CHKERRQ(ierr);
1040298fd71SBarry Smith     ierr = MatGetLocalSize(Atmp,&mA,NULL);CHKERRQ(ierr);
105a12302e2SJed Brown 
106a12302e2SJed Brown     maxnc = 0;
107a12302e2SJed Brown     for (i=0; i<mA; i++) {
1080298fd71SBarry Smith       ierr  = MatGetRow(Atmp,rstart+i,&nc,NULL,NULL);CHKERRQ(ierr);
1094b4408bbSJed Brown       ierr  = MatRestoreRow(Atmp,rstart+i,NULL,NULL,NULL);CHKERRQ(ierr);
110a12302e2SJed Brown       maxnc = PetscMax(nc,maxnc);
111a12302e2SJed Brown     }
112a12302e2SJed Brown     ierr = PetscMalloc(maxnc*sizeof(PetscInt),&ccols);CHKERRQ(ierr);
113a12302e2SJed Brown     for (i=0; i<mA; i++) {
1140298fd71SBarry Smith       ierr = MatGetRow(Atmp,rstart+i,&nc,&cols,NULL);CHKERRQ(ierr);
115a12302e2SJed Brown       /* remap the columns taking into how much they are shifted on each process */
116a12302e2SJed Brown       for (j=0; j<nc; j++) {
117a12302e2SJed Brown         proc = 0;
118a12302e2SJed Brown         while (cols[j] >= rstarts[proc+1]) proc++;
119a12302e2SJed Brown         ccols[j] = cols[j] + next->grstarts[proc] - rstarts[proc];
120a12302e2SJed Brown       }
121a12302e2SJed Brown       ierr = MatPreallocateSet(com->rstart+next->rstart+i,nc,ccols,dnz,onz);CHKERRQ(ierr);
1220298fd71SBarry Smith       ierr = MatRestoreRow(Atmp,rstart+i,&nc,&cols,NULL);CHKERRQ(ierr);
123a12302e2SJed Brown     }
124a12302e2SJed Brown     ierr = PetscFree(ccols);CHKERRQ(ierr);
125fcfd50ebSBarry Smith     ierr = MatDestroy(&Atmp);CHKERRQ(ierr);
126a12302e2SJed Brown     next = next->next;
127a12302e2SJed Brown   }
128a12302e2SJed Brown   if (com->FormCoupleLocations) {
1290298fd71SBarry Smith     ierr = (*com->FormCoupleLocations)(dm,NULL,dnz,onz,__rstart,__nrows,__start,__end);CHKERRQ(ierr);
130a12302e2SJed Brown   }
131a12302e2SJed Brown   ierr = MatMPIAIJSetPreallocation(*J,0,dnz,0,onz);CHKERRQ(ierr);
132a12302e2SJed Brown   ierr = MatSeqAIJSetPreallocation(*J,0,dnz);CHKERRQ(ierr);
133a12302e2SJed Brown   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
134a12302e2SJed Brown 
135fcfd50ebSBarry Smith   if (dm->prealloc_only) PetscFunctionReturn(0);
136ff6157d0SJed Brown 
137a12302e2SJed Brown   next = com->next;
138a12302e2SJed Brown   while (next) {
139a12302e2SJed Brown     PetscInt          nc,rstart,row,maxnc,*ccols;
140a12302e2SJed Brown     const PetscInt    *cols,*rstarts;
141a12302e2SJed Brown     const PetscScalar *values;
142a12302e2SJed Brown     PetscMPIInt       proc;
143a12302e2SJed Brown 
144b412c318SBarry Smith     ierr  = DMCreateMatrix(next->dm,&Atmp);CHKERRQ(ierr);
1450298fd71SBarry Smith     ierr  = MatGetOwnershipRange(Atmp,&rstart,NULL);CHKERRQ(ierr);
146a12302e2SJed Brown     ierr  = MatGetOwnershipRanges(Atmp,&rstarts);CHKERRQ(ierr);
1470298fd71SBarry Smith     ierr  = MatGetLocalSize(Atmp,&mA,NULL);CHKERRQ(ierr);
148a12302e2SJed Brown     maxnc = 0;
149a12302e2SJed Brown     for (i=0; i<mA; i++) {
1500298fd71SBarry Smith       ierr  = MatGetRow(Atmp,rstart+i,&nc,NULL,NULL);CHKERRQ(ierr);
151a12302e2SJed Brown       maxnc = PetscMax(nc,maxnc);
15218ff422cSJed Brown       ierr  = MatRestoreRow(Atmp,rstart+i,&nc,NULL,NULL);CHKERRQ(ierr);
153a12302e2SJed Brown     }
154a12302e2SJed Brown     ierr = PetscMalloc(maxnc*sizeof(PetscInt),&ccols);CHKERRQ(ierr);
155a12302e2SJed Brown     for (i=0; i<mA; i++) {
156a12302e2SJed Brown       ierr = MatGetRow(Atmp,rstart+i,&nc,(const PetscInt**)&cols,&values);CHKERRQ(ierr);
157a12302e2SJed Brown       for (j=0; j<nc; j++) {
158a12302e2SJed Brown         proc = 0;
159a12302e2SJed Brown         while (cols[j] >= rstarts[proc+1]) proc++;
160a12302e2SJed Brown         ccols[j] = cols[j] + next->grstarts[proc] - rstarts[proc];
161a12302e2SJed Brown       }
162a12302e2SJed Brown       row  = com->rstart+next->rstart+i;
163a12302e2SJed Brown       ierr = MatSetValues(*J,1,&row,nc,ccols,values,INSERT_VALUES);CHKERRQ(ierr);
164a12302e2SJed Brown       ierr = MatRestoreRow(Atmp,rstart+i,&nc,(const PetscInt**)&cols,&values);CHKERRQ(ierr);
165a12302e2SJed Brown     }
166a12302e2SJed Brown     ierr = PetscFree(ccols);CHKERRQ(ierr);
167fcfd50ebSBarry Smith     ierr = MatDestroy(&Atmp);CHKERRQ(ierr);
168a12302e2SJed Brown     next = next->next;
169a12302e2SJed Brown   }
170a12302e2SJed Brown   if (com->FormCoupleLocations) {
171a12302e2SJed Brown     PetscInt __rstart;
1720298fd71SBarry Smith     ierr = MatGetOwnershipRange(*J,&__rstart,NULL);CHKERRQ(ierr);
1730298fd71SBarry Smith     ierr = (*com->FormCoupleLocations)(dm,*J,NULL,NULL,__rstart,0,0,0);CHKERRQ(ierr);
174a12302e2SJed Brown   }
175a12302e2SJed Brown   ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
176a12302e2SJed Brown   ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
177a12302e2SJed Brown   PetscFunctionReturn(0);
178a12302e2SJed Brown }
179a12302e2SJed Brown 
180a12302e2SJed Brown #undef __FUNCT__
181950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_Composite"
182b412c318SBarry Smith PetscErrorCode DMCreateMatrix_Composite(DM dm,Mat *J)
183a12302e2SJed Brown {
184a12302e2SJed Brown   PetscErrorCode         ierr;
185a12302e2SJed Brown   PetscBool              usenest;
186a12302e2SJed Brown   ISLocalToGlobalMapping ltogmap,ltogmapb;
187a12302e2SJed Brown 
188a12302e2SJed Brown   PetscFunctionBegin;
189f692024eSJed Brown   ierr = DMSetUp(dm);CHKERRQ(ierr);
190b412c318SBarry Smith   ierr = PetscStrcmp(dm->mattype,MATNEST,&usenest);CHKERRQ(ierr);
191a12302e2SJed Brown   if (usenest) {
192b412c318SBarry Smith     ierr = DMCreateMatrix_Composite_Nest(dm,J);CHKERRQ(ierr);
193a12302e2SJed Brown   } else {
194b412c318SBarry Smith     ierr = DMCreateMatrix_Composite_AIJ(dm,J);CHKERRQ(ierr);
195a12302e2SJed Brown   }
196a12302e2SJed Brown 
197a12302e2SJed Brown   ierr = DMGetLocalToGlobalMapping(dm,&ltogmap);CHKERRQ(ierr);
198a12302e2SJed Brown   ierr = DMGetLocalToGlobalMappingBlock(dm,&ltogmapb);CHKERRQ(ierr);
199a12302e2SJed Brown   ierr = MatSetLocalToGlobalMapping(*J,ltogmap,ltogmap);CHKERRQ(ierr);
200a12302e2SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(*J,ltogmapb,ltogmapb);CHKERRQ(ierr);
201a12302e2SJed Brown   PetscFunctionReturn(0);
202a12302e2SJed Brown }
203