xref: /petsc/src/dm/dt/space/impls/point/spacepoint.c (revision 20cf1dd8cd656e27510885a71ea4be028bf48813)
1*20cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
2*20cf1dd8SToby Isaac #include <petsc/private/dtimpl.h> /*I "petscdt.h" I*/
3*20cf1dd8SToby Isaac 
4*20cf1dd8SToby Isaac PetscErrorCode PetscSpacePointView_Ascii(PetscSpace sp, PetscViewer viewer)
5*20cf1dd8SToby Isaac {
6*20cf1dd8SToby Isaac   PetscSpace_Point *pt = (PetscSpace_Point *) sp->data;
7*20cf1dd8SToby Isaac   PetscViewerFormat format;
8*20cf1dd8SToby Isaac   PetscErrorCode    ierr;
9*20cf1dd8SToby Isaac 
10*20cf1dd8SToby Isaac   PetscFunctionBegin;
11*20cf1dd8SToby Isaac   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
12*20cf1dd8SToby Isaac   if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
13*20cf1dd8SToby Isaac     ierr = PetscViewerASCIIPrintf(viewer, "Point space in dimension %d:\n", sp->Nv);CHKERRQ(ierr);
14*20cf1dd8SToby Isaac     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
15*20cf1dd8SToby Isaac     ierr = PetscQuadratureView(pt->quad, viewer);CHKERRQ(ierr);
16*20cf1dd8SToby Isaac     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
17*20cf1dd8SToby Isaac   } else {
18*20cf1dd8SToby Isaac     ierr = PetscViewerASCIIPrintf(viewer, "Point space in dimension %d on %d points\n", sp->Nv, pt->quad->numPoints);CHKERRQ(ierr);
19*20cf1dd8SToby Isaac   }
20*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
21*20cf1dd8SToby Isaac }
22*20cf1dd8SToby Isaac 
23*20cf1dd8SToby Isaac PetscErrorCode PetscSpaceView_Point(PetscSpace sp, PetscViewer viewer)
24*20cf1dd8SToby Isaac {
25*20cf1dd8SToby Isaac   PetscBool      iascii;
26*20cf1dd8SToby Isaac   PetscErrorCode ierr;
27*20cf1dd8SToby Isaac 
28*20cf1dd8SToby Isaac   PetscFunctionBegin;
29*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
30*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
31*20cf1dd8SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
32*20cf1dd8SToby Isaac   if (iascii) {ierr = PetscSpacePointView_Ascii(sp, viewer);CHKERRQ(ierr);}
33*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
34*20cf1dd8SToby Isaac }
35*20cf1dd8SToby Isaac 
36*20cf1dd8SToby Isaac PetscErrorCode PetscSpaceSetUp_Point(PetscSpace sp)
37*20cf1dd8SToby Isaac {
38*20cf1dd8SToby Isaac   PetscSpace_Point *pt = (PetscSpace_Point *) sp->data;
39*20cf1dd8SToby Isaac   PetscErrorCode    ierr;
40*20cf1dd8SToby Isaac 
41*20cf1dd8SToby Isaac   PetscFunctionBegin;
42*20cf1dd8SToby Isaac   if (!pt->quad->points && sp->degree >= 0) {
43*20cf1dd8SToby Isaac     ierr = PetscQuadratureDestroy(&pt->quad);CHKERRQ(ierr);
44*20cf1dd8SToby Isaac     ierr = PetscDTGaussJacobiQuadrature(sp->Nv, sp->Nc, PetscMax(sp->degree + 1, 1), -1.0, 1.0, &pt->quad);CHKERRQ(ierr);
45*20cf1dd8SToby Isaac   }
46*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
47*20cf1dd8SToby Isaac }
48*20cf1dd8SToby Isaac 
49*20cf1dd8SToby Isaac PetscErrorCode PetscSpaceDestroy_Point(PetscSpace sp)
50*20cf1dd8SToby Isaac {
51*20cf1dd8SToby Isaac   PetscSpace_Point *pt = (PetscSpace_Point *) sp->data;
52*20cf1dd8SToby Isaac   PetscErrorCode    ierr;
53*20cf1dd8SToby Isaac 
54*20cf1dd8SToby Isaac   PetscFunctionBegin;
55*20cf1dd8SToby Isaac   ierr = PetscQuadratureDestroy(&pt->quad);CHKERRQ(ierr);
56*20cf1dd8SToby Isaac   ierr = PetscFree(pt);CHKERRQ(ierr);
57*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
58*20cf1dd8SToby Isaac }
59*20cf1dd8SToby Isaac 
60*20cf1dd8SToby Isaac PetscErrorCode PetscSpaceGetDimension_Point(PetscSpace sp, PetscInt *dim)
61*20cf1dd8SToby Isaac {
62*20cf1dd8SToby Isaac   PetscSpace_Point *pt = (PetscSpace_Point *) sp->data;
63*20cf1dd8SToby Isaac 
64*20cf1dd8SToby Isaac   PetscFunctionBegin;
65*20cf1dd8SToby Isaac   *dim = pt->quad->numPoints;
66*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
67*20cf1dd8SToby Isaac }
68*20cf1dd8SToby Isaac 
69*20cf1dd8SToby Isaac PetscErrorCode PetscSpaceEvaluate_Point(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[])
70*20cf1dd8SToby Isaac {
71*20cf1dd8SToby Isaac   PetscSpace_Point *pt  = (PetscSpace_Point *) sp->data;
72*20cf1dd8SToby Isaac   PetscInt          dim = sp->Nv, pdim = pt->quad->numPoints, d, p, i, c;
73*20cf1dd8SToby Isaac   PetscErrorCode    ierr;
74*20cf1dd8SToby Isaac 
75*20cf1dd8SToby Isaac   PetscFunctionBegin;
76*20cf1dd8SToby Isaac   if (npoints != pt->quad->numPoints) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot evaluate Point space on %d points != %d size", npoints, pt->quad->numPoints);
77*20cf1dd8SToby Isaac   ierr = PetscMemzero(B, npoints*pdim * sizeof(PetscReal));CHKERRQ(ierr);
78*20cf1dd8SToby Isaac   for (p = 0; p < npoints; ++p) {
79*20cf1dd8SToby Isaac     for (i = 0; i < pdim; ++i) {
80*20cf1dd8SToby Isaac       for (d = 0; d < dim; ++d) {
81*20cf1dd8SToby Isaac         if (PetscAbsReal(points[p*dim+d] - pt->quad->points[p*dim+d]) > 1.0e-10) break;
82*20cf1dd8SToby Isaac       }
83*20cf1dd8SToby Isaac       if (d >= dim) {B[p*pdim+i] = 1.0; break;}
84*20cf1dd8SToby Isaac     }
85*20cf1dd8SToby Isaac   }
86*20cf1dd8SToby Isaac   /* Replicate for other components */
87*20cf1dd8SToby Isaac   for (c = 1; c < sp->Nc; ++c) {
88*20cf1dd8SToby Isaac     for (p = 0; p < npoints; ++p) {
89*20cf1dd8SToby Isaac       for (i = 0; i < pdim; ++i) {
90*20cf1dd8SToby Isaac         B[(c*npoints + p)*pdim + i] = B[p*pdim + i];
91*20cf1dd8SToby Isaac       }
92*20cf1dd8SToby Isaac     }
93*20cf1dd8SToby Isaac   }
94*20cf1dd8SToby Isaac   if (D) {ierr = PetscMemzero(D, npoints*pdim*dim * sizeof(PetscReal));CHKERRQ(ierr);}
95*20cf1dd8SToby Isaac   if (H) {ierr = PetscMemzero(H, npoints*pdim*dim*dim * sizeof(PetscReal));CHKERRQ(ierr);}
96*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
97*20cf1dd8SToby Isaac }
98*20cf1dd8SToby Isaac 
99*20cf1dd8SToby Isaac PetscErrorCode PetscSpaceInitialize_Point(PetscSpace sp)
100*20cf1dd8SToby Isaac {
101*20cf1dd8SToby Isaac   PetscFunctionBegin;
102*20cf1dd8SToby Isaac   sp->ops->setfromoptions = NULL;
103*20cf1dd8SToby Isaac   sp->ops->setup          = PetscSpaceSetUp_Point;
104*20cf1dd8SToby Isaac   sp->ops->view           = PetscSpaceView_Point;
105*20cf1dd8SToby Isaac   sp->ops->destroy        = PetscSpaceDestroy_Point;
106*20cf1dd8SToby Isaac   sp->ops->getdimension   = PetscSpaceGetDimension_Point;
107*20cf1dd8SToby Isaac   sp->ops->evaluate       = PetscSpaceEvaluate_Point;
108*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
109*20cf1dd8SToby Isaac }
110*20cf1dd8SToby Isaac 
111*20cf1dd8SToby Isaac /*MC
112*20cf1dd8SToby Isaac   PETSCSPACEPOINT = "point" - A PetscSpace object that encapsulates functions defined on a set of quadrature points.
113*20cf1dd8SToby Isaac 
114*20cf1dd8SToby Isaac   Level: intermediate
115*20cf1dd8SToby Isaac 
116*20cf1dd8SToby Isaac .seealso: PetscSpaceType, PetscSpaceCreate(), PetscSpaceSetType()
117*20cf1dd8SToby Isaac M*/
118*20cf1dd8SToby Isaac 
119*20cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscSpaceCreate_Point(PetscSpace sp)
120*20cf1dd8SToby Isaac {
121*20cf1dd8SToby Isaac   PetscSpace_Point *pt;
122*20cf1dd8SToby Isaac   PetscErrorCode    ierr;
123*20cf1dd8SToby Isaac 
124*20cf1dd8SToby Isaac   PetscFunctionBegin;
125*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
126*20cf1dd8SToby Isaac   ierr     = PetscNewLog(sp,&pt);CHKERRQ(ierr);
127*20cf1dd8SToby Isaac   sp->data = pt;
128*20cf1dd8SToby Isaac 
129*20cf1dd8SToby Isaac   sp->Nv = 0;
130*20cf1dd8SToby Isaac   sp->maxDegree = PETSC_MAX_INT;
131*20cf1dd8SToby Isaac   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &pt->quad);CHKERRQ(ierr);
132*20cf1dd8SToby Isaac   ierr = PetscQuadratureSetData(pt->quad, 0, 1, 0, NULL, NULL);CHKERRQ(ierr);
133*20cf1dd8SToby Isaac 
134*20cf1dd8SToby Isaac   ierr = PetscSpaceInitialize_Point(sp);CHKERRQ(ierr);
135*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
136*20cf1dd8SToby Isaac }
137*20cf1dd8SToby Isaac 
138*20cf1dd8SToby Isaac /*@
139*20cf1dd8SToby Isaac   PetscSpacePointSetPoints - Sets the evaluation points for the space to coincide with the points of a quadrature rule
140*20cf1dd8SToby Isaac 
141*20cf1dd8SToby Isaac   Logically collective
142*20cf1dd8SToby Isaac 
143*20cf1dd8SToby Isaac   Input Parameters:
144*20cf1dd8SToby Isaac + sp - The PetscSpace
145*20cf1dd8SToby Isaac - q  - The PetscQuadrature defining the points
146*20cf1dd8SToby Isaac 
147*20cf1dd8SToby Isaac   Level: intermediate
148*20cf1dd8SToby Isaac 
149*20cf1dd8SToby Isaac .keywords: PetscSpacePoint
150*20cf1dd8SToby Isaac .seealso: PetscSpaceCreate(), PetscSpaceSetType()
151*20cf1dd8SToby Isaac @*/
152*20cf1dd8SToby Isaac PetscErrorCode PetscSpacePointSetPoints(PetscSpace sp, PetscQuadrature q)
153*20cf1dd8SToby Isaac {
154*20cf1dd8SToby Isaac   PetscSpace_Point *pt = (PetscSpace_Point *) sp->data;
155*20cf1dd8SToby Isaac   PetscErrorCode    ierr;
156*20cf1dd8SToby Isaac 
157*20cf1dd8SToby Isaac   PetscFunctionBegin;
158*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
159*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 2);
160*20cf1dd8SToby Isaac   ierr = PetscQuadratureDestroy(&pt->quad);CHKERRQ(ierr);
161*20cf1dd8SToby Isaac   ierr = PetscQuadratureDuplicate(q, &pt->quad);CHKERRQ(ierr);
162*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
163*20cf1dd8SToby Isaac }
164*20cf1dd8SToby Isaac 
165*20cf1dd8SToby Isaac /*@
166*20cf1dd8SToby Isaac   PetscSpacePointGetPoints - Gets the evaluation points for the space as the points of a quadrature rule
167*20cf1dd8SToby Isaac 
168*20cf1dd8SToby Isaac   Logically collective
169*20cf1dd8SToby Isaac 
170*20cf1dd8SToby Isaac   Input Parameter:
171*20cf1dd8SToby Isaac . sp - The PetscSpace
172*20cf1dd8SToby Isaac 
173*20cf1dd8SToby Isaac   Output Parameter:
174*20cf1dd8SToby Isaac . q  - The PetscQuadrature defining the points
175*20cf1dd8SToby Isaac 
176*20cf1dd8SToby Isaac   Level: intermediate
177*20cf1dd8SToby Isaac 
178*20cf1dd8SToby Isaac .keywords: PetscSpacePoint
179*20cf1dd8SToby Isaac .seealso: PetscSpaceCreate(), PetscSpaceSetType()
180*20cf1dd8SToby Isaac @*/
181*20cf1dd8SToby Isaac PetscErrorCode PetscSpacePointGetPoints(PetscSpace sp, PetscQuadrature *q)
182*20cf1dd8SToby Isaac {
183*20cf1dd8SToby Isaac   PetscSpace_Point *pt = (PetscSpace_Point *) sp->data;
184*20cf1dd8SToby Isaac 
185*20cf1dd8SToby Isaac   PetscFunctionBegin;
186*20cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
187*20cf1dd8SToby Isaac   PetscValidPointer(q, 2);
188*20cf1dd8SToby Isaac   *q = pt->quad;
189*20cf1dd8SToby Isaac   PetscFunctionReturn(0);
190*20cf1dd8SToby Isaac }
191*20cf1dd8SToby Isaac 
192