Actual source code: space.c

  1: #include <petsc/private/petscfeimpl.h>
  2: #include <petscdmshell.h>

  4: PetscClassId PETSCSPACE_CLASSID = 0;

  6: PetscFunctionList PetscSpaceList              = NULL;
  7: PetscBool         PetscSpaceRegisterAllCalled = PETSC_FALSE;

  9: /*@C
 10:   PetscSpaceRegister - Adds a new PetscSpace implementation

 12:   Not Collective

 14:   Input Parameters:
 15: + name        - The name of a new user-defined creation routine
 16: - create_func - The creation routine for the implementation type

 18:   Notes:
 19:   PetscSpaceRegister() may be called multiple times to add several user-defined types of PetscSpaces.  The creation function is called
 20:   when the type is set to 'name'.

 22:   Sample usage:
 23: .vb
 24:     PetscSpaceRegister("my_space", MyPetscSpaceCreate);
 25: .ve

 27:   Then, your PetscSpace type can be chosen with the procedural interface via
 28: .vb
 29:     PetscSpaceCreate(MPI_Comm, PetscSpace *);
 30:     PetscSpaceSetType(PetscSpace, "my_space");
 31: .ve
 32:    or at runtime via the option
 33: .vb
 34:     -petscspace_type my_space
 35: .ve

 37:   Level: advanced

 39: .seealso: PetscSpaceRegisterAll(), PetscSpaceRegisterDestroy()

 41: @*/
 42: PetscErrorCode PetscSpaceRegister(const char sname[], PetscErrorCode (*function)(PetscSpace))
 43: {

 47:   PetscFunctionListAdd(&PetscSpaceList, sname, function);
 48:   return(0);
 49: }

 51: /*@C
 52:   PetscSpaceSetType - Builds a particular PetscSpace

 54:   Collective on sp

 56:   Input Parameters:
 57: + sp   - The PetscSpace object
 58: - name - The kind of space

 60:   Options Database Key:
 61: . -petscspace_type <type> - Sets the PetscSpace type; use -help for a list of available types

 63:   Level: intermediate

 65: .seealso: PetscSpaceGetType(), PetscSpaceCreate()
 66: @*/
 67: PetscErrorCode PetscSpaceSetType(PetscSpace sp, PetscSpaceType name)
 68: {
 69:   PetscErrorCode (*r)(PetscSpace);
 70:   PetscBool      match;

 75:   PetscObjectTypeCompare((PetscObject) sp, name, &match);
 76:   if (match) return(0);

 78:   PetscSpaceRegisterAll();
 79:   PetscFunctionListFind(PetscSpaceList, name, &r);
 80:   if (!r) SETERRQ1(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscSpace type: %s", name);

 82:   if (sp->ops->destroy) {
 83:     (*sp->ops->destroy)(sp);
 84:     sp->ops->destroy = NULL;
 85:   }
 86:   sp->dim = PETSC_DETERMINE;
 87:   (*r)(sp);
 88:   PetscObjectChangeTypeName((PetscObject) sp, name);
 89:   return(0);
 90: }

 92: /*@C
 93:   PetscSpaceGetType - Gets the PetscSpace type name (as a string) from the object.

 95:   Not Collective

 97:   Input Parameter:
 98: . sp  - The PetscSpace

100:   Output Parameter:
101: . name - The PetscSpace type name

103:   Level: intermediate

105: .seealso: PetscSpaceSetType(), PetscSpaceCreate()
106: @*/
107: PetscErrorCode PetscSpaceGetType(PetscSpace sp, PetscSpaceType *name)
108: {

114:   if (!PetscSpaceRegisterAllCalled) {
115:     PetscSpaceRegisterAll();
116:   }
117:   *name = ((PetscObject) sp)->type_name;
118:   return(0);
119: }

121: /*@C
122:    PetscSpaceViewFromOptions - View from Options

124:    Collective on PetscSpace

126:    Input Parameters:
127: +  A - the PetscSpace object
128: .  obj - Optional object
129: -  name - command line option

131:    Level: intermediate
132: .seealso:  PetscSpace, PetscSpaceView, PetscObjectViewFromOptions(), PetscSpaceCreate()
133: @*/
134: PetscErrorCode  PetscSpaceViewFromOptions(PetscSpace A,PetscObject obj,const char name[])
135: {

140:   PetscObjectViewFromOptions((PetscObject)A,obj,name);
141:   return(0);
142: }

144: /*@C
145:   PetscSpaceView - Views a PetscSpace

147:   Collective on sp

149:   Input Parameter:
150: + sp - the PetscSpace object to view
151: - v  - the viewer

153:   Level: beginner

155: .seealso PetscSpaceDestroy()
156: @*/
157: PetscErrorCode PetscSpaceView(PetscSpace sp, PetscViewer v)
158: {
159:   PetscInt       pdim;
160:   PetscBool      iascii;

166:   if (!v) {PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) sp), &v);}
167:   PetscSpaceGetDimension(sp, &pdim);
168:   PetscObjectPrintClassNamePrefixType((PetscObject)sp,v);
169:   PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii);
170:   PetscViewerASCIIPushTab(v);
171:   if (iascii) {PetscViewerASCIIPrintf(v, "Space in %D variables with %D components, size %D\n", sp->Nv, sp->Nc, pdim);}
172:   if (sp->ops->view) {(*sp->ops->view)(sp, v);}
173:   PetscViewerASCIIPopTab(v);
174:   return(0);
175: }

