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 */ 49*3f1c6e43SVijay 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; 158*3f1c6e43SVijay Mahadevan PetscReal factor; 159*3f1c6e43SVijay Mahadevan PetscBool eonbnd; 160b117cd09SVijay Mahadevan std::vector<int> bndrows; 161*3f1c6e43SVijay Mahadevan std::vector<PetscBool> dbdry; 162b117cd09SVijay Mahadevan 163b117cd09SVijay Mahadevan PetscFunctionBegin; 164b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm1,DM_CLASSID,1); 165b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm2,DM_CLASSID,2); 166b117cd09SVijay Mahadevan dmb1 = (DM_Moab*)(dm1)->data; 167b117cd09SVijay Mahadevan dmb2 = (DM_Moab*)(dm2)->data; 168b117cd09SVijay Mahadevan 169*3f1c6e43SVijay 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); 170*3f1c6e43SVijay Mahadevan 171b117cd09SVijay Mahadevan ierr = MatCreate(PetscObjectComm((PetscObject)dm1), interpl);CHKERRQ(ierr); 172e882eb38SVijay Mahadevan ierr = MatSetType(*interpl, dm1->mattype);CHKERRQ(ierr); 173b117cd09SVijay Mahadevan ierr = MatSetSizes(*interpl, dmb2->nloc, dmb1->nloc, dmb2->n, dmb1->n);CHKERRQ(ierr); 174b117cd09SVijay Mahadevan 175b117cd09SVijay Mahadevan /* TODO: This is a hack for the rectangular system - decipher NNZ pattern better */ 176b117cd09SVijay Mahadevan ierr = MatSetFromOptions(*interpl);CHKERRQ(ierr); 177b117cd09SVijay Mahadevan 178b117cd09SVijay Mahadevan ierr = MatSeqAIJSetPreallocation(*interpl,dmb2->nloc,0);CHKERRQ(ierr); 179b117cd09SVijay Mahadevan ierr = MatMPIAIJSetPreallocation(*interpl,dmb2->nloc,0,dmb2->nghost,0);CHKERRQ(ierr); 180b117cd09SVijay Mahadevan 181b117cd09SVijay Mahadevan /* set up internal matrix data-structures */ 182b117cd09SVijay Mahadevan ierr = MatSetUp(*interpl);CHKERRQ(ierr); 183b117cd09SVijay Mahadevan ierr = MatZeroEntries(*interpl);CHKERRQ(ierr); 184b117cd09SVijay Mahadevan 185b117cd09SVijay Mahadevan ierr = DMGetDimension(dm1, &dim);CHKERRQ(ierr); 186b117cd09SVijay Mahadevan 187*3f1c6e43SVijay Mahadevan factor = std::pow(2.0 /*degree_P_for_refinement*/,(dmb2->hlevel-dmb1->hlevel)*dmb1->dim*1.0); 188b117cd09SVijay Mahadevan 189e882eb38SVijay Mahadevan /* Loop through the remaining vertices. These vertices appear only on the current refined_level. */ 190b117cd09SVijay Mahadevan for(moab::Range::iterator iter = dmb1->elocal->begin(); iter!= dmb1->elocal->end(); iter++) { 191b117cd09SVijay Mahadevan 192b117cd09SVijay Mahadevan const moab::EntityHandle ehandle = *iter; 193b117cd09SVijay Mahadevan std::vector<moab::EntityHandle> children; 194b117cd09SVijay Mahadevan std::vector<moab::EntityHandle> connp, connc; 195b117cd09SVijay Mahadevan 196b117cd09SVijay Mahadevan /* Get the relation between the current (coarse) parent and its corresponding (finer) children elements */ 197b117cd09SVijay Mahadevan merr = dmb1->hierarchy->parent_to_child(ehandle, dmb1->hlevel, dmb2->hlevel, children);MBERRNM(merr); 198b117cd09SVijay Mahadevan 199b117cd09SVijay Mahadevan /* Get connectivity and coordinates of the parent vertices */ 200b117cd09SVijay Mahadevan merr = dmb1->hierarchy->get_connectivity(ehandle, dmb1->hlevel, connp);MBERRNM(merr); 201e882eb38SVijay Mahadevan for (unsigned ic=0; ic < children.size(); ic++) { 202b117cd09SVijay Mahadevan std::vector<moab::EntityHandle> tconnc; 203b117cd09SVijay Mahadevan /* Get coordinates of the parent vertices in canonical order */ 204b117cd09SVijay Mahadevan merr = dmb1->hierarchy->get_connectivity(children[ic], dmb2->hlevel, tconnc);MBERRNM(merr); 205e882eb38SVijay Mahadevan for (unsigned tc=0; tc<tconnc.size(); tc++) { 206b117cd09SVijay Mahadevan connc.push_back(tconnc[tc]); 207b117cd09SVijay Mahadevan } 208b117cd09SVijay Mahadevan } 209b117cd09SVijay Mahadevan 210b117cd09SVijay Mahadevan std::vector<double> pcoords(connp.size()*3), ccoords(connc.size()*3), values_phi(connc.size()); 211b117cd09SVijay Mahadevan /* Get coordinates for connectivity entities in canonical order for both coarse and finer levels */ 212b117cd09SVijay Mahadevan merr = dmb1->hierarchy->get_coordinates(&connp[0], connp.size(), dmb1->hlevel, &pcoords[0]);MBERRNM(merr); 213b117cd09SVijay Mahadevan merr = dmb2->hierarchy->get_coordinates(&connc[0], connc.size(), dmb2->hlevel, &ccoords[0]);MBERRNM(merr); 214b117cd09SVijay Mahadevan 215b117cd09SVijay Mahadevan std::vector<int> dofsp(connp.size()), dofsc(connc.size()); 216b117cd09SVijay Mahadevan /* TODO: specific to scalar system - use GetDofs */ 217b117cd09SVijay Mahadevan ierr = DMMoabGetFieldDofs(dm1, connp.size(), &connp[0], 0, &dofsp[0]);CHKERRQ(ierr); 218b117cd09SVijay Mahadevan ierr = DMMoabGetFieldDofs(dm2, connc.size(), &connc[0], 0, &dofsc[0]);CHKERRQ(ierr); 219b117cd09SVijay Mahadevan 220e882eb38SVijay Mahadevan /* Compute the interpolation weights by determining distance of 1-ring 221e882eb38SVijay Mahadevan neighbor vertices from current vertex */ 222e882eb38SVijay Mahadevan for (unsigned i=0;i<connp.size(); i++) { 223b117cd09SVijay Mahadevan double normsum=0.0; 224e882eb38SVijay Mahadevan for (unsigned j=0;j<connc.size(); j++) { 225*3f1c6e43SVijay Mahadevan values_phi[j] = 0.0; 226e882eb38SVijay Mahadevan for (unsigned k=0;k<3; k++) 227*3f1c6e43SVijay Mahadevan values_phi[j] += (pcoords[i*3+k]-ccoords[k+j*3])*(pcoords[i*3+k]-ccoords[k+j*3]); 228*3f1c6e43SVijay Mahadevan if (values_phi[j] < 1e-12) { 229*3f1c6e43SVijay Mahadevan values_phi[j] = 1e12; 230b117cd09SVijay Mahadevan } 231b117cd09SVijay Mahadevan else { 232*3f1c6e43SVijay Mahadevan values_phi[j] = 1.0/(values_phi[j]); 233*3f1c6e43SVijay Mahadevan normsum += values_phi[j]; 234b117cd09SVijay Mahadevan } 235b117cd09SVijay Mahadevan } 236e882eb38SVijay Mahadevan for (unsigned j=0;j<connc.size(); j++) { 237*3f1c6e43SVijay Mahadevan if (values_phi[j] > 1e11) 238*3f1c6e43SVijay Mahadevan values_phi[j] = factor*0.5/connc.size(); 239b117cd09SVijay Mahadevan else 240*3f1c6e43SVijay Mahadevan values_phi[j] = factor*values_phi[j]*0.5/(connc.size()*normsum); 241b117cd09SVijay Mahadevan } 242b117cd09SVijay Mahadevan ierr = MatSetValues(*interpl, connc.size(), &dofsc[0], 1, &dofsp[i], &values_phi[0], ADD_VALUES);CHKERRQ(ierr); 243b117cd09SVijay Mahadevan } 244b117cd09SVijay Mahadevan 245b117cd09SVijay Mahadevan /* check if element is on the boundary */ 246b117cd09SVijay Mahadevan //ierr = DMMoabIsEntityOnBoundary(dm1,ehandle,&eonbnd);CHKERRQ(ierr); 247*3f1c6e43SVijay Mahadevan dbdry.resize(connc.size()); 248*3f1c6e43SVijay Mahadevan ierr = DMMoabCheckBoundaryVertices(dm2,connc.size(),&connc[0],dbdry.data());CHKERRQ(ierr); 249b117cd09SVijay Mahadevan eonbnd=PETSC_FALSE; 250e882eb38SVijay Mahadevan for (unsigned i=0; i< connc.size(); ++i) 251b117cd09SVijay Mahadevan if (dbdry[i]) eonbnd=PETSC_TRUE; 252b117cd09SVijay Mahadevan 253b117cd09SVijay Mahadevan values_phi.clear(); 254b117cd09SVijay Mahadevan values_phi.resize(connp.size()); 255b117cd09SVijay Mahadevan /* apply dirichlet boundary conditions */ 256b117cd09SVijay Mahadevan if (eonbnd) { 257b117cd09SVijay Mahadevan 258b117cd09SVijay Mahadevan ierr = MatAssemblyBegin(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr); 259b117cd09SVijay Mahadevan ierr = MatAssemblyEnd(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr); 260b117cd09SVijay Mahadevan /* get the list of nodes on boundary so that we can enforce dirichlet conditions strongly */ 261b117cd09SVijay Mahadevan //ierr = DMMoabCheckBoundaryVertices(dm2,connc.size(),&connc[0],dbdry);CHKERRQ(ierr); 262e882eb38SVijay Mahadevan for (unsigned i=0; i < connc.size(); i++) { 263*3f1c6e43SVijay Mahadevan if (dbdry[i]) { /* dirichlet node */ 264b117cd09SVijay Mahadevan /* think about strongly imposing dirichlet */ 265b117cd09SVijay Mahadevan //bndrows.push_back(dofsc[i]); 266b117cd09SVijay Mahadevan 267b117cd09SVijay Mahadevan ierr = MatSetValues(*interpl, 1, &dofsc[i], connp.size(), &dofsp[0], &values_phi[0], INSERT_VALUES);CHKERRQ(ierr); 268b117cd09SVijay Mahadevan //values_phi[0]=1.0; 269b117cd09SVijay Mahadevan //ierr = MatSetValues(*interpl, 1, &dofsc[i], 1, &dofsc[i], &values_phi[0], INSERT_VALUES);CHKERRQ(ierr); 270b117cd09SVijay Mahadevan } 271b117cd09SVijay Mahadevan } 272b117cd09SVijay Mahadevan 273b117cd09SVijay Mahadevan ierr = MatAssemblyBegin(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr); 274b117cd09SVijay Mahadevan ierr = MatAssemblyEnd(*interpl,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr); 275b117cd09SVijay Mahadevan } 276b117cd09SVijay Mahadevan 277b117cd09SVijay Mahadevan //get interpolation weights 278b117cd09SVijay Mahadevan //ierr = Compute_Quad4_Basis(pcoords, 1, coord, values_phi);CHKERRQ(ierr); 279e882eb38SVijay Mahadevan // for (int j=0;j<dofs_per_element; j++) 280e882eb38SVijay Mahadevan // std::cout<<"values "<<values_phi[j]<<std::endl; 281b117cd09SVijay Mahadevan 282b117cd09SVijay Mahadevan //get row and column indices, zero weights are ignored 283b117cd09SVijay Mahadevan /* 284b117cd09SVijay Mahadevan int nz_ind = 0; 285b117cd09SVijay Mahadevan idx = dmb2->vowned->index(vhandle); 286b117cd09SVijay Mahadevan for (int j=0;j<dofs_per_element; j++){ 287b117cd09SVijay Mahadevan idy[nz_ind] = dmb1->vowned->index(connectivity[j]); 288b117cd09SVijay 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]); 289b117cd09SVijay Mahadevan //values_phi[nz_ind] = values_phi[j]; 290b117cd09SVijay Mahadevan nz_ind = nz_ind+1; 291b117cd09SVijay Mahadevan } 292b117cd09SVijay Mahadevan */ 293b117cd09SVijay Mahadevan 294b117cd09SVijay Mahadevan //ierr = MatSetValues(*interpl, nz_ind, idy, 1, &idx, values_phi, INSERT_VALUES);CHKERRQ(ierr); 295b117cd09SVijay Mahadevan //ierr = MatSetValues(*interpl, connp.size(), dofsp, connc.size(), dofsc, &values_phi[0], INSERT_VALUES);CHKERRQ(ierr); 296b117cd09SVijay Mahadevan } 297b117cd09SVijay Mahadevan 298b117cd09SVijay 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]); 299b117cd09SVijay Mahadevan //ierr = MatZeroRows(*interpl, bndrows.size(), &bndrows[0], 1.0, 0, 0);CHKERRQ(ierr); 300b117cd09SVijay Mahadevan 301b117cd09SVijay Mahadevan ierr = MatAssemblyBegin(*interpl,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 302b117cd09SVijay Mahadevan ierr = MatAssemblyEnd(*interpl,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 303b117cd09SVijay Mahadevan PetscFunctionReturn(0); 304b117cd09SVijay Mahadevan } 305b117cd09SVijay Mahadevan 306b117cd09SVijay Mahadevan #undef __FUNCT__ 307b117cd09SVijay Mahadevan #define __FUNCT__ "DMCreateInjection_Moab" 308b117cd09SVijay Mahadevan /*@ 309b117cd09SVijay Mahadevan DMCreateInjection_Moab - Generate a multi-level uniform refinement hierarchy 310b117cd09SVijay Mahadevan by succesively refining a coarse mesh, already defined in the DM object 311b117cd09SVijay Mahadevan provided by the user. 312b117cd09SVijay Mahadevan 313b117cd09SVijay Mahadevan Collective on MPI_Comm 314b117cd09SVijay Mahadevan 315b117cd09SVijay Mahadevan Input Parameter: 316b117cd09SVijay Mahadevan . dmb - The DMMoab object 317b117cd09SVijay Mahadevan 318b117cd09SVijay Mahadevan Output Parameter: 319b117cd09SVijay Mahadevan . nlevels - The number of levels of refinement needed to generate the hierarchy 320b117cd09SVijay Mahadevan + ldegrees - The degree of refinement at each level in the hierarchy 321b117cd09SVijay Mahadevan 322b117cd09SVijay Mahadevan Level: beginner 323b117cd09SVijay Mahadevan 324b117cd09SVijay Mahadevan .keywords: DMMoab, create, refinement 325b117cd09SVijay Mahadevan @*/ 326b117cd09SVijay Mahadevan PetscErrorCode DMCreateInjection_Moab(DM dm1,DM dm2,VecScatter* ctx) 327b117cd09SVijay Mahadevan { 328e882eb38SVijay Mahadevan //DM_Moab *dmmoab; 329b117cd09SVijay Mahadevan 330b117cd09SVijay Mahadevan PetscFunctionBegin; 331b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm1,DM_CLASSID,1); 332b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm2,DM_CLASSID,2); 333e882eb38SVijay Mahadevan //dmmoab = (DM_Moab*)(dm1)->data; 334b117cd09SVijay Mahadevan 335b117cd09SVijay Mahadevan PetscPrintf(PETSC_COMM_WORLD, "[DMCreateInjection_Moab] :: Placeholder\n"); 336b117cd09SVijay Mahadevan PetscFunctionReturn(0); 337b117cd09SVijay Mahadevan } 338b117cd09SVijay Mahadevan 339b117cd09SVijay Mahadevan 340b117cd09SVijay Mahadevan #undef __FUNCT__ 341b117cd09SVijay Mahadevan #define __FUNCT__ "DM_UMR_Moab_Private" 342b117cd09SVijay Mahadevan PetscErrorCode DM_UMR_Moab_Private(DM dm,MPI_Comm comm,PetscBool refine,DM *dmref) 343b117cd09SVijay Mahadevan { 344b117cd09SVijay Mahadevan PetscErrorCode ierr; 345e882eb38SVijay Mahadevan PetscInt i,dim; 346b117cd09SVijay Mahadevan DM dm2; 347b117cd09SVijay Mahadevan moab::ErrorCode merr; 348b117cd09SVijay Mahadevan DM_Moab *dmb = (DM_Moab*)dm->data,*dd2; 349b117cd09SVijay Mahadevan 350b117cd09SVijay Mahadevan PetscFunctionBegin; 351b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 352b117cd09SVijay Mahadevan PetscValidPointer(dmref,4); 353b117cd09SVijay Mahadevan 354b117cd09SVijay Mahadevan if ( (dmb->hlevel == dmb->nhlevels && refine) || (dmb->hlevel == 0 && !refine) ) { 355e882eb38SVijay 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); 356e882eb38SVijay Mahadevan if (dmb->hlevel-1 < 0 && !refine) PetscInfo1(NULL,"Invalid multigrid coarsen hierarchy level specified (%D). Creating a NULL object.\n",dmb->hlevel-1); 357e882eb38SVijay Mahadevan *dmref = PETSC_NULL; 358b117cd09SVijay Mahadevan PetscFunctionReturn(0); 359b117cd09SVijay Mahadevan } 360b117cd09SVijay Mahadevan 361b117cd09SVijay Mahadevan ierr = DMMoabCreate(PetscObjectComm((PetscObject)dm), &dm2);CHKERRQ(ierr); 362b117cd09SVijay Mahadevan dd2 = (DM_Moab*)dm2->data; 363b117cd09SVijay Mahadevan 364b117cd09SVijay Mahadevan dd2->mbiface = dmb->mbiface; 365b117cd09SVijay Mahadevan dd2->pcomm = dmb->pcomm; 366b117cd09SVijay Mahadevan dd2->icreatedinstance = PETSC_FALSE; 367b117cd09SVijay Mahadevan 368b117cd09SVijay Mahadevan /* set the new level based on refinement/coarsening */ 369b117cd09SVijay Mahadevan if (refine) { 370b117cd09SVijay Mahadevan dd2->hlevel=dmb->hlevel+1; 371b117cd09SVijay Mahadevan } 372b117cd09SVijay Mahadevan else { 373b117cd09SVijay Mahadevan dd2->hlevel=dmb->hlevel-1; 374b117cd09SVijay Mahadevan } 375b117cd09SVijay Mahadevan 376b117cd09SVijay Mahadevan /* Copy the multilevel hierarchy pointers in MOAB */ 377b117cd09SVijay Mahadevan dd2->hierarchy = dmb->hierarchy; 378b117cd09SVijay Mahadevan dd2->nhlevels = dmb->nhlevels; 379b117cd09SVijay Mahadevan ierr = PetscMalloc1(dd2->nhlevels+1,&dd2->hsets);CHKERRQ(ierr); 380b117cd09SVijay Mahadevan for (i=0; i<=dd2->nhlevels; i++) { 381b117cd09SVijay Mahadevan dd2->hsets[i]=dmb->hsets[i]; 382b117cd09SVijay Mahadevan } 383b117cd09SVijay Mahadevan dd2->fileset = dd2->hsets[dd2->hlevel]; 384b117cd09SVijay Mahadevan 385b117cd09SVijay Mahadevan /* do the remaining initializations for DMMoab */ 386b117cd09SVijay Mahadevan dd2->bs = dmb->bs; 387b117cd09SVijay Mahadevan dd2->numFields = dmb->numFields; 388b117cd09SVijay Mahadevan dd2->rw_dbglevel = dmb->rw_dbglevel; 389b117cd09SVijay Mahadevan dd2->partition_by_rank = dmb->partition_by_rank; 390b117cd09SVijay Mahadevan ierr = PetscStrcpy(dd2->extra_read_options, dmb->extra_read_options);CHKERRQ(ierr); 391b117cd09SVijay Mahadevan ierr = PetscStrcpy(dd2->extra_write_options, dmb->extra_write_options);CHKERRQ(ierr); 392b117cd09SVijay Mahadevan dd2->read_mode = dmb->read_mode; 393b117cd09SVijay Mahadevan dd2->write_mode = dmb->write_mode; 394b117cd09SVijay Mahadevan 395b117cd09SVijay Mahadevan /* set global ID tag handle */ 396b117cd09SVijay Mahadevan ierr = DMMoabSetLocalToGlobalTag(dm2, dmb->ltog_tag);CHKERRQ(ierr); 397b117cd09SVijay Mahadevan 398b117cd09SVijay Mahadevan merr = dd2->mbiface->tag_get_handle(MATERIAL_SET_TAG_NAME, dd2->material_tag);MBERRNM(merr); 399b117cd09SVijay Mahadevan 400b117cd09SVijay Mahadevan ierr = DMSetOptionsPrefix(dm2,((PetscObject)dm)->prefix);CHKERRQ(ierr); 401b117cd09SVijay Mahadevan ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 402b117cd09SVijay Mahadevan ierr = DMSetDimension(dm2,dim);CHKERRQ(ierr); 403b117cd09SVijay Mahadevan 404b117cd09SVijay Mahadevan /* allow overloaded (user replaced) operations to be inherited by refinement clones */ 405b117cd09SVijay Mahadevan dm2->ops->creatematrix = dm->ops->creatematrix; 406b117cd09SVijay Mahadevan 407b117cd09SVijay Mahadevan /* copy fill information if given */ 408b117cd09SVijay Mahadevan ierr = DMMoabSetBlockFills(dm2, dmb->dfill, dmb->ofill);CHKERRQ(ierr); 409b117cd09SVijay Mahadevan 410b117cd09SVijay Mahadevan /* copy vector type information */ 411b117cd09SVijay Mahadevan ierr = DMSetMatType(dm2,dm->mattype);CHKERRQ(ierr); 412b117cd09SVijay Mahadevan ierr = DMSetVecType(dm2,dm->vectype);CHKERRQ(ierr); 413*3f1c6e43SVijay Mahadevan dd2->numFields = dmb->numFields; 414*3f1c6e43SVijay Mahadevan if (dmb->numFields) { 415b117cd09SVijay Mahadevan ierr = DMMoabSetFieldNames(dm2,dmb->numFields,dmb->fieldNames);CHKERRQ(ierr); 416*3f1c6e43SVijay Mahadevan } 417b117cd09SVijay Mahadevan 418b117cd09SVijay Mahadevan ierr = DMSetFromOptions(dm2);CHKERRQ(ierr); 419b117cd09SVijay Mahadevan 420b117cd09SVijay Mahadevan /* recreate Dof numbering for the refined DM and make sure the distribution is correctly populated */ 421b117cd09SVijay Mahadevan ierr = DMSetUp(dm2);CHKERRQ(ierr); 422b117cd09SVijay Mahadevan 423b117cd09SVijay Mahadevan *dmref = dm2; 424b117cd09SVijay Mahadevan PetscFunctionReturn(0); 425b117cd09SVijay Mahadevan } 426b117cd09SVijay Mahadevan 427b117cd09SVijay Mahadevan 428b117cd09SVijay Mahadevan #undef __FUNCT__ 429b117cd09SVijay Mahadevan #define __FUNCT__ "DMRefine_Moab" 430b117cd09SVijay Mahadevan /*@ 431b117cd09SVijay Mahadevan DMRefine_Moab - Generate a multi-level uniform refinement hierarchy 432b117cd09SVijay Mahadevan by succesively refining a coarse mesh, already defined in the DM object 433b117cd09SVijay Mahadevan provided by the user. 434b117cd09SVijay Mahadevan 435e882eb38SVijay Mahadevan Collective on DM 436b117cd09SVijay Mahadevan 437b117cd09SVijay Mahadevan Input Parameter: 438e882eb38SVijay Mahadevan + dm - The DMMoab object 439e882eb38SVijay Mahadevan - comm - the communicator to contain the new DM object (or MPI_COMM_NULL) 440b117cd09SVijay Mahadevan 441b117cd09SVijay Mahadevan Output Parameter: 442e882eb38SVijay Mahadevan . dmf - the refined DM, or NULL 443b117cd09SVijay Mahadevan 444e882eb38SVijay Mahadevan Note: If no refinement was done, the return value is NULL 445e882eb38SVijay Mahadevan 446e882eb38SVijay Mahadevan Level: developer 447b117cd09SVijay Mahadevan 448b117cd09SVijay Mahadevan .keywords: DMMoab, create, refinement 449b117cd09SVijay Mahadevan @*/ 450b117cd09SVijay Mahadevan PetscErrorCode DMRefine_Moab(DM dm,MPI_Comm comm,DM* dmf) 451b117cd09SVijay Mahadevan { 452b117cd09SVijay Mahadevan PetscErrorCode ierr; 453b117cd09SVijay Mahadevan 454b117cd09SVijay Mahadevan PetscFunctionBegin; 455b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 456b117cd09SVijay Mahadevan 457b117cd09SVijay Mahadevan ierr = DM_UMR_Moab_Private(dm,comm,PETSC_TRUE,dmf);CHKERRQ(ierr); 458b117cd09SVijay Mahadevan PetscFunctionReturn(0); 459b117cd09SVijay Mahadevan } 460b117cd09SVijay Mahadevan 461b117cd09SVijay Mahadevan #undef __FUNCT__ 462b117cd09SVijay Mahadevan #define __FUNCT__ "DMCoarsen_Moab" 463b117cd09SVijay Mahadevan /*@ 464b117cd09SVijay Mahadevan DMCoarsen_Moab - Generate a multi-level uniform refinement hierarchy 465b117cd09SVijay Mahadevan by succesively refining a coarse mesh, already defined in the DM object 466b117cd09SVijay Mahadevan provided by the user. 467b117cd09SVijay Mahadevan 468e882eb38SVijay Mahadevan Collective on DM 469b117cd09SVijay Mahadevan 470b117cd09SVijay Mahadevan Input Parameter: 471e882eb38SVijay Mahadevan . dm - The DMMoab object 472e882eb38SVijay Mahadevan - comm - the communicator to contain the new DM object (or MPI_COMM_NULL) 473b117cd09SVijay Mahadevan 474b117cd09SVijay Mahadevan Output Parameter: 475e882eb38SVijay Mahadevan . dmf - the coarsened DM, or NULL 476b117cd09SVijay Mahadevan 477e882eb38SVijay Mahadevan Note: If no coarsening was done, the return value is NULL 478b117cd09SVijay Mahadevan 479e882eb38SVijay Mahadevan Level: developer 480e882eb38SVijay Mahadevan 481e882eb38SVijay Mahadevan .keywords: DMMoab, create, coarsening 482b117cd09SVijay Mahadevan @*/ 483b117cd09SVijay Mahadevan PetscErrorCode DMCoarsen_Moab(DM dm,MPI_Comm comm,DM* dmc) 484b117cd09SVijay Mahadevan { 485b117cd09SVijay Mahadevan PetscErrorCode ierr; 486b117cd09SVijay Mahadevan 487b117cd09SVijay Mahadevan PetscFunctionBegin; 488b117cd09SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 489b117cd09SVijay Mahadevan 490b117cd09SVijay Mahadevan ierr = DM_UMR_Moab_Private(dm,comm,PETSC_FALSE,dmc);CHKERRQ(ierr); 491b117cd09SVijay Mahadevan PetscFunctionReturn(0); 492b117cd09SVijay Mahadevan } 493