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; 160b117cd09SVijay Mahadevan std::vector<int> bndrows; 1613f1c6e43SVijay Mahadevan std::vector<PetscBool> dbdry; 162*941e0cffSVijay Mahadevan PetscInt innz, *nnz, ionz, *onz; 163*941e0cffSVijay Mahadevan PetscInt n_nnz, n_onz, nlsiz1, nlsiz2, ngsiz1, ngsiz2; 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; 170*941e0cffSVijay Mahadevan nlsiz1 = dmb1->nloc*dmb2->numFields; 171*941e0cffSVijay Mahadevan nlsiz2 = dmb2->nloc*dmb2->numFields; 172*941e0cffSVijay Mahadevan ngsiz1 = dmb1->n*dmb2->numFields; 173*941e0cffSVijay 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); 179*941e0cffSVijay 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 184*941e0cffSVijay Mahadevan /* allocate the nnz, onz arrays based on block size and local nodes */ 185*941e0cffSVijay Mahadevan ierr = PetscCalloc2(nlsiz2,&nnz,nlsiz2,&onz);CHKERRQ(ierr); 186*941e0cffSVijay Mahadevan 187*941e0cffSVijay Mahadevan /* Loop through the local elements and compute the relation between the current parent and the refined_level. */ 188*941e0cffSVijay Mahadevan for(moab::Range::iterator iter = dmb1->elocal->begin(); iter!= dmb1->elocal->end(); iter++) { 189*941e0cffSVijay Mahadevan 190*941e0cffSVijay Mahadevan const moab::EntityHandle ehandle = *iter; 191*941e0cffSVijay Mahadevan std::vector<moab::EntityHandle> children; 192*941e0cffSVijay Mahadevan std::vector<moab::EntityHandle> connp, connc; 193*941e0cffSVijay Mahadevan 194*941e0cffSVijay Mahadevan /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */ 195*941e0cffSVijay Mahadevan merr = dmb1->hierarchy->parent_to_child(ehandle, dmb1->hlevel, dmb2->hlevel, children);MBERRNM(merr); 196*941e0cffSVijay Mahadevan 197*941e0cffSVijay Mahadevan /* Get connectivity and coordinates of the parent vertices */ 198*941e0cffSVijay Mahadevan merr = dmb1->hierarchy->get_connectivity(ehandle, dmb1->hlevel, connp);MBERRNM(merr); 199*941e0cffSVijay Mahadevan for (unsigned ic=0; ic < children.size(); ic++) { 200*941e0cffSVijay Mahadevan std::vector<moab::EntityHandle> tconnc; 201*941e0cffSVijay Mahadevan /* Get coordinates of the parent vertices in canonical order */ 202*941e0cffSVijay Mahadevan merr = dmb2->hierarchy->get_connectivity(children[ic], dmb2->hlevel, tconnc);MBERRNM(merr); 203*941e0cffSVijay Mahadevan for (unsigned tc=0; tc<tconnc.size(); tc++) { 204*941e0cffSVijay Mahadevan if (std::find(connc.begin(), connc.end(), tconnc[tc]) == connc.end()) 205*941e0cffSVijay Mahadevan connc.push_back(tconnc[tc]); 206*941e0cffSVijay Mahadevan } 207*941e0cffSVijay Mahadevan } 208*941e0cffSVijay Mahadevan 209*941e0cffSVijay Mahadevan std::vector<int> dofsp(connp.size()), dofsc(connc.size()); 210*941e0cffSVijay Mahadevan /* TODO: specific to scalar system - use GetDofs */ 211*941e0cffSVijay Mahadevan //ierr = DMMoabGetFieldDofs(dm1, connp.size(), &connp[0], 0, &dofsp[0]);CHKERRQ(ierr); 212*941e0cffSVijay Mahadevan ierr = DMMoabGetFieldDofs(dm2, connc.size(), &connc[0], 0, &dofsc[0]);CHKERRQ(ierr); 213*941e0cffSVijay Mahadevan 214*941e0cffSVijay Mahadevan for (unsigned tp=0;tp<connp.size(); tp++) { 215*941e0cffSVijay Mahadevan // ierr = MatSetValues(*interpl, connc.size(), &dofsc[0], 1, &dofsp[i], &values_phi[0], ADD_VALUES);CHKERRQ(ierr); 216*941e0cffSVijay Mahadevan if (dmb1->vowned->find(connp[tp]) != dmb1->vowned->end()) { 217*941e0cffSVijay Mahadevan for (unsigned tc=0;tc<connc.size(); tc++) nnz[dofsc[tc]]++; 218*941e0cffSVijay Mahadevan } 219*941e0cffSVijay Mahadevan else if (dmb1->vghost->find(connp[tp]) != dmb1->vghost->end()) { 220*941e0cffSVijay Mahadevan for (unsigned tc=0;tc<connc.size(); tc++) onz[dofsc[tc]]++; 221*941e0cffSVijay Mahadevan } 222*941e0cffSVijay Mahadevan else SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Invalid entity in parent level %D\n", connc[tp]); 223*941e0cffSVijay Mahadevan } 224*941e0cffSVijay Mahadevan /* 225*941e0cffSVijay Mahadevan for(int tc = 0; tc < connc.size(); tc++) { 226*941e0cffSVijay Mahadevan if (dmb2->vowned->find(connc[tc]) != dmb2->vowned->end()) nnz[dofsc[tc]]++; 227*941e0cffSVijay Mahadevan else if (dmb2->vghost->find(connc[tc]) != dmb2->vghost->end()) onz[dofsc[tc]]++; 228*941e0cffSVijay Mahadevan else SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Invalid entity in child level %D\n", connc[tc]); 229*941e0cffSVijay Mahadevan }*/ 230*941e0cffSVijay Mahadevan } 231*941e0cffSVijay Mahadevan 232*941e0cffSVijay Mahadevan int i=0; 233*941e0cffSVijay Mahadevan std::vector<moab::EntityHandle> adjs; 234*941e0cffSVijay Mahadevan for(moab::Range::iterator iter = dmb1->vowned->begin(); iter!= dmb1->vowned->end(); iter++, i++) { 235*941e0cffSVijay Mahadevan merr = dmb1->hierarchy->get_adjacencies(*iter, 0, adjs);MBERRNM(merr); 236*941e0cffSVijay Mahadevan nnz[i] -= adjs.size(); 237*941e0cffSVijay Mahadevan adjs.clear(); 238*941e0cffSVijay Mahadevan } 239*941e0cffSVijay Mahadevan for(moab::Range::iterator iter = dmb1->vowned->begin(); iter!= dmb1->vowned->end(); iter++, i++) { 240*941e0cffSVijay Mahadevan //merr = dmb1->hierarchy->get_adjacencies(&(*iter), 1, 0, false, adjs, moab::Interface::UNION);MBERRNM(merr); 241*941e0cffSVijay Mahadevan merr = dmb1->hierarchy->get_adjacencies(*iter, 0, adjs);MBERRNM(merr); 242*941e0cffSVijay Mahadevan onz[i] -= adjs.size(); 243*941e0cffSVijay Mahadevan adjs.clear(); 244*941e0cffSVijay Mahadevan } 245*941e0cffSVijay Mahadevan 246*941e0cffSVijay Mahadevan ionz=onz[0]; 247*941e0cffSVijay Mahadevan innz=nnz[0]; 248*941e0cffSVijay Mahadevan for (int tc=0; tc < nlsiz2; tc++) { 249*941e0cffSVijay Mahadevan //nnz[tc]+=1; 250*941e0cffSVijay Mahadevan // check for maximum allowed sparsity = fully dense 251*941e0cffSVijay Mahadevan nnz[tc] = std::min(nlsiz1,nnz[tc]); 252*941e0cffSVijay Mahadevan onz[tc] = std::min(nlsiz1,onz[tc]); 253*941e0cffSVijay Mahadevan 254*941e0cffSVijay Mahadevan innz = (innz < nnz[tc] ? nnz[tc] : innz); 255*941e0cffSVijay Mahadevan ionz = (ionz < onz[tc] ? onz[tc] : ionz); 256*941e0cffSVijay Mahadevan 257*941e0cffSVijay Mahadevan //PetscPrintf(PETSC_COMM_WORLD, "[%D] NNZ = %D, ONZ = %D\n", tc, nnz[tc], onz[tc]); 258*941e0cffSVijay Mahadevan } 259*941e0cffSVijay Mahadevan //PetscPrintf(PETSC_COMM_WORLD, "Final: INNZ = %D, IONZ = %D\n", innz, ionz); 260*941e0cffSVijay Mahadevan 261*941e0cffSVijay Mahadevan ierr = MatSeqAIJSetPreallocation(*interpl,innz,nnz);CHKERRQ(ierr); 262*941e0cffSVijay Mahadevan ierr = MatMPIAIJSetPreallocation(*interpl,innz,nnz,ionz,onz);CHKERRQ(ierr); 263*941e0cffSVijay Mahadevan 264*941e0cffSVijay Mahadevan /* clean up temporary memory */ 265*941e0cffSVijay Mahadevan ierr = PetscFree2(nnz,onz);CHKERRQ(ierr); 266b117cd09SVijay Mahadevan 267b117cd09SVijay Mahadevan /* set up internal matrix data-structures */ 268b117cd09SVijay Mahadevan ierr = MatSetUp(*interpl);CHKERRQ(ierr); 269b117cd09SVijay Mahadevan ierr = MatZeroEntries(*interpl);CHKERRQ(ierr); 270b117cd09SVijay Mahadevan 271b117cd09SVijay Mahadevan ierr = DMGetDimension(dm1, &dim);CHKERRQ(ierr); 272b117cd09SVijay Mahadevan 2733f1c6e43SVijay Mahadevan factor = std::pow(2.0 /*degree_P_for_refinement*/,(dmb2->hlevel-dmb1->hlevel)*dmb1->dim*1.0); 274b117cd09SVijay Mahadevan 275e882eb38SVijay Mahadevan /* Loop through the remaining vertices. These vertices appear only on the current refined_level. */ 276b117cd09SVijay Mahadevan for(moab::Range::iterator iter = dmb1->elocal->begin(); iter!= dmb1->elocal->end(); iter++) { 277b117cd09SVijay Mahadevan 278b117cd09SVijay Mahadevan const moab::EntityHandle ehandle = *iter; 279b117cd09SVijay Mahadevan std::vector<moab::EntityHandle> children; 280b117cd09SVijay Mahadevan std::vector<moab::EntityHandle> connp, connc; 281b117cd09SVijay Mahadevan 282b117cd09SVijay Mahadevan /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */ 283b117cd09SVijay Mahadevan merr = dmb1->hierarchy->parent_to_child(ehandle, dmb1->hlevel, dmb2->hlevel, children);MBERRNM(merr); 284b117cd09SVijay Mahadevan 285b117cd09SVijay Mahadevan /* Get connectivity and coordinates of the parent vertices */ 286b117cd09SVijay Mahadevan merr = dmb1->hierarchy->get_connectivity(ehandle, dmb1->hlevel, connp);MBERRNM(merr); 287e882eb38SVijay Mahadevan for (unsigned ic=0; ic < children.size(); ic++) { 288b117cd09SVijay Mahadevan std::vector<moab::EntityHandle> tconnc; 289b117cd09SVijay Mahadevan /* Get coordinates of the parent vertices in canonical order */ 290b117cd09SVijay Mahadevan merr = dmb1->hierarchy->get_connectivity(children[ic], dmb2->hlevel, tconnc);MBERRNM(merr); 291e882eb38SVijay Mahadevan for (unsigned tc=0; tc<tconnc.size(); tc++) { 292b117cd09SVijay Mahadevan connc.push_back(tconnc[tc]); 293b117cd09SVijay Mahadevan } 294b117cd09SVijay Mahadevan } 295b117cd09SVijay Mahadevan 296b117cd09SVijay Mahadevan std::vector<double> pcoords(connp.size()*3), ccoords(connc.size()*3), values_phi(connc.size()); 297b117cd09SVijay Mahadevan /* Get coordinates for connectivity entities in canonical order for both coarse and finer levels */ 298b117cd09SVijay Mahadevan merr = dmb1->hierarchy->get_coordinates(&connp[0], connp.size(), dmb1->hlevel, &pcoords[0]);MBERRNM(merr); 299b117cd09SVijay Mahadevan merr = dmb2->hierarchy->get_coordinates(&connc[0], connc.size(), dmb2->hlevel, &ccoords[0]);MBERRNM(merr); 300b117cd09SVijay Mahadevan 301b117cd09SVijay Mahadevan std::vector<int> dofsp(connp.size()), dofsc(connc.size()); 302b117cd09SVijay Mahadevan /* TODO: specific to scalar system - use GetDofs */ 303b117cd09SVijay Mahadevan ierr = DMMoabGetFieldDofs(dm1, connp.size(), &connp[0], 0, &dofsp[0]);CHKERRQ(ierr); 304b117cd09SVijay Mahadevan ierr = DMMoabGetFieldDofs(dm2, connc.size(), &connc[0], 0, &dofsc[0]);CHKERRQ(ierr); 305b117cd09SVijay Mahadevan 306e882eb38SVijay Mahadevan /* Compute the interpolation weights by determining distance of 1-ring 307e882eb38SVijay Mahadevan neighbor vertices from current vertex */ 308e882eb38SVijay Mahadevan for (unsigned i=0;i<connp.size(); i++) { 309b117cd09SVijay Mahadevan double normsum=0.0; 310e882eb38SVijay Mahadevan for (unsigned j=0;j<connc.size(); j++) { 3113f1c6e43SVijay Mahadevan values_phi[j] = 0.0; 312e882eb38SVijay Mahadevan for (unsigned k=0;k<3; k++) 313*941e0cffSVijay Mahadevan values_phi[j] += std::pow(pcoords[i*3+k]-ccoords[k+j*3], dim); 3143f1c6e43SVijay Mahadevan if (values_phi[j] < 1e-12) { 3153f1c6e43SVijay Mahadevan values_phi[j] = 1e12; 316b117cd09SVijay Mahadevan } 317b117cd09SVijay Mahadevan else { 318*941e0cffSVijay Mahadevan //values_phi[j] = std::pow(values_phi[j], -1.0/dim); 319*941e0cffSVijay Mahadevan values_phi[j] = std::pow(values_phi[j], -1.0); 3203f1c6e43SVijay Mahadevan normsum += values_phi[j]; 321b117cd09SVijay Mahadevan } 322b117cd09SVijay Mahadevan } 323e882eb38SVijay Mahadevan for (unsigned j=0;j<connc.size(); j++) { 3243f1c6e43SVijay Mahadevan if (values_phi[j] > 1e11) 3253f1c6e43SVijay Mahadevan values_phi[j] = factor*0.5/connc.size(); 326b117cd09SVijay Mahadevan else 3273f1c6e43SVijay Mahadevan values_phi[j] = factor*values_phi[j]*0.5/(connc.size()*normsum); 328b117cd09SVijay Mahadevan } 329b117cd09SVijay Mahadevan ierr = MatSetValues(*interpl, connc.size(), &dofsc[0], 1, &dofsp[i], &values_phi[0], ADD_VALUES);CHKERRQ(ierr); 330b117cd09SVijay Mahadevan } 331b117cd09SVijay Mahadevan 332b117cd09SVijay Mahadevan /* check if element is on the boundary */ 333b117cd09SVijay Mahadevan //ierr = DMMoabIsEntityOnBoundary(dm1,ehandle,&eonbnd);CHKERRQ(ierr); 3343f1c6e43SVijay Mahadevan dbdry.resize(connc.size()); 3353f1c6e43SVijay Mahadevan ierr = DMMoabCheckBoundaryVertices(dm2,connc.size(),&connc[0],dbdry.data());CHKERRQ(ierr); 336b117cd09SVijay Mahadevan eonbnd=PETSC_FALSE; 337e882eb38SVijay Mahadevan for (unsigned i=0; i< connc.size(); ++i) 338b117cd09SVijay Mahadevan if (dbdry[i]) eonbnd=PETSC_TRUE; 339b117cd09SVijay Mahadevan 340b117cd09SVijay Mahadevan values_phi.clear(); 341b117cd09SVijay Mahadevan values_phi.resize(connp.size()); 342b117cd09SVijay Mahadevan /* apply dirichlet boundary conditions */ 343b117cd09SVijay Mahadevan if (eonbnd) { 344b117cd09SVijay Mahadevan 345b117cd09SVijay Mahadevan ierr = MatAssemblyBegin(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr); 346b117cd09SVijay Mahadevan ierr = MatAssemblyEnd(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr); 347b117cd09SVijay Mahadevan /* get the list of nodes on boundary so that we can enforce dirichlet conditions strongly */ 348b117cd09SVijay Mahadevan //ierr = DMMoabCheckBoundaryVertices(dm2,connc.size(),&connc[0],dbdry);CHKERRQ(ierr); 349e882eb38SVijay Mahadevan for (unsigned i=0; i < connc.size(); i++) { 3503f1c6e43SVijay Mahadevan if (dbdry[i]) { /* dirichlet node */ 351b117cd09SVijay Mahadevan /* think about strongly imposing dirichlet */ 352b117cd09SVijay Mahadevan //bndrows.push_back(dofsc[i]); 353b117cd09SVijay Mahadevan 354b117cd09SVijay Mahadevan ierr = MatSetValues(*interpl, 1, &dofsc[i], connp.size(), &dofsp[0], &values_phi[0], INSERT_VALUES);CHKERRQ(ierr); 355b117cd09SVijay Mahadevan //values_phi[0]=1.0; 356b117cd09SVijay Mahadevan //ierr = MatSetValues(*interpl, 1, &dofsc[i], 1, &dofsc[i], &values_phi[0], INSERT_VALUES);CHKERRQ(ierr); 357b117cd09SVijay Mahadevan } 358b117cd09SVijay Mahadevan } 359b117cd09SVijay Mahadevan 360b117cd09SVijay Mahadevan ierr = MatAssemblyBegin(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr); 361b117cd09SVijay Mahadevan ierr = MatAssemblyEnd(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr); 362b117cd09SVijay Mahadevan } 363b117cd09SVijay Mahadevan 364b117cd09SVijay Mahadevan //get interpolation weights 365b117cd09SVijay Mahadevan //ierr = Compute_Quad4_Basis(pcoords, 1, coord, values_phi);CHKERRQ(ierr); 366e882eb38SVijay Mahadevan // for (int j=0;j<dofs_per_element; j++) 367e882eb38SVijay Mahadevan // std::cout<<"values "<<values_phi[j]<<std::endl; 368b117cd09SVijay Mahadevan 369b117cd09SVijay Mahadevan //get row and column indices, zero weights are ignored 370b117cd09SVijay Mahadevan /* 371b117cd09SVijay Mahadevan int nz_ind = 0; 372b117cd09SVijay Mahadevan idx = dmb2->vowned->index(vhandle); 373b117cd09SVijay Mahadevan for (int j=0;j<dofs_per_element; j++){ 374b117cd09SVijay Mahadevan idy[nz_ind] = dmb1->vowned->index(connectivity[j]); 375b117cd09SVijay 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]); 376b117cd09SVijay Mahadevan //values_phi[nz_ind] = values_phi[j]; 377b117cd09SVijay Mahadevan nz_ind = nz_ind+1; 378b117cd09SVijay Mahadevan } 379b117cd09SVijay Mahadevan */ 380b117cd09SVijay Mahadevan 381b117cd09SVijay Mahadevan //ierr = MatSetValues(*interpl, nz_ind, idy, 1, &idx, values_phi, INSERT_VALUES);CHKERRQ(ierr); 382b117cd09SVijay Mahadevan //ierr = MatSetValues(*interpl, connp.size(), dofsp, connc.size(), dofsc, &values_phi[0], INSERT_VALUES);CHKERRQ(ierr); 383b117cd09SVijay Mahadevan } 384b117cd09SVijay Mahadevan 385b117cd09SVijay 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]); 386b117cd09SVijay Mahadevan //ierr = MatZeroRows(*interpl, bndrows.size(), &bndrows[0], 1.0, 0, 0);CHKERRQ(ierr); 387b117cd09SVijay Mahadevan 388b117cd09SVijay Mahadevan ierr = MatAssemblyBegin(*interpl,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 389b117cd09SVijay Mahadevan ierr = MatAssemblyEnd(*interpl,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 390b117cd09SVijay Mahadevan PetscFunctionReturn(0); 391b117cd09SVijay Mahadevan } 392b117cd09SVijay Mahadevan 393b117cd09SVijay Mahadevan #undef __FUNCT__ 394b117cd09SVijay Mahadevan #define __FUNCT__ "DMCreateInjection_Moab" 395b117cd09SVijay Mahadevan /*@ 396b117cd09SVijay Mahadevan DMCreateInjection_Moab - Generate a multi-level uniform refinement hierarchy 397b117cd09SVijay Mahadevan by succesively refining a coarse mesh, already defined in the DM object 398b117cd09SVijay Mahadevan provided by the user. 399b117cd09SVijay Mahadevan 400b117cd09SVijay Mahadevan Collective on MPI_Comm 401b117cd09SVijay Mahadevan 402b117cd09SVijay Mahadevan Input Parameter: 403b117cd09SVijay Mahadevan . dmb - The DMMoab object 404b117cd09SVijay Mahadevan 405b117cd09SVijay Mahadevan Output Parameter: 406b117cd09SVijay Mahadevan . nlevels - The number of levels of refinement needed to generate the hierarchy 407b117cd09SVijay Mahadevan + ldegrees - The degree of refinement at each level in the hierarchy 408b117cd09SVijay Mahadevan 409b117cd09SVijay Mahadevan Level: beginner 410b117cd09SVijay Mahadevan 411b117cd09SVijay Mahadevan .keywords: DMMoab, create, refinement 412b117cd09SVijay Mahadevan @*/ 413b117cd09SVijay Mahadevan PetscErrorCode DMCreateInjection_Moab(DM dm1,DM dm2,VecScatter* ctx) 414b117cd09SVijay Mahadevan { 415e882eb38SVijay Mahadevan //DM_Moab *dmmoab; 416b117cd09SVijay Mahadevan 417b117cd09SVijay Mahadevan PetscFunctionBegin; 418b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm1,DM_CLASSID,1); 419b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm2,DM_CLASSID,2); 420e882eb38SVijay Mahadevan //dmmoab = (DM_Moab*)(dm1)->data; 421b117cd09SVijay Mahadevan 422b117cd09SVijay Mahadevan PetscPrintf(PETSC_COMM_WORLD, "[DMCreateInjection_Moab] :: Placeholder\n"); 423b117cd09SVijay Mahadevan PetscFunctionReturn(0); 424b117cd09SVijay Mahadevan } 425b117cd09SVijay Mahadevan 426b117cd09SVijay Mahadevan 427b117cd09SVijay Mahadevan #undef __FUNCT__ 428b117cd09SVijay Mahadevan #define __FUNCT__ "DM_UMR_Moab_Private" 429b117cd09SVijay Mahadevan PetscErrorCode DM_UMR_Moab_Private(DM dm,MPI_Comm comm,PetscBool refine,DM *dmref) 430b117cd09SVijay Mahadevan { 431b117cd09SVijay Mahadevan PetscErrorCode ierr; 432e882eb38SVijay Mahadevan PetscInt i,dim; 433b117cd09SVijay Mahadevan DM dm2; 434b117cd09SVijay Mahadevan moab::ErrorCode merr; 435b117cd09SVijay Mahadevan DM_Moab *dmb = (DM_Moab*)dm->data,*dd2; 436b117cd09SVijay Mahadevan 437b117cd09SVijay Mahadevan PetscFunctionBegin; 438b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 439b117cd09SVijay Mahadevan PetscValidPointer(dmref,4); 440b117cd09SVijay Mahadevan 441b117cd09SVijay Mahadevan if ( (dmb->hlevel == dmb->nhlevels && refine) || (dmb->hlevel == 0 && !refine) ) { 442e882eb38SVijay 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); 443e882eb38SVijay Mahadevan if (dmb->hlevel-1 < 0 && !refine) PetscInfo1(NULL,"Invalid multigrid coarsen hierarchy level specified (%D). Creating a NULL object.\n",dmb->hlevel-1); 444e882eb38SVijay Mahadevan *dmref = PETSC_NULL; 445b117cd09SVijay Mahadevan PetscFunctionReturn(0); 446b117cd09SVijay Mahadevan } 447b117cd09SVijay Mahadevan 448b117cd09SVijay Mahadevan ierr = DMMoabCreate(PetscObjectComm((PetscObject)dm), &dm2);CHKERRQ(ierr); 449b117cd09SVijay Mahadevan dd2 = (DM_Moab*)dm2->data; 450b117cd09SVijay Mahadevan 451b117cd09SVijay Mahadevan dd2->mbiface = dmb->mbiface; 452b117cd09SVijay Mahadevan dd2->pcomm = dmb->pcomm; 453b117cd09SVijay Mahadevan dd2->icreatedinstance = PETSC_FALSE; 454b117cd09SVijay Mahadevan 455b117cd09SVijay Mahadevan /* set the new level based on refinement/coarsening */ 456b117cd09SVijay Mahadevan if (refine) { 457b117cd09SVijay Mahadevan dd2->hlevel=dmb->hlevel+1; 458b117cd09SVijay Mahadevan } 459b117cd09SVijay Mahadevan else { 460b117cd09SVijay Mahadevan dd2->hlevel=dmb->hlevel-1; 461b117cd09SVijay Mahadevan } 462b117cd09SVijay Mahadevan 463b117cd09SVijay Mahadevan /* Copy the multilevel hierarchy pointers in MOAB */ 464b117cd09SVijay Mahadevan dd2->hierarchy = dmb->hierarchy; 465b117cd09SVijay Mahadevan dd2->nhlevels = dmb->nhlevels; 466b117cd09SVijay Mahadevan ierr = PetscMalloc1(dd2->nhlevels+1,&dd2->hsets);CHKERRQ(ierr); 467b117cd09SVijay Mahadevan for (i=0; i<=dd2->nhlevels; i++) { 468b117cd09SVijay Mahadevan dd2->hsets[i]=dmb->hsets[i]; 469b117cd09SVijay Mahadevan } 470b117cd09SVijay Mahadevan dd2->fileset = dd2->hsets[dd2->hlevel]; 471b117cd09SVijay Mahadevan 472b117cd09SVijay Mahadevan /* do the remaining initializations for DMMoab */ 473b117cd09SVijay Mahadevan dd2->bs = dmb->bs; 474b117cd09SVijay Mahadevan dd2->numFields = dmb->numFields; 475b117cd09SVijay Mahadevan dd2->rw_dbglevel = dmb->rw_dbglevel; 476b117cd09SVijay Mahadevan dd2->partition_by_rank = dmb->partition_by_rank; 477b117cd09SVijay Mahadevan ierr = PetscStrcpy(dd2->extra_read_options, dmb->extra_read_options);CHKERRQ(ierr); 478b117cd09SVijay Mahadevan ierr = PetscStrcpy(dd2->extra_write_options, dmb->extra_write_options);CHKERRQ(ierr); 479b117cd09SVijay Mahadevan dd2->read_mode = dmb->read_mode; 480b117cd09SVijay Mahadevan dd2->write_mode = dmb->write_mode; 481b117cd09SVijay Mahadevan 482b117cd09SVijay Mahadevan /* set global ID tag handle */ 483b117cd09SVijay Mahadevan ierr = DMMoabSetLocalToGlobalTag(dm2, dmb->ltog_tag);CHKERRQ(ierr); 484b117cd09SVijay Mahadevan 485b117cd09SVijay Mahadevan merr = dd2->mbiface->tag_get_handle(MATERIAL_SET_TAG_NAME, dd2->material_tag);MBERRNM(merr); 486b117cd09SVijay Mahadevan 487b117cd09SVijay Mahadevan ierr = DMSetOptionsPrefix(dm2,((PetscObject)dm)->prefix);CHKERRQ(ierr); 488b117cd09SVijay Mahadevan ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 489b117cd09SVijay Mahadevan ierr = DMSetDimension(dm2,dim);CHKERRQ(ierr); 490b117cd09SVijay Mahadevan 491b117cd09SVijay Mahadevan /* allow overloaded (user replaced) operations to be inherited by refinement clones */ 492b117cd09SVijay Mahadevan dm2->ops->creatematrix = dm->ops->creatematrix; 493b117cd09SVijay Mahadevan 494b117cd09SVijay Mahadevan /* copy fill information if given */ 495b117cd09SVijay Mahadevan ierr = DMMoabSetBlockFills(dm2, dmb->dfill, dmb->ofill);CHKERRQ(ierr); 496b117cd09SVijay Mahadevan 497b117cd09SVijay Mahadevan /* copy vector type information */ 498b117cd09SVijay Mahadevan ierr = DMSetMatType(dm2,dm->mattype);CHKERRQ(ierr); 499b117cd09SVijay Mahadevan ierr = DMSetVecType(dm2,dm->vectype);CHKERRQ(ierr); 5003f1c6e43SVijay Mahadevan dd2->numFields = dmb->numFields; 5013f1c6e43SVijay Mahadevan if (dmb->numFields) { 502b117cd09SVijay Mahadevan ierr = DMMoabSetFieldNames(dm2,dmb->numFields,dmb->fieldNames);CHKERRQ(ierr); 5033f1c6e43SVijay Mahadevan } 504b117cd09SVijay Mahadevan 505b117cd09SVijay Mahadevan ierr = DMSetFromOptions(dm2);CHKERRQ(ierr); 506b117cd09SVijay Mahadevan 507b117cd09SVijay Mahadevan /* recreate Dof numbering for the refined DM and make sure the distribution is correctly populated */ 508b117cd09SVijay Mahadevan ierr = DMSetUp(dm2);CHKERRQ(ierr); 509b117cd09SVijay Mahadevan 510b117cd09SVijay Mahadevan *dmref = dm2; 511b117cd09SVijay Mahadevan PetscFunctionReturn(0); 512b117cd09SVijay Mahadevan } 513b117cd09SVijay Mahadevan 514b117cd09SVijay Mahadevan 515b117cd09SVijay Mahadevan #undef __FUNCT__ 516b117cd09SVijay Mahadevan #define __FUNCT__ "DMRefine_Moab" 517b117cd09SVijay Mahadevan /*@ 518b117cd09SVijay Mahadevan DMRefine_Moab - Generate a multi-level uniform refinement hierarchy 519b117cd09SVijay Mahadevan by succesively refining a coarse mesh, already defined in the DM object 520b117cd09SVijay Mahadevan provided by the user. 521b117cd09SVijay Mahadevan 522e882eb38SVijay Mahadevan Collective on DM 523b117cd09SVijay Mahadevan 524b117cd09SVijay Mahadevan Input Parameter: 525e882eb38SVijay Mahadevan + dm - The DMMoab object 526e882eb38SVijay Mahadevan - comm - the communicator to contain the new DM object (or MPI_COMM_NULL) 527b117cd09SVijay Mahadevan 528b117cd09SVijay Mahadevan Output Parameter: 529e882eb38SVijay Mahadevan . dmf - the refined DM, or NULL 530b117cd09SVijay Mahadevan 531e882eb38SVijay Mahadevan Note: If no refinement was done, the return value is NULL 532e882eb38SVijay Mahadevan 533e882eb38SVijay Mahadevan Level: developer 534b117cd09SVijay Mahadevan 535b117cd09SVijay Mahadevan .keywords: DMMoab, create, refinement 536b117cd09SVijay Mahadevan @*/ 537b117cd09SVijay Mahadevan PetscErrorCode DMRefine_Moab(DM dm,MPI_Comm comm,DM* dmf) 538b117cd09SVijay Mahadevan { 539b117cd09SVijay Mahadevan PetscErrorCode ierr; 540b117cd09SVijay Mahadevan 541b117cd09SVijay Mahadevan PetscFunctionBegin; 542b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 543b117cd09SVijay Mahadevan 544b117cd09SVijay Mahadevan ierr = DM_UMR_Moab_Private(dm,comm,PETSC_TRUE,dmf);CHKERRQ(ierr); 545b117cd09SVijay Mahadevan PetscFunctionReturn(0); 546b117cd09SVijay Mahadevan } 547b117cd09SVijay Mahadevan 548b117cd09SVijay Mahadevan #undef __FUNCT__ 549b117cd09SVijay Mahadevan #define __FUNCT__ "DMCoarsen_Moab" 550b117cd09SVijay Mahadevan /*@ 551b117cd09SVijay Mahadevan DMCoarsen_Moab - Generate a multi-level uniform refinement hierarchy 552b117cd09SVijay Mahadevan by succesively refining a coarse mesh, already defined in the DM object 553b117cd09SVijay Mahadevan provided by the user. 554b117cd09SVijay Mahadevan 555e882eb38SVijay Mahadevan Collective on DM 556b117cd09SVijay Mahadevan 557b117cd09SVijay Mahadevan Input Parameter: 558e882eb38SVijay Mahadevan . dm - The DMMoab object 559e882eb38SVijay Mahadevan - comm - the communicator to contain the new DM object (or MPI_COMM_NULL) 560b117cd09SVijay Mahadevan 561b117cd09SVijay Mahadevan Output Parameter: 562e882eb38SVijay Mahadevan . dmf - the coarsened DM, or NULL 563b117cd09SVijay Mahadevan 564e882eb38SVijay Mahadevan Note: If no coarsening was done, the return value is NULL 565b117cd09SVijay Mahadevan 566e882eb38SVijay Mahadevan Level: developer 567e882eb38SVijay Mahadevan 568e882eb38SVijay Mahadevan .keywords: DMMoab, create, coarsening 569b117cd09SVijay Mahadevan @*/ 570b117cd09SVijay Mahadevan PetscErrorCode DMCoarsen_Moab(DM dm,MPI_Comm comm,DM* dmc) 571b117cd09SVijay Mahadevan { 572b117cd09SVijay Mahadevan PetscErrorCode ierr; 573b117cd09SVijay Mahadevan 574b117cd09SVijay Mahadevan PetscFunctionBegin; 575b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 576b117cd09SVijay Mahadevan 577b117cd09SVijay Mahadevan ierr = DM_UMR_Moab_Private(dm,comm,PETSC_FALSE,dmc);CHKERRQ(ierr); 578b117cd09SVijay Mahadevan PetscFunctionReturn(0); 579b117cd09SVijay Mahadevan } 580