<div dir="ltr"><div dir="ltr">On Wed, Feb 19, 2025 at 11:18 AM Preda Silvia via petsc-users <<a href="mailto:petsc-users@mcs.anl.gov">petsc-users@mcs.anl.gov</a>> wrote:</div><div class="gmail_quote gmail_quote_container"><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div class="msg2285903082881261691">





<div lang="EN-US" style="overflow-wrap: break-word;">
<div class="m_2285903082881261691WordSection1">
<p class="MsoNormal"><span lang="IT">Hi,<u></u><u></u></span></p>
<p class="MsoNormal"><span lang="IT"><u></u> <u></u></span></p>
<p class="MsoNormal">I’m using the toycode below to manage adaptivity using a DMForest, based on p4est.<u></u><u></u></p>
<p class="MsoNormal"><u></u> <u></u></p>
<p class="MsoNormal">The code simply performs three adaptivity steps on a square grid [0,1]x[0,1], uniformly refined at level 2 at the beginning. The minimum and the maximum level of refinement are set to 1 and 5, respectively. During the adaptivity procedure,
 a quadrant is refined if its centroid is inside the circle of radius 0.25, centred in (0.5,0.5).
<u></u><u></u></p>
<p class="MsoNormal"><u></u> <u></u></p>
<p class="MsoNormal">I’m facing two main issues:<u></u><u></u></p>
<ul style="margin-top:0cm" type="disc">
<li class="m_2285903082881261691MsoListParagraph" style="margin-left:0cm">As far as I understood, when the label asks to refine a quadrant which is already at the maximum level, the adaption procedure errors out instead of ignoring the request. Thus I need
 to check the quadrant level before setting the adaptivity label. How can I get access to the quadrant-level of a cell?<u></u><u></u></li><li class="m_2285903082881261691MsoListParagraph" style="margin-left:0cm">Once the mesh is adapted, data need to be projected on the new mesh and in my method I’ll need to write an elaborate ad-hoc procedure for that. How can I access the map that provides
 the correspondence between the indexes of outcoming quadrants and the indexes of incoming ones?</li></ul></div></div></div></blockquote><div><br></div><div>Toby understands this better than me. Toby, two questions:</div><div><br></div><div>  1) Can we catch that error and just return?</div><div><br></div><div>  2) Will DMProjectField() work in the case of multiple refinements like this?</div><div><br></div><div>  Thanks,</div><div><br></div><div>      Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><div class="msg2285903082881261691"><div lang="EN-US" style="overflow-wrap: break-word;"><div class="m_2285903082881261691WordSection1">
