xref: /petsc/src/dm/field/impls/shell/dmfieldshell.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
151a454edSToby Isaac #include <petsc/private/dmfieldimpl.h> /*I "petscdmfield.h" I*/
251a454edSToby Isaac 
39371c9d4SSatish Balay typedef struct _n_DMField_Shell {
4028944aaSToby Isaac   void *ctx;
5028944aaSToby Isaac   PetscErrorCode (*destroy)(DMField);
69371c9d4SSatish Balay } DMField_Shell;
7028944aaSToby Isaac 
89371c9d4SSatish Balay PetscErrorCode DMFieldShellGetContext(DMField field, void *ctx) {
9028944aaSToby Isaac   PetscBool flg;
10028944aaSToby Isaac 
1151a454edSToby Isaac   PetscFunctionBegin;
12028944aaSToby Isaac   PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1);
13028944aaSToby Isaac   PetscValidPointer(ctx, 2);
149566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)field, DMFIELDSHELL, &flg));
15028944aaSToby Isaac   if (flg) *(void **)ctx = ((DMField_Shell *)(field->data))->ctx;
16028944aaSToby Isaac   else SETERRQ(PetscObjectComm((PetscObject)field), PETSC_ERR_SUP, "Cannot get context from non-shell shield");
1751a454edSToby Isaac   PetscFunctionReturn(0);
1851a454edSToby Isaac }
1951a454edSToby Isaac 
209371c9d4SSatish Balay static PetscErrorCode DMFieldDestroy_Shell(DMField field) {
21028944aaSToby Isaac   DMField_Shell *shell = (DMField_Shell *)field->data;
2251a454edSToby Isaac 
23028944aaSToby Isaac   PetscFunctionBegin;
249566063dSJacob Faibussowitsch   if (shell->destroy) PetscCall((*(shell->destroy))(field));
259566063dSJacob Faibussowitsch   PetscCall(PetscFree(field->data));
26028944aaSToby Isaac   PetscFunctionReturn(0);
27028944aaSToby Isaac }
28028944aaSToby Isaac 
299371c9d4SSatish Balay PetscErrorCode DMFieldShellEvaluateFEDefault(DMField field, IS pointIS, PetscQuadrature quad, PetscDataType type, void *B, void *D, void *H) {
30028944aaSToby Isaac   DM           dm = field->dm;
31028944aaSToby Isaac   DMField      coordField;
32028944aaSToby Isaac   PetscFEGeom *geom;
33028944aaSToby Isaac   Vec          pushforward;
34028944aaSToby Isaac   PetscInt     dimC, dim, numPoints, Nq, p, Nc;
35028944aaSToby Isaac   PetscScalar *pfArray;
36028944aaSToby Isaac 
37028944aaSToby Isaac   PetscFunctionBegin;
38028944aaSToby Isaac   Nc = field->numComponents;
399566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateField(dm, &coordField));
409566063dSJacob Faibussowitsch   PetscCall(DMFieldCreateFEGeom(coordField, pointIS, quad, PETSC_FALSE, &geom));
419566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &dimC));
429566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, &dim, NULL, &Nq, NULL, NULL));
439566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(pointIS, &numPoints));
449566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dimC * Nq * numPoints, &pfArray));
45028944aaSToby Isaac   for (p = 0; p < numPoints * dimC * Nq; p++) pfArray[p] = (PetscScalar)geom->v[p];
469566063dSJacob Faibussowitsch   PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)pointIS), dimC, dimC * Nq * numPoints, PETSC_DETERMINE, pfArray, &pushforward));
479566063dSJacob Faibussowitsch   PetscCall(DMFieldEvaluate(field, pushforward, type, B, D, H));
48028944aaSToby Isaac   /* TODO: handle covariant/contravariant pullbacks */
49028944aaSToby Isaac   if (D) {
50028944aaSToby Isaac     if (type == PETSC_SCALAR) {
51028944aaSToby Isaac       PetscScalar *sD = (PetscScalar *)D;
52028944aaSToby Isaac       PetscInt     q;
53028944aaSToby Isaac 
54028944aaSToby Isaac       for (p = 0; p < numPoints * Nq; p++) {
55028944aaSToby Isaac         for (q = 0; q < Nc; q++) {
56028944aaSToby Isaac           PetscScalar d[3];
57028944aaSToby Isaac 
58028944aaSToby Isaac           PetscInt i, j;
59028944aaSToby Isaac 
60028944aaSToby Isaac           for (i = 0; i < dimC; i++) d[i] = sD[(p * Nc + q) * dimC + i];
61028944aaSToby Isaac           for (i = 0; i < dimC; i++) sD[(p * Nc + q) * dimC + i] = 0.;
62028944aaSToby Isaac           for (i = 0; i < dimC; i++) {
63ad540459SPierre Jolivet             for (j = 0; j < dimC; j++) sD[(p * Nc + q) * dimC + i] += geom->J[(p * dimC + j) * dimC + i] * d[j];
64028944aaSToby Isaac           }
65028944aaSToby Isaac         }
66028944aaSToby Isaac       }
67028944aaSToby Isaac     } else {
68028944aaSToby Isaac       PetscReal *rD = (PetscReal *)D;
69028944aaSToby Isaac       PetscInt   q;
70028944aaSToby Isaac 
71028944aaSToby Isaac       for (p = 0; p < numPoints * Nq; p++) {
72028944aaSToby Isaac         for (q = 0; q < Nc; q++) {
73028944aaSToby Isaac           PetscReal d[3];
74028944aaSToby Isaac 
75028944aaSToby Isaac           PetscInt i, j;
76028944aaSToby Isaac 
77028944aaSToby Isaac           for (i = 0; i < dimC; i++) d[i] = rD[(p * Nc + q) * dimC + i];
78028944aaSToby Isaac           for (i = 0; i < dimC; i++) rD[(p * Nc + q) * dimC + i] = 0.;
79028944aaSToby Isaac           for (i = 0; i < dimC; i++) {
80ad540459SPierre Jolivet             for (j = 0; j < dimC; j++) rD[(p * Nc + q) * dimC + i] += geom->J[(p * dimC + j) * dimC + i] * d[j];
81028944aaSToby Isaac           }
82028944aaSToby Isaac         }
83028944aaSToby Isaac       }
84028944aaSToby Isaac     }
85028944aaSToby Isaac   }
86028944aaSToby Isaac   if (H) {
87028944aaSToby Isaac     if (type == PETSC_SCALAR) {
88028944aaSToby Isaac       PetscScalar *sH = (PetscScalar *)H;
89028944aaSToby Isaac       PetscInt     q;
90028944aaSToby Isaac 
91028944aaSToby Isaac       for (p = 0; p < numPoints * Nq; p++) {
92028944aaSToby Isaac         for (q = 0; q < Nc; q++) {
93028944aaSToby Isaac           PetscScalar d[3][3];
94028944aaSToby Isaac 
95028944aaSToby Isaac           PetscInt i, j, k, l;
96028944aaSToby Isaac 
979371c9d4SSatish Balay           for (i = 0; i < dimC; i++)
989371c9d4SSatish Balay             for (j = 0; j < dimC; j++) d[i][j] = sH[((p * Nc + q) * dimC + i) * dimC + j];
999371c9d4SSatish Balay           for (i = 0; i < dimC; i++)
1009371c9d4SSatish Balay             for (j = 0; j < dimC; j++) sH[((p * Nc + q) * dimC + i) * dimC + j] = 0.;
101028944aaSToby Isaac           for (i = 0; i < dimC; i++) {
102028944aaSToby Isaac             for (j = 0; j < dimC; j++) {
103028944aaSToby Isaac               for (k = 0; k < dimC; k++) {
104ad540459SPierre Jolivet                 for (l = 0; l < dimC; l++) sH[((p * Nc + q) * dimC + i) * dimC + j] += geom->J[(p * dimC + k) * dimC + i] * geom->J[(p * dimC + l) * dimC + j] * d[k][l];
105028944aaSToby Isaac               }
106028944aaSToby Isaac             }
107028944aaSToby Isaac           }
108028944aaSToby Isaac         }
109028944aaSToby Isaac       }
110028944aaSToby Isaac     } else {
111028944aaSToby Isaac       PetscReal *rH = (PetscReal *)H;
112028944aaSToby Isaac       PetscInt   q;
113028944aaSToby Isaac 
114028944aaSToby Isaac       for (p = 0; p < numPoints * Nq; p++) {
115028944aaSToby Isaac         for (q = 0; q < Nc; q++) {
116028944aaSToby Isaac           PetscReal d[3][3];
117028944aaSToby Isaac 
118028944aaSToby Isaac           PetscInt i, j, k, l;
119028944aaSToby Isaac 
1209371c9d4SSatish Balay           for (i = 0; i < dimC; i++)
1219371c9d4SSatish Balay             for (j = 0; j < dimC; j++) d[i][j] = rH[((p * Nc + q) * dimC + i) * dimC + j];
1229371c9d4SSatish Balay           for (i = 0; i < dimC; i++)
1239371c9d4SSatish Balay             for (j = 0; j < dimC; j++) rH[((p * Nc + q) * dimC + i) * dimC + j] = 0.;
124028944aaSToby Isaac           for (i = 0; i < dimC; i++) {
125028944aaSToby Isaac             for (j = 0; j < dimC; j++) {
126028944aaSToby Isaac               for (k = 0; k < dimC; k++) {
127ad540459SPierre Jolivet                 for (l = 0; l < dimC; l++) rH[((p * Nc + q) * dimC + i) * dimC + j] += geom->J[(p * dimC + k) * dimC + i] * geom->J[(p * dimC + l) * dimC + j] * d[k][l];
128028944aaSToby Isaac               }
129028944aaSToby Isaac             }
130028944aaSToby Isaac           }
131028944aaSToby Isaac         }
132028944aaSToby Isaac       }
133028944aaSToby Isaac     }
134028944aaSToby Isaac   }
1359566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pushforward));
1369566063dSJacob Faibussowitsch   PetscCall(PetscFree(pfArray));
1379566063dSJacob Faibussowitsch   PetscCall(PetscFEGeomDestroy(&geom));
138028944aaSToby Isaac   PetscFunctionReturn(0);
139028944aaSToby Isaac }
140028944aaSToby Isaac 
1419371c9d4SSatish Balay PetscErrorCode DMFieldShellEvaluateFVDefault(DMField field, IS pointIS, PetscDataType type, void *B, void *D, void *H) {
142028944aaSToby Isaac   DM              dm = field->dm;
143028944aaSToby Isaac   DMField         coordField;
144028944aaSToby Isaac   PetscFEGeom    *geom;
145028944aaSToby Isaac   Vec             pushforward;
146028944aaSToby Isaac   PetscInt        dimC, dim, numPoints, Nq, p;
147028944aaSToby Isaac   PetscScalar    *pfArray;
148028944aaSToby Isaac   PetscQuadrature quad;
149a2aba86cSMatthew G. Knepley   MPI_Comm        comm;
150028944aaSToby Isaac 
151028944aaSToby Isaac   PetscFunctionBegin;
1529566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)field, &comm));
1539566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
1549566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &dimC));
1559566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateField(dm, &coordField));
1569566063dSJacob Faibussowitsch   PetscCall(DMFieldGetFVQuadrature_Internal(coordField, pointIS, &quad));
157a2aba86cSMatthew G. Knepley   PetscCheck(quad, comm, PETSC_ERR_ARG_WRONGSTATE, "coordinate field must have default quadrature for FV computation");
1589566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &Nq, NULL, NULL));
159a2aba86cSMatthew G. Knepley   PetscCheck(Nq == 1, comm, PETSC_ERR_ARG_WRONGSTATE, "quadrature must have only one point");
1609566063dSJacob Faibussowitsch   PetscCall(DMFieldCreateFEGeom(coordField, pointIS, quad, PETSC_FALSE, &geom));
1619566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(pointIS, &numPoints));
1629566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(dimC * numPoints, &pfArray));
163028944aaSToby Isaac   for (p = 0; p < numPoints * dimC; p++) pfArray[p] = (PetscScalar)geom->v[p];
1649566063dSJacob Faibussowitsch   PetscCall(VecCreateMPIWithArray(PetscObjectComm((PetscObject)pointIS), dimC, dimC * numPoints, PETSC_DETERMINE, pfArray, &pushforward));
1659566063dSJacob Faibussowitsch   PetscCall(DMFieldEvaluate(field, pushforward, type, B, D, H));
1669566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&quad));
1679566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&pushforward));
1689566063dSJacob Faibussowitsch   PetscCall(PetscFree(pfArray));
1699566063dSJacob Faibussowitsch   PetscCall(PetscFEGeomDestroy(&geom));
170028944aaSToby Isaac   PetscFunctionReturn(0);
171028944aaSToby Isaac }
172028944aaSToby Isaac 
1739371c9d4SSatish Balay PetscErrorCode DMFieldShellSetDestroy(DMField field, PetscErrorCode (*destroy)(DMField)) {
174028944aaSToby Isaac   DMField_Shell *shell = (DMField_Shell *)field->data;
175028944aaSToby Isaac 
176028944aaSToby Isaac   PetscFunctionBegin;
177028944aaSToby Isaac   PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1);
178028944aaSToby Isaac   shell->destroy = destroy;
179028944aaSToby Isaac   PetscFunctionReturn(0);
180028944aaSToby Isaac }
181028944aaSToby Isaac 
1829371c9d4SSatish Balay PetscErrorCode DMFieldShellSetEvaluate(DMField field, PetscErrorCode (*evaluate)(DMField, Vec, PetscDataType, void *, void *, void *)) {
183028944aaSToby Isaac   PetscFunctionBegin;
184028944aaSToby Isaac   PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1);
185028944aaSToby Isaac   field->ops->evaluate = evaluate;
186028944aaSToby Isaac   PetscFunctionReturn(0);
187028944aaSToby Isaac }
188028944aaSToby Isaac 
1899371c9d4SSatish Balay PetscErrorCode DMFieldShellSetEvaluateFE(DMField field, PetscErrorCode (*evaluateFE)(DMField, IS, PetscQuadrature, PetscDataType, void *, void *, void *)) {
190028944aaSToby Isaac   PetscFunctionBegin;
191028944aaSToby Isaac   PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1);
192028944aaSToby Isaac   field->ops->evaluateFE = evaluateFE;
193028944aaSToby Isaac   PetscFunctionReturn(0);
194028944aaSToby Isaac }
195028944aaSToby Isaac 
1969371c9d4SSatish Balay PetscErrorCode DMFieldShellSetEvaluateFV(DMField field, PetscErrorCode (*evaluateFV)(DMField, IS, PetscDataType, void *, void *, void *)) {
197028944aaSToby Isaac   PetscFunctionBegin;
198028944aaSToby Isaac   PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1);
199028944aaSToby Isaac   field->ops->evaluateFV = evaluateFV;
200028944aaSToby Isaac   PetscFunctionReturn(0);
201028944aaSToby Isaac }
202028944aaSToby Isaac 
2039371c9d4SSatish Balay PetscErrorCode DMFieldShellSetGetDegree(DMField field, PetscErrorCode (*getDegree)(DMField, IS, PetscInt *, PetscInt *)) {
204028944aaSToby Isaac   PetscFunctionBegin;
205028944aaSToby Isaac   PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1);
206b7260050SToby Isaac   field->ops->getDegree = getDegree;
207028944aaSToby Isaac   PetscFunctionReturn(0);
208028944aaSToby Isaac }
209028944aaSToby Isaac 
2109371c9d4SSatish Balay PetscErrorCode DMFieldShellSetCreateDefaultQuadrature(DMField field, PetscErrorCode (*createDefaultQuadrature)(DMField, IS, PetscQuadrature *)) {
211028944aaSToby Isaac   PetscFunctionBegin;
212028944aaSToby Isaac   PetscValidHeaderSpecific(field, DMFIELD_CLASSID, 1);
213028944aaSToby Isaac   field->ops->createDefaultQuadrature = createDefaultQuadrature;
214028944aaSToby Isaac   PetscFunctionReturn(0);
215028944aaSToby Isaac }
216028944aaSToby Isaac 
2179371c9d4SSatish Balay static PetscErrorCode DMFieldInitialize_Shell(DMField field) {
218028944aaSToby Isaac   PetscFunctionBegin;
219028944aaSToby Isaac   field->ops->destroy                 = DMFieldDestroy_Shell;
220028944aaSToby Isaac   field->ops->evaluate                = NULL;
221028944aaSToby Isaac   field->ops->evaluateFE              = DMFieldShellEvaluateFEDefault;
222028944aaSToby Isaac   field->ops->evaluateFV              = DMFieldShellEvaluateFVDefault;
223b7260050SToby Isaac   field->ops->getDegree               = NULL;
224028944aaSToby Isaac   field->ops->createDefaultQuadrature = NULL;
225028944aaSToby Isaac   field->ops->view                    = NULL;
226028944aaSToby Isaac   PetscFunctionReturn(0);
227028944aaSToby Isaac }
228028944aaSToby Isaac 
2299371c9d4SSatish Balay PETSC_INTERN PetscErrorCode DMFieldCreate_Shell(DMField field) {
230028944aaSToby Isaac   DMField_Shell *shell;
231028944aaSToby Isaac 
232028944aaSToby Isaac   PetscFunctionBegin;
233*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&shell));
234028944aaSToby Isaac   field->data = shell;
2359566063dSJacob Faibussowitsch   PetscCall(DMFieldInitialize_Shell(field));
236028944aaSToby Isaac   PetscFunctionReturn(0);
237028944aaSToby Isaac }
238028944aaSToby Isaac 
2399371c9d4SSatish Balay PetscErrorCode DMFieldCreateShell(DM dm, PetscInt numComponents, DMFieldContinuity continuity, void *ctx, DMField *field) {
240028944aaSToby Isaac   DMField        b;
241028944aaSToby Isaac   DMField_Shell *shell;
242028944aaSToby Isaac 
243028944aaSToby Isaac   PetscFunctionBegin;
244028944aaSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
245028944aaSToby Isaac   if (ctx) PetscValidPointer(ctx, 4);
246028944aaSToby Isaac   PetscValidPointer(field, 5);
2479566063dSJacob Faibussowitsch   PetscCall(DMFieldCreate(dm, numComponents, continuity, &b));
2489566063dSJacob Faibussowitsch   PetscCall(DMFieldSetType(b, DMFIELDSHELL));
249028944aaSToby Isaac   shell      = (DMField_Shell *)b->data;
250028944aaSToby Isaac   shell->ctx = ctx;
251028944aaSToby Isaac   *field     = b;
252028944aaSToby Isaac   PetscFunctionReturn(0);
253028944aaSToby Isaac }
254