1b117cd09SVijay Mahadevan #include <petsc-private/dmmbimpl.h> /*I "petscdmmoab.h" I*/ 2b117cd09SVijay Mahadevan 3b117cd09SVijay Mahadevan #include <petscdmmoab.h> 4b117cd09SVijay Mahadevan #include <MBTagConventions.hpp> 5b117cd09SVijay Mahadevan #include <moab/NestedRefine.hpp> 6b117cd09SVijay Mahadevan 7b117cd09SVijay Mahadevan #undef __FUNCT__ 8b117cd09SVijay Mahadevan #define __FUNCT__ "DMMoabGenerateHierarchy" 9b117cd09SVijay Mahadevan /*@ 10b117cd09SVijay Mahadevan DMMoabGenerateHierarchy - Generate a multi-level uniform refinement hierarchy 11b117cd09SVijay Mahadevan by succesively refining a coarse mesh, already defined in the DM object 12b117cd09SVijay Mahadevan provided by the user. 13b117cd09SVijay Mahadevan 14b117cd09SVijay Mahadevan Collective on MPI_Comm 15b117cd09SVijay Mahadevan 16b117cd09SVijay Mahadevan Input Parameter: 17b117cd09SVijay Mahadevan + dmb - The DMMoab object 18b117cd09SVijay Mahadevan 19b117cd09SVijay Mahadevan Output Parameter: 20b117cd09SVijay Mahadevan + nlevels - The number of levels of refinement needed to generate the hierarchy 21b117cd09SVijay Mahadevan . ldegrees - The degree of refinement at each level in the hierarchy 22b117cd09SVijay Mahadevan 23b117cd09SVijay Mahadevan Level: beginner 24b117cd09SVijay Mahadevan 25b117cd09SVijay Mahadevan .keywords: DMMoab, create, refinement 26b117cd09SVijay Mahadevan @*/ 27b117cd09SVijay Mahadevan PetscErrorCode DMMoabGenerateHierarchy(DM dm,PetscInt nlevels,PetscInt *ldegrees) 28b117cd09SVijay Mahadevan { 29b117cd09SVijay Mahadevan DM_Moab *dmmoab; 30b117cd09SVijay Mahadevan PetscErrorCode ierr; 31b117cd09SVijay Mahadevan moab::ErrorCode merr; 32b117cd09SVijay Mahadevan PetscInt *pdegrees,i; 33e882eb38SVijay Mahadevan std::vector<moab::EntityHandle> hsets; 34b117cd09SVijay Mahadevan 35b117cd09SVijay Mahadevan PetscFunctionBegin; 36b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 37b117cd09SVijay Mahadevan dmmoab = (DM_Moab*)(dm)->data; 38b117cd09SVijay Mahadevan 39b117cd09SVijay Mahadevan if (!ldegrees) { 40b117cd09SVijay Mahadevan ierr = PetscMalloc1(nlevels,&pdegrees);CHKERRQ(ierr); 41b117cd09SVijay Mahadevan for (i=0; i<nlevels; i++) pdegrees[i]=2; /* default = Degree 2 refinement */ 42b117cd09SVijay Mahadevan } 43b117cd09SVijay Mahadevan else pdegrees=ldegrees; 44b117cd09SVijay Mahadevan 45b117cd09SVijay Mahadevan /* initialize set level refinement data for hierarchy */ 46b117cd09SVijay Mahadevan dmmoab->nhlevels=nlevels; 47b117cd09SVijay Mahadevan 48b117cd09SVijay Mahadevan /* Instantiate the nested refinement class */ 493f1c6e43SVijay Mahadevan dmmoab->hierarchy = new moab::NestedRefine(dynamic_cast<moab::Core*>(dmmoab->mbiface), dmmoab->pcomm, dmmoab->fileset); 50b117cd09SVijay Mahadevan 51e882eb38SVijay Mahadevan ierr = PetscMalloc1(nlevels+1,&dmmoab->hsets);CHKERRQ(ierr); 52e882eb38SVijay Mahadevan 53b117cd09SVijay Mahadevan /* generate the mesh hierarchy */ 54e882eb38SVijay Mahadevan merr = dmmoab->hierarchy->generate_mesh_hierarchy(pdegrees, nlevels, hsets);MBERRNM(merr); 55e882eb38SVijay Mahadevan 56e882eb38SVijay Mahadevan for (i=0; i<=nlevels; i++) dmmoab->hsets[i]=hsets[i]; 57b117cd09SVijay Mahadevan 58b117cd09SVijay Mahadevan if (!ldegrees) { 59b117cd09SVijay Mahadevan ierr = PetscFree(pdegrees);CHKERRQ(ierr); 60b117cd09SVijay Mahadevan } 61b117cd09SVijay Mahadevan PetscFunctionReturn(0); 62b117cd09SVijay Mahadevan } 63b117cd09SVijay Mahadevan 64b117cd09SVijay Mahadevan #undef __FUNCT__ 65b117cd09SVijay Mahadevan #define __FUNCT__ "DMRefineHierarchy_Moab" 66b117cd09SVijay Mahadevan /*@ 67e882eb38SVijay Mahadevan DMRefineHierarchy_Moab - Generate a multi-level DM hierarchy 68e882eb38SVijay Mahadevan by succesively refining a coarse mesh. 69b117cd09SVijay Mahadevan 70b117cd09SVijay Mahadevan Collective on MPI_Comm 71b117cd09SVijay Mahadevan 72b117cd09SVijay Mahadevan Input Parameter: 73e882eb38SVijay Mahadevan + dm - The DMMoab object 74b117cd09SVijay Mahadevan 75b117cd09SVijay Mahadevan Output Parameter: 76e882eb38SVijay Mahadevan + nlevels - The number of levels of refinement needed to generate the hierarchy 77e882eb38SVijay Mahadevan . dmf - The DM objects after successive refinement of the hierarchy 78b117cd09SVijay Mahadevan 79b117cd09SVijay Mahadevan Level: beginner 80b117cd09SVijay Mahadevan 81e882eb38SVijay Mahadevan .keywords: DMMoab, generate, refinement 82b117cd09SVijay Mahadevan @*/ 83b117cd09SVijay Mahadevan PetscErrorCode DMRefineHierarchy_Moab(DM dm,PetscInt nlevels,DM dmf[]) 84b117cd09SVijay Mahadevan { 85b117cd09SVijay Mahadevan PetscErrorCode ierr; 86b117cd09SVijay Mahadevan PetscInt i; 87b117cd09SVijay Mahadevan 88b117cd09SVijay Mahadevan PetscFunctionBegin; 89b117cd09SVijay Mahadevan 90e882eb38SVijay Mahadevan ierr = DMRefine(dm,PetscObjectComm((PetscObject)dm),&dmf[0]);CHKERRQ(ierr); 91e882eb38SVijay Mahadevan for (i=1; i<nlevels; i++) { 92e882eb38SVijay Mahadevan ierr = DMRefine(dmf[i-1],PetscObjectComm((PetscObject)dm),&dmf[i]);CHKERRQ(ierr); 93b117cd09SVijay Mahadevan } 94b117cd09SVijay Mahadevan PetscFunctionReturn(0); 95b117cd09SVijay Mahadevan } 96b117cd09SVijay Mahadevan 97b117cd09SVijay Mahadevan #undef __FUNCT__ 98b117cd09SVijay Mahadevan #define __FUNCT__ "DMCoarsenHierarchy_Moab" 99b117cd09SVijay Mahadevan /*@ 100e882eb38SVijay Mahadevan DMCoarsenHierarchy_Moab - Generate a multi-level DM hierarchy 101e882eb38SVijay Mahadevan by succesively coarsening a refined mesh. 102b117cd09SVijay Mahadevan 103b117cd09SVijay Mahadevan Collective on MPI_Comm 104b117cd09SVijay Mahadevan 105b117cd09SVijay Mahadevan Input Parameter: 106e882eb38SVijay Mahadevan + dm - The DMMoab object 107b117cd09SVijay Mahadevan 108b117cd09SVijay Mahadevan Output Parameter: 109e882eb38SVijay Mahadevan + nlevels - The number of levels of refinement needed to generate the hierarchy 110e882eb38SVijay Mahadevan . dmc - The DM objects after successive coarsening of the hierarchy 111b117cd09SVijay Mahadevan 112b117cd09SVijay Mahadevan Level: beginner 113b117cd09SVijay Mahadevan 114e882eb38SVijay Mahadevan .keywords: DMMoab, generate, coarsening 115b117cd09SVijay Mahadevan @*/ 116b117cd09SVijay Mahadevan PetscErrorCode DMCoarsenHierarchy_Moab(DM dm,PetscInt nlevels,DM dmc[]) 117b117cd09SVijay Mahadevan { 118b117cd09SVijay Mahadevan PetscErrorCode ierr; 119b117cd09SVijay Mahadevan PetscInt i; 120b117cd09SVijay Mahadevan 121b117cd09SVijay Mahadevan PetscFunctionBegin; 122b117cd09SVijay Mahadevan 123e882eb38SVijay Mahadevan ierr = DMCoarsen(dm,PetscObjectComm((PetscObject)dm),&dmc[0]);CHKERRQ(ierr); 124e882eb38SVijay Mahadevan for (i=1; i<nlevels; i++) { 125e882eb38SVijay Mahadevan ierr = DMCoarsen(dmc[i-1],PetscObjectComm((PetscObject)dm),&dmc[i]);CHKERRQ(ierr); 126b117cd09SVijay Mahadevan } 127b117cd09SVijay Mahadevan PetscFunctionReturn(0); 128b117cd09SVijay Mahadevan } 129b117cd09SVijay Mahadevan 130b117cd09SVijay Mahadevan 131b117cd09SVijay Mahadevan #undef __FUNCT__ 132b117cd09SVijay Mahadevan #define __FUNCT__ "DMCreateInterpolation_Moab" 133b117cd09SVijay Mahadevan /*@ 134e882eb38SVijay Mahadevan DMCreateInterpolation_Moab - Generate the interpolation operators to transform 135e882eb38SVijay Mahadevan operators (matrices, vectors) from parent level to child level as defined by 136e882eb38SVijay Mahadevan the DM inputs provided by the user. 137b117cd09SVijay Mahadevan 138b117cd09SVijay Mahadevan Collective on MPI_Comm 139b117cd09SVijay Mahadevan 140b117cd09SVijay Mahadevan Input Parameter: 141e882eb38SVijay Mahadevan + dm1 - The DMMoab object 142e882eb38SVijay Mahadevan - dm2 - the second, finer DMMoab object 143b117cd09SVijay Mahadevan 144b117cd09SVijay Mahadevan Output Parameter: 145e882eb38SVijay Mahadevan + interpl - The interpolation operator for transferring data between the levels 146e882eb38SVijay Mahadevan - vec - The scaling vector (optional) 147b117cd09SVijay Mahadevan 148e882eb38SVijay Mahadevan Level: developer 149b117cd09SVijay Mahadevan 150b117cd09SVijay Mahadevan .keywords: DMMoab, create, refinement 151b117cd09SVijay Mahadevan @*/ 152b117cd09SVijay Mahadevan PetscErrorCode DMCreateInterpolation_Moab(DM dm1,DM dm2,Mat* interpl,Vec* vec) 153b117cd09SVijay Mahadevan { 154b117cd09SVijay Mahadevan DM_Moab *dmb1, *dmb2; 155b117cd09SVijay Mahadevan PetscErrorCode ierr; 156b117cd09SVijay Mahadevan moab::ErrorCode merr; 157e882eb38SVijay Mahadevan PetscInt dim; 1583f1c6e43SVijay Mahadevan PetscReal factor; 1593f1c6e43SVijay Mahadevan PetscBool eonbnd; 160*ce27a4eeSVijay Mahadevan PetscInt innz, *nnz, ionz, *onz; 161*ce27a4eeSVijay Mahadevan PetscInt nlsiz1, nlsiz2, ngsiz1, ngsiz2; 162b117cd09SVijay Mahadevan std::vector<int> bndrows; 1633f1c6e43SVijay Mahadevan std::vector<PetscBool> dbdry; 164b117cd09SVijay Mahadevan 165b117cd09SVijay Mahadevan PetscFunctionBegin; 166b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm1,DM_CLASSID,1); 167b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm2,DM_CLASSID,2); 168b117cd09SVijay Mahadevan dmb1 = (DM_Moab*)(dm1)->data; 169b117cd09SVijay Mahadevan dmb2 = (DM_Moab*)(dm2)->data; 170941e0cffSVijay Mahadevan nlsiz1 = dmb1->nloc*dmb2->numFields; 171941e0cffSVijay Mahadevan nlsiz2 = dmb2->nloc*dmb2->numFields; 172941e0cffSVijay Mahadevan ngsiz1 = dmb1->n*dmb2->numFields; 173941e0cffSVijay Mahadevan ngsiz2 = dmb2->n*dmb2->numFields; 174b117cd09SVijay Mahadevan 1753f1c6e43SVijay Mahadevan PetscInfo4(dm1,"Creating interpolation matrix %D X %D to apply transformation between levels %D -> %D.\n",dmb1->nloc,dmb2->nloc,dmb1->hlevel,dmb2->hlevel); 1763f1c6e43SVijay Mahadevan 177b117cd09SVijay Mahadevan ierr = MatCreate(PetscObjectComm((PetscObject)dm1), interpl);CHKERRQ(ierr); 178e882eb38SVijay Mahadevan ierr = MatSetType(*interpl, dm1->mattype);CHKERRQ(ierr); 179941e0cffSVijay Mahadevan ierr = MatSetSizes(*interpl, nlsiz2, nlsiz1, ngsiz2, ngsiz1);CHKERRQ(ierr); 180b117cd09SVijay Mahadevan 181b117cd09SVijay Mahadevan /* TODO: This is a hack for the rectangular system - decipher NNZ pattern better */ 182b117cd09SVijay Mahadevan ierr = MatSetFromOptions(*interpl);CHKERRQ(ierr); 183b117cd09SVijay Mahadevan 184941e0cffSVijay Mahadevan /* allocate the nnz, onz arrays based on block size and local nodes */ 185941e0cffSVijay Mahadevan ierr = PetscCalloc2(nlsiz2,&nnz,nlsiz2,&onz);CHKERRQ(ierr); 186941e0cffSVijay Mahadevan 187941e0cffSVijay Mahadevan /* Loop through the local elements and compute the relation between the current parent and the refined_level. */ 188941e0cffSVijay Mahadevan for(moab::Range::iterator iter = dmb1->elocal->begin(); iter!= dmb1->elocal->end(); iter++) { 189941e0cffSVijay Mahadevan 190941e0cffSVijay Mahadevan const moab::EntityHandle ehandle = *iter; 191941e0cffSVijay Mahadevan std::vector<moab::EntityHandle> children; 192941e0cffSVijay Mahadevan std::vector<moab::EntityHandle> connp, connc; 193941e0cffSVijay Mahadevan 194941e0cffSVijay Mahadevan /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */ 195941e0cffSVijay Mahadevan merr = dmb1->hierarchy->parent_to_child(ehandle, dmb1->hlevel, dmb2->hlevel, children);MBERRNM(merr); 196941e0cffSVijay Mahadevan 197941e0cffSVijay Mahadevan /* Get connectivity and coordinates of the parent vertices */ 198941e0cffSVijay Mahadevan merr = dmb1->hierarchy->get_connectivity(ehandle, dmb1->hlevel, connp);MBERRNM(merr); 199941e0cffSVijay Mahadevan for (unsigned ic=0; ic < children.size(); ic++) { 200941e0cffSVijay Mahadevan std::vector<moab::EntityHandle> tconnc; 201941e0cffSVijay Mahadevan /* Get coordinates of the parent vertices in canonical order */ 202941e0cffSVijay Mahadevan merr = dmb2->hierarchy->get_connectivity(children[ic], dmb2->hlevel, tconnc);MBERRNM(merr); 203941e0cffSVijay Mahadevan for (unsigned tc=0; tc<tconnc.size(); tc++) { 204941e0cffSVijay Mahadevan if (std::find(connc.begin(), connc.end(), tconnc[tc]) == connc.end()) 205941e0cffSVijay Mahadevan connc.push_back(tconnc[tc]); 206941e0cffSVijay Mahadevan } 207941e0cffSVijay Mahadevan } 208941e0cffSVijay Mahadevan 209941e0cffSVijay Mahadevan std::vector<int> dofsp(connp.size()), dofsc(connc.size()); 210941e0cffSVijay Mahadevan /* TODO: specific to scalar system - use GetDofs */ 211941e0cffSVijay Mahadevan //ierr = DMMoabGetFieldDofs(dm1, connp.size(), &connp[0], 0, &dofsp[0]);CHKERRQ(ierr); 212941e0cffSVijay Mahadevan ierr = DMMoabGetFieldDofs(dm2, connc.size(), &connc[0], 0, &dofsc[0]);CHKERRQ(ierr); 213941e0cffSVijay Mahadevan 214941e0cffSVijay Mahadevan for (unsigned tp=0;tp<connp.size(); tp++) { 215941e0cffSVijay Mahadevan // ierr = MatSetValues(*interpl, connc.size(), &dofsc[0], 1, &dofsp[i], &values_phi[0], ADD_VALUES);CHKERRQ(ierr); 216941e0cffSVijay Mahadevan if (dmb1->vowned->find(connp[tp]) != dmb1->vowned->end()) { 217941e0cffSVijay Mahadevan for (unsigned tc=0;tc<connc.size(); tc++) nnz[dofsc[tc]]++; 218941e0cffSVijay Mahadevan } 219941e0cffSVijay Mahadevan else if (dmb1->vghost->find(connp[tp]) != dmb1->vghost->end()) { 220941e0cffSVijay Mahadevan for (unsigned tc=0;tc<connc.size(); tc++) onz[dofsc[tc]]++; 221941e0cffSVijay Mahadevan } 222941e0cffSVijay Mahadevan else SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Invalid entity in parent level %D\n", connc[tp]); 223941e0cffSVijay Mahadevan } 224941e0cffSVijay Mahadevan /* 225941e0cffSVijay Mahadevan for(int tc = 0; tc < connc.size(); tc++) { 226941e0cffSVijay Mahadevan if (dmb2->vowned->find(connc[tc]) != dmb2->vowned->end()) nnz[dofsc[tc]]++; 227941e0cffSVijay Mahadevan else if (dmb2->vghost->find(connc[tc]) != dmb2->vghost->end()) onz[dofsc[tc]]++; 228941e0cffSVijay Mahadevan else SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Invalid entity in child level %D\n", connc[tc]); 229941e0cffSVijay Mahadevan }*/ 230941e0cffSVijay Mahadevan } 231941e0cffSVijay Mahadevan 232941e0cffSVijay Mahadevan int i=0; 233941e0cffSVijay Mahadevan std::vector<moab::EntityHandle> adjs; 234941e0cffSVijay Mahadevan for(moab::Range::iterator iter = dmb1->vowned->begin(); iter!= dmb1->vowned->end(); iter++, i++) { 235941e0cffSVijay Mahadevan merr = dmb1->hierarchy->get_adjacencies(*iter, 0, adjs);MBERRNM(merr); 236941e0cffSVijay Mahadevan nnz[i] -= adjs.size(); 237941e0cffSVijay Mahadevan adjs.clear(); 238941e0cffSVijay Mahadevan } 239941e0cffSVijay Mahadevan for(moab::Range::iterator iter = dmb1->vowned->begin(); iter!= dmb1->vowned->end(); iter++, i++) { 240941e0cffSVijay Mahadevan //merr = dmb1->hierarchy->get_adjacencies(&(*iter), 1, 0, false, adjs, moab::Interface::UNION);MBERRNM(merr); 241941e0cffSVijay Mahadevan merr = dmb1->hierarchy->get_adjacencies(*iter, 0, adjs);MBERRNM(merr); 242941e0cffSVijay Mahadevan onz[i] -= adjs.size(); 243941e0cffSVijay Mahadevan adjs.clear(); 244941e0cffSVijay Mahadevan } 245941e0cffSVijay Mahadevan 246941e0cffSVijay Mahadevan ionz=onz[0]; 247941e0cffSVijay Mahadevan innz=nnz[0]; 248941e0cffSVijay Mahadevan for (int tc=0; tc < nlsiz2; tc++) { 249941e0cffSVijay Mahadevan // check for maximum allowed sparsity = fully dense 250941e0cffSVijay Mahadevan nnz[tc] = std::min(nlsiz1,nnz[tc]); 251941e0cffSVijay Mahadevan onz[tc] = std::min(nlsiz1,onz[tc]); 252941e0cffSVijay Mahadevan 253941e0cffSVijay Mahadevan innz = (innz < nnz[tc] ? nnz[tc] : innz); 254941e0cffSVijay Mahadevan ionz = (ionz < onz[tc] ? onz[tc] : ionz); 255*ce27a4eeSVijay Mahadevan //PetscPrintf(PETSC_COMM_SELF, "[%D] NNZ = %D, ONZ = %D\n", tc, nnz[tc], onz[tc]); 256941e0cffSVijay Mahadevan } 257*ce27a4eeSVijay Mahadevan //PetscPrintf(PETSC_COMM_SELF, "Final: INNZ = %D, IONZ = %D\n", innz, ionz); 258941e0cffSVijay Mahadevan 259941e0cffSVijay Mahadevan ierr = MatSeqAIJSetPreallocation(*interpl,innz,nnz);CHKERRQ(ierr); 260941e0cffSVijay Mahadevan ierr = MatMPIAIJSetPreallocation(*interpl,innz,nnz,ionz,onz);CHKERRQ(ierr); 261941e0cffSVijay Mahadevan 262941e0cffSVijay Mahadevan /* clean up temporary memory */ 263941e0cffSVijay Mahadevan ierr = PetscFree2(nnz,onz);CHKERRQ(ierr); 264b117cd09SVijay Mahadevan 265b117cd09SVijay Mahadevan /* set up internal matrix data-structures */ 266b117cd09SVijay Mahadevan ierr = MatSetUp(*interpl);CHKERRQ(ierr); 267b117cd09SVijay Mahadevan ierr = MatZeroEntries(*interpl);CHKERRQ(ierr); 268b117cd09SVijay Mahadevan 269b117cd09SVijay Mahadevan ierr = DMGetDimension(dm1, &dim);CHKERRQ(ierr); 270b117cd09SVijay Mahadevan 2713f1c6e43SVijay Mahadevan factor = std::pow(2.0 /*degree_P_for_refinement*/,(dmb2->hlevel-dmb1->hlevel)*dmb1->dim*1.0); 272b117cd09SVijay Mahadevan 273e882eb38SVijay Mahadevan /* Loop through the remaining vertices. These vertices appear only on the current refined_level. */ 274b117cd09SVijay Mahadevan for(moab::Range::iterator iter = dmb1->elocal->begin(); iter!= dmb1->elocal->end(); iter++) { 275b117cd09SVijay Mahadevan 276b117cd09SVijay Mahadevan const moab::EntityHandle ehandle = *iter; 277b117cd09SVijay Mahadevan std::vector<moab::EntityHandle> children; 278b117cd09SVijay Mahadevan std::vector<moab::EntityHandle> connp, connc; 279b117cd09SVijay Mahadevan 280b117cd09SVijay Mahadevan /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */ 281b117cd09SVijay Mahadevan merr = dmb1->hierarchy->parent_to_child(ehandle, dmb1->hlevel, dmb2->hlevel, children);MBERRNM(merr); 282b117cd09SVijay Mahadevan 283b117cd09SVijay Mahadevan /* Get connectivity and coordinates of the parent vertices */ 284b117cd09SVijay Mahadevan merr = dmb1->hierarchy->get_connectivity(ehandle, dmb1->hlevel, connp);MBERRNM(merr); 285e882eb38SVijay Mahadevan for (unsigned ic=0; ic < children.size(); ic++) { 286b117cd09SVijay Mahadevan std::vector<moab::EntityHandle> tconnc; 287b117cd09SVijay Mahadevan /* Get coordinates of the parent vertices in canonical order */ 288b117cd09SVijay Mahadevan merr = dmb1->hierarchy->get_connectivity(children[ic], dmb2->hlevel, tconnc);MBERRNM(merr); 289e882eb38SVijay Mahadevan for (unsigned tc=0; tc<tconnc.size(); tc++) { 290b117cd09SVijay Mahadevan connc.push_back(tconnc[tc]); 291b117cd09SVijay Mahadevan } 292b117cd09SVijay Mahadevan } 293b117cd09SVijay Mahadevan 294b117cd09SVijay Mahadevan std::vector<double> pcoords(connp.size()*3), ccoords(connc.size()*3), values_phi(connc.size()); 295b117cd09SVijay Mahadevan /* Get coordinates for connectivity entities in canonical order for both coarse and finer levels */ 296b117cd09SVijay Mahadevan merr = dmb1->hierarchy->get_coordinates(&connp[0], connp.size(), dmb1->hlevel, &pcoords[0]);MBERRNM(merr); 297b117cd09SVijay Mahadevan merr = dmb2->hierarchy->get_coordinates(&connc[0], connc.size(), dmb2->hlevel, &ccoords[0]);MBERRNM(merr); 298b117cd09SVijay Mahadevan 299b117cd09SVijay Mahadevan std::vector<int> dofsp(connp.size()), dofsc(connc.size()); 300b117cd09SVijay Mahadevan /* TODO: specific to scalar system - use GetDofs */ 301b117cd09SVijay Mahadevan ierr = DMMoabGetFieldDofs(dm1, connp.size(), &connp[0], 0, &dofsp[0]);CHKERRQ(ierr); 302b117cd09SVijay Mahadevan ierr = DMMoabGetFieldDofs(dm2, connc.size(), &connc[0], 0, &dofsc[0]);CHKERRQ(ierr); 303b117cd09SVijay Mahadevan 304e882eb38SVijay Mahadevan /* Compute the interpolation weights by determining distance of 1-ring 305e882eb38SVijay Mahadevan neighbor vertices from current vertex */ 306e882eb38SVijay Mahadevan for (unsigned i=0;i<connp.size(); i++) { 307b117cd09SVijay Mahadevan double normsum=0.0; 308e882eb38SVijay Mahadevan for (unsigned j=0;j<connc.size(); j++) { 3093f1c6e43SVijay Mahadevan values_phi[j] = 0.0; 310e882eb38SVijay Mahadevan for (unsigned k=0;k<3; k++) 311941e0cffSVijay Mahadevan values_phi[j] += std::pow(pcoords[i*3+k]-ccoords[k+j*3], dim); 3123f1c6e43SVijay Mahadevan if (values_phi[j] < 1e-12) { 3133f1c6e43SVijay Mahadevan values_phi[j] = 1e12; 314b117cd09SVijay Mahadevan } 315b117cd09SVijay Mahadevan else { 316941e0cffSVijay Mahadevan //values_phi[j] = std::pow(values_phi[j], -1.0/dim); 317941e0cffSVijay Mahadevan values_phi[j] = std::pow(values_phi[j], -1.0); 3183f1c6e43SVijay Mahadevan normsum += values_phi[j]; 319b117cd09SVijay Mahadevan } 320b117cd09SVijay Mahadevan } 321e882eb38SVijay Mahadevan for (unsigned j=0;j<connc.size(); j++) { 3223f1c6e43SVijay Mahadevan if (values_phi[j] > 1e11) 3233f1c6e43SVijay Mahadevan values_phi[j] = factor*0.5/connc.size(); 324b117cd09SVijay Mahadevan else 3253f1c6e43SVijay Mahadevan values_phi[j] = factor*values_phi[j]*0.5/(connc.size()*normsum); 326b117cd09SVijay Mahadevan } 327b117cd09SVijay Mahadevan ierr = MatSetValues(*interpl, connc.size(), &dofsc[0], 1, &dofsp[i], &values_phi[0], ADD_VALUES);CHKERRQ(ierr); 328b117cd09SVijay Mahadevan } 329b117cd09SVijay Mahadevan 330b117cd09SVijay Mahadevan /* check if element is on the boundary */ 331b117cd09SVijay Mahadevan //ierr = DMMoabIsEntityOnBoundary(dm1,ehandle,&eonbnd);CHKERRQ(ierr); 3323f1c6e43SVijay Mahadevan dbdry.resize(connc.size()); 3333f1c6e43SVijay Mahadevan ierr = DMMoabCheckBoundaryVertices(dm2,connc.size(),&connc[0],dbdry.data());CHKERRQ(ierr); 334b117cd09SVijay Mahadevan eonbnd=PETSC_FALSE; 335e882eb38SVijay Mahadevan for (unsigned i=0; i< connc.size(); ++i) 336b117cd09SVijay Mahadevan if (dbdry[i]) eonbnd=PETSC_TRUE; 337b117cd09SVijay Mahadevan 338b117cd09SVijay Mahadevan values_phi.clear(); 339b117cd09SVijay Mahadevan values_phi.resize(connp.size()); 340b117cd09SVijay Mahadevan /* apply dirichlet boundary conditions */ 341b117cd09SVijay Mahadevan if (eonbnd) { 342b117cd09SVijay Mahadevan 343b117cd09SVijay Mahadevan ierr = MatAssemblyBegin(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr); 344b117cd09SVijay Mahadevan ierr = MatAssemblyEnd(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr); 345b117cd09SVijay Mahadevan /* get the list of nodes on boundary so that we can enforce dirichlet conditions strongly */ 346b117cd09SVijay Mahadevan //ierr = DMMoabCheckBoundaryVertices(dm2,connc.size(),&connc[0],dbdry);CHKERRQ(ierr); 347e882eb38SVijay Mahadevan for (unsigned i=0; i < connc.size(); i++) { 3483f1c6e43SVijay Mahadevan if (dbdry[i]) { /* dirichlet node */ 349b117cd09SVijay Mahadevan /* think about strongly imposing dirichlet */ 350b117cd09SVijay Mahadevan //bndrows.push_back(dofsc[i]); 351b117cd09SVijay Mahadevan 352b117cd09SVijay Mahadevan ierr = MatSetValues(*interpl, 1, &dofsc[i], connp.size(), &dofsp[0], &values_phi[0], INSERT_VALUES);CHKERRQ(ierr); 353b117cd09SVijay Mahadevan //values_phi[0]=1.0; 354b117cd09SVijay Mahadevan //ierr = MatSetValues(*interpl, 1, &dofsc[i], 1, &dofsc[i], &values_phi[0], INSERT_VALUES);CHKERRQ(ierr); 355b117cd09SVijay Mahadevan } 356b117cd09SVijay Mahadevan } 357b117cd09SVijay Mahadevan 358b117cd09SVijay Mahadevan ierr = MatAssemblyBegin(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr); 359b117cd09SVijay Mahadevan ierr = MatAssemblyEnd(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr); 360b117cd09SVijay Mahadevan } 361b117cd09SVijay Mahadevan 362b117cd09SVijay Mahadevan //get interpolation weights 363b117cd09SVijay Mahadevan //ierr = Compute_Quad4_Basis(pcoords, 1, coord, values_phi);CHKERRQ(ierr); 364e882eb38SVijay Mahadevan // for (int j=0;j<dofs_per_element; j++) 365e882eb38SVijay Mahadevan // std::cout<<"values "<<values_phi[j]<<std::endl; 366b117cd09SVijay Mahadevan 367b117cd09SVijay Mahadevan //get row and column indices, zero weights are ignored 368b117cd09SVijay Mahadevan /* 369b117cd09SVijay Mahadevan int nz_ind = 0; 370b117cd09SVijay Mahadevan idx = dmb2->vowned->index(vhandle); 371b117cd09SVijay Mahadevan for (int j=0;j<dofs_per_element; j++){ 372b117cd09SVijay Mahadevan idy[nz_ind] = dmb1->vowned->index(connectivity[j]); 373b117cd09SVijay Mahadevan PetscPrintf(PETSC_COMM_WORLD, "Finding coarse connectivity vertex %D associated with [%D, %D] - set to %D\n", connectivity[j], parent.size(), vhandle, idy[nz_ind]); 374b117cd09SVijay Mahadevan //values_phi[nz_ind] = values_phi[j]; 375b117cd09SVijay Mahadevan nz_ind = nz_ind+1; 376b117cd09SVijay Mahadevan } 377b117cd09SVijay Mahadevan */ 378b117cd09SVijay Mahadevan 379b117cd09SVijay Mahadevan //ierr = MatSetValues(*interpl, nz_ind, idy, 1, &idx, values_phi, INSERT_VALUES);CHKERRQ(ierr); 380b117cd09SVijay Mahadevan //ierr = MatSetValues(*interpl, connp.size(), dofsp, connc.size(), dofsc, &values_phi[0], INSERT_VALUES);CHKERRQ(ierr); 381b117cd09SVijay Mahadevan } 382b117cd09SVijay Mahadevan 383b117cd09SVijay Mahadevan //PetscPrintf(PETSC_COMM_WORLD, "[Boundary vertices = %D] :: A few: %D %D %D %D \n", bndrows.size(), bndrows[0], bndrows[1], bndrows[2], bndrows[3]); 384b117cd09SVijay Mahadevan //ierr = MatZeroRows(*interpl, bndrows.size(), &bndrows[0], 1.0, 0, 0);CHKERRQ(ierr); 385b117cd09SVijay Mahadevan 386b117cd09SVijay Mahadevan ierr = MatAssemblyBegin(*interpl,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 387b117cd09SVijay Mahadevan ierr = MatAssemblyEnd(*interpl,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 388b117cd09SVijay Mahadevan PetscFunctionReturn(0); 389b117cd09SVijay Mahadevan } 390b117cd09SVijay Mahadevan 391b117cd09SVijay Mahadevan #undef __FUNCT__ 392b117cd09SVijay Mahadevan #define __FUNCT__ "DMCreateInjection_Moab" 393b117cd09SVijay Mahadevan /*@ 394b117cd09SVijay Mahadevan DMCreateInjection_Moab - Generate a multi-level uniform refinement hierarchy 395b117cd09SVijay Mahadevan by succesively refining a coarse mesh, already defined in the DM object 396b117cd09SVijay Mahadevan provided by the user. 397b117cd09SVijay Mahadevan 398b117cd09SVijay Mahadevan Collective on MPI_Comm 399b117cd09SVijay Mahadevan 400b117cd09SVijay Mahadevan Input Parameter: 401b117cd09SVijay Mahadevan . dmb - The DMMoab object 402b117cd09SVijay Mahadevan 403b117cd09SVijay Mahadevan Output Parameter: 404b117cd09SVijay Mahadevan . nlevels - The number of levels of refinement needed to generate the hierarchy 405b117cd09SVijay Mahadevan + ldegrees - The degree of refinement at each level in the hierarchy 406b117cd09SVijay Mahadevan 407b117cd09SVijay Mahadevan Level: beginner 408b117cd09SVijay Mahadevan 409b117cd09SVijay Mahadevan .keywords: DMMoab, create, refinement 410b117cd09SVijay Mahadevan @*/ 411b117cd09SVijay Mahadevan PetscErrorCode DMCreateInjection_Moab(DM dm1,DM dm2,VecScatter* ctx) 412b117cd09SVijay Mahadevan { 413e882eb38SVijay Mahadevan //DM_Moab *dmmoab; 414b117cd09SVijay Mahadevan 415b117cd09SVijay Mahadevan PetscFunctionBegin; 416b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm1,DM_CLASSID,1); 417b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm2,DM_CLASSID,2); 418e882eb38SVijay Mahadevan //dmmoab = (DM_Moab*)(dm1)->data; 419b117cd09SVijay Mahadevan 420b117cd09SVijay Mahadevan PetscPrintf(PETSC_COMM_WORLD, "[DMCreateInjection_Moab] :: Placeholder\n"); 421b117cd09SVijay Mahadevan PetscFunctionReturn(0); 422b117cd09SVijay Mahadevan } 423b117cd09SVijay Mahadevan 424b117cd09SVijay Mahadevan 425b117cd09SVijay Mahadevan #undef __FUNCT__ 426b117cd09SVijay Mahadevan #define __FUNCT__ "DM_UMR_Moab_Private" 427b117cd09SVijay Mahadevan PetscErrorCode DM_UMR_Moab_Private(DM dm,MPI_Comm comm,PetscBool refine,DM *dmref) 428b117cd09SVijay Mahadevan { 429b117cd09SVijay Mahadevan PetscErrorCode ierr; 430e882eb38SVijay Mahadevan PetscInt i,dim; 431b117cd09SVijay Mahadevan DM dm2; 432b117cd09SVijay Mahadevan moab::ErrorCode merr; 433b117cd09SVijay Mahadevan DM_Moab *dmb = (DM_Moab*)dm->data,*dd2; 434b117cd09SVijay Mahadevan 435b117cd09SVijay Mahadevan PetscFunctionBegin; 436b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 437b117cd09SVijay Mahadevan PetscValidPointer(dmref,4); 438b117cd09SVijay Mahadevan 439b117cd09SVijay Mahadevan if ( (dmb->hlevel == dmb->nhlevels && refine) || (dmb->hlevel == 0 && !refine) ) { 440e882eb38SVijay 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); 441e882eb38SVijay Mahadevan if (dmb->hlevel-1 < 0 && !refine) PetscInfo1(NULL,"Invalid multigrid coarsen hierarchy level specified (%D). Creating a NULL object.\n",dmb->hlevel-1); 442e882eb38SVijay Mahadevan *dmref = PETSC_NULL; 443b117cd09SVijay Mahadevan PetscFunctionReturn(0); 444b117cd09SVijay Mahadevan } 445b117cd09SVijay Mahadevan 446b117cd09SVijay Mahadevan ierr = DMMoabCreate(PetscObjectComm((PetscObject)dm), &dm2);CHKERRQ(ierr); 447b117cd09SVijay Mahadevan dd2 = (DM_Moab*)dm2->data; 448b117cd09SVijay Mahadevan 449b117cd09SVijay Mahadevan dd2->mbiface = dmb->mbiface; 450b117cd09SVijay Mahadevan dd2->pcomm = dmb->pcomm; 451b117cd09SVijay Mahadevan dd2->icreatedinstance = PETSC_FALSE; 452b117cd09SVijay Mahadevan 453b117cd09SVijay Mahadevan /* set the new level based on refinement/coarsening */ 454b117cd09SVijay Mahadevan if (refine) { 455b117cd09SVijay Mahadevan dd2->hlevel=dmb->hlevel+1; 456b117cd09SVijay Mahadevan } 457b117cd09SVijay Mahadevan else { 458b117cd09SVijay Mahadevan dd2->hlevel=dmb->hlevel-1; 459b117cd09SVijay Mahadevan } 460b117cd09SVijay Mahadevan 461b117cd09SVijay Mahadevan /* Copy the multilevel hierarchy pointers in MOAB */ 462b117cd09SVijay Mahadevan dd2->hierarchy = dmb->hierarchy; 463b117cd09SVijay Mahadevan dd2->nhlevels = dmb->nhlevels; 464b117cd09SVijay Mahadevan ierr = PetscMalloc1(dd2->nhlevels+1,&dd2->hsets);CHKERRQ(ierr); 465b117cd09SVijay Mahadevan for (i=0; i<=dd2->nhlevels; i++) { 466b117cd09SVijay Mahadevan dd2->hsets[i]=dmb->hsets[i]; 467b117cd09SVijay Mahadevan } 468b117cd09SVijay Mahadevan dd2->fileset = dd2->hsets[dd2->hlevel]; 469b117cd09SVijay Mahadevan 470b117cd09SVijay Mahadevan /* do the remaining initializations for DMMoab */ 471b117cd09SVijay Mahadevan dd2->bs = dmb->bs; 472b117cd09SVijay Mahadevan dd2->numFields = dmb->numFields; 473b117cd09SVijay Mahadevan dd2->rw_dbglevel = dmb->rw_dbglevel; 474b117cd09SVijay Mahadevan dd2->partition_by_rank = dmb->partition_by_rank; 475b117cd09SVijay Mahadevan ierr = PetscStrcpy(dd2->extra_read_options, dmb->extra_read_options);CHKERRQ(ierr); 476b117cd09SVijay Mahadevan ierr = PetscStrcpy(dd2->extra_write_options, dmb->extra_write_options);CHKERRQ(ierr); 477b117cd09SVijay Mahadevan dd2->read_mode = dmb->read_mode; 478b117cd09SVijay Mahadevan dd2->write_mode = dmb->write_mode; 479b117cd09SVijay Mahadevan 480b117cd09SVijay Mahadevan /* set global ID tag handle */ 481b117cd09SVijay Mahadevan ierr = DMMoabSetLocalToGlobalTag(dm2, dmb->ltog_tag);CHKERRQ(ierr); 482b117cd09SVijay Mahadevan 483b117cd09SVijay Mahadevan merr = dd2->mbiface->tag_get_handle(MATERIAL_SET_TAG_NAME, dd2->material_tag);MBERRNM(merr); 484b117cd09SVijay Mahadevan 485b117cd09SVijay Mahadevan ierr = DMSetOptionsPrefix(dm2,((PetscObject)dm)->prefix);CHKERRQ(ierr); 486b117cd09SVijay Mahadevan ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 487b117cd09SVijay Mahadevan ierr = DMSetDimension(dm2,dim);CHKERRQ(ierr); 488b117cd09SVijay Mahadevan 489b117cd09SVijay Mahadevan /* allow overloaded (user replaced) operations to be inherited by refinement clones */ 490b117cd09SVijay Mahadevan dm2->ops->creatematrix = dm->ops->creatematrix; 491b117cd09SVijay Mahadevan 492b117cd09SVijay Mahadevan /* copy fill information if given */ 493b117cd09SVijay Mahadevan ierr = DMMoabSetBlockFills(dm2, dmb->dfill, dmb->ofill);CHKERRQ(ierr); 494b117cd09SVijay Mahadevan 495b117cd09SVijay Mahadevan /* copy vector type information */ 496b117cd09SVijay Mahadevan ierr = DMSetMatType(dm2,dm->mattype);CHKERRQ(ierr); 497b117cd09SVijay Mahadevan ierr = DMSetVecType(dm2,dm->vectype);CHKERRQ(ierr); 4983f1c6e43SVijay Mahadevan dd2->numFields = dmb->numFields; 4993f1c6e43SVijay Mahadevan if (dmb->numFields) { 500b117cd09SVijay Mahadevan ierr = DMMoabSetFieldNames(dm2,dmb->numFields,dmb->fieldNames);CHKERRQ(ierr); 5013f1c6e43SVijay Mahadevan } 502b117cd09SVijay Mahadevan 503b117cd09SVijay Mahadevan ierr = DMSetFromOptions(dm2);CHKERRQ(ierr); 504b117cd09SVijay Mahadevan 505b117cd09SVijay Mahadevan /* recreate Dof numbering for the refined DM and make sure the distribution is correctly populated */ 506b117cd09SVijay Mahadevan ierr = DMSetUp(dm2);CHKERRQ(ierr); 507b117cd09SVijay Mahadevan 508b117cd09SVijay Mahadevan *dmref = dm2; 509b117cd09SVijay Mahadevan PetscFunctionReturn(0); 510b117cd09SVijay Mahadevan } 511b117cd09SVijay Mahadevan 512b117cd09SVijay Mahadevan 513b117cd09SVijay Mahadevan #undef __FUNCT__ 514b117cd09SVijay Mahadevan #define __FUNCT__ "DMRefine_Moab" 515b117cd09SVijay Mahadevan /*@ 516b117cd09SVijay Mahadevan DMRefine_Moab - Generate a multi-level uniform refinement hierarchy 517b117cd09SVijay Mahadevan by succesively refining a coarse mesh, already defined in the DM object 518b117cd09SVijay Mahadevan provided by the user. 519b117cd09SVijay Mahadevan 520e882eb38SVijay Mahadevan Collective on DM 521b117cd09SVijay Mahadevan 522b117cd09SVijay Mahadevan Input Parameter: 523e882eb38SVijay Mahadevan + dm - The DMMoab object 524e882eb38SVijay Mahadevan - comm - the communicator to contain the new DM object (or MPI_COMM_NULL) 525b117cd09SVijay Mahadevan 526b117cd09SVijay Mahadevan Output Parameter: 527e882eb38SVijay Mahadevan . dmf - the refined DM, or NULL 528b117cd09SVijay Mahadevan 529e882eb38SVijay Mahadevan Note: If no refinement was done, the return value is NULL 530e882eb38SVijay Mahadevan 531e882eb38SVijay Mahadevan Level: developer 532b117cd09SVijay Mahadevan 533b117cd09SVijay Mahadevan .keywords: DMMoab, create, refinement 534b117cd09SVijay Mahadevan @*/ 535b117cd09SVijay Mahadevan PetscErrorCode DMRefine_Moab(DM dm,MPI_Comm comm,DM* dmf) 536b117cd09SVijay Mahadevan { 537b117cd09SVijay Mahadevan PetscErrorCode ierr; 538b117cd09SVijay Mahadevan 539b117cd09SVijay Mahadevan PetscFunctionBegin; 540b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 541b117cd09SVijay Mahadevan 542b117cd09SVijay Mahadevan ierr = DM_UMR_Moab_Private(dm,comm,PETSC_TRUE,dmf);CHKERRQ(ierr); 543b117cd09SVijay Mahadevan PetscFunctionReturn(0); 544b117cd09SVijay Mahadevan } 545b117cd09SVijay Mahadevan 546b117cd09SVijay Mahadevan #undef __FUNCT__ 547b117cd09SVijay Mahadevan #define __FUNCT__ "DMCoarsen_Moab" 548b117cd09SVijay Mahadevan /*@ 549b117cd09SVijay Mahadevan DMCoarsen_Moab - Generate a multi-level uniform refinement hierarchy 550b117cd09SVijay Mahadevan by succesively refining a coarse mesh, already defined in the DM object 551b117cd09SVijay Mahadevan provided by the user. 552b117cd09SVijay Mahadevan 553e882eb38SVijay Mahadevan Collective on DM 554b117cd09SVijay Mahadevan 555b117cd09SVijay Mahadevan Input Parameter: 556e882eb38SVijay Mahadevan . dm - The DMMoab object 557e882eb38SVijay Mahadevan - comm - the communicator to contain the new DM object (or MPI_COMM_NULL) 558b117cd09SVijay Mahadevan 559b117cd09SVijay Mahadevan Output Parameter: 560e882eb38SVijay Mahadevan . dmf - the coarsened DM, or NULL 561b117cd09SVijay Mahadevan 562e882eb38SVijay Mahadevan Note: If no coarsening was done, the return value is NULL 563b117cd09SVijay Mahadevan 564e882eb38SVijay Mahadevan Level: developer 565e882eb38SVijay Mahadevan 566e882eb38SVijay Mahadevan .keywords: DMMoab, create, coarsening 567b117cd09SVijay Mahadevan @*/ 568b117cd09SVijay Mahadevan PetscErrorCode DMCoarsen_Moab(DM dm,MPI_Comm comm,DM* dmc) 569b117cd09SVijay Mahadevan { 570b117cd09SVijay Mahadevan PetscErrorCode ierr; 571b117cd09SVijay Mahadevan 572b117cd09SVijay Mahadevan PetscFunctionBegin; 573b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 574b117cd09SVijay Mahadevan 575b117cd09SVijay Mahadevan ierr = DM_UMR_Moab_Private(dm,comm,PETSC_FALSE,dmc);CHKERRQ(ierr); 576b117cd09SVijay Mahadevan PetscFunctionReturn(0); 577b117cd09SVijay Mahadevan } 578