Actual source code: feceed.c

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

  3: #ifdef PETSC_HAVE_LIBCEED
  4: #include <petscfeceed.h>

  6: /*@C
  7:   PetscFESetCeed - Set the `Ceed` object to a `PetscFE`

  9:   Not Collective

 11:   Input Parameters:
 12: + fe   - The `PetscFE`
 13: - ceed - The `Ceed` object

 15:   Level: intermediate

 17: .seealso: `PetscFE`, `PetscFEGetCeedBasis()`, `DMGetCeed()`
 18: @*/
 19: PetscErrorCode PetscFESetCeed(PetscFE fe, Ceed ceed)
 20: {
 21:   PetscFunctionBegin;
 23:   if (fe->ceed == ceed) PetscFunctionReturn(PETSC_SUCCESS);
 24:   PetscCallCEED(CeedReferenceCopy(ceed, &fe->ceed));
 25:   PetscFunctionReturn(PETSC_SUCCESS);
 26: }

 28: /*@C
 29:   PetscFEGetCeedBasis - Get the `Ceed` object mirroring this `PetscFE`

 31:   Not Collective

 33:   Input Parameter:
 34: . fe - The `PetscFE`

 36:   Output Parameter:
 37: . basis - The `CeedBasis`

 39:   Level: intermediate

 41:   Note:
 42:   This is a borrowed reference, so it is not freed.

 44: .seealso: `PetscFE`, `PetscFESetCeed()`, `DMGetCeed()`
 45: @*/
 46: PetscErrorCode PetscFEGetCeedBasis(PetscFE fe, CeedBasis *basis)
 47: {
 48:   PetscSpace      sp;
 49:   PetscQuadrature q;
 50:   PetscInt        dim, Nc, deg, ord;

 52:   PetscFunctionBegin;
 54:   PetscAssertPointer(basis, 2);
 55:   if (!fe->ceedBasis && fe->ceed) {
 56:     PetscCall(PetscFEGetSpatialDimension(fe, &dim));
 57:     PetscCall(PetscFEGetNumComponents(fe, &Nc));
 58:     PetscCall(PetscFEGetBasisSpace(fe, &sp));
 59:     PetscCall(PetscSpaceGetDegree(sp, &deg, NULL));
 60:     PetscCall(PetscFEGetQuadrature(fe, &q));
 61:     PetscCall(PetscQuadratureGetOrder(q, &ord));
 62:     PetscCallCEED(CeedBasisCreateTensorH1Lagrange(fe->ceed, dim, Nc, deg + 1, (ord + 1) / 2, CEED_GAUSS, &fe->ceedBasis));
 63:   }
 64:   *basis = fe->ceedBasis;
 65:   PetscFunctionReturn(PETSC_SUCCESS);
 66: }

 68: #endif