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
DMCreateMatrix_Moab(DM dm,Mat * J)10d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreateMatrix_Moab(DM dm, Mat *J)
11d71ae5a4SJacob Faibussowitsch {
12c2612ac9SVijay Mahadevan PetscInt innz = 0, ionz = 0, nlsiz;
13032b8ab6SVijay Mahadevan DM_Moab *dmmoab = (DM_Moab *)dm->data;
14da6192ceSVijay Mahadevan PetscInt *nnz = 0, *onz = 0;
15032b8ab6SVijay Mahadevan char *tmp = 0;
162494be97SVijay Mahadevan Mat A;
17baf0d1e0SVijay Mahadevan MatType mtype;
18032b8ab6SVijay Mahadevan
19032b8ab6SVijay Mahadevan PetscFunctionBegin;
20032b8ab6SVijay Mahadevan PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
214f572ea9SToby Isaac PetscAssertPointer(J, 3);
22032b8ab6SVijay Mahadevan
23db66d124SVijay Mahadevan /* next, need to allocate the non-zero arrays to enable pre-allocation */
24baf0d1e0SVijay Mahadevan mtype = dm->mattype;
259566063dSJacob Faibussowitsch PetscCall(PetscStrstr(mtype, MATBAIJ, &tmp));
2682dfd14aSVijay Mahadevan nlsiz = (tmp ? dmmoab->nloc : dmmoab->nloc * dmmoab->numFields);
27032b8ab6SVijay Mahadevan
28032b8ab6SVijay Mahadevan /* allocate the nnz, onz arrays based on block size and local nodes */
299566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(nlsiz, &nnz, nlsiz, &onz));
30032b8ab6SVijay Mahadevan
31032b8ab6SVijay Mahadevan /* compute the nonzero pattern based on MOAB connectivity data for local elements */
3257508eceSPierre Jolivet PetscCall(DMMoab_Compute_NNZ_From_Connectivity(dm, &innz, nnz, &ionz, onz, tmp ? PETSC_TRUE : PETSC_FALSE));
33032b8ab6SVijay Mahadevan
34da6192ceSVijay Mahadevan /* create the Matrix and set its type as specified by user */
3557508eceSPierre Jolivet PetscCall(MatCreate(((PetscObject)dm)->comm, &A));
369566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, dmmoab->nloc * dmmoab->numFields, dmmoab->nloc * dmmoab->numFields, PETSC_DETERMINE, PETSC_DETERMINE));
379566063dSJacob Faibussowitsch PetscCall(MatSetType(A, mtype));
389566063dSJacob Faibussowitsch PetscCall(MatSetBlockSize(A, dmmoab->bs));
399566063dSJacob Faibussowitsch PetscCall(MatSetDM(A, dm)); /* set DM reference */
409566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A));
41da6192ceSVijay Mahadevan
4257508eceSPierre Jolivet PetscCheck(dmmoab->ltog_map, ((PetscObject)dm)->comm, PETSC_ERR_ORDER, "Cannot create a DMMoab Mat without calling DMSetUp first.");
439566063dSJacob Faibussowitsch PetscCall(MatSetLocalToGlobalMapping(A, dmmoab->ltog_map, dmmoab->ltog_map));
44032b8ab6SVijay Mahadevan
45032b8ab6SVijay Mahadevan /* set preallocation based on different supported Mat types */
469566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(A, innz, nnz));
479566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(A, innz, nnz, ionz, onz));
489566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetPreallocation(A, dmmoab->bs, innz, nnz));
499566063dSJacob Faibussowitsch PetscCall(MatMPIBAIJSetPreallocation(A, dmmoab->bs, innz, nnz, ionz, onz));
50032b8ab6SVijay Mahadevan
51da8c5b52SVijay Mahadevan /* clean up temporary memory */
529566063dSJacob Faibussowitsch PetscCall(PetscFree2(nnz, onz));
53da8c5b52SVijay Mahadevan
54da8c5b52SVijay Mahadevan /* set up internal matrix data-structures */
559566063dSJacob Faibussowitsch PetscCall(MatSetUp(A));
56da8c5b52SVijay Mahadevan
57e92d1c7cSVijay Mahadevan /* MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE); */
589c368985SVijay Mahadevan
592494be97SVijay Mahadevan *J = A;
603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
61032b8ab6SVijay Mahadevan }
62032b8ab6SVijay Mahadevan
DMMoab_Compute_NNZ_From_Connectivity(DM dm,PetscInt * innz,PetscInt * nnz,PetscInt * ionz,PetscInt * onz,PetscBool isbaij)63d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMMoab_Compute_NNZ_From_Connectivity(DM dm, PetscInt *innz, PetscInt *nnz, PetscInt *ionz, PetscInt *onz, PetscBool isbaij)
64d71ae5a4SJacob Faibussowitsch {
65b117cd09SVijay 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;
689c368985SVijay Mahadevan moab::Range found;
699c368985SVijay Mahadevan std::vector<moab::EntityHandle> adjs, storage;
705905e1eaSVijay Mahadevan PetscBool isinterlaced;
71da6192ceSVijay Mahadevan moab::EntityHandle vtx;
72032b8ab6SVijay Mahadevan moab::ErrorCode merr;
73032b8ab6SVijay Mahadevan
74032b8ab6SVijay Mahadevan PetscFunctionBegin;
75032b8ab6SVijay Mahadevan bs = dmmoab->bs;
76032b8ab6SVijay Mahadevan nloc = dmmoab->nloc;
775905e1eaSVijay Mahadevan nfields = dmmoab->numFields;
7882dfd14aSVijay Mahadevan isinterlaced = (isbaij || bs == nfields ? PETSC_TRUE : PETSC_FALSE);
79c2612ac9SVijay Mahadevan nlsiz = (isinterlaced ? nloc : nloc * nfields);
80032b8ab6SVijay Mahadevan
81f6829af0SVijay Mahadevan /* loop over the locally owned vertices and figure out the NNZ pattern using connectivity information */
8264e1c140SVijay Mahadevan for (moab::Range::const_iterator iter = dmmoab->vowned->begin(); iter != dmmoab->vowned->end(); iter++, ivtx++) {
83da6192ceSVijay Mahadevan vtx = *iter;
84da6192ceSVijay Mahadevan /* Get adjacency information for current vertex - i.e., all elements of dimension (dim) that connects
85da6192ceSVijay Mahadevan to the current vertex. We can then decipher if a vertex is ghosted or not and compute the
86da6192ceSVijay Mahadevan non-zero pattern accordingly. */
8764e1c140SVijay Mahadevan adjs.clear();
88e92d1c7cSVijay Mahadevan if (dmmoab->hlevel && (dmmoab->pcomm->size() == 1)) {
899371c9d4SSatish Balay merr = dmmoab->hierarchy->get_adjacencies(vtx, dmmoab->dim, adjs);
909371c9d4SSatish Balay MBERRNM(merr);
911baa6e33SBarry Smith } else {
929371c9d4SSatish Balay merr = dmmoab->mbiface->get_adjacencies(&vtx, 1, dmmoab->dim, true, adjs, moab::Interface::UNION);
939371c9d4SSatish Balay MBERRNM(merr);
949c368985SVijay Mahadevan }
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 */
1019c368985SVijay Mahadevan for (unsigned jter = 0; jter < adjs.size(); ++jter) {
10215de014fSVijay Mahadevan /* Get connectivity information in canonical ordering for the local element */
1039c368985SVijay Mahadevan const moab::EntityHandle *connect;
1049c368985SVijay Mahadevan std::vector<moab::EntityHandle> cconnect;
1059371c9d4SSatish Balay merr = dmmoab->mbiface->get_connectivity(adjs[jter], connect, vpere, false, &storage);
1069371c9d4SSatish Balay 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 */
121*ac530a7eSPierre Jolivet if (onz) onz[ivtx] = n_onz; /* add ghost non-owned nodes */
1221baa6e33SBarry Smith } else { /* AIJ matrices */
1235905e1eaSVijay Mahadevan if (!isinterlaced) {
1245905e1eaSVijay Mahadevan for (f = 0; f < nfields; f++) {
1255905e1eaSVijay Mahadevan nnz[f * nloc + ivtx] = n_nnz; /* leave out self to avoid repeats -> node shared by multiple elements */
1269371c9d4SSatish Balay if (onz) onz[f * nloc + ivtx] = n_onz; /* add ghost non-owned nodes */
1275905e1eaSVijay Mahadevan }
1281baa6e33SBarry Smith } else {
1295905e1eaSVijay Mahadevan for (f = 0; f < nfields; f++) {
1305905e1eaSVijay Mahadevan nnz[nfields * ivtx + f] = n_nnz; /* leave out self to avoid repeats -> node shared by multiple elements */
1319371c9d4SSatish Balay if (onz) onz[nfields * ivtx + f] = n_onz; /* add ghost non-owned nodes */
1325905e1eaSVijay Mahadevan }
1330df6e276SVijay Mahadevan }
1340df6e276SVijay Mahadevan }
1350df6e276SVijay Mahadevan }
1360df6e276SVijay Mahadevan
1379371c9d4SSatish Balay for (i = 0; i < nlsiz; i++) nnz[i] += 1; /* self count the node */
1385905e1eaSVijay Mahadevan
13982dfd14aSVijay Mahadevan for (ivtx = 0; ivtx < nloc; ivtx++) {
1405905e1eaSVijay Mahadevan if (!isbaij) {
1415905e1eaSVijay Mahadevan for (ibs = 0; ibs < nfields; ibs++) {
14282dfd14aSVijay Mahadevan if (dmmoab->dfill) { /* first address the diagonal block */
1435905e1eaSVijay Mahadevan /* just add up the ints -- easier/faster rather than branching based on "1" */
1449371c9d4SSatish Balay for (jbs = 0, inbsize = 0; jbs < nfields; jbs++) inbsize += dmmoab->dfill[ibs * nfields + jbs];
1459371c9d4SSatish Balay } else inbsize = nfields; /* dense coupling since user didn't specify the component fill explicitly */
14682dfd14aSVijay Mahadevan if (isinterlaced) nnz[ivtx * nfields + ibs] *= inbsize;
14782dfd14aSVijay Mahadevan else nnz[ibs * nloc + ivtx] *= inbsize;
1485905e1eaSVijay Mahadevan
1495905e1eaSVijay Mahadevan if (onz) {
15082dfd14aSVijay Mahadevan if (dmmoab->ofill) { /* next address the off-diagonal block */
1515905e1eaSVijay Mahadevan /* just add up the ints -- easier/faster rather than branching based on "1" */
1529371c9d4SSatish Balay for (jbs = 0, iobsize = 0; jbs < nfields; jbs++) iobsize += dmmoab->dfill[ibs * nfields + jbs];
1539371c9d4SSatish Balay } else iobsize = nfields; /* dense coupling since user didn't specify the component fill explicitly */
15482dfd14aSVijay Mahadevan if (isinterlaced) onz[ivtx * nfields + ibs] *= iobsize;
15582dfd14aSVijay Mahadevan else onz[ibs * nloc + ivtx] *= iobsize;
1565905e1eaSVijay Mahadevan }
1575905e1eaSVijay Mahadevan }
1581baa6e33SBarry Smith } else {
15982dfd14aSVijay Mahadevan /* check if we got overzealous in our nnz and onz computations */
16082dfd14aSVijay Mahadevan nnz[ivtx] = (nnz[ivtx] > dmmoab->nloc ? dmmoab->nloc : nnz[ivtx]);
16182dfd14aSVijay Mahadevan if (onz) onz[ivtx] = (onz[ivtx] > dmmoab->nloc ? dmmoab->nloc : onz[ivtx]);
1625905e1eaSVijay Mahadevan }
1635905e1eaSVijay Mahadevan }
1645905e1eaSVijay Mahadevan /* update innz and ionz based on local maxima */
1655905e1eaSVijay Mahadevan if (innz || ionz) {
1660df6e276SVijay Mahadevan if (innz) *innz = 0;
1670df6e276SVijay Mahadevan if (ionz) *ionz = 0;
168c2612ac9SVijay Mahadevan for (i = 0; i < nlsiz; i++) {
169b6555650SPierre Jolivet if (innz && nnz[i] > *innz) *innz = nnz[i];
170b6555650SPierre Jolivet if (ionz && onz && onz[i] > *ionz) *ionz = onz[i];
171032b8ab6SVijay Mahadevan }
1725905e1eaSVijay Mahadevan }
1733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
174032b8ab6SVijay Mahadevan }
175032b8ab6SVijay Mahadevan
DMMoabSetBlockFills_Private(PetscInt w,const PetscInt * fill,PetscInt ** rfill)176d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMMoabSetBlockFills_Private(PetscInt w, const PetscInt *fill, PetscInt **rfill)
177d71ae5a4SJacob Faibussowitsch {
1785905e1eaSVijay Mahadevan PetscInt i, j, *ifill;
1795905e1eaSVijay Mahadevan
1805905e1eaSVijay Mahadevan PetscFunctionBegin;
1813ba16761SJacob Faibussowitsch if (!fill) PetscFunctionReturn(PETSC_SUCCESS);
1829566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(w * w, &ifill));
1835905e1eaSVijay Mahadevan for (i = 0; i < w; i++) {
1849371c9d4SSatish Balay for (j = 0; j < w; j++) ifill[i * w + j] = fill[i * w + j];
1855905e1eaSVijay Mahadevan }
1865905e1eaSVijay Mahadevan
1875905e1eaSVijay Mahadevan *rfill = ifill;
1883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1895905e1eaSVijay Mahadevan }
1905905e1eaSVijay Mahadevan
191cab5ea25SPierre Jolivet /*@C
1925905e1eaSVijay Mahadevan DMMoabSetBlockFills - Sets the fill pattern in each block for a multi-component problem
19320f4b53cSBarry Smith of the matrix returned by `DMCreateMatrix()`.
1945905e1eaSVijay Mahadevan
19520f4b53cSBarry Smith Logically Collective
1965905e1eaSVijay Mahadevan
197d8d19677SJose E. Roman Input Parameters:
19820f4b53cSBarry Smith + dm - the `DMMOAB` object
19920f4b53cSBarry Smith . dfill - the fill pattern in the diagonal block (may be `NULL`, means use dense block)
2005905e1eaSVijay Mahadevan - ofill - the fill pattern in the off-diagonal blocks
2015905e1eaSVijay Mahadevan
2025905e1eaSVijay Mahadevan Level: developer
2035905e1eaSVijay Mahadevan
20495452b02SPatrick Sanan Notes:
20595452b02SPatrick Sanan This only makes sense when you are doing multicomponent problems but using the
2065905e1eaSVijay Mahadevan MPIAIJ matrix format
2075905e1eaSVijay Mahadevan
2085905e1eaSVijay Mahadevan The format for dfill and ofill is a 2 dimensional dof by dof matrix with 1 entries
2095905e1eaSVijay Mahadevan representing coupling and 0 entries for missing coupling. For example
21020f4b53cSBarry Smith .vb
21120f4b53cSBarry Smith dfill[9] = {1, 0, 0,
21220f4b53cSBarry Smith 1, 1, 0,
21320f4b53cSBarry Smith 0, 1, 1}
21420f4b53cSBarry Smith .ve
2155905e1eaSVijay Mahadevan means that row 0 is coupled with only itself in the diagonal block, row 1 is coupled with
2165905e1eaSVijay Mahadevan itself and row 0 (in the diagonal block) and row 2 is coupled with itself and row 1 (in the
2175905e1eaSVijay Mahadevan diagonal block).
2185905e1eaSVijay Mahadevan
21920f4b53cSBarry Smith `DMDASetGetMatrix()` allows you to provide general code for those more complicated nonzero patterns then
2205905e1eaSVijay Mahadevan can be represented in the dfill, ofill format
2215905e1eaSVijay Mahadevan
2221d27aa22SBarry Smith Contributed by\:
2231d27aa22SBarry Smith Glenn Hammond
2245905e1eaSVijay Mahadevan
225a4e35b19SJacob Faibussowitsch .seealso: `DMCreateMatrix()`, `DMDASetGetMatrix()`, `DMSetMatrixPreallocateOnly()`
2265905e1eaSVijay Mahadevan @*/
DMMoabSetBlockFills(DM dm,const PetscInt * dfill,const PetscInt * ofill)227d71ae5a4SJacob Faibussowitsch PetscErrorCode DMMoabSetBlockFills(DM dm, const PetscInt *dfill, const PetscInt *ofill)
228d71ae5a4SJacob Faibussowitsch {
2295905e1eaSVijay Mahadevan DM_Moab *dmmoab = (DM_Moab *)dm->data;
2305905e1eaSVijay Mahadevan
2315905e1eaSVijay Mahadevan PetscFunctionBegin;
2325905e1eaSVijay Mahadevan PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2339566063dSJacob Faibussowitsch PetscCall(DMMoabSetBlockFills_Private(dmmoab->numFields, dfill, &dmmoab->dfill));
2349566063dSJacob Faibussowitsch PetscCall(DMMoabSetBlockFills_Private(dmmoab->numFields, ofill, &dmmoab->ofill));
2353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2365905e1eaSVijay Mahadevan }
237