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 8b117cd09SVijay Mahadevan 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: 14a2b725a8SWilliam Gropp . 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 @*/ 23b117cd09SVijay Mahadevan PetscErrorCode DMMoabGenerateHierarchy(DM dm, PetscInt nlevels, PetscInt *ldegrees) 24b117cd09SVijay Mahadevan { 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 */ 37b117cd09SVijay Mahadevan } 38b117cd09SVijay Mahadevan else pdegrees = ldegrees; 39b117cd09SVijay Mahadevan 40b117cd09SVijay Mahadevan /* initialize set level refinement data for hierarchy */ 41b117cd09SVijay Mahadevan dmmoab->nhlevels = nlevels; 42b117cd09SVijay Mahadevan 43b117cd09SVijay Mahadevan /* Instantiate the nested refinement class */ 449daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 453f1c6e43SVijay Mahadevan dmmoab->hierarchy = new moab::NestedRefine(dynamic_cast<moab::Core*>(dmmoab->mbiface), dmmoab->pcomm, dmmoab->fileset); 469daf19fdSVijay Mahadevan #else 479daf19fdSVijay Mahadevan dmmoab->hierarchy = new moab::NestedRefine(dynamic_cast<moab::Core*>(dmmoab->mbiface), NULL, dmmoab->fileset); 489daf19fdSVijay Mahadevan #endif 49b117cd09SVijay Mahadevan 509566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nlevels + 1, &dmmoab->hsets)); 51e882eb38SVijay Mahadevan 52b117cd09SVijay Mahadevan /* generate the mesh hierarchy */ 539c368985SVijay Mahadevan merr = dmmoab->hierarchy->generate_mesh_hierarchy(nlevels, pdegrees, hsets, false);MBERRNM(merr); 54e882eb38SVijay Mahadevan 559daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 56755f3dfbSVijay Mahadevan if (dmmoab->pcomm->size() > 1) { 5764e1c140SVijay Mahadevan merr = dmmoab->hierarchy->exchange_ghosts(hsets, dmmoab->nghostrings);MBERRNM(merr); 58755f3dfbSVijay Mahadevan } 599daf19fdSVijay Mahadevan #endif 6049d66b22SVijay Mahadevan 6164e1c140SVijay Mahadevan /* copy the mesh sets for nested refinement hierarchy */ 62e92d1c7cSVijay Mahadevan dmmoab->hsets[0] = hsets[0]; 63e92d1c7cSVijay Mahadevan for (ilevel = 1; ilevel <= nlevels; ilevel++) 642417220eSVijay Mahadevan { 652417220eSVijay Mahadevan dmmoab->hsets[ilevel] = hsets[ilevel]; 66e92d1c7cSVijay Mahadevan 679c368985SVijay Mahadevan #ifdef MOAB_HAVE_MPI 68e92d1c7cSVijay Mahadevan merr = dmmoab->pcomm->assign_global_ids(hsets[ilevel], dmmoab->dim, 0, false, true, false);MBERRNM(merr); 699c368985SVijay Mahadevan #endif 70e92d1c7cSVijay Mahadevan 71e92d1c7cSVijay Mahadevan /* Update material and other geometric tags from parent to child sets */ 72e92d1c7cSVijay Mahadevan merr = dmmoab->hierarchy->update_special_tags(ilevel, hsets[ilevel]);MBERRNM(merr); 732417220eSVijay Mahadevan } 7464e1c140SVijay Mahadevan 7564e1c140SVijay Mahadevan hsets.clear(); 76b117cd09SVijay Mahadevan if (!ldegrees) { 779566063dSJacob Faibussowitsch PetscCall(PetscFree(pdegrees)); 78b117cd09SVijay Mahadevan } 79b117cd09SVijay Mahadevan PetscFunctionReturn(0); 80b117cd09SVijay Mahadevan } 81b117cd09SVijay Mahadevan 82cab5ea25SPierre Jolivet /*@C 83e882eb38SVijay Mahadevan 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: 89a2b725a8SWilliam Gropp . 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 @*/ 98304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMRefineHierarchy_Moab(DM dm, PetscInt nlevels, DM dmf[]) 99b117cd09SVijay Mahadevan { 100b117cd09SVijay Mahadevan PetscInt i; 101b117cd09SVijay Mahadevan 102b117cd09SVijay Mahadevan PetscFunctionBegin; 103b117cd09SVijay Mahadevan 1049566063dSJacob Faibussowitsch PetscCall(DMRefine(dm, PetscObjectComm((PetscObject)dm), &dmf[0])); 105e882eb38SVijay Mahadevan for (i = 1; i < nlevels; i++) { 1069566063dSJacob Faibussowitsch PetscCall(DMRefine(dmf[i - 1], PetscObjectComm((PetscObject)dm), &dmf[i])); 107b117cd09SVijay Mahadevan } 108b117cd09SVijay Mahadevan PetscFunctionReturn(0); 109b117cd09SVijay Mahadevan } 110b117cd09SVijay Mahadevan 111cab5ea25SPierre Jolivet /*@C 112e882eb38SVijay Mahadevan DMCoarsenHierarchy_Moab - Generate a multi-level DM hierarchy 113e882eb38SVijay Mahadevan by succesively coarsening a refined mesh. 114b117cd09SVijay Mahadevan 115d083f849SBarry Smith Collective 116b117cd09SVijay Mahadevan 117b117cd09SVijay Mahadevan Input Parameter: 118a2b725a8SWilliam Gropp . dm - The DMMoab object 119b117cd09SVijay Mahadevan 120d8d19677SJose E. Roman Output Parameters: 121e882eb38SVijay Mahadevan + nlevels - The number of levels of refinement needed to generate the hierarchy 122a2b725a8SWilliam Gropp - dmc - The DM objects after successive coarsening of the hierarchy 123b117cd09SVijay Mahadevan 124b117cd09SVijay Mahadevan Level: beginner 125b117cd09SVijay Mahadevan 126b117cd09SVijay Mahadevan @*/ 127304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMCoarsenHierarchy_Moab(DM dm, PetscInt nlevels, DM dmc[]) 128b117cd09SVijay Mahadevan { 129b117cd09SVijay Mahadevan PetscInt i; 130b117cd09SVijay Mahadevan 131b117cd09SVijay Mahadevan PetscFunctionBegin; 132b117cd09SVijay Mahadevan 1339566063dSJacob Faibussowitsch PetscCall(DMCoarsen(dm, PetscObjectComm((PetscObject)dm), &dmc[0])); 134e882eb38SVijay Mahadevan for (i = 1; i < nlevels; i++) { 1359566063dSJacob Faibussowitsch PetscCall(DMCoarsen(dmc[i - 1], PetscObjectComm((PetscObject)dm), &dmc[i])); 136b117cd09SVijay Mahadevan } 137b117cd09SVijay Mahadevan PetscFunctionReturn(0); 138b117cd09SVijay Mahadevan } 139b117cd09SVijay Mahadevan 1402417220eSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoab_Compute_NNZ_From_Connectivity(DM, PetscInt*, PetscInt*, PetscInt*, PetscInt*, PetscBool); 141b117cd09SVijay Mahadevan 142cab5ea25SPierre Jolivet /*@C 143e882eb38SVijay Mahadevan DMCreateInterpolation_Moab - Generate the interpolation operators to transform 144e882eb38SVijay Mahadevan operators (matrices, vectors) from parent level to child level as defined by 145e882eb38SVijay Mahadevan the DM inputs provided by the user. 146b117cd09SVijay Mahadevan 147d083f849SBarry Smith Collective 148b117cd09SVijay Mahadevan 149d8d19677SJose E. Roman Input Parameters: 150e882eb38SVijay Mahadevan + dm1 - The DMMoab object 151e882eb38SVijay Mahadevan - dm2 - the second, finer DMMoab object 152b117cd09SVijay Mahadevan 153d8d19677SJose E. Roman Output Parameters: 154e882eb38SVijay Mahadevan + interpl - The interpolation operator for transferring data between the levels 155e882eb38SVijay Mahadevan - vec - The scaling vector (optional) 156b117cd09SVijay Mahadevan 157e882eb38SVijay Mahadevan Level: developer 158b117cd09SVijay Mahadevan 159b117cd09SVijay Mahadevan @*/ 160304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMCreateInterpolation_Moab(DM dmp, DM dmc, Mat* interpl, Vec* vec) 161b117cd09SVijay Mahadevan { 162755f3dfbSVijay Mahadevan DM_Moab *dmbp, *dmbc; 163b117cd09SVijay Mahadevan moab::ErrorCode merr; 164e882eb38SVijay Mahadevan PetscInt dim; 1653f1c6e43SVijay Mahadevan PetscReal factor; 166ce27a4eeSVijay Mahadevan PetscInt innz, *nnz, ionz, *onz; 167755f3dfbSVijay Mahadevan PetscInt nlsizp, nlsizc, nlghsizp, ngsizp, ngsizc; 168a86ed394SVijay Mahadevan const PetscBool use_consistent_bases=PETSC_TRUE; 169b117cd09SVijay Mahadevan 170b117cd09SVijay Mahadevan PetscFunctionBegin; 171755f3dfbSVijay Mahadevan PetscValidHeaderSpecific(dmp, DM_CLASSID, 1); 172755f3dfbSVijay Mahadevan PetscValidHeaderSpecific(dmc, DM_CLASSID, 2); 173755f3dfbSVijay Mahadevan dmbp = (DM_Moab*)(dmp)->data; 174755f3dfbSVijay Mahadevan dmbc = (DM_Moab*)(dmc)->data; 175755f3dfbSVijay Mahadevan nlsizp = dmbp->nloc;// *dmb1->numFields; 176755f3dfbSVijay Mahadevan nlsizc = dmbc->nloc;// *dmb2->numFields; 177755f3dfbSVijay Mahadevan ngsizp = dmbp->n; // *dmb1->numFields; 178755f3dfbSVijay Mahadevan ngsizc = dmbc->n; // *dmb2->numFields; 179755f3dfbSVijay Mahadevan nlghsizp = (dmbp->nloc + dmbp->nghost); // *dmb1->numFields; 180b117cd09SVijay Mahadevan 1812417220eSVijay Mahadevan // Columns = Parent DoFs ; Rows = Child DoFs 182755f3dfbSVijay Mahadevan // Interpolation matrix: \sum_{i=1}^P Owned(Child) * (Owned(Parent) + Ghosted(Parent)) 183755f3dfbSVijay Mahadevan // Size: nlsizc * nlghsizp 184*63a3b9bcSJacob Faibussowitsch 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); 185b117cd09SVijay Mahadevan 1869566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dmp, &dim)); 187a86ed394SVijay Mahadevan 188941e0cffSVijay Mahadevan /* allocate the nnz, onz arrays based on block size and local nodes */ 1899566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(nlsizc, &nnz, nlsizc, &onz)); 190941e0cffSVijay Mahadevan 191941e0cffSVijay Mahadevan /* Loop through the local elements and compute the relation between the current parent and the refined_level. */ 1922417220eSVijay Mahadevan for (moab::Range::iterator iter = dmbc->vowned->begin(); iter != dmbc->vowned->end(); iter++) { 193941e0cffSVijay Mahadevan 1942417220eSVijay Mahadevan const moab::EntityHandle vhandle = *iter; 1952417220eSVijay Mahadevan /* define local variables */ 1962417220eSVijay Mahadevan moab::EntityHandle parent; 1972417220eSVijay Mahadevan std::vector<moab::EntityHandle> adjs; 1982417220eSVijay Mahadevan moab::Range found; 1992417220eSVijay Mahadevan 2002417220eSVijay Mahadevan /* store the vertex DoF number */ 2012417220eSVijay Mahadevan const int ldof = dmbc->lidmap[vhandle - dmbc->seqstart]; 2022417220eSVijay Mahadevan 2032417220eSVijay Mahadevan /* Get adjacency information for current vertex - i.e., all elements of dimension (dim) that connects 2042417220eSVijay Mahadevan to the current vertex. We can then decipher if a vertex is ghosted or not and compute the 2052417220eSVijay Mahadevan non-zero pattern accordingly. */ 2062417220eSVijay Mahadevan merr = dmbc->hierarchy->get_adjacencies(vhandle, dmbc->dim, adjs);MBERRNM(merr); 2072417220eSVijay Mahadevan 2082417220eSVijay Mahadevan /* loop over vertices and update the number of connectivity */ 2092417220eSVijay Mahadevan for (unsigned jter = 0; jter < adjs.size(); jter++) { 2102417220eSVijay Mahadevan 2112417220eSVijay Mahadevan const moab::EntityHandle jhandle = adjs[jter]; 212941e0cffSVijay Mahadevan 213941e0cffSVijay Mahadevan /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */ 2142417220eSVijay Mahadevan merr = dmbc->hierarchy->child_to_parent(jhandle, dmbc->hlevel, dmbp->hlevel, &parent);MBERRNM(merr); 215941e0cffSVijay Mahadevan 2162417220eSVijay Mahadevan /* Get connectivity information in canonical ordering for the local element */ 2172417220eSVijay Mahadevan std::vector<moab::EntityHandle> connp; 2182417220eSVijay Mahadevan merr = dmbp->hierarchy->get_connectivity(parent, dmbp->hlevel, connp);MBERRNM(merr); 219a86ed394SVijay Mahadevan 2202417220eSVijay Mahadevan for (unsigned ic=0; ic < connp.size(); ++ic) { 2212417220eSVijay Mahadevan 2222417220eSVijay Mahadevan /* loop over each element connected to the adjacent vertex and update as needed */ 2232417220eSVijay Mahadevan /* find the truly user-expected layer of ghosted entities to decipher NNZ pattern */ 2242417220eSVijay Mahadevan if (found.find(connp[ic]) != found.end()) continue; /* make sure we don't double count shared vertices */ 2252417220eSVijay Mahadevan if (dmbp->vghost->find(connp[ic]) != dmbp->vghost->end()) onz[ldof]++; /* update out-of-proc onz */ 2262417220eSVijay Mahadevan else nnz[ldof]++; /* else local vertex */ 2272417220eSVijay Mahadevan found.insert(connp[ic]); 2282417220eSVijay Mahadevan } 2292417220eSVijay Mahadevan } 230941e0cffSVijay Mahadevan } 231941e0cffSVijay Mahadevan 2322417220eSVijay Mahadevan for (int i = 0; i < nlsizc; i++) 2332417220eSVijay Mahadevan nnz[i] += 1; /* self count the node */ 234941e0cffSVijay Mahadevan 235941e0cffSVijay Mahadevan ionz = onz[0]; 236941e0cffSVijay Mahadevan innz = nnz[0]; 237755f3dfbSVijay Mahadevan for (int tc = 0; tc < nlsizc; tc++) { 238941e0cffSVijay Mahadevan // check for maximum allowed sparsity = fully dense 239755f3dfbSVijay Mahadevan nnz[tc] = std::min(nlsizp, nnz[tc]); 2402417220eSVijay Mahadevan onz[tc] = std::min(ngsizp - nlsizp, onz[tc]); 2412417220eSVijay Mahadevan 2427d3de750SJacob Faibussowitsch PetscInfo(NULL, " %d: NNZ = %d, ONZ = %d\n", tc, nnz[tc], onz[tc]); 243941e0cffSVijay Mahadevan 244941e0cffSVijay Mahadevan innz = (innz < nnz[tc] ? nnz[tc] : innz); 245941e0cffSVijay Mahadevan ionz = (ionz < onz[tc] ? onz[tc] : ionz); 246941e0cffSVijay Mahadevan } 24764e1c140SVijay Mahadevan 24864e1c140SVijay Mahadevan /* create interpolation matrix */ 2499566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)dmc), interpl)); 2509566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*interpl, nlsizc, nlsizp, ngsizc, ngsizp)); 2519566063dSJacob Faibussowitsch PetscCall(MatSetType(*interpl, MATAIJ)); 2529566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(*interpl)); 25364e1c140SVijay Mahadevan 2549566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(*interpl, innz, nnz)); 2559566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(*interpl, innz, nnz, ionz, onz)); 256941e0cffSVijay Mahadevan 257941e0cffSVijay Mahadevan /* clean up temporary memory */ 2589566063dSJacob Faibussowitsch PetscCall(PetscFree2(nnz, onz)); 259b117cd09SVijay Mahadevan 260b117cd09SVijay Mahadevan /* set up internal matrix data-structures */ 2619566063dSJacob Faibussowitsch PetscCall(MatSetUp(*interpl)); 262b117cd09SVijay Mahadevan 263a86ed394SVijay Mahadevan /* Define variables for assembly */ 264a86ed394SVijay Mahadevan std::vector<moab::EntityHandle> children; 265a86ed394SVijay Mahadevan std::vector<moab::EntityHandle> connp, connc; 266a86ed394SVijay Mahadevan std::vector<PetscReal> pcoords, ccoords, values_phi; 267a86ed394SVijay Mahadevan 268a86ed394SVijay Mahadevan if (use_consistent_bases) { 2692417220eSVijay Mahadevan const moab::EntityHandle ehandle = dmbp->elocal->front(); 270a86ed394SVijay Mahadevan 271a86ed394SVijay Mahadevan merr = dmbp->hierarchy->parent_to_child(ehandle, dmbp->hlevel, dmbc->hlevel, children);MBERRNM(merr); 272a86ed394SVijay Mahadevan 273a86ed394SVijay Mahadevan /* Get connectivity and coordinates of the parent vertices */ 274a86ed394SVijay Mahadevan merr = dmbp->hierarchy->get_connectivity(ehandle, dmbp->hlevel, connp);MBERRNM(merr); 2752417220eSVijay Mahadevan merr = dmbc->mbiface->get_connectivity(&children[0], children.size(), connc);MBERRNM(merr); 276a86ed394SVijay Mahadevan 277a86ed394SVijay Mahadevan std::vector<PetscReal> natparam(3*connc.size(), 0.0); 278a86ed394SVijay Mahadevan pcoords.resize(connp.size() * 3); 279a86ed394SVijay Mahadevan ccoords.resize(connc.size() * 3); 280a86ed394SVijay Mahadevan values_phi.resize(connp.size()*connc.size()); 281a86ed394SVijay Mahadevan /* Get coordinates for connectivity entities in canonical order for both coarse and finer levels */ 282a86ed394SVijay Mahadevan merr = dmbp->hierarchy->get_coordinates(&connp[0], connp.size(), dmbp->hlevel, &pcoords[0]);MBERRNM(merr); 283a86ed394SVijay Mahadevan merr = dmbc->hierarchy->get_coordinates(&connc[0], connc.size(), dmbc->hlevel, &ccoords[0]);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 } 292a86ed394SVijay Mahadevan } 293a86ed394SVijay Mahadevan else { 294a86ed394SVijay Mahadevan factor = std::pow(2.0 /*degree_P_for_refinement*/, (dmbc->hlevel - dmbp->hlevel) * dmbp->dim * 1.0); 295a86ed394SVijay Mahadevan } 296a86ed394SVijay Mahadevan 297755f3dfbSVijay Mahadevan /* TODO: Decipher the correct non-zero pattern. There is still some issue with onz allocation */ 2989566063dSJacob Faibussowitsch PetscCall(MatSetOption(*interpl, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 299b117cd09SVijay Mahadevan 300e882eb38SVijay Mahadevan /* Loop through the remaining vertices. These vertices appear only on the current refined_level. */ 3012417220eSVijay Mahadevan for (moab::Range::iterator iter = dmbp->elocal->begin(); iter != dmbp->elocal->end(); iter++) { 302b117cd09SVijay Mahadevan 303b117cd09SVijay Mahadevan const moab::EntityHandle ehandle = *iter; 304b117cd09SVijay Mahadevan 305b117cd09SVijay Mahadevan /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */ 306a86ed394SVijay Mahadevan children.clear(); 307a86ed394SVijay Mahadevan connc.clear(); 308755f3dfbSVijay Mahadevan merr = dmbp->hierarchy->parent_to_child(ehandle, dmbp->hlevel, dmbc->hlevel, children);MBERRNM(merr); 309b117cd09SVijay Mahadevan 310b117cd09SVijay Mahadevan /* Get connectivity and coordinates of the parent vertices */ 311755f3dfbSVijay Mahadevan merr = dmbp->hierarchy->get_connectivity(ehandle, dmbp->hlevel, connp);MBERRNM(merr); 3122417220eSVijay Mahadevan merr = dmbc->mbiface->get_connectivity(&children[0], children.size(), connc);MBERRNM(merr); 313b117cd09SVijay Mahadevan 314a86ed394SVijay Mahadevan pcoords.resize(connp.size() * 3); 315a86ed394SVijay Mahadevan ccoords.resize(connc.size() * 3); 316b117cd09SVijay Mahadevan /* Get coordinates for connectivity entities in canonical order for both coarse and finer levels */ 317755f3dfbSVijay Mahadevan merr = dmbp->hierarchy->get_coordinates(&connp[0], connp.size(), dmbp->hlevel, &pcoords[0]);MBERRNM(merr); 318755f3dfbSVijay Mahadevan merr = dmbc->hierarchy->get_coordinates(&connc[0], connc.size(), dmbc->hlevel, &ccoords[0]);MBERRNM(merr); 319b117cd09SVijay Mahadevan 320b117cd09SVijay Mahadevan std::vector<int> dofsp(connp.size()), dofsc(connc.size()); 321b117cd09SVijay Mahadevan /* TODO: specific to scalar system - use GetDofs */ 3229566063dSJacob Faibussowitsch PetscCall(DMMoabGetDofsBlocked(dmp, connp.size(), &connp[0], &dofsp[0])); 3239566063dSJacob Faibussowitsch PetscCall(DMMoabGetDofsBlocked(dmc, connc.size(), &connc[0], &dofsc[0])); 324b117cd09SVijay Mahadevan 325a86ed394SVijay Mahadevan /* Compute the actual interpolation weights when projecting solution/residual between levels */ 326a86ed394SVijay Mahadevan if (use_consistent_bases) { 327a86ed394SVijay Mahadevan 328a86ed394SVijay Mahadevan /* Use the cached values of natural parameteric coordinates and basis pre-evaluated. 329a86ed394SVijay Mahadevan We are making an assumption here that UMR used in GMG to generate the hierarchy uses 330a86ed394SVijay Mahadevan the same template for all elements; This will fail for mixed element meshes (TRI/QUAD). 331a86ed394SVijay Mahadevan 3322417220eSVijay Mahadevan TODO: Fix the above assumption by caching data for families (especially for Tets and mixed meshes) 333a86ed394SVijay Mahadevan */ 334a86ed394SVijay Mahadevan 335a86ed394SVijay Mahadevan /* Set values: For each DOF in coarse grid cell, set the contribution or PHI evaluated at each fine grid DOF point */ 336a86ed394SVijay Mahadevan for (unsigned tc = 0; tc < connc.size(); tc++) { 337a86ed394SVijay Mahadevan /* TODO: Check if we should be using INSERT_VALUES instead */ 3389566063dSJacob Faibussowitsch PetscCall(MatSetValues(*interpl, 1, &dofsc[tc], connp.size(), &dofsp[0], &values_phi[connp.size()*tc], ADD_VALUES)); 339a86ed394SVijay Mahadevan } 340a86ed394SVijay Mahadevan } 341a86ed394SVijay Mahadevan else { 342e882eb38SVijay Mahadevan /* Compute the interpolation weights by determining distance of 1-ring 343755f3dfbSVijay Mahadevan neighbor vertices from current vertex 344755f3dfbSVijay Mahadevan 345a86ed394SVijay Mahadevan This should be used only when FEM basis is not used for the discretization. 346a86ed394SVijay Mahadevan Else, the consistent interface to compute the basis function for interpolation 347a86ed394SVijay Mahadevan between the levels should be evaluated correctly to preserve convergence of GMG. 3482417220eSVijay Mahadevan Shephard's basis will be terrible for any unsmooth problems. 349755f3dfbSVijay Mahadevan */ 350755f3dfbSVijay Mahadevan values_phi.resize(connp.size()); 351755f3dfbSVijay Mahadevan for (unsigned tc = 0; tc < connc.size(); tc++) { 352755f3dfbSVijay Mahadevan 353a86ed394SVijay Mahadevan PetscReal normsum = 0.0; 354755f3dfbSVijay Mahadevan for (unsigned tp = 0; tp < connp.size(); tp++) { 355755f3dfbSVijay Mahadevan values_phi[tp] = 0.0; 356e882eb38SVijay Mahadevan for (unsigned k = 0; k < 3; k++) 357755f3dfbSVijay Mahadevan values_phi[tp] += std::pow(pcoords[tp * 3 + k] - ccoords[k + tc * 3], dim); 358755f3dfbSVijay Mahadevan if (values_phi[tp] < 1e-12) { 359755f3dfbSVijay Mahadevan values_phi[tp] = 1e12; 360b117cd09SVijay Mahadevan } 361b117cd09SVijay Mahadevan else { 362755f3dfbSVijay Mahadevan //values_phi[tp] = std::pow(values_phi[tp], -1.0/dim); 363755f3dfbSVijay Mahadevan values_phi[tp] = std::pow(values_phi[tp], -1.0); 364755f3dfbSVijay Mahadevan normsum += values_phi[tp]; 365b117cd09SVijay Mahadevan } 366b117cd09SVijay Mahadevan } 367755f3dfbSVijay Mahadevan for (unsigned tp = 0; tp < connp.size(); tp++) { 368755f3dfbSVijay Mahadevan if (values_phi[tp] > 1e11) 369755f3dfbSVijay Mahadevan values_phi[tp] = factor * 0.5 / connp.size(); 370b117cd09SVijay Mahadevan else 371755f3dfbSVijay Mahadevan values_phi[tp] = factor * values_phi[tp] * 0.5 / (connp.size() * normsum); 372b117cd09SVijay Mahadevan } 3739566063dSJacob Faibussowitsch PetscCall(MatSetValues(*interpl, 1, &dofsc[tc], connp.size(), &dofsp[0], &values_phi[0], ADD_VALUES)); 374b117cd09SVijay Mahadevan } 375a86ed394SVijay Mahadevan } 376b117cd09SVijay Mahadevan } 377304006b3SVijay Mahadevan if (vec) *vec = NULL; 3789566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*interpl, MAT_FINAL_ASSEMBLY)); 3799566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*interpl, MAT_FINAL_ASSEMBLY)); 380b117cd09SVijay Mahadevan PetscFunctionReturn(0); 381b117cd09SVijay Mahadevan } 382b117cd09SVijay Mahadevan 383cab5ea25SPierre Jolivet /*@C 384b117cd09SVijay Mahadevan DMCreateInjection_Moab - Generate a multi-level uniform refinement hierarchy 385b117cd09SVijay Mahadevan by succesively refining a coarse mesh, already defined in the DM object 386b117cd09SVijay Mahadevan provided by the user. 387b117cd09SVijay Mahadevan 388d083f849SBarry Smith Collective 389b117cd09SVijay Mahadevan 390b117cd09SVijay Mahadevan Input Parameter: 391b117cd09SVijay Mahadevan . dmb - The DMMoab object 392b117cd09SVijay Mahadevan 393d8d19677SJose E. Roman Output Parameters: 394a2b725a8SWilliam Gropp + nlevels - The number of levels of refinement needed to generate the hierarchy 395a2b725a8SWilliam Gropp - ldegrees - The degree of refinement at each level in the hierarchy 396b117cd09SVijay Mahadevan 397b117cd09SVijay Mahadevan Level: beginner 398b117cd09SVijay Mahadevan 399b117cd09SVijay Mahadevan @*/ 400304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMCreateInjection_Moab(DM dm1, DM dm2, VecScatter* ctx) 401b117cd09SVijay Mahadevan { 402e882eb38SVijay Mahadevan //DM_Moab *dmmoab; 403b117cd09SVijay Mahadevan 404b117cd09SVijay Mahadevan PetscFunctionBegin; 405b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm1, DM_CLASSID, 1); 406b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm2, DM_CLASSID, 2); 407e882eb38SVijay Mahadevan //dmmoab = (DM_Moab*)(dm1)->data; 408b117cd09SVijay Mahadevan 409b117cd09SVijay Mahadevan PetscPrintf(PETSC_COMM_WORLD, "[DMCreateInjection_Moab] :: Placeholder\n"); 410b117cd09SVijay Mahadevan PetscFunctionReturn(0); 411b117cd09SVijay Mahadevan } 412b117cd09SVijay Mahadevan 413470880abSPatrick Sanan static PetscErrorCode DMMoab_UMR_Private(DM dm, MPI_Comm comm, PetscBool refine, DM *dmref) 414b117cd09SVijay Mahadevan { 415e882eb38SVijay Mahadevan PetscInt i, dim; 416b117cd09SVijay Mahadevan DM dm2; 417b117cd09SVijay Mahadevan moab::ErrorCode merr; 418b117cd09SVijay Mahadevan DM_Moab *dmb = (DM_Moab*)dm->data, *dd2; 419b117cd09SVijay Mahadevan 420b117cd09SVijay Mahadevan PetscFunctionBegin; 421b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 422b117cd09SVijay Mahadevan PetscValidPointer(dmref, 4); 423b117cd09SVijay Mahadevan 424b117cd09SVijay Mahadevan if ((dmb->hlevel == dmb->nhlevels && refine) || (dmb->hlevel == 0 && !refine)) { 425*63a3b9bcSJacob Faibussowitsch if (dmb->hlevel + 1 > dmb->nhlevels && refine) 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); 426*63a3b9bcSJacob Faibussowitsch if (dmb->hlevel - 1 < 0 && !refine) PetscInfo(NULL, "Invalid multigrid coarsen hierarchy level specified (%" PetscInt_FMT "). Creating a NULL object.\n", dmb->hlevel - 1); 427e882eb38SVijay Mahadevan *dmref = PETSC_NULL; 428b117cd09SVijay Mahadevan PetscFunctionReturn(0); 429b117cd09SVijay Mahadevan } 430b117cd09SVijay Mahadevan 4319566063dSJacob Faibussowitsch PetscCall(DMMoabCreate(PetscObjectComm((PetscObject)dm), &dm2)); 432b117cd09SVijay Mahadevan dd2 = (DM_Moab*)dm2->data; 433b117cd09SVijay Mahadevan 434b117cd09SVijay Mahadevan dd2->mbiface = dmb->mbiface; 4359daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 436b117cd09SVijay Mahadevan dd2->pcomm = dmb->pcomm; 4379daf19fdSVijay Mahadevan #endif 438b117cd09SVijay Mahadevan dd2->icreatedinstance = PETSC_FALSE; 43964e1c140SVijay Mahadevan dd2->nghostrings = dmb->nghostrings; 440b117cd09SVijay Mahadevan 441b117cd09SVijay Mahadevan /* set the new level based on refinement/coarsening */ 442b117cd09SVijay Mahadevan if (refine) { 443b117cd09SVijay Mahadevan dd2->hlevel = dmb->hlevel + 1; 444b117cd09SVijay Mahadevan } 445b117cd09SVijay Mahadevan else { 446b117cd09SVijay Mahadevan dd2->hlevel = dmb->hlevel - 1; 447b117cd09SVijay Mahadevan } 448b117cd09SVijay Mahadevan 449b117cd09SVijay Mahadevan /* Copy the multilevel hierarchy pointers in MOAB */ 450b117cd09SVijay Mahadevan dd2->hierarchy = dmb->hierarchy; 451b117cd09SVijay Mahadevan dd2->nhlevels = dmb->nhlevels; 4529566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd2->nhlevels + 1, &dd2->hsets)); 453b117cd09SVijay Mahadevan for (i = 0; i <= dd2->nhlevels; i++) { 454b117cd09SVijay Mahadevan dd2->hsets[i] = dmb->hsets[i]; 455b117cd09SVijay Mahadevan } 456b117cd09SVijay Mahadevan dd2->fileset = dd2->hsets[dd2->hlevel]; 457b117cd09SVijay Mahadevan 458b117cd09SVijay Mahadevan /* do the remaining initializations for DMMoab */ 459b117cd09SVijay Mahadevan dd2->bs = dmb->bs; 460b117cd09SVijay Mahadevan dd2->numFields = dmb->numFields; 461b117cd09SVijay Mahadevan dd2->rw_dbglevel = dmb->rw_dbglevel; 462b117cd09SVijay Mahadevan dd2->partition_by_rank = dmb->partition_by_rank; 4639566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(dd2->extra_read_options, dmb->extra_read_options)); 4649566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(dd2->extra_write_options, dmb->extra_write_options)); 465b117cd09SVijay Mahadevan dd2->read_mode = dmb->read_mode; 466b117cd09SVijay Mahadevan dd2->write_mode = dmb->write_mode; 467b117cd09SVijay Mahadevan 468b117cd09SVijay Mahadevan /* set global ID tag handle */ 4699566063dSJacob Faibussowitsch PetscCall(DMMoabSetLocalToGlobalTag(dm2, dmb->ltog_tag)); 470b117cd09SVijay Mahadevan 471b117cd09SVijay Mahadevan merr = dd2->mbiface->tag_get_handle(MATERIAL_SET_TAG_NAME, dd2->material_tag);MBERRNM(merr); 472b117cd09SVijay Mahadevan 4739566063dSJacob Faibussowitsch PetscCall(DMSetOptionsPrefix(dm2, ((PetscObject)dm)->prefix)); 4749566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 4759566063dSJacob Faibussowitsch PetscCall(DMSetDimension(dm2, dim)); 476b117cd09SVijay Mahadevan 477b117cd09SVijay Mahadevan /* allow overloaded (user replaced) operations to be inherited by refinement clones */ 478b117cd09SVijay Mahadevan dm2->ops->creatematrix = dm->ops->creatematrix; 479b117cd09SVijay Mahadevan 480b117cd09SVijay Mahadevan /* copy fill information if given */ 4819566063dSJacob Faibussowitsch PetscCall(DMMoabSetBlockFills(dm2, dmb->dfill, dmb->ofill)); 482b117cd09SVijay Mahadevan 483b117cd09SVijay Mahadevan /* copy vector type information */ 4849566063dSJacob Faibussowitsch PetscCall(DMSetMatType(dm2, dm->mattype)); 4859566063dSJacob Faibussowitsch PetscCall(DMSetVecType(dm2, dm->vectype)); 4863f1c6e43SVijay Mahadevan dd2->numFields = dmb->numFields; 4873f1c6e43SVijay Mahadevan if (dmb->numFields) { 4889566063dSJacob Faibussowitsch PetscCall(DMMoabSetFieldNames(dm2, dmb->numFields, dmb->fieldNames)); 4893f1c6e43SVijay Mahadevan } 490b117cd09SVijay Mahadevan 4919566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm2)); 492b117cd09SVijay Mahadevan 493b117cd09SVijay Mahadevan /* recreate Dof numbering for the refined DM and make sure the distribution is correctly populated */ 4949566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm2)); 495b117cd09SVijay Mahadevan 496b117cd09SVijay Mahadevan *dmref = dm2; 497b117cd09SVijay Mahadevan PetscFunctionReturn(0); 498b117cd09SVijay Mahadevan } 499b117cd09SVijay Mahadevan 500cab5ea25SPierre Jolivet /*@C 501b117cd09SVijay Mahadevan DMRefine_Moab - Generate a multi-level uniform refinement hierarchy 502b117cd09SVijay Mahadevan by succesively refining a coarse mesh, already defined in the DM object 503b117cd09SVijay Mahadevan provided by the user. 504b117cd09SVijay Mahadevan 505d083f849SBarry Smith Collective on dm 506b117cd09SVijay Mahadevan 507d8d19677SJose E. Roman Input Parameters: 508e882eb38SVijay Mahadevan + dm - The DMMoab object 509e882eb38SVijay Mahadevan - comm - the communicator to contain the new DM object (or MPI_COMM_NULL) 510b117cd09SVijay Mahadevan 511b117cd09SVijay Mahadevan Output Parameter: 512e882eb38SVijay Mahadevan . dmf - the refined DM, or NULL 513b117cd09SVijay Mahadevan 514e882eb38SVijay Mahadevan Note: If no refinement was done, the return value is NULL 515e882eb38SVijay Mahadevan 516e882eb38SVijay Mahadevan Level: developer 517b117cd09SVijay Mahadevan 518b117cd09SVijay Mahadevan @*/ 519304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMRefine_Moab(DM dm, MPI_Comm comm, DM* dmf) 520b117cd09SVijay Mahadevan { 521b117cd09SVijay Mahadevan PetscFunctionBegin; 522b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 523b117cd09SVijay Mahadevan 5249566063dSJacob Faibussowitsch PetscCall(DMMoab_UMR_Private(dm, comm, PETSC_TRUE, dmf)); 525b117cd09SVijay Mahadevan PetscFunctionReturn(0); 526b117cd09SVijay Mahadevan } 527b117cd09SVijay Mahadevan 528cab5ea25SPierre Jolivet /*@C 529b117cd09SVijay Mahadevan DMCoarsen_Moab - Generate a multi-level uniform refinement hierarchy 530b117cd09SVijay Mahadevan by succesively refining a coarse mesh, already defined in the DM object 531b117cd09SVijay Mahadevan provided by the user. 532b117cd09SVijay Mahadevan 533d083f849SBarry Smith Collective on dm 534b117cd09SVijay Mahadevan 535d8d19677SJose E. Roman Input Parameters: 536a2b725a8SWilliam Gropp + dm - The DMMoab object 537e882eb38SVijay Mahadevan - comm - the communicator to contain the new DM object (or MPI_COMM_NULL) 538b117cd09SVijay Mahadevan 539b117cd09SVijay Mahadevan Output Parameter: 540e882eb38SVijay Mahadevan . dmf - the coarsened DM, or NULL 541b117cd09SVijay Mahadevan 542e882eb38SVijay Mahadevan Note: If no coarsening was done, the return value is NULL 543b117cd09SVijay Mahadevan 544e882eb38SVijay Mahadevan Level: developer 545e882eb38SVijay Mahadevan 546b117cd09SVijay Mahadevan @*/ 547304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMCoarsen_Moab(DM dm, MPI_Comm comm, DM* dmc) 548b117cd09SVijay Mahadevan { 549b117cd09SVijay Mahadevan PetscFunctionBegin; 550b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5519566063dSJacob Faibussowitsch PetscCall(DMMoab_UMR_Private(dm, comm, PETSC_FALSE, dmc)); 552b117cd09SVijay Mahadevan PetscFunctionReturn(0); 553b117cd09SVijay Mahadevan } 554