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