1*c6f61ee2SBarry Smith #include <petscdmadaptor.h> /*I "petscdmadaptor.h" I*/ 2*c6f61ee2SBarry Smith #include <petscdmplex.h> 3*c6f61ee2SBarry Smith #include <petscdmforest.h> 4*c6f61ee2SBarry Smith #include <petscds.h> 5*c6f61ee2SBarry Smith #include <petscblaslapack.h> 6*c6f61ee2SBarry Smith 7*c6f61ee2SBarry Smith #include <petsc/private/dmadaptorimpl.h> 8*c6f61ee2SBarry Smith #include <petsc/private/dmpleximpl.h> 9*c6f61ee2SBarry Smith #include <petsc/private/petscfeimpl.h> 10*c6f61ee2SBarry Smith 11*c6f61ee2SBarry Smith static PetscErrorCode DMAdaptorSimpleErrorIndicator_Private(DMAdaptor, PetscInt, PetscInt, const PetscScalar *, const PetscScalar *, const PetscFVCellGeom *, PetscReal *, void *); 12*c6f61ee2SBarry Smith 13*c6f61ee2SBarry Smith static PetscErrorCode DMAdaptorTransferSolution_Exact_Private(DMAdaptor adaptor, DM dm, Vec u, DM adm, Vec au, void *ctx) 14*c6f61ee2SBarry Smith { 15*c6f61ee2SBarry Smith PetscFunctionBeginUser; 16*c6f61ee2SBarry Smith PetscCall(DMProjectFunction(adm, 0.0, adaptor->exactSol, adaptor->exactCtx, INSERT_ALL_VALUES, au)); 17*c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 18*c6f61ee2SBarry Smith } 19*c6f61ee2SBarry Smith 20*c6f61ee2SBarry Smith /*@ 21*c6f61ee2SBarry Smith DMAdaptorCreate - Create a `DMAdaptor` object. Its purpose is to construct a adaptation `DMLabel` or metric `Vec` that can be used to modify the `DM`. 22*c6f61ee2SBarry Smith 23*c6f61ee2SBarry Smith Collective 24*c6f61ee2SBarry Smith 25*c6f61ee2SBarry Smith Input Parameter: 26*c6f61ee2SBarry Smith . comm - The communicator for the `DMAdaptor` object 27*c6f61ee2SBarry Smith 28*c6f61ee2SBarry Smith Output Parameter: 29*c6f61ee2SBarry Smith . adaptor - The `DMAdaptor` object 30*c6f61ee2SBarry Smith 31*c6f61ee2SBarry Smith Level: beginner 32*c6f61ee2SBarry Smith 33*c6f61ee2SBarry Smith .seealso: `DM`, `DMAdaptor`, `DMAdaptorDestroy()`, `DMAdaptorAdapt()`, `PetscConvEst`, `PetscConvEstCreate()` 34*c6f61ee2SBarry Smith @*/ 35*c6f61ee2SBarry Smith PetscErrorCode DMAdaptorCreate(MPI_Comm comm, DMAdaptor *adaptor) 36*c6f61ee2SBarry Smith { 37*c6f61ee2SBarry Smith VecTaggerBox refineBox, coarsenBox; 38*c6f61ee2SBarry Smith 39*c6f61ee2SBarry Smith PetscFunctionBegin; 40*c6f61ee2SBarry Smith PetscAssertPointer(adaptor, 2); 41*c6f61ee2SBarry Smith PetscCall(PetscSysInitializePackage()); 42*c6f61ee2SBarry Smith PetscCall(PetscHeaderCreate(*adaptor, DM_CLASSID, "DMAdaptor", "DM Adaptor", "SNES", comm, DMAdaptorDestroy, DMAdaptorView)); 43*c6f61ee2SBarry Smith 44*c6f61ee2SBarry Smith (*adaptor)->monitor = PETSC_FALSE; 45*c6f61ee2SBarry Smith (*adaptor)->adaptCriterion = DM_ADAPTATION_NONE; 46*c6f61ee2SBarry Smith (*adaptor)->numSeq = 1; 47*c6f61ee2SBarry Smith (*adaptor)->Nadapt = -1; 48*c6f61ee2SBarry Smith (*adaptor)->refinementFactor = 2.0; 49*c6f61ee2SBarry Smith (*adaptor)->ops->computeerrorindicator = DMAdaptorSimpleErrorIndicator_Private; 50*c6f61ee2SBarry Smith refineBox.min = refineBox.max = PETSC_MAX_REAL; 51*c6f61ee2SBarry Smith PetscCall(VecTaggerCreate(PetscObjectComm((PetscObject)*adaptor), &(*adaptor)->refineTag)); 52*c6f61ee2SBarry Smith PetscCall(PetscObjectSetOptionsPrefix((PetscObject)(*adaptor)->refineTag, "refine_")); 53*c6f61ee2SBarry Smith PetscCall(VecTaggerSetType((*adaptor)->refineTag, VECTAGGERABSOLUTE)); 54*c6f61ee2SBarry Smith PetscCall(VecTaggerAbsoluteSetBox((*adaptor)->refineTag, &refineBox)); 55*c6f61ee2SBarry Smith coarsenBox.min = coarsenBox.max = PETSC_MAX_REAL; 56*c6f61ee2SBarry Smith PetscCall(VecTaggerCreate(PetscObjectComm((PetscObject)*adaptor), &(*adaptor)->coarsenTag)); 57*c6f61ee2SBarry Smith PetscCall(PetscObjectSetOptionsPrefix((PetscObject)(*adaptor)->coarsenTag, "coarsen_")); 58*c6f61ee2SBarry Smith PetscCall(VecTaggerSetType((*adaptor)->coarsenTag, VECTAGGERABSOLUTE)); 59*c6f61ee2SBarry Smith PetscCall(VecTaggerAbsoluteSetBox((*adaptor)->coarsenTag, &coarsenBox)); 60*c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 61*c6f61ee2SBarry Smith } 62*c6f61ee2SBarry Smith 63*c6f61ee2SBarry Smith /*@ 64*c6f61ee2SBarry Smith DMAdaptorDestroy - Destroys a `DMAdaptor` object 65*c6f61ee2SBarry Smith 66*c6f61ee2SBarry Smith Collective 67*c6f61ee2SBarry Smith 68*c6f61ee2SBarry Smith Input Parameter: 69*c6f61ee2SBarry Smith . adaptor - The `DMAdaptor` object 70*c6f61ee2SBarry Smith 71*c6f61ee2SBarry Smith Level: beginner 72*c6f61ee2SBarry Smith 73*c6f61ee2SBarry Smith .seealso: `DM`, `DMAdaptor`, `DMAdaptorCreate()`, `DMAdaptorAdapt()` 74*c6f61ee2SBarry Smith @*/ 75*c6f61ee2SBarry Smith PetscErrorCode DMAdaptorDestroy(DMAdaptor *adaptor) 76*c6f61ee2SBarry Smith { 77*c6f61ee2SBarry Smith PetscFunctionBegin; 78*c6f61ee2SBarry Smith if (!*adaptor) PetscFunctionReturn(PETSC_SUCCESS); 79*c6f61ee2SBarry Smith PetscValidHeaderSpecific((*adaptor), DM_CLASSID, 1); 80*c6f61ee2SBarry Smith if (--((PetscObject)(*adaptor))->refct > 0) { 81*c6f61ee2SBarry Smith *adaptor = NULL; 82*c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 83*c6f61ee2SBarry Smith } 84*c6f61ee2SBarry Smith PetscCall(VecTaggerDestroy(&(*adaptor)->refineTag)); 85*c6f61ee2SBarry Smith PetscCall(VecTaggerDestroy(&(*adaptor)->coarsenTag)); 86*c6f61ee2SBarry Smith PetscCall(PetscFree2((*adaptor)->exactSol, (*adaptor)->exactCtx)); 87*c6f61ee2SBarry Smith PetscCall(PetscHeaderDestroy(adaptor)); 88*c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 89*c6f61ee2SBarry Smith } 90*c6f61ee2SBarry Smith 91*c6f61ee2SBarry Smith /*@ 92*c6f61ee2SBarry Smith DMAdaptorSetFromOptions - Sets properties of a `DMAdaptor` object from values in the options database 93*c6f61ee2SBarry Smith 94*c6f61ee2SBarry Smith Collective 95*c6f61ee2SBarry Smith 96*c6f61ee2SBarry Smith Input Parameter: 97*c6f61ee2SBarry Smith . adaptor - The `DMAdaptor` object 98*c6f61ee2SBarry Smith 99*c6f61ee2SBarry Smith Options Database Keys: 100*c6f61ee2SBarry Smith + -adaptor_monitor <bool> - Monitor the adaptation process 101*c6f61ee2SBarry Smith . -adaptor_sequence_num <num> - Number of adaptations to generate an optimal grid 102*c6f61ee2SBarry Smith . -adaptor_target_num <num> - Set the target number of vertices N_adapt, -1 for automatic determination 103*c6f61ee2SBarry Smith - -adaptor_refinement_factor <r> - Set r such that N_adapt = r^dim N_orig 104*c6f61ee2SBarry Smith 105*c6f61ee2SBarry Smith Level: beginner 106*c6f61ee2SBarry Smith 107*c6f61ee2SBarry Smith .seealso: `DM`, `DMAdaptor`, `DMAdaptorCreate()`, `DMAdaptorAdapt()` 108*c6f61ee2SBarry Smith @*/ 109*c6f61ee2SBarry Smith PetscErrorCode DMAdaptorSetFromOptions(DMAdaptor adaptor) 110*c6f61ee2SBarry Smith { 111*c6f61ee2SBarry Smith PetscFunctionBegin; 112*c6f61ee2SBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)adaptor), "", "DM Adaptor Options", "DMAdaptor"); 113*c6f61ee2SBarry Smith PetscCall(PetscOptionsBool("-adaptor_monitor", "Monitor the adaptation process", "DMAdaptorMonitor", adaptor->monitor, &adaptor->monitor, NULL)); 114*c6f61ee2SBarry Smith PetscCall(PetscOptionsInt("-adaptor_sequence_num", "Number of adaptations to generate an optimal grid", "DMAdaptorSetSequenceLength", adaptor->numSeq, &adaptor->numSeq, NULL)); 115*c6f61ee2SBarry Smith PetscCall(PetscOptionsInt("-adaptor_target_num", "Set the target number of vertices N_adapt, -1 for automatic determination", "DMAdaptor", adaptor->Nadapt, &adaptor->Nadapt, NULL)); 116*c6f61ee2SBarry Smith PetscCall(PetscOptionsReal("-adaptor_refinement_factor", "Set r such that N_adapt = r^dim N_orig", "DMAdaptor", adaptor->refinementFactor, &adaptor->refinementFactor, NULL)); 117*c6f61ee2SBarry Smith PetscOptionsEnd(); 118*c6f61ee2SBarry Smith PetscCall(VecTaggerSetFromOptions(adaptor->refineTag)); 119*c6f61ee2SBarry Smith PetscCall(VecTaggerSetFromOptions(adaptor->coarsenTag)); 120*c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 121*c6f61ee2SBarry Smith } 122*c6f61ee2SBarry Smith 123*c6f61ee2SBarry Smith /*@ 124*c6f61ee2SBarry Smith DMAdaptorView - Views a `DMAdaptor` object 125*c6f61ee2SBarry Smith 126*c6f61ee2SBarry Smith Collective 127*c6f61ee2SBarry Smith 128*c6f61ee2SBarry Smith Input Parameters: 129*c6f61ee2SBarry Smith + adaptor - The `DMAdaptor` object 130*c6f61ee2SBarry Smith - viewer - The `PetscViewer` object 131*c6f61ee2SBarry Smith 132*c6f61ee2SBarry Smith Level: beginner 133*c6f61ee2SBarry Smith 134*c6f61ee2SBarry Smith .seealso: `DM`, `DMAdaptor`, `DMAdaptorCreate()`, `DMAdaptorAdapt()` 135*c6f61ee2SBarry Smith @*/ 136*c6f61ee2SBarry Smith PetscErrorCode DMAdaptorView(DMAdaptor adaptor, PetscViewer viewer) 137*c6f61ee2SBarry Smith { 138*c6f61ee2SBarry Smith PetscFunctionBegin; 139*c6f61ee2SBarry Smith PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)adaptor, viewer)); 140*c6f61ee2SBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, "DM Adaptor\n")); 141*c6f61ee2SBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, " sequence length: %" PetscInt_FMT "\n", adaptor->numSeq)); 142*c6f61ee2SBarry Smith PetscCall(VecTaggerView(adaptor->refineTag, viewer)); 143*c6f61ee2SBarry Smith PetscCall(VecTaggerView(adaptor->coarsenTag, viewer)); 144*c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 145*c6f61ee2SBarry Smith } 146*c6f61ee2SBarry Smith 147*c6f61ee2SBarry Smith /*@ 148*c6f61ee2SBarry Smith DMAdaptorGetSolver - Gets the solver used to produce discrete solutions 149*c6f61ee2SBarry Smith 150*c6f61ee2SBarry Smith Not Collective 151*c6f61ee2SBarry Smith 152*c6f61ee2SBarry Smith Input Parameter: 153*c6f61ee2SBarry Smith . adaptor - The `DMAdaptor` object 154*c6f61ee2SBarry Smith 155*c6f61ee2SBarry Smith Output Parameter: 156*c6f61ee2SBarry Smith . snes - The solver 157*c6f61ee2SBarry Smith 158*c6f61ee2SBarry Smith Level: intermediate 159*c6f61ee2SBarry Smith 160*c6f61ee2SBarry Smith .seealso: `DM`, `DMAdaptor`, `DMAdaptorSetSolver()`, `DMAdaptorCreate()`, `DMAdaptorAdapt()` 161*c6f61ee2SBarry Smith @*/ 162*c6f61ee2SBarry Smith PetscErrorCode DMAdaptorGetSolver(DMAdaptor adaptor, SNES *snes) 163*c6f61ee2SBarry Smith { 164*c6f61ee2SBarry Smith PetscFunctionBegin; 165*c6f61ee2SBarry Smith PetscValidHeaderSpecific(adaptor, DM_CLASSID, 1); 166*c6f61ee2SBarry Smith PetscAssertPointer(snes, 2); 167*c6f61ee2SBarry Smith *snes = adaptor->snes; 168*c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 169*c6f61ee2SBarry Smith } 170*c6f61ee2SBarry Smith 171*c6f61ee2SBarry Smith /*@ 172*c6f61ee2SBarry Smith DMAdaptorSetSolver - Sets the solver used to produce discrete solutions 173*c6f61ee2SBarry Smith 174*c6f61ee2SBarry Smith Not Collective 175*c6f61ee2SBarry Smith 176*c6f61ee2SBarry Smith Input Parameters: 177*c6f61ee2SBarry Smith + adaptor - The `DMAdaptor` object 178*c6f61ee2SBarry Smith - snes - The solver, this MUST have an attached `DM`/`PetscDS`, so that the exact solution can be computed 179*c6f61ee2SBarry Smith 180*c6f61ee2SBarry Smith Level: intermediate 181*c6f61ee2SBarry Smith 182*c6f61ee2SBarry Smith .seealso: `DMAdaptor`, `DMAdaptorGetSolver()`, `DMAdaptorCreate()`, `DMAdaptorAdapt()` 183*c6f61ee2SBarry Smith @*/ 184*c6f61ee2SBarry Smith PetscErrorCode DMAdaptorSetSolver(DMAdaptor adaptor, SNES snes) 185*c6f61ee2SBarry Smith { 186*c6f61ee2SBarry Smith PetscFunctionBegin; 187*c6f61ee2SBarry Smith PetscValidHeaderSpecific(adaptor, DM_CLASSID, 1); 188*c6f61ee2SBarry Smith PetscValidHeaderSpecific(snes, SNES_CLASSID, 2); 189*c6f61ee2SBarry Smith adaptor->snes = snes; 190*c6f61ee2SBarry Smith PetscCall(SNESGetDM(adaptor->snes, &adaptor->idm)); 191*c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 192*c6f61ee2SBarry Smith } 193*c6f61ee2SBarry Smith 194*c6f61ee2SBarry Smith /*@ 195*c6f61ee2SBarry Smith DMAdaptorGetSequenceLength - Gets the number of sequential adaptations used by an adapter 196*c6f61ee2SBarry Smith 197*c6f61ee2SBarry Smith Not Collective 198*c6f61ee2SBarry Smith 199*c6f61ee2SBarry Smith Input Parameter: 200*c6f61ee2SBarry Smith . adaptor - The `DMAdaptor` object 201*c6f61ee2SBarry Smith 202*c6f61ee2SBarry Smith Output Parameter: 203*c6f61ee2SBarry Smith . num - The number of adaptations 204*c6f61ee2SBarry Smith 205*c6f61ee2SBarry Smith Level: intermediate 206*c6f61ee2SBarry Smith 207*c6f61ee2SBarry Smith .seealso: `DMAdaptor`, `DMAdaptorSetSequenceLength()`, `DMAdaptorCreate()`, `DMAdaptorAdapt()` 208*c6f61ee2SBarry Smith @*/ 209*c6f61ee2SBarry Smith PetscErrorCode DMAdaptorGetSequenceLength(DMAdaptor adaptor, PetscInt *num) 210*c6f61ee2SBarry Smith { 211*c6f61ee2SBarry Smith PetscFunctionBegin; 212*c6f61ee2SBarry Smith PetscValidHeaderSpecific(adaptor, DM_CLASSID, 1); 213*c6f61ee2SBarry Smith PetscAssertPointer(num, 2); 214*c6f61ee2SBarry Smith *num = adaptor->numSeq; 215*c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 216*c6f61ee2SBarry Smith } 217*c6f61ee2SBarry Smith 218*c6f61ee2SBarry Smith /*@ 219*c6f61ee2SBarry Smith DMAdaptorSetSequenceLength - Sets the number of sequential adaptations 220*c6f61ee2SBarry Smith 221*c6f61ee2SBarry Smith Not Collective 222*c6f61ee2SBarry Smith 223*c6f61ee2SBarry Smith Input Parameters: 224*c6f61ee2SBarry Smith + adaptor - The `DMAdaptor` object 225*c6f61ee2SBarry Smith - num - The number of adaptations 226*c6f61ee2SBarry Smith 227*c6f61ee2SBarry Smith Level: intermediate 228*c6f61ee2SBarry Smith 229*c6f61ee2SBarry Smith .seealso: `DMAdaptorGetSequenceLength()`, `DMAdaptorCreate()`, `DMAdaptorAdapt()` 230*c6f61ee2SBarry Smith @*/ 231*c6f61ee2SBarry Smith PetscErrorCode DMAdaptorSetSequenceLength(DMAdaptor adaptor, PetscInt num) 232*c6f61ee2SBarry Smith { 233*c6f61ee2SBarry Smith PetscFunctionBegin; 234*c6f61ee2SBarry Smith PetscValidHeaderSpecific(adaptor, DM_CLASSID, 1); 235*c6f61ee2SBarry Smith adaptor->numSeq = num; 236*c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 237*c6f61ee2SBarry Smith } 238*c6f61ee2SBarry Smith 239*c6f61ee2SBarry Smith /*@ 240*c6f61ee2SBarry Smith DMAdaptorSetUp - After the solver is specified, creates data structures for controlling adaptivity 241*c6f61ee2SBarry Smith 242*c6f61ee2SBarry Smith Collective 243*c6f61ee2SBarry Smith 244*c6f61ee2SBarry Smith Input Parameter: 245*c6f61ee2SBarry Smith . adaptor - The `DMAdaptor` object 246*c6f61ee2SBarry Smith 247*c6f61ee2SBarry Smith Level: beginner 248*c6f61ee2SBarry Smith 249*c6f61ee2SBarry Smith .seealso: `DMAdaptor`, `DMAdaptorCreate()`, `DMAdaptorAdapt()` 250*c6f61ee2SBarry Smith @*/ 251*c6f61ee2SBarry Smith PetscErrorCode DMAdaptorSetUp(DMAdaptor adaptor) 252*c6f61ee2SBarry Smith { 253*c6f61ee2SBarry Smith PetscDS prob; 254*c6f61ee2SBarry Smith PetscInt Nf, f; 255*c6f61ee2SBarry Smith 256*c6f61ee2SBarry Smith PetscFunctionBegin; 257*c6f61ee2SBarry Smith PetscCall(DMGetDS(adaptor->idm, &prob)); 258*c6f61ee2SBarry Smith PetscCall(VecTaggerSetUp(adaptor->refineTag)); 259*c6f61ee2SBarry Smith PetscCall(VecTaggerSetUp(adaptor->coarsenTag)); 260*c6f61ee2SBarry Smith PetscCall(PetscDSGetNumFields(prob, &Nf)); 261*c6f61ee2SBarry Smith PetscCall(PetscMalloc2(Nf, &adaptor->exactSol, Nf, &adaptor->exactCtx)); 262*c6f61ee2SBarry Smith for (f = 0; f < Nf; ++f) { 263*c6f61ee2SBarry Smith PetscCall(PetscDSGetExactSolution(prob, f, &adaptor->exactSol[f], &adaptor->exactCtx[f])); 264*c6f61ee2SBarry Smith /* TODO Have a flag that forces projection rather than using the exact solution */ 265*c6f61ee2SBarry Smith if (adaptor->exactSol[0]) PetscCall(DMAdaptorSetTransferFunction(adaptor, DMAdaptorTransferSolution_Exact_Private)); 266*c6f61ee2SBarry Smith } 267*c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 268*c6f61ee2SBarry Smith } 269*c6f61ee2SBarry Smith 270*c6f61ee2SBarry Smith PetscErrorCode DMAdaptorGetTransferFunction(DMAdaptor adaptor, PetscErrorCode (**tfunc)(DMAdaptor, DM, Vec, DM, Vec, void *)) 271*c6f61ee2SBarry Smith { 272*c6f61ee2SBarry Smith PetscFunctionBegin; 273*c6f61ee2SBarry Smith *tfunc = adaptor->ops->transfersolution; 274*c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 275*c6f61ee2SBarry Smith } 276*c6f61ee2SBarry Smith 277*c6f61ee2SBarry Smith PetscErrorCode DMAdaptorSetTransferFunction(DMAdaptor adaptor, PetscErrorCode (*tfunc)(DMAdaptor, DM, Vec, DM, Vec, void *)) 278*c6f61ee2SBarry Smith { 279*c6f61ee2SBarry Smith PetscFunctionBegin; 280*c6f61ee2SBarry Smith adaptor->ops->transfersolution = tfunc; 281*c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 282*c6f61ee2SBarry Smith } 283*c6f61ee2SBarry Smith 284*c6f61ee2SBarry Smith static PetscErrorCode DMAdaptorPreAdapt(DMAdaptor adaptor, Vec locX) 285*c6f61ee2SBarry Smith { 286*c6f61ee2SBarry Smith DM plex; 287*c6f61ee2SBarry Smith PetscDS prob; 288*c6f61ee2SBarry Smith PetscObject obj; 289*c6f61ee2SBarry Smith PetscClassId id; 290*c6f61ee2SBarry Smith PetscBool isForest; 291*c6f61ee2SBarry Smith 292*c6f61ee2SBarry Smith PetscFunctionBegin; 293*c6f61ee2SBarry Smith PetscCall(DMConvert(adaptor->idm, DMPLEX, &plex)); 294*c6f61ee2SBarry Smith PetscCall(DMGetDS(adaptor->idm, &prob)); 295*c6f61ee2SBarry Smith PetscCall(PetscDSGetDiscretization(prob, 0, &obj)); 296*c6f61ee2SBarry Smith PetscCall(PetscObjectGetClassId(obj, &id)); 297*c6f61ee2SBarry Smith PetscCall(DMIsForest(adaptor->idm, &isForest)); 298*c6f61ee2SBarry Smith if (adaptor->adaptCriterion == DM_ADAPTATION_NONE) { 299*c6f61ee2SBarry Smith if (isForest) { 300*c6f61ee2SBarry Smith adaptor->adaptCriterion = DM_ADAPTATION_LABEL; 301*c6f61ee2SBarry Smith } 302*c6f61ee2SBarry Smith #if defined(PETSC_HAVE_PRAGMATIC) 303*c6f61ee2SBarry Smith else { 304*c6f61ee2SBarry Smith adaptor->adaptCriterion = DM_ADAPTATION_METRIC; 305*c6f61ee2SBarry Smith } 306*c6f61ee2SBarry Smith #elif defined(PETSC_HAVE_MMG) 307*c6f61ee2SBarry Smith else { 308*c6f61ee2SBarry Smith adaptor->adaptCriterion = DM_ADAPTATION_METRIC; 309*c6f61ee2SBarry Smith } 310*c6f61ee2SBarry Smith #elif defined(PETSC_HAVE_PARMMG) 311*c6f61ee2SBarry Smith else { 312*c6f61ee2SBarry Smith adaptor->adaptCriterion = DM_ADAPTATION_METRIC; 313*c6f61ee2SBarry Smith } 314*c6f61ee2SBarry Smith #else 315*c6f61ee2SBarry Smith else { 316*c6f61ee2SBarry Smith adaptor->adaptCriterion = DM_ADAPTATION_REFINE; 317*c6f61ee2SBarry Smith } 318*c6f61ee2SBarry Smith #endif 319*c6f61ee2SBarry Smith } 320*c6f61ee2SBarry Smith if (id == PETSCFV_CLASSID) { 321*c6f61ee2SBarry Smith adaptor->femType = PETSC_FALSE; 322*c6f61ee2SBarry Smith } else { 323*c6f61ee2SBarry Smith adaptor->femType = PETSC_TRUE; 324*c6f61ee2SBarry Smith } 325*c6f61ee2SBarry Smith if (adaptor->femType) { 326*c6f61ee2SBarry Smith /* Compute local solution bc */ 327*c6f61ee2SBarry Smith PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, 0.0, adaptor->faceGeom, adaptor->cellGeom, NULL)); 328*c6f61ee2SBarry Smith } else { 329*c6f61ee2SBarry Smith PetscFV fvm = (PetscFV)obj; 330*c6f61ee2SBarry Smith PetscLimiter noneLimiter; 331*c6f61ee2SBarry Smith Vec grad; 332*c6f61ee2SBarry Smith 333*c6f61ee2SBarry Smith PetscCall(PetscFVGetComputeGradients(fvm, &adaptor->computeGradient)); 334*c6f61ee2SBarry Smith PetscCall(PetscFVSetComputeGradients(fvm, PETSC_TRUE)); 335*c6f61ee2SBarry Smith /* Use no limiting when reconstructing gradients for adaptivity */ 336*c6f61ee2SBarry Smith PetscCall(PetscFVGetLimiter(fvm, &adaptor->limiter)); 337*c6f61ee2SBarry Smith PetscCall(PetscObjectReference((PetscObject)adaptor->limiter)); 338*c6f61ee2SBarry Smith PetscCall(PetscLimiterCreate(PetscObjectComm((PetscObject)fvm), &noneLimiter)); 339*c6f61ee2SBarry Smith PetscCall(PetscLimiterSetType(noneLimiter, PETSCLIMITERNONE)); 340*c6f61ee2SBarry Smith PetscCall(PetscFVSetLimiter(fvm, noneLimiter)); 341*c6f61ee2SBarry Smith /* Get FVM data */ 342*c6f61ee2SBarry Smith PetscCall(DMPlexGetDataFVM(plex, fvm, &adaptor->cellGeom, &adaptor->faceGeom, &adaptor->gradDM)); 343*c6f61ee2SBarry Smith PetscCall(VecGetDM(adaptor->cellGeom, &adaptor->cellDM)); 344*c6f61ee2SBarry Smith PetscCall(VecGetArrayRead(adaptor->cellGeom, &adaptor->cellGeomArray)); 345*c6f61ee2SBarry Smith /* Compute local solution bc */ 346*c6f61ee2SBarry Smith PetscCall(DMPlexInsertBoundaryValues(plex, PETSC_TRUE, locX, 0.0, adaptor->faceGeom, adaptor->cellGeom, NULL)); 347*c6f61ee2SBarry Smith /* Compute gradients */ 348*c6f61ee2SBarry Smith PetscCall(DMCreateGlobalVector(adaptor->gradDM, &grad)); 349*c6f61ee2SBarry Smith PetscCall(DMPlexReconstructGradientsFVM(plex, locX, grad)); 350*c6f61ee2SBarry Smith PetscCall(DMGetLocalVector(adaptor->gradDM, &adaptor->cellGrad)); 351*c6f61ee2SBarry Smith PetscCall(DMGlobalToLocalBegin(adaptor->gradDM, grad, INSERT_VALUES, adaptor->cellGrad)); 352*c6f61ee2SBarry Smith PetscCall(DMGlobalToLocalEnd(adaptor->gradDM, grad, INSERT_VALUES, adaptor->cellGrad)); 353*c6f61ee2SBarry Smith PetscCall(VecDestroy(&grad)); 354*c6f61ee2SBarry Smith PetscCall(VecGetArrayRead(adaptor->cellGrad, &adaptor->cellGradArray)); 355*c6f61ee2SBarry Smith } 356*c6f61ee2SBarry Smith PetscCall(DMDestroy(&plex)); 357*c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 358*c6f61ee2SBarry Smith } 359*c6f61ee2SBarry Smith 360*c6f61ee2SBarry Smith static PetscErrorCode DMAdaptorTransferSolution(DMAdaptor adaptor, DM dm, Vec x, DM adm, Vec ax) 361*c6f61ee2SBarry Smith { 362*c6f61ee2SBarry Smith PetscReal time = 0.0; 363*c6f61ee2SBarry Smith Mat interp; 364*c6f61ee2SBarry Smith void *ctx; 365*c6f61ee2SBarry Smith 366*c6f61ee2SBarry Smith PetscFunctionBegin; 367*c6f61ee2SBarry Smith PetscCall(DMGetApplicationContext(dm, &ctx)); 368*c6f61ee2SBarry Smith if (adaptor->ops->transfersolution) PetscUseTypeMethod(adaptor, transfersolution, dm, x, adm, ax, ctx); 369*c6f61ee2SBarry Smith else { 370*c6f61ee2SBarry Smith switch (adaptor->adaptCriterion) { 371*c6f61ee2SBarry Smith case DM_ADAPTATION_LABEL: 372*c6f61ee2SBarry Smith PetscCall(DMForestTransferVec(dm, x, adm, ax, PETSC_TRUE, time)); 373*c6f61ee2SBarry Smith break; 374*c6f61ee2SBarry Smith case DM_ADAPTATION_REFINE: 375*c6f61ee2SBarry Smith case DM_ADAPTATION_METRIC: 376*c6f61ee2SBarry Smith PetscCall(DMCreateInterpolation(dm, adm, &interp, NULL)); 377*c6f61ee2SBarry Smith PetscCall(MatInterpolate(interp, x, ax)); 378*c6f61ee2SBarry Smith PetscCall(DMInterpolate(dm, interp, adm)); 379*c6f61ee2SBarry Smith PetscCall(MatDestroy(&interp)); 380*c6f61ee2SBarry Smith break; 381*c6f61ee2SBarry Smith default: 382*c6f61ee2SBarry Smith SETERRQ(PetscObjectComm((PetscObject)adaptor), PETSC_ERR_SUP, "No built-in projection for this adaptation criterion: %d", adaptor->adaptCriterion); 383*c6f61ee2SBarry Smith } 384*c6f61ee2SBarry Smith } 385*c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 386*c6f61ee2SBarry Smith } 387*c6f61ee2SBarry Smith 388*c6f61ee2SBarry Smith static PetscErrorCode DMAdaptorPostAdapt(DMAdaptor adaptor) 389*c6f61ee2SBarry Smith { 390*c6f61ee2SBarry Smith PetscDS prob; 391*c6f61ee2SBarry Smith PetscObject obj; 392*c6f61ee2SBarry Smith PetscClassId id; 393*c6f61ee2SBarry Smith 394*c6f61ee2SBarry Smith PetscFunctionBegin; 395*c6f61ee2SBarry Smith PetscCall(DMGetDS(adaptor->idm, &prob)); 396*c6f61ee2SBarry Smith PetscCall(PetscDSGetDiscretization(prob, 0, &obj)); 397*c6f61ee2SBarry Smith PetscCall(PetscObjectGetClassId(obj, &id)); 398*c6f61ee2SBarry Smith if (id == PETSCFV_CLASSID) { 399*c6f61ee2SBarry Smith PetscFV fvm = (PetscFV)obj; 400*c6f61ee2SBarry Smith 401*c6f61ee2SBarry Smith PetscCall(PetscFVSetComputeGradients(fvm, adaptor->computeGradient)); 402*c6f61ee2SBarry Smith /* Restore original limiter */ 403*c6f61ee2SBarry Smith PetscCall(PetscFVSetLimiter(fvm, adaptor->limiter)); 404*c6f61ee2SBarry Smith 405*c6f61ee2SBarry Smith PetscCall(VecRestoreArrayRead(adaptor->cellGeom, &adaptor->cellGeomArray)); 406*c6f61ee2SBarry Smith PetscCall(VecRestoreArrayRead(adaptor->cellGrad, &adaptor->cellGradArray)); 407*c6f61ee2SBarry Smith PetscCall(DMRestoreLocalVector(adaptor->gradDM, &adaptor->cellGrad)); 408*c6f61ee2SBarry Smith } 409*c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 410*c6f61ee2SBarry Smith } 411*c6f61ee2SBarry Smith 412*c6f61ee2SBarry Smith /* 413*c6f61ee2SBarry Smith DMAdaptorSimpleErrorIndicator - Use the integrated gradient as an error indicator in the `DMAdaptor` 414*c6f61ee2SBarry Smith 415*c6f61ee2SBarry Smith Input Parameters: 416*c6f61ee2SBarry Smith + adaptor - The `DMAdaptor` object 417*c6f61ee2SBarry Smith . dim - The topological dimension 418*c6f61ee2SBarry Smith . cell - The cell 419*c6f61ee2SBarry Smith . field - The field integrated over the cell 420*c6f61ee2SBarry Smith . gradient - The gradient integrated over the cell 421*c6f61ee2SBarry Smith . cg - A `PetscFVCellGeom` struct 422*c6f61ee2SBarry Smith - ctx - A user context 423*c6f61ee2SBarry Smith 424*c6f61ee2SBarry Smith Output Parameter: 425*c6f61ee2SBarry Smith . errInd - The error indicator 426*c6f61ee2SBarry Smith 427*c6f61ee2SBarry Smith Developer Note: 428*c6f61ee2SBarry Smith Some of the input arguments are absurdly specialized to special situations, it is not clear this is a good general API 429*c6f61ee2SBarry Smith 430*c6f61ee2SBarry Smith .seealso: `DMAdaptor`, `DMAdaptorComputeErrorIndicator()` 431*c6f61ee2SBarry Smith */ 432*c6f61ee2SBarry Smith static PetscErrorCode DMAdaptorSimpleErrorIndicator_Private(DMAdaptor adaptor, PetscInt dim, PetscInt Nc, const PetscScalar *field, const PetscScalar *gradient, const PetscFVCellGeom *cg, PetscReal *errInd, void *ctx) 433*c6f61ee2SBarry Smith { 434*c6f61ee2SBarry Smith PetscReal err = 0.; 435*c6f61ee2SBarry Smith PetscInt c, d; 436*c6f61ee2SBarry Smith 437*c6f61ee2SBarry Smith PetscFunctionBeginHot; 438*c6f61ee2SBarry Smith for (c = 0; c < Nc; c++) { 439*c6f61ee2SBarry Smith for (d = 0; d < dim; ++d) err += PetscSqr(PetscRealPart(gradient[c * dim + d])); 440*c6f61ee2SBarry Smith } 441*c6f61ee2SBarry Smith *errInd = cg->volume * err; 442*c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 443*c6f61ee2SBarry Smith } 444*c6f61ee2SBarry Smith 445*c6f61ee2SBarry Smith static PetscErrorCode DMAdaptorComputeErrorIndicator_Private(DMAdaptor adaptor, DM plex, PetscInt cell, Vec locX, PetscReal *errInd) 446*c6f61ee2SBarry Smith { 447*c6f61ee2SBarry Smith PetscDS prob; 448*c6f61ee2SBarry Smith PetscObject obj; 449*c6f61ee2SBarry Smith PetscClassId id; 450*c6f61ee2SBarry Smith void *ctx; 451*c6f61ee2SBarry Smith PetscQuadrature quad; 452*c6f61ee2SBarry Smith PetscInt dim, d, cdim, Nc; 453*c6f61ee2SBarry Smith 454*c6f61ee2SBarry Smith PetscFunctionBegin; 455*c6f61ee2SBarry Smith *errInd = 0.; 456*c6f61ee2SBarry Smith PetscCall(DMGetDimension(plex, &dim)); 457*c6f61ee2SBarry Smith PetscCall(DMGetCoordinateDim(plex, &cdim)); 458*c6f61ee2SBarry Smith PetscCall(DMGetApplicationContext(plex, &ctx)); 459*c6f61ee2SBarry Smith PetscCall(DMGetDS(plex, &prob)); 460*c6f61ee2SBarry Smith PetscCall(PetscDSGetDiscretization(prob, 0, &obj)); 461*c6f61ee2SBarry Smith PetscCall(PetscObjectGetClassId(obj, &id)); 462*c6f61ee2SBarry Smith if (id == PETSCFV_CLASSID) { 463*c6f61ee2SBarry Smith const PetscScalar *pointSols; 464*c6f61ee2SBarry Smith const PetscScalar *pointSol; 465*c6f61ee2SBarry Smith const PetscScalar *pointGrad; 466*c6f61ee2SBarry Smith PetscFVCellGeom *cg; 467*c6f61ee2SBarry Smith 468*c6f61ee2SBarry Smith PetscCall(PetscFVGetNumComponents((PetscFV)obj, &Nc)); 469*c6f61ee2SBarry Smith PetscCall(VecGetArrayRead(locX, &pointSols)); 470*c6f61ee2SBarry Smith PetscCall(DMPlexPointLocalRead(plex, cell, pointSols, (void *)&pointSol)); 471*c6f61ee2SBarry Smith PetscCall(DMPlexPointLocalRead(adaptor->gradDM, cell, adaptor->cellGradArray, (void *)&pointGrad)); 472*c6f61ee2SBarry Smith PetscCall(DMPlexPointLocalRead(adaptor->cellDM, cell, adaptor->cellGeomArray, &cg)); 473*c6f61ee2SBarry Smith PetscUseTypeMethod(adaptor, computeerrorindicator, dim, Nc, pointSol, pointGrad, cg, errInd, ctx); 474*c6f61ee2SBarry Smith PetscCall(VecRestoreArrayRead(locX, &pointSols)); 475*c6f61ee2SBarry Smith } else { 476*c6f61ee2SBarry Smith PetscScalar *x = NULL, *field, *gradient, *interpolant, *interpolantGrad; 477*c6f61ee2SBarry Smith PetscFVCellGeom cg; 478*c6f61ee2SBarry Smith PetscFEGeom fegeom; 479*c6f61ee2SBarry Smith const PetscReal *quadWeights; 480*c6f61ee2SBarry Smith PetscReal *coords; 481*c6f61ee2SBarry Smith PetscInt Nb, fc, Nq, qNc, Nf, f, fieldOffset; 482*c6f61ee2SBarry Smith 483*c6f61ee2SBarry Smith fegeom.dim = dim; 484*c6f61ee2SBarry Smith fegeom.dimEmbed = cdim; 485*c6f61ee2SBarry Smith PetscCall(PetscDSGetNumFields(prob, &Nf)); 486*c6f61ee2SBarry Smith PetscCall(PetscFEGetQuadrature((PetscFE)obj, &quad)); 487*c6f61ee2SBarry Smith PetscCall(DMPlexVecGetClosure(plex, NULL, locX, cell, NULL, &x)); 488*c6f61ee2SBarry Smith PetscCall(PetscFEGetDimension((PetscFE)obj, &Nb)); 489*c6f61ee2SBarry Smith PetscCall(PetscFEGetNumComponents((PetscFE)obj, &Nc)); 490*c6f61ee2SBarry Smith PetscCall(PetscQuadratureGetData(quad, NULL, &qNc, &Nq, NULL, &quadWeights)); 491*c6f61ee2SBarry Smith PetscCall(PetscCalloc6(Nc, &field, cdim * Nc, &gradient, cdim * Nq, &coords, Nq, &fegeom.detJ, cdim * cdim * Nq, &fegeom.J, cdim * cdim * Nq, &fegeom.invJ)); 492*c6f61ee2SBarry Smith PetscCall(PetscMalloc2(Nc, &interpolant, cdim * Nc, &interpolantGrad)); 493*c6f61ee2SBarry Smith PetscCall(DMPlexComputeCellGeometryFEM(plex, cell, quad, coords, fegeom.J, fegeom.invJ, fegeom.detJ)); 494*c6f61ee2SBarry Smith PetscCall(DMPlexComputeCellGeometryFVM(plex, cell, &cg.volume, NULL, NULL)); 495*c6f61ee2SBarry Smith PetscCall(PetscArrayzero(gradient, cdim * Nc)); 496*c6f61ee2SBarry Smith for (f = 0, fieldOffset = 0; f < Nf; ++f) { 497*c6f61ee2SBarry Smith PetscInt qc = 0, q; 498*c6f61ee2SBarry Smith 499*c6f61ee2SBarry Smith PetscCall(PetscDSGetDiscretization(prob, f, &obj)); 500*c6f61ee2SBarry Smith PetscCall(PetscArrayzero(interpolant, Nc)); 501*c6f61ee2SBarry Smith PetscCall(PetscArrayzero(interpolantGrad, cdim * Nc)); 502*c6f61ee2SBarry Smith for (q = 0; q < Nq; ++q) { 503*c6f61ee2SBarry Smith PetscCall(PetscFEInterpolateFieldAndGradient_Static((PetscFE)obj, 1, x, &fegeom, q, interpolant, interpolantGrad)); 504*c6f61ee2SBarry Smith for (fc = 0; fc < Nc; ++fc) { 505*c6f61ee2SBarry Smith const PetscReal wt = quadWeights[q * qNc + qc + fc]; 506*c6f61ee2SBarry Smith 507*c6f61ee2SBarry Smith field[fc] += interpolant[fc] * wt * fegeom.detJ[q]; 508*c6f61ee2SBarry Smith for (d = 0; d < cdim; ++d) gradient[fc * cdim + d] += interpolantGrad[fc * dim + d] * wt * fegeom.detJ[q]; 509*c6f61ee2SBarry Smith } 510*c6f61ee2SBarry Smith } 511*c6f61ee2SBarry Smith fieldOffset += Nb; 512*c6f61ee2SBarry Smith qc += Nc; 513*c6f61ee2SBarry Smith } 514*c6f61ee2SBarry Smith PetscCall(PetscFree2(interpolant, interpolantGrad)); 515*c6f61ee2SBarry Smith PetscCall(DMPlexVecRestoreClosure(plex, NULL, locX, cell, NULL, &x)); 516*c6f61ee2SBarry Smith for (fc = 0; fc < Nc; ++fc) { 517*c6f61ee2SBarry Smith field[fc] /= cg.volume; 518*c6f61ee2SBarry Smith for (d = 0; d < cdim; ++d) gradient[fc * cdim + d] /= cg.volume; 519*c6f61ee2SBarry Smith } 520*c6f61ee2SBarry Smith PetscUseTypeMethod(adaptor, computeerrorindicator, dim, Nc, field, gradient, &cg, errInd, ctx); 521*c6f61ee2SBarry Smith PetscCall(PetscFree6(field, gradient, coords, fegeom.detJ, fegeom.J, fegeom.invJ)); 522*c6f61ee2SBarry Smith } 523*c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 524*c6f61ee2SBarry Smith } 525*c6f61ee2SBarry Smith 526*c6f61ee2SBarry Smith static void identityFunc(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f[]) 527*c6f61ee2SBarry Smith { 528*c6f61ee2SBarry Smith PetscInt i, j; 529*c6f61ee2SBarry Smith 530*c6f61ee2SBarry Smith for (i = 0; i < dim; ++i) { 531*c6f61ee2SBarry Smith for (j = 0; j < dim; ++j) f[i + dim * j] = u[i + dim * j]; 532*c6f61ee2SBarry Smith } 533*c6f61ee2SBarry Smith } 534*c6f61ee2SBarry Smith 535*c6f61ee2SBarry Smith static PetscErrorCode DMAdaptorAdapt_Sequence_Private(DMAdaptor adaptor, Vec inx, PetscBool doSolve, DM *adm, Vec *ax) 536*c6f61ee2SBarry Smith { 537*c6f61ee2SBarry Smith PetscDS prob; 538*c6f61ee2SBarry Smith void *ctx; 539*c6f61ee2SBarry Smith MPI_Comm comm; 540*c6f61ee2SBarry Smith PetscInt numAdapt = adaptor->numSeq, adaptIter; 541*c6f61ee2SBarry Smith PetscInt dim, coordDim, numFields, cStart, cEnd, c; 542*c6f61ee2SBarry Smith 543*c6f61ee2SBarry Smith PetscFunctionBegin; 544*c6f61ee2SBarry Smith PetscCall(DMViewFromOptions(adaptor->idm, NULL, "-dm_adapt_pre_view")); 545*c6f61ee2SBarry Smith PetscCall(VecViewFromOptions(inx, NULL, "-sol_adapt_pre_view")); 546*c6f61ee2SBarry Smith PetscCall(PetscObjectGetComm((PetscObject)adaptor, &comm)); 547*c6f61ee2SBarry Smith PetscCall(DMGetDimension(adaptor->idm, &dim)); 548*c6f61ee2SBarry Smith PetscCall(DMGetCoordinateDim(adaptor->idm, &coordDim)); 549*c6f61ee2SBarry Smith PetscCall(DMGetApplicationContext(adaptor->idm, &ctx)); 550*c6f61ee2SBarry Smith PetscCall(DMGetDS(adaptor->idm, &prob)); 551*c6f61ee2SBarry Smith PetscCall(PetscDSGetNumFields(prob, &numFields)); 552*c6f61ee2SBarry Smith PetscCheck(numFields != 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of fields is zero!"); 553*c6f61ee2SBarry Smith 554*c6f61ee2SBarry Smith /* Adapt until nothing changes */ 555*c6f61ee2SBarry Smith /* Adapt for a specified number of iterates */ 556*c6f61ee2SBarry Smith for (adaptIter = 0; adaptIter < numAdapt - 1; ++adaptIter) PetscCall(PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_(comm))); 557*c6f61ee2SBarry Smith for (adaptIter = 0; adaptIter < numAdapt; ++adaptIter) { 558*c6f61ee2SBarry Smith PetscBool adapted = PETSC_FALSE; 559*c6f61ee2SBarry Smith DM dm = adaptIter ? *adm : adaptor->idm, odm; 560*c6f61ee2SBarry Smith Vec x = adaptIter ? *ax : inx, locX, ox; 561*c6f61ee2SBarry Smith 562*c6f61ee2SBarry Smith PetscCall(DMGetLocalVector(dm, &locX)); 563*c6f61ee2SBarry Smith PetscCall(DMGlobalToLocalBegin(dm, adaptIter ? *ax : x, INSERT_VALUES, locX)); 564*c6f61ee2SBarry Smith PetscCall(DMGlobalToLocalEnd(dm, adaptIter ? *ax : x, INSERT_VALUES, locX)); 565*c6f61ee2SBarry Smith PetscCall(DMAdaptorPreAdapt(adaptor, locX)); 566*c6f61ee2SBarry Smith if (doSolve) { 567*c6f61ee2SBarry Smith SNES snes; 568*c6f61ee2SBarry Smith 569*c6f61ee2SBarry Smith PetscCall(DMAdaptorGetSolver(adaptor, &snes)); 570*c6f61ee2SBarry Smith PetscCall(SNESSolve(snes, NULL, adaptIter ? *ax : x)); 571*c6f61ee2SBarry Smith } 572*c6f61ee2SBarry Smith /* PetscCall(DMAdaptorMonitor(adaptor)); 573*c6f61ee2SBarry Smith Print iterate, memory used, DM, solution */ 574*c6f61ee2SBarry Smith switch (adaptor->adaptCriterion) { 575*c6f61ee2SBarry Smith case DM_ADAPTATION_REFINE: 576*c6f61ee2SBarry Smith PetscCall(DMRefine(dm, comm, &odm)); 577*c6f61ee2SBarry Smith PetscCheck(odm, comm, PETSC_ERR_ARG_INCOMP, "DMRefine() did not perform any refinement, cannot continue grid sequencing"); 578*c6f61ee2SBarry Smith adapted = PETSC_TRUE; 579*c6f61ee2SBarry Smith break; 580*c6f61ee2SBarry Smith case DM_ADAPTATION_LABEL: { 581*c6f61ee2SBarry Smith /* Adapt DM 582*c6f61ee2SBarry Smith Create local solution 583*c6f61ee2SBarry Smith Reconstruct gradients (FVM) or solve adjoint equation (FEM) 584*c6f61ee2SBarry Smith Produce cellwise error indicator */ 585*c6f61ee2SBarry Smith DM plex; 586*c6f61ee2SBarry Smith DMLabel adaptLabel; 587*c6f61ee2SBarry Smith IS refineIS, coarsenIS; 588*c6f61ee2SBarry Smith Vec errVec; 589*c6f61ee2SBarry Smith PetscScalar *errArray; 590*c6f61ee2SBarry Smith const PetscScalar *pointSols; 591*c6f61ee2SBarry Smith PetscReal minMaxInd[2] = {PETSC_MAX_REAL, PETSC_MIN_REAL}, minMaxIndGlobal[2]; 592*c6f61ee2SBarry Smith PetscInt nRefine, nCoarsen; 593*c6f61ee2SBarry Smith 594*c6f61ee2SBarry Smith PetscCall(DMConvert(dm, DMPLEX, &plex)); 595*c6f61ee2SBarry Smith PetscCall(DMLabelCreate(PETSC_COMM_SELF, "adapt", &adaptLabel)); 596*c6f61ee2SBarry Smith PetscCall(DMPlexGetSimplexOrBoxCells(plex, 0, &cStart, &cEnd)); 597*c6f61ee2SBarry Smith 598*c6f61ee2SBarry Smith PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)adaptor), cEnd - cStart, PETSC_DETERMINE, &errVec)); 599*c6f61ee2SBarry Smith PetscCall(VecSetUp(errVec)); 600*c6f61ee2SBarry Smith PetscCall(VecGetArray(errVec, &errArray)); 601*c6f61ee2SBarry Smith PetscCall(VecGetArrayRead(locX, &pointSols)); 602*c6f61ee2SBarry Smith for (c = cStart; c < cEnd; ++c) { 603*c6f61ee2SBarry Smith PetscReal errInd; 604*c6f61ee2SBarry Smith 605*c6f61ee2SBarry Smith PetscCall(DMAdaptorComputeErrorIndicator_Private(adaptor, plex, c, locX, &errInd)); 606*c6f61ee2SBarry Smith errArray[c - cStart] = errInd; 607*c6f61ee2SBarry Smith minMaxInd[0] = PetscMin(minMaxInd[0], errInd); 608*c6f61ee2SBarry Smith minMaxInd[1] = PetscMax(minMaxInd[1], errInd); 609*c6f61ee2SBarry Smith } 610*c6f61ee2SBarry Smith PetscCall(VecRestoreArrayRead(locX, &pointSols)); 611*c6f61ee2SBarry Smith PetscCall(VecRestoreArray(errVec, &errArray)); 612*c6f61ee2SBarry Smith PetscCall(PetscGlobalMinMaxReal(PetscObjectComm((PetscObject)adaptor), minMaxInd, minMaxIndGlobal)); 613*c6f61ee2SBarry Smith PetscCall(PetscInfo(adaptor, "DMAdaptor: error indicator range (%g, %g)\n", (double)minMaxIndGlobal[0], (double)minMaxIndGlobal[1])); 614*c6f61ee2SBarry Smith /* Compute IS from VecTagger */ 615*c6f61ee2SBarry Smith PetscCall(VecTaggerComputeIS(adaptor->refineTag, errVec, &refineIS, NULL)); 616*c6f61ee2SBarry Smith PetscCall(VecTaggerComputeIS(adaptor->coarsenTag, errVec, &coarsenIS, NULL)); 617*c6f61ee2SBarry Smith PetscCall(ISGetSize(refineIS, &nRefine)); 618*c6f61ee2SBarry Smith PetscCall(ISGetSize(coarsenIS, &nCoarsen)); 619*c6f61ee2SBarry Smith PetscCall(PetscInfo(adaptor, "DMAdaptor: numRefine %" PetscInt_FMT ", numCoarsen %" PetscInt_FMT "\n", nRefine, nCoarsen)); 620*c6f61ee2SBarry Smith if (nRefine) PetscCall(DMLabelSetStratumIS(adaptLabel, DM_ADAPT_REFINE, refineIS)); 621*c6f61ee2SBarry Smith if (nCoarsen) PetscCall(DMLabelSetStratumIS(adaptLabel, DM_ADAPT_COARSEN, coarsenIS)); 622*c6f61ee2SBarry Smith PetscCall(ISDestroy(&coarsenIS)); 623*c6f61ee2SBarry Smith PetscCall(ISDestroy(&refineIS)); 624*c6f61ee2SBarry Smith PetscCall(VecDestroy(&errVec)); 625*c6f61ee2SBarry Smith /* Adapt DM from label */ 626*c6f61ee2SBarry Smith if (nRefine || nCoarsen) { 627*c6f61ee2SBarry Smith PetscCall(DMAdaptLabel(dm, adaptLabel, &odm)); 628*c6f61ee2SBarry Smith adapted = PETSC_TRUE; 629*c6f61ee2SBarry Smith } 630*c6f61ee2SBarry Smith PetscCall(DMLabelDestroy(&adaptLabel)); 631*c6f61ee2SBarry Smith PetscCall(DMDestroy(&plex)); 632*c6f61ee2SBarry Smith } break; 633*c6f61ee2SBarry Smith case DM_ADAPTATION_METRIC: { 634*c6f61ee2SBarry Smith DM dmGrad, dmHess, dmMetric, dmDet; 635*c6f61ee2SBarry Smith Vec xGrad, xHess, metric, determinant; 636*c6f61ee2SBarry Smith PetscReal N; 637*c6f61ee2SBarry Smith DMLabel bdLabel = NULL, rgLabel = NULL; 638*c6f61ee2SBarry Smith PetscBool higherOrder = PETSC_FALSE; 639*c6f61ee2SBarry Smith PetscInt Nd = coordDim * coordDim, f, vStart, vEnd; 640*c6f61ee2SBarry Smith void (**funcs)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]); 641*c6f61ee2SBarry Smith 642*c6f61ee2SBarry Smith PetscCall(PetscMalloc(1, &funcs)); 643*c6f61ee2SBarry Smith funcs[0] = identityFunc; 644*c6f61ee2SBarry Smith 645*c6f61ee2SBarry Smith /* Setup finite element spaces */ 646*c6f61ee2SBarry Smith PetscCall(DMClone(dm, &dmGrad)); 647*c6f61ee2SBarry Smith PetscCall(DMClone(dm, &dmHess)); 648*c6f61ee2SBarry Smith PetscCheck(numFields <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Adaptation with multiple fields not yet considered"); // TODO 649*c6f61ee2SBarry Smith for (f = 0; f < numFields; ++f) { 650*c6f61ee2SBarry Smith PetscFE fe, feGrad, feHess; 651*c6f61ee2SBarry Smith PetscDualSpace Q; 652*c6f61ee2SBarry Smith PetscSpace space; 653*c6f61ee2SBarry Smith DM K; 654*c6f61ee2SBarry Smith PetscQuadrature q; 655*c6f61ee2SBarry Smith PetscInt Nc, qorder, p; 656*c6f61ee2SBarry Smith const char *prefix; 657*c6f61ee2SBarry Smith 658*c6f61ee2SBarry Smith PetscCall(PetscDSGetDiscretization(prob, f, (PetscObject *)&fe)); 659*c6f61ee2SBarry Smith PetscCall(PetscFEGetNumComponents(fe, &Nc)); 660*c6f61ee2SBarry Smith PetscCheck(Nc <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Adaptation with multiple components not yet considered"); // TODO 661*c6f61ee2SBarry Smith PetscCall(PetscFEGetBasisSpace(fe, &space)); 662*c6f61ee2SBarry Smith PetscCall(PetscSpaceGetDegree(space, NULL, &p)); 663*c6f61ee2SBarry Smith if (p > 1) higherOrder = PETSC_TRUE; 664*c6f61ee2SBarry Smith PetscCall(PetscFEGetDualSpace(fe, &Q)); 665*c6f61ee2SBarry Smith PetscCall(PetscDualSpaceGetDM(Q, &K)); 666*c6f61ee2SBarry Smith PetscCall(DMPlexGetDepthStratum(K, 0, &vStart, &vEnd)); 667*c6f61ee2SBarry Smith PetscCall(PetscFEGetQuadrature(fe, &q)); 668*c6f61ee2SBarry Smith PetscCall(PetscQuadratureGetOrder(q, &qorder)); 669*c6f61ee2SBarry Smith PetscCall(PetscObjectGetOptionsPrefix((PetscObject)fe, &prefix)); 670*c6f61ee2SBarry Smith PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dmGrad), dim, Nc * coordDim, PETSC_TRUE, prefix, qorder, &feGrad)); 671*c6f61ee2SBarry Smith PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dmHess), dim, Nc * Nd, PETSC_TRUE, prefix, qorder, &feHess)); 672*c6f61ee2SBarry Smith PetscCall(DMSetField(dmGrad, f, NULL, (PetscObject)feGrad)); 673*c6f61ee2SBarry Smith PetscCall(DMSetField(dmHess, f, NULL, (PetscObject)feHess)); 674*c6f61ee2SBarry Smith PetscCall(DMCreateDS(dmGrad)); 675*c6f61ee2SBarry Smith PetscCall(DMCreateDS(dmHess)); 676*c6f61ee2SBarry Smith PetscCall(PetscFEDestroy(&feGrad)); 677*c6f61ee2SBarry Smith PetscCall(PetscFEDestroy(&feHess)); 678*c6f61ee2SBarry Smith } 679*c6f61ee2SBarry Smith /* Compute vertexwise gradients from cellwise gradients */ 680*c6f61ee2SBarry Smith PetscCall(DMCreateLocalVector(dmGrad, &xGrad)); 681*c6f61ee2SBarry Smith PetscCall(VecViewFromOptions(locX, NULL, "-sol_adapt_loc_pre_view")); 682*c6f61ee2SBarry Smith PetscCall(DMPlexComputeGradientClementInterpolant(dm, locX, xGrad)); 683*c6f61ee2SBarry Smith PetscCall(VecViewFromOptions(xGrad, NULL, "-adapt_gradient_view")); 684*c6f61ee2SBarry Smith /* Compute vertexwise Hessians from cellwise Hessians */ 685*c6f61ee2SBarry Smith PetscCall(DMCreateLocalVector(dmHess, &xHess)); 686*c6f61ee2SBarry Smith PetscCall(DMPlexComputeGradientClementInterpolant(dmGrad, xGrad, xHess)); 687*c6f61ee2SBarry Smith PetscCall(VecViewFromOptions(xHess, NULL, "-adapt_hessian_view")); 688*c6f61ee2SBarry Smith PetscCall(VecDestroy(&xGrad)); 689*c6f61ee2SBarry Smith PetscCall(DMDestroy(&dmGrad)); 690*c6f61ee2SBarry Smith /* Compute L-p normalized metric */ 691*c6f61ee2SBarry Smith PetscCall(DMClone(dm, &dmMetric)); 692*c6f61ee2SBarry Smith N = adaptor->Nadapt >= 0 ? adaptor->Nadapt : PetscPowRealInt(adaptor->refinementFactor, dim) * ((PetscReal)(vEnd - vStart)); 693*c6f61ee2SBarry Smith if (adaptor->monitor) { 694*c6f61ee2SBarry Smith PetscMPIInt rank, size; 695*c6f61ee2SBarry Smith PetscCallMPI(MPI_Comm_rank(comm, &size)); 696*c6f61ee2SBarry Smith PetscCallMPI(MPI_Comm_rank(comm, &rank)); 697*c6f61ee2SBarry Smith PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] N_orig: %" PetscInt_FMT " N_adapt: %g\n", rank, vEnd - vStart, (double)N)); 698*c6f61ee2SBarry Smith } 699*c6f61ee2SBarry Smith PetscCall(DMPlexMetricSetTargetComplexity(dmMetric, (PetscReal)N)); 700*c6f61ee2SBarry Smith if (higherOrder) { 701*c6f61ee2SBarry Smith /* Project Hessian into P1 space, if required */ 702*c6f61ee2SBarry Smith PetscCall(DMPlexMetricCreate(dmMetric, 0, &metric)); 703*c6f61ee2SBarry Smith PetscCall(DMProjectFieldLocal(dmMetric, 0.0, xHess, funcs, INSERT_ALL_VALUES, metric)); 704*c6f61ee2SBarry Smith PetscCall(VecDestroy(&xHess)); 705*c6f61ee2SBarry Smith xHess = metric; 706*c6f61ee2SBarry Smith } 707*c6f61ee2SBarry Smith PetscCall(PetscFree(funcs)); 708*c6f61ee2SBarry Smith PetscCall(DMPlexMetricCreate(dmMetric, 0, &metric)); 709*c6f61ee2SBarry Smith PetscCall(DMPlexMetricDeterminantCreate(dmMetric, 0, &determinant, &dmDet)); 710*c6f61ee2SBarry Smith PetscCall(DMPlexMetricNormalize(dmMetric, xHess, PETSC_TRUE, PETSC_TRUE, metric, determinant)); 711*c6f61ee2SBarry Smith PetscCall(VecDestroy(&determinant)); 712*c6f61ee2SBarry Smith PetscCall(DMDestroy(&dmDet)); 713*c6f61ee2SBarry Smith PetscCall(VecDestroy(&xHess)); 714*c6f61ee2SBarry Smith PetscCall(DMDestroy(&dmHess)); 715*c6f61ee2SBarry Smith /* Adapt DM from metric */ 716*c6f61ee2SBarry Smith PetscCall(DMGetLabel(dm, "marker", &bdLabel)); 717*c6f61ee2SBarry Smith PetscCall(DMAdaptMetric(dm, metric, bdLabel, rgLabel, &odm)); 718*c6f61ee2SBarry Smith adapted = PETSC_TRUE; 719*c6f61ee2SBarry Smith /* Cleanup */ 720*c6f61ee2SBarry Smith PetscCall(VecDestroy(&metric)); 721*c6f61ee2SBarry Smith PetscCall(DMDestroy(&dmMetric)); 722*c6f61ee2SBarry Smith } break; 723*c6f61ee2SBarry Smith default: 724*c6f61ee2SBarry Smith SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Invalid adaptation type: %d", adaptor->adaptCriterion); 725*c6f61ee2SBarry Smith } 726*c6f61ee2SBarry Smith PetscCall(DMAdaptorPostAdapt(adaptor)); 727*c6f61ee2SBarry Smith PetscCall(DMRestoreLocalVector(dm, &locX)); 728*c6f61ee2SBarry Smith /* If DM was adapted, replace objects and recreate solution */ 729*c6f61ee2SBarry Smith if (adapted) { 730*c6f61ee2SBarry Smith const char *name; 731*c6f61ee2SBarry Smith 732*c6f61ee2SBarry Smith PetscCall(PetscObjectGetName((PetscObject)dm, &name)); 733*c6f61ee2SBarry Smith PetscCall(PetscObjectSetName((PetscObject)odm, name)); 734*c6f61ee2SBarry Smith /* Reconfigure solver */ 735*c6f61ee2SBarry Smith PetscCall(SNESReset(adaptor->snes)); 736*c6f61ee2SBarry Smith PetscCall(SNESSetDM(adaptor->snes, odm)); 737*c6f61ee2SBarry Smith PetscCall(DMAdaptorSetSolver(adaptor, adaptor->snes)); 738*c6f61ee2SBarry Smith PetscCall(DMPlexSetSNESLocalFEM(odm, PETSC_FALSE, ctx)); 739*c6f61ee2SBarry Smith PetscCall(SNESSetFromOptions(adaptor->snes)); 740*c6f61ee2SBarry Smith /* Transfer system */ 741*c6f61ee2SBarry Smith PetscCall(DMCopyDisc(dm, odm)); 742*c6f61ee2SBarry Smith /* Transfer solution */ 743*c6f61ee2SBarry Smith PetscCall(DMCreateGlobalVector(odm, &ox)); 744*c6f61ee2SBarry Smith PetscCall(PetscObjectGetName((PetscObject)x, &name)); 745*c6f61ee2SBarry Smith PetscCall(PetscObjectSetName((PetscObject)ox, name)); 746*c6f61ee2SBarry Smith PetscCall(DMAdaptorTransferSolution(adaptor, dm, x, odm, ox)); 747*c6f61ee2SBarry Smith /* Cleanup adaptivity info */ 748*c6f61ee2SBarry Smith if (adaptIter > 0) PetscCall(PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_(comm))); 749*c6f61ee2SBarry Smith PetscCall(DMForestSetAdaptivityForest(dm, NULL)); /* clear internal references to the previous dm */ 750*c6f61ee2SBarry Smith PetscCall(DMDestroy(&dm)); 751*c6f61ee2SBarry Smith PetscCall(VecDestroy(&x)); 752*c6f61ee2SBarry Smith *adm = odm; 753*c6f61ee2SBarry Smith *ax = ox; 754*c6f61ee2SBarry Smith } else { 755*c6f61ee2SBarry Smith *adm = dm; 756*c6f61ee2SBarry Smith *ax = x; 757*c6f61ee2SBarry Smith adaptIter = numAdapt; 758*c6f61ee2SBarry Smith } 759*c6f61ee2SBarry Smith if (adaptIter < numAdapt - 1) { 760*c6f61ee2SBarry Smith PetscCall(DMViewFromOptions(odm, NULL, "-dm_adapt_iter_view")); 761*c6f61ee2SBarry Smith PetscCall(VecViewFromOptions(ox, NULL, "-sol_adapt_iter_view")); 762*c6f61ee2SBarry Smith } 763*c6f61ee2SBarry Smith } 764*c6f61ee2SBarry Smith PetscCall(DMViewFromOptions(*adm, NULL, "-dm_adapt_view")); 765*c6f61ee2SBarry Smith PetscCall(VecViewFromOptions(*ax, NULL, "-sol_adapt_view")); 766*c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 767*c6f61ee2SBarry Smith } 768*c6f61ee2SBarry Smith 769*c6f61ee2SBarry Smith /*@ 770*c6f61ee2SBarry Smith DMAdaptorAdapt - Creates a new `DM` that is adapted to the problem 771*c6f61ee2SBarry Smith 772*c6f61ee2SBarry Smith Not Collective 773*c6f61ee2SBarry Smith 774*c6f61ee2SBarry Smith Input Parameters: 775*c6f61ee2SBarry Smith + adaptor - The `DMAdaptor` object 776*c6f61ee2SBarry Smith . x - The global approximate solution 777*c6f61ee2SBarry Smith - strategy - The adaptation strategy, see `DMAdaptationStrategy` 778*c6f61ee2SBarry Smith 779*c6f61ee2SBarry Smith Output Parameters: 780*c6f61ee2SBarry Smith + adm - The adapted `DM` 781*c6f61ee2SBarry Smith - ax - The adapted solution 782*c6f61ee2SBarry Smith 783*c6f61ee2SBarry Smith Options Database Keys: 784*c6f61ee2SBarry Smith + -snes_adapt <strategy> - initial, sequential, multigrid 785*c6f61ee2SBarry Smith . -adapt_gradient_view - View the Clement interpolant of the solution gradient 786*c6f61ee2SBarry Smith . -adapt_hessian_view - View the Clement interpolant of the solution Hessian 787*c6f61ee2SBarry Smith - -adapt_metric_view - View the metric tensor for adaptive mesh refinement 788*c6f61ee2SBarry Smith 789*c6f61ee2SBarry Smith Level: intermediate 790*c6f61ee2SBarry Smith 791*c6f61ee2SBarry Smith .seealso: `DMAdaptor`, `DMAdaptationStrategy`, `DMAdaptorSetSolver()`, `DMAdaptorCreate()` 792*c6f61ee2SBarry Smith @*/ 793*c6f61ee2SBarry Smith PetscErrorCode DMAdaptorAdapt(DMAdaptor adaptor, Vec x, DMAdaptationStrategy strategy, DM *adm, Vec *ax) 794*c6f61ee2SBarry Smith { 795*c6f61ee2SBarry Smith PetscFunctionBegin; 796*c6f61ee2SBarry Smith switch (strategy) { 797*c6f61ee2SBarry Smith case DM_ADAPTATION_INITIAL: 798*c6f61ee2SBarry Smith PetscCall(DMAdaptorAdapt_Sequence_Private(adaptor, x, PETSC_FALSE, adm, ax)); 799*c6f61ee2SBarry Smith break; 800*c6f61ee2SBarry Smith case DM_ADAPTATION_SEQUENTIAL: 801*c6f61ee2SBarry Smith PetscCall(DMAdaptorAdapt_Sequence_Private(adaptor, x, PETSC_TRUE, adm, ax)); 802*c6f61ee2SBarry Smith break; 803*c6f61ee2SBarry Smith default: 804*c6f61ee2SBarry Smith SETERRQ(PetscObjectComm((PetscObject)adaptor), PETSC_ERR_ARG_WRONG, "Unrecognized adaptation strategy %d", strategy); 805*c6f61ee2SBarry Smith } 806*c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 807*c6f61ee2SBarry Smith } 808