1755f3dfbSVijay Mahadevan #include <petsc/private/dmmbimpl.h> 2b117cd09SVijay Mahadevan #include <petscdmmoab.h> 3b117cd09SVijay Mahadevan #include <MBTagConventions.hpp> 4b117cd09SVijay Mahadevan #include <moab/NestedRefine.hpp> 5b117cd09SVijay Mahadevan 6cab5ea25SPierre Jolivet /*@C 7b117cd09SVijay Mahadevan DMMoabGenerateHierarchy - Generate a multi-level uniform refinement hierarchy 8*20f4b53cSBarry Smith by succesively refining a coarse mesh, already defined in the `DM` object 9b117cd09SVijay Mahadevan provided by the user. 10b117cd09SVijay Mahadevan 11d083f849SBarry Smith Collective 12b117cd09SVijay Mahadevan 13b117cd09SVijay Mahadevan Input Parameter: 14*20f4b53cSBarry Smith . dmb - The `DMMOAB` object 15b117cd09SVijay Mahadevan 16d8d19677SJose E. Roman Output Parameters: 17b117cd09SVijay Mahadevan + nlevels - The number of levels of refinement needed to generate the hierarchy 18a2b725a8SWilliam Gropp - ldegrees - The degree of refinement at each level in the hierarchy 19b117cd09SVijay Mahadevan 20b117cd09SVijay Mahadevan Level: beginner 21b117cd09SVijay Mahadevan 22b117cd09SVijay Mahadevan @*/ 23d71ae5a4SJacob Faibussowitsch PetscErrorCode DMMoabGenerateHierarchy(DM dm, PetscInt nlevels, PetscInt *ldegrees) 24d71ae5a4SJacob Faibussowitsch { 25b117cd09SVijay Mahadevan DM_Moab *dmmoab; 26b117cd09SVijay Mahadevan moab::ErrorCode merr; 272417220eSVijay Mahadevan PetscInt *pdegrees, ilevel; 28e882eb38SVijay Mahadevan std::vector<moab::EntityHandle> hsets; 29b117cd09SVijay Mahadevan 30b117cd09SVijay Mahadevan PetscFunctionBegin; 31b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 32b117cd09SVijay Mahadevan dmmoab = (DM_Moab *)(dm)->data; 33b117cd09SVijay Mahadevan 34b117cd09SVijay Mahadevan if (!ldegrees) { 359566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nlevels, &pdegrees)); 362417220eSVijay Mahadevan for (ilevel = 0; ilevel < nlevels; ilevel++) pdegrees[ilevel] = 2; /* default = Degree 2 refinement */ 379371c9d4SSatish Balay } else pdegrees = ldegrees; 38b117cd09SVijay Mahadevan 39b117cd09SVijay Mahadevan /* initialize set level refinement data for hierarchy */ 40b117cd09SVijay Mahadevan dmmoab->nhlevels = nlevels; 41b117cd09SVijay Mahadevan 42b117cd09SVijay Mahadevan /* Instantiate the nested refinement class */ 439daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 443f1c6e43SVijay Mahadevan dmmoab->hierarchy = new moab::NestedRefine(dynamic_cast<moab::Core *>(dmmoab->mbiface), dmmoab->pcomm, dmmoab->fileset); 459daf19fdSVijay Mahadevan #else 469daf19fdSVijay Mahadevan dmmoab->hierarchy = new moab::NestedRefine(dynamic_cast<moab::Core *>(dmmoab->mbiface), NULL, dmmoab->fileset); 479daf19fdSVijay Mahadevan #endif 48b117cd09SVijay Mahadevan 499566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nlevels + 1, &dmmoab->hsets)); 50e882eb38SVijay Mahadevan 51b117cd09SVijay Mahadevan /* generate the mesh hierarchy */ 529371c9d4SSatish Balay merr = dmmoab->hierarchy->generate_mesh_hierarchy(nlevels, pdegrees, hsets, false); 539371c9d4SSatish Balay MBERRNM(merr); 54e882eb38SVijay Mahadevan 559daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 56755f3dfbSVijay Mahadevan if (dmmoab->pcomm->size() > 1) { 579371c9d4SSatish Balay merr = dmmoab->hierarchy->exchange_ghosts(hsets, dmmoab->nghostrings); 589371c9d4SSatish Balay MBERRNM(merr); 59755f3dfbSVijay Mahadevan } 609daf19fdSVijay Mahadevan #endif 6149d66b22SVijay Mahadevan 6264e1c140SVijay Mahadevan /* copy the mesh sets for nested refinement hierarchy */ 63e92d1c7cSVijay Mahadevan dmmoab->hsets[0] = hsets[0]; 649371c9d4SSatish Balay for (ilevel = 1; ilevel <= nlevels; ilevel++) { 652417220eSVijay Mahadevan dmmoab->hsets[ilevel] = hsets[ilevel]; 66e92d1c7cSVijay Mahadevan 679c368985SVijay Mahadevan #ifdef MOAB_HAVE_MPI 689371c9d4SSatish Balay merr = dmmoab->pcomm->assign_global_ids(hsets[ilevel], dmmoab->dim, 0, false, true, false); 699371c9d4SSatish Balay MBERRNM(merr); 709c368985SVijay Mahadevan #endif 71e92d1c7cSVijay Mahadevan 72e92d1c7cSVijay Mahadevan /* Update material and other geometric tags from parent to child sets */ 739371c9d4SSatish Balay merr = dmmoab->hierarchy->update_special_tags(ilevel, hsets[ilevel]); 749371c9d4SSatish Balay MBERRNM(merr); 752417220eSVijay Mahadevan } 7664e1c140SVijay Mahadevan 7764e1c140SVijay Mahadevan hsets.clear(); 7848a46eb9SPierre Jolivet if (!ldegrees) PetscCall(PetscFree(pdegrees)); 793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 80b117cd09SVijay Mahadevan } 81b117cd09SVijay Mahadevan 82cab5ea25SPierre Jolivet /*@C 83*20f4b53cSBarry Smith DMRefineHierarchy_Moab - Generate a multi-level `DM` hierarchy 84e882eb38SVijay Mahadevan by succesively refining a coarse mesh. 85b117cd09SVijay Mahadevan 86d083f849SBarry Smith Collective 87b117cd09SVijay Mahadevan 88b117cd09SVijay Mahadevan Input Parameter: 89*20f4b53cSBarry Smith . dm - The `DMMOAB` object 90b117cd09SVijay Mahadevan 91d8d19677SJose E. Roman Output Parameters: 92e882eb38SVijay Mahadevan + nlevels - The number of levels of refinement needed to generate the hierarchy 93a2b725a8SWilliam Gropp - dmf - The DM objects after successive refinement of the hierarchy 94b117cd09SVijay Mahadevan 95b117cd09SVijay Mahadevan Level: beginner 96b117cd09SVijay Mahadevan 97b117cd09SVijay Mahadevan @*/ 98d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMRefineHierarchy_Moab(DM dm, PetscInt nlevels, DM dmf[]) 99d71ae5a4SJacob Faibussowitsch { 100b117cd09SVijay Mahadevan PetscInt i; 101b117cd09SVijay Mahadevan 102b117cd09SVijay Mahadevan PetscFunctionBegin; 103b117cd09SVijay Mahadevan 1049566063dSJacob Faibussowitsch PetscCall(DMRefine(dm, PetscObjectComm((PetscObject)dm), &dmf[0])); 10548a46eb9SPierre Jolivet for (i = 1; i < nlevels; i++) PetscCall(DMRefine(dmf[i - 1], PetscObjectComm((PetscObject)dm), &dmf[i])); 1063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 107b117cd09SVijay Mahadevan } 108b117cd09SVijay Mahadevan 109cab5ea25SPierre Jolivet /*@C 110*20f4b53cSBarry Smith DMCoarsenHierarchy_Moab - Generate a multi-level `DM` hierarchy 111e882eb38SVijay Mahadevan by succesively coarsening a refined mesh. 112b117cd09SVijay Mahadevan 113d083f849SBarry Smith Collective 114b117cd09SVijay Mahadevan 115b117cd09SVijay Mahadevan Input Parameter: 116*20f4b53cSBarry Smith . dm - The `DMMOAB` object 117b117cd09SVijay Mahadevan 118d8d19677SJose E. Roman Output Parameters: 119e882eb38SVijay Mahadevan + nlevels - The number of levels of refinement needed to generate the hierarchy 120*20f4b53cSBarry Smith - dmc - The `DM` objects after successive coarsening of the hierarchy 121b117cd09SVijay Mahadevan 122b117cd09SVijay Mahadevan Level: beginner 123b117cd09SVijay Mahadevan 124b117cd09SVijay Mahadevan @*/ 125d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCoarsenHierarchy_Moab(DM dm, PetscInt nlevels, DM dmc[]) 126d71ae5a4SJacob Faibussowitsch { 127b117cd09SVijay Mahadevan PetscInt i; 128b117cd09SVijay Mahadevan 129b117cd09SVijay Mahadevan PetscFunctionBegin; 130b117cd09SVijay Mahadevan 1319566063dSJacob Faibussowitsch PetscCall(DMCoarsen(dm, PetscObjectComm((PetscObject)dm), &dmc[0])); 13248a46eb9SPierre Jolivet for (i = 1; i < nlevels; i++) PetscCall(DMCoarsen(dmc[i - 1], PetscObjectComm((PetscObject)dm), &dmc[i])); 1333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 134b117cd09SVijay Mahadevan } 135b117cd09SVijay Mahadevan 1362417220eSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoab_Compute_NNZ_From_Connectivity(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscBool); 137b117cd09SVijay Mahadevan 138cab5ea25SPierre Jolivet /*@C 139e882eb38SVijay Mahadevan DMCreateInterpolation_Moab - Generate the interpolation operators to transform 140e882eb38SVijay Mahadevan operators (matrices, vectors) from parent level to child level as defined by 141*20f4b53cSBarry Smith the `DM` inputs provided by the user. 142b117cd09SVijay Mahadevan 143d083f849SBarry Smith Collective 144b117cd09SVijay Mahadevan 145d8d19677SJose E. Roman Input Parameters: 146*20f4b53cSBarry Smith + dm1 - The `DMMOAB` object 147*20f4b53cSBarry Smith - dm2 - the second, finer `DMMOAB` object 148b117cd09SVijay Mahadevan 149d8d19677SJose E. Roman Output Parameters: 150e882eb38SVijay Mahadevan + interpl - The interpolation operator for transferring data between the levels 151e882eb38SVijay Mahadevan - vec - The scaling vector (optional) 152b117cd09SVijay Mahadevan 153e882eb38SVijay Mahadevan Level: developer 154b117cd09SVijay Mahadevan 155b117cd09SVijay Mahadevan @*/ 156d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreateInterpolation_Moab(DM dmp, DM dmc, Mat *interpl, Vec *vec) 157d71ae5a4SJacob Faibussowitsch { 158755f3dfbSVijay Mahadevan DM_Moab *dmbp, *dmbc; 159b117cd09SVijay Mahadevan moab::ErrorCode merr; 160e882eb38SVijay Mahadevan PetscInt dim; 1613f1c6e43SVijay Mahadevan PetscReal factor; 162ce27a4eeSVijay Mahadevan PetscInt innz, *nnz, ionz, *onz; 163755f3dfbSVijay Mahadevan PetscInt nlsizp, nlsizc, nlghsizp, ngsizp, ngsizc; 164a86ed394SVijay Mahadevan const PetscBool use_consistent_bases = PETSC_TRUE; 165b117cd09SVijay Mahadevan 166b117cd09SVijay Mahadevan PetscFunctionBegin; 167755f3dfbSVijay Mahadevan PetscValidHeaderSpecific(dmp, DM_CLASSID, 1); 168755f3dfbSVijay Mahadevan PetscValidHeaderSpecific(dmc, DM_CLASSID, 2); 169755f3dfbSVijay Mahadevan dmbp = (DM_Moab *)(dmp)->data; 170755f3dfbSVijay Mahadevan dmbc = (DM_Moab *)(dmc)->data; 171755f3dfbSVijay Mahadevan nlsizp = dmbp->nloc; // *dmb1->numFields; 172755f3dfbSVijay Mahadevan nlsizc = dmbc->nloc; // *dmb2->numFields; 173755f3dfbSVijay Mahadevan ngsizp = dmbp->n; // *dmb1->numFields; 174755f3dfbSVijay Mahadevan ngsizc = dmbc->n; // *dmb2->numFields; 175755f3dfbSVijay Mahadevan nlghsizp = (dmbp->nloc + dmbp->nghost); // *dmb1->numFields; 176b117cd09SVijay Mahadevan 1772417220eSVijay Mahadevan // Columns = Parent DoFs ; Rows = Child DoFs 178755f3dfbSVijay Mahadevan // Interpolation matrix: \sum_{i=1}^P Owned(Child) * (Owned(Parent) + Ghosted(Parent)) 179755f3dfbSVijay Mahadevan // Size: nlsizc * nlghsizp 1803ba16761SJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Creating interpolation matrix %" PetscInt_FMT " X %" PetscInt_FMT " to apply transformation between levels %" PetscInt_FMT " -> %" PetscInt_FMT ".\n", ngsizc, nlghsizp, dmbp->hlevel, dmbc->hlevel)); 181b117cd09SVijay Mahadevan 1829566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dmp, &dim)); 183a86ed394SVijay Mahadevan 184941e0cffSVijay Mahadevan /* allocate the nnz, onz arrays based on block size and local nodes */ 1859566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(nlsizc, &nnz, nlsizc, &onz)); 186941e0cffSVijay Mahadevan 187941e0cffSVijay Mahadevan /* Loop through the local elements and compute the relation between the current parent and the refined_level. */ 1882417220eSVijay Mahadevan for (moab::Range::iterator iter = dmbc->vowned->begin(); iter != dmbc->vowned->end(); iter++) { 1892417220eSVijay Mahadevan const moab::EntityHandle vhandle = *iter; 1902417220eSVijay Mahadevan /* define local variables */ 1912417220eSVijay Mahadevan moab::EntityHandle parent; 1922417220eSVijay Mahadevan std::vector<moab::EntityHandle> adjs; 1932417220eSVijay Mahadevan moab::Range found; 1942417220eSVijay Mahadevan 1952417220eSVijay Mahadevan /* store the vertex DoF number */ 1962417220eSVijay Mahadevan const int ldof = dmbc->lidmap[vhandle - dmbc->seqstart]; 1972417220eSVijay Mahadevan 1982417220eSVijay Mahadevan /* Get adjacency information for current vertex - i.e., all elements of dimension (dim) that connects 1992417220eSVijay Mahadevan to the current vertex. We can then decipher if a vertex is ghosted or not and compute the 2002417220eSVijay Mahadevan non-zero pattern accordingly. */ 2019371c9d4SSatish Balay merr = dmbc->hierarchy->get_adjacencies(vhandle, dmbc->dim, adjs); 2029371c9d4SSatish Balay MBERRNM(merr); 2032417220eSVijay Mahadevan 2042417220eSVijay Mahadevan /* loop over vertices and update the number of connectivity */ 2052417220eSVijay Mahadevan for (unsigned jter = 0; jter < adjs.size(); jter++) { 2062417220eSVijay Mahadevan const moab::EntityHandle jhandle = adjs[jter]; 207941e0cffSVijay Mahadevan 208941e0cffSVijay Mahadevan /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */ 2099371c9d4SSatish Balay merr = dmbc->hierarchy->child_to_parent(jhandle, dmbc->hlevel, dmbp->hlevel, &parent); 2109371c9d4SSatish Balay MBERRNM(merr); 211941e0cffSVijay Mahadevan 2122417220eSVijay Mahadevan /* Get connectivity information in canonical ordering for the local element */ 2132417220eSVijay Mahadevan std::vector<moab::EntityHandle> connp; 2149371c9d4SSatish Balay merr = dmbp->hierarchy->get_connectivity(parent, dmbp->hlevel, connp); 2159371c9d4SSatish Balay MBERRNM(merr); 216a86ed394SVijay Mahadevan 2172417220eSVijay Mahadevan for (unsigned ic = 0; ic < connp.size(); ++ic) { 2182417220eSVijay Mahadevan /* loop over each element connected to the adjacent vertex and update as needed */ 2192417220eSVijay Mahadevan /* find the truly user-expected layer of ghosted entities to decipher NNZ pattern */ 2202417220eSVijay Mahadevan if (found.find(connp[ic]) != found.end()) continue; /* make sure we don't double count shared vertices */ 2212417220eSVijay Mahadevan if (dmbp->vghost->find(connp[ic]) != dmbp->vghost->end()) onz[ldof]++; /* update out-of-proc onz */ 2222417220eSVijay Mahadevan else nnz[ldof]++; /* else local vertex */ 2232417220eSVijay Mahadevan found.insert(connp[ic]); 2242417220eSVijay Mahadevan } 2252417220eSVijay Mahadevan } 226941e0cffSVijay Mahadevan } 227941e0cffSVijay Mahadevan 2289371c9d4SSatish Balay for (int i = 0; i < nlsizc; i++) nnz[i] += 1; /* self count the node */ 229941e0cffSVijay Mahadevan 230941e0cffSVijay Mahadevan ionz = onz[0]; 231941e0cffSVijay Mahadevan innz = nnz[0]; 232755f3dfbSVijay Mahadevan for (int tc = 0; tc < nlsizc; tc++) { 233941e0cffSVijay Mahadevan // check for maximum allowed sparsity = fully dense 234755f3dfbSVijay Mahadevan nnz[tc] = std::min(nlsizp, nnz[tc]); 2352417220eSVijay Mahadevan onz[tc] = std::min(ngsizp - nlsizp, onz[tc]); 2362417220eSVijay Mahadevan 2373ba16761SJacob Faibussowitsch PetscCall(PetscInfo(NULL, " %d: NNZ = %d, ONZ = %d\n", tc, nnz[tc], onz[tc])); 238941e0cffSVijay Mahadevan 239941e0cffSVijay Mahadevan innz = (innz < nnz[tc] ? nnz[tc] : innz); 240941e0cffSVijay Mahadevan ionz = (ionz < onz[tc] ? onz[tc] : ionz); 241941e0cffSVijay Mahadevan } 24264e1c140SVijay Mahadevan 24364e1c140SVijay Mahadevan /* create interpolation matrix */ 2449566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)dmc), interpl)); 2459566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*interpl, nlsizc, nlsizp, ngsizc, ngsizp)); 2469566063dSJacob Faibussowitsch PetscCall(MatSetType(*interpl, MATAIJ)); 2479566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(*interpl)); 24864e1c140SVijay Mahadevan 2499566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(*interpl, innz, nnz)); 2509566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(*interpl, innz, nnz, ionz, onz)); 251941e0cffSVijay Mahadevan 252941e0cffSVijay Mahadevan /* clean up temporary memory */ 2539566063dSJacob Faibussowitsch PetscCall(PetscFree2(nnz, onz)); 254b117cd09SVijay Mahadevan 255b117cd09SVijay Mahadevan /* set up internal matrix data-structures */ 2569566063dSJacob Faibussowitsch PetscCall(MatSetUp(*interpl)); 257b117cd09SVijay Mahadevan 258a86ed394SVijay Mahadevan /* Define variables for assembly */ 259a86ed394SVijay Mahadevan std::vector<moab::EntityHandle> children; 260a86ed394SVijay Mahadevan std::vector<moab::EntityHandle> connp, connc; 261a86ed394SVijay Mahadevan std::vector<PetscReal> pcoords, ccoords, values_phi; 262a86ed394SVijay Mahadevan 263a86ed394SVijay Mahadevan if (use_consistent_bases) { 2642417220eSVijay Mahadevan const moab::EntityHandle ehandle = dmbp->elocal->front(); 265a86ed394SVijay Mahadevan 2669371c9d4SSatish Balay merr = dmbp->hierarchy->parent_to_child(ehandle, dmbp->hlevel, dmbc->hlevel, children); 2679371c9d4SSatish Balay MBERRNM(merr); 268a86ed394SVijay Mahadevan 269a86ed394SVijay Mahadevan /* Get connectivity and coordinates of the parent vertices */ 2709371c9d4SSatish Balay merr = dmbp->hierarchy->get_connectivity(ehandle, dmbp->hlevel, connp); 2719371c9d4SSatish Balay MBERRNM(merr); 2729371c9d4SSatish Balay merr = dmbc->mbiface->get_connectivity(&children[0], children.size(), connc); 2739371c9d4SSatish Balay MBERRNM(merr); 274a86ed394SVijay Mahadevan 275a86ed394SVijay Mahadevan std::vector<PetscReal> natparam(3 * connc.size(), 0.0); 276a86ed394SVijay Mahadevan pcoords.resize(connp.size() * 3); 277a86ed394SVijay Mahadevan ccoords.resize(connc.size() * 3); 278a86ed394SVijay Mahadevan values_phi.resize(connp.size() * connc.size()); 279a86ed394SVijay Mahadevan /* Get coordinates for connectivity entities in canonical order for both coarse and finer levels */ 2809371c9d4SSatish Balay merr = dmbp->hierarchy->get_coordinates(&connp[0], connp.size(), dmbp->hlevel, &pcoords[0]); 2819371c9d4SSatish Balay MBERRNM(merr); 2829371c9d4SSatish Balay merr = dmbc->hierarchy->get_coordinates(&connc[0], connc.size(), dmbc->hlevel, &ccoords[0]); 2839371c9d4SSatish Balay MBERRNM(merr); 284a86ed394SVijay Mahadevan 285a86ed394SVijay Mahadevan /* Set values: For each DOF in coarse grid cell, set the contribution or PHI evaluated at each fine grid DOF point */ 286a86ed394SVijay Mahadevan for (unsigned tc = 0; tc < connc.size(); tc++) { 287a86ed394SVijay Mahadevan const PetscInt offset = tc * 3; 288a86ed394SVijay Mahadevan 289a86ed394SVijay Mahadevan /* Scale ccoords relative to pcoords */ 2909566063dSJacob Faibussowitsch PetscCall(DMMoabPToRMapping(dim, connp.size(), &pcoords[0], &ccoords[offset], &natparam[offset], &values_phi[connp.size() * tc])); 291a86ed394SVijay Mahadevan } 2921baa6e33SBarry Smith } else { 293a86ed394SVijay Mahadevan factor = std::pow(2.0 /*degree_P_for_refinement*/, (dmbc->hlevel - dmbp->hlevel) * dmbp->dim * 1.0); 294a86ed394SVijay Mahadevan } 295a86ed394SVijay Mahadevan 296755f3dfbSVijay Mahadevan /* TODO: Decipher the correct non-zero pattern. There is still some issue with onz allocation */ 2979566063dSJacob Faibussowitsch PetscCall(MatSetOption(*interpl, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 298b117cd09SVijay Mahadevan 299e882eb38SVijay Mahadevan /* Loop through the remaining vertices. These vertices appear only on the current refined_level. */ 3002417220eSVijay Mahadevan for (moab::Range::iterator iter = dmbp->elocal->begin(); iter != dmbp->elocal->end(); iter++) { 301b117cd09SVijay Mahadevan const moab::EntityHandle ehandle = *iter; 302b117cd09SVijay Mahadevan 303b117cd09SVijay Mahadevan /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */ 304a86ed394SVijay Mahadevan children.clear(); 305a86ed394SVijay Mahadevan connc.clear(); 3069371c9d4SSatish Balay merr = dmbp->hierarchy->parent_to_child(ehandle, dmbp->hlevel, dmbc->hlevel, children); 3079371c9d4SSatish Balay MBERRNM(merr); 308b117cd09SVijay Mahadevan 309b117cd09SVijay Mahadevan /* Get connectivity and coordinates of the parent vertices */ 3109371c9d4SSatish Balay merr = dmbp->hierarchy->get_connectivity(ehandle, dmbp->hlevel, connp); 3119371c9d4SSatish Balay MBERRNM(merr); 3129371c9d4SSatish Balay merr = dmbc->mbiface->get_connectivity(&children[0], children.size(), connc); 3139371c9d4SSatish Balay MBERRNM(merr); 314b117cd09SVijay Mahadevan 315a86ed394SVijay Mahadevan pcoords.resize(connp.size() * 3); 316a86ed394SVijay Mahadevan ccoords.resize(connc.size() * 3); 317b117cd09SVijay Mahadevan /* Get coordinates for connectivity entities in canonical order for both coarse and finer levels */ 3189371c9d4SSatish Balay merr = dmbp->hierarchy->get_coordinates(&connp[0], connp.size(), dmbp->hlevel, &pcoords[0]); 3199371c9d4SSatish Balay MBERRNM(merr); 3209371c9d4SSatish Balay merr = dmbc->hierarchy->get_coordinates(&connc[0], connc.size(), dmbc->hlevel, &ccoords[0]); 3219371c9d4SSatish Balay MBERRNM(merr); 322b117cd09SVijay Mahadevan 323b117cd09SVijay Mahadevan std::vector<int> dofsp(connp.size()), dofsc(connc.size()); 324b117cd09SVijay Mahadevan /* TODO: specific to scalar system - use GetDofs */ 3259566063dSJacob Faibussowitsch PetscCall(DMMoabGetDofsBlocked(dmp, connp.size(), &connp[0], &dofsp[0])); 3269566063dSJacob Faibussowitsch PetscCall(DMMoabGetDofsBlocked(dmc, connc.size(), &connc[0], &dofsc[0])); 327b117cd09SVijay Mahadevan 328a86ed394SVijay Mahadevan /* Compute the actual interpolation weights when projecting solution/residual between levels */ 329a86ed394SVijay Mahadevan if (use_consistent_bases) { 330a86ed394SVijay Mahadevan /* Use the cached values of natural parameteric coordinates and basis pre-evaluated. 331a86ed394SVijay Mahadevan We are making an assumption here that UMR used in GMG to generate the hierarchy uses 332a86ed394SVijay Mahadevan the same template for all elements; This will fail for mixed element meshes (TRI/QUAD). 333a86ed394SVijay Mahadevan 3342417220eSVijay Mahadevan TODO: Fix the above assumption by caching data for families (especially for Tets and mixed meshes) 335a86ed394SVijay Mahadevan */ 336a86ed394SVijay Mahadevan 337a86ed394SVijay Mahadevan /* Set values: For each DOF in coarse grid cell, set the contribution or PHI evaluated at each fine grid DOF point */ 338a86ed394SVijay Mahadevan for (unsigned tc = 0; tc < connc.size(); tc++) { 339a86ed394SVijay Mahadevan /* TODO: Check if we should be using INSERT_VALUES instead */ 3409566063dSJacob Faibussowitsch PetscCall(MatSetValues(*interpl, 1, &dofsc[tc], connp.size(), &dofsp[0], &values_phi[connp.size() * tc], ADD_VALUES)); 341a86ed394SVijay Mahadevan } 3421baa6e33SBarry Smith } else { 343e882eb38SVijay Mahadevan /* Compute the interpolation weights by determining distance of 1-ring 344755f3dfbSVijay Mahadevan neighbor vertices from current vertex 345755f3dfbSVijay Mahadevan 346a86ed394SVijay Mahadevan This should be used only when FEM basis is not used for the discretization. 347a86ed394SVijay Mahadevan Else, the consistent interface to compute the basis function for interpolation 348a86ed394SVijay Mahadevan between the levels should be evaluated correctly to preserve convergence of GMG. 3492417220eSVijay Mahadevan Shephard's basis will be terrible for any unsmooth problems. 350755f3dfbSVijay Mahadevan */ 351755f3dfbSVijay Mahadevan values_phi.resize(connp.size()); 352755f3dfbSVijay Mahadevan for (unsigned tc = 0; tc < connc.size(); tc++) { 353a86ed394SVijay Mahadevan PetscReal normsum = 0.0; 354755f3dfbSVijay Mahadevan for (unsigned tp = 0; tp < connp.size(); tp++) { 355755f3dfbSVijay Mahadevan values_phi[tp] = 0.0; 3569371c9d4SSatish Balay for (unsigned k = 0; k < 3; k++) values_phi[tp] += std::pow(pcoords[tp * 3 + k] - ccoords[k + tc * 3], dim); 357755f3dfbSVijay Mahadevan if (values_phi[tp] < 1e-12) { 358755f3dfbSVijay Mahadevan values_phi[tp] = 1e12; 3591baa6e33SBarry Smith } else { 360755f3dfbSVijay Mahadevan //values_phi[tp] = std::pow(values_phi[tp], -1.0/dim); 361755f3dfbSVijay Mahadevan values_phi[tp] = std::pow(values_phi[tp], -1.0); 362755f3dfbSVijay Mahadevan normsum += values_phi[tp]; 363b117cd09SVijay Mahadevan } 364b117cd09SVijay Mahadevan } 365755f3dfbSVijay Mahadevan for (unsigned tp = 0; tp < connp.size(); tp++) { 3669371c9d4SSatish Balay if (values_phi[tp] > 1e11) values_phi[tp] = factor * 0.5 / connp.size(); 3679371c9d4SSatish Balay else values_phi[tp] = factor * values_phi[tp] * 0.5 / (connp.size() * normsum); 368b117cd09SVijay Mahadevan } 3699566063dSJacob Faibussowitsch PetscCall(MatSetValues(*interpl, 1, &dofsc[tc], connp.size(), &dofsp[0], &values_phi[0], ADD_VALUES)); 370b117cd09SVijay Mahadevan } 371a86ed394SVijay Mahadevan } 372b117cd09SVijay Mahadevan } 373304006b3SVijay Mahadevan if (vec) *vec = NULL; 3749566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*interpl, MAT_FINAL_ASSEMBLY)); 3759566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*interpl, MAT_FINAL_ASSEMBLY)); 3763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 377b117cd09SVijay Mahadevan } 378b117cd09SVijay Mahadevan 379cab5ea25SPierre Jolivet /*@C 380b117cd09SVijay Mahadevan DMCreateInjection_Moab - Generate a multi-level uniform refinement hierarchy 381*20f4b53cSBarry Smith by succesively refining a coarse mesh, already defined in the `DM` object 382b117cd09SVijay Mahadevan provided by the user. 383b117cd09SVijay Mahadevan 384d083f849SBarry Smith Collective 385b117cd09SVijay Mahadevan 386b117cd09SVijay Mahadevan Input Parameter: 387*20f4b53cSBarry Smith . dmb - The `DMMOAB` object 388b117cd09SVijay Mahadevan 389d8d19677SJose E. Roman Output Parameters: 390a2b725a8SWilliam Gropp + nlevels - The number of levels of refinement needed to generate the hierarchy 391a2b725a8SWilliam Gropp - ldegrees - The degree of refinement at each level in the hierarchy 392b117cd09SVijay Mahadevan 393b117cd09SVijay Mahadevan Level: beginner 394b117cd09SVijay Mahadevan 395b117cd09SVijay Mahadevan @*/ 396d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreateInjection_Moab(DM dm1, DM dm2, VecScatter *ctx) 397d71ae5a4SJacob Faibussowitsch { 398e882eb38SVijay Mahadevan //DM_Moab *dmmoab; 399b117cd09SVijay Mahadevan 400b117cd09SVijay Mahadevan PetscFunctionBegin; 401b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm1, DM_CLASSID, 1); 402b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm2, DM_CLASSID, 2); 403e882eb38SVijay Mahadevan //dmmoab = (DM_Moab*)(dm1)->data; 404b117cd09SVijay Mahadevan 4053ba16761SJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[DMCreateInjection_Moab] :: Placeholder\n")); 4063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 407b117cd09SVijay Mahadevan } 408b117cd09SVijay Mahadevan 409d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMMoab_UMR_Private(DM dm, MPI_Comm comm, PetscBool refine, DM *dmref) 410d71ae5a4SJacob Faibussowitsch { 411e882eb38SVijay Mahadevan PetscInt i, dim; 412b117cd09SVijay Mahadevan DM dm2; 413b117cd09SVijay Mahadevan moab::ErrorCode merr; 414b117cd09SVijay Mahadevan DM_Moab *dmb = (DM_Moab *)dm->data, *dd2; 415b117cd09SVijay Mahadevan 416b117cd09SVijay Mahadevan PetscFunctionBegin; 417b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 418b117cd09SVijay Mahadevan PetscValidPointer(dmref, 4); 419b117cd09SVijay Mahadevan 420b117cd09SVijay Mahadevan if ((dmb->hlevel == dmb->nhlevels && refine) || (dmb->hlevel == 0 && !refine)) { 4213ba16761SJacob Faibussowitsch if (dmb->hlevel + 1 > dmb->nhlevels && refine) { 4223ba16761SJacob Faibussowitsch PetscCall(PetscInfo(NULL, "Invalid multigrid refinement hierarchy level specified (%" PetscInt_FMT "). MOAB UMR max levels = %" PetscInt_FMT ". Creating a NULL object.\n", dmb->hlevel + 1, dmb->nhlevels)); 4233ba16761SJacob Faibussowitsch } 4243ba16761SJacob Faibussowitsch if (dmb->hlevel - 1 < 0 && !refine) PetscCall(PetscInfo(NULL, "Invalid multigrid coarsen hierarchy level specified (%" PetscInt_FMT "). Creating a NULL object.\n", dmb->hlevel - 1)); 425f3fa974cSJacob Faibussowitsch *dmref = NULL; 4263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 427b117cd09SVijay Mahadevan } 428b117cd09SVijay Mahadevan 4299566063dSJacob Faibussowitsch PetscCall(DMMoabCreate(PetscObjectComm((PetscObject)dm), &dm2)); 430b117cd09SVijay Mahadevan dd2 = (DM_Moab *)dm2->data; 431b117cd09SVijay Mahadevan 432b117cd09SVijay Mahadevan dd2->mbiface = dmb->mbiface; 4339daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 434b117cd09SVijay Mahadevan dd2->pcomm = dmb->pcomm; 4359daf19fdSVijay Mahadevan #endif 436b117cd09SVijay Mahadevan dd2->icreatedinstance = PETSC_FALSE; 43764e1c140SVijay Mahadevan dd2->nghostrings = dmb->nghostrings; 438b117cd09SVijay Mahadevan 439b117cd09SVijay Mahadevan /* set the new level based on refinement/coarsening */ 440b117cd09SVijay Mahadevan if (refine) { 441b117cd09SVijay Mahadevan dd2->hlevel = dmb->hlevel + 1; 4421baa6e33SBarry Smith } else { 443b117cd09SVijay Mahadevan dd2->hlevel = dmb->hlevel - 1; 444b117cd09SVijay Mahadevan } 445b117cd09SVijay Mahadevan 446b117cd09SVijay Mahadevan /* Copy the multilevel hierarchy pointers in MOAB */ 447b117cd09SVijay Mahadevan dd2->hierarchy = dmb->hierarchy; 448b117cd09SVijay Mahadevan dd2->nhlevels = dmb->nhlevels; 4499566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd2->nhlevels + 1, &dd2->hsets)); 450ad540459SPierre Jolivet for (i = 0; i <= dd2->nhlevels; i++) dd2->hsets[i] = dmb->hsets[i]; 451b117cd09SVijay Mahadevan dd2->fileset = dd2->hsets[dd2->hlevel]; 452b117cd09SVijay Mahadevan 453b117cd09SVijay Mahadevan /* do the remaining initializations for DMMoab */ 454b117cd09SVijay Mahadevan dd2->bs = dmb->bs; 455b117cd09SVijay Mahadevan dd2->numFields = dmb->numFields; 456b117cd09SVijay Mahadevan dd2->rw_dbglevel = dmb->rw_dbglevel; 457b117cd09SVijay Mahadevan dd2->partition_by_rank = dmb->partition_by_rank; 458c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(dd2->extra_read_options, dmb->extra_read_options, sizeof(dd2->extra_read_options))); 459c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(dd2->extra_write_options, dmb->extra_write_options, sizeof(dd2->extra_write_options))); 460b117cd09SVijay Mahadevan dd2->read_mode = dmb->read_mode; 461b117cd09SVijay Mahadevan dd2->write_mode = dmb->write_mode; 462b117cd09SVijay Mahadevan 463b117cd09SVijay Mahadevan /* set global ID tag handle */ 4649566063dSJacob Faibussowitsch PetscCall(DMMoabSetLocalToGlobalTag(dm2, dmb->ltog_tag)); 465b117cd09SVijay Mahadevan 4669371c9d4SSatish Balay merr = dd2->mbiface->tag_get_handle(MATERIAL_SET_TAG_NAME, dd2->material_tag); 4679371c9d4SSatish Balay MBERRNM(merr); 468b117cd09SVijay Mahadevan 4699566063dSJacob Faibussowitsch PetscCall(DMSetOptionsPrefix(dm2, ((PetscObject)dm)->prefix)); 4709566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 4719566063dSJacob Faibussowitsch PetscCall(DMSetDimension(dm2, dim)); 472b117cd09SVijay Mahadevan 473b117cd09SVijay Mahadevan /* allow overloaded (user replaced) operations to be inherited by refinement clones */ 474b117cd09SVijay Mahadevan dm2->ops->creatematrix = dm->ops->creatematrix; 475b117cd09SVijay Mahadevan 476b117cd09SVijay Mahadevan /* copy fill information if given */ 4779566063dSJacob Faibussowitsch PetscCall(DMMoabSetBlockFills(dm2, dmb->dfill, dmb->ofill)); 478b117cd09SVijay Mahadevan 479b117cd09SVijay Mahadevan /* copy vector type information */ 4809566063dSJacob Faibussowitsch PetscCall(DMSetMatType(dm2, dm->mattype)); 4819566063dSJacob Faibussowitsch PetscCall(DMSetVecType(dm2, dm->vectype)); 4823f1c6e43SVijay Mahadevan dd2->numFields = dmb->numFields; 48348a46eb9SPierre Jolivet if (dmb->numFields) PetscCall(DMMoabSetFieldNames(dm2, dmb->numFields, dmb->fieldNames)); 484b117cd09SVijay Mahadevan 4859566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm2)); 486b117cd09SVijay Mahadevan 487b117cd09SVijay Mahadevan /* recreate Dof numbering for the refined DM and make sure the distribution is correctly populated */ 4889566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm2)); 489b117cd09SVijay Mahadevan 490b117cd09SVijay Mahadevan *dmref = dm2; 4913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 492b117cd09SVijay Mahadevan } 493b117cd09SVijay Mahadevan 494cab5ea25SPierre Jolivet /*@C 495b117cd09SVijay Mahadevan DMRefine_Moab - Generate a multi-level uniform refinement hierarchy 496*20f4b53cSBarry Smith by succesively refining a coarse mesh, already defined in the `DM` object 497b117cd09SVijay Mahadevan provided by the user. 498b117cd09SVijay Mahadevan 499*20f4b53cSBarry Smith Collective 500b117cd09SVijay Mahadevan 501d8d19677SJose E. Roman Input Parameters: 502*20f4b53cSBarry Smith + dm - The `DMMOAB` object 503*20f4b53cSBarry Smith - comm - the communicator to contain the new DM object (or `MPI_COMM_NULL`) 504b117cd09SVijay Mahadevan 505b117cd09SVijay Mahadevan Output Parameter: 506*20f4b53cSBarry Smith . dmf - the refined `DM`, or `NULL` 507e882eb38SVijay Mahadevan 508e882eb38SVijay Mahadevan Level: developer 509b117cd09SVijay Mahadevan 510*20f4b53cSBarry Smith Note: 511*20f4b53cSBarry Smith If no refinement was done, the return value is `NULL` 512*20f4b53cSBarry Smith 513b117cd09SVijay Mahadevan @*/ 514d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMRefine_Moab(DM dm, MPI_Comm comm, DM *dmf) 515d71ae5a4SJacob Faibussowitsch { 516b117cd09SVijay Mahadevan PetscFunctionBegin; 517b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 518b117cd09SVijay Mahadevan 5199566063dSJacob Faibussowitsch PetscCall(DMMoab_UMR_Private(dm, comm, PETSC_TRUE, dmf)); 5203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 521b117cd09SVijay Mahadevan } 522b117cd09SVijay Mahadevan 523cab5ea25SPierre Jolivet /*@C 524b117cd09SVijay Mahadevan DMCoarsen_Moab - Generate a multi-level uniform refinement hierarchy 525*20f4b53cSBarry Smith by succesively refining a coarse mesh, already defined in the `DM` object 526b117cd09SVijay Mahadevan provided by the user. 527b117cd09SVijay Mahadevan 528*20f4b53cSBarry Smith Collective 529b117cd09SVijay Mahadevan 530d8d19677SJose E. Roman Input Parameters: 531*20f4b53cSBarry Smith + dm - The `DMMOAB` object 532*20f4b53cSBarry Smith - comm - the communicator to contain the new `DM` object (or `MPI_COMM_NULL`) 533b117cd09SVijay Mahadevan 534b117cd09SVijay Mahadevan Output Parameter: 535*20f4b53cSBarry Smith . dmf - the coarsened `DM`, or `NULL` 536b117cd09SVijay Mahadevan 537e882eb38SVijay Mahadevan Level: developer 538e882eb38SVijay Mahadevan 539*20f4b53cSBarry Smith Note: 540*20f4b53cSBarry Smith If no coarsening was done, the return value is `NULL` 541*20f4b53cSBarry Smith 542b117cd09SVijay Mahadevan @*/ 543d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCoarsen_Moab(DM dm, MPI_Comm comm, DM *dmc) 544d71ae5a4SJacob Faibussowitsch { 545b117cd09SVijay Mahadevan PetscFunctionBegin; 546b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5479566063dSJacob Faibussowitsch PetscCall(DMMoab_UMR_Private(dm, comm, PETSC_FALSE, dmc)); 5483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 549b117cd09SVijay Mahadevan } 550