xref: /petsc/src/dm/dt/space/impls/point/spacepoint.c (revision 29b5c1159092c0c4e109657cced02aef37644e8e)
120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
220cf1dd8SToby Isaac #include <petsc/private/dtimpl.h> /*I "petscdt.h" I*/
320cf1dd8SToby Isaac 
4*29b5c115SMatthew G. Knepley static PetscErrorCode PetscSpacePointView_Ascii(PetscSpace sp, PetscViewer viewer)
520cf1dd8SToby Isaac {
620cf1dd8SToby Isaac   PetscSpace_Point *pt = (PetscSpace_Point *) sp->data;
720cf1dd8SToby Isaac   PetscViewerFormat format;
820cf1dd8SToby Isaac   PetscErrorCode    ierr;
920cf1dd8SToby Isaac 
1020cf1dd8SToby Isaac   PetscFunctionBegin;
1120cf1dd8SToby Isaac   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
1220cf1dd8SToby Isaac   if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
1320cf1dd8SToby Isaac     ierr = PetscViewerASCIIPrintf(viewer, "Point space in dimension %d:\n", sp->Nv);CHKERRQ(ierr);
1420cf1dd8SToby Isaac     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1520cf1dd8SToby Isaac     ierr = PetscQuadratureView(pt->quad, viewer);CHKERRQ(ierr);
1620cf1dd8SToby Isaac     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
1720cf1dd8SToby Isaac   } else {
1820cf1dd8SToby Isaac     ierr = PetscViewerASCIIPrintf(viewer, "Point space in dimension %d on %d points\n", sp->Nv, pt->quad->numPoints);CHKERRQ(ierr);
1920cf1dd8SToby Isaac   }
2020cf1dd8SToby Isaac   PetscFunctionReturn(0);
2120cf1dd8SToby Isaac }
2220cf1dd8SToby Isaac 
23*29b5c115SMatthew G. Knepley static PetscErrorCode PetscSpaceView_Point(PetscSpace sp, PetscViewer viewer)
2420cf1dd8SToby Isaac {
2520cf1dd8SToby Isaac   PetscBool      iascii;
2620cf1dd8SToby Isaac   PetscErrorCode ierr;
2720cf1dd8SToby Isaac 
2820cf1dd8SToby Isaac   PetscFunctionBegin;
2920cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
3020cf1dd8SToby Isaac   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
3120cf1dd8SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
3220cf1dd8SToby Isaac   if (iascii) {ierr = PetscSpacePointView_Ascii(sp, viewer);CHKERRQ(ierr);}
3320cf1dd8SToby Isaac   PetscFunctionReturn(0);
3420cf1dd8SToby Isaac }
3520cf1dd8SToby Isaac 
36*29b5c115SMatthew G. Knepley static PetscErrorCode PetscSpaceSetUp_Point(PetscSpace sp)
3720cf1dd8SToby Isaac {
3820cf1dd8SToby Isaac   PetscSpace_Point *pt = (PetscSpace_Point *) sp->data;
3920cf1dd8SToby Isaac   PetscErrorCode    ierr;
4020cf1dd8SToby Isaac 
4120cf1dd8SToby Isaac   PetscFunctionBegin;
4220cf1dd8SToby Isaac   if (!pt->quad->points && sp->degree >= 0) {
4320cf1dd8SToby Isaac     ierr = PetscQuadratureDestroy(&pt->quad);CHKERRQ(ierr);
4420cf1dd8SToby Isaac     ierr = PetscDTGaussJacobiQuadrature(sp->Nv, sp->Nc, PetscMax(sp->degree + 1, 1), -1.0, 1.0, &pt->quad);CHKERRQ(ierr);
4520cf1dd8SToby Isaac   }
4620cf1dd8SToby Isaac   PetscFunctionReturn(0);
4720cf1dd8SToby Isaac }
4820cf1dd8SToby Isaac 
49*29b5c115SMatthew G. Knepley static PetscErrorCode PetscSpaceDestroy_Point(PetscSpace sp)
5020cf1dd8SToby Isaac {
5120cf1dd8SToby Isaac   PetscSpace_Point *pt = (PetscSpace_Point *) sp->data;
5220cf1dd8SToby Isaac   PetscErrorCode    ierr;
5320cf1dd8SToby Isaac 
5420cf1dd8SToby Isaac   PetscFunctionBegin;
5520cf1dd8SToby Isaac   ierr = PetscQuadratureDestroy(&pt->quad);CHKERRQ(ierr);
5620cf1dd8SToby Isaac   ierr = PetscFree(pt);CHKERRQ(ierr);
5720cf1dd8SToby Isaac   PetscFunctionReturn(0);
5820cf1dd8SToby Isaac }
5920cf1dd8SToby Isaac 
60*29b5c115SMatthew G. Knepley static PetscErrorCode PetscSpaceGetDimension_Point(PetscSpace sp, PetscInt *dim)
6120cf1dd8SToby Isaac {
6220cf1dd8SToby Isaac   PetscSpace_Point *pt = (PetscSpace_Point *) sp->data;
6320cf1dd8SToby Isaac 
6420cf1dd8SToby Isaac   PetscFunctionBegin;
6520cf1dd8SToby Isaac   *dim = pt->quad->numPoints;
6620cf1dd8SToby Isaac   PetscFunctionReturn(0);
6720cf1dd8SToby Isaac }
6820cf1dd8SToby Isaac 
69*29b5c115SMatthew G. Knepley static PetscErrorCode PetscSpaceEvaluate_Point(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[])
7020cf1dd8SToby Isaac {
7120cf1dd8SToby Isaac   PetscSpace_Point *pt  = (PetscSpace_Point *) sp->data;
7220cf1dd8SToby Isaac   PetscInt          dim = sp->Nv, pdim = pt->quad->numPoints, d, p, i, c;
7320cf1dd8SToby Isaac   PetscErrorCode    ierr;
7420cf1dd8SToby Isaac 
7520cf1dd8SToby Isaac   PetscFunctionBegin;
7620cf1dd8SToby 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);
77580bdb30SBarry Smith   ierr = PetscArrayzero(B, npoints*pdim);CHKERRQ(ierr);
7820cf1dd8SToby Isaac   for (p = 0; p < npoints; ++p) {
7920cf1dd8SToby Isaac     for (i = 0; i < pdim; ++i) {
8020cf1dd8SToby Isaac       for (d = 0; d < dim; ++d) {
8120cf1dd8SToby Isaac         if (PetscAbsReal(points[p*dim+d] - pt->quad->points[p*dim+d]) > 1.0e-10) break;
8220cf1dd8SToby Isaac       }
8320cf1dd8SToby Isaac       if (d >= dim) {B[p*pdim+i] = 1.0; break;}
8420cf1dd8SToby Isaac     }
8520cf1dd8SToby Isaac   }
8620cf1dd8SToby Isaac   /* Replicate for other components */
8720cf1dd8SToby Isaac   for (c = 1; c < sp->Nc; ++c) {
8820cf1dd8SToby Isaac     for (p = 0; p < npoints; ++p) {
8920cf1dd8SToby Isaac       for (i = 0; i < pdim; ++i) {
9020cf1dd8SToby Isaac         B[(c*npoints + p)*pdim + i] = B[p*pdim + i];
9120cf1dd8SToby Isaac       }
9220cf1dd8SToby Isaac     }
9320cf1dd8SToby Isaac   }
94580bdb30SBarry Smith   if (D) {ierr = PetscArrayzero(D, npoints*pdim*dim);CHKERRQ(ierr);}
95580bdb30SBarry Smith   if (H) {ierr = PetscArrayzero(H, npoints*pdim*dim*dim);CHKERRQ(ierr);}
9620cf1dd8SToby Isaac   PetscFunctionReturn(0);
9720cf1dd8SToby Isaac }
9820cf1dd8SToby Isaac 
99*29b5c115SMatthew G. Knepley static PetscErrorCode PetscSpaceInitialize_Point(PetscSpace sp)
10020cf1dd8SToby Isaac {
10120cf1dd8SToby Isaac   PetscFunctionBegin;
10220cf1dd8SToby Isaac   sp->ops->setfromoptions = NULL;
10320cf1dd8SToby Isaac   sp->ops->setup          = PetscSpaceSetUp_Point;
10420cf1dd8SToby Isaac   sp->ops->view           = PetscSpaceView_Point;
10520cf1dd8SToby Isaac   sp->ops->destroy        = PetscSpaceDestroy_Point;
10620cf1dd8SToby Isaac   sp->ops->getdimension   = PetscSpaceGetDimension_Point;
10720cf1dd8SToby Isaac   sp->ops->evaluate       = PetscSpaceEvaluate_Point;
10820cf1dd8SToby Isaac   PetscFunctionReturn(0);
10920cf1dd8SToby Isaac }
11020cf1dd8SToby Isaac 
11120cf1dd8SToby Isaac /*MC
11220cf1dd8SToby Isaac   PETSCSPACEPOINT = "point" - A PetscSpace object that encapsulates functions defined on a set of quadrature points.
11320cf1dd8SToby Isaac 
11420cf1dd8SToby Isaac   Level: intermediate
11520cf1dd8SToby Isaac 
11620cf1dd8SToby Isaac .seealso: PetscSpaceType, PetscSpaceCreate(), PetscSpaceSetType()
11720cf1dd8SToby Isaac M*/
11820cf1dd8SToby Isaac 
11920cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscSpaceCreate_Point(PetscSpace sp)
12020cf1dd8SToby Isaac {
12120cf1dd8SToby Isaac   PetscSpace_Point *pt;
12220cf1dd8SToby Isaac   PetscErrorCode    ierr;
12320cf1dd8SToby Isaac 
12420cf1dd8SToby Isaac   PetscFunctionBegin;
12520cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
12620cf1dd8SToby Isaac   ierr     = PetscNewLog(sp,&pt);CHKERRQ(ierr);
12720cf1dd8SToby Isaac   sp->data = pt;
12820cf1dd8SToby Isaac 
12920cf1dd8SToby Isaac   sp->Nv = 0;
13020cf1dd8SToby Isaac   sp->maxDegree = PETSC_MAX_INT;
13120cf1dd8SToby Isaac   ierr = PetscQuadratureCreate(PETSC_COMM_SELF, &pt->quad);CHKERRQ(ierr);
13220cf1dd8SToby Isaac   ierr = PetscQuadratureSetData(pt->quad, 0, 1, 0, NULL, NULL);CHKERRQ(ierr);
13320cf1dd8SToby Isaac 
13420cf1dd8SToby Isaac   ierr = PetscSpaceInitialize_Point(sp);CHKERRQ(ierr);
13520cf1dd8SToby Isaac   PetscFunctionReturn(0);
13620cf1dd8SToby Isaac }
13720cf1dd8SToby Isaac 
13820cf1dd8SToby Isaac /*@
13920cf1dd8SToby Isaac   PetscSpacePointSetPoints - Sets the evaluation points for the space to coincide with the points of a quadrature rule
14020cf1dd8SToby Isaac 
14120cf1dd8SToby Isaac   Logically collective
14220cf1dd8SToby Isaac 
14320cf1dd8SToby Isaac   Input Parameters:
14420cf1dd8SToby Isaac + sp - The PetscSpace
14520cf1dd8SToby Isaac - q  - The PetscQuadrature defining the points
14620cf1dd8SToby Isaac 
14720cf1dd8SToby Isaac   Level: intermediate
14820cf1dd8SToby Isaac 
14920cf1dd8SToby Isaac .seealso: PetscSpaceCreate(), PetscSpaceSetType()
15020cf1dd8SToby Isaac @*/
15120cf1dd8SToby Isaac PetscErrorCode PetscSpacePointSetPoints(PetscSpace sp, PetscQuadrature q)
15220cf1dd8SToby Isaac {
15320cf1dd8SToby Isaac   PetscSpace_Point *pt = (PetscSpace_Point *) sp->data;
15420cf1dd8SToby Isaac   PetscErrorCode    ierr;
15520cf1dd8SToby Isaac 
15620cf1dd8SToby Isaac   PetscFunctionBegin;
15720cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
15820cf1dd8SToby Isaac   PetscValidHeaderSpecific(q, PETSC_OBJECT_CLASSID, 2);
15920cf1dd8SToby Isaac   ierr = PetscQuadratureDestroy(&pt->quad);CHKERRQ(ierr);
16020cf1dd8SToby Isaac   ierr = PetscQuadratureDuplicate(q, &pt->quad);CHKERRQ(ierr);
16120cf1dd8SToby Isaac   PetscFunctionReturn(0);
16220cf1dd8SToby Isaac }
16320cf1dd8SToby Isaac 
16420cf1dd8SToby Isaac /*@
16520cf1dd8SToby Isaac   PetscSpacePointGetPoints - Gets the evaluation points for the space as the points of a quadrature rule
16620cf1dd8SToby Isaac 
16720cf1dd8SToby Isaac   Logically collective
16820cf1dd8SToby Isaac 
16920cf1dd8SToby Isaac   Input Parameter:
17020cf1dd8SToby Isaac . sp - The PetscSpace
17120cf1dd8SToby Isaac 
17220cf1dd8SToby Isaac   Output Parameter:
17320cf1dd8SToby Isaac . q  - The PetscQuadrature defining the points
17420cf1dd8SToby Isaac 
17520cf1dd8SToby Isaac   Level: intermediate
17620cf1dd8SToby Isaac 
17720cf1dd8SToby Isaac .seealso: PetscSpaceCreate(), PetscSpaceSetType()
17820cf1dd8SToby Isaac @*/
17920cf1dd8SToby Isaac PetscErrorCode PetscSpacePointGetPoints(PetscSpace sp, PetscQuadrature *q)
18020cf1dd8SToby Isaac {
18120cf1dd8SToby Isaac   PetscSpace_Point *pt = (PetscSpace_Point *) sp->data;
18220cf1dd8SToby Isaac 
18320cf1dd8SToby Isaac   PetscFunctionBegin;
18420cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
18520cf1dd8SToby Isaac   PetscValidPointer(q, 2);
18620cf1dd8SToby Isaac   *q = pt->quad;
18720cf1dd8SToby Isaac   PetscFunctionReturn(0);
18820cf1dd8SToby Isaac }
189