177: /*@
178:   PetscSpaceSetFromOptions - sets parameters in a PetscSpace from the options database

180:   Collective on sp

182:   Input Parameter:
183: . sp - the PetscSpace object to set options for

185:   Options Database:
186: + -petscspace_degree <deg> - the approximation order of the space
187: . -petscspace_variables <n> - the number of different variables, e.g. x and y
188: - -petscspace_components <c> - the number of components, say d for a vector field

190:   Level: intermediate

192: .seealso PetscSpaceView()
193: @*/
194: PetscErrorCode PetscSpaceSetFromOptions(PetscSpace sp)
195: {
196:   const char    *defaultType;
197:   char           name[256];
198:   PetscBool      flg;

203:   if (!((PetscObject) sp)->type_name) {
204:     defaultType = PETSCSPACEPOLYNOMIAL;
205:   } else {
206:     defaultType = ((PetscObject) sp)->type_name;
207:   }
208:   if (!PetscSpaceRegisterAllCalled) {PetscSpaceRegisterAll();}

210:   PetscObjectOptionsBegin((PetscObject) sp);
211:   PetscOptionsFList("-petscspace_type", "Linear space", "PetscSpaceSetType", PetscSpaceList, defaultType, name, 256, &flg);
212:   if (flg) {
213:     PetscSpaceSetType(sp, name);
214:   } else if (!((PetscObject) sp)->type_name) {
215:     PetscSpaceSetType(sp, defaultType);
216:   }
217:   {
218:     PetscOptionsDeprecated("-petscspace_order","-petscspace_degree","3.11",NULL);
219:     PetscOptionsBoundedInt("-petscspace_order", "DEPRECATED: The approximation order", "PetscSpaceSetDegree", sp->degree, &sp->degree, NULL,0);
220:   }
221:   PetscOptionsBoundedInt("-petscspace_degree", "The (maximally included) polynomial degree", "PetscSpaceSetDegree", sp->degree, &sp->degree, NULL,0);
222:   PetscOptionsBoundedInt("-petscspace_variables", "The number of different variables, e.g. x and y", "PetscSpaceSetNumVariables", sp->Nv, &sp->Nv, NULL,0);
223:   PetscOptionsBoundedInt("-petscspace_components", "The number of components", "PetscSpaceSetNumComponents", sp->Nc, &sp->Nc, NULL,0);
224:   if (sp->ops->setfromoptions) {
225:     (*sp->ops->setfromoptions)(PetscOptionsObject,sp);
226:   }
227:   /* process any options handlers added with PetscObjectAddOptionsHandler() */
228:   PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) sp);
229:   PetscOptionsEnd();
230:   PetscSpaceViewFromOptions(sp, NULL, "-petscspace_view");
231:   return(0);
232: }

234: /*@C
235:   PetscSpaceSetUp - Construct data structures for the PetscSpace

237:   Collective on sp

239:   Input Parameter:
240: . sp - the PetscSpace object to setup

242:   Level: intermediate

244: .seealso PetscSpaceView(), PetscSpaceDestroy()
245: @*/
246: PetscErrorCode PetscSpaceSetUp(PetscSpace sp)
247: {

252:   if (sp->ops->setup) {(*sp->ops->setup)(sp);}
253:   return(0);
254: }

