15ffaa15fSToby Isaac const char help[] = "Construct and set a Lagrange dual space from options, then view it to\n"
25ffaa15fSToby Isaac "understand the effects of different parameters.";
35ffaa15fSToby Isaac
45ffaa15fSToby Isaac #include <petscfe.h>
55ffaa15fSToby Isaac #include <petscdmplex.h>
65ffaa15fSToby Isaac
main(int argc,char ** argv)7*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
8*d71ae5a4SJacob Faibussowitsch {
95ffaa15fSToby Isaac PetscInt dim;
105ffaa15fSToby Isaac PetscBool tensorCell;
115ffaa15fSToby Isaac DM K;
125ffaa15fSToby Isaac PetscDualSpace dsp;
135ffaa15fSToby Isaac
14327415f7SBarry Smith PetscFunctionBeginUser;
159566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
165ffaa15fSToby Isaac dim = 2;
175ffaa15fSToby Isaac tensorCell = PETSC_FALSE;
18d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, "", "Options for PETSCDUALSPACELAGRANGE test", "none");
199566063dSJacob Faibussowitsch PetscCall(PetscOptionsRangeInt("-dim", "The spatial dimension", "ex1.c", dim, &dim, NULL, 0, 3));
209566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-tensor", "Whether the cell is a tensor product cell or a simplex", "ex1.c", tensorCell, &tensorCell, NULL));
21d0609cedSBarry Smith PetscOptionsEnd();
225ffaa15fSToby Isaac
239566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceCreate(PETSC_COMM_WORLD, &dsp));
249566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)dsp, "Lagrange dual space"));
259566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetType(dsp, PETSCDUALSPACELAGRANGE));
265ffaa15fSToby Isaac /* While Lagrange nodes don't require the existence of a reference cell to
275ffaa15fSToby Isaac * be refined, when we construct finite element dual spaces we have to be
285ffaa15fSToby Isaac * careful about what kind of continuity is maintained when cells are glued
295ffaa15fSToby Isaac * together to make a mesh. The PetscDualSpace object is responsible for
305ffaa15fSToby Isaac * conveying continuity requirements to a finite element assembly routines,
315ffaa15fSToby Isaac * so a PetscDualSpace needs a reference element: a single element mesh,
325ffaa15fSToby Isaac * whose boundary points are the interstitial points in a mesh */
339566063dSJacob Faibussowitsch PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_WORLD, DMPolytopeTypeSimpleShape(dim, (PetscBool)!tensorCell), &K));
349566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetDM(dsp, K));
355ffaa15fSToby Isaac /* This gives us the opportunity to change the parameters of the dual space
365ffaa15fSToby Isaac * from the command line, as we do in the tests below. When
375ffaa15fSToby Isaac * PetscDualSpaceSetFromOptions() is called, it also enables other optional
385ffaa15fSToby Isaac * behavior (see the next step) */
399566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetFromOptions(dsp));
405ffaa15fSToby Isaac /* This step parses the parameters of the dual space into
415ffaa15fSToby Isaac * sets of functionals that are assigned to each of the mesh points in K.
425ffaa15fSToby Isaac *
435ffaa15fSToby Isaac * The functionals can be accessed individually by
445ffaa15fSToby Isaac * PetscDualSpaceGetFunctional(), or more efficiently all at once by
455ffaa15fSToby Isaac * PetscDualSpaceGetAllData(), which returns a set of quadrature points
465ffaa15fSToby Isaac * at which to evaluate a function, and a matrix that takes those
475ffaa15fSToby Isaac * evaluations and turns them into the evaluation of the dual space's
485ffaa15fSToby Isaac * functionals on the function.
495ffaa15fSToby Isaac *
505ffaa15fSToby Isaac * (TODO: tutorial for PetscDualSpaceGetAllData() and
515ffaa15fSToby Isaac * PetscDualSpaceGetInteriorData().)
525ffaa15fSToby Isaac *
535ffaa15fSToby Isaac * Because we called PetscDualSpaceSetFromOptions(), we have the opportunity
545ffaa15fSToby Isaac * to inspect the results of PetscDualSpaceSetUp() from the command line
555ffaa15fSToby Isaac * with "-petscdualspace_view", followed by an optional description of how
565ffaa15fSToby Isaac * we would like to see the dual space (see examples in the tests below).
575ffaa15fSToby Isaac * */
589566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetUp(dsp));
599566063dSJacob Faibussowitsch PetscCall(DMDestroy(&K));
609566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceDestroy(&dsp));
619566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
62b122ec5aSJacob Faibussowitsch return 0;
635ffaa15fSToby Isaac }
645ffaa15fSToby Isaac
655ffaa15fSToby Isaac /*TEST
665ffaa15fSToby Isaac
675ffaa15fSToby Isaac # quadratic nodes on the triangle
685ffaa15fSToby Isaac test:
695ffaa15fSToby Isaac suffix: 0
706ff15688SToby Isaac filter: sed -E "s/\(\+0, \+0\)/(+0., +0.)/g"
715ffaa15fSToby Isaac args: -dim 2 -tensor 0 -petscdualspace_order 2 -petscdualspace_view ascii::ascii_info_detail
725ffaa15fSToby Isaac
735ffaa15fSToby Isaac # linear nodes on the quadrilateral
745ffaa15fSToby Isaac test:
755ffaa15fSToby Isaac suffix: 1
766ff15688SToby Isaac args: -dim 2 -tensor 1 -petscdualspace_order 1 -petscdualspace_lagrange_tensor 1 -petscdualspace_view ascii::ascii_info_detail
775ffaa15fSToby Isaac
785ffaa15fSToby Isaac # lowest order Raviart-Thomas / Nedelec edge nodes on the hexahedron
795ffaa15fSToby Isaac test:
805ffaa15fSToby Isaac suffix: 2
815ffaa15fSToby Isaac args: -dim 3 -tensor 1 -petscdualspace_order 1 -petscdualspace_components 3 -petscdualspace_form_degree 1 -petscdualspace_lagrange_trimmed 1 -petscdualspace_lagrange_tensor 1 -petscdualspace_view ascii::ascii_info_detail
825ffaa15fSToby Isaac
835ffaa15fSToby Isaac # first order Nedelec second type face nodes on the tetrahedron
845ffaa15fSToby Isaac test:
855ffaa15fSToby Isaac suffix: 3
865ffaa15fSToby Isaac args: -dim 3 -tensor 0 -petscdualspace_order 1 -petscdualspace_components 3 -petscdualspace_form_degree -2 -petscdualspace_view ascii::ascii_info_detail
875ffaa15fSToby Isaac
885ffaa15fSToby Isaac ## Comparing different node types
895ffaa15fSToby Isaac
905ffaa15fSToby Isaac test:
915ffaa15fSToby Isaac suffix: 4
925ffaa15fSToby Isaac args: -dim 2 -tensor 0 -petscdualspace_order 3 -petscdualspace_lagrange_continuity 0 -petscdualspace_lagrange_node_type equispaced -petscdualspace_lagrange_node_endpoints 0 -petscdualspace_view ascii::ascii_info_detail
935ffaa15fSToby Isaac
945ffaa15fSToby Isaac test:
955ffaa15fSToby Isaac suffix: 5
965ffaa15fSToby Isaac args: -dim 2 -tensor 0 -petscdualspace_order 3 -petscdualspace_lagrange_continuity 0 -petscdualspace_lagrange_node_type equispaced -petscdualspace_lagrange_node_endpoints 1 -petscdualspace_view ascii::ascii_info_detail
975ffaa15fSToby Isaac
985ffaa15fSToby Isaac test:
995ffaa15fSToby Isaac suffix: 6
1005ffaa15fSToby Isaac args: -dim 2 -tensor 0 -petscdualspace_order 3 -petscdualspace_lagrange_continuity 0 -petscdualspace_lagrange_node_type gaussjacobi -petscdualspace_lagrange_node_endpoints 0 -petscdualspace_view ascii::ascii_info_detail
1015ffaa15fSToby Isaac
1025ffaa15fSToby Isaac test:
1035ffaa15fSToby Isaac suffix: 7
1045ffaa15fSToby Isaac args: -dim 2 -tensor 0 -petscdualspace_order 3 -petscdualspace_lagrange_continuity 0 -petscdualspace_lagrange_node_type gaussjacobi -petscdualspace_lagrange_node_endpoints 1 -petscdualspace_view ascii::ascii_info_detail
1055ffaa15fSToby Isaac
1065ffaa15fSToby Isaac TEST*/
107