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 10304006b3SVijay Mahadevan PETSC_EXTERN 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; 2652c5f739Sprj- ierr = PetscStrstr(mtype, MATBAIJ, &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 */ 369daf19fdSVijay Mahadevan ierr = MatCreate((((PetscObject)dm)->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 439daf19fdSVijay Mahadevan if (!dmmoab->ltog_map) SETERRQ((((PetscObject)dm)->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 58e92d1c7cSVijay Mahadevan /* MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE); */ 599c368985SVijay Mahadevan 602494be97SVijay Mahadevan *J = A; 61032b8ab6SVijay Mahadevan PetscFunctionReturn(0); 62032b8ab6SVijay Mahadevan } 63032b8ab6SVijay Mahadevan 64032b8ab6SVijay Mahadevan 65a86ed394SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoab_Compute_NNZ_From_Connectivity(DM dm, PetscInt* innz, PetscInt* nnz, PetscInt* ionz, PetscInt* onz, PetscBool isbaij) 66032b8ab6SVijay Mahadevan { 67b117cd09SVijay Mahadevan PetscInt i, f, nloc, vpere, bs, n_nnz, n_onz, ivtx = 0; 68c2612ac9SVijay Mahadevan PetscInt ibs, jbs, inbsize, iobsize, nfields, nlsiz; 69032b8ab6SVijay Mahadevan DM_Moab *dmmoab = (DM_Moab*)dm->data; 709c368985SVijay Mahadevan moab::Range found; 719c368985SVijay Mahadevan std::vector<moab::EntityHandle> adjs, storage; 725905e1eaSVijay Mahadevan PetscBool isinterlaced; 73da6192ceSVijay Mahadevan moab::EntityHandle vtx; 74032b8ab6SVijay Mahadevan moab::ErrorCode merr; 75032b8ab6SVijay Mahadevan 76032b8ab6SVijay Mahadevan PetscFunctionBegin; 77032b8ab6SVijay Mahadevan bs = dmmoab->bs; 78032b8ab6SVijay Mahadevan nloc = dmmoab->nloc; 795905e1eaSVijay Mahadevan nfields = dmmoab->numFields; 8082dfd14aSVijay Mahadevan isinterlaced = (isbaij || bs == nfields ? PETSC_TRUE : PETSC_FALSE); 81c2612ac9SVijay Mahadevan nlsiz = (isinterlaced ? nloc : nloc * nfields); 82032b8ab6SVijay Mahadevan 83f6829af0SVijay Mahadevan /* loop over the locally owned vertices and figure out the NNZ pattern using connectivity information */ 8464e1c140SVijay Mahadevan for (moab::Range::const_iterator iter = dmmoab->vowned->begin(); iter != dmmoab->vowned->end(); iter++, ivtx++) { 85032b8ab6SVijay Mahadevan 86da6192ceSVijay Mahadevan vtx = *iter; 87da6192ceSVijay Mahadevan /* Get adjacency information for current vertex - i.e., all elements of dimension (dim) that connects 88da6192ceSVijay Mahadevan to the current vertex. We can then decipher if a vertex is ghosted or not and compute the 89da6192ceSVijay Mahadevan non-zero pattern accordingly. */ 9064e1c140SVijay Mahadevan adjs.clear(); 91e92d1c7cSVijay Mahadevan if (dmmoab->hlevel && (dmmoab->pcomm->size() == 1)) { 929c368985SVijay Mahadevan merr = dmmoab->hierarchy->get_adjacencies(vtx, dmmoab->dim, adjs); MBERRNM(merr); 939c368985SVijay Mahadevan } 949c368985SVijay Mahadevan else { 95e92d1c7cSVijay Mahadevan merr = dmmoab->mbiface->get_adjacencies(&vtx, 1, dmmoab->dim, true, adjs, moab::Interface::UNION); MBERRNM(merr); 969c368985SVijay Mahadevan } 970df6e276SVijay Mahadevan 98e427d9c9SVijay Mahadevan /* reset counters */ 99e427d9c9SVijay Mahadevan n_nnz = n_onz = 0; 100e427d9c9SVijay Mahadevan found.clear(); 101e427d9c9SVijay Mahadevan 10215de014fSVijay Mahadevan /* loop over vertices and update the number of connectivity */ 1039c368985SVijay Mahadevan for (unsigned jter = 0; jter < adjs.size(); ++jter) { 1040df6e276SVijay Mahadevan 10515de014fSVijay Mahadevan /* Get connectivity information in canonical ordering for the local element */ 1069c368985SVijay Mahadevan const moab::EntityHandle *connect; 1079c368985SVijay Mahadevan std::vector<moab::EntityHandle> cconnect; 1089c368985SVijay Mahadevan merr = dmmoab->mbiface->get_connectivity(adjs[jter], connect, vpere, false, &storage); MBERRNM(merr); 1090df6e276SVijay Mahadevan 110f6829af0SVijay Mahadevan /* loop over each element connected to the adjacent vertex and update as needed */ 11115de014fSVijay Mahadevan for (i = 0; i < vpere; ++i) { 11264e1c140SVijay Mahadevan /* find the truly user-expected layer of ghosted entities to decipher NNZ pattern */ 113f6829af0SVijay Mahadevan if (connect[i] == vtx || found.find(connect[i]) != found.end()) continue; /* make sure we don't double count shared vertices */ 11464e1c140SVijay Mahadevan if (dmmoab->vghost->find(connect[i]) != dmmoab->vghost->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 } 1199daf19fdSVijay Mahadevan storage.clear(); 12015de014fSVijay Mahadevan 12115de014fSVijay Mahadevan if (isbaij) { 12215de014fSVijay Mahadevan nnz[ivtx] = n_nnz; /* leave out self to avoid repeats -> node shared by multiple elements */ 1235905e1eaSVijay Mahadevan if (onz) { 1245905e1eaSVijay Mahadevan onz[ivtx] = n_onz; /* add ghost non-owned nodes */ 1255905e1eaSVijay Mahadevan } 12615de014fSVijay Mahadevan } 127f6829af0SVijay Mahadevan else { /* AIJ matrices */ 1285905e1eaSVijay Mahadevan if (!isinterlaced) { 1295905e1eaSVijay Mahadevan for (f = 0; f < nfields; f++) { 1305905e1eaSVijay Mahadevan nnz[f * nloc + ivtx] = n_nnz; /* leave out self to avoid repeats -> node shared by multiple elements */ 1315905e1eaSVijay Mahadevan if (onz) 1325905e1eaSVijay Mahadevan onz[f * nloc + ivtx] = n_onz; /* add ghost non-owned nodes */ 1335905e1eaSVijay Mahadevan } 1345905e1eaSVijay Mahadevan } 1355905e1eaSVijay Mahadevan else { 1365905e1eaSVijay Mahadevan for (f = 0; f < nfields; f++) { 1375905e1eaSVijay Mahadevan nnz[nfields * ivtx + f] = n_nnz; /* leave out self to avoid repeats -> node shared by multiple elements */ 1385905e1eaSVijay Mahadevan if (onz) 1395905e1eaSVijay Mahadevan onz[nfields * ivtx + f] = n_onz; /* add ghost non-owned nodes */ 1405905e1eaSVijay Mahadevan } 1410df6e276SVijay Mahadevan } 1420df6e276SVijay Mahadevan } 1430df6e276SVijay Mahadevan } 1440df6e276SVijay Mahadevan 145c2612ac9SVijay Mahadevan for (i = 0; i < nlsiz; i++) 1465905e1eaSVijay Mahadevan nnz[i] += 1; /* self count the node */ 1475905e1eaSVijay Mahadevan 14882dfd14aSVijay Mahadevan for (ivtx = 0; ivtx < nloc; ivtx++) { 1495905e1eaSVijay Mahadevan if (!isbaij) { 1505905e1eaSVijay Mahadevan for (ibs = 0; ibs < nfields; ibs++) { 15182dfd14aSVijay Mahadevan if (dmmoab->dfill) { /* first address the diagonal block */ 1525905e1eaSVijay Mahadevan /* just add up the ints -- easier/faster rather than branching based on "1" */ 1535905e1eaSVijay Mahadevan for (jbs = 0, inbsize = 0; jbs < nfields; jbs++) 1545905e1eaSVijay Mahadevan inbsize += dmmoab->dfill[ibs * nfields + jbs]; 1555905e1eaSVijay Mahadevan } 15682dfd14aSVijay Mahadevan else inbsize = nfields; /* dense coupling since user didn't specify the component fill explicitly */ 15782dfd14aSVijay Mahadevan if (isinterlaced) nnz[ivtx * nfields + ibs] *= inbsize; 15882dfd14aSVijay Mahadevan else nnz[ibs * nloc + ivtx] *= inbsize; 1595905e1eaSVijay Mahadevan 1605905e1eaSVijay Mahadevan if (onz) { 16182dfd14aSVijay Mahadevan if (dmmoab->ofill) { /* next address the off-diagonal block */ 1625905e1eaSVijay Mahadevan /* just add up the ints -- easier/faster rather than branching based on "1" */ 1635905e1eaSVijay Mahadevan for (jbs = 0, iobsize = 0; jbs < nfields; jbs++) 1645905e1eaSVijay Mahadevan iobsize += dmmoab->dfill[ibs * nfields + jbs]; 1655905e1eaSVijay Mahadevan } 16682dfd14aSVijay Mahadevan else iobsize = nfields; /* dense coupling since user didn't specify the component fill explicitly */ 16782dfd14aSVijay Mahadevan if (isinterlaced) onz[ivtx * nfields + ibs] *= iobsize; 16882dfd14aSVijay Mahadevan else onz[ibs * nloc + ivtx] *= iobsize; 1695905e1eaSVijay Mahadevan } 1705905e1eaSVijay Mahadevan } 1715905e1eaSVijay Mahadevan } 1725905e1eaSVijay Mahadevan else { 17382dfd14aSVijay Mahadevan /* check if we got overzealous in our nnz and onz computations */ 17482dfd14aSVijay Mahadevan nnz[ivtx] = (nnz[ivtx] > dmmoab->nloc ? dmmoab->nloc : nnz[ivtx]); 17582dfd14aSVijay Mahadevan if (onz) onz[ivtx] = (onz[ivtx] > dmmoab->nloc ? dmmoab->nloc : onz[ivtx]); 1765905e1eaSVijay Mahadevan } 1775905e1eaSVijay Mahadevan } 1785905e1eaSVijay Mahadevan /* update innz and ionz based on local maxima */ 1795905e1eaSVijay Mahadevan if (innz || ionz) { 1800df6e276SVijay Mahadevan if (innz) *innz = 0; 1810df6e276SVijay Mahadevan if (ionz) *ionz = 0; 182c2612ac9SVijay Mahadevan for (i = 0; i < nlsiz; i++) { 183da6192ceSVijay Mahadevan if (innz && (nnz[i] > *innz)) *innz = nnz[i]; 184da6192ceSVijay Mahadevan if ((ionz && onz) && (onz[i] > *ionz)) *ionz = onz[i]; 185032b8ab6SVijay Mahadevan } 1865905e1eaSVijay Mahadevan } 1870df6e276SVijay Mahadevan PetscFunctionReturn(0); 188032b8ab6SVijay Mahadevan } 189032b8ab6SVijay Mahadevan 190da8c5b52SVijay Mahadevan 1915905e1eaSVijay Mahadevan static PetscErrorCode DMMoabSetBlockFills_Private(PetscInt w, const PetscInt *fill, PetscInt **rfill) 1925905e1eaSVijay Mahadevan { 1935905e1eaSVijay Mahadevan PetscErrorCode ierr; 1945905e1eaSVijay Mahadevan PetscInt i, j, *ifill; 1955905e1eaSVijay Mahadevan 1965905e1eaSVijay Mahadevan PetscFunctionBegin; 1975905e1eaSVijay Mahadevan if (!fill) PetscFunctionReturn(0); 1985905e1eaSVijay Mahadevan ierr = PetscMalloc1(w * w, &ifill);CHKERRQ(ierr); 1995905e1eaSVijay Mahadevan for (i = 0; i < w; i++) { 2005905e1eaSVijay Mahadevan for (j = 0; j < w; j++) 2015905e1eaSVijay Mahadevan ifill[i * w + j] = fill[i * w + j]; 2025905e1eaSVijay Mahadevan } 2035905e1eaSVijay Mahadevan 2045905e1eaSVijay Mahadevan *rfill = ifill; 2055905e1eaSVijay Mahadevan PetscFunctionReturn(0); 2065905e1eaSVijay Mahadevan } 2075905e1eaSVijay Mahadevan 2085905e1eaSVijay Mahadevan 209*cab5ea25SPierre Jolivet /*@C 2105905e1eaSVijay Mahadevan DMMoabSetBlockFills - Sets the fill pattern in each block for a multi-component problem 2115905e1eaSVijay Mahadevan of the matrix returned by DMCreateMatrix(). 2125905e1eaSVijay Mahadevan 213d083f849SBarry Smith Logically Collective on da 2145905e1eaSVijay Mahadevan 2155905e1eaSVijay Mahadevan Input Parameter: 2165905e1eaSVijay Mahadevan + dm - the DMMoab object 2175905e1eaSVijay Mahadevan . dfill - the fill pattern in the diagonal block (may be NULL, means use dense block) 2185905e1eaSVijay Mahadevan - ofill - the fill pattern in the off-diagonal blocks 2195905e1eaSVijay Mahadevan 2205905e1eaSVijay Mahadevan Level: developer 2215905e1eaSVijay Mahadevan 22295452b02SPatrick Sanan Notes: 22395452b02SPatrick Sanan This only makes sense when you are doing multicomponent problems but using the 2245905e1eaSVijay Mahadevan MPIAIJ matrix format 2255905e1eaSVijay Mahadevan 2265905e1eaSVijay Mahadevan The format for dfill and ofill is a 2 dimensional dof by dof matrix with 1 entries 2275905e1eaSVijay Mahadevan representing coupling and 0 entries for missing coupling. For example 2285905e1eaSVijay Mahadevan $ dfill[9] = {1, 0, 0, 2295905e1eaSVijay Mahadevan $ 1, 1, 0, 2305905e1eaSVijay Mahadevan $ 0, 1, 1} 2315905e1eaSVijay Mahadevan means that row 0 is coupled with only itself in the diagonal block, row 1 is coupled with 2325905e1eaSVijay Mahadevan itself and row 0 (in the diagonal block) and row 2 is coupled with itself and row 1 (in the 2335905e1eaSVijay Mahadevan diagonal block). 2345905e1eaSVijay Mahadevan 2355905e1eaSVijay Mahadevan DMDASetGetMatrix() allows you to provide general code for those more complicated nonzero patterns then 2365905e1eaSVijay Mahadevan can be represented in the dfill, ofill format 2375905e1eaSVijay Mahadevan 2385905e1eaSVijay Mahadevan Contributed by Glenn Hammond 2395905e1eaSVijay Mahadevan 2405905e1eaSVijay Mahadevan .seealso DMCreateMatrix(), DMDASetGetMatrix(), DMSetMatrixPreallocateOnly() 2415905e1eaSVijay Mahadevan 2425905e1eaSVijay Mahadevan @*/ 2435905e1eaSVijay Mahadevan PetscErrorCode DMMoabSetBlockFills(DM dm, const PetscInt *dfill, const PetscInt *ofill) 2445905e1eaSVijay Mahadevan { 2455905e1eaSVijay Mahadevan DM_Moab *dmmoab = (DM_Moab*)dm->data; 2465905e1eaSVijay Mahadevan PetscErrorCode ierr; 2475905e1eaSVijay Mahadevan 2485905e1eaSVijay Mahadevan PetscFunctionBegin; 2495905e1eaSVijay Mahadevan PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2505905e1eaSVijay Mahadevan ierr = DMMoabSetBlockFills_Private(dmmoab->numFields, dfill, &dmmoab->dfill);CHKERRQ(ierr); 2515905e1eaSVijay Mahadevan ierr = DMMoabSetBlockFills_Private(dmmoab->numFields, ofill, &dmmoab->ofill);CHKERRQ(ierr); 2525905e1eaSVijay Mahadevan PetscFunctionReturn(0); 2535905e1eaSVijay Mahadevan } 254