256: /*@
257:   PetscSpaceDestroy - Destroys a PetscSpace object

259:   Collective on sp

261:   Input Parameter:
262: . sp - the PetscSpace object to destroy

264:   Level: beginner

266: .seealso PetscSpaceView()
267: @*/
268: PetscErrorCode PetscSpaceDestroy(PetscSpace *sp)
269: {

273:   if (!*sp) return(0);

276:   if (--((PetscObject)(*sp))->refct > 0) {*sp = NULL; return(0);}
277:   ((PetscObject) (*sp))->refct = 0;
278:   DMDestroy(&(*sp)->dm);

280:   (*(*sp)->ops->destroy)(*sp);
281:   PetscHeaderDestroy(sp);
282:   return(0);
283: }

285: /*@
286:   PetscSpaceCreate - Creates an empty PetscSpace object. The type can then be set with PetscSpaceSetType().

288:   Collective

290:   Input Parameter:
291: . comm - The communicator for the PetscSpace object

293:   Output Parameter:
294: . sp - The PetscSpace object

296:   Level: beginner

298: .seealso: PetscSpaceSetType(), PETSCSPACEPOLYNOMIAL
299: @*/
300: PetscErrorCode PetscSpaceCreate(MPI_Comm comm, PetscSpace *sp)
301: {
302:   PetscSpace     s;

307:   PetscCitationsRegister(FECitation,&FEcite);
308:   *sp  = NULL;
309:   PetscFEInitializePackage();

311:   PetscHeaderCreate(s, PETSCSPACE_CLASSID, "PetscSpace", "Linear Space", "PetscSpace", comm, PetscSpaceDestroy, PetscSpaceView);

313:   s->degree    = 0;
314:   s->maxDegree = PETSC_DETERMINE;
315:   s->Nc        = 1;
316:   s->Nv        = 0;
317:   s->dim       = PETSC_DETERMINE;
318:   DMShellCreate(comm, &s->dm);
319:   PetscSpaceSetType(s, PETSCSPACEPOLYNOMIAL);

321:   *sp = s;
322:   return(0);
323: }

325: /*@
326:   PetscSpaceGetDimension - Return the dimension of this space, i.e. the number of basis vectors

328:   Input Parameter:
329: . sp - The PetscSpace

331:   Output Parameter:
332: . dim - The dimension

334:   Level: intermediate

336: .seealso: PetscSpaceGetDegree(), PetscSpaceCreate(), PetscSpace
337: @*/
338: PetscErrorCode PetscSpaceGetDimension(PetscSpace sp, PetscInt *dim)
339: {

345:   if (sp->dim == PETSC_DETERMINE) {
346:     if (sp->ops->getdimension) {(*sp->ops->getdimension)(sp, &sp->dim);}
347:   }
348:   *dim = sp->dim;
349:   return(0);
350: }

352: /*@
353:   PetscSpaceGetDegree - Return the polynomial degrees that characterize this space

355:   Input Parameter:
356: . sp - The PetscSpace

358:   Output Parameter:
359: + minDegree - The degree of the largest polynomial space contained in the space
360: - maxDegree - The degree of the smallest polynomial space containing the space


363:   Level: intermediate

365: .seealso: PetscSpaceSetDegree(), PetscSpaceGetDimension(), PetscSpaceCreate(), PetscSpace
366: @*/
367: PetscErrorCode PetscSpaceGetDegree(PetscSpace sp, PetscInt *minDegree, PetscInt *maxDegree)
368: {
373:   if (minDegree) *minDegree = sp->degree;
374:   if (maxDegree) *maxDegree = sp->maxDegree;
375:   return(0);
376: }

378: /*@
379:   PetscSpaceSetDegree - Set the degree of approximation for this space.

381:   Input Parameters:
382: + sp - The PetscSpace
383: . degree - The degree of the largest polynomial space contained in the space
384: - maxDegree - The degree of the largest polynomial space containing the space.  One of degree and maxDegree can be PETSC_DETERMINE.

386:   Level: intermediate

388: .seealso: PetscSpaceGetDegree(), PetscSpaceCreate(), PetscSpace
389: @*/
390: PetscErrorCode PetscSpaceSetDegree(PetscSpace sp, PetscInt degree, PetscInt maxDegree)
391: {
394:   sp->degree = degree;
395:   sp->maxDegree = maxDegree;
396:   return(0);
397: }

