120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
220cf1dd8SToby Isaac #include <petsc/private/dtimpl.h> /*I "petscdt.h" I*/
320cf1dd8SToby Isaac
PetscSpacePointView_Ascii(PetscSpace sp,PetscViewer viewer)4d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpacePointView_Ascii(PetscSpace sp, PetscViewer viewer)
5d71ae5a4SJacob Faibussowitsch {
620cf1dd8SToby Isaac PetscSpace_Point *pt = (PetscSpace_Point *)sp->data;
720cf1dd8SToby Isaac PetscViewerFormat format;
820cf1dd8SToby Isaac
920cf1dd8SToby Isaac PetscFunctionBegin;
109566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format));
1120cf1dd8SToby Isaac if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
1263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Point space in dimension %" PetscInt_FMT ":\n", sp->Nv));
139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer));
149566063dSJacob Faibussowitsch PetscCall(PetscQuadratureView(pt->quad, viewer));
159566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer));
161baa6e33SBarry Smith } else PetscCall(PetscViewerASCIIPrintf(viewer, "Point space in dimension %" PetscInt_FMT " on %" PetscInt_FMT " points\n", sp->Nv, pt->quad->numPoints));
173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1820cf1dd8SToby Isaac }
1920cf1dd8SToby Isaac
PetscSpaceView_Point(PetscSpace sp,PetscViewer viewer)20d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpaceView_Point(PetscSpace sp, PetscViewer viewer)
21d71ae5a4SJacob Faibussowitsch {
22*9f196a02SMartin Diehl PetscBool isascii;
2320cf1dd8SToby Isaac
2420cf1dd8SToby Isaac PetscFunctionBegin;
2520cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
2620cf1dd8SToby Isaac PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
27*9f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
28*9f196a02SMartin Diehl if (isascii) PetscCall(PetscSpacePointView_Ascii(sp, viewer));
293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3020cf1dd8SToby Isaac }
3120cf1dd8SToby Isaac
PetscSpaceSetUp_Point(PetscSpace sp)32d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpaceSetUp_Point(PetscSpace sp)
33d71ae5a4SJacob Faibussowitsch {
3420cf1dd8SToby Isaac PetscSpace_Point *pt = (PetscSpace_Point *)sp->data;
3520cf1dd8SToby Isaac
3620cf1dd8SToby Isaac PetscFunctionBegin;
3720cf1dd8SToby Isaac if (!pt->quad->points && sp->degree >= 0) {
389566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&pt->quad));
399566063dSJacob Faibussowitsch PetscCall(PetscDTStroudConicalQuadrature(sp->Nv, sp->Nc, PetscMax(sp->degree + 1, 1), -1.0, 1.0, &pt->quad));
4020cf1dd8SToby Isaac }
413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4220cf1dd8SToby Isaac }
4320cf1dd8SToby Isaac
PetscSpaceDestroy_Point(PetscSpace sp)44d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpaceDestroy_Point(PetscSpace sp)
45d71ae5a4SJacob Faibussowitsch {
4620cf1dd8SToby Isaac PetscSpace_Point *pt = (PetscSpace_Point *)sp->data;
4720cf1dd8SToby Isaac
4820cf1dd8SToby Isaac PetscFunctionBegin;
499566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&pt->quad));
509566063dSJacob Faibussowitsch PetscCall(PetscFree(pt));
513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5220cf1dd8SToby Isaac }
5320cf1dd8SToby Isaac
PetscSpaceGetDimension_Point(PetscSpace sp,PetscInt * dim)54d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpaceGetDimension_Point(PetscSpace sp, PetscInt *dim)
55d71ae5a4SJacob Faibussowitsch {
5620cf1dd8SToby Isaac PetscSpace_Point *pt = (PetscSpace_Point *)sp->data;
5720cf1dd8SToby Isaac
5820cf1dd8SToby Isaac PetscFunctionBegin;
5920cf1dd8SToby Isaac *dim = pt->quad->numPoints;
603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
6120cf1dd8SToby Isaac }
6220cf1dd8SToby Isaac
PetscSpaceEvaluate_Point(PetscSpace sp,PetscInt npoints,const PetscReal points[],PetscReal B[],PetscReal D[],PetscReal H[])63d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpaceEvaluate_Point(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[])
64d71ae5a4SJacob Faibussowitsch {
6520cf1dd8SToby Isaac PetscSpace_Point *pt = (PetscSpace_Point *)sp->data;
6620cf1dd8SToby Isaac PetscInt dim = sp->Nv, pdim = pt->quad->numPoints, d, p, i, c;
6720cf1dd8SToby Isaac
6820cf1dd8SToby Isaac PetscFunctionBegin;
6963a3b9bcSJacob Faibussowitsch PetscCheck(npoints == pt->quad->numPoints, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot evaluate Point space on %" PetscInt_FMT " points != %" PetscInt_FMT " size", npoints, pt->quad->numPoints);
709566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(B, npoints * pdim));
7120cf1dd8SToby Isaac for (p = 0; p < npoints; ++p) {
7220cf1dd8SToby Isaac for (i = 0; i < pdim; ++i) {
7320cf1dd8SToby Isaac for (d = 0; d < dim; ++d) {
7420cf1dd8SToby Isaac if (PetscAbsReal(points[p * dim + d] - pt->quad->points[p * dim + d]) > 1.0e-10) break;
7520cf1dd8SToby Isaac }
769371c9d4SSatish Balay if (d >= dim) {
779371c9d4SSatish Balay B[p * pdim + i] = 1.0;
789371c9d4SSatish Balay break;
799371c9d4SSatish Balay }
8020cf1dd8SToby Isaac }
8120cf1dd8SToby Isaac }
8220cf1dd8SToby Isaac /* Replicate for other components */
8320cf1dd8SToby Isaac for (c = 1; c < sp->Nc; ++c) {
8420cf1dd8SToby Isaac for (p = 0; p < npoints; ++p) {
85ad540459SPierre Jolivet for (i = 0; i < pdim; ++i) B[(c * npoints + p) * pdim + i] = B[p * pdim + i];
8620cf1dd8SToby Isaac }
8720cf1dd8SToby Isaac }
889566063dSJacob Faibussowitsch if (D) PetscCall(PetscArrayzero(D, npoints * pdim * dim));
899566063dSJacob Faibussowitsch if (H) PetscCall(PetscArrayzero(H, npoints * pdim * dim * dim));
903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
9120cf1dd8SToby Isaac }
9220cf1dd8SToby Isaac
PetscSpaceInitialize_Point(PetscSpace sp)93d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpaceInitialize_Point(PetscSpace sp)
94d71ae5a4SJacob Faibussowitsch {
9520cf1dd8SToby Isaac PetscFunctionBegin;
9620cf1dd8SToby Isaac sp->ops->setfromoptions = NULL;
9720cf1dd8SToby Isaac sp->ops->setup = PetscSpaceSetUp_Point;
9820cf1dd8SToby Isaac sp->ops->view = PetscSpaceView_Point;
9920cf1dd8SToby Isaac sp->ops->destroy = PetscSpaceDestroy_Point;
10020cf1dd8SToby Isaac sp->ops->getdimension = PetscSpaceGetDimension_Point;
10120cf1dd8SToby Isaac sp->ops->evaluate = PetscSpaceEvaluate_Point;
1023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
10320cf1dd8SToby Isaac }
10420cf1dd8SToby Isaac
10520cf1dd8SToby Isaac /*MC
106dce8aebaSBarry Smith PETSCSPACEPOINT = "point" - A `PetscSpace` object that encapsulates functions defined on a set of quadrature points.
10720cf1dd8SToby Isaac
10820cf1dd8SToby Isaac Level: intermediate
10920cf1dd8SToby Isaac
110dce8aebaSBarry Smith .seealso: `PetscSpace`, `PetscSpaceType`, `PetscSpaceCreate()`, `PetscSpaceSetType()`
11120cf1dd8SToby Isaac M*/
11220cf1dd8SToby Isaac
PetscSpaceCreate_Point(PetscSpace sp)113d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscSpaceCreate_Point(PetscSpace sp)
114d71ae5a4SJacob Faibussowitsch {
11520cf1dd8SToby Isaac PetscSpace_Point *pt;
11620cf1dd8SToby Isaac
11720cf1dd8SToby Isaac PetscFunctionBegin;
11820cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
1194dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&pt));
12020cf1dd8SToby Isaac sp->data = pt;
12120cf1dd8SToby Isaac
12220cf1dd8SToby Isaac sp->Nv = 0;
1231690c2aeSBarry Smith sp->maxDegree = PETSC_INT_MAX;
1249566063dSJacob Faibussowitsch PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &pt->quad));
1259566063dSJacob Faibussowitsch PetscCall(PetscQuadratureSetData(pt->quad, 0, 1, 0, NULL, NULL));
12620cf1dd8SToby Isaac
1279566063dSJacob Faibussowitsch PetscCall(PetscSpaceInitialize_Point(sp));
1283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
12920cf1dd8SToby Isaac }
13020cf1dd8SToby Isaac
13120cf1dd8SToby Isaac /*@
13220cf1dd8SToby Isaac PetscSpacePointSetPoints - Sets the evaluation points for the space to coincide with the points of a quadrature rule
13320cf1dd8SToby Isaac
13420f4b53cSBarry Smith Logically Collective
13520cf1dd8SToby Isaac
13620cf1dd8SToby Isaac Input Parameters:
137dce8aebaSBarry Smith + sp - The `PetscSpace`
138dce8aebaSBarry Smith - q - The `PetscQuadrature` defining the points
13920cf1dd8SToby Isaac
14020cf1dd8SToby Isaac Level: intermediate
14120cf1dd8SToby Isaac
142dce8aebaSBarry Smith .seealso: `PetscSpace`, `PetscQuadrature`, `PetscSpaceCreate()`, `PetscSpaceSetType()`
14320cf1dd8SToby Isaac @*/
PetscSpacePointSetPoints(PetscSpace sp,PetscQuadrature q)144d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSpacePointSetPoints(PetscSpace sp, PetscQuadrature q)
145d71ae5a4SJacob Faibussowitsch {
14620cf1dd8SToby Isaac PetscSpace_Point *pt = (PetscSpace_Point *)sp->data;
14720cf1dd8SToby Isaac
14820cf1dd8SToby Isaac PetscFunctionBegin;
14920cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
150064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(q, PETSCQUADRATURE_CLASSID, 2);
1519566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&pt->quad));
1529566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDuplicate(q, &pt->quad));
1533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
15420cf1dd8SToby Isaac }
15520cf1dd8SToby Isaac
15620cf1dd8SToby Isaac /*@
15720cf1dd8SToby Isaac PetscSpacePointGetPoints - Gets the evaluation points for the space as the points of a quadrature rule
15820cf1dd8SToby Isaac
15920f4b53cSBarry Smith Logically Collective
16020cf1dd8SToby Isaac
16120cf1dd8SToby Isaac Input Parameter:
162dce8aebaSBarry Smith . sp - The `PetscSpace`
16320cf1dd8SToby Isaac
16420cf1dd8SToby Isaac Output Parameter:
165dce8aebaSBarry Smith . q - The `PetscQuadrature` defining the points
16620cf1dd8SToby Isaac
16720cf1dd8SToby Isaac Level: intermediate
16820cf1dd8SToby Isaac
169dce8aebaSBarry Smith .seealso: `PetscSpace`, `PetscQuadrature`, `PetscSpaceCreate()`, `PetscSpaceSetType()`
17020cf1dd8SToby Isaac @*/
PetscSpacePointGetPoints(PetscSpace sp,PetscQuadrature * q)171d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSpacePointGetPoints(PetscSpace sp, PetscQuadrature *q)
172d71ae5a4SJacob Faibussowitsch {
17320cf1dd8SToby Isaac PetscSpace_Point *pt = (PetscSpace_Point *)sp->data;
17420cf1dd8SToby Isaac
17520cf1dd8SToby Isaac PetscFunctionBegin;
17620cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
1774f572ea9SToby Isaac PetscAssertPointer(q, 2);
17820cf1dd8SToby Isaac *q = pt->quad;
1793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
18020cf1dd8SToby Isaac }
181