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 6*ca4445c7SIlya Fursov // A helper function to convert Real vector to Scalar vector (required by MatSetValues) 7*ca4445c7SIlya Fursov static inline std::vector<PetscScalar> VecReal_to_VecScalar(const std::vector<PetscReal> &v) 8*ca4445c7SIlya Fursov { 9*ca4445c7SIlya Fursov std::vector<PetscScalar> res(v.size()); 10*ca4445c7SIlya Fursov for (size_t i = 0; i < res.size(); i++) res[i] = v[i]; 11*ca4445c7SIlya Fursov return res; 12*ca4445c7SIlya Fursov } 13*ca4445c7SIlya Fursov 14cab5ea25SPierre Jolivet /*@C 15b117cd09SVijay Mahadevan DMMoabGenerateHierarchy - Generate a multi-level uniform refinement hierarchy 1620f4b53cSBarry Smith by succesively refining a coarse mesh, already defined in the `DM` object 17b117cd09SVijay Mahadevan provided by the user. 18b117cd09SVijay Mahadevan 19d083f849SBarry Smith Collective 20b117cd09SVijay Mahadevan 21b117cd09SVijay Mahadevan Input Parameter: 2260225df5SJacob Faibussowitsch . dm - The `DMMOAB` object 23b117cd09SVijay Mahadevan 24d8d19677SJose E. Roman Output Parameters: 25b117cd09SVijay Mahadevan + nlevels - The number of levels of refinement needed to generate the hierarchy 26a2b725a8SWilliam Gropp - ldegrees - The degree of refinement at each level in the hierarchy 27b117cd09SVijay Mahadevan 28b117cd09SVijay Mahadevan Level: beginner 29b117cd09SVijay Mahadevan 30a4e35b19SJacob Faibussowitsch .seealso: `DMMoabCreate()` 31b117cd09SVijay Mahadevan @*/ 32d71ae5a4SJacob Faibussowitsch PetscErrorCode DMMoabGenerateHierarchy(DM dm, PetscInt nlevels, PetscInt *ldegrees) 33d71ae5a4SJacob Faibussowitsch { 34b117cd09SVijay Mahadevan DM_Moab *dmmoab; 35b117cd09SVijay Mahadevan moab::ErrorCode merr; 362417220eSVijay Mahadevan PetscInt *pdegrees, ilevel; 37e882eb38SVijay Mahadevan std::vector<moab::EntityHandle> hsets; 38b117cd09SVijay Mahadevan 39b117cd09SVijay Mahadevan PetscFunctionBegin; 40b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 41b117cd09SVijay Mahadevan dmmoab = (DM_Moab *)(dm)->data; 42b117cd09SVijay Mahadevan 43b117cd09SVijay Mahadevan if (!ldegrees) { 449566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nlevels, &pdegrees)); 452417220eSVijay Mahadevan for (ilevel = 0; ilevel < nlevels; ilevel++) pdegrees[ilevel] = 2; /* default = Degree 2 refinement */ 469371c9d4SSatish Balay } else pdegrees = ldegrees; 47b117cd09SVijay Mahadevan 48b117cd09SVijay Mahadevan /* initialize set level refinement data for hierarchy */ 49b117cd09SVijay Mahadevan dmmoab->nhlevels = nlevels; 50b117cd09SVijay Mahadevan 51b117cd09SVijay Mahadevan /* Instantiate the nested refinement class */ 529daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 533f1c6e43SVijay Mahadevan dmmoab->hierarchy = new moab::NestedRefine(dynamic_cast<moab::Core *>(dmmoab->mbiface), dmmoab->pcomm, dmmoab->fileset); 549daf19fdSVijay Mahadevan #else 559daf19fdSVijay Mahadevan dmmoab->hierarchy = new moab::NestedRefine(dynamic_cast<moab::Core *>(dmmoab->mbiface), NULL, dmmoab->fileset); 569daf19fdSVijay Mahadevan #endif 57b117cd09SVijay Mahadevan 589566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nlevels + 1, &dmmoab->hsets)); 59e882eb38SVijay Mahadevan 60b117cd09SVijay Mahadevan /* generate the mesh hierarchy */ 619371c9d4SSatish Balay merr = dmmoab->hierarchy->generate_mesh_hierarchy(nlevels, pdegrees, hsets, false); 629371c9d4SSatish Balay MBERRNM(merr); 63e882eb38SVijay Mahadevan 649daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 65755f3dfbSVijay Mahadevan if (dmmoab->pcomm->size() > 1) { 669371c9d4SSatish Balay merr = dmmoab->hierarchy->exchange_ghosts(hsets, dmmoab->nghostrings); 679371c9d4SSatish Balay MBERRNM(merr); 68755f3dfbSVijay Mahadevan } 699daf19fdSVijay Mahadevan #endif 7049d66b22SVijay Mahadevan 7164e1c140SVijay Mahadevan /* copy the mesh sets for nested refinement hierarchy */ 72e92d1c7cSVijay Mahadevan dmmoab->hsets[0] = hsets[0]; 739371c9d4SSatish Balay for (ilevel = 1; ilevel <= nlevels; ilevel++) { 742417220eSVijay Mahadevan dmmoab->hsets[ilevel] = hsets[ilevel]; 75e92d1c7cSVijay Mahadevan 769c368985SVijay Mahadevan #ifdef MOAB_HAVE_MPI 779371c9d4SSatish Balay merr = dmmoab->pcomm->assign_global_ids(hsets[ilevel], dmmoab->dim, 0, false, true, false); 789371c9d4SSatish Balay MBERRNM(merr); 799c368985SVijay Mahadevan #endif 80e92d1c7cSVijay Mahadevan 81e92d1c7cSVijay Mahadevan /* Update material and other geometric tags from parent to child sets */ 829371c9d4SSatish Balay merr = dmmoab->hierarchy->update_special_tags(ilevel, hsets[ilevel]); 839371c9d4SSatish Balay MBERRNM(merr); 842417220eSVijay Mahadevan } 8564e1c140SVijay Mahadevan 8664e1c140SVijay Mahadevan hsets.clear(); 8748a46eb9SPierre Jolivet if (!ldegrees) PetscCall(PetscFree(pdegrees)); 883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 89b117cd09SVijay Mahadevan } 90b117cd09SVijay Mahadevan 91a4e35b19SJacob Faibussowitsch // PetscClangLinter pragma ignore: -fdoc-* 92a4e35b19SJacob Faibussowitsch /* 9320f4b53cSBarry Smith DMRefineHierarchy_Moab - Generate a multi-level `DM` hierarchy 94e882eb38SVijay Mahadevan by succesively refining a coarse mesh. 95b117cd09SVijay Mahadevan 96d083f849SBarry Smith Collective 97b117cd09SVijay Mahadevan 98b117cd09SVijay Mahadevan Input Parameter: 9920f4b53cSBarry Smith . dm - The `DMMOAB` object 100b117cd09SVijay Mahadevan 101d8d19677SJose E. Roman Output Parameters: 102e882eb38SVijay Mahadevan + nlevels - The number of levels of refinement needed to generate the hierarchy 103a2b725a8SWilliam Gropp - dmf - The DM objects after successive refinement of the hierarchy 104b117cd09SVijay Mahadevan 105b117cd09SVijay Mahadevan Level: beginner 106a4e35b19SJacob Faibussowitsch */ 107d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMRefineHierarchy_Moab(DM dm, PetscInt nlevels, DM dmf[]) 108d71ae5a4SJacob Faibussowitsch { 109b117cd09SVijay Mahadevan PetscInt i; 110b117cd09SVijay Mahadevan 111b117cd09SVijay Mahadevan PetscFunctionBegin; 112b117cd09SVijay Mahadevan 1139566063dSJacob Faibussowitsch PetscCall(DMRefine(dm, PetscObjectComm((PetscObject)dm), &dmf[0])); 11448a46eb9SPierre Jolivet for (i = 1; i < nlevels; i++) PetscCall(DMRefine(dmf[i - 1], PetscObjectComm((PetscObject)dm), &dmf[i])); 1153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 116b117cd09SVijay Mahadevan } 117b117cd09SVijay Mahadevan 118a4e35b19SJacob Faibussowitsch // PetscClangLinter pragma ignore: -fdoc-* 119a4e35b19SJacob Faibussowitsch /* 12020f4b53cSBarry Smith DMCoarsenHierarchy_Moab - Generate a multi-level `DM` hierarchy 121e882eb38SVijay Mahadevan by succesively coarsening a refined mesh. 122b117cd09SVijay Mahadevan 123d083f849SBarry Smith Collective 124b117cd09SVijay Mahadevan 125b117cd09SVijay Mahadevan Input Parameter: 12620f4b53cSBarry Smith . dm - The `DMMOAB` object 127b117cd09SVijay Mahadevan 128d8d19677SJose E. Roman Output Parameters: 129e882eb38SVijay Mahadevan + nlevels - The number of levels of refinement needed to generate the hierarchy 13020f4b53cSBarry Smith - dmc - The `DM` objects after successive coarsening of the hierarchy 131b117cd09SVijay Mahadevan 132b117cd09SVijay Mahadevan Level: beginner 133a4e35b19SJacob Faibussowitsch */ 134d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCoarsenHierarchy_Moab(DM dm, PetscInt nlevels, DM dmc[]) 135d71ae5a4SJacob Faibussowitsch { 136b117cd09SVijay Mahadevan PetscInt i; 137b117cd09SVijay Mahadevan 138b117cd09SVijay Mahadevan PetscFunctionBegin; 139b117cd09SVijay Mahadevan 1409566063dSJacob Faibussowitsch PetscCall(DMCoarsen(dm, PetscObjectComm((PetscObject)dm), &dmc[0])); 14148a46eb9SPierre Jolivet for (i = 1; i < nlevels; i++) PetscCall(DMCoarsen(dmc[i - 1], PetscObjectComm((PetscObject)dm), &dmc[i])); 1423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 143b117cd09SVijay Mahadevan } 144b117cd09SVijay Mahadevan 1452417220eSVijay Mahadevan PETSC_EXTERN PetscErrorCode DMMoab_Compute_NNZ_From_Connectivity(DM, PetscInt *, PetscInt *, PetscInt *, PetscInt *, PetscBool); 146b117cd09SVijay Mahadevan 147a4e35b19SJacob Faibussowitsch // PetscClangLinter pragma ignore: -fdoc-* 148a4e35b19SJacob Faibussowitsch /* 149e882eb38SVijay Mahadevan DMCreateInterpolation_Moab - Generate the interpolation operators to transform 150e882eb38SVijay Mahadevan operators (matrices, vectors) from parent level to child level as defined by 15120f4b53cSBarry Smith the `DM` inputs provided by the user. 152b117cd09SVijay Mahadevan 153d083f849SBarry Smith Collective 154b117cd09SVijay Mahadevan 155d8d19677SJose E. Roman Input Parameters: 15660225df5SJacob Faibussowitsch + dmp - The `DMMOAB` object 15760225df5SJacob Faibussowitsch - dmc - the second, finer `DMMOAB` object 158b117cd09SVijay Mahadevan 159d8d19677SJose E. Roman Output Parameters: 160e882eb38SVijay Mahadevan + interpl - The interpolation operator for transferring data between the levels 161e882eb38SVijay Mahadevan - vec - The scaling vector (optional) 162b117cd09SVijay Mahadevan 163e882eb38SVijay Mahadevan Level: developer 164a4e35b19SJacob Faibussowitsch */ 165d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreateInterpolation_Moab(DM dmp, DM dmc, Mat *interpl, Vec *vec) 166d71ae5a4SJacob Faibussowitsch { 167755f3dfbSVijay Mahadevan DM_Moab *dmbp, *dmbc; 168b117cd09SVijay Mahadevan moab::ErrorCode merr; 169e882eb38SVijay Mahadevan PetscInt dim; 1703f1c6e43SVijay Mahadevan PetscReal factor; 171ce27a4eeSVijay Mahadevan PetscInt innz, *nnz, ionz, *onz; 172755f3dfbSVijay Mahadevan PetscInt nlsizp, nlsizc, nlghsizp, ngsizp, ngsizc; 173a86ed394SVijay Mahadevan const PetscBool use_consistent_bases = PETSC_TRUE; 174b117cd09SVijay Mahadevan 175b117cd09SVijay Mahadevan PetscFunctionBegin; 176755f3dfbSVijay Mahadevan PetscValidHeaderSpecific(dmp, DM_CLASSID, 1); 177755f3dfbSVijay Mahadevan PetscValidHeaderSpecific(dmc, DM_CLASSID, 2); 178755f3dfbSVijay Mahadevan dmbp = (DM_Moab *)(dmp)->data; 179755f3dfbSVijay Mahadevan dmbc = (DM_Moab *)(dmc)->data; 180755f3dfbSVijay Mahadevan nlsizp = dmbp->nloc; // *dmb1->numFields; 181755f3dfbSVijay Mahadevan nlsizc = dmbc->nloc; // *dmb2->numFields; 182755f3dfbSVijay Mahadevan ngsizp = dmbp->n; // *dmb1->numFields; 183755f3dfbSVijay Mahadevan ngsizc = dmbc->n; // *dmb2->numFields; 184755f3dfbSVijay Mahadevan nlghsizp = (dmbp->nloc + dmbp->nghost); // *dmb1->numFields; 185b117cd09SVijay Mahadevan 1862417220eSVijay Mahadevan // Columns = Parent DoFs ; Rows = Child DoFs 187755f3dfbSVijay Mahadevan // Interpolation matrix: \sum_{i=1}^P Owned(Child) * (Owned(Parent) + Ghosted(Parent)) 188755f3dfbSVijay Mahadevan // Size: nlsizc * nlghsizp 1893ba16761SJacob 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)); 190b117cd09SVijay Mahadevan 1919566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dmp, &dim)); 192a86ed394SVijay Mahadevan 193941e0cffSVijay Mahadevan /* allocate the nnz, onz arrays based on block size and local nodes */ 1949566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(nlsizc, &nnz, nlsizc, &onz)); 195941e0cffSVijay Mahadevan 196941e0cffSVijay Mahadevan /* Loop through the local elements and compute the relation between the current parent and the refined_level. */ 1972417220eSVijay Mahadevan for (moab::Range::iterator iter = dmbc->vowned->begin(); iter != dmbc->vowned->end(); iter++) { 1982417220eSVijay Mahadevan const moab::EntityHandle vhandle = *iter; 1992417220eSVijay Mahadevan /* define local variables */ 2002417220eSVijay Mahadevan moab::EntityHandle parent; 2012417220eSVijay Mahadevan std::vector<moab::EntityHandle> adjs; 2022417220eSVijay Mahadevan moab::Range found; 2032417220eSVijay Mahadevan 2042417220eSVijay Mahadevan /* store the vertex DoF number */ 2052417220eSVijay Mahadevan const int ldof = dmbc->lidmap[vhandle - dmbc->seqstart]; 2062417220eSVijay Mahadevan 2072417220eSVijay Mahadevan /* Get adjacency information for current vertex - i.e., all elements of dimension (dim) that connects 2082417220eSVijay Mahadevan to the current vertex. We can then decipher if a vertex is ghosted or not and compute the 2092417220eSVijay Mahadevan non-zero pattern accordingly. */ 2109371c9d4SSatish Balay merr = dmbc->hierarchy->get_adjacencies(vhandle, dmbc->dim, adjs); 2119371c9d4SSatish Balay MBERRNM(merr); 2122417220eSVijay Mahadevan 2132417220eSVijay Mahadevan /* loop over vertices and update the number of connectivity */ 2142417220eSVijay Mahadevan for (unsigned jter = 0; jter < adjs.size(); jter++) { 2152417220eSVijay Mahadevan const moab::EntityHandle jhandle = adjs[jter]; 216941e0cffSVijay Mahadevan 217941e0cffSVijay Mahadevan /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */ 2189371c9d4SSatish Balay merr = dmbc->hierarchy->child_to_parent(jhandle, dmbc->hlevel, dmbp->hlevel, &parent); 2199371c9d4SSatish Balay MBERRNM(merr); 220941e0cffSVijay Mahadevan 2212417220eSVijay Mahadevan /* Get connectivity information in canonical ordering for the local element */ 2222417220eSVijay Mahadevan std::vector<moab::EntityHandle> connp; 2239371c9d4SSatish Balay merr = dmbp->hierarchy->get_connectivity(parent, dmbp->hlevel, connp); 2249371c9d4SSatish Balay MBERRNM(merr); 225a86ed394SVijay Mahadevan 2262417220eSVijay Mahadevan for (unsigned ic = 0; ic < connp.size(); ++ic) { 2272417220eSVijay Mahadevan /* loop over each element connected to the adjacent vertex and update as needed */ 2282417220eSVijay Mahadevan /* find the truly user-expected layer of ghosted entities to decipher NNZ pattern */ 2292417220eSVijay Mahadevan if (found.find(connp[ic]) != found.end()) continue; /* make sure we don't double count shared vertices */ 2302417220eSVijay Mahadevan if (dmbp->vghost->find(connp[ic]) != dmbp->vghost->end()) onz[ldof]++; /* update out-of-proc onz */ 2312417220eSVijay Mahadevan else nnz[ldof]++; /* else local vertex */ 2322417220eSVijay Mahadevan found.insert(connp[ic]); 2332417220eSVijay Mahadevan } 2342417220eSVijay Mahadevan } 235941e0cffSVijay Mahadevan } 236941e0cffSVijay Mahadevan 2379371c9d4SSatish Balay for (int i = 0; i < nlsizc; i++) nnz[i] += 1; /* self count the node */ 238941e0cffSVijay Mahadevan 239941e0cffSVijay Mahadevan ionz = onz[0]; 240941e0cffSVijay Mahadevan innz = nnz[0]; 241755f3dfbSVijay Mahadevan for (int tc = 0; tc < nlsizc; tc++) { 242941e0cffSVijay Mahadevan // check for maximum allowed sparsity = fully dense 243755f3dfbSVijay Mahadevan nnz[tc] = std::min(nlsizp, nnz[tc]); 2442417220eSVijay Mahadevan onz[tc] = std::min(ngsizp - nlsizp, onz[tc]); 2452417220eSVijay Mahadevan 2463ba16761SJacob Faibussowitsch PetscCall(PetscInfo(NULL, " %d: NNZ = %d, ONZ = %d\n", tc, nnz[tc], onz[tc])); 247941e0cffSVijay Mahadevan 248941e0cffSVijay Mahadevan innz = (innz < nnz[tc] ? nnz[tc] : innz); 249941e0cffSVijay Mahadevan ionz = (ionz < onz[tc] ? onz[tc] : ionz); 250941e0cffSVijay Mahadevan } 25164e1c140SVijay Mahadevan 25264e1c140SVijay Mahadevan /* create interpolation matrix */ 2539566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)dmc), interpl)); 2549566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*interpl, nlsizc, nlsizp, ngsizc, ngsizp)); 2559566063dSJacob Faibussowitsch PetscCall(MatSetType(*interpl, MATAIJ)); 2569566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(*interpl)); 25764e1c140SVijay Mahadevan 2589566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(*interpl, innz, nnz)); 2599566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(*interpl, innz, nnz, ionz, onz)); 260941e0cffSVijay Mahadevan 261941e0cffSVijay Mahadevan /* clean up temporary memory */ 2629566063dSJacob Faibussowitsch PetscCall(PetscFree2(nnz, onz)); 263b117cd09SVijay Mahadevan 264b117cd09SVijay Mahadevan /* set up internal matrix data-structures */ 2659566063dSJacob Faibussowitsch PetscCall(MatSetUp(*interpl)); 266b117cd09SVijay Mahadevan 267a86ed394SVijay Mahadevan /* Define variables for assembly */ 268a86ed394SVijay Mahadevan std::vector<moab::EntityHandle> children; 269a86ed394SVijay Mahadevan std::vector<moab::EntityHandle> connp, connc; 270a86ed394SVijay Mahadevan std::vector<PetscReal> pcoords, ccoords, values_phi; 271*ca4445c7SIlya Fursov std::vector<PetscScalar> values_phi_scalar; 272a86ed394SVijay Mahadevan 273a86ed394SVijay Mahadevan if (use_consistent_bases) { 2742417220eSVijay Mahadevan const moab::EntityHandle ehandle = dmbp->elocal->front(); 275a86ed394SVijay Mahadevan 2769371c9d4SSatish Balay merr = dmbp->hierarchy->parent_to_child(ehandle, dmbp->hlevel, dmbc->hlevel, children); 2779371c9d4SSatish Balay MBERRNM(merr); 278a86ed394SVijay Mahadevan 279a86ed394SVijay Mahadevan /* Get connectivity and coordinates of the parent vertices */ 2809371c9d4SSatish Balay merr = dmbp->hierarchy->get_connectivity(ehandle, dmbp->hlevel, connp); 2819371c9d4SSatish Balay MBERRNM(merr); 2829371c9d4SSatish Balay merr = dmbc->mbiface->get_connectivity(&children[0], children.size(), connc); 2839371c9d4SSatish Balay MBERRNM(merr); 284a86ed394SVijay Mahadevan 285a86ed394SVijay Mahadevan std::vector<PetscReal> natparam(3 * connc.size(), 0.0); 286a86ed394SVijay Mahadevan pcoords.resize(connp.size() * 3); 287a86ed394SVijay Mahadevan ccoords.resize(connc.size() * 3); 288a86ed394SVijay Mahadevan values_phi.resize(connp.size() * connc.size()); 289a86ed394SVijay Mahadevan /* Get coordinates for connectivity entities in canonical order for both coarse and finer levels */ 2909371c9d4SSatish Balay merr = dmbp->hierarchy->get_coordinates(&connp[0], connp.size(), dmbp->hlevel, &pcoords[0]); 2919371c9d4SSatish Balay MBERRNM(merr); 2929371c9d4SSatish Balay merr = dmbc->hierarchy->get_coordinates(&connc[0], connc.size(), dmbc->hlevel, &ccoords[0]); 2939371c9d4SSatish Balay MBERRNM(merr); 294a86ed394SVijay Mahadevan 295a86ed394SVijay Mahadevan /* Set values: For each DOF in coarse grid cell, set the contribution or PHI evaluated at each fine grid DOF point */ 296a86ed394SVijay Mahadevan for (unsigned tc = 0; tc < connc.size(); tc++) { 297a86ed394SVijay Mahadevan const PetscInt offset = tc * 3; 298a86ed394SVijay Mahadevan 299a86ed394SVijay Mahadevan /* Scale ccoords relative to pcoords */ 3009566063dSJacob Faibussowitsch PetscCall(DMMoabPToRMapping(dim, connp.size(), &pcoords[0], &ccoords[offset], &natparam[offset], &values_phi[connp.size() * tc])); 301a86ed394SVijay Mahadevan } 3021baa6e33SBarry Smith } else { 303a86ed394SVijay Mahadevan factor = std::pow(2.0 /*degree_P_for_refinement*/, (dmbc->hlevel - dmbp->hlevel) * dmbp->dim * 1.0); 304a86ed394SVijay Mahadevan } 305a86ed394SVijay Mahadevan 306755f3dfbSVijay Mahadevan /* TODO: Decipher the correct non-zero pattern. There is still some issue with onz allocation */ 3079566063dSJacob Faibussowitsch PetscCall(MatSetOption(*interpl, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 308b117cd09SVijay Mahadevan 309e882eb38SVijay Mahadevan /* Loop through the remaining vertices. These vertices appear only on the current refined_level. */ 310*ca4445c7SIlya Fursov values_phi_scalar = VecReal_to_VecScalar(values_phi); 3112417220eSVijay Mahadevan for (moab::Range::iterator iter = dmbp->elocal->begin(); iter != dmbp->elocal->end(); iter++) { 312b117cd09SVijay Mahadevan const moab::EntityHandle ehandle = *iter; 313b117cd09SVijay Mahadevan 314b117cd09SVijay Mahadevan /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */ 315a86ed394SVijay Mahadevan children.clear(); 316a86ed394SVijay Mahadevan connc.clear(); 3179371c9d4SSatish Balay merr = dmbp->hierarchy->parent_to_child(ehandle, dmbp->hlevel, dmbc->hlevel, children); 3189371c9d4SSatish Balay MBERRNM(merr); 319b117cd09SVijay Mahadevan 320b117cd09SVijay Mahadevan /* Get connectivity and coordinates of the parent vertices */ 3219371c9d4SSatish Balay merr = dmbp->hierarchy->get_connectivity(ehandle, dmbp->hlevel, connp); 3229371c9d4SSatish Balay MBERRNM(merr); 3239371c9d4SSatish Balay merr = dmbc->mbiface->get_connectivity(&children[0], children.size(), connc); 3249371c9d4SSatish Balay MBERRNM(merr); 325b117cd09SVijay Mahadevan 326a86ed394SVijay Mahadevan pcoords.resize(connp.size() * 3); 327a86ed394SVijay Mahadevan ccoords.resize(connc.size() * 3); 328b117cd09SVijay Mahadevan /* Get coordinates for connectivity entities in canonical order for both coarse and finer levels */ 3299371c9d4SSatish Balay merr = dmbp->hierarchy->get_coordinates(&connp[0], connp.size(), dmbp->hlevel, &pcoords[0]); 3309371c9d4SSatish Balay MBERRNM(merr); 3319371c9d4SSatish Balay merr = dmbc->hierarchy->get_coordinates(&connc[0], connc.size(), dmbc->hlevel, &ccoords[0]); 3329371c9d4SSatish Balay MBERRNM(merr); 333b117cd09SVijay Mahadevan 334b117cd09SVijay Mahadevan std::vector<int> dofsp(connp.size()), dofsc(connc.size()); 335b117cd09SVijay Mahadevan /* TODO: specific to scalar system - use GetDofs */ 3369566063dSJacob Faibussowitsch PetscCall(DMMoabGetDofsBlocked(dmp, connp.size(), &connp[0], &dofsp[0])); 3379566063dSJacob Faibussowitsch PetscCall(DMMoabGetDofsBlocked(dmc, connc.size(), &connc[0], &dofsc[0])); 338b117cd09SVijay Mahadevan 339a86ed394SVijay Mahadevan /* Compute the actual interpolation weights when projecting solution/residual between levels */ 340a86ed394SVijay Mahadevan if (use_consistent_bases) { 341145b44c9SPierre Jolivet /* Use the cached values of natural parametric coordinates and basis pre-evaluated. 342a86ed394SVijay Mahadevan We are making an assumption here that UMR used in GMG to generate the hierarchy uses 343a86ed394SVijay Mahadevan the same template for all elements; This will fail for mixed element meshes (TRI/QUAD). 344a86ed394SVijay Mahadevan 3452417220eSVijay Mahadevan TODO: Fix the above assumption by caching data for families (especially for Tets and mixed meshes) 346a86ed394SVijay Mahadevan */ 347a86ed394SVijay Mahadevan 348a86ed394SVijay Mahadevan /* Set values: For each DOF in coarse grid cell, set the contribution or PHI evaluated at each fine grid DOF point */ 349a86ed394SVijay Mahadevan for (unsigned tc = 0; tc < connc.size(); tc++) { 350a86ed394SVijay Mahadevan /* TODO: Check if we should be using INSERT_VALUES instead */ 351*ca4445c7SIlya Fursov PetscCall(MatSetValues(*interpl, 1, &dofsc[tc], connp.size(), &dofsp[0], &values_phi_scalar[connp.size() * tc], ADD_VALUES)); 352a86ed394SVijay Mahadevan } 3531baa6e33SBarry Smith } else { 354e882eb38SVijay Mahadevan /* Compute the interpolation weights by determining distance of 1-ring 355755f3dfbSVijay Mahadevan neighbor vertices from current vertex 356755f3dfbSVijay Mahadevan 357a86ed394SVijay Mahadevan This should be used only when FEM basis is not used for the discretization. 358a86ed394SVijay Mahadevan Else, the consistent interface to compute the basis function for interpolation 359a86ed394SVijay Mahadevan between the levels should be evaluated correctly to preserve convergence of GMG. 3602417220eSVijay Mahadevan Shephard's basis will be terrible for any unsmooth problems. 361755f3dfbSVijay Mahadevan */ 362755f3dfbSVijay Mahadevan values_phi.resize(connp.size()); 363755f3dfbSVijay Mahadevan for (unsigned tc = 0; tc < connc.size(); tc++) { 364a86ed394SVijay Mahadevan PetscReal normsum = 0.0; 365*ca4445c7SIlya Fursov std::vector<PetscScalar> values_phi_scalar2; 366755f3dfbSVijay Mahadevan for (unsigned tp = 0; tp < connp.size(); tp++) { 367755f3dfbSVijay Mahadevan values_phi[tp] = 0.0; 3689371c9d4SSatish Balay for (unsigned k = 0; k < 3; k++) values_phi[tp] += std::pow(pcoords[tp * 3 + k] - ccoords[k + tc * 3], dim); 369755f3dfbSVijay Mahadevan if (values_phi[tp] < 1e-12) { 370755f3dfbSVijay Mahadevan values_phi[tp] = 1e12; 3711baa6e33SBarry Smith } else { 372755f3dfbSVijay Mahadevan //values_phi[tp] = std::pow(values_phi[tp], -1.0/dim); 373755f3dfbSVijay Mahadevan values_phi[tp] = std::pow(values_phi[tp], -1.0); 374755f3dfbSVijay Mahadevan normsum += values_phi[tp]; 375b117cd09SVijay Mahadevan } 376b117cd09SVijay Mahadevan } 377755f3dfbSVijay Mahadevan for (unsigned tp = 0; tp < connp.size(); tp++) { 3789371c9d4SSatish Balay if (values_phi[tp] > 1e11) values_phi[tp] = factor * 0.5 / connp.size(); 3799371c9d4SSatish Balay else values_phi[tp] = factor * values_phi[tp] * 0.5 / (connp.size() * normsum); 380b117cd09SVijay Mahadevan } 381*ca4445c7SIlya Fursov values_phi_scalar2 = VecReal_to_VecScalar(values_phi); 382*ca4445c7SIlya Fursov PetscCall(MatSetValues(*interpl, 1, &dofsc[tc], connp.size(), &dofsp[0], &values_phi_scalar2[0], ADD_VALUES)); 383b117cd09SVijay Mahadevan } 384a86ed394SVijay Mahadevan } 385b117cd09SVijay Mahadevan } 386*ca4445c7SIlya Fursov if (vec) *vec = NULL; // TODO: <-- is it safe/appropriate? 3879566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*interpl, MAT_FINAL_ASSEMBLY)); 3889566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*interpl, MAT_FINAL_ASSEMBLY)); 3893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 390b117cd09SVijay Mahadevan } 391b117cd09SVijay Mahadevan 392a4e35b19SJacob Faibussowitsch // PetscClangLinter pragma ignore: -fdoc-* 393a4e35b19SJacob Faibussowitsch /* 394b117cd09SVijay Mahadevan DMCreateInjection_Moab - Generate a multi-level uniform refinement hierarchy 39520f4b53cSBarry Smith by succesively refining a coarse mesh, already defined in the `DM` object 396b117cd09SVijay Mahadevan provided by the user. 397b117cd09SVijay Mahadevan 398d083f849SBarry Smith Collective 399b117cd09SVijay Mahadevan 400b117cd09SVijay Mahadevan Input Parameter: 40120f4b53cSBarry Smith . dmb - The `DMMOAB` object 402b117cd09SVijay Mahadevan 403d8d19677SJose E. Roman Output Parameters: 404a2b725a8SWilliam Gropp + nlevels - The number of levels of refinement needed to generate the hierarchy 405a2b725a8SWilliam Gropp - ldegrees - The degree of refinement at each level in the hierarchy 406b117cd09SVijay Mahadevan 407b117cd09SVijay Mahadevan Level: beginner 408a4e35b19SJacob Faibussowitsch */ 409d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreateInjection_Moab(DM dm1, DM dm2, VecScatter *ctx) 410d71ae5a4SJacob Faibussowitsch { 411e882eb38SVijay Mahadevan //DM_Moab *dmmoab; 412b117cd09SVijay Mahadevan 413b117cd09SVijay Mahadevan PetscFunctionBegin; 414b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm1, DM_CLASSID, 1); 415b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm2, DM_CLASSID, 2); 416e882eb38SVijay Mahadevan //dmmoab = (DM_Moab*)(dm1)->data; 417b117cd09SVijay Mahadevan 4183ba16761SJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[DMCreateInjection_Moab] :: Placeholder\n")); 4193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 420b117cd09SVijay Mahadevan } 421b117cd09SVijay Mahadevan 422d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMMoab_UMR_Private(DM dm, MPI_Comm comm, PetscBool refine, DM *dmref) 423d71ae5a4SJacob Faibussowitsch { 424e882eb38SVijay Mahadevan PetscInt i, dim; 425b117cd09SVijay Mahadevan DM dm2; 426b117cd09SVijay Mahadevan moab::ErrorCode merr; 427b117cd09SVijay Mahadevan DM_Moab *dmb = (DM_Moab *)dm->data, *dd2; 428b117cd09SVijay Mahadevan 429b117cd09SVijay Mahadevan PetscFunctionBegin; 430b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4314f572ea9SToby Isaac PetscAssertPointer(dmref, 4); 432b117cd09SVijay Mahadevan 433b117cd09SVijay Mahadevan if ((dmb->hlevel == dmb->nhlevels && refine) || (dmb->hlevel == 0 && !refine)) { 4343ba16761SJacob Faibussowitsch if (dmb->hlevel + 1 > dmb->nhlevels && refine) { 4353ba16761SJacob 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)); 4363ba16761SJacob Faibussowitsch } 4373ba16761SJacob 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)); 438f3fa974cSJacob Faibussowitsch *dmref = NULL; 4393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 440b117cd09SVijay Mahadevan } 441b117cd09SVijay Mahadevan 4429566063dSJacob Faibussowitsch PetscCall(DMMoabCreate(PetscObjectComm((PetscObject)dm), &dm2)); 443b117cd09SVijay Mahadevan dd2 = (DM_Moab *)dm2->data; 444b117cd09SVijay Mahadevan 445b117cd09SVijay Mahadevan dd2->mbiface = dmb->mbiface; 4469daf19fdSVijay Mahadevan #ifdef MOAB_HAVE_MPI 447b117cd09SVijay Mahadevan dd2->pcomm = dmb->pcomm; 4489daf19fdSVijay Mahadevan #endif 449b117cd09SVijay Mahadevan dd2->icreatedinstance = PETSC_FALSE; 45064e1c140SVijay Mahadevan dd2->nghostrings = dmb->nghostrings; 451b117cd09SVijay Mahadevan 452b117cd09SVijay Mahadevan /* set the new level based on refinement/coarsening */ 453b117cd09SVijay Mahadevan if (refine) { 454b117cd09SVijay Mahadevan dd2->hlevel = dmb->hlevel + 1; 4551baa6e33SBarry Smith } else { 456b117cd09SVijay Mahadevan dd2->hlevel = dmb->hlevel - 1; 457b117cd09SVijay Mahadevan } 458b117cd09SVijay Mahadevan 459b117cd09SVijay Mahadevan /* Copy the multilevel hierarchy pointers in MOAB */ 460b117cd09SVijay Mahadevan dd2->hierarchy = dmb->hierarchy; 461b117cd09SVijay Mahadevan dd2->nhlevels = dmb->nhlevels; 4629566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dd2->nhlevels + 1, &dd2->hsets)); 463ad540459SPierre Jolivet for (i = 0; i <= dd2->nhlevels; i++) dd2->hsets[i] = dmb->hsets[i]; 464b117cd09SVijay Mahadevan dd2->fileset = dd2->hsets[dd2->hlevel]; 465b117cd09SVijay Mahadevan 466b117cd09SVijay Mahadevan /* do the remaining initializations for DMMoab */ 467b117cd09SVijay Mahadevan dd2->bs = dmb->bs; 468b117cd09SVijay Mahadevan dd2->numFields = dmb->numFields; 469b117cd09SVijay Mahadevan dd2->rw_dbglevel = dmb->rw_dbglevel; 470b117cd09SVijay Mahadevan dd2->partition_by_rank = dmb->partition_by_rank; 471c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(dd2->extra_read_options, dmb->extra_read_options, sizeof(dd2->extra_read_options))); 472c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(dd2->extra_write_options, dmb->extra_write_options, sizeof(dd2->extra_write_options))); 473b117cd09SVijay Mahadevan dd2->read_mode = dmb->read_mode; 474b117cd09SVijay Mahadevan dd2->write_mode = dmb->write_mode; 475b117cd09SVijay Mahadevan 476b117cd09SVijay Mahadevan /* set global ID tag handle */ 4779566063dSJacob Faibussowitsch PetscCall(DMMoabSetLocalToGlobalTag(dm2, dmb->ltog_tag)); 478b117cd09SVijay Mahadevan 4799371c9d4SSatish Balay merr = dd2->mbiface->tag_get_handle(MATERIAL_SET_TAG_NAME, dd2->material_tag); 4809371c9d4SSatish Balay MBERRNM(merr); 481b117cd09SVijay Mahadevan 4829566063dSJacob Faibussowitsch PetscCall(DMSetOptionsPrefix(dm2, ((PetscObject)dm)->prefix)); 4839566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 4849566063dSJacob Faibussowitsch PetscCall(DMSetDimension(dm2, dim)); 485b117cd09SVijay Mahadevan 486b117cd09SVijay Mahadevan /* allow overloaded (user replaced) operations to be inherited by refinement clones */ 487b117cd09SVijay Mahadevan dm2->ops->creatematrix = dm->ops->creatematrix; 488b117cd09SVijay Mahadevan 489b117cd09SVijay Mahadevan /* copy fill information if given */ 4909566063dSJacob Faibussowitsch PetscCall(DMMoabSetBlockFills(dm2, dmb->dfill, dmb->ofill)); 491b117cd09SVijay Mahadevan 492b117cd09SVijay Mahadevan /* copy vector type information */ 4939566063dSJacob Faibussowitsch PetscCall(DMSetMatType(dm2, dm->mattype)); 4949566063dSJacob Faibussowitsch PetscCall(DMSetVecType(dm2, dm->vectype)); 4953f1c6e43SVijay Mahadevan dd2->numFields = dmb->numFields; 49648a46eb9SPierre Jolivet if (dmb->numFields) PetscCall(DMMoabSetFieldNames(dm2, dmb->numFields, dmb->fieldNames)); 497b117cd09SVijay Mahadevan 4989566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm2)); 499b117cd09SVijay Mahadevan 500b117cd09SVijay Mahadevan /* recreate Dof numbering for the refined DM and make sure the distribution is correctly populated */ 5019566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm2)); 502b117cd09SVijay Mahadevan 503b117cd09SVijay Mahadevan *dmref = dm2; 5043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 505b117cd09SVijay Mahadevan } 506b117cd09SVijay Mahadevan 507a4e35b19SJacob Faibussowitsch // PetscClangLinter pragma ignore: -fdoc-* 508a4e35b19SJacob Faibussowitsch /* 509b117cd09SVijay Mahadevan DMRefine_Moab - Generate a multi-level uniform refinement hierarchy 51020f4b53cSBarry Smith by succesively refining a coarse mesh, already defined in the `DM` object 511b117cd09SVijay Mahadevan provided by the user. 512b117cd09SVijay Mahadevan 51320f4b53cSBarry Smith Collective 514b117cd09SVijay Mahadevan 515d8d19677SJose E. Roman Input Parameters: 51620f4b53cSBarry Smith + dm - The `DMMOAB` object 51720f4b53cSBarry Smith - comm - the communicator to contain the new DM object (or `MPI_COMM_NULL`) 518b117cd09SVijay Mahadevan 519b117cd09SVijay Mahadevan Output Parameter: 52020f4b53cSBarry Smith . dmf - the refined `DM`, or `NULL` 521e882eb38SVijay Mahadevan 522e882eb38SVijay Mahadevan Level: developer 523b117cd09SVijay Mahadevan 52420f4b53cSBarry Smith Note: 52520f4b53cSBarry Smith If no refinement was done, the return value is `NULL` 526a4e35b19SJacob Faibussowitsch */ 527d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMRefine_Moab(DM dm, MPI_Comm comm, DM *dmf) 528d71ae5a4SJacob Faibussowitsch { 529b117cd09SVijay Mahadevan PetscFunctionBegin; 530b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 531b117cd09SVijay Mahadevan 5329566063dSJacob Faibussowitsch PetscCall(DMMoab_UMR_Private(dm, comm, PETSC_TRUE, dmf)); 5333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 534b117cd09SVijay Mahadevan } 535b117cd09SVijay Mahadevan 536a4e35b19SJacob Faibussowitsch // PetscClangLinter pragma ignore: -fdoc-* 537a4e35b19SJacob Faibussowitsch /* 538b117cd09SVijay Mahadevan DMCoarsen_Moab - Generate a multi-level uniform refinement hierarchy 53920f4b53cSBarry Smith by succesively refining a coarse mesh, already defined in the `DM` object 540b117cd09SVijay Mahadevan provided by the user. 541b117cd09SVijay Mahadevan 54220f4b53cSBarry Smith Collective 543b117cd09SVijay Mahadevan 544d8d19677SJose E. Roman Input Parameters: 54520f4b53cSBarry Smith + dm - The `DMMOAB` object 54620f4b53cSBarry Smith - comm - the communicator to contain the new `DM` object (or `MPI_COMM_NULL`) 547b117cd09SVijay Mahadevan 548b117cd09SVijay Mahadevan Output Parameter: 54960225df5SJacob Faibussowitsch . dmc - the coarsened `DM`, or `NULL` 550b117cd09SVijay Mahadevan 551e882eb38SVijay Mahadevan Level: developer 552e882eb38SVijay Mahadevan 55320f4b53cSBarry Smith Note: 55420f4b53cSBarry Smith If no coarsening was done, the return value is `NULL` 555a4e35b19SJacob Faibussowitsch */ 556d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCoarsen_Moab(DM dm, MPI_Comm comm, DM *dmc) 557d71ae5a4SJacob Faibussowitsch { 558b117cd09SVijay Mahadevan PetscFunctionBegin; 559b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5609566063dSJacob Faibussowitsch PetscCall(DMMoab_UMR_Private(dm, comm, PETSC_FALSE, dmc)); 5613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 562b117cd09SVijay Mahadevan } 563