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; 19*2494be97SVijay Mahadevan Mat A; 20baf0d1e0SVijay Mahadevan MatType mtype; 21032b8ab6SVijay Mahadevan 22032b8ab6SVijay Mahadevan PetscFunctionBegin; 23032b8ab6SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 24032b8ab6SVijay Mahadevan PetscValidPointer(J,3); 25032b8ab6SVijay Mahadevan 26db66d124SVijay Mahadevan /* next, need to allocate the non-zero arrays to enable pre-allocation */ 27baf0d1e0SVijay Mahadevan mtype = dm->mattype; 28baf0d1e0SVijay Mahadevan ierr = PetscStrstr(mtype, "baij", &tmp);CHKERRQ(ierr); 2982dfd14aSVijay Mahadevan nlsiz = (tmp ? dmmoab->nloc:dmmoab->nloc*dmmoab->numFields); 30032b8ab6SVijay Mahadevan 31032b8ab6SVijay Mahadevan /* allocate the nnz, onz arrays based on block size and local nodes */ 32da6192ceSVijay Mahadevan ierr = PetscMalloc((nlsiz)*sizeof(PetscInt),&nnz);CHKERRQ(ierr); 33da6192ceSVijay Mahadevan ierr = PetscMemzero(nnz,sizeof(PetscInt)*(nlsiz));CHKERRQ(ierr); 34da6192ceSVijay Mahadevan ierr = PetscMalloc(nlsiz*sizeof(PetscInt),&onz);CHKERRQ(ierr); 35da6192ceSVijay Mahadevan ierr = PetscMemzero(onz,sizeof(PetscInt)*nlsiz);CHKERRQ(ierr); 36032b8ab6SVijay Mahadevan 37032b8ab6SVijay Mahadevan /* compute the nonzero pattern based on MOAB connectivity data for local elements */ 38032b8ab6SVijay Mahadevan ierr = DMMoab_Compute_NNZ_From_Connectivity(dm,&innz,nnz,&ionz,onz,(tmp?PETSC_TRUE:PETSC_FALSE));CHKERRQ(ierr); 39032b8ab6SVijay Mahadevan 40da6192ceSVijay Mahadevan /* create the Matrix and set its type as specified by user */ 41*2494be97SVijay Mahadevan ierr = MatCreate(dmmoab->pcomm->comm(), &A);CHKERRQ(ierr); 42*2494be97SVijay Mahadevan ierr = MatSetSizes(A, dmmoab->nloc*dmmoab->numFields, dmmoab->nloc*dmmoab->numFields, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 43*2494be97SVijay Mahadevan ierr = MatSetType(A, mtype);CHKERRQ(ierr); 44*2494be97SVijay Mahadevan ierr = MatSetBlockSize(A, dmmoab->bs);CHKERRQ(ierr); 45*2494be97SVijay Mahadevan ierr = MatSetDM(A, dm);CHKERRQ(ierr); /* set DM reference */ 46*2494be97SVijay Mahadevan ierr = MatSetFromOptions(A);CHKERRQ(ierr); 47da6192ceSVijay Mahadevan 48da6192ceSVijay Mahadevan if (!dmmoab->ltog_map) SETERRQ(dmmoab->pcomm->comm(), PETSC_ERR_ORDER, "Cannot create a DMMoab Mat without calling DMSetUp first."); 49*2494be97SVijay Mahadevan ierr = MatSetLocalToGlobalMapping(A,dmmoab->ltog_map,dmmoab->ltog_map);CHKERRQ(ierr); 50032b8ab6SVijay Mahadevan 51032b8ab6SVijay Mahadevan /* set preallocation based on different supported Mat types */ 52*2494be97SVijay Mahadevan ierr = MatSeqAIJSetPreallocation(A, innz, nnz);CHKERRQ(ierr); 53*2494be97SVijay Mahadevan ierr = MatMPIAIJSetPreallocation(A, innz, nnz, ionz, onz);CHKERRQ(ierr); 54*2494be97SVijay Mahadevan ierr = MatSeqBAIJSetPreallocation(A, dmmoab->bs, innz, nnz);CHKERRQ(ierr); 55*2494be97SVijay Mahadevan ierr = MatMPIBAIJSetPreallocation(A, dmmoab->bs, innz, nnz, ionz, onz);CHKERRQ(ierr); 56032b8ab6SVijay Mahadevan 57da8c5b52SVijay Mahadevan /* clean up temporary memory */ 58032b8ab6SVijay Mahadevan ierr = PetscFree(nnz);CHKERRQ(ierr); 59032b8ab6SVijay Mahadevan ierr = PetscFree(onz);CHKERRQ(ierr); 60da8c5b52SVijay Mahadevan 61da8c5b52SVijay Mahadevan /* set up internal matrix data-structures */ 62*2494be97SVijay Mahadevan ierr = MatSetUp(A);CHKERRQ(ierr); 63da8c5b52SVijay Mahadevan 64da8c5b52SVijay Mahadevan /* set the correct NNZ pattern by setting matrix entries - make the matrix ready to use */ 65*2494be97SVijay Mahadevan ierr = DMMoab_MatFillMatrixEntries_Private(dm,A);CHKERRQ(ierr); 66*2494be97SVijay Mahadevan 67*2494be97SVijay Mahadevan *J = A; 68032b8ab6SVijay Mahadevan PetscFunctionReturn(0); 69032b8ab6SVijay Mahadevan } 70032b8ab6SVijay Mahadevan 71032b8ab6SVijay Mahadevan 72032b8ab6SVijay Mahadevan #undef __FUNCT__ 73032b8ab6SVijay Mahadevan #define __FUNCT__ "DMMoab_Compute_NNZ_From_Connectivity" 74032b8ab6SVijay Mahadevan PetscErrorCode DMMoab_Compute_NNZ_From_Connectivity(DM dm,PetscInt* innz,PetscInt* nnz,PetscInt* ionz,PetscInt* onz,PetscBool isbaij) 75032b8ab6SVijay Mahadevan { 765905e1eaSVijay Mahadevan PetscInt i,f,nloc,vpere,bs,ivtx,n_nnz,n_onz; 775905e1eaSVijay Mahadevan PetscInt ibs,jbs,inbsize,iobsize,nfields; 78032b8ab6SVijay Mahadevan DM_Moab *dmmoab = (DM_Moab*)dm->data; 79032b8ab6SVijay Mahadevan const moab::EntityHandle *connect; 80da6192ceSVijay Mahadevan moab::Range adjs,found,allvlocal,allvghost; 81da6192ceSVijay Mahadevan moab::Range::iterator iter,jter; 82da6192ceSVijay Mahadevan std::vector<moab::EntityHandle> storage; 835905e1eaSVijay Mahadevan PetscBool isinterlaced; 84da6192ceSVijay Mahadevan moab::EntityHandle vtx; 85032b8ab6SVijay Mahadevan moab::ErrorCode merr; 86032b8ab6SVijay Mahadevan 87032b8ab6SVijay Mahadevan PetscFunctionBegin; 88032b8ab6SVijay Mahadevan bs = dmmoab->bs; 89032b8ab6SVijay Mahadevan nloc = dmmoab->nloc; 905905e1eaSVijay Mahadevan nfields = dmmoab->numFields; 9182dfd14aSVijay Mahadevan isinterlaced=(isbaij || bs==nfields ? PETSC_TRUE : PETSC_FALSE); 92032b8ab6SVijay Mahadevan 93da6192ceSVijay Mahadevan /* find the truly user-expected layer of ghosted entities to decipher NNZ pattern */ 94da6192ceSVijay Mahadevan merr = dmmoab->mbiface->get_entities_by_type(dmmoab->fileset,moab::MBVERTEX,allvlocal,true);MBERRNM(merr); 95da6192ceSVijay Mahadevan merr = dmmoab->pcomm->filter_pstatus(allvlocal,PSTATUS_NOT_OWNED,PSTATUS_NOT,-1,&adjs);MBERRNM(merr); 96da6192ceSVijay Mahadevan allvghost = moab::subtract(allvlocal, adjs); 97032b8ab6SVijay Mahadevan 98f6829af0SVijay Mahadevan /* loop over the locally owned vertices and figure out the NNZ pattern using connectivity information */ 99e427d9c9SVijay Mahadevan for(iter = dmmoab->vowned->begin(),ivtx=0; iter != dmmoab->vowned->end(); iter++,ivtx++) { 100032b8ab6SVijay Mahadevan 101da6192ceSVijay Mahadevan vtx = *iter; 102da6192ceSVijay Mahadevan adjs.clear(); 103da6192ceSVijay Mahadevan /* Get adjacency information for current vertex - i.e., all elements of dimension (dim) that connects 104da6192ceSVijay Mahadevan to the current vertex. We can then decipher if a vertex is ghosted or not and compute the 105da6192ceSVijay Mahadevan non-zero pattern accordingly. */ 106da6192ceSVijay Mahadevan merr = dmmoab->mbiface->get_adjacencies(&vtx,1,dmmoab->dim,false,adjs,moab::Interface::INTERSECT); 1070df6e276SVijay Mahadevan 108e427d9c9SVijay Mahadevan /* reset counters */ 109e427d9c9SVijay Mahadevan n_nnz=n_onz=0; 110e427d9c9SVijay Mahadevan found.clear(); 111e427d9c9SVijay Mahadevan 11215de014fSVijay Mahadevan /* loop over vertices and update the number of connectivity */ 11315de014fSVijay Mahadevan for(jter = adjs.begin(); jter != adjs.end(); jter++) { 1140df6e276SVijay Mahadevan 11515de014fSVijay Mahadevan /* Get connectivity information in canonical ordering for the local element */ 11615de014fSVijay Mahadevan merr = dmmoab->mbiface->get_connectivity(*jter,connect,vpere,false,&storage);MBERRNM(merr); 1170df6e276SVijay Mahadevan 118f6829af0SVijay Mahadevan /* loop over each element connected to the adjacent vertex and update as needed */ 11915de014fSVijay Mahadevan for (i=0; i<vpere; ++i) { 120f6829af0SVijay Mahadevan if (connect[i] == vtx || found.find(connect[i]) != found.end()) continue; /* make sure we don't double count shared vertices */ 121f6829af0SVijay Mahadevan if (allvghost.find(connect[i]) != allvghost.end()) n_onz++; /* update out-of-proc onz */ 122f6829af0SVijay Mahadevan else n_nnz++; /* else local vertex */ 12315de014fSVijay Mahadevan found.insert(connect[i]); 12415de014fSVijay Mahadevan } 12515de014fSVijay Mahadevan } 12615de014fSVijay Mahadevan 12715de014fSVijay Mahadevan if (isbaij) { 12815de014fSVijay Mahadevan nnz[ivtx]=n_nnz; /* leave out self to avoid repeats -> node shared by multiple elements */ 1295905e1eaSVijay Mahadevan if (onz) { 1305905e1eaSVijay Mahadevan onz[ivtx]=n_onz; /* add ghost non-owned nodes */ 1315905e1eaSVijay Mahadevan } 13215de014fSVijay Mahadevan } 133f6829af0SVijay Mahadevan else { /* AIJ matrices */ 1345905e1eaSVijay Mahadevan if (!isinterlaced) { 1355905e1eaSVijay Mahadevan for (f=0;f<nfields;f++) { 1365905e1eaSVijay Mahadevan nnz[f*nloc+ivtx]=n_nnz; /* leave out self to avoid repeats -> node shared by multiple elements */ 1375905e1eaSVijay Mahadevan if (onz) 1385905e1eaSVijay Mahadevan onz[f*nloc+ivtx]=n_onz; /* add ghost non-owned nodes */ 1395905e1eaSVijay Mahadevan } 1405905e1eaSVijay Mahadevan } 1415905e1eaSVijay Mahadevan else { 1425905e1eaSVijay Mahadevan for (f=0;f<nfields;f++) { 1435905e1eaSVijay Mahadevan nnz[nfields*ivtx+f]=n_nnz; /* leave out self to avoid repeats -> node shared by multiple elements */ 1445905e1eaSVijay Mahadevan if (onz) 1455905e1eaSVijay Mahadevan onz[nfields*ivtx+f]=n_onz; /* add ghost non-owned nodes */ 1465905e1eaSVijay Mahadevan } 1470df6e276SVijay Mahadevan } 1480df6e276SVijay Mahadevan } 1490df6e276SVijay Mahadevan } 1500df6e276SVijay Mahadevan 15182dfd14aSVijay Mahadevan for (i=0;i<nloc*(isbaij?1:nfields);i++) 1525905e1eaSVijay Mahadevan nnz[i]+=1; /* self count the node */ 1535905e1eaSVijay Mahadevan 15482dfd14aSVijay Mahadevan for (ivtx=0;ivtx<nloc;ivtx++) { 1555905e1eaSVijay Mahadevan if (!isbaij) { 1565905e1eaSVijay Mahadevan for (ibs=0; ibs<nfields; ibs++) { 15782dfd14aSVijay Mahadevan if (dmmoab->dfill) { /* first address the diagonal block */ 1585905e1eaSVijay Mahadevan /* just add up the ints -- easier/faster rather than branching based on "1" */ 1595905e1eaSVijay Mahadevan for (jbs=0,inbsize=0; jbs<nfields; jbs++) 1605905e1eaSVijay Mahadevan inbsize += dmmoab->dfill[ibs*nfields+jbs]; 1615905e1eaSVijay Mahadevan } 16282dfd14aSVijay Mahadevan else inbsize=nfields; /* dense coupling since user didn't specify the component fill explicitly */ 16382dfd14aSVijay Mahadevan if (isinterlaced) nnz[ivtx*nfields+ibs]*=inbsize; 16482dfd14aSVijay Mahadevan else nnz[ibs*nloc+ivtx]*=inbsize; 1655905e1eaSVijay Mahadevan 1665905e1eaSVijay Mahadevan if (onz) { 16782dfd14aSVijay Mahadevan if (dmmoab->ofill) { /* next address the off-diagonal block */ 1685905e1eaSVijay Mahadevan /* just add up the ints -- easier/faster rather than branching based on "1" */ 1695905e1eaSVijay Mahadevan for (jbs=0,iobsize=0; jbs<nfields; jbs++) 1705905e1eaSVijay Mahadevan iobsize += dmmoab->dfill[ibs*nfields+jbs]; 1715905e1eaSVijay Mahadevan } 17282dfd14aSVijay Mahadevan else iobsize=nfields; /* dense coupling since user didn't specify the component fill explicitly */ 17382dfd14aSVijay Mahadevan if (isinterlaced) onz[ivtx*nfields+ibs]*=iobsize; 17482dfd14aSVijay Mahadevan else onz[ibs*nloc+ivtx]*=iobsize; 1755905e1eaSVijay Mahadevan } 1765905e1eaSVijay Mahadevan } 1775905e1eaSVijay Mahadevan } 1785905e1eaSVijay Mahadevan else { 17982dfd14aSVijay Mahadevan /* check if we got overzealous in our nnz and onz computations */ 18082dfd14aSVijay Mahadevan nnz[ivtx]=(nnz[ivtx]>dmmoab->nloc?dmmoab->nloc:nnz[ivtx]); 18182dfd14aSVijay Mahadevan if (onz) onz[ivtx]=(onz[ivtx]>dmmoab->nloc?dmmoab->nloc:onz[ivtx]); 1825905e1eaSVijay Mahadevan } 1835905e1eaSVijay Mahadevan } 1845905e1eaSVijay Mahadevan /* update innz and ionz based on local maxima */ 1855905e1eaSVijay Mahadevan if (innz || ionz) { 1860df6e276SVijay Mahadevan if (innz) *innz=0; 1870df6e276SVijay Mahadevan if (ionz) *ionz=0; 1885905e1eaSVijay Mahadevan for (i=0;i<nloc*nfields;i++) { 189da6192ceSVijay Mahadevan if (innz && (nnz[i]>*innz)) *innz=nnz[i]; 190da6192ceSVijay Mahadevan if ((ionz && onz) && (onz[i]>*ionz)) *ionz=onz[i]; 191032b8ab6SVijay Mahadevan } 1925905e1eaSVijay Mahadevan } 1930df6e276SVijay Mahadevan PetscFunctionReturn(0); 194032b8ab6SVijay Mahadevan } 195032b8ab6SVijay Mahadevan 196da8c5b52SVijay Mahadevan 197da8c5b52SVijay Mahadevan #undef __FUNCT__ 198da8c5b52SVijay Mahadevan #define __FUNCT__ "DMMoab_MatFillMatrixEntries_Private" 199da8c5b52SVijay Mahadevan PetscErrorCode DMMoab_MatFillMatrixEntries_Private(DM dm, Mat A) 200da8c5b52SVijay Mahadevan { 201da8c5b52SVijay Mahadevan DM_Moab *dmmoab = (DM_Moab*)dm->data; 20282dfd14aSVijay Mahadevan PetscInt nconn = 0,prev_nconn = 0,bs,nloc,nfields; 203da8c5b52SVijay Mahadevan const moab::EntityHandle *connect; 204da8c5b52SVijay Mahadevan PetscScalar *locala=NULL; 205da8c5b52SVijay Mahadevan PetscInt *dof_indices=NULL; 20682dfd14aSVijay Mahadevan PetscBool isinterlaced; 20782dfd14aSVijay Mahadevan char* tmp=0; 20882dfd14aSVijay Mahadevan MatType mtype; 209da8c5b52SVijay Mahadevan PetscErrorCode ierr; 210da8c5b52SVijay Mahadevan 211da8c5b52SVijay Mahadevan PetscFunctionBegin; 21282dfd14aSVijay Mahadevan bs = dmmoab->bs; 21382dfd14aSVijay Mahadevan nloc = dmmoab->nloc; 21482dfd14aSVijay Mahadevan nfields = dmmoab->numFields; 21582dfd14aSVijay Mahadevan 21682dfd14aSVijay Mahadevan /* check whether we are updating BAIJ or AIJ matrix */ 21782dfd14aSVijay Mahadevan ierr = MatGetType(A, &mtype);CHKERRQ(ierr); 21882dfd14aSVijay Mahadevan ierr = PetscStrstr(mtype, "baij", &tmp);CHKERRQ(ierr); 21982dfd14aSVijay Mahadevan isinterlaced=(tmp || bs==nfields ? PETSC_TRUE : PETSC_FALSE); 22082dfd14aSVijay Mahadevan 221da8c5b52SVijay Mahadevan /* loop over local elements */ 222da8c5b52SVijay Mahadevan for(moab::Range::iterator iter = dmmoab->elocal->begin(); iter != dmmoab->elocal->end(); iter++) { 223da8c5b52SVijay Mahadevan const moab::EntityHandle ehandle = *iter; 224da8c5b52SVijay Mahadevan 225da8c5b52SVijay Mahadevan /* Get connectivity information: */ 226da8c5b52SVijay Mahadevan ierr = DMMoabGetElementConnectivity(dm, ehandle, &nconn, &connect);CHKERRQ(ierr); 227da8c5b52SVijay Mahadevan 228da8c5b52SVijay Mahadevan /* if we have mixed elements or arrays have not been initialized - Allocate now */ 229da8c5b52SVijay Mahadevan if (prev_nconn != nconn) { 230da8c5b52SVijay Mahadevan if (locala) { 231da8c5b52SVijay Mahadevan ierr = PetscFree(locala);CHKERRQ(ierr); 232da8c5b52SVijay Mahadevan ierr = PetscFree(dof_indices);CHKERRQ(ierr); 233da8c5b52SVijay Mahadevan } 23482dfd14aSVijay Mahadevan ierr = PetscMalloc(sizeof(PetscScalar)*nconn*nconn*nfields*nfields,&locala);CHKERRQ(ierr); 23582dfd14aSVijay Mahadevan ierr = PetscMemzero(locala,sizeof(PetscScalar)*nconn*nconn*nfields*nfields);CHKERRQ(ierr); 23682dfd14aSVijay Mahadevan ierr = PetscMalloc(sizeof(PetscInt)*nconn*(isinterlaced?1:nfields),&dof_indices);CHKERRQ(ierr); 237da8c5b52SVijay Mahadevan prev_nconn=nconn; 238da8c5b52SVijay Mahadevan } 239da8c5b52SVijay Mahadevan 24082dfd14aSVijay Mahadevan if (isinterlaced) { 241da8c5b52SVijay Mahadevan /* get the global DOF number to appropriately set the element contribution in the RHS vector */ 242da8c5b52SVijay Mahadevan ierr = DMMoabGetDofsBlockedLocal(dm, nconn, connect, dof_indices);CHKERRQ(ierr); 243da8c5b52SVijay Mahadevan 244da8c5b52SVijay Mahadevan /* set the values directly into appropriate locations. Can alternately use VecSetValues */ 245da8c5b52SVijay Mahadevan ierr = MatSetValuesBlockedLocal(A, nconn, dof_indices, nconn, dof_indices, locala, INSERT_VALUES);CHKERRQ(ierr); 246da8c5b52SVijay Mahadevan } 24782dfd14aSVijay Mahadevan else { 24882dfd14aSVijay Mahadevan /* get the global DOF number to appropriately set the element contribution in the RHS vector */ 24982dfd14aSVijay Mahadevan ierr = DMMoabGetDofsLocal(dm, nconn, connect, dof_indices);CHKERRQ(ierr); 25082dfd14aSVijay Mahadevan 25182dfd14aSVijay Mahadevan /* set the values directly into appropriate locations. Can alternately use VecSetValues */ 25282dfd14aSVijay Mahadevan ierr = MatSetValuesLocal(A, nconn, dof_indices, nconn, dof_indices, locala, INSERT_VALUES);CHKERRQ(ierr); 25382dfd14aSVijay Mahadevan } 25482dfd14aSVijay Mahadevan } 255da8c5b52SVijay Mahadevan 256da8c5b52SVijay Mahadevan /* clean up memory */ 257da8c5b52SVijay Mahadevan ierr = PetscFree(locala);CHKERRQ(ierr); 258da8c5b52SVijay Mahadevan ierr = PetscFree(dof_indices);CHKERRQ(ierr); 259da8c5b52SVijay Mahadevan 260da8c5b52SVijay Mahadevan /* finish assembly */ 261da8c5b52SVijay Mahadevan ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 262da8c5b52SVijay Mahadevan ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 263da8c5b52SVijay Mahadevan PetscFunctionReturn(0); 264da8c5b52SVijay Mahadevan } 265da8c5b52SVijay Mahadevan 2665905e1eaSVijay Mahadevan 2675905e1eaSVijay Mahadevan #undef __FUNCT__ 2685905e1eaSVijay Mahadevan #define __FUNCT__ "DMMoabSetBlockFills_Private" 2695905e1eaSVijay Mahadevan static PetscErrorCode DMMoabSetBlockFills_Private(PetscInt w,const PetscInt *fill,PetscInt **rfill) 2705905e1eaSVijay Mahadevan { 2715905e1eaSVijay Mahadevan PetscErrorCode ierr; 2725905e1eaSVijay Mahadevan PetscInt i,j,*ifill; 2735905e1eaSVijay Mahadevan 2745905e1eaSVijay Mahadevan PetscFunctionBegin; 2755905e1eaSVijay Mahadevan if (!fill) PetscFunctionReturn(0); 2765905e1eaSVijay Mahadevan ierr = PetscMalloc1(w*w,&ifill);CHKERRQ(ierr); 2775905e1eaSVijay Mahadevan for (i=0;i<w;i++) { 2785905e1eaSVijay Mahadevan for (j=0; j<w; j++) 2795905e1eaSVijay Mahadevan ifill[i*w+j]=fill[i*w+j]; 2805905e1eaSVijay Mahadevan } 2815905e1eaSVijay Mahadevan 2825905e1eaSVijay Mahadevan *rfill = ifill; 2835905e1eaSVijay Mahadevan PetscFunctionReturn(0); 2845905e1eaSVijay Mahadevan } 2855905e1eaSVijay Mahadevan 2865905e1eaSVijay Mahadevan 2875905e1eaSVijay Mahadevan #undef __FUNCT__ 2885905e1eaSVijay Mahadevan #define __FUNCT__ "DMMoabSetBlockFills" 2895905e1eaSVijay Mahadevan /*@ 2905905e1eaSVijay Mahadevan DMMoabSetBlockFills - Sets the fill pattern in each block for a multi-component problem 2915905e1eaSVijay Mahadevan of the matrix returned by DMCreateMatrix(). 2925905e1eaSVijay Mahadevan 2935905e1eaSVijay Mahadevan Logically Collective on DMDA 2945905e1eaSVijay Mahadevan 2955905e1eaSVijay Mahadevan Input Parameter: 2965905e1eaSVijay Mahadevan + dm - the DMMoab object 2975905e1eaSVijay Mahadevan . dfill - the fill pattern in the diagonal block (may be NULL, means use dense block) 2985905e1eaSVijay Mahadevan - ofill - the fill pattern in the off-diagonal blocks 2995905e1eaSVijay Mahadevan 3005905e1eaSVijay Mahadevan Level: developer 3015905e1eaSVijay Mahadevan 3025905e1eaSVijay Mahadevan Notes: This only makes sense when you are doing multicomponent problems but using the 3035905e1eaSVijay Mahadevan MPIAIJ matrix format 3045905e1eaSVijay Mahadevan 3055905e1eaSVijay Mahadevan The format for dfill and ofill is a 2 dimensional dof by dof matrix with 1 entries 3065905e1eaSVijay Mahadevan representing coupling and 0 entries for missing coupling. For example 3075905e1eaSVijay Mahadevan $ dfill[9] = {1, 0, 0, 3085905e1eaSVijay Mahadevan $ 1, 1, 0, 3095905e1eaSVijay Mahadevan $ 0, 1, 1} 3105905e1eaSVijay Mahadevan means that row 0 is coupled with only itself in the diagonal block, row 1 is coupled with 3115905e1eaSVijay Mahadevan itself and row 0 (in the diagonal block) and row 2 is coupled with itself and row 1 (in the 3125905e1eaSVijay Mahadevan diagonal block). 3135905e1eaSVijay Mahadevan 3145905e1eaSVijay Mahadevan DMDASetGetMatrix() allows you to provide general code for those more complicated nonzero patterns then 3155905e1eaSVijay Mahadevan can be represented in the dfill, ofill format 3165905e1eaSVijay Mahadevan 3175905e1eaSVijay Mahadevan Contributed by Glenn Hammond 3185905e1eaSVijay Mahadevan 3195905e1eaSVijay Mahadevan .seealso DMCreateMatrix(), DMDASetGetMatrix(), DMSetMatrixPreallocateOnly() 3205905e1eaSVijay Mahadevan 3215905e1eaSVijay Mahadevan @*/ 3225905e1eaSVijay Mahadevan PetscErrorCode DMMoabSetBlockFills(DM dm,const PetscInt *dfill,const PetscInt *ofill) 3235905e1eaSVijay Mahadevan { 3245905e1eaSVijay Mahadevan DM_Moab *dmmoab=(DM_Moab*)dm->data; 3255905e1eaSVijay Mahadevan PetscErrorCode ierr; 3265905e1eaSVijay Mahadevan 3275905e1eaSVijay Mahadevan PetscFunctionBegin; 3285905e1eaSVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3295905e1eaSVijay Mahadevan ierr = DMMoabSetBlockFills_Private(dmmoab->numFields,dfill,&dmmoab->dfill);CHKERRQ(ierr); 3305905e1eaSVijay Mahadevan ierr = DMMoabSetBlockFills_Private(dmmoab->numFields,ofill,&dmmoab->ofill);CHKERRQ(ierr); 3315905e1eaSVijay Mahadevan PetscFunctionReturn(0); 3325905e1eaSVijay Mahadevan } 333