Actual source code: plexgenerate.c

petsc-3.12.5 2020-03-29
Report Typos and Errors
  1:  #include <petsc/private/dmpleximpl.h>

  3: PetscErrorCode DMPlexInvertCell_Internal(PetscInt dim, PetscInt numCorners, PetscInt cone[])
  4: {
  5:   int tmpc;

  8:   if (dim != 3) return(0);
  9:   switch (numCorners) {
 10:   case 4:
 11:     tmpc    = cone[0];
 12:     cone[0] = cone[1];
 13:     cone[1] = tmpc;
 14:     break;
 15:   case 8:
 16:     tmpc    = cone[1];
 17:     cone[1] = cone[3];
 18:     cone[3] = tmpc;
 19:     break;
 20:   default: break;
 21:   }
 22:   return(0);
 23: }

 25: /*@C
 26:   DMPlexInvertCell - This flips tetrahedron and hexahedron orientation since Plex stores them internally with outward normals. Other cells are left untouched.

 28:   Input Parameters:
 29: + numCorners - The number of vertices in a cell
 30: - cone - The incoming cone

 32:   Output Parameter:
 33: . cone - The inverted cone (in-place)

 35:   Level: developer

 37: .seealso: DMPlexGenerate()
 38: @*/
 39: PetscErrorCode DMPlexInvertCell(PetscInt dim, PetscInt numCorners, int cone[])
 40: {
 41:   int tmpc;

 44:   if (dim != 3) return(0);
 45:   switch (numCorners) {
 46:   case 4:
 47:     tmpc    = cone[0];
 48:     cone[0] = cone[1];
 49:     cone[1] = tmpc;
 50:     break;
 51:   case 8:
 52:     tmpc    = cone[1];
 53:     cone[1] = cone[3];
 54:     cone[3] = tmpc;
 55:     break;
 56:   default: break;
 57:   }
 58:   return(0);
 59: }


 62: /*@C
 63:   DMPlexTriangleSetOptions - Set the options used for the Triangle mesh generator

 65:   Not Collective

 67:   Inputs Parameters:
 68: + dm - The DMPlex object
 69: - opts - The command line options

 71:   Level: developer

 73: .seealso: DMPlexTetgenSetOptions(), DMPlexGenerate()
 74: @*/
 75: PetscErrorCode DMPlexTriangleSetOptions(DM dm, const char *opts)
 76: {
 77:   DM_Plex       *mesh = (DM_Plex*) dm->data;

 83:   PetscFree(mesh->triangleOpts);
 84:   PetscStrallocpy(opts, &mesh->triangleOpts);
 85:   return(0);
 86: }

 88: /*@C
 89:   DMPlexTetgenSetOptions - Set the options used for the Tetgen mesh generator

 91:   Not Collective

 93:   Inputs Parameters:
 94: + dm - The DMPlex object
 95: - opts - The command line options

 97:   Level: developer

 99: .seealso: DMPlexTriangleSetOptions(), DMPlexGenerate()
100: @*/
101: PetscErrorCode DMPlexTetgenSetOptions(DM dm, const char *opts)
102: {
103:   DM_Plex       *mesh = (DM_Plex*) dm->data;

109:   PetscFree(mesh->tetgenOpts);
110:   PetscStrallocpy(opts, &mesh->tetgenOpts);
111:   return(0);
112: }

114: /*
115:    Contains the list of registered DMPlexGenerators routines
116: */
117: PetscFunctionList DMPlexGenerateList = NULL;

119: struct _n_PetscFunctionList {
120:   PetscErrorCode    (*generate)(DM, PetscBool, DM*);
121:   PetscErrorCode    (*refine)(DM,double*, DM*);
122:   char              *name;               /* string to identify routine */
123:   PetscInt          dim;
124:   PetscFunctionList next;                /* next pointer */
125: };

