xref: /petsc/src/dm/impls/moab/dmmbmat.cxx (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
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 
8a86ed394SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoab_Compute_NNZ_From_Connectivity(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscBool);
9032b8ab6SVijay Mahadevan 
10*9371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMCreateMatrix_Moab(DM dm, Mat *J) {
11c2612ac9SVijay Mahadevan   PetscInt  innz = 0, ionz = 0, nlsiz;
12032b8ab6SVijay Mahadevan   DM_Moab  *dmmoab = (DM_Moab *)dm->data;
13da6192ceSVijay Mahadevan   PetscInt *nnz = 0, *onz = 0;
14032b8ab6SVijay Mahadevan   char     *tmp = 0;
152494be97SVijay Mahadevan   Mat       A;
16baf0d1e0SVijay Mahadevan   MatType   mtype;
17032b8ab6SVijay Mahadevan 
18032b8ab6SVijay Mahadevan   PetscFunctionBegin;
19032b8ab6SVijay Mahadevan   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
20032b8ab6SVijay Mahadevan   PetscValidPointer(J, 3);
21032b8ab6SVijay Mahadevan 
22db66d124SVijay Mahadevan   /* next, need to allocate the non-zero arrays to enable pre-allocation */
23baf0d1e0SVijay Mahadevan   mtype = dm->mattype;
249566063dSJacob Faibussowitsch   PetscCall(PetscStrstr(mtype, MATBAIJ, &tmp));
2582dfd14aSVijay Mahadevan   nlsiz = (tmp ? dmmoab->nloc : dmmoab->nloc * dmmoab->numFields);
26032b8ab6SVijay Mahadevan 
27032b8ab6SVijay Mahadevan   /* allocate the nnz, onz arrays based on block size and local nodes */
289566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(nlsiz, &nnz, nlsiz, &onz));
29032b8ab6SVijay Mahadevan 
30032b8ab6SVijay Mahadevan   /* compute the nonzero pattern based on MOAB connectivity data for local elements */
319566063dSJacob Faibussowitsch   PetscCall(DMMoab_Compute_NNZ_From_Connectivity(dm, &innz, nnz, &ionz, onz, (tmp ? PETSC_TRUE : PETSC_FALSE)));
32032b8ab6SVijay Mahadevan 
33da6192ceSVijay Mahadevan   /* create the Matrix and set its type as specified by user */
349566063dSJacob Faibussowitsch   PetscCall(MatCreate((((PetscObject)dm)->comm), &A));
359566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, dmmoab->nloc * dmmoab->numFields, dmmoab->nloc * dmmoab->numFields, PETSC_DETERMINE, PETSC_DETERMINE));
369566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, mtype));
379566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSize(A, dmmoab->bs));
389566063dSJacob Faibussowitsch   PetscCall(MatSetDM(A, dm)); /* set DM reference */
399566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
40da6192ceSVijay Mahadevan 
417a8be351SBarry Smith   PetscCheck(dmmoab->ltog_map, (((PetscObject)dm)->comm), PETSC_ERR_ORDER, "Cannot create a DMMoab Mat without calling DMSetUp first.");
429566063dSJacob Faibussowitsch   PetscCall(MatSetLocalToGlobalMapping(A, dmmoab->ltog_map, dmmoab->ltog_map));
43032b8ab6SVijay Mahadevan 
44032b8ab6SVijay Mahadevan   /* set preallocation based on different supported Mat types */
459566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(A, innz, nnz));
469566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(A, innz, nnz, ionz, onz));
479566063dSJacob Faibussowitsch   PetscCall(MatSeqBAIJSetPreallocation(A, dmmoab->bs, innz, nnz));
489566063dSJacob Faibussowitsch   PetscCall(MatMPIBAIJSetPreallocation(A, dmmoab->bs, innz, nnz, ionz, onz));
49032b8ab6SVijay Mahadevan 
50da8c5b52SVijay Mahadevan   /* clean up temporary memory */
519566063dSJacob Faibussowitsch   PetscCall(PetscFree2(nnz, onz));
52da8c5b52SVijay Mahadevan 
53da8c5b52SVijay Mahadevan   /* set up internal matrix data-structures */
549566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
55da8c5b52SVijay Mahadevan 
56e92d1c7cSVijay Mahadevan   /* MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE); */
579c368985SVijay Mahadevan 
582494be97SVijay Mahadevan   *J = A;
59032b8ab6SVijay Mahadevan   PetscFunctionReturn(0);
60032b8ab6SVijay Mahadevan }
61032b8ab6SVijay Mahadevan 
62*9371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMMoab_Compute_NNZ_From_Connectivity(DM dm, PetscInt *innz, PetscInt *nnz, PetscInt *ionz, PetscInt *onz, PetscBool isbaij) {
63b117cd09SVijay Mahadevan   PetscInt                        i, f, nloc, vpere, bs, n_nnz, n_onz, ivtx = 0;
64c2612ac9SVijay Mahadevan   PetscInt                        ibs, jbs, inbsize, iobsize, nfields, nlsiz;
65032b8ab6SVijay Mahadevan   DM_Moab                        *dmmoab = (DM_Moab *)dm->data;
669c368985SVijay Mahadevan   moab::Range                     found;
679c368985SVijay Mahadevan   std::vector<moab::EntityHandle> adjs, storage;
685905e1eaSVijay Mahadevan   PetscBool                       isinterlaced;
69da6192ceSVijay Mahadevan   moab::EntityHandle              vtx;
70032b8ab6SVijay Mahadevan   moab::ErrorCode                 merr;
71032b8ab6SVijay Mahadevan 
72032b8ab6SVijay Mahadevan   PetscFunctionBegin;
73032b8ab6SVijay Mahadevan   bs           = dmmoab->bs;
74032b8ab6SVijay Mahadevan   nloc         = dmmoab->nloc;
755905e1eaSVijay Mahadevan   nfields      = dmmoab->numFields;
7682dfd14aSVijay Mahadevan   isinterlaced = (isbaij || bs == nfields ? PETSC_TRUE : PETSC_FALSE);
77c2612ac9SVijay Mahadevan   nlsiz        = (isinterlaced ? nloc : nloc * nfields);
78032b8ab6SVijay Mahadevan 
79f6829af0SVijay Mahadevan   /* loop over the locally owned vertices and figure out the NNZ pattern using connectivity information */
8064e1c140SVijay Mahadevan   for (moab::Range::const_iterator iter = dmmoab->vowned->begin(); iter != dmmoab->vowned->end(); iter++, ivtx++) {
81da6192ceSVijay Mahadevan     vtx = *iter;
82da6192ceSVijay Mahadevan     /* Get adjacency information for current vertex - i.e., all elements of dimension (dim) that connects
83da6192ceSVijay Mahadevan        to the current vertex. We can then decipher if a vertex is ghosted or not and compute the
84da6192ceSVijay Mahadevan        non-zero pattern accordingly. */
8564e1c140SVijay Mahadevan     adjs.clear();
86e92d1c7cSVijay Mahadevan     if (dmmoab->hlevel && (dmmoab->pcomm->size() == 1)) {
87*9371c9d4SSatish Balay       merr = dmmoab->hierarchy->get_adjacencies(vtx, dmmoab->dim, adjs);
88*9371c9d4SSatish Balay       MBERRNM(merr);
891baa6e33SBarry Smith     } else {
90*9371c9d4SSatish Balay       merr = dmmoab->mbiface->get_adjacencies(&vtx, 1, dmmoab->dim, true, adjs, moab::Interface::UNION);
91*9371c9d4SSatish Balay       MBERRNM(merr);
929c368985SVijay Mahadevan     }
930df6e276SVijay Mahadevan 
94e427d9c9SVijay Mahadevan     /* reset counters */
95e427d9c9SVijay Mahadevan     n_nnz = n_onz = 0;
96e427d9c9SVijay Mahadevan     found.clear();
97e427d9c9SVijay Mahadevan 
9815de014fSVijay Mahadevan     /* loop over vertices and update the number of connectivity */
999c368985SVijay Mahadevan     for (unsigned jter = 0; jter < adjs.size(); ++jter) {
10015de014fSVijay Mahadevan       /* Get connectivity information in canonical ordering for the local element */
1019c368985SVijay Mahadevan       const moab::EntityHandle       *connect;
1029c368985SVijay Mahadevan       std::vector<moab::EntityHandle> cconnect;
103*9371c9d4SSatish Balay       merr = dmmoab->mbiface->get_connectivity(adjs[jter], connect, vpere, false, &storage);
104*9371c9d4SSatish Balay       MBERRNM(merr);
1050df6e276SVijay Mahadevan 
106f6829af0SVijay Mahadevan       /* loop over each element connected to the adjacent vertex and update as needed */
10715de014fSVijay Mahadevan       for (i = 0; i < vpere; ++i) {
10864e1c140SVijay Mahadevan         /* find the truly user-expected layer of ghosted entities to decipher NNZ pattern */
109f6829af0SVijay Mahadevan         if (connect[i] == vtx || found.find(connect[i]) != found.end()) continue; /* make sure we don't double count shared vertices */
11064e1c140SVijay Mahadevan         if (dmmoab->vghost->find(connect[i]) != dmmoab->vghost->end()) n_onz++;   /* update out-of-proc onz */
111f6829af0SVijay Mahadevan         else n_nnz++;                                                             /* else local vertex */
11215de014fSVijay Mahadevan         found.insert(connect[i]);
11315de014fSVijay Mahadevan       }
11415de014fSVijay Mahadevan     }
1159daf19fdSVijay Mahadevan     storage.clear();
11615de014fSVijay Mahadevan 
11715de014fSVijay Mahadevan     if (isbaij) {
11815de014fSVijay Mahadevan       nnz[ivtx] = n_nnz; /* leave out self to avoid repeats -> node shared by multiple elements */
119*9371c9d4SSatish Balay       if (onz) { onz[ivtx] = n_onz; /* add ghost non-owned nodes */ }
1201baa6e33SBarry Smith     } else { /* AIJ matrices */
1215905e1eaSVijay Mahadevan       if (!isinterlaced) {
1225905e1eaSVijay Mahadevan         for (f = 0; f < nfields; f++) {
1235905e1eaSVijay Mahadevan           nnz[f * nloc + ivtx] = n_nnz;          /* leave out self to avoid repeats -> node shared by multiple elements */
124*9371c9d4SSatish Balay           if (onz) onz[f * nloc + ivtx] = n_onz; /* add ghost non-owned nodes */
1255905e1eaSVijay Mahadevan         }
1261baa6e33SBarry Smith       } else {
1275905e1eaSVijay Mahadevan         for (f = 0; f < nfields; f++) {
1285905e1eaSVijay Mahadevan           nnz[nfields * ivtx + f] = n_nnz;          /* leave out self to avoid repeats -> node shared by multiple elements */
129*9371c9d4SSatish Balay           if (onz) onz[nfields * ivtx + f] = n_onz; /* add ghost non-owned nodes */
1305905e1eaSVijay Mahadevan         }
1310df6e276SVijay Mahadevan       }
1320df6e276SVijay Mahadevan     }
1330df6e276SVijay Mahadevan   }
1340df6e276SVijay Mahadevan 
135*9371c9d4SSatish Balay   for (i = 0; i < nlsiz; i++) nnz[i] += 1; /* self count the node */
1365905e1eaSVijay Mahadevan 
13782dfd14aSVijay Mahadevan   for (ivtx = 0; ivtx < nloc; ivtx++) {
1385905e1eaSVijay Mahadevan     if (!isbaij) {
1395905e1eaSVijay Mahadevan       for (ibs = 0; ibs < nfields; ibs++) {
14082dfd14aSVijay Mahadevan         if (dmmoab->dfill) { /* first address the diagonal block */
1415905e1eaSVijay Mahadevan           /* just add up the ints -- easier/faster rather than branching based on "1" */
142*9371c9d4SSatish Balay           for (jbs = 0, inbsize = 0; jbs < nfields; jbs++) inbsize += dmmoab->dfill[ibs * nfields + jbs];
143*9371c9d4SSatish Balay         } else inbsize = nfields; /* dense coupling since user didn't specify the component fill explicitly */
14482dfd14aSVijay Mahadevan         if (isinterlaced) nnz[ivtx * nfields + ibs] *= inbsize;
14582dfd14aSVijay Mahadevan         else nnz[ibs * nloc + ivtx] *= inbsize;
1465905e1eaSVijay Mahadevan 
1475905e1eaSVijay Mahadevan         if (onz) {
14882dfd14aSVijay Mahadevan           if (dmmoab->ofill) { /* next address the off-diagonal block */
1495905e1eaSVijay Mahadevan             /* just add up the ints -- easier/faster rather than branching based on "1" */
150*9371c9d4SSatish Balay             for (jbs = 0, iobsize = 0; jbs < nfields; jbs++) iobsize += dmmoab->dfill[ibs * nfields + jbs];
151*9371c9d4SSatish Balay           } else iobsize = nfields; /* dense coupling since user didn't specify the component fill explicitly */
15282dfd14aSVijay Mahadevan           if (isinterlaced) onz[ivtx * nfields + ibs] *= iobsize;
15382dfd14aSVijay Mahadevan           else onz[ibs * nloc + ivtx] *= iobsize;
1545905e1eaSVijay Mahadevan         }
1555905e1eaSVijay Mahadevan       }
1561baa6e33SBarry Smith     } else {
15782dfd14aSVijay Mahadevan       /* check if we got overzealous in our nnz and onz computations */
15882dfd14aSVijay Mahadevan       nnz[ivtx] = (nnz[ivtx] > dmmoab->nloc ? dmmoab->nloc : nnz[ivtx]);
15982dfd14aSVijay Mahadevan       if (onz) onz[ivtx] = (onz[ivtx] > dmmoab->nloc ? dmmoab->nloc : onz[ivtx]);
1605905e1eaSVijay Mahadevan     }
1615905e1eaSVijay Mahadevan   }
1625905e1eaSVijay Mahadevan   /* update innz and ionz based on local maxima */
1635905e1eaSVijay Mahadevan   if (innz || ionz) {
1640df6e276SVijay Mahadevan     if (innz) *innz = 0;
1650df6e276SVijay Mahadevan     if (ionz) *ionz = 0;
166c2612ac9SVijay Mahadevan     for (i = 0; i < nlsiz; i++) {
167da6192ceSVijay Mahadevan       if (innz && (nnz[i] > *innz)) *innz = nnz[i];
168da6192ceSVijay Mahadevan       if ((ionz && onz) && (onz[i] > *ionz)) *ionz = onz[i];
169032b8ab6SVijay Mahadevan     }
1705905e1eaSVijay Mahadevan   }
1710df6e276SVijay Mahadevan   PetscFunctionReturn(0);
172032b8ab6SVijay Mahadevan }
173032b8ab6SVijay Mahadevan 
174*9371c9d4SSatish Balay static PetscErrorCode DMMoabSetBlockFills_Private(PetscInt w, const PetscInt *fill, PetscInt **rfill) {
1755905e1eaSVijay Mahadevan   PetscInt i, j, *ifill;
1765905e1eaSVijay Mahadevan 
1775905e1eaSVijay Mahadevan   PetscFunctionBegin;
1785905e1eaSVijay Mahadevan   if (!fill) PetscFunctionReturn(0);
1799566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(w * w, &ifill));
1805905e1eaSVijay Mahadevan   for (i = 0; i < w; i++) {
181*9371c9d4SSatish Balay     for (j = 0; j < w; j++) ifill[i * w + j] = fill[i * w + j];
1825905e1eaSVijay Mahadevan   }
1835905e1eaSVijay Mahadevan 
1845905e1eaSVijay Mahadevan   *rfill = ifill;
1855905e1eaSVijay Mahadevan   PetscFunctionReturn(0);
1865905e1eaSVijay Mahadevan }
1875905e1eaSVijay Mahadevan 
188cab5ea25SPierre Jolivet /*@C
1895905e1eaSVijay Mahadevan     DMMoabSetBlockFills - Sets the fill pattern in each block for a multi-component problem
1905905e1eaSVijay Mahadevan     of the matrix returned by DMCreateMatrix().
1915905e1eaSVijay Mahadevan 
192d083f849SBarry Smith     Logically Collective on da
1935905e1eaSVijay Mahadevan 
194d8d19677SJose E. Roman     Input Parameters:
1955905e1eaSVijay Mahadevan +   dm - the DMMoab object
1965905e1eaSVijay Mahadevan .   dfill - the fill pattern in the diagonal block (may be NULL, means use dense block)
1975905e1eaSVijay Mahadevan -   ofill - the fill pattern in the off-diagonal blocks
1985905e1eaSVijay Mahadevan 
1995905e1eaSVijay Mahadevan     Level: developer
2005905e1eaSVijay Mahadevan 
20195452b02SPatrick Sanan     Notes:
20295452b02SPatrick Sanan     This only makes sense when you are doing multicomponent problems but using the
2035905e1eaSVijay Mahadevan        MPIAIJ matrix format
2045905e1eaSVijay Mahadevan 
2055905e1eaSVijay Mahadevan            The format for dfill and ofill is a 2 dimensional dof by dof matrix with 1 entries
2065905e1eaSVijay Mahadevan        representing coupling and 0 entries for missing coupling. For example
2075905e1eaSVijay Mahadevan $             dfill[9] = {1, 0, 0,
2085905e1eaSVijay Mahadevan $                         1, 1, 0,
2095905e1eaSVijay Mahadevan $                         0, 1, 1}
2105905e1eaSVijay Mahadevan        means that row 0 is coupled with only itself in the diagonal block, row 1 is coupled with
2115905e1eaSVijay Mahadevan        itself and row 0 (in the diagonal block) and row 2 is coupled with itself and row 1 (in the
2125905e1eaSVijay Mahadevan        diagonal block).
2135905e1eaSVijay Mahadevan 
2145905e1eaSVijay Mahadevan      DMDASetGetMatrix() allows you to provide general code for those more complicated nonzero patterns then
2155905e1eaSVijay Mahadevan      can be represented in the dfill, ofill format
2165905e1eaSVijay Mahadevan 
2175905e1eaSVijay Mahadevan    Contributed by Glenn Hammond
2185905e1eaSVijay Mahadevan 
219db781477SPatrick Sanan .seealso `DMCreateMatrix()`, `DMDASetGetMatrix()`, `DMSetMatrixPreallocateOnly()`
2205905e1eaSVijay Mahadevan 
2215905e1eaSVijay Mahadevan @*/
222*9371c9d4SSatish Balay PetscErrorCode DMMoabSetBlockFills(DM dm, const PetscInt *dfill, const PetscInt *ofill) {
2235905e1eaSVijay Mahadevan   DM_Moab *dmmoab = (DM_Moab *)dm->data;
2245905e1eaSVijay Mahadevan 
2255905e1eaSVijay Mahadevan   PetscFunctionBegin;
2265905e1eaSVijay Mahadevan   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2279566063dSJacob Faibussowitsch   PetscCall(DMMoabSetBlockFills_Private(dmmoab->numFields, dfill, &dmmoab->dfill));
2289566063dSJacob Faibussowitsch   PetscCall(DMMoabSetBlockFills_Private(dmmoab->numFields, ofill, &dmmoab->ofill));
2295905e1eaSVijay Mahadevan   PetscFunctionReturn(0);
2305905e1eaSVijay Mahadevan }
231