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 6b117cd09SVijay Mahadevan #undef __FUNCT__ 7b117cd09SVijay Mahadevan #define __FUNCT__ "DMMoabGenerateHierarchy" 8b117cd09SVijay Mahadevan /*@ 9b117cd09SVijay Mahadevan DMMoabGenerateHierarchy - Generate a multi-level uniform refinement hierarchy 10b117cd09SVijay Mahadevan by succesively refining a coarse mesh, already defined in the DM object 11b117cd09SVijay Mahadevan provided by the user. 12b117cd09SVijay Mahadevan 13b117cd09SVijay Mahadevan Collective on MPI_Comm 14b117cd09SVijay Mahadevan 15b117cd09SVijay Mahadevan Input Parameter: 16b117cd09SVijay Mahadevan + dmb - The DMMoab object 17b117cd09SVijay Mahadevan 18b117cd09SVijay Mahadevan Output Parameter: 19b117cd09SVijay Mahadevan + nlevels - The number of levels of refinement needed to generate the hierarchy 20b117cd09SVijay Mahadevan . ldegrees - The degree of refinement at each level in the hierarchy 21b117cd09SVijay Mahadevan 22b117cd09SVijay Mahadevan Level: beginner 23b117cd09SVijay Mahadevan 24b117cd09SVijay Mahadevan .keywords: DMMoab, create, refinement 25b117cd09SVijay Mahadevan @*/ 26b117cd09SVijay Mahadevan PetscErrorCode DMMoabGenerateHierarchy(DM dm, PetscInt nlevels, PetscInt *ldegrees) 27b117cd09SVijay Mahadevan { 28b117cd09SVijay Mahadevan DM_Moab *dmmoab; 29b117cd09SVijay Mahadevan PetscErrorCode ierr; 30b117cd09SVijay Mahadevan moab::ErrorCode merr; 312417220eSVijay Mahadevan PetscInt *pdegrees, ilevel; 32e882eb38SVijay Mahadevan std::vector<moab::EntityHandle> hsets; 33b117cd09SVijay Mahadevan 34b117cd09SVijay Mahadevan PetscFunctionBegin; 35b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 36b117cd09SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 37b117cd09SVijay Mahadevan 38b117cd09SVijay Mahadevan if (!ldegrees) { 39b117cd09SVijay Mahadevan ierr = PetscMalloc1(nlevels, &pdegrees);CHKERRQ(ierr); 402417220eSVijay Mahadevan for (ilevel = 0; ilevel < nlevels; ilevel++) pdegrees[ilevel] = 2; /* default = Degree 2 refinement */ 41b117cd09SVijay Mahadevan } 42b117cd09SVijay Mahadevan else pdegrees = ldegrees; 43b117cd09SVijay Mahadevan 44b117cd09SVijay Mahadevan /* initialize set level refinement data for hierarchy */ 45b117cd09SVijay Mahadevan dmmoab->nhlevels = nlevels; 46b117cd09SVijay Mahadevan 47b117cd09SVijay Mahadevan /* Instantiate the nested refinement class */ 489daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 493f1c6e43SVijay Mahadevan dmmoab->hierarchy = new moab::NestedRefine(dynamic_cast<moab::Core*>(dmmoab->mbiface), dmmoab->pcomm, dmmoab->fileset); 509daf19fdSVijay Mahadevan #else 519daf19fdSVijay Mahadevan dmmoab->hierarchy = new moab::NestedRefine(dynamic_cast<moab::Core*>(dmmoab->mbiface), NULL, dmmoab->fileset); 529daf19fdSVijay Mahadevan #endif 53b117cd09SVijay Mahadevan 54e882eb38SVijay Mahadevan ierr = PetscMalloc1(nlevels + 1, &dmmoab->hsets);CHKERRQ(ierr); 55e882eb38SVijay Mahadevan 56b117cd09SVijay Mahadevan /* generate the mesh hierarchy */ 579c368985SVijay Mahadevan merr = dmmoab->hierarchy->generate_mesh_hierarchy(nlevels, pdegrees, hsets, false);MBERRNM(merr); 58e882eb38SVijay Mahadevan 599daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 60755f3dfbSVijay Mahadevan if (dmmoab->pcomm->size() > 1) { 6164e1c140SVijay Mahadevan merr = dmmoab->hierarchy->exchange_ghosts(hsets, dmmoab->nghostrings);MBERRNM(merr); 62755f3dfbSVijay Mahadevan } 639daf19fdSVijay Mahadevan #endif 6449d66b22SVijay Mahadevan 6564e1c140SVijay Mahadevan /* copy the mesh sets for nested refinement hierarchy */ 66*e92d1c7cSVijay Mahadevan dmmoab->hsets[0] = hsets[0]; 67*e92d1c7cSVijay Mahadevan for (ilevel = 1; ilevel <= nlevels; ilevel++) 682417220eSVijay Mahadevan { 692417220eSVijay Mahadevan dmmoab->hsets[ilevel] = hsets[ilevel]; 70*e92d1c7cSVijay Mahadevan 719c368985SVijay Mahadevan #ifdef MOAB_HAVE_MPI 72*e92d1c7cSVijay Mahadevan merr = dmmoab->pcomm->assign_global_ids(hsets[ilevel], dmmoab->dim, 0, false, true, false);MBERRNM(merr); 739c368985SVijay Mahadevan #endif 74*e92d1c7cSVijay Mahadevan 75*e92d1c7cSVijay Mahadevan /* Update material and other geometric tags from parent to child sets */ 76*e92d1c7cSVijay Mahadevan merr = dmmoab->hierarchy->update_special_tags(ilevel, hsets[ilevel]);MBERRNM(merr); 772417220eSVijay Mahadevan } 7864e1c140SVijay Mahadevan 7964e1c140SVijay Mahadevan hsets.clear(); 80b117cd09SVijay Mahadevan if (!ldegrees) { 81b117cd09SVijay Mahadevan ierr = PetscFree(pdegrees);CHKERRQ(ierr); 82b117cd09SVijay Mahadevan } 83b117cd09SVijay Mahadevan PetscFunctionReturn(0); 84b117cd09SVijay Mahadevan } 85b117cd09SVijay Mahadevan 86b117cd09SVijay Mahadevan #undef __FUNCT__ 87b117cd09SVijay Mahadevan #define __FUNCT__ "DMRefineHierarchy_Moab" 88b117cd09SVijay Mahadevan /*@ 89e882eb38SVijay Mahadevan DMRefineHierarchy_Moab - Generate a multi-level DM hierarchy 90e882eb38SVijay Mahadevan by succesively refining a coarse mesh. 91b117cd09SVijay Mahadevan 92b117cd09SVijay Mahadevan Collective on MPI_Comm 93b117cd09SVijay Mahadevan 94b117cd09SVijay Mahadevan Input Parameter: 95e882eb38SVijay Mahadevan + dm - The DMMoab object 96b117cd09SVijay Mahadevan 97b117cd09SVijay Mahadevan Output Parameter: 98e882eb38SVijay Mahadevan + nlevels - The number of levels of refinement needed to generate the hierarchy 99e882eb38SVijay Mahadevan . dmf - The DM objects after successive refinement of the hierarchy 100b117cd09SVijay Mahadevan 101b117cd09SVijay Mahadevan Level: beginner 102b117cd09SVijay Mahadevan 103e882eb38SVijay Mahadevan .keywords: DMMoab, generate, refinement 104b117cd09SVijay Mahadevan @*/ 105304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMRefineHierarchy_Moab(DM dm, PetscInt nlevels, DM dmf[]) 106b117cd09SVijay Mahadevan { 107b117cd09SVijay Mahadevan PetscErrorCode ierr; 108b117cd09SVijay Mahadevan PetscInt i; 109b117cd09SVijay Mahadevan 110b117cd09SVijay Mahadevan PetscFunctionBegin; 111b117cd09SVijay Mahadevan 112e882eb38SVijay Mahadevan ierr = DMRefine(dm, PetscObjectComm((PetscObject)dm), &dmf[0]);CHKERRQ(ierr); 113e882eb38SVijay Mahadevan for (i = 1; i < nlevels; i++) { 114e882eb38SVijay Mahadevan ierr = DMRefine(dmf[i - 1], PetscObjectComm((PetscObject)dm), &dmf[i]);CHKERRQ(ierr); 115b117cd09SVijay Mahadevan } 116b117cd09SVijay Mahadevan PetscFunctionReturn(0); 117b117cd09SVijay Mahadevan } 118b117cd09SVijay Mahadevan 119b117cd09SVijay Mahadevan #undef __FUNCT__ 120b117cd09SVijay Mahadevan #define __FUNCT__ "DMCoarsenHierarchy_Moab" 121b117cd09SVijay Mahadevan /*@ 122e882eb38SVijay Mahadevan DMCoarsenHierarchy_Moab - Generate a multi-level DM hierarchy 123e882eb38SVijay Mahadevan by succesively coarsening a refined mesh. 124b117cd09SVijay Mahadevan 125b117cd09SVijay Mahadevan Collective on MPI_Comm 126b117cd09SVijay Mahadevan 127b117cd09SVijay Mahadevan Input Parameter: 128e882eb38SVijay Mahadevan + dm - The DMMoab object 129b117cd09SVijay Mahadevan 130b117cd09SVijay Mahadevan Output Parameter: 131e882eb38SVijay Mahadevan + nlevels - The number of levels of refinement needed to generate the hierarchy 132e882eb38SVijay Mahadevan . dmc - The DM objects after successive coarsening of the hierarchy 133b117cd09SVijay Mahadevan 134b117cd09SVijay Mahadevan Level: beginner 135b117cd09SVijay Mahadevan 136e882eb38SVijay Mahadevan .keywords: DMMoab, generate, coarsening 137b117cd09SVijay Mahadevan @*/ 138304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMCoarsenHierarchy_Moab(DM dm, PetscInt nlevels, DM dmc[]) 139b117cd09SVijay Mahadevan { 140b117cd09SVijay Mahadevan PetscErrorCode ierr; 141b117cd09SVijay Mahadevan PetscInt i; 142b117cd09SVijay Mahadevan 143b117cd09SVijay Mahadevan PetscFunctionBegin; 144b117cd09SVijay Mahadevan 145e882eb38SVijay Mahadevan ierr = DMCoarsen(dm, PetscObjectComm((PetscObject)dm), &dmc[0]);CHKERRQ(ierr); 146e882eb38SVijay Mahadevan for (i = 1; i < nlevels; i++) { 147e882eb38SVijay Mahadevan ierr = DMCoarsen(dmc[i - 1], PetscObjectComm((PetscObject)dm), &dmc[i]);CHKERRQ(ierr); 148b117cd09SVijay Mahadevan } 149b117cd09SVijay Mahadevan PetscFunctionReturn(0); 150b117cd09SVijay Mahadevan } 151b117cd09SVijay Mahadevan 1522417220eSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoab_Compute_NNZ_From_Connectivity(DM, PetscInt*, PetscInt*, PetscInt*, PetscInt*, PetscBool); 153b117cd09SVijay Mahadevan 154b117cd09SVijay Mahadevan #undef __FUNCT__ 155b117cd09SVijay Mahadevan #define __FUNCT__ "DMCreateInterpolation_Moab" 156b117cd09SVijay Mahadevan /*@ 157e882eb38SVijay Mahadevan DMCreateInterpolation_Moab - Generate the interpolation operators to transform 158e882eb38SVijay Mahadevan operators (matrices, vectors) from parent level to child level as defined by 159e882eb38SVijay Mahadevan the DM inputs provided by the user. 160b117cd09SVijay Mahadevan 161b117cd09SVijay Mahadevan Collective on MPI_Comm 162b117cd09SVijay Mahadevan 163b117cd09SVijay Mahadevan Input Parameter: 164e882eb38SVijay Mahadevan + dm1 - The DMMoab object 165e882eb38SVijay Mahadevan - dm2 - the second, finer DMMoab object 166b117cd09SVijay Mahadevan 167b117cd09SVijay Mahadevan Output Parameter: 168e882eb38SVijay Mahadevan + interpl - The interpolation operator for transferring data between the levels 169e882eb38SVijay Mahadevan - vec - The scaling vector (optional) 170b117cd09SVijay Mahadevan 171e882eb38SVijay Mahadevan Level: developer 172b117cd09SVijay Mahadevan 173b117cd09SVijay Mahadevan .keywords: DMMoab, create, refinement 174b117cd09SVijay Mahadevan @*/ 175304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMCreateInterpolation_Moab(DM dmp, DM dmc, Mat* interpl, Vec* vec) 176b117cd09SVijay Mahadevan { 177755f3dfbSVijay Mahadevan DM_Moab *dmbp, *dmbc; 178b117cd09SVijay Mahadevan PetscErrorCode ierr; 179b117cd09SVijay Mahadevan moab::ErrorCode merr; 180e882eb38SVijay Mahadevan PetscInt dim; 1813f1c6e43SVijay Mahadevan PetscReal factor; 182ce27a4eeSVijay Mahadevan PetscInt innz, *nnz, ionz, *onz; 183755f3dfbSVijay Mahadevan PetscInt nlsizp, nlsizc, nlghsizp, ngsizp, ngsizc; 184a86ed394SVijay Mahadevan const PetscBool use_consistent_bases=PETSC_TRUE; 185b117cd09SVijay Mahadevan 186b117cd09SVijay Mahadevan PetscFunctionBegin; 187755f3dfbSVijay Mahadevan PetscValidHeaderSpecific(dmp, DM_CLASSID, 1); 188755f3dfbSVijay Mahadevan PetscValidHeaderSpecific(dmc, DM_CLASSID, 2); 189755f3dfbSVijay Mahadevan dmbp = (DM_Moab*)(dmp)->data; 190755f3dfbSVijay Mahadevan dmbc = (DM_Moab*)(dmc)->data; 191755f3dfbSVijay Mahadevan nlsizp = dmbp->nloc;// *dmb1->numFields; 192755f3dfbSVijay Mahadevan nlsizc = dmbc->nloc;// *dmb2->numFields; 193755f3dfbSVijay Mahadevan ngsizp = dmbp->n; // *dmb1->numFields; 194755f3dfbSVijay Mahadevan ngsizc = dmbc->n; // *dmb2->numFields; 195755f3dfbSVijay Mahadevan nlghsizp = (dmbp->nloc + dmbp->nghost); // *dmb1->numFields; 196b117cd09SVijay Mahadevan 1972417220eSVijay Mahadevan // Columns = Parent DoFs ; Rows = Child DoFs 198755f3dfbSVijay Mahadevan // Interpolation matrix: \sum_{i=1}^P Owned(Child) * (Owned(Parent) + Ghosted(Parent)) 199755f3dfbSVijay Mahadevan // Size: nlsizc * nlghsizp 200755f3dfbSVijay Mahadevan PetscInfo4(NULL, "Creating interpolation matrix %D X %D to apply transformation between levels %D -> %D.\n", ngsizc, nlghsizp, dmbp->hlevel, dmbc->hlevel); 201b117cd09SVijay Mahadevan 202a86ed394SVijay Mahadevan ierr = DMGetDimension(dmp, &dim);CHKERRQ(ierr); 203a86ed394SVijay Mahadevan 204941e0cffSVijay Mahadevan /* allocate the nnz, onz arrays based on block size and local nodes */ 205755f3dfbSVijay Mahadevan ierr = PetscCalloc2(nlsizc, &nnz, nlsizc, &onz);CHKERRQ(ierr); 206941e0cffSVijay Mahadevan 207941e0cffSVijay Mahadevan /* Loop through the local elements and compute the relation between the current parent and the refined_level. */ 2082417220eSVijay Mahadevan for (moab::Range::iterator iter = dmbc->vowned->begin(); iter != dmbc->vowned->end(); iter++) { 209941e0cffSVijay Mahadevan 2102417220eSVijay Mahadevan const moab::EntityHandle vhandle = *iter; 2112417220eSVijay Mahadevan /* define local variables */ 2122417220eSVijay Mahadevan moab::EntityHandle parent; 2132417220eSVijay Mahadevan std::vector<moab::EntityHandle> adjs; 2142417220eSVijay Mahadevan moab::Range found; 2152417220eSVijay Mahadevan 2162417220eSVijay Mahadevan /* store the vertex DoF number */ 2172417220eSVijay Mahadevan const int ldof = dmbc->lidmap[vhandle - dmbc->seqstart]; 2182417220eSVijay Mahadevan 2192417220eSVijay Mahadevan /* Get adjacency information for current vertex - i.e., all elements of dimension (dim) that connects 2202417220eSVijay Mahadevan to the current vertex. We can then decipher if a vertex is ghosted or not and compute the 2212417220eSVijay Mahadevan non-zero pattern accordingly. */ 2222417220eSVijay Mahadevan merr = dmbc->hierarchy->get_adjacencies(vhandle, dmbc->dim, adjs);MBERRNM(merr); 2232417220eSVijay Mahadevan 2242417220eSVijay Mahadevan /* loop over vertices and update the number of connectivity */ 2252417220eSVijay Mahadevan for (unsigned jter = 0; jter < adjs.size(); jter++) { 2262417220eSVijay Mahadevan 2272417220eSVijay Mahadevan const moab::EntityHandle jhandle = adjs[jter]; 228941e0cffSVijay Mahadevan 229941e0cffSVijay Mahadevan /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */ 2302417220eSVijay Mahadevan merr = dmbc->hierarchy->child_to_parent(jhandle, dmbc->hlevel, dmbp->hlevel, &parent);MBERRNM(merr); 231941e0cffSVijay Mahadevan 2322417220eSVijay Mahadevan /* Get connectivity information in canonical ordering for the local element */ 2332417220eSVijay Mahadevan std::vector<moab::EntityHandle> connp; 2342417220eSVijay Mahadevan merr = dmbp->hierarchy->get_connectivity(parent, dmbp->hlevel, connp);MBERRNM(merr); 235a86ed394SVijay Mahadevan 2362417220eSVijay Mahadevan for (unsigned ic=0; ic < connp.size(); ++ic) { 2372417220eSVijay Mahadevan 2382417220eSVijay Mahadevan /* loop over each element connected to the adjacent vertex and update as needed */ 2392417220eSVijay Mahadevan /* find the truly user-expected layer of ghosted entities to decipher NNZ pattern */ 2402417220eSVijay Mahadevan if (found.find(connp[ic]) != found.end()) continue; /* make sure we don't double count shared vertices */ 2412417220eSVijay Mahadevan if (dmbp->vghost->find(connp[ic]) != dmbp->vghost->end()) onz[ldof]++; /* update out-of-proc onz */ 2422417220eSVijay Mahadevan else nnz[ldof]++; /* else local vertex */ 2432417220eSVijay Mahadevan found.insert(connp[ic]); 2442417220eSVijay Mahadevan } 2452417220eSVijay Mahadevan } 246941e0cffSVijay Mahadevan } 247941e0cffSVijay Mahadevan 2482417220eSVijay Mahadevan for (int i = 0; i < nlsizc; i++) 2492417220eSVijay Mahadevan nnz[i] += 1; /* self count the node */ 250941e0cffSVijay Mahadevan 251941e0cffSVijay Mahadevan ionz = onz[0]; 252941e0cffSVijay Mahadevan innz = nnz[0]; 253755f3dfbSVijay Mahadevan for (int tc = 0; tc < nlsizc; tc++) { 254941e0cffSVijay Mahadevan // check for maximum allowed sparsity = fully dense 255755f3dfbSVijay Mahadevan nnz[tc] = std::min(nlsizp, nnz[tc]); 2562417220eSVijay Mahadevan onz[tc] = std::min(ngsizp - nlsizp, onz[tc]); 2572417220eSVijay Mahadevan 2582417220eSVijay Mahadevan PetscInfo3(NULL, " %d: NNZ = %d, ONZ = %d\n", tc, nnz[tc], onz[tc]); 259941e0cffSVijay Mahadevan 260941e0cffSVijay Mahadevan innz = (innz < nnz[tc] ? nnz[tc] : innz); 261941e0cffSVijay Mahadevan ionz = (ionz < onz[tc] ? onz[tc] : ionz); 262941e0cffSVijay Mahadevan } 26364e1c140SVijay Mahadevan 26464e1c140SVijay Mahadevan /* create interpolation matrix */ 265755f3dfbSVijay Mahadevan ierr = MatCreate(PetscObjectComm((PetscObject)dmc), interpl);CHKERRQ(ierr); 266755f3dfbSVijay Mahadevan ierr = MatSetSizes(*interpl, nlsizc, nlsizp, ngsizc, ngsizp);CHKERRQ(ierr); 26764e1c140SVijay Mahadevan ierr = MatSetType(*interpl, MATAIJ);CHKERRQ(ierr); 26864e1c140SVijay Mahadevan ierr = MatSetFromOptions(*interpl);CHKERRQ(ierr); 26964e1c140SVijay Mahadevan 270941e0cffSVijay Mahadevan ierr = MatSeqAIJSetPreallocation(*interpl, innz, nnz);CHKERRQ(ierr); 271941e0cffSVijay Mahadevan ierr = MatMPIAIJSetPreallocation(*interpl, innz, nnz, ionz, onz);CHKERRQ(ierr); 272941e0cffSVijay Mahadevan 273941e0cffSVijay Mahadevan /* clean up temporary memory */ 274941e0cffSVijay Mahadevan ierr = PetscFree2(nnz, onz);CHKERRQ(ierr); 275b117cd09SVijay Mahadevan 276b117cd09SVijay Mahadevan /* set up internal matrix data-structures */ 277b117cd09SVijay Mahadevan ierr = MatSetUp(*interpl);CHKERRQ(ierr); 278b117cd09SVijay Mahadevan 279a86ed394SVijay Mahadevan /* Define variables for assembly */ 280a86ed394SVijay Mahadevan std::vector<moab::EntityHandle> children; 281a86ed394SVijay Mahadevan std::vector<moab::EntityHandle> connp, connc; 282a86ed394SVijay Mahadevan std::vector<PetscReal> pcoords, ccoords, values_phi; 283a86ed394SVijay Mahadevan 284a86ed394SVijay Mahadevan if (use_consistent_bases) { 2852417220eSVijay Mahadevan const moab::EntityHandle ehandle = dmbp->elocal->front(); 286a86ed394SVijay Mahadevan 287a86ed394SVijay Mahadevan merr = dmbp->hierarchy->parent_to_child(ehandle, dmbp->hlevel, dmbc->hlevel, children);MBERRNM(merr); 288a86ed394SVijay Mahadevan 289a86ed394SVijay Mahadevan /* Get connectivity and coordinates of the parent vertices */ 290a86ed394SVijay Mahadevan merr = dmbp->hierarchy->get_connectivity(ehandle, dmbp->hlevel, connp);MBERRNM(merr); 2912417220eSVijay Mahadevan merr = dmbc->mbiface->get_connectivity(&children[0], children.size(), connc);MBERRNM(merr); 292a86ed394SVijay Mahadevan 293a86ed394SVijay Mahadevan std::vector<PetscReal> natparam(3*connc.size(), 0.0); 294a86ed394SVijay Mahadevan pcoords.resize(connp.size() * 3); 295a86ed394SVijay Mahadevan ccoords.resize(connc.size() * 3); 296a86ed394SVijay Mahadevan values_phi.resize(connp.size()*connc.size()); 297a86ed394SVijay Mahadevan /* Get coordinates for connectivity entities in canonical order for both coarse and finer levels */ 298a86ed394SVijay Mahadevan merr = dmbp->hierarchy->get_coordinates(&connp[0], connp.size(), dmbp->hlevel, &pcoords[0]);MBERRNM(merr); 299a86ed394SVijay Mahadevan merr = dmbc->hierarchy->get_coordinates(&connc[0], connc.size(), dmbc->hlevel, &ccoords[0]);MBERRNM(merr); 300a86ed394SVijay Mahadevan 301a86ed394SVijay Mahadevan /* Set values: For each DOF in coarse grid cell, set the contribution or PHI evaluated at each fine grid DOF point */ 302a86ed394SVijay Mahadevan for (unsigned tc = 0; tc < connc.size(); tc++) { 303a86ed394SVijay Mahadevan const PetscInt offset = tc * 3; 304a86ed394SVijay Mahadevan 305a86ed394SVijay Mahadevan /* Scale ccoords relative to pcoords */ 306a86ed394SVijay Mahadevan ierr = DMMoabPToRMapping(dim, connp.size(), &pcoords[0], &ccoords[offset], &natparam[offset], &values_phi[connp.size()*tc]);CHKERRQ(ierr); 307a86ed394SVijay Mahadevan } 308a86ed394SVijay Mahadevan } 309a86ed394SVijay Mahadevan else { 310a86ed394SVijay Mahadevan factor = std::pow(2.0 /*degree_P_for_refinement*/, (dmbc->hlevel - dmbp->hlevel) * dmbp->dim * 1.0); 311a86ed394SVijay Mahadevan } 312a86ed394SVijay Mahadevan 313755f3dfbSVijay Mahadevan /* TODO: Decipher the correct non-zero pattern. There is still some issue with onz allocation */ 3142417220eSVijay Mahadevan ierr = MatSetOption(*interpl, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE);CHKERRQ(ierr); 315b117cd09SVijay Mahadevan 316e882eb38SVijay Mahadevan /* Loop through the remaining vertices. These vertices appear only on the current refined_level. */ 3172417220eSVijay Mahadevan for (moab::Range::iterator iter = dmbp->elocal->begin(); iter != dmbp->elocal->end(); iter++) { 318b117cd09SVijay Mahadevan 319b117cd09SVijay Mahadevan const moab::EntityHandle ehandle = *iter; 320b117cd09SVijay Mahadevan 321b117cd09SVijay Mahadevan /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */ 322a86ed394SVijay Mahadevan children.clear(); 323a86ed394SVijay Mahadevan connc.clear(); 324755f3dfbSVijay Mahadevan merr = dmbp->hierarchy->parent_to_child(ehandle, dmbp->hlevel, dmbc->hlevel, children);MBERRNM(merr); 325b117cd09SVijay Mahadevan 326b117cd09SVijay Mahadevan /* Get connectivity and coordinates of the parent vertices */ 327755f3dfbSVijay Mahadevan merr = dmbp->hierarchy->get_connectivity(ehandle, dmbp->hlevel, connp);MBERRNM(merr); 3282417220eSVijay Mahadevan merr = dmbc->mbiface->get_connectivity(&children[0], children.size(), connc);MBERRNM(merr); 329b117cd09SVijay Mahadevan 330a86ed394SVijay Mahadevan pcoords.resize(connp.size() * 3); 331a86ed394SVijay Mahadevan ccoords.resize(connc.size() * 3); 332b117cd09SVijay Mahadevan /* Get coordinates for connectivity entities in canonical order for both coarse and finer levels */ 333755f3dfbSVijay Mahadevan merr = dmbp->hierarchy->get_coordinates(&connp[0], connp.size(), dmbp->hlevel, &pcoords[0]);MBERRNM(merr); 334755f3dfbSVijay Mahadevan merr = dmbc->hierarchy->get_coordinates(&connc[0], connc.size(), dmbc->hlevel, &ccoords[0]);MBERRNM(merr); 335b117cd09SVijay Mahadevan 336b117cd09SVijay Mahadevan std::vector<int> dofsp(connp.size()), dofsc(connc.size()); 337b117cd09SVijay Mahadevan /* TODO: specific to scalar system - use GetDofs */ 338755f3dfbSVijay Mahadevan ierr = DMMoabGetDofsBlocked(dmp, connp.size(), &connp[0], &dofsp[0]);CHKERRQ(ierr); 339755f3dfbSVijay Mahadevan ierr = DMMoabGetDofsBlocked(dmc, connc.size(), &connc[0], &dofsc[0]);CHKERRQ(ierr); 340b117cd09SVijay Mahadevan 341a86ed394SVijay Mahadevan /* Compute the actual interpolation weights when projecting solution/residual between levels */ 342a86ed394SVijay Mahadevan if (use_consistent_bases) { 343a86ed394SVijay Mahadevan 344a86ed394SVijay Mahadevan /* Use the cached values of natural parameteric coordinates and basis pre-evaluated. 345a86ed394SVijay Mahadevan We are making an assumption here that UMR used in GMG to generate the hierarchy uses 346a86ed394SVijay Mahadevan the same template for all elements; This will fail for mixed element meshes (TRI/QUAD). 347a86ed394SVijay Mahadevan 3482417220eSVijay Mahadevan TODO: Fix the above assumption by caching data for families (especially for Tets and mixed meshes) 349a86ed394SVijay Mahadevan */ 350a86ed394SVijay Mahadevan 351a86ed394SVijay Mahadevan /* Set values: For each DOF in coarse grid cell, set the contribution or PHI evaluated at each fine grid DOF point */ 352a86ed394SVijay Mahadevan for (unsigned tc = 0; tc < connc.size(); tc++) { 353a86ed394SVijay Mahadevan /* TODO: Check if we should be using INSERT_VALUES instead */ 354a86ed394SVijay Mahadevan ierr = MatSetValues(*interpl, 1, &dofsc[tc], connp.size(), &dofsp[0], &values_phi[connp.size()*tc], ADD_VALUES);CHKERRQ(ierr); 355a86ed394SVijay Mahadevan } 356a86ed394SVijay Mahadevan } 357a86ed394SVijay Mahadevan else { 358e882eb38SVijay Mahadevan /* Compute the interpolation weights by determining distance of 1-ring 359755f3dfbSVijay Mahadevan neighbor vertices from current vertex 360755f3dfbSVijay Mahadevan 361a86ed394SVijay Mahadevan This should be used only when FEM basis is not used for the discretization. 362a86ed394SVijay Mahadevan Else, the consistent interface to compute the basis function for interpolation 363a86ed394SVijay Mahadevan between the levels should be evaluated correctly to preserve convergence of GMG. 3642417220eSVijay Mahadevan Shephard's basis will be terrible for any unsmooth problems. 365755f3dfbSVijay Mahadevan */ 366755f3dfbSVijay Mahadevan values_phi.resize(connp.size()); 367755f3dfbSVijay Mahadevan for (unsigned tc = 0; tc < connc.size(); tc++) { 368755f3dfbSVijay Mahadevan 369a86ed394SVijay Mahadevan PetscReal normsum = 0.0; 370755f3dfbSVijay Mahadevan for (unsigned tp = 0; tp < connp.size(); tp++) { 371755f3dfbSVijay Mahadevan values_phi[tp] = 0.0; 372e882eb38SVijay Mahadevan for (unsigned k = 0; k < 3; k++) 373755f3dfbSVijay Mahadevan values_phi[tp] += std::pow(pcoords[tp * 3 + k] - ccoords[k + tc * 3], dim); 374755f3dfbSVijay Mahadevan if (values_phi[tp] < 1e-12) { 375755f3dfbSVijay Mahadevan values_phi[tp] = 1e12; 376b117cd09SVijay Mahadevan } 377b117cd09SVijay Mahadevan else { 378755f3dfbSVijay Mahadevan //values_phi[tp] = std::pow(values_phi[tp], -1.0/dim); 379755f3dfbSVijay Mahadevan values_phi[tp] = std::pow(values_phi[tp], -1.0); 380755f3dfbSVijay Mahadevan normsum += values_phi[tp]; 381b117cd09SVijay Mahadevan } 382b117cd09SVijay Mahadevan } 383755f3dfbSVijay Mahadevan for (unsigned tp = 0; tp < connp.size(); tp++) { 384755f3dfbSVijay Mahadevan if (values_phi[tp] > 1e11) 385755f3dfbSVijay Mahadevan values_phi[tp] = factor * 0.5 / connp.size(); 386b117cd09SVijay Mahadevan else 387755f3dfbSVijay Mahadevan values_phi[tp] = factor * values_phi[tp] * 0.5 / (connp.size() * normsum); 388b117cd09SVijay Mahadevan } 389755f3dfbSVijay Mahadevan ierr = MatSetValues(*interpl, 1, &dofsc[tc], connp.size(), &dofsp[0], &values_phi[0], ADD_VALUES);CHKERRQ(ierr); 390b117cd09SVijay Mahadevan } 391a86ed394SVijay Mahadevan } 392b117cd09SVijay Mahadevan } 393304006b3SVijay Mahadevan if (vec) *vec = NULL; 394b117cd09SVijay Mahadevan ierr = MatAssemblyBegin(*interpl, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 395b117cd09SVijay Mahadevan ierr = MatAssemblyEnd(*interpl, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 396b117cd09SVijay Mahadevan PetscFunctionReturn(0); 397b117cd09SVijay Mahadevan } 398b117cd09SVijay Mahadevan 399b117cd09SVijay Mahadevan #undef __FUNCT__ 400b117cd09SVijay Mahadevan #define __FUNCT__ "DMCreateInjection_Moab" 401b117cd09SVijay Mahadevan /*@ 402b117cd09SVijay Mahadevan DMCreateInjection_Moab - Generate a multi-level uniform refinement hierarchy 403b117cd09SVijay Mahadevan by succesively refining a coarse mesh, already defined in the DM object 404b117cd09SVijay Mahadevan provided by the user. 405b117cd09SVijay Mahadevan 406b117cd09SVijay Mahadevan Collective on MPI_Comm 407b117cd09SVijay Mahadevan 408b117cd09SVijay Mahadevan Input Parameter: 409b117cd09SVijay Mahadevan . dmb - The DMMoab object 410b117cd09SVijay Mahadevan 411b117cd09SVijay Mahadevan Output Parameter: 412b117cd09SVijay Mahadevan . nlevels - The number of levels of refinement needed to generate the hierarchy 413b117cd09SVijay Mahadevan + ldegrees - The degree of refinement at each level in the hierarchy 414b117cd09SVijay Mahadevan 415b117cd09SVijay Mahadevan Level: beginner 416b117cd09SVijay Mahadevan 417b117cd09SVijay Mahadevan .keywords: DMMoab, create, refinement 418b117cd09SVijay Mahadevan @*/ 419304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMCreateInjection_Moab(DM dm1, DM dm2, VecScatter* ctx) 420b117cd09SVijay Mahadevan { 421e882eb38SVijay Mahadevan //DM_Moab *dmmoab; 422b117cd09SVijay Mahadevan 423b117cd09SVijay Mahadevan PetscFunctionBegin; 424b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm1, DM_CLASSID, 1); 425b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm2, DM_CLASSID, 2); 426e882eb38SVijay Mahadevan //dmmoab = (DM_Moab*)(dm1)->data; 427b117cd09SVijay Mahadevan 428b117cd09SVijay Mahadevan PetscPrintf(PETSC_COMM_WORLD, "[DMCreateInjection_Moab] :: Placeholder\n"); 429b117cd09SVijay Mahadevan PetscFunctionReturn(0); 430b117cd09SVijay Mahadevan } 431b117cd09SVijay Mahadevan 432b117cd09SVijay Mahadevan 433b117cd09SVijay Mahadevan #undef __FUNCT__ 434b117cd09SVijay Mahadevan #define __FUNCT__ "DM_UMR_Moab_Private" 435b117cd09SVijay Mahadevan PetscErrorCode DM_UMR_Moab_Private(DM dm, MPI_Comm comm, PetscBool refine, DM *dmref) 436b117cd09SVijay Mahadevan { 437b117cd09SVijay Mahadevan PetscErrorCode ierr; 438e882eb38SVijay Mahadevan PetscInt i, dim; 439b117cd09SVijay Mahadevan DM dm2; 440b117cd09SVijay Mahadevan moab::ErrorCode merr; 441b117cd09SVijay Mahadevan DM_Moab *dmb = (DM_Moab*)dm->data, *dd2; 442b117cd09SVijay Mahadevan 443b117cd09SVijay Mahadevan PetscFunctionBegin; 444b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 445b117cd09SVijay Mahadevan PetscValidPointer(dmref, 4); 446b117cd09SVijay Mahadevan 447b117cd09SVijay Mahadevan if ( (dmb->hlevel == dmb->nhlevels && refine) || (dmb->hlevel == 0 && !refine) ) { 448e882eb38SVijay Mahadevan if (dmb->hlevel + 1 > dmb->nhlevels && refine) PetscInfo2(NULL, "Invalid multigrid refinement hierarchy level specified (%D). MOAB UMR max levels = %D. Creating a NULL object.\n", dmb->hlevel + 1, dmb->nhlevels); 449e882eb38SVijay Mahadevan if (dmb->hlevel - 1 < 0 && !refine) PetscInfo1(NULL, "Invalid multigrid coarsen hierarchy level specified (%D). Creating a NULL object.\n", dmb->hlevel - 1); 450e882eb38SVijay Mahadevan *dmref = PETSC_NULL; 451b117cd09SVijay Mahadevan PetscFunctionReturn(0); 452b117cd09SVijay Mahadevan } 453b117cd09SVijay Mahadevan 454b117cd09SVijay Mahadevan ierr = DMMoabCreate(PetscObjectComm((PetscObject)dm), &dm2);CHKERRQ(ierr); 455b117cd09SVijay Mahadevan dd2 = (DM_Moab*)dm2->data; 456b117cd09SVijay Mahadevan 457b117cd09SVijay Mahadevan dd2->mbiface = dmb->mbiface; 4589daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 459b117cd09SVijay Mahadevan dd2->pcomm = dmb->pcomm; 4609daf19fdSVijay Mahadevan #endif 461b117cd09SVijay Mahadevan dd2->icreatedinstance = PETSC_FALSE; 46264e1c140SVijay Mahadevan dd2->nghostrings = dmb->nghostrings; 463b117cd09SVijay Mahadevan 464b117cd09SVijay Mahadevan /* set the new level based on refinement/coarsening */ 465b117cd09SVijay Mahadevan if (refine) { 466b117cd09SVijay Mahadevan dd2->hlevel = dmb->hlevel + 1; 467b117cd09SVijay Mahadevan } 468b117cd09SVijay Mahadevan else { 469b117cd09SVijay Mahadevan dd2->hlevel = dmb->hlevel - 1; 470b117cd09SVijay Mahadevan } 471b117cd09SVijay Mahadevan 472b117cd09SVijay Mahadevan /* Copy the multilevel hierarchy pointers in MOAB */ 473b117cd09SVijay Mahadevan dd2->hierarchy = dmb->hierarchy; 474b117cd09SVijay Mahadevan dd2->nhlevels = dmb->nhlevels; 475b117cd09SVijay Mahadevan ierr = PetscMalloc1(dd2->nhlevels + 1, &dd2->hsets);CHKERRQ(ierr); 476b117cd09SVijay Mahadevan for (i = 0; i <= dd2->nhlevels; i++) { 477b117cd09SVijay Mahadevan dd2->hsets[i] = dmb->hsets[i]; 478b117cd09SVijay Mahadevan } 479b117cd09SVijay Mahadevan dd2->fileset = dd2->hsets[dd2->hlevel]; 480b117cd09SVijay Mahadevan 481b117cd09SVijay Mahadevan /* do the remaining initializations for DMMoab */ 482b117cd09SVijay Mahadevan dd2->bs = dmb->bs; 483b117cd09SVijay Mahadevan dd2->numFields = dmb->numFields; 484b117cd09SVijay Mahadevan dd2->rw_dbglevel = dmb->rw_dbglevel; 485b117cd09SVijay Mahadevan dd2->partition_by_rank = dmb->partition_by_rank; 486b117cd09SVijay Mahadevan ierr = PetscStrcpy(dd2->extra_read_options, dmb->extra_read_options);CHKERRQ(ierr); 487b117cd09SVijay Mahadevan ierr = PetscStrcpy(dd2->extra_write_options, dmb->extra_write_options);CHKERRQ(ierr); 488b117cd09SVijay Mahadevan dd2->read_mode = dmb->read_mode; 489b117cd09SVijay Mahadevan dd2->write_mode = dmb->write_mode; 490b117cd09SVijay Mahadevan 491b117cd09SVijay Mahadevan /* set global ID tag handle */ 492b117cd09SVijay Mahadevan ierr = DMMoabSetLocalToGlobalTag(dm2, dmb->ltog_tag);CHKERRQ(ierr); 493b117cd09SVijay Mahadevan 494b117cd09SVijay Mahadevan merr = dd2->mbiface->tag_get_handle(MATERIAL_SET_TAG_NAME, dd2->material_tag);MBERRNM(merr); 495b117cd09SVijay Mahadevan 496b117cd09SVijay Mahadevan ierr = DMSetOptionsPrefix(dm2, ((PetscObject)dm)->prefix);CHKERRQ(ierr); 497b117cd09SVijay Mahadevan ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 498b117cd09SVijay Mahadevan ierr = DMSetDimension(dm2, dim);CHKERRQ(ierr); 499b117cd09SVijay Mahadevan 500b117cd09SVijay Mahadevan /* allow overloaded (user replaced) operations to be inherited by refinement clones */ 501b117cd09SVijay Mahadevan dm2->ops->creatematrix = dm->ops->creatematrix; 502b117cd09SVijay Mahadevan 503b117cd09SVijay Mahadevan /* copy fill information if given */ 504b117cd09SVijay Mahadevan ierr = DMMoabSetBlockFills(dm2, dmb->dfill, dmb->ofill);CHKERRQ(ierr); 505b117cd09SVijay Mahadevan 506b117cd09SVijay Mahadevan /* copy vector type information */ 507b117cd09SVijay Mahadevan ierr = DMSetMatType(dm2, dm->mattype);CHKERRQ(ierr); 508b117cd09SVijay Mahadevan ierr = DMSetVecType(dm2, dm->vectype);CHKERRQ(ierr); 5093f1c6e43SVijay Mahadevan dd2->numFields = dmb->numFields; 5103f1c6e43SVijay Mahadevan if (dmb->numFields) { 511b117cd09SVijay Mahadevan ierr = DMMoabSetFieldNames(dm2, dmb->numFields, dmb->fieldNames);CHKERRQ(ierr); 5123f1c6e43SVijay Mahadevan } 513b117cd09SVijay Mahadevan 514b117cd09SVijay Mahadevan ierr = DMSetFromOptions(dm2);CHKERRQ(ierr); 515b117cd09SVijay Mahadevan 516b117cd09SVijay Mahadevan /* recreate Dof numbering for the refined DM and make sure the distribution is correctly populated */ 517b117cd09SVijay Mahadevan ierr = DMSetUp(dm2);CHKERRQ(ierr); 518b117cd09SVijay Mahadevan 519b117cd09SVijay Mahadevan *dmref = dm2; 520b117cd09SVijay Mahadevan PetscFunctionReturn(0); 521b117cd09SVijay Mahadevan } 522b117cd09SVijay Mahadevan 523b117cd09SVijay Mahadevan 524b117cd09SVijay Mahadevan #undef __FUNCT__ 525b117cd09SVijay Mahadevan #define __FUNCT__ "DMRefine_Moab" 526b117cd09SVijay Mahadevan /*@ 527b117cd09SVijay Mahadevan DMRefine_Moab - Generate a multi-level uniform refinement hierarchy 528b117cd09SVijay Mahadevan by succesively refining a coarse mesh, already defined in the DM object 529b117cd09SVijay Mahadevan provided by the user. 530b117cd09SVijay Mahadevan 531e882eb38SVijay Mahadevan Collective on DM 532b117cd09SVijay Mahadevan 533b117cd09SVijay Mahadevan Input Parameter: 534e882eb38SVijay Mahadevan + dm - The DMMoab object 535e882eb38SVijay Mahadevan - comm - the communicator to contain the new DM object (or MPI_COMM_NULL) 536b117cd09SVijay Mahadevan 537b117cd09SVijay Mahadevan Output Parameter: 538e882eb38SVijay Mahadevan . dmf - the refined DM, or NULL 539b117cd09SVijay Mahadevan 540e882eb38SVijay Mahadevan Note: If no refinement was done, the return value is NULL 541e882eb38SVijay Mahadevan 542e882eb38SVijay Mahadevan Level: developer 543b117cd09SVijay Mahadevan 544b117cd09SVijay Mahadevan .keywords: DMMoab, create, refinement 545b117cd09SVijay Mahadevan @*/ 546304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMRefine_Moab(DM dm, MPI_Comm comm, DM* dmf) 547b117cd09SVijay Mahadevan { 548b117cd09SVijay Mahadevan PetscErrorCode ierr; 549b117cd09SVijay Mahadevan 550b117cd09SVijay Mahadevan PetscFunctionBegin; 551b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 552b117cd09SVijay Mahadevan 553b117cd09SVijay Mahadevan ierr = DM_UMR_Moab_Private(dm, comm, PETSC_TRUE, dmf);CHKERRQ(ierr); 554b117cd09SVijay Mahadevan PetscFunctionReturn(0); 555b117cd09SVijay Mahadevan } 556b117cd09SVijay Mahadevan 557b117cd09SVijay Mahadevan #undef __FUNCT__ 558b117cd09SVijay Mahadevan #define __FUNCT__ "DMCoarsen_Moab" 559b117cd09SVijay Mahadevan /*@ 560b117cd09SVijay Mahadevan DMCoarsen_Moab - Generate a multi-level uniform refinement hierarchy 561b117cd09SVijay Mahadevan by succesively refining a coarse mesh, already defined in the DM object 562b117cd09SVijay Mahadevan provided by the user. 563b117cd09SVijay Mahadevan 564e882eb38SVijay Mahadevan Collective on DM 565b117cd09SVijay Mahadevan 566b117cd09SVijay Mahadevan Input Parameter: 567e882eb38SVijay Mahadevan . dm - The DMMoab object 568e882eb38SVijay Mahadevan - comm - the communicator to contain the new DM object (or MPI_COMM_NULL) 569b117cd09SVijay Mahadevan 570b117cd09SVijay Mahadevan Output Parameter: 571e882eb38SVijay Mahadevan . dmf - the coarsened DM, or NULL 572b117cd09SVijay Mahadevan 573e882eb38SVijay Mahadevan Note: If no coarsening was done, the return value is NULL 574b117cd09SVijay Mahadevan 575e882eb38SVijay Mahadevan Level: developer 576e882eb38SVijay Mahadevan 577e882eb38SVijay Mahadevan .keywords: DMMoab, create, coarsening 578b117cd09SVijay Mahadevan @*/ 579304006b3SVijay Mahadevan PETSC_EXTERN PetscErrorCode DMCoarsen_Moab(DM dm, MPI_Comm comm, DM* dmc) 580b117cd09SVijay Mahadevan { 581b117cd09SVijay Mahadevan PetscErrorCode ierr; 582b117cd09SVijay Mahadevan 583b117cd09SVijay Mahadevan PetscFunctionBegin; 584b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 585b117cd09SVijay Mahadevan 586b117cd09SVijay Mahadevan ierr = DM_UMR_Moab_Private(dm, comm, PETSC_FALSE, dmc);CHKERRQ(ierr); 587b117cd09SVijay Mahadevan PetscFunctionReturn(0); 588b117cd09SVijay Mahadevan } 589