xref: /petsc/src/dm/impls/moab/dmmbmat.cxx (revision a86ed394d0a6627fd6ef9fd16cfa3d0db6a8bc8a)
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>
6b117cd09SVijay Mahadevan #include <moab/NestedRefine.hpp>
7032b8ab6SVijay Mahadevan 
8*a86ed394SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoab_Compute_NNZ_From_Connectivity(DM, PetscInt*, PetscInt*, PetscInt*, PetscInt*, PetscBool);
9032b8ab6SVijay Mahadevan 
10304006b3SVijay Mahadevan #undef __FUNCT__
11304006b3SVijay Mahadevan #define __FUNCT__ "DMCreateMatrix_Moab"
12304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMCreateMatrix_Moab(DM dm, Mat *J)
13032b8ab6SVijay Mahadevan {
14032b8ab6SVijay Mahadevan   PetscErrorCode  ierr;
15c2612ac9SVijay Mahadevan   PetscInt        innz = 0, ionz = 0, nlsiz;
16032b8ab6SVijay Mahadevan   DM_Moab         *dmmoab = (DM_Moab*)dm->data;
17da6192ceSVijay Mahadevan   PetscInt        *nnz = 0, *onz = 0;
18032b8ab6SVijay Mahadevan   char            *tmp = 0;
192494be97SVijay 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 */
327ae5e5b6SVijay Mahadevan   ierr = PetscCalloc2(nlsiz, &nnz, nlsiz, &onz);CHKERRQ(ierr);
33032b8ab6SVijay Mahadevan 
34032b8ab6SVijay Mahadevan   /* compute the nonzero pattern based on MOAB connectivity data for local elements */
35032b8ab6SVijay Mahadevan   ierr = DMMoab_Compute_NNZ_From_Connectivity(dm, &innz, nnz, &ionz, onz, (tmp ? PETSC_TRUE : PETSC_FALSE));CHKERRQ(ierr);
36032b8ab6SVijay Mahadevan 
37da6192ceSVijay Mahadevan   /* create the Matrix and set its type as specified by user */
389daf19fdSVijay Mahadevan   ierr = MatCreate((((PetscObject)dm)->comm), &A);CHKERRQ(ierr);
392494be97SVijay Mahadevan   ierr = MatSetSizes(A, dmmoab->nloc * dmmoab->numFields, dmmoab->nloc * dmmoab->numFields, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
402494be97SVijay Mahadevan   ierr = MatSetType(A, mtype);CHKERRQ(ierr);
412494be97SVijay Mahadevan   ierr = MatSetBlockSize(A, dmmoab->bs);CHKERRQ(ierr);
422494be97SVijay Mahadevan   ierr = MatSetDM(A, dm);CHKERRQ(ierr); /* set DM reference */
432494be97SVijay Mahadevan   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
44da6192ceSVijay Mahadevan 
459daf19fdSVijay Mahadevan   if (!dmmoab->ltog_map) SETERRQ((((PetscObject)dm)->comm), PETSC_ERR_ORDER, "Cannot create a DMMoab Mat without calling DMSetUp first.");
462494be97SVijay Mahadevan   ierr = MatSetLocalToGlobalMapping(A, dmmoab->ltog_map, dmmoab->ltog_map);CHKERRQ(ierr);
47032b8ab6SVijay Mahadevan 
48032b8ab6SVijay Mahadevan   /* set preallocation based on different supported Mat types */
492494be97SVijay Mahadevan   ierr = MatSeqAIJSetPreallocation(A, innz, nnz);CHKERRQ(ierr);
502494be97SVijay Mahadevan   ierr = MatMPIAIJSetPreallocation(A, innz, nnz, ionz, onz);CHKERRQ(ierr);
512494be97SVijay Mahadevan   ierr = MatSeqBAIJSetPreallocation(A, dmmoab->bs, innz, nnz);CHKERRQ(ierr);
522494be97SVijay Mahadevan   ierr = MatMPIBAIJSetPreallocation(A, dmmoab->bs, innz, nnz, ionz, onz);CHKERRQ(ierr);
53032b8ab6SVijay Mahadevan 
54da8c5b52SVijay Mahadevan   /* clean up temporary memory */
557ae5e5b6SVijay Mahadevan   ierr = PetscFree2(nnz, onz);CHKERRQ(ierr);
56da8c5b52SVijay Mahadevan 
57da8c5b52SVijay Mahadevan   /* set up internal matrix data-structures */
582494be97SVijay Mahadevan   ierr = MatSetUp(A);CHKERRQ(ierr);
59da8c5b52SVijay Mahadevan 
602494be97SVijay Mahadevan   *J = A;
61032b8ab6SVijay Mahadevan   PetscFunctionReturn(0);
62032b8ab6SVijay Mahadevan }
63032b8ab6SVijay Mahadevan 
64032b8ab6SVijay Mahadevan 
6563d025dbSVijay Mahadevan #undef __FUNCT__
6663d025dbSVijay Mahadevan #define __FUNCT__ "DMMoab_Compute_NNZ_From_Connectivity"
67*a86ed394SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoab_Compute_NNZ_From_Connectivity(DM dm, PetscInt* innz, PetscInt* nnz, PetscInt* ionz, PetscInt* onz, PetscBool isbaij)
68032b8ab6SVijay Mahadevan {
69b117cd09SVijay Mahadevan   PetscInt        i, f, nloc, vpere, bs, n_nnz, n_onz, ivtx = 0;
70c2612ac9SVijay Mahadevan   PetscInt        ibs, jbs, inbsize, iobsize, nfields, nlsiz;
71032b8ab6SVijay Mahadevan   DM_Moab         *dmmoab = (DM_Moab*)dm->data;
72032b8ab6SVijay Mahadevan   const moab::EntityHandle *connect;
7364e1c140SVijay Mahadevan   moab::Range     adjs, found;
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 
86f6829af0SVijay Mahadevan   /* loop over the locally owned vertices and figure out the NNZ pattern using connectivity information */
8764e1c140SVijay Mahadevan   for (moab::Range::const_iterator iter = dmmoab->vowned->begin(); iter != dmmoab->vowned->end(); iter++, ivtx++) {
88032b8ab6SVijay Mahadevan 
89da6192ceSVijay Mahadevan     vtx = *iter;
90da6192ceSVijay Mahadevan     /* Get adjacency information for current vertex - i.e., all elements of dimension (dim) that connects
91da6192ceSVijay Mahadevan        to the current vertex. We can then decipher if a vertex is ghosted or not and compute the
92da6192ceSVijay Mahadevan        non-zero pattern accordingly. */
9364e1c140SVijay Mahadevan     adjs.clear();
94b117cd09SVijay Mahadevan     merr = dmmoab->mbiface->get_adjacencies(&vtx, 1, dmmoab->dim, true, adjs, moab::Interface::INTERSECT); MBERRNM(merr);
950df6e276SVijay Mahadevan 
96e427d9c9SVijay Mahadevan     /* reset counters */
97e427d9c9SVijay Mahadevan     n_nnz = n_onz = 0;
98e427d9c9SVijay Mahadevan     found.clear();
99e427d9c9SVijay Mahadevan 
10015de014fSVijay Mahadevan     /* loop over vertices and update the number of connectivity */
10164e1c140SVijay Mahadevan     for (moab::Range::const_iterator jter = adjs.begin(); jter != adjs.end(); jter++) {
102b117cd09SVijay Mahadevan 
10364e1c140SVijay Mahadevan       const moab::EntityHandle jhandle = *jter;
1040df6e276SVijay Mahadevan 
10515de014fSVijay Mahadevan       /* Get connectivity information in canonical ordering for the local element */
106b117cd09SVijay Mahadevan       merr = dmmoab->mbiface->get_connectivity(jhandle, connect, vpere, false, &storage); MBERRNM(merr);
1070df6e276SVijay Mahadevan 
108f6829af0SVijay Mahadevan       /* loop over each element connected to the adjacent vertex and update as needed */
10915de014fSVijay Mahadevan       for (i = 0; i < vpere; ++i) {
11064e1c140SVijay Mahadevan         /* find the truly user-expected layer of ghosted entities to decipher NNZ pattern */
111f6829af0SVijay Mahadevan         if (connect[i] == vtx || found.find(connect[i]) != found.end()) continue; /* make sure we don't double count shared vertices */
11264e1c140SVijay Mahadevan         if (dmmoab->vghost->find(connect[i]) != dmmoab->vghost->end()) n_onz++; /* update out-of-proc onz */
113f6829af0SVijay Mahadevan         else n_nnz++; /* else local vertex */
11415de014fSVijay Mahadevan         found.insert(connect[i]);
11515de014fSVijay Mahadevan       }
11615de014fSVijay Mahadevan     }
1179daf19fdSVijay Mahadevan     storage.clear();
11815de014fSVijay Mahadevan 
11915de014fSVijay Mahadevan     if (isbaij) {
12015de014fSVijay Mahadevan       nnz[ivtx] = n_nnz;  /* leave out self to avoid repeats -> node shared by multiple elements */
1215905e1eaSVijay Mahadevan       if (onz) {
1225905e1eaSVijay Mahadevan         onz[ivtx] = n_onz; /* add ghost non-owned nodes */
1235905e1eaSVijay Mahadevan       }
12415de014fSVijay Mahadevan     }
125f6829af0SVijay Mahadevan     else { /* AIJ matrices */
1265905e1eaSVijay Mahadevan       if (!isinterlaced) {
1275905e1eaSVijay Mahadevan         for (f = 0; f < nfields; f++) {
1285905e1eaSVijay Mahadevan           nnz[f * nloc + ivtx] = n_nnz; /* leave out self to avoid repeats -> node shared by multiple elements */
1295905e1eaSVijay Mahadevan           if (onz)
1305905e1eaSVijay Mahadevan             onz[f * nloc + ivtx] = n_onz; /* add ghost non-owned nodes */
1315905e1eaSVijay Mahadevan         }
1325905e1eaSVijay Mahadevan       }
1335905e1eaSVijay Mahadevan       else {
1345905e1eaSVijay Mahadevan         for (f = 0; f < nfields; f++) {
1355905e1eaSVijay Mahadevan           nnz[nfields * ivtx + f] = n_nnz; /* leave out self to avoid repeats -> node shared by multiple elements */
1365905e1eaSVijay Mahadevan           if (onz)
1375905e1eaSVijay Mahadevan             onz[nfields * ivtx + f] = n_onz; /* add ghost non-owned nodes */
1385905e1eaSVijay Mahadevan         }
1390df6e276SVijay Mahadevan       }
1400df6e276SVijay Mahadevan     }
1410df6e276SVijay Mahadevan   }
1420df6e276SVijay Mahadevan 
143c2612ac9SVijay Mahadevan   for (i = 0; i < nlsiz; i++)
1445905e1eaSVijay Mahadevan     nnz[i] += 1; /* self count the node */
1455905e1eaSVijay Mahadevan 
14682dfd14aSVijay Mahadevan   for (ivtx = 0; ivtx < nloc; ivtx++) {
1475905e1eaSVijay Mahadevan     if (!isbaij) {
1485905e1eaSVijay Mahadevan       for (ibs = 0; ibs < nfields; ibs++) {
14982dfd14aSVijay Mahadevan         if (dmmoab->dfill) {  /* first address the diagonal block */
1505905e1eaSVijay Mahadevan           /* just add up the ints -- easier/faster rather than branching based on "1" */
1515905e1eaSVijay Mahadevan           for (jbs = 0, inbsize = 0; jbs < nfields; jbs++)
1525905e1eaSVijay Mahadevan             inbsize += dmmoab->dfill[ibs * nfields + jbs];
1535905e1eaSVijay Mahadevan         }
15482dfd14aSVijay Mahadevan         else inbsize = nfields; /* dense coupling since user didn't specify the component fill explicitly */
15582dfd14aSVijay Mahadevan         if (isinterlaced) nnz[ivtx * nfields + ibs] *= inbsize;
15682dfd14aSVijay Mahadevan         else nnz[ibs * nloc + ivtx] *= inbsize;
1575905e1eaSVijay Mahadevan 
1585905e1eaSVijay Mahadevan         if (onz) {
15982dfd14aSVijay Mahadevan           if (dmmoab->ofill) {  /* next address the off-diagonal block */
1605905e1eaSVijay Mahadevan             /* just add up the ints -- easier/faster rather than branching based on "1" */
1615905e1eaSVijay Mahadevan             for (jbs = 0, iobsize = 0; jbs < nfields; jbs++)
1625905e1eaSVijay Mahadevan               iobsize += dmmoab->dfill[ibs * nfields + jbs];
1635905e1eaSVijay Mahadevan           }
16482dfd14aSVijay Mahadevan           else iobsize = nfields; /* dense coupling since user didn't specify the component fill explicitly */
16582dfd14aSVijay Mahadevan           if (isinterlaced) onz[ivtx * nfields + ibs] *= iobsize;
16682dfd14aSVijay Mahadevan           else onz[ibs * nloc + ivtx] *= iobsize;
1675905e1eaSVijay Mahadevan         }
1685905e1eaSVijay Mahadevan       }
1695905e1eaSVijay Mahadevan     }
1705905e1eaSVijay Mahadevan     else {
17182dfd14aSVijay Mahadevan       /* check if we got overzealous in our nnz and onz computations */
17282dfd14aSVijay Mahadevan       nnz[ivtx] = (nnz[ivtx] > dmmoab->nloc ? dmmoab->nloc : nnz[ivtx]);
17382dfd14aSVijay Mahadevan       if (onz) onz[ivtx] = (onz[ivtx] > dmmoab->nloc ? dmmoab->nloc : onz[ivtx]);
1745905e1eaSVijay Mahadevan     }
1755905e1eaSVijay Mahadevan   }
1765905e1eaSVijay Mahadevan   /* update innz and ionz based on local maxima */
1775905e1eaSVijay Mahadevan   if (innz || ionz) {
1780df6e276SVijay Mahadevan     if (innz) *innz = 0;
1790df6e276SVijay Mahadevan     if (ionz) *ionz = 0;
1809daf19fdSVijay Mahadevan     if (isbaij) {
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     }
1869daf19fdSVijay Mahadevan     else {
1879daf19fdSVijay Mahadevan       for (i = 0; i < nlsiz * nfields; i++) {
1889daf19fdSVijay Mahadevan         if (innz && (nnz[i] > *innz)) *innz = nnz[i];
1899daf19fdSVijay Mahadevan         if ((ionz && onz) && (onz[i] > *ionz)) *ionz = onz[i];
1909daf19fdSVijay Mahadevan       }
1919daf19fdSVijay Mahadevan     }
1929daf19fdSVijay Mahadevan   }
1930df6e276SVijay Mahadevan   PetscFunctionReturn(0);
194032b8ab6SVijay Mahadevan }
195032b8ab6SVijay Mahadevan 
196da8c5b52SVijay Mahadevan 
19763d025dbSVijay Mahadevan #undef __FUNCT__
19863d025dbSVijay Mahadevan #define __FUNCT__ "DMMoabSetBlockFills_Private"
1995905e1eaSVijay Mahadevan static PetscErrorCode DMMoabSetBlockFills_Private(PetscInt w, const PetscInt *fill, PetscInt **rfill)
2005905e1eaSVijay Mahadevan {
2015905e1eaSVijay Mahadevan   PetscErrorCode ierr;
2025905e1eaSVijay Mahadevan   PetscInt       i, j, *ifill;
2035905e1eaSVijay Mahadevan 
2045905e1eaSVijay Mahadevan   PetscFunctionBegin;
2055905e1eaSVijay Mahadevan   if (!fill) PetscFunctionReturn(0);
2065905e1eaSVijay Mahadevan   ierr  = PetscMalloc1(w * w, &ifill);CHKERRQ(ierr);
2075905e1eaSVijay Mahadevan   for (i = 0; i < w; i++) {
2085905e1eaSVijay Mahadevan     for (j = 0; j < w; j++)
2095905e1eaSVijay Mahadevan       ifill[i * w + j] = fill[i * w + j];
2105905e1eaSVijay Mahadevan   }
2115905e1eaSVijay Mahadevan 
2125905e1eaSVijay Mahadevan   *rfill = ifill;
2135905e1eaSVijay Mahadevan   PetscFunctionReturn(0);
2145905e1eaSVijay Mahadevan }
2155905e1eaSVijay Mahadevan 
2165905e1eaSVijay Mahadevan 
2175905e1eaSVijay Mahadevan /*@
2185905e1eaSVijay Mahadevan     DMMoabSetBlockFills - Sets the fill pattern in each block for a multi-component problem
2195905e1eaSVijay Mahadevan     of the matrix returned by DMCreateMatrix().
2205905e1eaSVijay Mahadevan 
2215905e1eaSVijay Mahadevan     Logically Collective on DMDA
2225905e1eaSVijay Mahadevan 
2235905e1eaSVijay Mahadevan     Input Parameter:
2245905e1eaSVijay Mahadevan +   dm - the DMMoab object
2255905e1eaSVijay Mahadevan .   dfill - the fill pattern in the diagonal block (may be NULL, means use dense block)
2265905e1eaSVijay Mahadevan -   ofill - the fill pattern in the off-diagonal blocks
2275905e1eaSVijay Mahadevan 
2285905e1eaSVijay Mahadevan     Level: developer
2295905e1eaSVijay Mahadevan 
2305905e1eaSVijay Mahadevan     Notes: This only makes sense when you are doing multicomponent problems but using the
2315905e1eaSVijay Mahadevan        MPIAIJ matrix format
2325905e1eaSVijay Mahadevan 
2335905e1eaSVijay Mahadevan            The format for dfill and ofill is a 2 dimensional dof by dof matrix with 1 entries
2345905e1eaSVijay Mahadevan        representing coupling and 0 entries for missing coupling. For example
2355905e1eaSVijay Mahadevan $             dfill[9] = {1, 0, 0,
2365905e1eaSVijay Mahadevan $                         1, 1, 0,
2375905e1eaSVijay Mahadevan $                         0, 1, 1}
2385905e1eaSVijay Mahadevan        means that row 0 is coupled with only itself in the diagonal block, row 1 is coupled with
2395905e1eaSVijay Mahadevan        itself and row 0 (in the diagonal block) and row 2 is coupled with itself and row 1 (in the
2405905e1eaSVijay Mahadevan        diagonal block).
2415905e1eaSVijay Mahadevan 
2425905e1eaSVijay Mahadevan      DMDASetGetMatrix() allows you to provide general code for those more complicated nonzero patterns then
2435905e1eaSVijay Mahadevan      can be represented in the dfill, ofill format
2445905e1eaSVijay Mahadevan 
2455905e1eaSVijay Mahadevan    Contributed by Glenn Hammond
2465905e1eaSVijay Mahadevan 
2475905e1eaSVijay Mahadevan .seealso DMCreateMatrix(), DMDASetGetMatrix(), DMSetMatrixPreallocateOnly()
2485905e1eaSVijay Mahadevan 
2495905e1eaSVijay Mahadevan @*/
2505905e1eaSVijay Mahadevan PetscErrorCode  DMMoabSetBlockFills(DM dm, const PetscInt *dfill, const PetscInt *ofill)
2515905e1eaSVijay Mahadevan {
2525905e1eaSVijay Mahadevan   DM_Moab       *dmmoab = (DM_Moab*)dm->data;
2535905e1eaSVijay Mahadevan   PetscErrorCode ierr;
2545905e1eaSVijay Mahadevan 
2555905e1eaSVijay Mahadevan   PetscFunctionBegin;
2565905e1eaSVijay Mahadevan   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2575905e1eaSVijay Mahadevan   ierr = DMMoabSetBlockFills_Private(dmmoab->numFields, dfill, &dmmoab->dfill);CHKERRQ(ierr);
2585905e1eaSVijay Mahadevan   ierr = DMMoabSetBlockFills_Private(dmmoab->numFields, ofill, &dmmoab->ofill);CHKERRQ(ierr);
2595905e1eaSVijay Mahadevan   PetscFunctionReturn(0);
2605905e1eaSVijay Mahadevan }
261