<p class="MsoNormal">Below is the code, to be run with these options: <u></u><u></u></p>
<p class="MsoNormal">-dm_type p4est <u></u><u></u></p>
<p class="MsoNormal">-dm_forest_topology brick <u></u><u></u></p>
<p class="MsoNormal">-dm_p4est_brick_size 1,1 <u></u><u></u></p>
<p class="MsoNormal">-dm_view vtk:brick.vtu<u></u><u></u></p>
<p class="MsoNormal">-dm_forest_initial_refinement 2 <u></u><u></u></p>
<p class="MsoNormal">-dm_forest_minimum_refinement 1 <u></u><u></u></p>
<p class="MsoNormal">-dm_forest_maximum_refinement 5<u></u><u></u></p>
<div style="border-top:none;border-right:none;border-left:none;border-bottom:2.25pt double windowtext;padding:0cm 0cm 1pt">
<p class="MsoNormal" style="border:none;padding:0cm"><u></u><u></u></p>
</div>
<p class="MsoNormal"><u></u> <u></u></p>
<div style="border-top:none;border-right:none;border-left:none;border-bottom:2.25pt double windowtext;padding:0cm 0cm 1pt">
<p class="MsoNormal" style="border:none;padding:0cm">static char help[] = "Create and view a forest mesh\n\n";<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm"><u></u> <u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">#include <petscdmforest.h><u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">#include <petscdmplex.h><u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">#include <petscoptions.h><u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm"><u></u> <u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">static PetscErrorCode CreateAdaptLabel(DM dm, DM dmConv, PetscInt *nAdaptLoc)
<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">{<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  DMLabel  adaptLabel;<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  PetscInt cStart, cEnd, c;<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm"><u></u> <u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  PetscFunctionBeginUser;<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  PetscCall(DMGetCoordinatesLocalSetUp(dm));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  PetscCall(DMCreateLabel(dm, "adaptLabel"));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  PetscCall(DMGetLabel(dm, "adaptLabel", &adaptLabel));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  PetscCall(DMForestGetCellChart(dm, &cStart, &cEnd));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  for (c = cStart; c < cEnd; ++c) {<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">    PetscReal centroid[3], volume, x, y;<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">    PetscCall(DMPlexComputeCellGeometryFVM(dmConv, c, &volume, centroid, NULL));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">    <span lang="FR">x = centroid[0];<u></u><u></u></span></p>
<p class="MsoNormal" style="border:none;padding:0cm"><span lang="FR">    y = centroid[1];<u></u><u></u></span></p>
<p class="MsoNormal" style="border:none;padding:0cm"><span lang="FR">    if (std::sqrt((x-0.5)*(x-0.5)+(y-0.5)*(y-0.5))<0.25) {<u></u><u></u></span></p>
<p class="MsoNormal" style="border:none;padding:0cm"><span lang="FR">      PetscCall(DMLabelSetValue(adaptLabel, c, DM_ADAPT_REFINE));<u></u><u></u></span></p>
<p class="MsoNormal" style="border:none;padding:0cm"><span lang="FR">      ++nAdaptLoc[0];<u></u><u></u></span></p>
<p class="MsoNormal" style="border:none;padding:0cm"><span lang="FR">    } else {<u></u><u></u></span></p>
<p class="MsoNormal" style="border:none;padding:0cm"><span lang="FR">      PetscCall(DMLabelSetValue(adaptLabel, c, DM_ADAPT_KEEP));<u></u><u></u></span></p>
<p class="MsoNormal" style="border:none;padding:0cm"><span lang="FR">      ++nAdaptLoc[1];<u></u><u></u></span></p>
<p class="MsoNormal" style="border:none;padding:0cm"><span lang="FR">    }<u></u><u></u></span></p>
<p class="MsoNormal" style="border:none;padding:0cm"><span lang="FR">  }<u></u><u></u></span></p>
<p class="MsoNormal" style="border:none;padding:0cm"><span lang="FR">  PetscFunctionReturn(PETSC_SUCCESS);<u></u><u></u></span></p>
<p class="MsoNormal" style="border:none;padding:0cm"><span lang="FR">}<u></u><u></u></span></p>
<p class="MsoNormal" style="border:none;padding:0cm"><span lang="FR"><u></u> <u></u></span></p>
<p class="MsoNormal" style="border:none;padding:0cm"><span lang="FR">static PetscErrorCode ForestToPlex(DM *dm, DM *dmConv)<u></u><u></u></span></p>
<p class="MsoNormal" style="border:none;padding:0cm"><span lang="FR">{<u></u><u></u></span></p>
<p class="MsoNormal" style="border:none;padding:0cm"><span lang="FR">  PetscFunctionBeginUser;<u></u><u></u></span></p>
<p class="MsoNormal" style="border:none;padding:0cm"><span lang="FR">  PetscCall(DMConvert(*dm, DMPLEX, dmConv));<u></u><u></u></span></p>
<p class="MsoNormal" style="border:none;padding:0cm"><span lang="FR">  PetscCall(DMLocalizeCoordinates(*dmConv));<u></u><u></u></span></p>
<p class="MsoNormal" style="border:none;padding:0cm"><span lang="FR">  </span>PetscCall(DMViewFromOptions(*dmConv, NULL, "-dm_conv_view"));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  PetscCall(DMPlexCheckCellShape(*dmConv, PETSC_FALSE, PETSC_DETERMINE));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  PetscFunctionReturn(PETSC_SUCCESS);<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">}<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm"><u></u> <u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">static PetscErrorCode AdaptMesh(DM *dm)<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">{<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  DM              dmCur = *dm;<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  PetscBool       hasLabel=PETSC_FALSE, adapt=PETSC_TRUE;<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  PetscInt        adaptIter=0, maxAdaptIter=3;<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  <u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  PetscFunctionBeginUser;<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  while (adapt) {<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">    DM       dmAdapt;<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">    DMLabel  adaptLabel;<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">    PetscInt nAdaptLoc[2]={0,0}, nAdapt[2]={0,0};<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">    <u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">    ++adaptIter;<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">    PetscCall(PetscPrintf(PETSC_COMM_SELF,"\nADAPT ITER %d\n",adaptIter));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">    <u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">    DM dmConv;<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">    PetscCall(ForestToPlex(&dmCur,&dmConv));
<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">    PetscCall(CreateAdaptLabel(dmCur,dmConv,nAdaptLoc));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">    PetscCallMPI(MPIU_Allreduce(&nAdaptLoc, &nAdapt, 2, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dmCur)));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">    PetscCall(DMGetLabel(dmCur, "adaptLabel", &adaptLabel));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">    PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Cell to refine = %d\n",nAdapt[0]));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">    PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Cell to keep = %d\n",nAdapt[1]));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">    <u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">    if (nAdapt[0]) {<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">      PetscCall(DMAdaptLabel(dmCur, adaptLabel, &dmAdapt));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">      PetscCall(DMHasLabel(dmAdapt, "adaptLabel", &hasLabel));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">      PetscCall(DMDestroy(&dmCur));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">      PetscCall(DMViewFromOptions(dmAdapt, NULL, "-adapt_dm_view"));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">      dmCur = dmAdapt;<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">    }<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">    //PetscCall(DMLabelDestroy(&adaptLabel));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">    PetscCall(DMDestroy(&dmConv));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">    if (adaptIter==maxAdaptIter) adapt=PETSC_FALSE;<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  }<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  *dm = dmCur;<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  PetscFunctionReturn(PETSC_SUCCESS);<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm"><span lang="DE">}<u></u><u></u></span></p>
<p class="MsoNormal" style="border:none;padding:0cm"><span lang="DE"><u></u> <u></u></span></p>
<p class="MsoNormal" style="border:none;padding:0cm"><span lang="DE">int main(int argc, char **argv)<u></u><u></u></span></p>
<p class="MsoNormal" style="border:none;padding:0cm"><span lang="DE">{<u></u><u></u></span></p>
<p class="MsoNormal" style="border:none;padding:0cm"><span lang="DE">  DM          dm;<u></u><u></u></span></p>
<p class="MsoNormal" style="border:none;padding:0cm"><span lang="DE">  </span>char        typeString[256] = {'\0'};<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  PetscViewer viewer          = NULL;<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  <u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  PetscFunctionBeginUser;<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  PetscCall(PetscInitialize(&argc, &argv, NULL, help));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  PetscCall(PetscStrncpy(typeString, DMFOREST, 256));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "DM Forest example options", NULL);<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  PetscCall(PetscOptionsString("-dm_type", "The type of the dm", NULL, DMFOREST, typeString, sizeof(typeString), NULL));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  PetscOptionsEnd();<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  <u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  PetscCall(PetscPrintf(PETSC_COMM_SELF,"\n ==== TOY CODE DMFOREST WITH AMR ====\n"));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  <u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  PetscCall(DMSetType(dm, (DMType)typeString));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  PetscCall(DMSetFromOptions(dm));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  PetscCall(DMSetUp(dm));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">    <u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  /* Adapt */<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  PetscCall(PetscPrintf(PETSC_COMM_SELF,"\nADAPTIVITY PHASE STARTED\n"));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  PetscCall(AdaptMesh(&dm));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  PetscCall(PetscPrintf(PETSC_COMM_SELF,"\nADAPTIVITY PHASE ENDED\n\n"));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  <u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  PetscCall(PetscViewerDestroy(&viewer));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm"><u></u> <u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  PetscCall(DMDestroy(&dm));<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  PetscCall(PetscFinalize());<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">  return 0;<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm">}<u></u><u></u></p>
<p class="MsoNormal" style="border:none;padding:0cm"><u></u> <u></u></p>
</div>
<p class="MsoNormal"><u></u> <u></u></p>
<p class="MsoNormal">Thanks a lot,<u></u><u></u></p>
<p class="MsoNormal"><u></u> <u></u></p>
<p class="MsoNormal">Silvia<u></u><u></u></p>
<p class="MsoNormal"><u></u> <u></u></p>
</div>
</div>

</div></blockquote></div><div><br clear="all"></div><div><br></div><span class="gmail_signature_prefix">-- </span><br><div dir="ltr" class="gmail_signature"><div dir="ltr"><div><div dir="ltr"><div><div dir="ltr"><div>What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br>-- Norbert Wiener</div><div><br></div><div><a href="https://urldefense.us/v3/__http://www.cse.buffalo.edu/*knepley/__;fg!!G_uCfscf7eWS!ZSCq5o18nPkr_fesXZsTE1dDngSGlwHv1_gJXKWrkgdaUqTlCeNxjF-VSyRK8bjgaHXIu9LJw-fH5Wo8ZisR$" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br></div></div></div></div></div></div></div></div>