127: /*@C
128:   DMPlexGenerate - Generates a mesh.

130:   Not Collective

132:   Input Parameters:
133: + boundary - The DMPlex boundary object
134: . name - The mesh generation package name
135: - interpolate - Flag to create intermediate mesh elements

137:   Output Parameter:
138: . mesh - The DMPlex object

140:   Level: intermediate

142: .seealso: DMPlexCreate(), DMRefine()
143: @*/
144: PetscErrorCode DMPlexGenerate(DM boundary, const char name[], PetscBool interpolate, DM *mesh)
145: {
146:   PetscInt          dim;
147:   char              genname[1024];
148:   PetscBool         flg;
149:   PetscErrorCode    ierr;
150:   PetscFunctionList fl;

155:   DMGetDimension(boundary, &dim);
156:   PetscOptionsGetString(((PetscObject) boundary)->options,((PetscObject) boundary)->prefix, "-dm_plex_generator", genname, 1024, &flg);
157:   if (flg) name = genname;

159:   fl = DMPlexGenerateList;
160:   if (name) {
161:     while (fl) {
162:       PetscStrcmp(fl->name,name,&flg);
163:       if (flg) {
164:         (*fl->generate)(boundary,interpolate,mesh);
165:         return(0);
166:       }
167:       fl = fl->next;
168:     }
169:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Grid generator %g not registered",name);
170:   } else {
171:     while (fl) {
172:       if (boundary->dim == fl->dim) {
173:         (*fl->generate)(boundary,interpolate,mesh);
174:         return(0);
175:       }
176:       fl = fl->next;
177:     }
178:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"No grid generator of dimension %D registered",boundary->dim+1);
179:   }
180:   return(0);
181: }

183: /*@C
184:   DMPlexGenerateRegister -  Adds a grid generator to DMPlex

186:    Not Collective

188:    Input Parameters:
189: +  name_solver - name of a new user-defined grid generator
190: .  fnc - generator function
191: .  rfnc - refinement function
192: -  dim - dimension of boundary of domain

194:    Notes:
195:    DMPlexGenerateRegister() may be called multiple times to add several user-defined solvers.

197:    Sample usage:
198: .vb
199:    DMPlexGenerateRegister("my_generator",MyGeneratorCreate,MyGeneratorRefiner,dim);
200: .ve

202:    Then, your generator can be chosen with the procedural interface via
203: $     DMPlexGenerate(dm,"my_generator",...)
204:    or at runtime via the option
205: $     -dm_plex_generator my_generator

207:    Level: advanced

209: .seealso: DMPlexGenerateRegisterAll(), DMPlexGenerate(), DMPlexGenerateRegisterDestroy()

211: @*/
212: PetscErrorCode  DMPlexGenerateRegister(const char sname[],PetscErrorCode (*fnc)(DM, PetscBool,DM*), PetscErrorCode (*rfnc)(DM, double*,DM*),PetscInt dim)
213: {
214:   PetscErrorCode    ierr;
215:   PetscFunctionList entry;

218:   PetscNew(&entry);
219:   PetscStrallocpy(sname,&entry->name);
220:   entry->generate = fnc;
221:   entry->refine   = rfnc;
222:   entry->dim      = dim;
223:   entry->next     = NULL;
224:   if (!DMPlexGenerateList) DMPlexGenerateList = entry;
225:   else {
226:     PetscFunctionList fl = DMPlexGenerateList;
227:     while (fl->next) fl = fl->next;
228:     fl->next = entry;
229:   }
230:   return(0);
231: }

233: extern PetscBool DMPlexGenerateRegisterAllCalled;

235: PetscErrorCode  DMPlexGenerateRegisterDestroy(void)
236: {
237:   PetscFunctionList next,fl;
238:   PetscErrorCode    ierr;

241:   next = fl =  DMPlexGenerateList;
242:     while (next) {
243:     next = fl ? fl->next : NULL;
244:     if (fl) {PetscFree(fl->name);}
245:     PetscFree(fl);
246:     fl = next;
247:   }
248:   DMPlexGenerateList              = NULL;
249:   DMPlexGenerateRegisterAllCalled = PETSC_FALSE;
250:   return(0);
251: }