xref: /petsc/src/dm/field/impls/shell/dmfieldshell.c (revision 028944aa1150740f314ac582cc1ff363c357a432)
151a454edSToby Isaac #include <petsc/private/dmfieldimpl.h> /*I "petscdmfield.h" I*/
251a454edSToby Isaac 
3*028944aaSToby Isaac typedef struct _n_DMField_Shell
451a454edSToby Isaac {
5*028944aaSToby Isaac   void *ctx;
6*028944aaSToby Isaac   PetscErrorCode (*destroy) (DMField);
7*028944aaSToby Isaac }
8*028944aaSToby Isaac DMField_Shell;
9*028944aaSToby Isaac 
10*028944aaSToby Isaac PetscErrorCode DMFieldShellGetContext(DMField field, void *ctx)
11*028944aaSToby Isaac {
12*028944aaSToby Isaac   PetscBool      flg;
13*028944aaSToby Isaac   PetscErrorCode ierr;
14*028944aaSToby Isaac 
1551a454edSToby Isaac   PetscFunctionBegin;
16*028944aaSToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
17*028944aaSToby Isaac   PetscValidPointer(ctx,2);
18*028944aaSToby Isaac   ierr = PetscObjectTypeCompare((PetscObject)field,DMFIELDSHELL,&flg);CHKERRQ(ierr);
19*028944aaSToby Isaac   if (flg) *(void**)ctx = ((DMField_Shell*)(field->data))->ctx;
20*028944aaSToby Isaac   else SETERRQ(PetscObjectComm((PetscObject)field),PETSC_ERR_SUP,"Cannot get context from non-shell shield");
2151a454edSToby Isaac   PetscFunctionReturn(0);
2251a454edSToby Isaac }
2351a454edSToby Isaac 
24*028944aaSToby Isaac static PetscErrorCode DMFieldDestroy_Shell(DMField field)
25*028944aaSToby Isaac {
26*028944aaSToby Isaac   DMField_Shell *shell = (DMField_Shell *) field->data;
27*028944aaSToby Isaac   PetscErrorCode ierr;
2851a454edSToby Isaac 
29*028944aaSToby Isaac   PetscFunctionBegin;
30*028944aaSToby Isaac   if (shell->destroy) {ierr = (*(shell->destroy)) (field);CHKERRQ(ierr);}
31*028944aaSToby Isaac   ierr = PetscFree(field->data);CHKERRQ(ierr);
32*028944aaSToby Isaac   PetscFunctionReturn(0);
33*028944aaSToby Isaac }
34*028944aaSToby Isaac 
35*028944aaSToby Isaac PetscErrorCode DMFieldShellEvaluateFEDefault(DMField field, IS pointIS, PetscQuadrature quad, PetscDataType type, void *B, void *D, void *H)
36*028944aaSToby Isaac {
37*028944aaSToby Isaac   DM              dm = field->dm;
38*028944aaSToby Isaac   DMField         coordField;
39*028944aaSToby Isaac   PetscFEGeom    *geom;
40*028944aaSToby Isaac   Vec             pushforward;
41*028944aaSToby Isaac   PetscInt        dimC, dim, numPoints, Nq, p, Nc;
42*028944aaSToby Isaac   PetscScalar    *pfArray;
43*028944aaSToby Isaac   PetscErrorCode  ierr;
44*028944aaSToby Isaac 
45*028944aaSToby Isaac   PetscFunctionBegin;
46*028944aaSToby Isaac   Nc   = field->numComponents;
47*028944aaSToby Isaac   ierr = DMGetCoordinateField(dm, &coordField);CHKERRQ(ierr);
48*028944aaSToby Isaac   ierr = DMFieldCreateFEGeom(coordField, pointIS, quad, PETSC_FALSE, &geom);CHKERRQ(ierr);
49*028944aaSToby Isaac   ierr = DMGetCoordinateDim(dm, &dimC);CHKERRQ(ierr);
50*028944aaSToby Isaac   ierr = PetscQuadratureGetData(quad, &dim, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);
51*028944aaSToby Isaac   ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr);
52*028944aaSToby Isaac   ierr = PetscMalloc1(dimC * Nq * numPoints, &pfArray);CHKERRQ(ierr);
53*028944aaSToby Isaac   for (p = 0; p < numPoints * dimC * Nq; p++) pfArray[p] = (PetscScalar) geom->v[p];
54*028944aaSToby Isaac   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)pointIS), dimC, dimC * Nq * numPoints, PETSC_DETERMINE, pfArray, &pushforward);CHKERRQ(ierr);
55*028944aaSToby Isaac   ierr = DMFieldEvaluate(field, pushforward, type, B, D, H);CHKERRQ(ierr);
56*028944aaSToby Isaac   /* TODO: handle covariant/contravariant pullbacks */
57*028944aaSToby Isaac   if (D) {
58*028944aaSToby Isaac     if (type == PETSC_SCALAR) {
59*028944aaSToby Isaac       PetscScalar *sD = (PetscScalar *) D;
60*028944aaSToby Isaac       PetscInt q;
61*028944aaSToby Isaac 
62*028944aaSToby Isaac       for (p = 0; p < numPoints * Nq; p++) {
63*028944aaSToby Isaac         for (q = 0; q < Nc; q++) {
64*028944aaSToby Isaac           PetscScalar d[3];
65*028944aaSToby Isaac 
66*028944aaSToby Isaac           PetscInt i, j;
67*028944aaSToby Isaac 
68*028944aaSToby Isaac           for (i = 0; i < dimC; i++) d[i] = sD[(p * Nc + q) * dimC + i];
69*028944aaSToby Isaac           for (i = 0; i < dimC; i++) sD[(p * Nc + q) * dimC + i] = 0.;
70*028944aaSToby Isaac           for (i = 0; i < dimC; i++) {
71*028944aaSToby Isaac             for (j = 0; j < dimC; j++) {
72*028944aaSToby Isaac               sD[(p * Nc + q) * dimC + i] += geom->J[(p * dimC + j) * dimC + i] * d[j];
73*028944aaSToby Isaac             }
74*028944aaSToby Isaac           }
75*028944aaSToby Isaac         }
76*028944aaSToby Isaac       }
77*028944aaSToby Isaac     } else {
78*028944aaSToby Isaac       PetscReal *rD = (PetscReal *) D;
79*028944aaSToby Isaac       PetscInt q;
80*028944aaSToby Isaac 
81*028944aaSToby Isaac       for (p = 0; p < numPoints * Nq; p++) {
82*028944aaSToby Isaac         for (q = 0; q < Nc; q++) {
83*028944aaSToby Isaac           PetscReal d[3];
84*028944aaSToby Isaac 
85*028944aaSToby Isaac           PetscInt i, j;
86*028944aaSToby Isaac 
87*028944aaSToby Isaac           for (i = 0; i < dimC; i++) d[i] = rD[(p * Nc + q) * dimC + i];
88*028944aaSToby Isaac           for (i = 0; i < dimC; i++) rD[(p * Nc + q) * dimC + i] = 0.;
89*028944aaSToby Isaac           for (i = 0; i < dimC; i++) {
90*028944aaSToby Isaac             for (j = 0; j < dimC; j++) {
91*028944aaSToby Isaac               rD[(p * Nc + q) * dimC + i] += geom->J[(p * dimC + j) * dimC + i] * d[j];
92*028944aaSToby Isaac             }
93*028944aaSToby Isaac           }
94*028944aaSToby Isaac         }
95*028944aaSToby Isaac       }
96*028944aaSToby Isaac     }
97*028944aaSToby Isaac   }
98*028944aaSToby Isaac   if (H) {
99*028944aaSToby Isaac     if (type == PETSC_SCALAR) {
100*028944aaSToby Isaac       PetscScalar *sH = (PetscScalar *) H;
101*028944aaSToby Isaac       PetscInt q;
102*028944aaSToby Isaac 
103*028944aaSToby Isaac       for (p = 0; p < numPoints * Nq; p++) {
104*028944aaSToby Isaac         for (q = 0; q < Nc; q++) {
105*028944aaSToby Isaac           PetscScalar d[3][3];
106*028944aaSToby Isaac 
107*028944aaSToby Isaac           PetscInt i, j, k, l;
108*028944aaSToby Isaac 
109*028944aaSToby Isaac           for (i = 0; i < dimC; i++) for (j = 0; j < dimC; j++) d[i][j] = sH[((p * Nc + q) * dimC + i) * dimC + j];
110*028944aaSToby Isaac           for (i = 0; i < dimC; i++) for (j = 0; j < dimC; j++) sH[((p * Nc + q) * dimC + i) * dimC + j] = 0.;
111*028944aaSToby Isaac           for (i = 0; i < dimC; i++) {
112*028944aaSToby Isaac             for (j = 0; j < dimC; j++) {
113*028944aaSToby Isaac               for (k = 0; k < dimC; k++) {
114*028944aaSToby Isaac                 for (l = 0; l < dimC; l++) {
115*028944aaSToby Isaac                   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];
116*028944aaSToby Isaac                 }
117*028944aaSToby Isaac               }
118*028944aaSToby Isaac             }
119*028944aaSToby Isaac           }
120*028944aaSToby Isaac         }
121*028944aaSToby Isaac       }
122*028944aaSToby Isaac     } else {
123*028944aaSToby Isaac       PetscReal *rH = (PetscReal *) H;
124*028944aaSToby Isaac       PetscInt q;
125*028944aaSToby Isaac 
126*028944aaSToby Isaac       for (p = 0; p < numPoints * Nq; p++) {
127*028944aaSToby Isaac         for (q = 0; q < Nc; q++) {
128*028944aaSToby Isaac           PetscReal d[3][3];
129*028944aaSToby Isaac 
130*028944aaSToby Isaac           PetscInt i, j, k, l;
131*028944aaSToby Isaac 
132*028944aaSToby Isaac           for (i = 0; i < dimC; i++) for (j = 0; j < dimC; j++) d[i][j] = rH[((p * Nc + q) * dimC + i) * dimC + j];
133*028944aaSToby Isaac           for (i = 0; i < dimC; i++) for (j = 0; j < dimC; j++) rH[((p * Nc + q) * dimC + i) * dimC + j] = 0.;
134*028944aaSToby Isaac           for (i = 0; i < dimC; i++) {
135*028944aaSToby Isaac             for (j = 0; j < dimC; j++) {
136*028944aaSToby Isaac               for (k = 0; k < dimC; k++) {
137*028944aaSToby Isaac                 for (l = 0; l < dimC; l++) {
138*028944aaSToby Isaac                   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];
139*028944aaSToby Isaac                 }
140*028944aaSToby Isaac               }
141*028944aaSToby Isaac             }
142*028944aaSToby Isaac           }
143*028944aaSToby Isaac         }
144*028944aaSToby Isaac       }
145*028944aaSToby Isaac     }
146*028944aaSToby Isaac   }
147*028944aaSToby Isaac   ierr = VecDestroy(&pushforward);CHKERRQ(ierr);
148*028944aaSToby Isaac   ierr = PetscFree(pfArray);CHKERRQ(ierr);
149*028944aaSToby Isaac   ierr = PetscFEGeomDestroy(&geom);CHKERRQ(ierr);
150*028944aaSToby Isaac   PetscFunctionReturn(0);
151*028944aaSToby Isaac }
152*028944aaSToby Isaac 
153*028944aaSToby Isaac PetscErrorCode DMFieldShellEvaluateFVDefault(DMField field, IS pointIS, PetscDataType type, void *B, void *D, void *H)
154*028944aaSToby Isaac {
155*028944aaSToby Isaac   DM              dm = field->dm;
156*028944aaSToby Isaac   DMField         coordField;
157*028944aaSToby Isaac   PetscFEGeom    *geom;
158*028944aaSToby Isaac   Vec             pushforward;
159*028944aaSToby Isaac   PetscInt        dimC, dim, numPoints, Nq, p;
160*028944aaSToby Isaac   PetscScalar    *pfArray;
161*028944aaSToby Isaac   PetscQuadrature quad;
162*028944aaSToby Isaac   PetscErrorCode  ierr;
163*028944aaSToby Isaac 
164*028944aaSToby Isaac   PetscFunctionBegin;
165*028944aaSToby Isaac   ierr = DMGetCoordinateField(dm, &coordField);CHKERRQ(ierr);
166*028944aaSToby Isaac   ierr = DMFieldCreateDefaultQuadrature(coordField, pointIS, &quad);CHKERRQ(ierr);
167*028944aaSToby Isaac   if (!quad) SETERRQ(PetscObjectComm((PetscObject) pointIS), PETSC_ERR_ARG_WRONGSTATE, "coordinate field must have default quadrature for FV computation");
168*028944aaSToby Isaac   ierr = DMFieldCreateFEGeom(coordField, pointIS, quad, PETSC_FALSE, &geom);CHKERRQ(ierr);
169*028944aaSToby Isaac   ierr = DMGetCoordinateDim(dm, &dimC);CHKERRQ(ierr);
170*028944aaSToby Isaac   ierr = PetscQuadratureGetData(quad, &dim, NULL, &Nq, NULL, NULL);CHKERRQ(ierr);
171*028944aaSToby Isaac   if (Nq != 1) SETERRQ(PetscObjectComm((PetscObject) quad), PETSC_ERR_ARG_WRONGSTATE, "quadrature must have only one point");
172*028944aaSToby Isaac   ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr);
173*028944aaSToby Isaac   ierr = PetscMalloc1(dimC * numPoints, &pfArray);CHKERRQ(ierr);
174*028944aaSToby Isaac   for (p = 0; p < numPoints * dimC; p++) pfArray[p] = (PetscScalar) geom->v[p];
175*028944aaSToby Isaac   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)pointIS), dimC, dimC * numPoints, PETSC_DETERMINE, pfArray, &pushforward);CHKERRQ(ierr);
176*028944aaSToby Isaac   ierr = DMFieldEvaluate(field, pushforward, type, B, D, H);CHKERRQ(ierr);
177*028944aaSToby Isaac   ierr = PetscQuadratureDestroy(&quad);CHKERRQ(ierr);
178*028944aaSToby Isaac   ierr = VecDestroy(&pushforward);CHKERRQ(ierr);
179*028944aaSToby Isaac   ierr = PetscFree(pfArray);CHKERRQ(ierr);
180*028944aaSToby Isaac   ierr = PetscFEGeomDestroy(&geom);CHKERRQ(ierr);
181*028944aaSToby Isaac   PetscFunctionReturn(0);
182*028944aaSToby Isaac }
183*028944aaSToby Isaac 
184*028944aaSToby Isaac PetscErrorCode DMFieldShellSetDestroy(DMField field, PetscErrorCode (*destroy)(DMField))
185*028944aaSToby Isaac {
186*028944aaSToby Isaac   DMField_Shell *shell = (DMField_Shell *) field->data;
187*028944aaSToby Isaac 
188*028944aaSToby Isaac   PetscFunctionBegin;
189*028944aaSToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
190*028944aaSToby Isaac   shell->destroy = destroy;
191*028944aaSToby Isaac   PetscFunctionReturn(0);
192*028944aaSToby Isaac }
193*028944aaSToby Isaac 
194*028944aaSToby Isaac PetscErrorCode DMFieldShellSetEvaluate(DMField field, PetscErrorCode (*evaluate)(DMField,Vec,PetscDataType,void*,void*,void*))
195*028944aaSToby Isaac {
196*028944aaSToby Isaac   PetscFunctionBegin;
197*028944aaSToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
198*028944aaSToby Isaac   field->ops->evaluate = evaluate;
199*028944aaSToby Isaac   PetscFunctionReturn(0);
200*028944aaSToby Isaac }
201*028944aaSToby Isaac 
202*028944aaSToby Isaac PetscErrorCode DMFieldShellSetEvaluateFE(DMField field, PetscErrorCode (*evaluateFE)(DMField,IS,PetscQuadrature,PetscDataType,void*,void*,void*))
203*028944aaSToby Isaac {
204*028944aaSToby Isaac   PetscFunctionBegin;
205*028944aaSToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
206*028944aaSToby Isaac   field->ops->evaluateFE = evaluateFE;
207*028944aaSToby Isaac   PetscFunctionReturn(0);
208*028944aaSToby Isaac }
209*028944aaSToby Isaac 
210*028944aaSToby Isaac PetscErrorCode DMFieldShellSetEvaluateFV(DMField field, PetscErrorCode (*evaluateFV)(DMField,IS,PetscDataType,void*,void*,void*))
211*028944aaSToby Isaac {
212*028944aaSToby Isaac   PetscFunctionBegin;
213*028944aaSToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
214*028944aaSToby Isaac   field->ops->evaluateFV = evaluateFV;
215*028944aaSToby Isaac   PetscFunctionReturn(0);
216*028944aaSToby Isaac }
217*028944aaSToby Isaac 
218*028944aaSToby Isaac PetscErrorCode DMFieldShellSetGetFEInvariance(DMField field, PetscErrorCode (*getFEInvariance)(DMField,IS,PetscBool*,PetscBool*,PetscBool*))
219*028944aaSToby Isaac {
220*028944aaSToby Isaac   PetscFunctionBegin;
221*028944aaSToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
222*028944aaSToby Isaac   field->ops->getFEInvariance = getFEInvariance;
223*028944aaSToby Isaac   PetscFunctionReturn(0);
224*028944aaSToby Isaac }
225*028944aaSToby Isaac 
226*028944aaSToby Isaac PetscErrorCode DMFieldShellSetCreateDefaultQuadrature(DMField field, PetscErrorCode (*createDefaultQuadrature)(DMField,IS,PetscQuadrature*))
227*028944aaSToby Isaac {
228*028944aaSToby Isaac   PetscFunctionBegin;
229*028944aaSToby Isaac   PetscValidHeaderSpecific(field,DMFIELD_CLASSID,1);
230*028944aaSToby Isaac   field->ops->createDefaultQuadrature = createDefaultQuadrature;
231*028944aaSToby Isaac   PetscFunctionReturn(0);
232*028944aaSToby Isaac }
233*028944aaSToby Isaac 
234*028944aaSToby Isaac static PetscErrorCode DMFieldInitialize_Shell(DMField field)
235*028944aaSToby Isaac {
236*028944aaSToby Isaac   PetscFunctionBegin;
237*028944aaSToby Isaac   field->ops->destroy                 = DMFieldDestroy_Shell;
238*028944aaSToby Isaac   field->ops->evaluate                = NULL;
239*028944aaSToby Isaac   field->ops->evaluateFE              = DMFieldShellEvaluateFEDefault;
240*028944aaSToby Isaac   field->ops->evaluateFV              = DMFieldShellEvaluateFVDefault;
241*028944aaSToby Isaac   field->ops->getFEInvariance         = NULL;
242*028944aaSToby Isaac   field->ops->createDefaultQuadrature = NULL;
243*028944aaSToby Isaac   field->ops->view                    = NULL;
244*028944aaSToby Isaac   PetscFunctionReturn(0);
245*028944aaSToby Isaac }
246*028944aaSToby Isaac 
247*028944aaSToby Isaac PETSC_INTERN PetscErrorCode DMFieldCreate_Shell(DMField field)
248*028944aaSToby Isaac {
249*028944aaSToby Isaac   DMField_Shell *shell;
250*028944aaSToby Isaac   PetscErrorCode ierr;
251*028944aaSToby Isaac 
252*028944aaSToby Isaac   PetscFunctionBegin;
253*028944aaSToby Isaac   ierr = PetscNewLog(field,&shell);CHKERRQ(ierr);
254*028944aaSToby Isaac   field->data = shell;
255*028944aaSToby Isaac   ierr = DMFieldInitialize_Shell(field);CHKERRQ(ierr);
256*028944aaSToby Isaac   PetscFunctionReturn(0);
257*028944aaSToby Isaac }
258*028944aaSToby Isaac 
259*028944aaSToby Isaac PetscErrorCode DMFieldCreateShell(DM dm, PetscInt numComponents, DMFieldContinuity continuity, void *ctx, DMField *field)
260*028944aaSToby Isaac {
261*028944aaSToby Isaac   DMField        b;
262*028944aaSToby Isaac   DMField_Shell  *shell;
263*028944aaSToby Isaac   PetscErrorCode ierr;
264*028944aaSToby Isaac 
265*028944aaSToby Isaac   PetscFunctionBegin;
266*028944aaSToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
267*028944aaSToby Isaac   if (ctx) PetscValidPointer(ctx, 4);
268*028944aaSToby Isaac   PetscValidPointer(field, 5);
269*028944aaSToby Isaac   ierr = DMFieldCreate(dm,numComponents,continuity,&b);CHKERRQ(ierr);
270*028944aaSToby Isaac   ierr = DMFieldSetType(b,DMFIELDSHELL);CHKERRQ(ierr);
271*028944aaSToby Isaac   shell = (DMField_Shell *) b->data;
272*028944aaSToby Isaac   shell->ctx = ctx;
273*028944aaSToby Isaac   *field = b;
274*028944aaSToby Isaac   PetscFunctionReturn(0);
275*028944aaSToby Isaac }
276