399: /*@
400:   PetscSpaceGetNumComponents - Return the number of components for this space

402:   Input Parameter:
403: . sp - The PetscSpace

405:   Output Parameter:
406: . Nc - The number of components

408:   Note: A vector space, for example, will have d components, where d is the spatial dimension

410:   Level: intermediate

412: .seealso: PetscSpaceSetNumComponents(), PetscSpaceGetNumVariables(), PetscSpaceGetDimension(), PetscSpaceCreate(), PetscSpace
413: @*/
414: PetscErrorCode PetscSpaceGetNumComponents(PetscSpace sp, PetscInt *Nc)
415: {
419:   *Nc = sp->Nc;
420:   return(0);
421: }

423: /*@
424:   PetscSpaceSetNumComponents - Set the number of components for this space

426:   Input Parameters:
427: + sp - The PetscSpace
428: - order - The number of components

430:   Level: intermediate

432: .seealso: PetscSpaceGetNumComponents(), PetscSpaceSetNumVariables(), PetscSpaceCreate(), PetscSpace
433: @*/
434: PetscErrorCode PetscSpaceSetNumComponents(PetscSpace sp, PetscInt Nc)
435: {
438:   sp->Nc = Nc;
439:   return(0);
440: }

442: /*@
443:   PetscSpaceSetNumVariables - Set the number of variables for this space

445:   Input Parameters:
446: + sp - The PetscSpace
447: - n - The number of variables, e.g. x, y, z...

449:   Level: intermediate

451: .seealso: PetscSpaceGetNumVariables(), PetscSpaceSetNumComponents(), PetscSpaceCreate(), PetscSpace
452: @*/
453: PetscErrorCode PetscSpaceSetNumVariables(PetscSpace sp, PetscInt n)
454: {
457:   sp->Nv = n;
458:   return(0);
459: }

461: /*@
462:   PetscSpaceGetNumVariables - Return the number of variables for this space

464:   Input Parameter:
465: . sp - The PetscSpace

467:   Output Parameter:
468: . Nc - The number of variables, e.g. x, y, z...

470:   Level: intermediate

472: .seealso: PetscSpaceSetNumVariables(), PetscSpaceGetNumComponents(), PetscSpaceGetDimension(), PetscSpaceCreate(), PetscSpace
473: @*/
474: PetscErrorCode PetscSpaceGetNumVariables(PetscSpace sp, PetscInt *n)
475: {
479:   *n = sp->Nv;
480:   return(0);
481: }

483: /*@C
484:   PetscSpaceEvaluate - Evaluate the basis functions and their derivatives (jet) at each point

486:   Input Parameters:
487: + sp      - The PetscSpace
488: . npoints - The number of evaluation points, in reference coordinates
489: - points  - The point coordinates

491:   Output Parameters:
492: + B - The function evaluations in a npoints x nfuncs array
493: . D - The derivative evaluations in a npoints x nfuncs x dim array
494: - H - The second derivative evaluations in a npoints x nfuncs x dim x dim array

496:   Note: Above nfuncs is the dimension of the space, and dim is the spatial dimension. The coordinates are given
497:   on the reference cell, not in real space.

499:   Level: beginner

501: .seealso: PetscFECreateTabulation(), PetscFEGetCellTabulation(), PetscSpaceCreate()
502: @*/
503: PetscErrorCode PetscSpaceEvaluate(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[])
504: {

508:   if (!npoints) return(0);
514:   if (sp->ops->evaluate) {(*sp->ops->evaluate)(sp, npoints, points, B, D, H);}
515:   return(0);
516: }

518: /*@
519:   PetscSpaceGetHeightSubspace - Get the subset of the primal space basis that is supported on a mesh point of a given height.

521:   If the space is not defined on mesh points of the given height (e.g. if the space is discontinuous and
522:   pointwise values are not defined on the element boundaries), or if the implementation of PetscSpace does not
523:   support extracting subspaces, then NULL is returned.

525:   This does not increment the reference count on the returned space, and the user should not destroy it.

527:   Not collective

529:   Input Parameters:
530: + sp - the PetscSpace object
531: - height - the height of the mesh point for which the subspace is desired

533:   Output Parameter:
534: . subsp - the subspace

536:   Level: advanced

538: .seealso: PetscDualSpaceGetHeightSubspace(), PetscSpace
539: @*/
540: PetscErrorCode PetscSpaceGetHeightSubspace(PetscSpace sp, PetscInt height, PetscSpace *subsp)
541: {

547:   *subsp = NULL;
548:   if (sp->ops->getheightsubspace) {
549:     (*sp->ops->getheightsubspace)(sp, height, subsp);
550:   }
551:   return(0);
552: }