xref: /petsc/src/dm/impls/moab/dmmbmat.cxx (revision b117cd09e766aba90b52b0bef2091bbb861995da)
1af0996ceSBarry Smith #include <petsc/private/dmmbimpl.h> /*I  "petscdmmoab.h"   I*/
2af0996ceSBarry Smith #include <petsc/private/vecimpl.h>
3032b8ab6SVijay Mahadevan 
4032b8ab6SVijay Mahadevan #include <petscdmmoab.h>
5da6192ceSVijay Mahadevan #include <MBTagConventions.hpp>
6*b117cd09SVijay Mahadevan #include <moab/NestedRefine.hpp>
7032b8ab6SVijay Mahadevan 
8032b8ab6SVijay Mahadevan static PetscErrorCode DMMoab_Compute_NNZ_From_Connectivity(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscBool);
9032b8ab6SVijay Mahadevan 
10baf0d1e0SVijay Mahadevan PetscErrorCode DMCreateMatrix_Moab(DM dm,Mat *J)
11032b8ab6SVijay Mahadevan {
12032b8ab6SVijay Mahadevan   PetscErrorCode  ierr;
13c2612ac9SVijay Mahadevan   PetscInt        innz=0,ionz=0,nlsiz;
14032b8ab6SVijay Mahadevan   DM_Moab         *dmmoab=(DM_Moab*)dm->data;
15da6192ceSVijay Mahadevan   PetscInt        *nnz=0,*onz=0;
16032b8ab6SVijay Mahadevan   char            *tmp=0;
172494be97SVijay Mahadevan   Mat             A;
18baf0d1e0SVijay Mahadevan   MatType         mtype;
19032b8ab6SVijay Mahadevan 
20032b8ab6SVijay Mahadevan   PetscFunctionBegin;
21032b8ab6SVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
22032b8ab6SVijay Mahadevan   PetscValidPointer(J,3);
23032b8ab6SVijay Mahadevan 
24db66d124SVijay Mahadevan   /* next, need to allocate the non-zero arrays to enable pre-allocation */
25baf0d1e0SVijay Mahadevan   mtype = dm->mattype;
26baf0d1e0SVijay Mahadevan   ierr = PetscStrstr(mtype, "baij", &tmp);CHKERRQ(ierr);
2782dfd14aSVijay Mahadevan   nlsiz = (tmp ? dmmoab->nloc:dmmoab->nloc*dmmoab->numFields);
28032b8ab6SVijay Mahadevan 
29032b8ab6SVijay Mahadevan   /* allocate the nnz, onz arrays based on block size and local nodes */
307ae5e5b6SVijay Mahadevan   ierr = PetscCalloc2(nlsiz,&nnz,nlsiz,&onz);CHKERRQ(ierr);
31032b8ab6SVijay Mahadevan 
32032b8ab6SVijay Mahadevan   /* compute the nonzero pattern based on MOAB connectivity data for local elements */
33032b8ab6SVijay Mahadevan   ierr = DMMoab_Compute_NNZ_From_Connectivity(dm,&innz,nnz,&ionz,onz,(tmp?PETSC_TRUE:PETSC_FALSE));CHKERRQ(ierr);
34032b8ab6SVijay Mahadevan 
35da6192ceSVijay Mahadevan   /* create the Matrix and set its type as specified by user */
362494be97SVijay Mahadevan   ierr = MatCreate(dmmoab->pcomm->comm(), &A);CHKERRQ(ierr);
372494be97SVijay Mahadevan   ierr = MatSetSizes(A, dmmoab->nloc*dmmoab->numFields, dmmoab->nloc*dmmoab->numFields, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
382494be97SVijay Mahadevan   ierr = MatSetType(A, mtype);CHKERRQ(ierr);
392494be97SVijay Mahadevan   ierr = MatSetBlockSize(A, dmmoab->bs);CHKERRQ(ierr);
402494be97SVijay Mahadevan   ierr = MatSetDM(A, dm);CHKERRQ(ierr);  /* set DM reference */
412494be97SVijay Mahadevan   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
42da6192ceSVijay Mahadevan 
43da6192ceSVijay Mahadevan   if (!dmmoab->ltog_map) SETERRQ(dmmoab->pcomm->comm(), PETSC_ERR_ORDER, "Cannot create a DMMoab Mat without calling DMSetUp first.");
442494be97SVijay Mahadevan   ierr = MatSetLocalToGlobalMapping(A,dmmoab->ltog_map,dmmoab->ltog_map);CHKERRQ(ierr);
45032b8ab6SVijay Mahadevan 
46032b8ab6SVijay Mahadevan   /* set preallocation based on different supported Mat types */
472494be97SVijay Mahadevan   ierr = MatSeqAIJSetPreallocation(A, innz, nnz);CHKERRQ(ierr);
482494be97SVijay Mahadevan   ierr = MatMPIAIJSetPreallocation(A, innz, nnz, ionz, onz);CHKERRQ(ierr);
492494be97SVijay Mahadevan   ierr = MatSeqBAIJSetPreallocation(A, dmmoab->bs, innz, nnz);CHKERRQ(ierr);
502494be97SVijay Mahadevan   ierr = MatMPIBAIJSetPreallocation(A, dmmoab->bs, innz, nnz, ionz, onz);CHKERRQ(ierr);
51032b8ab6SVijay Mahadevan 
52da8c5b52SVijay Mahadevan   /* clean up temporary memory */
537ae5e5b6SVijay Mahadevan   ierr = PetscFree2(nnz,onz);CHKERRQ(ierr);
54da8c5b52SVijay Mahadevan 
55da8c5b52SVijay Mahadevan   /* set up internal matrix data-structures */
562494be97SVijay Mahadevan   ierr = MatSetUp(A);CHKERRQ(ierr);
57da8c5b52SVijay Mahadevan 
582494be97SVijay Mahadevan   *J = A;
59032b8ab6SVijay Mahadevan   PetscFunctionReturn(0);
60032b8ab6SVijay Mahadevan }
61032b8ab6SVijay Mahadevan 
62032b8ab6SVijay Mahadevan 
63032b8ab6SVijay Mahadevan PetscErrorCode DMMoab_Compute_NNZ_From_Connectivity(DM dm,PetscInt* innz,PetscInt* nnz,PetscInt* ionz,PetscInt* onz,PetscBool isbaij)
64032b8ab6SVijay Mahadevan {
65*b117cd09SVijay Mahadevan   PetscInt        i,f,nloc,vpere,bs,n_nnz,n_onz,ivtx=0;
66c2612ac9SVijay Mahadevan   PetscInt        ibs,jbs,inbsize,iobsize,nfields,nlsiz;
67032b8ab6SVijay Mahadevan   DM_Moab         *dmmoab = (DM_Moab*)dm->data;
68032b8ab6SVijay Mahadevan   const moab::EntityHandle *connect;
69da6192ceSVijay Mahadevan   moab::Range     adjs,found,allvlocal,allvghost;
70da6192ceSVijay Mahadevan   std::vector<moab::EntityHandle> storage;
715905e1eaSVijay Mahadevan   PetscBool isinterlaced;
72da6192ceSVijay Mahadevan   moab::EntityHandle vtx;
73032b8ab6SVijay Mahadevan   moab::ErrorCode merr;
74032b8ab6SVijay Mahadevan 
75032b8ab6SVijay Mahadevan   PetscFunctionBegin;
76032b8ab6SVijay Mahadevan   bs = dmmoab->bs;
77032b8ab6SVijay Mahadevan   nloc = dmmoab->nloc;
785905e1eaSVijay Mahadevan   nfields = dmmoab->numFields;
7982dfd14aSVijay Mahadevan   isinterlaced=(isbaij || bs==nfields ? PETSC_TRUE : PETSC_FALSE);
80c2612ac9SVijay Mahadevan   nlsiz = (isinterlaced ? nloc:nloc*nfields);
81032b8ab6SVijay Mahadevan 
82da6192ceSVijay Mahadevan   /* find the truly user-expected layer of ghosted entities to decipher NNZ pattern */
83*b117cd09SVijay Mahadevan   merr = dmmoab->mbiface->get_entities_by_type(dmmoab->fileset,moab::MBVERTEX,allvlocal/*,true*/);MBERRNM(merr);
84da6192ceSVijay Mahadevan   merr = dmmoab->pcomm->filter_pstatus(allvlocal,PSTATUS_NOT_OWNED,PSTATUS_NOT,-1,&adjs);MBERRNM(merr);
85da6192ceSVijay Mahadevan   allvghost = moab::subtract(allvlocal, adjs);
86032b8ab6SVijay Mahadevan 
87f6829af0SVijay Mahadevan   /* loop over the locally owned vertices and figure out the NNZ pattern using connectivity information */
88*b117cd09SVijay Mahadevan   for(moab::Range::iterator iter = dmmoab->vowned->begin(); iter != dmmoab->vowned->end(); iter++,ivtx++) {
89032b8ab6SVijay Mahadevan 
90da6192ceSVijay Mahadevan     vtx = *iter;
91da6192ceSVijay Mahadevan     adjs.clear();
92da6192ceSVijay Mahadevan     /* Get adjacency information for current vertex - i.e., all elements of dimension (dim) that connects
93da6192ceSVijay Mahadevan        to the current vertex. We can then decipher if a vertex is ghosted or not and compute the
94da6192ceSVijay Mahadevan        non-zero pattern accordingly. */
95*b117cd09SVijay Mahadevan     std::vector<moab::EntityHandle> adjstl;
96*b117cd09SVijay Mahadevan     if (dmmoab->hierarchy) {
97*b117cd09SVijay Mahadevan       merr = dmmoab->hierarchy->get_adjacencies(vtx,dmmoab->dim,adjstl);MBERRNM(merr);
98*b117cd09SVijay Mahadevan     }
99*b117cd09SVijay Mahadevan     else {
100*b117cd09SVijay Mahadevan       merr = dmmoab->mbiface->get_adjacencies(&vtx,1,dmmoab->dim,true,adjs,moab::Interface::INTERSECT);MBERRNM(merr);
101*b117cd09SVijay Mahadevan     }
1020df6e276SVijay Mahadevan 
103e427d9c9SVijay Mahadevan     /* reset counters */
104e427d9c9SVijay Mahadevan     n_nnz=n_onz=0;
105e427d9c9SVijay Mahadevan     found.clear();
106e427d9c9SVijay Mahadevan 
10715de014fSVijay Mahadevan     /* loop over vertices and update the number of connectivity */
108*b117cd09SVijay Mahadevan     for(unsigned jter = 0; jter < (dmmoab->hierarchy ? adjstl.size(): adjs.size()); jter++) {
109*b117cd09SVijay Mahadevan 
110*b117cd09SVijay Mahadevan       moab::EntityHandle jhandle;
111*b117cd09SVijay Mahadevan       if (dmmoab->hierarchy) jhandle = adjstl[jter];
112*b117cd09SVijay Mahadevan       else jhandle = adjs[jter];
1130df6e276SVijay Mahadevan 
11415de014fSVijay Mahadevan       /* Get connectivity information in canonical ordering for the local element */
115*b117cd09SVijay Mahadevan       merr = dmmoab->mbiface->get_connectivity(jhandle,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 
150c2612ac9SVijay Mahadevan   for (i=0;i<nlsiz;i++)
1515905e1eaSVijay Mahadevan     nnz[i]+=1;  /* self count the node */
1525905e1eaSVijay Mahadevan 
15382dfd14aSVijay Mahadevan   for (ivtx=0;ivtx<nloc;ivtx++) {
1545905e1eaSVijay Mahadevan     if (!isbaij) {
1555905e1eaSVijay Mahadevan       for (ibs=0; ibs<nfields; ibs++) {
15682dfd14aSVijay 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         }
16182dfd14aSVijay Mahadevan         else inbsize=nfields; /* dense coupling since user didn't specify the component fill explicitly */
16282dfd14aSVijay Mahadevan         if (isinterlaced) nnz[ivtx*nfields+ibs]*=inbsize;
16382dfd14aSVijay Mahadevan         else nnz[ibs*nloc+ivtx]*=inbsize;
1645905e1eaSVijay Mahadevan 
1655905e1eaSVijay Mahadevan         if (onz) {
16682dfd14aSVijay 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           }
17182dfd14aSVijay Mahadevan           else iobsize=nfields; /* dense coupling since user didn't specify the component fill explicitly */
17282dfd14aSVijay Mahadevan           if (isinterlaced) onz[ivtx*nfields+ibs]*=iobsize;
17382dfd14aSVijay Mahadevan           else onz[ibs*nloc+ivtx]*=iobsize;
1745905e1eaSVijay Mahadevan         }
1755905e1eaSVijay Mahadevan       }
1765905e1eaSVijay Mahadevan     }
1775905e1eaSVijay Mahadevan     else {
17882dfd14aSVijay Mahadevan       /* check if we got overzealous in our nnz and onz computations */
17982dfd14aSVijay Mahadevan       nnz[ivtx]=(nnz[ivtx]>dmmoab->nloc?dmmoab->nloc:nnz[ivtx]);
18082dfd14aSVijay 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;
187c2612ac9SVijay Mahadevan     for (i=0;i<nlsiz;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 
1965905e1eaSVijay Mahadevan static PetscErrorCode DMMoabSetBlockFills_Private(PetscInt w,const PetscInt *fill,PetscInt **rfill)
1975905e1eaSVijay Mahadevan {
1985905e1eaSVijay Mahadevan   PetscErrorCode ierr;
1995905e1eaSVijay Mahadevan   PetscInt       i,j,*ifill;
2005905e1eaSVijay Mahadevan 
2015905e1eaSVijay Mahadevan   PetscFunctionBegin;
2025905e1eaSVijay Mahadevan   if (!fill) PetscFunctionReturn(0);
2035905e1eaSVijay Mahadevan   ierr  = PetscMalloc1(w*w,&ifill);CHKERRQ(ierr);
2045905e1eaSVijay Mahadevan   for (i=0;i<w;i++) {
2055905e1eaSVijay Mahadevan     for (j=0; j<w; j++)
2065905e1eaSVijay Mahadevan       ifill[i*w+j]=fill[i*w+j];
2075905e1eaSVijay Mahadevan   }
2085905e1eaSVijay Mahadevan 
2095905e1eaSVijay Mahadevan   *rfill = ifill;
2105905e1eaSVijay Mahadevan   PetscFunctionReturn(0);
2115905e1eaSVijay Mahadevan }
2125905e1eaSVijay Mahadevan 
2135905e1eaSVijay Mahadevan 
2145905e1eaSVijay Mahadevan /*@
2155905e1eaSVijay Mahadevan     DMMoabSetBlockFills - Sets the fill pattern in each block for a multi-component problem
2165905e1eaSVijay Mahadevan     of the matrix returned by DMCreateMatrix().
2175905e1eaSVijay Mahadevan 
2185905e1eaSVijay Mahadevan     Logically Collective on DMDA
2195905e1eaSVijay Mahadevan 
2205905e1eaSVijay Mahadevan     Input Parameter:
2215905e1eaSVijay Mahadevan +   dm - the DMMoab object
2225905e1eaSVijay Mahadevan .   dfill - the fill pattern in the diagonal block (may be NULL, means use dense block)
2235905e1eaSVijay Mahadevan -   ofill - the fill pattern in the off-diagonal blocks
2245905e1eaSVijay Mahadevan 
2255905e1eaSVijay Mahadevan     Level: developer
2265905e1eaSVijay Mahadevan 
2275905e1eaSVijay Mahadevan     Notes: This only makes sense when you are doing multicomponent problems but using the
2285905e1eaSVijay Mahadevan        MPIAIJ matrix format
2295905e1eaSVijay Mahadevan 
2305905e1eaSVijay Mahadevan            The format for dfill and ofill is a 2 dimensional dof by dof matrix with 1 entries
2315905e1eaSVijay Mahadevan        representing coupling and 0 entries for missing coupling. For example
2325905e1eaSVijay Mahadevan $             dfill[9] = {1, 0, 0,
2335905e1eaSVijay Mahadevan $                         1, 1, 0,
2345905e1eaSVijay Mahadevan $                         0, 1, 1}
2355905e1eaSVijay Mahadevan        means that row 0 is coupled with only itself in the diagonal block, row 1 is coupled with
2365905e1eaSVijay Mahadevan        itself and row 0 (in the diagonal block) and row 2 is coupled with itself and row 1 (in the
2375905e1eaSVijay Mahadevan        diagonal block).
2385905e1eaSVijay Mahadevan 
2395905e1eaSVijay Mahadevan      DMDASetGetMatrix() allows you to provide general code for those more complicated nonzero patterns then
2405905e1eaSVijay Mahadevan      can be represented in the dfill, ofill format
2415905e1eaSVijay Mahadevan 
2425905e1eaSVijay Mahadevan    Contributed by Glenn Hammond
2435905e1eaSVijay Mahadevan 
2445905e1eaSVijay Mahadevan .seealso DMCreateMatrix(), DMDASetGetMatrix(), DMSetMatrixPreallocateOnly()
2455905e1eaSVijay Mahadevan 
2465905e1eaSVijay Mahadevan @*/
2475905e1eaSVijay Mahadevan PetscErrorCode  DMMoabSetBlockFills(DM dm,const PetscInt *dfill,const PetscInt *ofill)
2485905e1eaSVijay Mahadevan {
2495905e1eaSVijay Mahadevan   DM_Moab       *dmmoab=(DM_Moab*)dm->data;
2505905e1eaSVijay Mahadevan   PetscErrorCode ierr;
2515905e1eaSVijay Mahadevan 
2525905e1eaSVijay Mahadevan   PetscFunctionBegin;
2535905e1eaSVijay Mahadevan   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2545905e1eaSVijay Mahadevan   ierr = DMMoabSetBlockFills_Private(dmmoab->numFields,dfill,&dmmoab->dfill);CHKERRQ(ierr);
2555905e1eaSVijay Mahadevan   ierr = DMMoabSetBlockFills_Private(dmmoab->numFields,ofill,&dmmoab->ofill);CHKERRQ(ierr);
2565905e1eaSVijay Mahadevan   PetscFunctionReturn(0);
2575905e1eaSVijay Mahadevan }
258