xref: /petsc/src/dm/dt/fe/interface/ceed/feceed.c (revision 5a236de668401f703e133ea891e0305529ed6f83)
1 #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
2 
3 #include <petscfeceed.h>
4 
5 /*@C
6   PetscFESetCeed - Set the `Ceed` object to a `PetscFE`
7 
8   Not Collective
9 
10   Input Parameters:
11 + fe   - The `PetscFE`
12 - ceed - The `Ceed` object
13 
14   Level: intermediate
15 
16 .seealso: `PetscFE`, `PetscFEGetCeedBasis()`, `DMGetCeed()`
17 @*/
PetscFESetCeed(PetscFE fe,Ceed ceed)18 PetscErrorCode PetscFESetCeed(PetscFE fe, Ceed ceed)
19 {
20   PetscFunctionBegin;
21   PetscValidHeaderSpecific(fe, PETSCFE_CLASSID, 1);
22   if (fe->ceed == ceed) PetscFunctionReturn(PETSC_SUCCESS);
23   PetscCallCEED(CeedReferenceCopy(ceed, &fe->ceed));
24   PetscFunctionReturn(PETSC_SUCCESS);
25 }
26 
27 /*@C
28   PetscFEGetCeedBasis - Get the `Ceed` object mirroring this `PetscFE`
29 
30   Not Collective
31 
32   Input Parameter:
33 . fe - The `PetscFE`
34 
35   Output Parameter:
36 . basis - The `CeedBasis`
37 
38   Level: intermediate
39 
40   Note:
41   This is a borrowed reference, so it is not freed.
42 
43 .seealso: `PetscFE`, `PetscFESetCeed()`, `DMGetCeed()`
44 @*/
PetscFEGetCeedBasis(PetscFE fe,CeedBasis * basis)45 PetscErrorCode PetscFEGetCeedBasis(PetscFE fe, CeedBasis *basis)
46 {
47   PetscSpace      sp;
48   PetscQuadrature q;
49   PetscInt        dim, Nc, deg, ord;
50 
51   PetscFunctionBegin;
52   PetscValidHeaderSpecific(fe, PETSCFE_CLASSID, 1);
53   PetscAssertPointer(basis, 2);
54   if (!fe->ceedBasis && fe->ceed) {
55     PetscCall(PetscFEGetSpatialDimension(fe, &dim));
56     PetscCall(PetscFEGetNumComponents(fe, &Nc));
57     PetscCall(PetscFEGetBasisSpace(fe, &sp));
58     PetscCall(PetscSpaceGetDegree(sp, &deg, NULL));
59     PetscCall(PetscFEGetQuadrature(fe, &q));
60     PetscCall(PetscQuadratureGetOrder(q, &ord));
61     PetscCallCEED(CeedBasisCreateTensorH1Lagrange(fe->ceed, dim, Nc, deg + 1, (ord + 1) / 2, CEED_GAUSS, &fe->ceedBasis));
62   }
63   *basis = fe->ceedBasis;
64   PetscFunctionReturn(PETSC_SUCCESS);
65 }
66