xref: /petsc/src/dm/impls/moab/dmmbmat.cxx (revision 7ae5e5b63e531d1144854fb663ae16bda2a50fcf)
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);
8032b8ab6SVijay Mahadevan 
9032b8ab6SVijay Mahadevan #undef __FUNCT__
10032b8ab6SVijay Mahadevan #define __FUNCT__ "DMCreateMatrix_Moab"
11baf0d1e0SVijay Mahadevan PetscErrorCode DMCreateMatrix_Moab(DM dm,Mat *J)
12032b8ab6SVijay Mahadevan {
13032b8ab6SVijay Mahadevan   PetscErrorCode  ierr;
14c2612ac9SVijay Mahadevan   PetscInt        innz=0,ionz=0,nlsiz;
15032b8ab6SVijay Mahadevan   DM_Moab         *dmmoab=(DM_Moab*)dm->data;
16da6192ceSVijay Mahadevan   PetscInt        *nnz=0,*onz=0;
17032b8ab6SVijay Mahadevan   char            *tmp=0;
182494be97SVijay Mahadevan   Mat             A;
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);
2882dfd14aSVijay 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 */
31*7ae5e5b6SVijay Mahadevan   ierr = PetscCalloc2(nlsiz,&nnz,nlsiz,&onz);CHKERRQ(ierr);
32032b8ab6SVijay Mahadevan 
33032b8ab6SVijay Mahadevan   /* compute the nonzero pattern based on MOAB connectivity data for local elements */
34032b8ab6SVijay Mahadevan   ierr = DMMoab_Compute_NNZ_From_Connectivity(dm,&innz,nnz,&ionz,onz,(tmp?PETSC_TRUE:PETSC_FALSE));CHKERRQ(ierr);
35032b8ab6SVijay Mahadevan 
36da6192ceSVijay Mahadevan   /* create the Matrix and set its type as specified by user */
372494be97SVijay Mahadevan   ierr = MatCreate(dmmoab->pcomm->comm(), &A);CHKERRQ(ierr);
382494be97SVijay Mahadevan   ierr = MatSetSizes(A, dmmoab->nloc*dmmoab->numFields, dmmoab->nloc*dmmoab->numFields, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
392494be97SVijay Mahadevan   ierr = MatSetType(A, mtype);CHKERRQ(ierr);
402494be97SVijay Mahadevan   ierr = MatSetBlockSize(A, dmmoab->bs);CHKERRQ(ierr);
412494be97SVijay Mahadevan   ierr = MatSetDM(A, dm);CHKERRQ(ierr);  /* set DM reference */
422494be97SVijay Mahadevan   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
43da6192ceSVijay Mahadevan 
44da6192ceSVijay Mahadevan   if (!dmmoab->ltog_map) SETERRQ(dmmoab->pcomm->comm(), PETSC_ERR_ORDER, "Cannot create a DMMoab Mat without calling DMSetUp first.");
452494be97SVijay Mahadevan   ierr = MatSetLocalToGlobalMapping(A,dmmoab->ltog_map,dmmoab->ltog_map);CHKERRQ(ierr);
46032b8ab6SVijay Mahadevan 
47032b8ab6SVijay Mahadevan   /* set preallocation based on different supported Mat types */
482494be97SVijay Mahadevan   ierr = MatSeqAIJSetPreallocation(A, innz, nnz);CHKERRQ(ierr);
492494be97SVijay Mahadevan   ierr = MatMPIAIJSetPreallocation(A, innz, nnz, ionz, onz);CHKERRQ(ierr);
502494be97SVijay Mahadevan   ierr = MatSeqBAIJSetPreallocation(A, dmmoab->bs, innz, nnz);CHKERRQ(ierr);
512494be97SVijay Mahadevan   ierr = MatMPIBAIJSetPreallocation(A, dmmoab->bs, innz, nnz, ionz, onz);CHKERRQ(ierr);
52032b8ab6SVijay Mahadevan 
53da8c5b52SVijay Mahadevan   /* clean up temporary memory */
54*7ae5e5b6SVijay Mahadevan   ierr = PetscFree2(nnz,onz);CHKERRQ(ierr);
55da8c5b52SVijay Mahadevan 
56da8c5b52SVijay Mahadevan   /* set up internal matrix data-structures */
572494be97SVijay Mahadevan   ierr = MatSetUp(A);CHKERRQ(ierr);
58da8c5b52SVijay Mahadevan 
592494be97SVijay Mahadevan   *J = A;
60032b8ab6SVijay Mahadevan   PetscFunctionReturn(0);
61032b8ab6SVijay Mahadevan }
62032b8ab6SVijay Mahadevan 
63032b8ab6SVijay Mahadevan 
64032b8ab6SVijay Mahadevan #undef __FUNCT__
65032b8ab6SVijay Mahadevan #define __FUNCT__ "DMMoab_Compute_NNZ_From_Connectivity"
66032b8ab6SVijay Mahadevan PetscErrorCode DMMoab_Compute_NNZ_From_Connectivity(DM dm,PetscInt* innz,PetscInt* nnz,PetscInt* ionz,PetscInt* onz,PetscBool isbaij)
67032b8ab6SVijay Mahadevan {
685905e1eaSVijay Mahadevan   PetscInt        i,f,nloc,vpere,bs,ivtx,n_nnz,n_onz;
69c2612ac9SVijay Mahadevan   PetscInt        ibs,jbs,inbsize,iobsize,nfields,nlsiz;
70032b8ab6SVijay Mahadevan   DM_Moab         *dmmoab = (DM_Moab*)dm->data;
71032b8ab6SVijay Mahadevan   const moab::EntityHandle *connect;
72da6192ceSVijay Mahadevan   moab::Range     adjs,found,allvlocal,allvghost;
73da6192ceSVijay Mahadevan   moab::Range::iterator iter,jter;
74da6192ceSVijay Mahadevan   std::vector<moab::EntityHandle> storage;
755905e1eaSVijay Mahadevan   PetscBool isinterlaced;
76da6192ceSVijay Mahadevan   moab::EntityHandle vtx;
77032b8ab6SVijay Mahadevan   moab::ErrorCode merr;
78032b8ab6SVijay Mahadevan 
79032b8ab6SVijay Mahadevan   PetscFunctionBegin;
80032b8ab6SVijay Mahadevan   bs = dmmoab->bs;
81032b8ab6SVijay Mahadevan   nloc = dmmoab->nloc;
825905e1eaSVijay Mahadevan   nfields = dmmoab->numFields;
8382dfd14aSVijay Mahadevan   isinterlaced=(isbaij || bs==nfields ? PETSC_TRUE : PETSC_FALSE);
84c2612ac9SVijay Mahadevan   nlsiz = (isinterlaced ? nloc:nloc*nfields);
85032b8ab6SVijay Mahadevan 
86da6192ceSVijay Mahadevan   /* find the truly user-expected layer of ghosted entities to decipher NNZ pattern */
87da6192ceSVijay Mahadevan   merr = dmmoab->mbiface->get_entities_by_type(dmmoab->fileset,moab::MBVERTEX,allvlocal,true);MBERRNM(merr);
88da6192ceSVijay Mahadevan   merr = dmmoab->pcomm->filter_pstatus(allvlocal,PSTATUS_NOT_OWNED,PSTATUS_NOT,-1,&adjs);MBERRNM(merr);
89da6192ceSVijay Mahadevan   allvghost = moab::subtract(allvlocal, adjs);
90032b8ab6SVijay Mahadevan 
91f6829af0SVijay Mahadevan   /* loop over the locally owned vertices and figure out the NNZ pattern using connectivity information */
92e427d9c9SVijay Mahadevan   for(iter = dmmoab->vowned->begin(),ivtx=0; iter != dmmoab->vowned->end(); iter++,ivtx++) {
93032b8ab6SVijay Mahadevan 
94da6192ceSVijay Mahadevan     vtx = *iter;
95da6192ceSVijay Mahadevan     adjs.clear();
96da6192ceSVijay Mahadevan     /* Get adjacency information for current vertex - i.e., all elements of dimension (dim) that connects
97da6192ceSVijay Mahadevan        to the current vertex. We can then decipher if a vertex is ghosted or not and compute the
98da6192ceSVijay Mahadevan        non-zero pattern accordingly. */
99da6192ceSVijay Mahadevan     merr = dmmoab->mbiface->get_adjacencies(&vtx,1,dmmoab->dim,false,adjs,moab::Interface::INTERSECT);
1000df6e276SVijay Mahadevan 
101e427d9c9SVijay Mahadevan     /* reset counters */
102e427d9c9SVijay Mahadevan     n_nnz=n_onz=0;
103e427d9c9SVijay Mahadevan     found.clear();
104e427d9c9SVijay Mahadevan 
10515de014fSVijay Mahadevan     /* loop over vertices and update the number of connectivity */
10615de014fSVijay Mahadevan     for(jter = adjs.begin(); jter != adjs.end(); jter++) {
1070df6e276SVijay Mahadevan 
10815de014fSVijay Mahadevan       /* Get connectivity information in canonical ordering for the local element */
10915de014fSVijay Mahadevan       merr = dmmoab->mbiface->get_connectivity(*jter,connect,vpere,false,&storage);MBERRNM(merr);
1100df6e276SVijay Mahadevan 
111f6829af0SVijay Mahadevan       /* loop over each element connected to the adjacent vertex and update as needed */
11215de014fSVijay Mahadevan       for (i=0; i<vpere; ++i) {
113f6829af0SVijay Mahadevan         if (connect[i] == vtx || found.find(connect[i]) != found.end()) continue; /* make sure we don't double count shared vertices */
114f6829af0SVijay Mahadevan         if (allvghost.find(connect[i]) != allvghost.end()) n_onz++; /* update out-of-proc onz */
115f6829af0SVijay Mahadevan         else n_nnz++; /* else local vertex */
11615de014fSVijay Mahadevan         found.insert(connect[i]);
11715de014fSVijay Mahadevan       }
11815de014fSVijay Mahadevan     }
11915de014fSVijay Mahadevan 
12015de014fSVijay Mahadevan     if (isbaij) {
12115de014fSVijay Mahadevan       nnz[ivtx]=n_nnz;    /* leave out self to avoid repeats -> node shared by multiple elements */
1225905e1eaSVijay Mahadevan       if (onz) {
1235905e1eaSVijay Mahadevan         onz[ivtx]=n_onz;  /* add ghost non-owned nodes */
1245905e1eaSVijay Mahadevan       }
12515de014fSVijay Mahadevan     }
126f6829af0SVijay Mahadevan     else { /* AIJ matrices */
1275905e1eaSVijay Mahadevan       if (!isinterlaced) {
1285905e1eaSVijay Mahadevan         for (f=0;f<nfields;f++) {
1295905e1eaSVijay Mahadevan           nnz[f*nloc+ivtx]=n_nnz;    /* leave out self to avoid repeats -> node shared by multiple elements */
1305905e1eaSVijay Mahadevan           if (onz)
1315905e1eaSVijay Mahadevan             onz[f*nloc+ivtx]=n_onz;  /* add ghost non-owned nodes */
1325905e1eaSVijay Mahadevan         }
1335905e1eaSVijay Mahadevan       }
1345905e1eaSVijay Mahadevan       else {
1355905e1eaSVijay Mahadevan         for (f=0;f<nfields;f++) {
1365905e1eaSVijay Mahadevan           nnz[nfields*ivtx+f]=n_nnz;    /* leave out self to avoid repeats -> node shared by multiple elements */
1375905e1eaSVijay Mahadevan           if (onz)
1385905e1eaSVijay Mahadevan             onz[nfields*ivtx+f]=n_onz;  /* add ghost non-owned nodes */
1395905e1eaSVijay Mahadevan         }
1400df6e276SVijay Mahadevan       }
1410df6e276SVijay Mahadevan     }
1420df6e276SVijay Mahadevan   }
1430df6e276SVijay Mahadevan 
144c2612ac9SVijay Mahadevan   for (i=0;i<nlsiz;i++)
1455905e1eaSVijay Mahadevan     nnz[i]+=1;  /* self count the node */
1465905e1eaSVijay Mahadevan 
14782dfd14aSVijay Mahadevan   for (ivtx=0;ivtx<nloc;ivtx++) {
1485905e1eaSVijay Mahadevan     if (!isbaij) {
1495905e1eaSVijay Mahadevan       for (ibs=0; ibs<nfields; ibs++) {
15082dfd14aSVijay Mahadevan         if (dmmoab->dfill) {  /* first address the diagonal block */
1515905e1eaSVijay Mahadevan           /* just add up the ints -- easier/faster rather than branching based on "1" */
1525905e1eaSVijay Mahadevan           for (jbs=0,inbsize=0; jbs<nfields; jbs++)
1535905e1eaSVijay Mahadevan             inbsize += dmmoab->dfill[ibs*nfields+jbs];
1545905e1eaSVijay Mahadevan         }
15582dfd14aSVijay Mahadevan         else inbsize=nfields; /* dense coupling since user didn't specify the component fill explicitly */
15682dfd14aSVijay Mahadevan         if (isinterlaced) nnz[ivtx*nfields+ibs]*=inbsize;
15782dfd14aSVijay Mahadevan         else nnz[ibs*nloc+ivtx]*=inbsize;
1585905e1eaSVijay Mahadevan 
1595905e1eaSVijay Mahadevan         if (onz) {
16082dfd14aSVijay Mahadevan           if (dmmoab->ofill) {  /* next address the off-diagonal block */
1615905e1eaSVijay Mahadevan             /* just add up the ints -- easier/faster rather than branching based on "1" */
1625905e1eaSVijay Mahadevan             for (jbs=0,iobsize=0; jbs<nfields; jbs++)
1635905e1eaSVijay Mahadevan               iobsize += dmmoab->dfill[ibs*nfields+jbs];
1645905e1eaSVijay Mahadevan           }
16582dfd14aSVijay Mahadevan           else iobsize=nfields; /* dense coupling since user didn't specify the component fill explicitly */
16682dfd14aSVijay Mahadevan           if (isinterlaced) onz[ivtx*nfields+ibs]*=iobsize;
16782dfd14aSVijay Mahadevan           else onz[ibs*nloc+ivtx]*=iobsize;
1685905e1eaSVijay Mahadevan         }
1695905e1eaSVijay Mahadevan       }
1705905e1eaSVijay Mahadevan     }
1715905e1eaSVijay Mahadevan     else {
17282dfd14aSVijay Mahadevan       /* check if we got overzealous in our nnz and onz computations */
17382dfd14aSVijay Mahadevan       nnz[ivtx]=(nnz[ivtx]>dmmoab->nloc?dmmoab->nloc:nnz[ivtx]);
17482dfd14aSVijay Mahadevan       if (onz) onz[ivtx]=(onz[ivtx]>dmmoab->nloc?dmmoab->nloc:onz[ivtx]);
1755905e1eaSVijay Mahadevan     }
1765905e1eaSVijay Mahadevan   }
1775905e1eaSVijay Mahadevan   /* update innz and ionz based on local maxima */
1785905e1eaSVijay Mahadevan   if (innz || ionz) {
1790df6e276SVijay Mahadevan     if (innz) *innz=0;
1800df6e276SVijay Mahadevan     if (ionz) *ionz=0;
181c2612ac9SVijay Mahadevan     for (i=0;i<nlsiz;i++) {
182da6192ceSVijay Mahadevan       if (innz && (nnz[i]>*innz)) *innz=nnz[i];
183da6192ceSVijay Mahadevan       if ((ionz && onz) && (onz[i]>*ionz)) *ionz=onz[i];
184032b8ab6SVijay Mahadevan     }
1855905e1eaSVijay Mahadevan   }
1860df6e276SVijay Mahadevan   PetscFunctionReturn(0);
187032b8ab6SVijay Mahadevan }
188032b8ab6SVijay Mahadevan 
189da8c5b52SVijay Mahadevan 
190da8c5b52SVijay Mahadevan #undef __FUNCT__
1915905e1eaSVijay Mahadevan #define __FUNCT__ "DMMoabSetBlockFills_Private"
1925905e1eaSVijay Mahadevan static PetscErrorCode DMMoabSetBlockFills_Private(PetscInt w,const PetscInt *fill,PetscInt **rfill)
1935905e1eaSVijay Mahadevan {
1945905e1eaSVijay Mahadevan   PetscErrorCode ierr;
1955905e1eaSVijay Mahadevan   PetscInt       i,j,*ifill;
1965905e1eaSVijay Mahadevan 
1975905e1eaSVijay Mahadevan   PetscFunctionBegin;
1985905e1eaSVijay Mahadevan   if (!fill) PetscFunctionReturn(0);
1995905e1eaSVijay Mahadevan   ierr  = PetscMalloc1(w*w,&ifill);CHKERRQ(ierr);
2005905e1eaSVijay Mahadevan   for (i=0;i<w;i++) {
2015905e1eaSVijay Mahadevan     for (j=0; j<w; j++)
2025905e1eaSVijay Mahadevan       ifill[i*w+j]=fill[i*w+j];
2035905e1eaSVijay Mahadevan   }
2045905e1eaSVijay Mahadevan 
2055905e1eaSVijay Mahadevan   *rfill = ifill;
2065905e1eaSVijay Mahadevan   PetscFunctionReturn(0);
2075905e1eaSVijay Mahadevan }
2085905e1eaSVijay Mahadevan 
2095905e1eaSVijay Mahadevan 
2105905e1eaSVijay Mahadevan #undef __FUNCT__
2115905e1eaSVijay Mahadevan #define __FUNCT__ "DMMoabSetBlockFills"
2125905e1eaSVijay Mahadevan /*@
2135905e1eaSVijay Mahadevan     DMMoabSetBlockFills - Sets the fill pattern in each block for a multi-component problem
2145905e1eaSVijay Mahadevan     of the matrix returned by DMCreateMatrix().
2155905e1eaSVijay Mahadevan 
2165905e1eaSVijay Mahadevan     Logically Collective on DMDA
2175905e1eaSVijay Mahadevan 
2185905e1eaSVijay Mahadevan     Input Parameter:
2195905e1eaSVijay Mahadevan +   dm - the DMMoab object
2205905e1eaSVijay Mahadevan .   dfill - the fill pattern in the diagonal block (may be NULL, means use dense block)
2215905e1eaSVijay Mahadevan -   ofill - the fill pattern in the off-diagonal blocks
2225905e1eaSVijay Mahadevan 
2235905e1eaSVijay Mahadevan     Level: developer
2245905e1eaSVijay Mahadevan 
2255905e1eaSVijay Mahadevan     Notes: This only makes sense when you are doing multicomponent problems but using the
2265905e1eaSVijay Mahadevan        MPIAIJ matrix format
2275905e1eaSVijay Mahadevan 
2285905e1eaSVijay Mahadevan            The format for dfill and ofill is a 2 dimensional dof by dof matrix with 1 entries
2295905e1eaSVijay Mahadevan        representing coupling and 0 entries for missing coupling. For example
2305905e1eaSVijay Mahadevan $             dfill[9] = {1, 0, 0,
2315905e1eaSVijay Mahadevan $                         1, 1, 0,
2325905e1eaSVijay Mahadevan $                         0, 1, 1}
2335905e1eaSVijay Mahadevan        means that row 0 is coupled with only itself in the diagonal block, row 1 is coupled with
2345905e1eaSVijay Mahadevan        itself and row 0 (in the diagonal block) and row 2 is coupled with itself and row 1 (in the
2355905e1eaSVijay Mahadevan        diagonal block).
2365905e1eaSVijay Mahadevan 
2375905e1eaSVijay Mahadevan      DMDASetGetMatrix() allows you to provide general code for those more complicated nonzero patterns then
2385905e1eaSVijay Mahadevan      can be represented in the dfill, ofill format
2395905e1eaSVijay Mahadevan 
2405905e1eaSVijay Mahadevan    Contributed by Glenn Hammond
2415905e1eaSVijay Mahadevan 
2425905e1eaSVijay Mahadevan .seealso DMCreateMatrix(), DMDASetGetMatrix(), DMSetMatrixPreallocateOnly()
2435905e1eaSVijay Mahadevan 
2445905e1eaSVijay Mahadevan @*/
2455905e1eaSVijay Mahadevan PetscErrorCode  DMMoabSetBlockFills(DM dm,const PetscInt *dfill,const PetscInt *ofill)
2465905e1eaSVijay Mahadevan {
2475905e1eaSVijay Mahadevan   DM_Moab       *dmmoab=(DM_Moab*)dm->data;
2485905e1eaSVijay Mahadevan   PetscErrorCode ierr;
2495905e1eaSVijay Mahadevan 
2505905e1eaSVijay Mahadevan   PetscFunctionBegin;
2515905e1eaSVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2525905e1eaSVijay Mahadevan   ierr = DMMoabSetBlockFills_Private(dmmoab->numFields,dfill,&dmmoab->dfill);CHKERRQ(ierr);
2535905e1eaSVijay Mahadevan   ierr = DMMoabSetBlockFills_Private(dmmoab->numFields,ofill,&dmmoab->ofill);CHKERRQ(ierr);
2545905e1eaSVijay Mahadevan   PetscFunctionReturn(0);
2555905e1eaSVijay Mahadevan }
256