PetscErrorCode SNESSetUp_FAS(SNES snes) { SNES_FAS *fas = (SNES_FAS*) snes->data; PetscErrorCode ierr; PetscInt dm_levels; Vec vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */ SNES next; PetscBool isFine; SNESLineSearch linesearch; SNESLineSearch slinesearch; void *lsprectx,*lspostctx; PetscErrorCode (*precheck)(SNESLineSearch,Vec,Vec,PetscBool*,void*); PetscErrorCode (*postcheck)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*); PetscFunctionBegin; ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); if (fas->usedmfornumberoflevels && isFine) { ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr); dm_levels++; if (dm_levels > fas->levels) { /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESFASSetLevels_FAS*/ vec_sol = snes->vec_sol; vec_func = snes->vec_func; vec_sol_update = snes->vec_sol_update; vec_rhs = snes->vec_rhs; snes->vec_sol = NULL; snes->vec_func = NULL; snes->vec_sol_update = NULL; snes->vec_rhs = NULL; /* reset the number of levels */ ierr = SNESFASSetLevels(snes,dm_levels,NULL);CHKERRQ(ierr); ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); snes->vec_sol = vec_sol; snes->vec_func = vec_func; snes->vec_rhs = vec_rhs; snes->vec_sol_update = vec_sol_update; } } ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); if (!isFine) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */ ierr = SNESSetWorkVecs(snes, 2);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ /* set up the smoothers if they haven't already been set up */ if (!fas->smoothd) { ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr); } if (snes->dm) { /* set the smoother DMs properly */ if (fas->smoothu) ierr = SNESSetDM(fas->smoothu, snes->dm);CHKERRQ(ierr); ierr = SNESSetDM(fas->smoothd, snes->dm);CHKERRQ(ierr); /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ if (next) { /* for now -- assume the DM and the evaluation functions have been set externally */ if (!next->dm) { ierr = DMCoarsen(snes->dm, PetscObjectComm((PetscObject)next), &next->dm);CHKERRQ(ierr); ierr = SNESSetDM(next, next->dm);CHKERRQ(ierr); } /* set the interpolation and restriction from the DM */ if (!fas->interpolate) { ierr = DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); if (!fas->restrct) { ierr = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr); fas->restrct = fas->interpolate; } } /* set the injection from the DM */ if (!fas->inject) { ierr = DMCreateInjection(next->dm, snes->dm, &fas->inject);CHKERRQ(ierr); } } } /*pass the smoother, function, and jacobian up to the next level if it's not user set already */ if (fas->galerkin) { if (next) { ierr = SNESSetFunction(next, NULL, SNESFASGalerkinDefaultFunction, next);CHKERRQ(ierr); } if (fas->smoothd && fas->level != fas->levels - 1) { ierr = SNESSetFunction(fas->smoothd, NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr); } if (fas->smoothu && fas->level != fas->levels - 1) { ierr = SNESSetFunction(fas->smoothu, NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr); } } /* sets the down (pre) smoother's default norm and sets it from options */ if (fas->smoothd) { if (fas->level == 0 && fas->levels != 1) { ierr = SNESSetNormSchedule(fas->smoothd, SNES_NORM_NONE);CHKERRQ(ierr); } else { ierr = SNESSetNormSchedule(fas->smoothd, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr); } ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)fas->smoothd);CHKERRQ(ierr); ierr = SNESSetFromOptions(fas->smoothd);CHKERRQ(ierr); ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr); ierr = SNESGetLineSearch(fas->smoothd,&slinesearch);CHKERRQ(ierr); ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr); ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr); ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr); ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr); ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr); fas->smoothd->vec_sol = snes->vec_sol; ierr = PetscObjectReference((PetscObject)snes->vec_sol);CHKERRQ(ierr); fas->smoothd->vec_sol_update = snes->vec_sol_update; ierr = PetscObjectReference((PetscObject)snes->vec_sol_update);CHKERRQ(ierr); fas->smoothd->vec_func = snes->vec_func; ierr = PetscObjectReference((PetscObject)snes->vec_func);CHKERRQ(ierr); if (fas->eventsmoothsetup) {ierr = PetscLogEventBegin(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} ierr = SNESSetUp(fas->smoothd);CHKERRQ(ierr); if (fas->eventsmoothsetup) {ierr = PetscLogEventEnd(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} } /* sets the up (post) smoother's default norm and sets it from options */ if (fas->smoothu) { if (fas->level != fas->levels - 1) { ierr = SNESSetNormSchedule(fas->smoothu, SNES_NORM_NONE);CHKERRQ(ierr); } else { ierr = SNESSetNormSchedule(fas->smoothu, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr); } ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)fas->smoothu);CHKERRQ(ierr); ierr = SNESSetFromOptions(fas->smoothu);CHKERRQ(ierr); ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr); ierr = SNESGetLineSearch(fas->smoothu,&slinesearch);CHKERRQ(ierr); ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr); ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr); ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr); ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr); ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr); fas->smoothu->vec_sol = snes->vec_sol; ierr = PetscObjectReference((PetscObject)snes->vec_sol);CHKERRQ(ierr); fas->smoothu->vec_sol_update = snes->vec_sol_update; ierr = PetscObjectReference((PetscObject)snes->vec_sol_update);CHKERRQ(ierr); fas->smoothu->vec_func = snes->vec_func; ierr = PetscObjectReference((PetscObject)snes->vec_func);CHKERRQ(ierr); if (fas->eventsmoothsetup) {ierr = PetscLogEventBegin(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} ierr = SNESSetUp(fas->smoothu);CHKERRQ(ierr); if (fas->eventsmoothsetup) {ierr = PetscLogEventEnd(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} } if (next) { /* gotta set up the solution vector for this to work */ if (!next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_sol);CHKERRQ(ierr);} if (!next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_rhs);CHKERRQ(ierr);} ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)next);CHKERRQ(ierr); ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr); ierr = SNESGetLineSearch(fas->next,&slinesearch);CHKERRQ(ierr); ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr); ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr); ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr); ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr); ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr); ierr = SNESSetUp(next);CHKERRQ(ierr); } /* setup FAS work vectors */ if (fas->galerkin) { ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr); ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr); } PetscFunctionReturn(0); }
PetscErrorCode DMCoarsen_DA(DM da, MPI_Comm comm,DM *daref) { PetscErrorCode ierr; PetscInt M,N,P,i; DM da2; DM_DA *dd = (DM_DA*)da->data,*dd2; PetscFunctionBegin; PetscValidHeaderSpecific(da,DM_CLASSID,1); PetscValidPointer(daref,3); if (dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){ M = dd->M / dd->coarsen_x; } else { M = 1 + (dd->M - 1) / dd->coarsen_x; } if (dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){ if (dd->dim > 1) { N = dd->N / dd->coarsen_y; } else { N = 1; } } else { N = 1 + (dd->N - 1) / dd->coarsen_y; } if (dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0){ if (dd->dim > 2) { P = dd->P / dd->coarsen_z; } else { P = 1; } } else { P = 1 + (dd->P - 1) / dd->coarsen_z; } ierr = DMDACreate(((PetscObject)da)->comm,&da2);CHKERRQ(ierr); ierr = DMSetOptionsPrefix(da2,((PetscObject)da)->prefix);CHKERRQ(ierr); ierr = DMDASetDim(da2,dd->dim);CHKERRQ(ierr); ierr = DMDASetSizes(da2,M,N,P);CHKERRQ(ierr); ierr = DMDASetNumProcs(da2,dd->m,dd->n,dd->p);CHKERRQ(ierr); ierr = DMDASetBoundaryType(da2,dd->bx,dd->by,dd->bz);CHKERRQ(ierr); ierr = DMDASetDof(da2,dd->w);CHKERRQ(ierr); ierr = DMDASetStencilType(da2,dd->stencil_type);CHKERRQ(ierr); ierr = DMDASetStencilWidth(da2,dd->s);CHKERRQ(ierr); if (dd->dim == 3) { PetscInt *lx,*ly,*lz; ierr = PetscMalloc3(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly,dd->p,PetscInt,&lz);CHKERRQ(ierr); ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr); ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_y,dd->n,dd->ly,ly);CHKERRQ(ierr); ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bz == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_z,dd->p,dd->lz,lz);CHKERRQ(ierr); ierr = DMDASetOwnershipRanges(da2,lx,ly,lz);CHKERRQ(ierr); ierr = PetscFree3(lx,ly,lz);CHKERRQ(ierr); } else if (dd->dim == 2) { PetscInt *lx,*ly; ierr = PetscMalloc2(dd->m,PetscInt,&lx,dd->n,PetscInt,&ly);CHKERRQ(ierr); ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr); ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->by == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_y,dd->n,dd->ly,ly);CHKERRQ(ierr); ierr = DMDASetOwnershipRanges(da2,lx,ly,PETSC_NULL);CHKERRQ(ierr); ierr = PetscFree2(lx,ly);CHKERRQ(ierr); } else if (dd->dim == 1) { PetscInt *lx; ierr = PetscMalloc(dd->m*sizeof(PetscInt),&lx);CHKERRQ(ierr); ierr = DMDACoarsenOwnershipRanges(da,(PetscBool)(dd->bx == DMDA_BOUNDARY_PERIODIC || dd->interptype == DMDA_Q0),dd->s,dd->coarsen_x,dd->m,dd->lx,lx);CHKERRQ(ierr); ierr = DMDASetOwnershipRanges(da2,lx,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); ierr = PetscFree(lx);CHKERRQ(ierr); } dd2 = (DM_DA*)da2->data; /* allow overloaded (user replaced) operations to be inherited by refinement clones; why are only some inherited and not all? */ /* da2->ops->createinterpolation = da->ops->createinterpolation; copying this one causes trouble for DMSetVI */ da2->ops->creatematrix = da->ops->creatematrix; da2->ops->getcoloring = da->ops->getcoloring; dd2->interptype = dd->interptype; /* copy fill information if given */ if (dd->dfill) { ierr = PetscMalloc((dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->dfill);CHKERRQ(ierr); ierr = PetscMemcpy(dd2->dfill,dd->dfill,(dd->dfill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); } if (dd->ofill) { ierr = PetscMalloc((dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt),&dd2->ofill);CHKERRQ(ierr); ierr = PetscMemcpy(dd2->ofill,dd->ofill,(dd->ofill[dd->w]+dd->w+1)*sizeof(PetscInt));CHKERRQ(ierr); } /* copy the refine information */ dd2->coarsen_x = dd2->refine_x = dd->coarsen_x; dd2->coarsen_y = dd2->refine_y = dd->coarsen_y; dd2->coarsen_z = dd2->refine_z = dd->coarsen_z; /* copy vector type information */ ierr = PetscFree(da2->vectype);CHKERRQ(ierr); ierr = PetscStrallocpy(da->vectype,(char**)&da2->vectype);CHKERRQ(ierr); dd2->lf = dd->lf; dd2->lj = dd->lj; da2->leveldown = da->leveldown + 1; da2->levelup = da->levelup; ierr = DMSetFromOptions(da2);CHKERRQ(ierr); ierr = DMSetUp(da2);CHKERRQ(ierr); ierr = DMViewFromOptions(da2,"-dm_view");CHKERRQ(ierr); /* inject coordinates if they are set on the fine grid */ if (da->coordinates) { DM cdaf,cdac; Vec coordsc,coordsf; VecScatter inject; ierr = DMGetCoordinateDM(da,&cdaf);CHKERRQ(ierr); ierr = DMGetCoordinates(da,&coordsf);CHKERRQ(ierr); ierr = DMGetCoordinateDM(da2,&cdac);CHKERRQ(ierr); /* force creation of the coordinate vector */ ierr = DMDASetUniformCoordinates(da2,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr); ierr = DMGetCoordinates(da2,&coordsc);CHKERRQ(ierr); ierr = DMCreateInjection(cdac,cdaf,&inject);CHKERRQ(ierr); ierr = VecScatterBegin(inject,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); ierr = VecScatterEnd(inject ,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); ierr = VecScatterDestroy(&inject);CHKERRQ(ierr); } for (i=0; i<da->bs; i++) { const char *fieldname; ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr); ierr = DMDASetFieldName(da2,i,fieldname);CHKERRQ(ierr); } *daref = da2; PetscFunctionReturn(0); }
PetscErrorCode test1_DAInjection3d(PetscInt mx, PetscInt my, PetscInt mz) { PetscErrorCode ierr; DM dac,daf; PetscViewer vv; Vec ac,af; PetscInt periodicity; DMBoundaryType bx,by,bz; PetscFunctionBeginUser; bx = DM_BOUNDARY_NONE; by = DM_BOUNDARY_NONE; bz = DM_BOUNDARY_NONE; periodicity = 0; ierr = PetscOptionsGetInt(NULL,NULL,"-periodic", &periodicity, NULL);CHKERRQ(ierr); if (periodicity==1) { bx = DM_BOUNDARY_PERIODIC; } else if (periodicity==2) { by = DM_BOUNDARY_PERIODIC; } else if (periodicity==3) { bz = DM_BOUNDARY_PERIODIC; } ierr = DMDACreate3d(PETSC_COMM_WORLD, bx,by,bz, DMDA_STENCIL_BOX,mx+1, my+1,mz+1,PETSC_DECIDE, PETSC_DECIDE,PETSC_DECIDE,1, /* 1 dof */ 1, /* stencil = 1 */NULL,NULL,NULL,&daf);CHKERRQ(ierr); ierr = DMSetFromOptions(daf);CHKERRQ(ierr); ierr = DMSetUp(daf);CHKERRQ(ierr); ierr = DMCoarsen(daf,MPI_COMM_NULL,&dac);CHKERRQ(ierr); ierr = DMDASetUniformCoordinates(dac, -1.0,1.0, -1.0,1.0, -1.0,1.0);CHKERRQ(ierr); ierr = DMDASetUniformCoordinates(daf, -1.0,1.0, -1.0,1.0, -1.0,1.0);CHKERRQ(ierr); { DM cdaf,cdac; Vec coordsc,coordsf,coordsf2; Mat inject; VecScatter vscat; Mat interp; PetscReal norm; ierr = DMGetCoordinateDM(dac,&cdac);CHKERRQ(ierr); ierr = DMGetCoordinateDM(daf,&cdaf);CHKERRQ(ierr); ierr = DMGetCoordinates(dac,&coordsc);CHKERRQ(ierr); ierr = DMGetCoordinates(daf,&coordsf);CHKERRQ(ierr); ierr = DMCreateInjection(cdac,cdaf,&inject);CHKERRQ(ierr); ierr = MatScatterGetVecScatter(inject,&vscat);CHKERRQ(ierr); ierr = VecScatterBegin(vscat,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); ierr = VecScatterEnd(vscat ,coordsf,coordsc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); ierr = MatDestroy(&inject);CHKERRQ(ierr); ierr = DMCreateInterpolation(cdac,cdaf,&interp,NULL);CHKERRQ(ierr); ierr = VecDuplicate(coordsf,&coordsf2);CHKERRQ(ierr); ierr = MatInterpolate(interp,coordsc,coordsf2);CHKERRQ(ierr); ierr = VecAXPY(coordsf2,-1.0,coordsf);CHKERRQ(ierr); ierr = VecNorm(coordsf2,NORM_MAX,&norm);CHKERRQ(ierr); /* The fine coordinates are only reproduced in certain cases */ if (!bx && !by && !bz && norm > PETSC_SQRT_MACHINE_EPSILON) {ierr = PetscPrintf(PETSC_COMM_WORLD,"Norm %g\n",(double)norm);CHKERRQ(ierr);} ierr = VecDestroy(&coordsf2);CHKERRQ(ierr); ierr = MatDestroy(&interp);CHKERRQ(ierr); } if (0) { ierr = DMCreateGlobalVector(dac,&ac);CHKERRQ(ierr); ierr = VecZeroEntries(ac);CHKERRQ(ierr); ierr = DMCreateGlobalVector(daf,&af);CHKERRQ(ierr); ierr = VecZeroEntries(af);CHKERRQ(ierr); ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD, "dac_7.vtk", &vv);CHKERRQ(ierr); ierr = PetscViewerPushFormat(vv, PETSC_VIEWER_ASCII_VTK);CHKERRQ(ierr); ierr = DMView(dac, vv);CHKERRQ(ierr); ierr = VecView(ac, vv);CHKERRQ(ierr); ierr = PetscViewerPopFormat(vv);CHKERRQ(ierr); ierr = PetscViewerDestroy(&vv);CHKERRQ(ierr); ierr = PetscViewerASCIIOpen(PETSC_COMM_WORLD, "daf_7.vtk", &vv);CHKERRQ(ierr); ierr = PetscViewerPushFormat(vv, PETSC_VIEWER_ASCII_VTK);CHKERRQ(ierr); ierr = DMView(daf, vv);CHKERRQ(ierr); ierr = VecView(af, vv);CHKERRQ(ierr); ierr = PetscViewerPopFormat(vv);CHKERRQ(ierr); ierr = PetscViewerDestroy(&vv);CHKERRQ(ierr); ierr = VecDestroy(&ac);CHKERRQ(ierr); ierr = VecDestroy(&af);CHKERRQ(ierr); } ierr = DMDestroy(&dac);CHKERRQ(ierr); ierr = DMDestroy(&daf);CHKERRQ(ierr); PetscFunctionReturn(0); }
/* DMCoarsen_SNESVI - Computes the regular coarsened DM then computes additional information about its inactive set */ PetscErrorCode DMCoarsen_SNESVI(DM dm1,MPI_Comm comm,DM *dm2) { PetscErrorCode ierr; PetscContainer isnes; DM_SNESVI *dmsnesvi1; Vec finemarked,coarsemarked; IS inactive; VecScatter inject; const PetscInt *index; PetscInt n,k,cnt = 0,rstart,*coarseindex; PetscScalar *marked; PetscFunctionBegin; ierr = PetscObjectQuery((PetscObject)dm1,"VI",(PetscObject *)&isnes);CHKERRQ(ierr); if (!isnes) SETERRQ(((PetscObject)dm1)->comm,PETSC_ERR_PLIB,"Composed VI data structure is missing"); ierr = PetscContainerGetPointer(isnes,(void**)&dmsnesvi1);CHKERRQ(ierr); /* get the original coarsen */ ierr = (*dmsnesvi1->coarsen)(dm1,comm,dm2);CHKERRQ(ierr); /* not sure why this extra reference is needed, but without the dm2 disappears too early */ ierr = PetscObjectReference((PetscObject)*dm2);CHKERRQ(ierr); /* need to set back global vectors in order to use the original injection */ ierr = DMClearGlobalVectors(dm1);CHKERRQ(ierr); dm1->ops->createglobalvector = dmsnesvi1->createglobalvector; ierr = DMCreateGlobalVector(dm1,&finemarked);CHKERRQ(ierr); ierr = DMCreateGlobalVector(*dm2,&coarsemarked);CHKERRQ(ierr); /* fill finemarked with locations of inactive points */ ierr = ISGetIndices(dmsnesvi1->inactive,&index);CHKERRQ(ierr); ierr = ISGetLocalSize(dmsnesvi1->inactive,&n);CHKERRQ(ierr); ierr = VecSet(finemarked,0.0);CHKERRQ(ierr); for (k=0;k<n;k++){ ierr = VecSetValue(finemarked,index[k],1.0,INSERT_VALUES);CHKERRQ(ierr); } ierr = VecAssemblyBegin(finemarked);CHKERRQ(ierr); ierr = VecAssemblyEnd(finemarked);CHKERRQ(ierr); ierr = DMCreateInjection(*dm2,dm1,&inject);CHKERRQ(ierr); ierr = VecScatterBegin(inject,finemarked,coarsemarked,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); ierr = VecScatterEnd(inject,finemarked,coarsemarked,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); ierr = VecScatterDestroy(&inject);CHKERRQ(ierr); /* create index set list of coarse inactive points from coarsemarked */ ierr = VecGetLocalSize(coarsemarked,&n);CHKERRQ(ierr); ierr = VecGetOwnershipRange(coarsemarked,&rstart,PETSC_NULL);CHKERRQ(ierr); ierr = VecGetArray(coarsemarked,&marked);CHKERRQ(ierr); for (k=0; k<n; k++) { if (marked[k] != 0.0) cnt++; } ierr = PetscMalloc(cnt*sizeof(PetscInt),&coarseindex);CHKERRQ(ierr); cnt = 0; for (k=0; k<n; k++) { if (marked[k] != 0.0) coarseindex[cnt++] = k + rstart; } ierr = VecRestoreArray(coarsemarked,&marked);CHKERRQ(ierr); ierr = ISCreateGeneral(((PetscObject)coarsemarked)->comm,cnt,coarseindex,PETSC_OWN_POINTER,&inactive);CHKERRQ(ierr); ierr = DMClearGlobalVectors(dm1);CHKERRQ(ierr); dm1->ops->createglobalvector = DMCreateGlobalVector_SNESVI; ierr = DMSetVI(*dm2,inactive);CHKERRQ(ierr); ierr = VecDestroy(&finemarked);CHKERRQ(ierr); ierr = VecDestroy(&coarsemarked);CHKERRQ(ierr); ierr = ISDestroy(&inactive);CHKERRQ(ierr); PetscFunctionReturn(0); }