1 #include <petsc/private/dmimpl.h> /*I "petscdm.h" I*/ 2 3 PETSC_EXTERN PetscErrorCode DMIsForest(DM, PetscBool *); 4 5 DMGeneratorFunctionList DMGenerateList = NULL; 6 PetscBool DMGenerateRegisterAllCalled = PETSC_FALSE; 7 8 #if defined(PETSC_HAVE_TRIANGLE) 9 PETSC_EXTERN PetscErrorCode DMPlexGenerate_Triangle(DM, PetscBool, DM *); 10 PETSC_EXTERN PetscErrorCode DMPlexRefine_Triangle(DM, PetscReal *, DM *); 11 #endif 12 #if defined(PETSC_HAVE_TETGEN) 13 PETSC_EXTERN PetscErrorCode DMPlexGenerate_Tetgen(DM, PetscBool, DM *); 14 PETSC_EXTERN PetscErrorCode DMPlexRefine_Tetgen(DM, double *, DM *); 15 #endif 16 #if defined(PETSC_HAVE_CTETGEN) 17 PETSC_EXTERN PetscErrorCode DMPlexGenerate_CTetgen(DM, PetscBool, DM *); 18 PETSC_EXTERN PetscErrorCode DMPlexRefine_CTetgen(DM, double *, DM *); 19 #endif 20 #if defined(PETSC_HAVE_PRAGMATIC) 21 PETSC_EXTERN PetscErrorCode DMAdaptMetric_Pragmatic_Plex(DM, Vec, DMLabel, DMLabel, DM *); 22 #endif 23 #if defined(PETSC_HAVE_MMG) 24 PETSC_EXTERN PetscErrorCode DMAdaptMetric_Mmg_Plex(DM, Vec, DMLabel, DMLabel, DM *); 25 #endif 26 #if defined(PETSC_HAVE_PARMMG) 27 PETSC_EXTERN PetscErrorCode DMAdaptMetric_ParMmg_Plex(DM, Vec, DMLabel, DMLabel, DM *); 28 #endif 29 PETSC_EXTERN PetscErrorCode DMPlexTransformAdaptLabel(DM, Vec, DMLabel, DMLabel, DM *); 30 PETSC_EXTERN PetscErrorCode DMAdaptLabel_Forest(DM, Vec, DMLabel, DMLabel, DM *); 31 32 /*@C 33 DMGenerateRegisterAll - Registers all of the mesh generation methods in the `DM` package. 34 35 Not Collective 36 37 Level: advanced 38 39 .seealso: `DM`, `DMGenerateRegisterDestroy()` 40 @*/ 41 PetscErrorCode DMGenerateRegisterAll(void) 42 { 43 PetscFunctionBegin; 44 if (DMGenerateRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS); 45 DMGenerateRegisterAllCalled = PETSC_TRUE; 46 #if defined(PETSC_HAVE_TRIANGLE) 47 PetscCall(DMGenerateRegister("triangle", DMPlexGenerate_Triangle, DMPlexRefine_Triangle, NULL, 1)); 48 #endif 49 #if defined(PETSC_HAVE_CTETGEN) 50 PetscCall(DMGenerateRegister("ctetgen", DMPlexGenerate_CTetgen, DMPlexRefine_CTetgen, NULL, 2)); 51 #endif 52 #if defined(PETSC_HAVE_TETGEN) 53 PetscCall(DMGenerateRegister("tetgen", DMPlexGenerate_Tetgen, DMPlexRefine_Tetgen, NULL, 2)); 54 #endif 55 #if defined(PETSC_HAVE_PRAGMATIC) 56 PetscCall(DMGenerateRegister("pragmatic", NULL, NULL, DMAdaptMetric_Pragmatic_Plex, -1)); 57 #endif 58 #if defined(PETSC_HAVE_MMG) 59 PetscCall(DMGenerateRegister("mmg", NULL, NULL, DMAdaptMetric_Mmg_Plex, -1)); 60 #endif 61 #if defined(PETSC_HAVE_PARMMG) 62 PetscCall(DMGenerateRegister("parmmg", NULL, NULL, DMAdaptMetric_ParMmg_Plex, -1)); 63 #endif 64 PetscCall(DMGenerateRegister("cellrefiner", NULL, NULL, DMPlexTransformAdaptLabel, -1)); 65 PetscCall(DMGenerateRegister("forest", NULL, NULL, DMAdaptLabel_Forest, -1)); 66 PetscFunctionReturn(PETSC_SUCCESS); 67 } 68 69 /*@C 70 DMGenerateRegister - Adds a grid generator to `DM` 71 72 Not Collective, No Fortran Support 73 74 Input Parameters: 75 + sname - name of a new user-defined grid generator 76 . fnc - generator function 77 . rfnc - refinement function 78 . alfnc - adapt by label function 79 - dim - dimension of boundary of domain 80 81 Example Usage: 82 .vb 83 DMGenerateRegister("my_generator", MyGeneratorCreate, MyGeneratorRefiner, MyGeneratorAdaptor, dim); 84 .ve 85 86 Then, your generator can be chosen with the procedural interface via 87 $ DMGenerate(dm, "my_generator",...) 88 or at runtime via the option 89 $ -dm_generator my_generator 90 91 Level: advanced 92 93 Note: 94 `DMGenerateRegister()` may be called multiple times to add several user-defined generators 95 96 .seealso: `DM`, `DMGenerateRegisterAll()`, `DMPlexGenerate()`, `DMGenerateRegisterDestroy()` 97 @*/ 98 PetscErrorCode DMGenerateRegister(const char sname[], PetscErrorCode (*fnc)(DM, PetscBool, DM *), PetscErrorCode (*rfnc)(DM, PetscReal *, DM *), PetscErrorCode (*alfnc)(DM, Vec, DMLabel, DMLabel, DM *), PetscInt dim) 99 { 100 DMGeneratorFunctionList entry; 101 102 PetscFunctionBegin; 103 PetscCall(PetscNew(&entry)); 104 PetscCall(PetscStrallocpy(sname, &entry->name)); 105 entry->generate = fnc; 106 entry->refine = rfnc; 107 entry->adapt = alfnc; 108 entry->dim = dim; 109 entry->next = NULL; 110 if (!DMGenerateList) DMGenerateList = entry; 111 else { 112 DMGeneratorFunctionList fl = DMGenerateList; 113 while (fl->next) fl = fl->next; 114 fl->next = entry; 115 } 116 PetscFunctionReturn(PETSC_SUCCESS); 117 } 118 119 PetscErrorCode DMGenerateRegisterDestroy(void) 120 { 121 DMGeneratorFunctionList next, fl; 122 123 PetscFunctionBegin; 124 next = fl = DMGenerateList; 125 while (next) { 126 next = fl ? fl->next : NULL; 127 if (fl) PetscCall(PetscFree(fl->name)); 128 PetscCall(PetscFree(fl)); 129 fl = next; 130 } 131 DMGenerateList = NULL; 132 DMGenerateRegisterAllCalled = PETSC_FALSE; 133 PetscFunctionReturn(PETSC_SUCCESS); 134 } 135 136 /*@ 137 DMAdaptLabel - Adapt a `DM` based on a `DMLabel` with values interpreted as coarsening and refining flags. Specific implementations of `DM` maybe have 138 specialized flags, but all implementations should accept flag values `DM_ADAPT_DETERMINE`, `DM_ADAPT_KEEP`, `DM_ADAPT_REFINE`, and, 139 `DM_ADAPT_COARSEN`. 140 141 Collective 142 143 Input Parameters: 144 + dm - the pre-adaptation `DM` object 145 - label - label with the flags 146 147 Output Parameters: 148 . dmAdapt - the adapted `DM` object: may be `NULL` if an adapted `DM` could not be produced. 149 150 Level: intermediate 151 152 .seealso: `DM`, `DMAdaptMetric()`, `DMCoarsen()`, `DMRefine()` 153 @*/ 154 PetscErrorCode DMAdaptLabel(DM dm, DMLabel label, DM *dmAdapt) 155 { 156 DMGeneratorFunctionList fl; 157 char adaptname[PETSC_MAX_PATH_LEN]; 158 const char *name; 159 PetscInt dim; 160 PetscBool flg, isForest, found = PETSC_FALSE; 161 162 PetscFunctionBegin; 163 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 164 if (label) PetscValidHeaderSpecific(label, DMLABEL_CLASSID, 2); 165 PetscAssertPointer(dmAdapt, 3); 166 *dmAdapt = NULL; 167 PetscCall(DMGetDimension(dm, &dim)); 168 PetscCall(DMIsForest(dm, &isForest)); 169 name = isForest ? "forest" : "cellrefiner"; 170 PetscCall(PetscOptionsGetString(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_adaptor", adaptname, sizeof(adaptname), &flg)); 171 if (flg) name = adaptname; 172 173 fl = DMGenerateList; 174 while (fl) { 175 PetscCall(PetscStrcmp(fl->name, name, &flg)); 176 if (flg) { 177 PetscCall((*fl->adapt)(dm, NULL, label, NULL, dmAdapt)); 178 found = PETSC_TRUE; 179 } 180 fl = fl->next; 181 } 182 PetscCheck(found, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Grid adaptor %s not registered; you may need to add --download-%s to your ./configure options", name, name); 183 if (*dmAdapt) { 184 void *ctx; 185 186 (*dmAdapt)->prealloc_only = dm->prealloc_only; /* maybe this should go .... */ 187 PetscCall(PetscFree((*dmAdapt)->vectype)); 188 PetscCall(PetscStrallocpy(dm->vectype, (char **)&(*dmAdapt)->vectype)); 189 PetscCall(PetscFree((*dmAdapt)->mattype)); 190 PetscCall(PetscStrallocpy(dm->mattype, (char **)&(*dmAdapt)->mattype)); 191 PetscCall(DMGetApplicationContext(dm, &ctx)); 192 PetscCall(DMSetApplicationContext(*dmAdapt, ctx)); 193 } 194 PetscFunctionReturn(PETSC_SUCCESS); 195 } 196 197 /*@ 198 DMAdaptMetric - Generates a mesh adapted to the specified metric field. 199 200 Input Parameters: 201 + dm - The DM object 202 . metric - The metric to which the mesh is adapted, defined vertex-wise. 203 . bdLabel - Label for boundary tags, which will be preserved in the output mesh. `bdLabel` should be `NULL` if there is no such label, and should be different from "_boundary_". 204 - rgLabel - Label for cell tags, which will be preserved in the output mesh. `rgLabel` should be `NULL` if there is no such label, and should be different from "_regions_". 205 206 Output Parameter: 207 . dmAdapt - Pointer to the `DM` object containing the adapted mesh 208 209 Note: 210 The label in the adapted mesh will be registered under the name of the input `DMLabel` object 211 212 Level: advanced 213 214 .seealso: `DMAdaptLabel()`, `DMCoarsen()`, `DMRefine()` 215 @*/ 216 PetscErrorCode DMAdaptMetric(DM dm, Vec metric, DMLabel bdLabel, DMLabel rgLabel, DM *dmAdapt) 217 { 218 DMGeneratorFunctionList fl; 219 char adaptname[PETSC_MAX_PATH_LEN]; 220 const char *name; 221 const char *const adaptors[3] = {"pragmatic", "mmg", "parmmg"}; 222 PetscInt dim; 223 PetscBool flg, found = PETSC_FALSE; 224 225 PetscFunctionBegin; 226 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 227 PetscValidHeaderSpecific(metric, VEC_CLASSID, 2); 228 if (bdLabel) PetscValidHeaderSpecific(bdLabel, DMLABEL_CLASSID, 3); 229 if (rgLabel) PetscValidHeaderSpecific(rgLabel, DMLABEL_CLASSID, 4); 230 PetscAssertPointer(dmAdapt, 5); 231 *dmAdapt = NULL; 232 PetscCall(DMGetDimension(dm, &dim)); 233 PetscCall(PetscOptionsGetString(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_adaptor", adaptname, sizeof(adaptname), &flg)); 234 235 /* Default to Mmg in serial and ParMmg in parallel */ 236 if (flg) name = adaptname; 237 else { 238 MPI_Comm comm; 239 PetscMPIInt size; 240 241 PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 242 PetscCallMPI(MPI_Comm_size(comm, &size)); 243 if (size == 1) name = adaptors[1]; 244 else name = adaptors[2]; 245 } 246 247 fl = DMGenerateList; 248 while (fl) { 249 PetscCall(PetscStrcmp(fl->name, name, &flg)); 250 if (flg) { 251 PetscCall((*fl->adapt)(dm, metric, bdLabel, rgLabel, dmAdapt)); 252 found = PETSC_TRUE; 253 } 254 fl = fl->next; 255 } 256 PetscCheck(found, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Grid adaptor %s not registered; you may need to add --download-%s to your ./configure options", name, name); 257 if (*dmAdapt) { 258 (*dmAdapt)->prealloc_only = dm->prealloc_only; /* maybe this should go .... */ 259 PetscCall(PetscFree((*dmAdapt)->vectype)); 260 PetscCall(PetscStrallocpy(dm->vectype, (char **)&(*dmAdapt)->vectype)); 261 PetscCall(PetscFree((*dmAdapt)->mattype)); 262 PetscCall(PetscStrallocpy(dm->mattype, (char **)&(*dmAdapt)->mattype)); 263 } 264 PetscFunctionReturn(PETSC_SUCCESS); 265 } 266