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