xref: /petsc/src/dm/impls/moab/dmmbmat.cxx (revision 82dfd14a11d3f504e0a9b937a8040e969f448de4)
1b8ecf6d3SVijay Mahadevan #include <petsc-private/dmmbimpl.h> /*I  "petscdmmoab.h"   I*/
2b8ecf6d3SVijay Mahadevan #include <petsc-private/vecimpl.h>
3032b8ab6SVijay Mahadevan 
4032b8ab6SVijay Mahadevan #include <petscdmmoab.h>
5da6192ceSVijay Mahadevan #include <MBTagConventions.hpp>
6032b8ab6SVijay Mahadevan 
7032b8ab6SVijay Mahadevan static PetscErrorCode DMMoab_Compute_NNZ_From_Connectivity(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscBool);
8da8c5b52SVijay Mahadevan static PetscErrorCode DMMoab_MatFillMatrixEntries_Private(DM,Mat);
9032b8ab6SVijay Mahadevan 
10032b8ab6SVijay Mahadevan #undef __FUNCT__
11032b8ab6SVijay Mahadevan #define __FUNCT__ "DMCreateMatrix_Moab"
12baf0d1e0SVijay Mahadevan PetscErrorCode DMCreateMatrix_Moab(DM dm,Mat *J)
13032b8ab6SVijay Mahadevan {
14032b8ab6SVijay Mahadevan   PetscErrorCode  ierr;
15da6192ceSVijay Mahadevan   PetscInt        innz,ionz,nlsiz;
16032b8ab6SVijay Mahadevan   DM_Moab         *dmmoab=(DM_Moab*)dm->data;
17da6192ceSVijay Mahadevan   PetscInt        *nnz=0,*onz=0;
18032b8ab6SVijay Mahadevan   char            *tmp=0;
19baf0d1e0SVijay Mahadevan   MatType         mtype;
20032b8ab6SVijay Mahadevan 
21032b8ab6SVijay Mahadevan   PetscFunctionBegin;
22032b8ab6SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
23032b8ab6SVijay Mahadevan   PetscValidPointer(J,3);
24032b8ab6SVijay Mahadevan 
25db66d124SVijay Mahadevan   /* next, need to allocate the non-zero arrays to enable pre-allocation */
26baf0d1e0SVijay Mahadevan   mtype = dm->mattype;
27baf0d1e0SVijay Mahadevan   ierr = PetscStrstr(mtype, "baij", &tmp);CHKERRQ(ierr);
28*82dfd14aSVijay Mahadevan   nlsiz = (tmp ? dmmoab->nloc:dmmoab->nloc*dmmoab->numFields);
29032b8ab6SVijay Mahadevan 
30032b8ab6SVijay Mahadevan   /* allocate the nnz, onz arrays based on block size and local nodes */
31da6192ceSVijay Mahadevan   ierr = PetscMalloc((nlsiz)*sizeof(PetscInt),&nnz);CHKERRQ(ierr);
32da6192ceSVijay Mahadevan   ierr = PetscMemzero(nnz,sizeof(PetscInt)*(nlsiz));CHKERRQ(ierr);
33da6192ceSVijay Mahadevan   ierr = PetscMalloc(nlsiz*sizeof(PetscInt),&onz);CHKERRQ(ierr);
34da6192ceSVijay Mahadevan   ierr = PetscMemzero(onz,sizeof(PetscInt)*nlsiz);CHKERRQ(ierr);
35032b8ab6SVijay Mahadevan 
36032b8ab6SVijay Mahadevan   /* compute the nonzero pattern based on MOAB connectivity data for local elements */
37032b8ab6SVijay Mahadevan   ierr = DMMoab_Compute_NNZ_From_Connectivity(dm,&innz,nnz,&ionz,onz,(tmp?PETSC_TRUE:PETSC_FALSE));CHKERRQ(ierr);
38032b8ab6SVijay Mahadevan 
39da6192ceSVijay Mahadevan   /* create the Matrix and set its type as specified by user */
40da6192ceSVijay Mahadevan   ierr = MatCreate(dmmoab->pcomm->comm(), J);CHKERRQ(ierr);
41addae81cSVijay Mahadevan   ierr = MatSetSizes(*J, dmmoab->nloc*dmmoab->numFields, dmmoab->nloc*dmmoab->numFields, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
42032b8ab6SVijay Mahadevan   ierr = MatSetBlockSize(*J, dmmoab->bs);CHKERRQ(ierr);
43da6192ceSVijay Mahadevan   ierr = MatSetType(*J, mtype);CHKERRQ(ierr);
44da6192ceSVijay Mahadevan   ierr = MatSetFromOptions(*J);CHKERRQ(ierr);
45da6192ceSVijay Mahadevan 
46da6192ceSVijay Mahadevan   if (!dmmoab->ltog_map) SETERRQ(dmmoab->pcomm->comm(), PETSC_ERR_ORDER, "Cannot create a DMMoab Mat without calling DMSetUp first.");
47da6192ceSVijay Mahadevan   ierr = MatSetLocalToGlobalMapping(*J,dmmoab->ltog_map,dmmoab->ltog_map);CHKERRQ(ierr);
48032b8ab6SVijay Mahadevan 
49032b8ab6SVijay Mahadevan   /* set preallocation based on different supported Mat types */
50e23c60ebSVijay Mahadevan   ierr = MatSeqAIJSetPreallocation(*J, innz, nnz);CHKERRQ(ierr);
51e23c60ebSVijay Mahadevan   ierr = MatMPIAIJSetPreallocation(*J, innz, nnz, ionz, onz);CHKERRQ(ierr);
52032b8ab6SVijay Mahadevan   ierr = MatSeqBAIJSetPreallocation(*J, dmmoab->bs, innz, nnz);CHKERRQ(ierr);
53032b8ab6SVijay Mahadevan   ierr = MatMPIBAIJSetPreallocation(*J, dmmoab->bs, innz, nnz, ionz, onz);CHKERRQ(ierr);
54032b8ab6SVijay Mahadevan 
55da8c5b52SVijay Mahadevan   /* clean up temporary memory */
56032b8ab6SVijay Mahadevan   ierr = PetscFree(nnz);CHKERRQ(ierr);
57032b8ab6SVijay Mahadevan   ierr = PetscFree(onz);CHKERRQ(ierr);
58da8c5b52SVijay Mahadevan 
59da8c5b52SVijay Mahadevan   /* set up internal matrix data-structures */
60da8c5b52SVijay Mahadevan   ierr = MatSetUp(*J);CHKERRQ(ierr);
61da8c5b52SVijay Mahadevan 
62da8c5b52SVijay Mahadevan   /* set DM reference */
63da8c5b52SVijay Mahadevan   ierr = MatSetDM(*J, dm);CHKERRQ(ierr);
64da8c5b52SVijay Mahadevan 
65da8c5b52SVijay Mahadevan   /* set the correct NNZ pattern by setting matrix entries - make the matrix ready to use */
66da8c5b52SVijay Mahadevan   ierr = DMMoab_MatFillMatrixEntries_Private(dm,*J);CHKERRQ(ierr);
67032b8ab6SVijay Mahadevan   PetscFunctionReturn(0);
68032b8ab6SVijay Mahadevan }
69032b8ab6SVijay Mahadevan 
70032b8ab6SVijay Mahadevan 
71032b8ab6SVijay Mahadevan #undef __FUNCT__
72032b8ab6SVijay Mahadevan #define __FUNCT__ "DMMoab_Compute_NNZ_From_Connectivity"
73032b8ab6SVijay Mahadevan PetscErrorCode DMMoab_Compute_NNZ_From_Connectivity(DM dm,PetscInt* innz,PetscInt* nnz,PetscInt* ionz,PetscInt* onz,PetscBool isbaij)
74032b8ab6SVijay Mahadevan {
755905e1eaSVijay Mahadevan   PetscInt        i,f,nloc,vpere,bs,ivtx,n_nnz,n_onz;
765905e1eaSVijay Mahadevan   PetscInt        ibs,jbs,inbsize,iobsize,nfields;
77032b8ab6SVijay Mahadevan   DM_Moab         *dmmoab = (DM_Moab*)dm->data;
78032b8ab6SVijay Mahadevan   const moab::EntityHandle *connect;
79da6192ceSVijay Mahadevan   moab::Range     adjs,found,allvlocal,allvghost;
80da6192ceSVijay Mahadevan   moab::Range::iterator iter,jter;
81da6192ceSVijay Mahadevan   std::vector<moab::EntityHandle> storage;
825905e1eaSVijay Mahadevan   PetscBool isinterlaced;
83da6192ceSVijay Mahadevan   moab::EntityHandle vtx;
84032b8ab6SVijay Mahadevan   moab::ErrorCode merr;
85032b8ab6SVijay Mahadevan 
86032b8ab6SVijay Mahadevan   PetscFunctionBegin;
87032b8ab6SVijay Mahadevan   bs = dmmoab->bs;
88032b8ab6SVijay Mahadevan   nloc = dmmoab->nloc;
895905e1eaSVijay Mahadevan   nfields = dmmoab->numFields;
90*82dfd14aSVijay Mahadevan   isinterlaced=(isbaij || bs==nfields ? PETSC_TRUE : PETSC_FALSE);
91032b8ab6SVijay Mahadevan 
92da6192ceSVijay Mahadevan   /* find the truly user-expected layer of ghosted entities to decipher NNZ pattern */
93da6192ceSVijay Mahadevan   merr = dmmoab->mbiface->get_entities_by_type(dmmoab->fileset,moab::MBVERTEX,allvlocal,true);MBERRNM(merr);
94da6192ceSVijay Mahadevan   merr = dmmoab->pcomm->filter_pstatus(allvlocal,PSTATUS_NOT_OWNED,PSTATUS_NOT,-1,&adjs);MBERRNM(merr);
95da6192ceSVijay Mahadevan   allvghost = moab::subtract(allvlocal, adjs);
96032b8ab6SVijay Mahadevan 
97f6829af0SVijay Mahadevan   /* loop over the locally owned vertices and figure out the NNZ pattern using connectivity information */
98e427d9c9SVijay Mahadevan   for(iter = dmmoab->vowned->begin(),ivtx=0; iter != dmmoab->vowned->end(); iter++,ivtx++) {
99032b8ab6SVijay Mahadevan 
100da6192ceSVijay Mahadevan     vtx = *iter;
101da6192ceSVijay Mahadevan     adjs.clear();
102da6192ceSVijay Mahadevan     /* Get adjacency information for current vertex - i.e., all elements of dimension (dim) that connects
103da6192ceSVijay Mahadevan        to the current vertex. We can then decipher if a vertex is ghosted or not and compute the
104da6192ceSVijay Mahadevan        non-zero pattern accordingly. */
105da6192ceSVijay Mahadevan     merr = dmmoab->mbiface->get_adjacencies(&vtx,1,dmmoab->dim,false,adjs,moab::Interface::INTERSECT);
1060df6e276SVijay Mahadevan 
107e427d9c9SVijay Mahadevan     /* reset counters */
108e427d9c9SVijay Mahadevan     n_nnz=n_onz=0;
109e427d9c9SVijay Mahadevan     found.clear();
110e427d9c9SVijay Mahadevan 
11115de014fSVijay Mahadevan     /* loop over vertices and update the number of connectivity */
11215de014fSVijay Mahadevan     for(jter = adjs.begin(); jter != adjs.end(); jter++) {
1130df6e276SVijay Mahadevan 
11415de014fSVijay Mahadevan       /* Get connectivity information in canonical ordering for the local element */
11515de014fSVijay Mahadevan       merr = dmmoab->mbiface->get_connectivity(*jter,connect,vpere,false,&storage);MBERRNM(merr);
1160df6e276SVijay Mahadevan 
117f6829af0SVijay Mahadevan       /* loop over each element connected to the adjacent vertex and update as needed */
11815de014fSVijay Mahadevan       for (i=0; i<vpere; ++i) {
119f6829af0SVijay Mahadevan         if (connect[i] == vtx || found.find(connect[i]) != found.end()) continue; /* make sure we don't double count shared vertices */
120f6829af0SVijay Mahadevan         if (allvghost.find(connect[i]) != allvghost.end()) n_onz++; /* update out-of-proc onz */
121f6829af0SVijay Mahadevan         else n_nnz++; /* else local vertex */
12215de014fSVijay Mahadevan         found.insert(connect[i]);
12315de014fSVijay Mahadevan       }
12415de014fSVijay Mahadevan     }
12515de014fSVijay Mahadevan 
12615de014fSVijay Mahadevan     if (isbaij) {
12715de014fSVijay Mahadevan       nnz[ivtx]=n_nnz;    /* leave out self to avoid repeats -> node shared by multiple elements */
1285905e1eaSVijay Mahadevan       if (onz) {
1295905e1eaSVijay Mahadevan         onz[ivtx]=n_onz;  /* add ghost non-owned nodes */
1305905e1eaSVijay Mahadevan       }
13115de014fSVijay Mahadevan     }
132f6829af0SVijay Mahadevan     else { /* AIJ matrices */
1335905e1eaSVijay Mahadevan       if (!isinterlaced) {
1345905e1eaSVijay Mahadevan         for (f=0;f<nfields;f++) {
1355905e1eaSVijay Mahadevan           nnz[f*nloc+ivtx]=n_nnz;    /* leave out self to avoid repeats -> node shared by multiple elements */
1365905e1eaSVijay Mahadevan           if (onz)
1375905e1eaSVijay Mahadevan             onz[f*nloc+ivtx]=n_onz;  /* add ghost non-owned nodes */
1385905e1eaSVijay Mahadevan         }
1395905e1eaSVijay Mahadevan       }
1405905e1eaSVijay Mahadevan       else {
1415905e1eaSVijay Mahadevan         for (f=0;f<nfields;f++) {
1425905e1eaSVijay Mahadevan           nnz[nfields*ivtx+f]=n_nnz;    /* leave out self to avoid repeats -> node shared by multiple elements */
1435905e1eaSVijay Mahadevan           if (onz)
1445905e1eaSVijay Mahadevan             onz[nfields*ivtx+f]=n_onz;  /* add ghost non-owned nodes */
1455905e1eaSVijay Mahadevan         }
1460df6e276SVijay Mahadevan       }
1470df6e276SVijay Mahadevan     }
1480df6e276SVijay Mahadevan   }
1490df6e276SVijay Mahadevan 
150*82dfd14aSVijay Mahadevan   for (i=0;i<nloc*(isbaij?1:nfields);i++)
1515905e1eaSVijay Mahadevan     nnz[i]+=1;  /* self count the node */
1525905e1eaSVijay Mahadevan 
153*82dfd14aSVijay Mahadevan   for (ivtx=0;ivtx<nloc;ivtx++) {
1545905e1eaSVijay Mahadevan     if (!isbaij) {
1555905e1eaSVijay Mahadevan       for (ibs=0; ibs<nfields; ibs++) {
156*82dfd14aSVijay Mahadevan         if (dmmoab->dfill) {  /* first address the diagonal block */
1575905e1eaSVijay Mahadevan           /* just add up the ints -- easier/faster rather than branching based on "1" */
1585905e1eaSVijay Mahadevan           for (jbs=0,inbsize=0; jbs<nfields; jbs++)
1595905e1eaSVijay Mahadevan             inbsize += dmmoab->dfill[ibs*nfields+jbs];
1605905e1eaSVijay Mahadevan         }
161*82dfd14aSVijay Mahadevan         else inbsize=nfields; /* dense coupling since user didn't specify the component fill explicitly */
162*82dfd14aSVijay Mahadevan         if (isinterlaced) nnz[ivtx*nfields+ibs]*=inbsize;
163*82dfd14aSVijay Mahadevan         else nnz[ibs*nloc+ivtx]*=inbsize;
1645905e1eaSVijay Mahadevan 
1655905e1eaSVijay Mahadevan         if (onz) {
166*82dfd14aSVijay Mahadevan           if (dmmoab->ofill) {  /* next address the off-diagonal block */
1675905e1eaSVijay Mahadevan             /* just add up the ints -- easier/faster rather than branching based on "1" */
1685905e1eaSVijay Mahadevan             for (jbs=0,iobsize=0; jbs<nfields; jbs++)
1695905e1eaSVijay Mahadevan               iobsize += dmmoab->dfill[ibs*nfields+jbs];
1705905e1eaSVijay Mahadevan           }
171*82dfd14aSVijay Mahadevan           else iobsize=nfields; /* dense coupling since user didn't specify the component fill explicitly */
172*82dfd14aSVijay Mahadevan           if (isinterlaced) onz[ivtx*nfields+ibs]*=iobsize;
173*82dfd14aSVijay Mahadevan           else onz[ibs*nloc+ivtx]*=iobsize;
1745905e1eaSVijay Mahadevan         }
1755905e1eaSVijay Mahadevan       }
1765905e1eaSVijay Mahadevan     }
1775905e1eaSVijay Mahadevan     else {
178*82dfd14aSVijay Mahadevan       /* check if we got overzealous in our nnz and onz computations */
179*82dfd14aSVijay Mahadevan       nnz[ivtx]=(nnz[ivtx]>dmmoab->nloc?dmmoab->nloc:nnz[ivtx]);
180*82dfd14aSVijay Mahadevan       if (onz) onz[ivtx]=(onz[ivtx]>dmmoab->nloc?dmmoab->nloc:onz[ivtx]);
1815905e1eaSVijay Mahadevan     }
1825905e1eaSVijay Mahadevan   }
1835905e1eaSVijay Mahadevan   /* update innz and ionz based on local maxima */
1845905e1eaSVijay Mahadevan   if (innz || ionz) {
1850df6e276SVijay Mahadevan     if (innz) *innz=0;
1860df6e276SVijay Mahadevan     if (ionz) *ionz=0;
1875905e1eaSVijay Mahadevan     for (i=0;i<nloc*nfields;i++) {
188da6192ceSVijay Mahadevan       if (innz && (nnz[i]>*innz)) *innz=nnz[i];
189da6192ceSVijay Mahadevan       if ((ionz && onz) && (onz[i]>*ionz)) *ionz=onz[i];
190032b8ab6SVijay Mahadevan     }
1915905e1eaSVijay Mahadevan   }
1920df6e276SVijay Mahadevan   PetscFunctionReturn(0);
193032b8ab6SVijay Mahadevan }
194032b8ab6SVijay Mahadevan 
195da8c5b52SVijay Mahadevan 
196da8c5b52SVijay Mahadevan #undef __FUNCT__
197da8c5b52SVijay Mahadevan #define __FUNCT__ "DMMoab_MatFillMatrixEntries_Private"
198da8c5b52SVijay Mahadevan PetscErrorCode DMMoab_MatFillMatrixEntries_Private(DM dm, Mat A)
199da8c5b52SVijay Mahadevan {
200da8c5b52SVijay Mahadevan   DM_Moab                   *dmmoab = (DM_Moab*)dm->data;
201*82dfd14aSVijay Mahadevan   PetscInt                  nconn = 0,prev_nconn = 0,bs,nloc,nfields;
202da8c5b52SVijay Mahadevan   const moab::EntityHandle  *connect;
203da8c5b52SVijay Mahadevan   PetscScalar               *locala=NULL;
204da8c5b52SVijay Mahadevan   PetscInt                  *dof_indices=NULL;
205*82dfd14aSVijay Mahadevan   PetscBool                 isinterlaced;
206*82dfd14aSVijay Mahadevan   char*                     tmp=0;
207*82dfd14aSVijay Mahadevan   MatType                   mtype;
208da8c5b52SVijay Mahadevan   PetscErrorCode            ierr;
209da8c5b52SVijay Mahadevan 
210da8c5b52SVijay Mahadevan   PetscFunctionBegin;
211*82dfd14aSVijay Mahadevan   bs = dmmoab->bs;
212*82dfd14aSVijay Mahadevan   nloc = dmmoab->nloc;
213*82dfd14aSVijay Mahadevan   nfields = dmmoab->numFields;
214*82dfd14aSVijay Mahadevan 
215*82dfd14aSVijay Mahadevan   /* check whether we are updating BAIJ or AIJ matrix */
216*82dfd14aSVijay Mahadevan   ierr = MatGetType(A, &mtype);CHKERRQ(ierr);
217*82dfd14aSVijay Mahadevan   ierr = PetscStrstr(mtype, "baij", &tmp);CHKERRQ(ierr);
218*82dfd14aSVijay Mahadevan   isinterlaced=(tmp || bs==nfields ? PETSC_TRUE : PETSC_FALSE);
219*82dfd14aSVijay Mahadevan 
220da8c5b52SVijay Mahadevan   /* loop over local elements */
221da8c5b52SVijay Mahadevan   for(moab::Range::iterator iter = dmmoab->elocal->begin(); iter != dmmoab->elocal->end(); iter++) {
222da8c5b52SVijay Mahadevan     const moab::EntityHandle ehandle = *iter;
223da8c5b52SVijay Mahadevan 
224da8c5b52SVijay Mahadevan     /* Get connectivity information: */
225da8c5b52SVijay Mahadevan     ierr = DMMoabGetElementConnectivity(dm, ehandle, &nconn, &connect);CHKERRQ(ierr);
226da8c5b52SVijay Mahadevan 
227da8c5b52SVijay Mahadevan     /* if we have mixed elements or arrays have not been initialized - Allocate now */
228da8c5b52SVijay Mahadevan     if (prev_nconn != nconn) {
229da8c5b52SVijay Mahadevan       if (locala) {
230da8c5b52SVijay Mahadevan         ierr = PetscFree(locala);CHKERRQ(ierr);
231da8c5b52SVijay Mahadevan         ierr = PetscFree(dof_indices);CHKERRQ(ierr);
232da8c5b52SVijay Mahadevan       }
233*82dfd14aSVijay Mahadevan       ierr = PetscMalloc(sizeof(PetscScalar)*nconn*nconn*nfields*nfields,&locala);CHKERRQ(ierr);
234*82dfd14aSVijay Mahadevan       ierr = PetscMemzero(locala,sizeof(PetscScalar)*nconn*nconn*nfields*nfields);CHKERRQ(ierr);
235*82dfd14aSVijay Mahadevan       ierr = PetscMalloc(sizeof(PetscInt)*nconn*(isinterlaced?1:nfields),&dof_indices);CHKERRQ(ierr);
236da8c5b52SVijay Mahadevan       prev_nconn=nconn;
237da8c5b52SVijay Mahadevan     }
238da8c5b52SVijay Mahadevan 
239*82dfd14aSVijay Mahadevan     if (isinterlaced) {
240da8c5b52SVijay Mahadevan       /* get the global DOF number to appropriately set the element contribution in the RHS vector */
241da8c5b52SVijay Mahadevan       ierr = DMMoabGetDofsBlockedLocal(dm, nconn, connect, dof_indices);CHKERRQ(ierr);
242da8c5b52SVijay Mahadevan 
243da8c5b52SVijay Mahadevan       /* set the values directly into appropriate locations. Can alternately use VecSetValues */
244da8c5b52SVijay Mahadevan       ierr = MatSetValuesBlockedLocal(A, nconn, dof_indices, nconn, dof_indices, locala, INSERT_VALUES);CHKERRQ(ierr);
245da8c5b52SVijay Mahadevan     }
246*82dfd14aSVijay Mahadevan     else {
247*82dfd14aSVijay Mahadevan       /* get the global DOF number to appropriately set the element contribution in the RHS vector */
248*82dfd14aSVijay Mahadevan       ierr = DMMoabGetDofsLocal(dm, nconn, connect, dof_indices);CHKERRQ(ierr);
249*82dfd14aSVijay Mahadevan 
250*82dfd14aSVijay Mahadevan       /* set the values directly into appropriate locations. Can alternately use VecSetValues */
251*82dfd14aSVijay Mahadevan       ierr = MatSetValuesLocal(A, nconn, dof_indices, nconn, dof_indices, locala, INSERT_VALUES);CHKERRQ(ierr);
252*82dfd14aSVijay Mahadevan     }
253*82dfd14aSVijay Mahadevan   }
254da8c5b52SVijay Mahadevan 
255da8c5b52SVijay Mahadevan   /* clean up memory */
256da8c5b52SVijay Mahadevan   ierr = PetscFree(locala);CHKERRQ(ierr);
257da8c5b52SVijay Mahadevan   ierr = PetscFree(dof_indices);CHKERRQ(ierr);
258da8c5b52SVijay Mahadevan 
259da8c5b52SVijay Mahadevan   /* finish assembly */
260da8c5b52SVijay Mahadevan   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
261da8c5b52SVijay Mahadevan   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
262da8c5b52SVijay Mahadevan   PetscFunctionReturn(0);
263da8c5b52SVijay Mahadevan }
264da8c5b52SVijay Mahadevan 
2655905e1eaSVijay Mahadevan 
2665905e1eaSVijay Mahadevan #undef __FUNCT__
2675905e1eaSVijay Mahadevan #define __FUNCT__ "DMMoabSetBlockFills_Private"
2685905e1eaSVijay Mahadevan static PetscErrorCode DMMoabSetBlockFills_Private(PetscInt w,const PetscInt *fill,PetscInt **rfill)
2695905e1eaSVijay Mahadevan {
2705905e1eaSVijay Mahadevan   PetscErrorCode ierr;
2715905e1eaSVijay Mahadevan   PetscInt       i,j,*ifill;
2725905e1eaSVijay Mahadevan 
2735905e1eaSVijay Mahadevan   PetscFunctionBegin;
2745905e1eaSVijay Mahadevan   if (!fill) PetscFunctionReturn(0);
2755905e1eaSVijay Mahadevan   ierr  = PetscMalloc1(w*w,&ifill);CHKERRQ(ierr);
2765905e1eaSVijay Mahadevan   for (i=0;i<w;i++) {
2775905e1eaSVijay Mahadevan     for (j=0; j<w; j++)
2785905e1eaSVijay Mahadevan       ifill[i*w+j]=fill[i*w+j];
2795905e1eaSVijay Mahadevan   }
2805905e1eaSVijay Mahadevan 
2815905e1eaSVijay Mahadevan   *rfill = ifill;
2825905e1eaSVijay Mahadevan   PetscFunctionReturn(0);
2835905e1eaSVijay Mahadevan }
2845905e1eaSVijay Mahadevan 
2855905e1eaSVijay Mahadevan 
2865905e1eaSVijay Mahadevan #undef __FUNCT__
2875905e1eaSVijay Mahadevan #define __FUNCT__ "DMMoabSetBlockFills"
2885905e1eaSVijay Mahadevan /*@
2895905e1eaSVijay Mahadevan     DMMoabSetBlockFills - Sets the fill pattern in each block for a multi-component problem
2905905e1eaSVijay Mahadevan     of the matrix returned by DMCreateMatrix().
2915905e1eaSVijay Mahadevan 
2925905e1eaSVijay Mahadevan     Logically Collective on DMDA
2935905e1eaSVijay Mahadevan 
2945905e1eaSVijay Mahadevan     Input Parameter:
2955905e1eaSVijay Mahadevan +   dm - the DMMoab object
2965905e1eaSVijay Mahadevan .   dfill - the fill pattern in the diagonal block (may be NULL, means use dense block)
2975905e1eaSVijay Mahadevan -   ofill - the fill pattern in the off-diagonal blocks
2985905e1eaSVijay Mahadevan 
2995905e1eaSVijay Mahadevan     Level: developer
3005905e1eaSVijay Mahadevan 
3015905e1eaSVijay Mahadevan     Notes: This only makes sense when you are doing multicomponent problems but using the
3025905e1eaSVijay Mahadevan        MPIAIJ matrix format
3035905e1eaSVijay Mahadevan 
3045905e1eaSVijay Mahadevan            The format for dfill and ofill is a 2 dimensional dof by dof matrix with 1 entries
3055905e1eaSVijay Mahadevan        representing coupling and 0 entries for missing coupling. For example
3065905e1eaSVijay Mahadevan $             dfill[9] = {1, 0, 0,
3075905e1eaSVijay Mahadevan $                         1, 1, 0,
3085905e1eaSVijay Mahadevan $                         0, 1, 1}
3095905e1eaSVijay Mahadevan        means that row 0 is coupled with only itself in the diagonal block, row 1 is coupled with
3105905e1eaSVijay Mahadevan        itself and row 0 (in the diagonal block) and row 2 is coupled with itself and row 1 (in the
3115905e1eaSVijay Mahadevan        diagonal block).
3125905e1eaSVijay Mahadevan 
3135905e1eaSVijay Mahadevan      DMDASetGetMatrix() allows you to provide general code for those more complicated nonzero patterns then
3145905e1eaSVijay Mahadevan      can be represented in the dfill, ofill format
3155905e1eaSVijay Mahadevan 
3165905e1eaSVijay Mahadevan    Contributed by Glenn Hammond
3175905e1eaSVijay Mahadevan 
3185905e1eaSVijay Mahadevan .seealso DMCreateMatrix(), DMDASetGetMatrix(), DMSetMatrixPreallocateOnly()
3195905e1eaSVijay Mahadevan 
3205905e1eaSVijay Mahadevan @*/
3215905e1eaSVijay Mahadevan PetscErrorCode  DMMoabSetBlockFills(DM dm,const PetscInt *dfill,const PetscInt *ofill)
3225905e1eaSVijay Mahadevan {
3235905e1eaSVijay Mahadevan   DM_Moab       *dmmoab=(DM_Moab*)dm->data;
3245905e1eaSVijay Mahadevan   PetscErrorCode ierr;
3255905e1eaSVijay Mahadevan 
3265905e1eaSVijay Mahadevan   PetscFunctionBegin;
3275905e1eaSVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
3285905e1eaSVijay Mahadevan   ierr = DMMoabSetBlockFills_Private(dmmoab->numFields,dfill,&dmmoab->dfill);CHKERRQ(ierr);
3295905e1eaSVijay Mahadevan   ierr = DMMoabSetBlockFills_Private(dmmoab->numFields,ofill,&dmmoab->ofill);CHKERRQ(ierr);
3305905e1eaSVijay Mahadevan   PetscFunctionReturn(0);
3315905e1eaSVijay Mahadevan }
332