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