xref: /petsc/src/dm/impls/plex/plexfem.c (revision d7ddef95b6c504adce762f22d2858ae4704940f7)
1cb1e1211SMatthew G Knepley #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2cb1e1211SMatthew G Knepley 
30f2d7e86SMatthew G. Knepley #include <petsc-private/petscfeimpl.h>
415496722SMatthew G. Knepley #include <petsc-private/petscfvimpl.h>
5a0845e3aSMatthew G. Knepley 
6cb1e1211SMatthew G Knepley #undef __FUNCT__
7cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexGetScale"
8cb1e1211SMatthew G Knepley PetscErrorCode DMPlexGetScale(DM dm, PetscUnit unit, PetscReal *scale)
9cb1e1211SMatthew G Knepley {
10cb1e1211SMatthew G Knepley   DM_Plex *mesh = (DM_Plex*) dm->data;
11cb1e1211SMatthew G Knepley 
12cb1e1211SMatthew G Knepley   PetscFunctionBegin;
13cb1e1211SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
14cb1e1211SMatthew G Knepley   PetscValidPointer(scale, 3);
15cb1e1211SMatthew G Knepley   *scale = mesh->scale[unit];
16cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
17cb1e1211SMatthew G Knepley }
18cb1e1211SMatthew G Knepley 
19cb1e1211SMatthew G Knepley #undef __FUNCT__
20cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexSetScale"
21cb1e1211SMatthew G Knepley PetscErrorCode DMPlexSetScale(DM dm, PetscUnit unit, PetscReal scale)
22cb1e1211SMatthew G Knepley {
23cb1e1211SMatthew G Knepley   DM_Plex *mesh = (DM_Plex*) dm->data;
24cb1e1211SMatthew G Knepley 
25cb1e1211SMatthew G Knepley   PetscFunctionBegin;
26cb1e1211SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
27cb1e1211SMatthew G Knepley   mesh->scale[unit] = scale;
28cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
29cb1e1211SMatthew G Knepley }
30cb1e1211SMatthew G Knepley 
31cb1e1211SMatthew G Knepley PETSC_STATIC_INLINE PetscInt epsilon(PetscInt i, PetscInt j, PetscInt k)
32cb1e1211SMatthew G Knepley {
33cb1e1211SMatthew G Knepley   switch (i) {
34cb1e1211SMatthew G Knepley   case 0:
35cb1e1211SMatthew G Knepley     switch (j) {
36cb1e1211SMatthew G Knepley     case 0: return 0;
37cb1e1211SMatthew G Knepley     case 1:
38cb1e1211SMatthew G Knepley       switch (k) {
39cb1e1211SMatthew G Knepley       case 0: return 0;
40cb1e1211SMatthew G Knepley       case 1: return 0;
41cb1e1211SMatthew G Knepley       case 2: return 1;
42cb1e1211SMatthew G Knepley       }
43cb1e1211SMatthew G Knepley     case 2:
44cb1e1211SMatthew G Knepley       switch (k) {
45cb1e1211SMatthew G Knepley       case 0: return 0;
46cb1e1211SMatthew G Knepley       case 1: return -1;
47cb1e1211SMatthew G Knepley       case 2: return 0;
48cb1e1211SMatthew G Knepley       }
49cb1e1211SMatthew G Knepley     }
50cb1e1211SMatthew G Knepley   case 1:
51cb1e1211SMatthew G Knepley     switch (j) {
52cb1e1211SMatthew G Knepley     case 0:
53cb1e1211SMatthew G Knepley       switch (k) {
54cb1e1211SMatthew G Knepley       case 0: return 0;
55cb1e1211SMatthew G Knepley       case 1: return 0;
56cb1e1211SMatthew G Knepley       case 2: return -1;
57cb1e1211SMatthew G Knepley       }
58cb1e1211SMatthew G Knepley     case 1: return 0;
59cb1e1211SMatthew G Knepley     case 2:
60cb1e1211SMatthew G Knepley       switch (k) {
61cb1e1211SMatthew G Knepley       case 0: return 1;
62cb1e1211SMatthew G Knepley       case 1: return 0;
63cb1e1211SMatthew G Knepley       case 2: return 0;
64cb1e1211SMatthew G Knepley       }
65cb1e1211SMatthew G Knepley     }
66cb1e1211SMatthew G Knepley   case 2:
67cb1e1211SMatthew G Knepley     switch (j) {
68cb1e1211SMatthew G Knepley     case 0:
69cb1e1211SMatthew G Knepley       switch (k) {
70cb1e1211SMatthew G Knepley       case 0: return 0;
71cb1e1211SMatthew G Knepley       case 1: return 1;
72cb1e1211SMatthew G Knepley       case 2: return 0;
73cb1e1211SMatthew G Knepley       }
74cb1e1211SMatthew G Knepley     case 1:
75cb1e1211SMatthew G Knepley       switch (k) {
76cb1e1211SMatthew G Knepley       case 0: return -1;
77cb1e1211SMatthew G Knepley       case 1: return 0;
78cb1e1211SMatthew G Knepley       case 2: return 0;
79cb1e1211SMatthew G Knepley       }
80cb1e1211SMatthew G Knepley     case 2: return 0;
81cb1e1211SMatthew G Knepley     }
82cb1e1211SMatthew G Knepley   }
83cb1e1211SMatthew G Knepley   return 0;
84cb1e1211SMatthew G Knepley }
85cb1e1211SMatthew G Knepley 
86cb1e1211SMatthew G Knepley #undef __FUNCT__
87cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexCreateRigidBody"
88cb1e1211SMatthew G Knepley /*@C
89cb1e1211SMatthew G Knepley   DMPlexCreateRigidBody - create rigid body modes from coordinates
90cb1e1211SMatthew G Knepley 
91cb1e1211SMatthew G Knepley   Collective on DM
92cb1e1211SMatthew G Knepley 
93cb1e1211SMatthew G Knepley   Input Arguments:
94cb1e1211SMatthew G Knepley + dm - the DM
95cb1e1211SMatthew G Knepley . section - the local section associated with the rigid field, or NULL for the default section
96cb1e1211SMatthew G Knepley - globalSection - the global section associated with the rigid field, or NULL for the default section
97cb1e1211SMatthew G Knepley 
98cb1e1211SMatthew G Knepley   Output Argument:
99cb1e1211SMatthew G Knepley . sp - the null space
100cb1e1211SMatthew G Knepley 
101cb1e1211SMatthew G Knepley   Note: This is necessary to take account of Dirichlet conditions on the displacements
102cb1e1211SMatthew G Knepley 
103cb1e1211SMatthew G Knepley   Level: advanced
104cb1e1211SMatthew G Knepley 
105cb1e1211SMatthew G Knepley .seealso: MatNullSpaceCreate()
106cb1e1211SMatthew G Knepley @*/
107cb1e1211SMatthew G Knepley PetscErrorCode DMPlexCreateRigidBody(DM dm, PetscSection section, PetscSection globalSection, MatNullSpace *sp)
108cb1e1211SMatthew G Knepley {
109cb1e1211SMatthew G Knepley   MPI_Comm       comm;
110cb1e1211SMatthew G Knepley   Vec            coordinates, localMode, mode[6];
111cb1e1211SMatthew G Knepley   PetscSection   coordSection;
112cb1e1211SMatthew G Knepley   PetscScalar   *coords;
113cb1e1211SMatthew G Knepley   PetscInt       dim, vStart, vEnd, v, n, m, d, i, j;
114cb1e1211SMatthew G Knepley   PetscErrorCode ierr;
115cb1e1211SMatthew G Knepley 
116cb1e1211SMatthew G Knepley   PetscFunctionBegin;
117cb1e1211SMatthew G Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
118c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
119cb1e1211SMatthew G Knepley   if (dim == 1) {
120cb1e1211SMatthew G Knepley     ierr = MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, sp);CHKERRQ(ierr);
121cb1e1211SMatthew G Knepley     PetscFunctionReturn(0);
122cb1e1211SMatthew G Knepley   }
123cb1e1211SMatthew G Knepley   if (!section)       {ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);}
124cb1e1211SMatthew G Knepley   if (!globalSection) {ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);}
125cb1e1211SMatthew G Knepley   ierr = PetscSectionGetConstrainedStorageSize(globalSection, &n);CHKERRQ(ierr);
126cb1e1211SMatthew G Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
12769d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr);
128cb1e1211SMatthew G Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
129cb1e1211SMatthew G Knepley   m    = (dim*(dim+1))/2;
130cb1e1211SMatthew G Knepley   ierr = VecCreate(comm, &mode[0]);CHKERRQ(ierr);
131cb1e1211SMatthew G Knepley   ierr = VecSetSizes(mode[0], n, PETSC_DETERMINE);CHKERRQ(ierr);
132cb1e1211SMatthew G Knepley   ierr = VecSetUp(mode[0]);CHKERRQ(ierr);
133cb1e1211SMatthew G Knepley   for (i = 1; i < m; ++i) {ierr = VecDuplicate(mode[0], &mode[i]);CHKERRQ(ierr);}
134cb1e1211SMatthew G Knepley   /* Assume P1 */
135cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localMode);CHKERRQ(ierr);
136cb1e1211SMatthew G Knepley   for (d = 0; d < dim; ++d) {
137cb1e1211SMatthew G Knepley     PetscScalar values[3] = {0.0, 0.0, 0.0};
138cb1e1211SMatthew G Knepley 
139cb1e1211SMatthew G Knepley     values[d] = 1.0;
140cb1e1211SMatthew G Knepley     ierr      = VecSet(localMode, 0.0);CHKERRQ(ierr);
141cb1e1211SMatthew G Knepley     for (v = vStart; v < vEnd; ++v) {
142cb1e1211SMatthew G Knepley       ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr);
143cb1e1211SMatthew G Knepley     }
144cb1e1211SMatthew G Knepley     ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
145cb1e1211SMatthew G Knepley     ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
146cb1e1211SMatthew G Knepley   }
147cb1e1211SMatthew G Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
148cb1e1211SMatthew G Knepley   for (d = dim; d < dim*(dim+1)/2; ++d) {
149cb1e1211SMatthew G Knepley     PetscInt i, j, k = dim > 2 ? d - dim : d;
150cb1e1211SMatthew G Knepley 
151cb1e1211SMatthew G Knepley     ierr = VecSet(localMode, 0.0);CHKERRQ(ierr);
152cb1e1211SMatthew G Knepley     for (v = vStart; v < vEnd; ++v) {
153cb1e1211SMatthew G Knepley       PetscScalar values[3] = {0.0, 0.0, 0.0};
154cb1e1211SMatthew G Knepley       PetscInt    off;
155cb1e1211SMatthew G Knepley 
156cb1e1211SMatthew G Knepley       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
157cb1e1211SMatthew G Knepley       for (i = 0; i < dim; ++i) {
158cb1e1211SMatthew G Knepley         for (j = 0; j < dim; ++j) {
159cb1e1211SMatthew G Knepley           values[j] += epsilon(i, j, k)*PetscRealPart(coords[off+i]);
160cb1e1211SMatthew G Knepley         }
161cb1e1211SMatthew G Knepley       }
162cb1e1211SMatthew G Knepley       ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr);
163cb1e1211SMatthew G Knepley     }
164cb1e1211SMatthew G Knepley     ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
165cb1e1211SMatthew G Knepley     ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr);
166cb1e1211SMatthew G Knepley   }
167cb1e1211SMatthew G Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
168cb1e1211SMatthew G Knepley   ierr = DMRestoreLocalVector(dm, &localMode);CHKERRQ(ierr);
169cb1e1211SMatthew G Knepley   for (i = 0; i < dim; ++i) {ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);}
170cb1e1211SMatthew G Knepley   /* Orthonormalize system */
171cb1e1211SMatthew G Knepley   for (i = dim; i < m; ++i) {
172cb1e1211SMatthew G Knepley     PetscScalar dots[6];
173cb1e1211SMatthew G Knepley 
174cb1e1211SMatthew G Knepley     ierr = VecMDot(mode[i], i, mode, dots);CHKERRQ(ierr);
175cb1e1211SMatthew G Knepley     for (j = 0; j < i; ++j) dots[j] *= -1.0;
176cb1e1211SMatthew G Knepley     ierr = VecMAXPY(mode[i], i, dots, mode);CHKERRQ(ierr);
177cb1e1211SMatthew G Knepley     ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);
178cb1e1211SMatthew G Knepley   }
179cb1e1211SMatthew G Knepley   ierr = MatNullSpaceCreate(comm, PETSC_FALSE, m, mode, sp);CHKERRQ(ierr);
180cb1e1211SMatthew G Knepley   for (i = 0; i< m; ++i) {ierr = VecDestroy(&mode[i]);CHKERRQ(ierr);}
181cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
182cb1e1211SMatthew G Knepley }
183cb1e1211SMatthew G Knepley 
184cb1e1211SMatthew G Knepley #undef __FUNCT__
185a18a7fb9SMatthew G. Knepley #define __FUNCT__ "DMPlexProjectFunctionLabelLocal"
186bf3434eeSMatthew G. Knepley PetscErrorCode DMPlexProjectFunctionLabelLocal(DM dm, DMLabel label, PetscInt numIds, const PetscInt ids[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
187a18a7fb9SMatthew G. Knepley {
188bf3434eeSMatthew G. Knepley   PetscFE         fe;
189a18a7fb9SMatthew G. Knepley   PetscDualSpace *sp;
190a18a7fb9SMatthew G. Knepley   PetscSection    section;
191a18a7fb9SMatthew G. Knepley   PetscScalar    *values;
192ad96f515SMatthew G. Knepley   PetscBool      *fieldActive;
193a18a7fb9SMatthew G. Knepley   PetscInt        numFields, numComp, dim, spDim, totDim = 0, numValues, cStart, cEnd, f, d, v, i, comp;
194a18a7fb9SMatthew G. Knepley   PetscErrorCode  ierr;
195a18a7fb9SMatthew G. Knepley 
196a18a7fb9SMatthew G. Knepley   PetscFunctionBegin;
197c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
198a18a7fb9SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
199a18a7fb9SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
200e1d0b1adSMatthew G. Knepley   ierr = PetscMalloc1(numFields,&sp);CHKERRQ(ierr);
201a18a7fb9SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {
202bf3434eeSMatthew G. Knepley     ierr = DMGetField(dm, f, (PetscObject *) &fe);CHKERRQ(ierr);
203bf3434eeSMatthew G. Knepley     ierr = PetscFEGetDualSpace(fe, &sp[f]);CHKERRQ(ierr);
204bf3434eeSMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr);
205a18a7fb9SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
206a18a7fb9SMatthew G. Knepley     totDim += spDim*numComp;
207a18a7fb9SMatthew G. Knepley   }
208a18a7fb9SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
209a18a7fb9SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
210a18a7fb9SMatthew G. Knepley   if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim);
211a18a7fb9SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
212ad96f515SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr);
213ad96f515SMatthew G. Knepley   for (f = 0; f < numFields; ++f) fieldActive[f] = funcs[f] ? PETSC_TRUE : PETSC_FALSE;
214a18a7fb9SMatthew G. Knepley   for (i = 0; i < numIds; ++i) {
215a18a7fb9SMatthew G. Knepley     IS              pointIS;
216a18a7fb9SMatthew G. Knepley     const PetscInt *points;
217a18a7fb9SMatthew G. Knepley     PetscInt        n, p;
218a18a7fb9SMatthew G. Knepley 
219a18a7fb9SMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr);
220a18a7fb9SMatthew G. Knepley     ierr = ISGetLocalSize(pointIS, &n);CHKERRQ(ierr);
221a18a7fb9SMatthew G. Knepley     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
222a18a7fb9SMatthew G. Knepley     for (p = 0; p < n; ++p) {
223a18a7fb9SMatthew G. Knepley       const PetscInt  point = points[p];
224e1d0b1adSMatthew G. Knepley       PetscFECellGeom geom;
225a18a7fb9SMatthew G. Knepley 
226a18a7fb9SMatthew G. Knepley       if ((point < cStart) || (point >= cEnd)) continue;
227e1d0b1adSMatthew G. Knepley       ierr = DMPlexComputeCellGeometryFEM(dm, point, NULL, geom.v0, geom.J, NULL, &geom.detJ);CHKERRQ(ierr);
228a18a7fb9SMatthew G. Knepley       for (f = 0, v = 0; f < numFields; ++f) {
229a18a7fb9SMatthew G. Knepley         void * const ctx = ctxs ? ctxs[f] : NULL;
230bf3434eeSMatthew G. Knepley 
231bf3434eeSMatthew G. Knepley         ierr = DMGetField(dm, f, (PetscObject *) &fe);CHKERRQ(ierr);
232bf3434eeSMatthew G. Knepley         ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr);
233a18a7fb9SMatthew G. Knepley         ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
234a18a7fb9SMatthew G. Knepley         for (d = 0; d < spDim; ++d) {
235a18a7fb9SMatthew G. Knepley           if (funcs[f]) {
236e1d0b1adSMatthew G. Knepley             ierr = PetscDualSpaceApply(sp[f], d, &geom, numComp, funcs[f], ctx, &values[v]);CHKERRQ(ierr);
237a18a7fb9SMatthew G. Knepley           } else {
238a18a7fb9SMatthew G. Knepley             for (comp = 0; comp < numComp; ++comp) values[v+comp] = 0.0;
239a18a7fb9SMatthew G. Knepley           }
240a18a7fb9SMatthew G. Knepley           v += numComp;
241a18a7fb9SMatthew G. Knepley         }
242a18a7fb9SMatthew G. Knepley       }
243ad96f515SMatthew G. Knepley       ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, point, values, mode);CHKERRQ(ierr);
244a18a7fb9SMatthew G. Knepley     }
245a18a7fb9SMatthew G. Knepley     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
246a18a7fb9SMatthew G. Knepley     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
247a18a7fb9SMatthew G. Knepley   }
248a18a7fb9SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
249ad96f515SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr);
250e1d0b1adSMatthew G. Knepley   ierr = PetscFree(sp);CHKERRQ(ierr);
251a18a7fb9SMatthew G. Knepley   PetscFunctionReturn(0);
252a18a7fb9SMatthew G. Knepley }
253a18a7fb9SMatthew G. Knepley 
254a18a7fb9SMatthew G. Knepley #undef __FUNCT__
255cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexProjectFunctionLocal"
2560f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexProjectFunctionLocal(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX)
257cb1e1211SMatthew G Knepley {
25872f94c41SMatthew G. Knepley   PetscDualSpace *sp;
25915496722SMatthew G. Knepley   PetscInt       *numComp;
26072f94c41SMatthew G. Knepley   PetscSection    section;
26172f94c41SMatthew G. Knepley   PetscScalar    *values;
26215496722SMatthew G. Knepley   PetscInt        numFields, dim, spDim, totDim = 0, numValues, cStart, cEnd, c, f, d, v, comp;
263cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
264cb1e1211SMatthew G Knepley 
265cb1e1211SMatthew G Knepley   PetscFunctionBegin;
266cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
26772f94c41SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
26815496722SMatthew G. Knepley   ierr = PetscMalloc2(numFields, &sp, numFields, &numComp);CHKERRQ(ierr);
26972f94c41SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {
27015496722SMatthew G. Knepley     PetscObject  obj;
27115496722SMatthew G. Knepley     PetscClassId id;
2720f2d7e86SMatthew G. Knepley 
27315496722SMatthew G. Knepley     ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr);
27415496722SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
27515496722SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
27615496722SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
27715496722SMatthew G. Knepley 
27815496722SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &numComp[f]);CHKERRQ(ierr);
2790f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDualSpace(fe, &sp[f]);CHKERRQ(ierr);
28015496722SMatthew G. Knepley       ierr = PetscObjectReference((PetscObject) sp[f]);CHKERRQ(ierr);
28115496722SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
28215496722SMatthew G. Knepley       PetscFV         fv = (PetscFV) obj;
28315496722SMatthew G. Knepley       PetscQuadrature q;
28415496722SMatthew G. Knepley 
28515496722SMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &numComp[f]);CHKERRQ(ierr);
28615496722SMatthew G. Knepley       ierr = PetscFVGetQuadrature(fv, &q);CHKERRQ(ierr);
28715496722SMatthew G. Knepley       ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject) fv), &sp[f]);CHKERRQ(ierr);
28815496722SMatthew G. Knepley       ierr = PetscDualSpaceSetDM(sp[f], dm);CHKERRQ(ierr);
28915496722SMatthew G. Knepley       ierr = PetscDualSpaceSetType(sp[f], PETSCDUALSPACESIMPLE);CHKERRQ(ierr);
29015496722SMatthew G. Knepley       ierr = PetscDualSpaceSimpleSetDimension(sp[f], 1);CHKERRQ(ierr);
29115496722SMatthew G. Knepley       ierr = PetscDualSpaceSimpleSetFunctional(sp[f], 0, q);CHKERRQ(ierr);
29215496722SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f);
29372f94c41SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
29415496722SMatthew G. Knepley     totDim += spDim*numComp[f];
295cb1e1211SMatthew G Knepley   }
296c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
29772f94c41SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
29872f94c41SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
29972f94c41SMatthew G. Knepley   if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim);
30072f94c41SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
30172f94c41SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
302e1d0b1adSMatthew G. Knepley     PetscFECellGeom geom;
303cb1e1211SMatthew G Knepley 
304e1d0b1adSMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, geom.v0, geom.J, NULL, &geom.detJ);CHKERRQ(ierr);
30572f94c41SMatthew G. Knepley     for (f = 0, v = 0; f < numFields; ++f) {
306c110b1eeSGeoffrey Irving       void *const ctx = ctxs ? ctxs[f] : NULL;
3070f2d7e86SMatthew G. Knepley 
30872f94c41SMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr);
30972f94c41SMatthew G. Knepley       for (d = 0; d < spDim; ++d) {
310120386c5SMatthew G. Knepley         if (funcs[f]) {
311e1d0b1adSMatthew G. Knepley           ierr = PetscDualSpaceApply(sp[f], d, &geom, numComp[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr);
312120386c5SMatthew G. Knepley         } else {
31315496722SMatthew G. Knepley           for (comp = 0; comp < numComp[f]; ++comp) values[v+comp] = 0.0;
314120386c5SMatthew G. Knepley         }
31515496722SMatthew G. Knepley         v += numComp[f];
316cb1e1211SMatthew G Knepley       }
317cb1e1211SMatthew G Knepley     }
31872f94c41SMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr);
319cb1e1211SMatthew G Knepley   }
32072f94c41SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
32115496722SMatthew G. Knepley   for (f = 0; f < numFields; ++f) {ierr = PetscDualSpaceDestroy(&sp[f]);CHKERRQ(ierr);}
32215496722SMatthew G. Knepley   ierr = PetscFree2(sp, numComp);CHKERRQ(ierr);
323cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
324cb1e1211SMatthew G Knepley }
325cb1e1211SMatthew G Knepley 
326cb1e1211SMatthew G Knepley #undef __FUNCT__
327cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexProjectFunction"
328cb1e1211SMatthew G Knepley /*@C
329cb1e1211SMatthew G Knepley   DMPlexProjectFunction - This projects the given function into the function space provided.
330cb1e1211SMatthew G Knepley 
331cb1e1211SMatthew G Knepley   Input Parameters:
332cb1e1211SMatthew G Knepley + dm      - The DM
33372f94c41SMatthew G. Knepley . funcs   - The coordinate functions to evaluate, one per field
334c110b1eeSGeoffrey Irving . ctxs    - Optional array of contexts to pass to each coordinate function.  ctxs itself may be null.
335cb1e1211SMatthew G Knepley - mode    - The insertion mode for values
336cb1e1211SMatthew G Knepley 
337cb1e1211SMatthew G Knepley   Output Parameter:
338cb1e1211SMatthew G Knepley . X - vector
339cb1e1211SMatthew G Knepley 
340cb1e1211SMatthew G Knepley   Level: developer
341cb1e1211SMatthew G Knepley 
342878cb397SSatish Balay .seealso: DMPlexComputeL2Diff()
343878cb397SSatish Balay @*/
3440f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexProjectFunction(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec X)
345cb1e1211SMatthew G Knepley {
346cb1e1211SMatthew G Knepley   Vec            localX;
347cb1e1211SMatthew G Knepley   PetscErrorCode ierr;
348cb1e1211SMatthew G Knepley 
349cb1e1211SMatthew G Knepley   PetscFunctionBegin;
3509a800dd8SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
351cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
3520f2d7e86SMatthew G. Knepley   ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, mode, localX);CHKERRQ(ierr);
353cb1e1211SMatthew G Knepley   ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr);
354cb1e1211SMatthew G Knepley   ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr);
355cb1e1211SMatthew G Knepley   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
356cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
357cb1e1211SMatthew G Knepley }
358cb1e1211SMatthew G Knepley 
35955f2e967SMatthew G. Knepley #undef __FUNCT__
3600f2d7e86SMatthew G. Knepley #define __FUNCT__ "DMPlexProjectFieldLocal"
3613bc3b0a0SMatthew G. Knepley PetscErrorCode DMPlexProjectFieldLocal(DM dm, Vec localU, void (**funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal [], PetscScalar []), InsertMode mode, Vec localX)
3620f2d7e86SMatthew G. Knepley {
3630f2d7e86SMatthew G. Knepley   DM                dmAux;
3642764a2aaSMatthew G. Knepley   PetscDS           prob, probAux;
3650f2d7e86SMatthew G. Knepley   Vec               A;
366326413afSMatthew G. Knepley   PetscSection      section, sectionAux;
367326413afSMatthew G. Knepley   PetscScalar      *values, *u, *u_x, *a, *a_x;
3680f2d7e86SMatthew G. Knepley   PetscReal        *x, *v0, *J, *invJ, detJ, **basisField, **basisFieldDer, **basisFieldAux, **basisFieldDerAux;
369326413afSMatthew G. Knepley   PetscInt          Nf, dim, spDim, totDim, numValues, cStart, cEnd, c, f, d, v, comp;
3700f2d7e86SMatthew G. Knepley   PetscErrorCode    ierr;
3710f2d7e86SMatthew G. Knepley 
3720f2d7e86SMatthew G. Knepley   PetscFunctionBegin;
3732764a2aaSMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
374c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
3750f2d7e86SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
3760f2d7e86SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
3770f2d7e86SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
3782764a2aaSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
3792764a2aaSMatthew G. Knepley   ierr = PetscDSGetTabulation(prob, &basisField, &basisFieldDer);CHKERRQ(ierr);
3802764a2aaSMatthew G. Knepley   ierr = PetscDSGetEvaluationArrays(prob, &u, NULL, &u_x);CHKERRQ(ierr);
3812764a2aaSMatthew G. Knepley   ierr = PetscDSGetRefCoordArrays(prob, &x, NULL);CHKERRQ(ierr);
3820f2d7e86SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
3830f2d7e86SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
3840f2d7e86SMatthew G. Knepley   if (dmAux) {
3852764a2aaSMatthew G. Knepley     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
386326413afSMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
3872764a2aaSMatthew G. Knepley     ierr = PetscDSGetTabulation(prob, &basisFieldAux, &basisFieldDerAux);CHKERRQ(ierr);
3882764a2aaSMatthew G. Knepley     ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);CHKERRQ(ierr);
3890f2d7e86SMatthew G. Knepley   }
390*d7ddef95SMatthew G. Knepley   ierr = DMPlexInsertBoundaryValues(dm, localU, 0.0, NULL, NULL, NULL);CHKERRQ(ierr);
3910f2d7e86SMatthew G. Knepley   ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr);
3920f2d7e86SMatthew G. Knepley   if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim);
3930f2d7e86SMatthew G. Knepley   ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
3940f2d7e86SMatthew G. Knepley   ierr = PetscMalloc3(dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
3950f2d7e86SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
396326413afSMatthew G. Knepley     PetscScalar *coefficients = NULL, *coefficientsAux = NULL;
397326413afSMatthew G. Knepley 
3988e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
399326413afSMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, localU, c, NULL, &coefficients);CHKERRQ(ierr);
400326413afSMatthew G. Knepley     if (dmAux) {ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &coefficientsAux);CHKERRQ(ierr);}
4010f2d7e86SMatthew G. Knepley     for (f = 0, v = 0; f < Nf; ++f) {
4023113607cSMatthew G. Knepley       PetscFE        fe;
4033113607cSMatthew G. Knepley       PetscDualSpace sp;
4043113607cSMatthew G. Knepley       PetscInt       Ncf;
4053113607cSMatthew G. Knepley 
4062764a2aaSMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
4073113607cSMatthew G. Knepley       ierr = PetscFEGetDualSpace(fe, &sp);CHKERRQ(ierr);
4083113607cSMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Ncf);CHKERRQ(ierr);
4093113607cSMatthew G. Knepley       ierr = PetscDualSpaceGetDimension(sp, &spDim);CHKERRQ(ierr);
4100f2d7e86SMatthew G. Knepley       for (d = 0; d < spDim; ++d) {
4110f2d7e86SMatthew G. Knepley         PetscQuadrature  quad;
4120f2d7e86SMatthew G. Knepley         const PetscReal *points, *weights;
4130f2d7e86SMatthew G. Knepley         PetscInt         numPoints, q;
4140f2d7e86SMatthew G. Knepley 
4150f2d7e86SMatthew G. Knepley         if (funcs[f]) {
4163113607cSMatthew G. Knepley           ierr = PetscDualSpaceGetFunctional(sp, d, &quad);CHKERRQ(ierr);
4170f2d7e86SMatthew G. Knepley           ierr = PetscQuadratureGetData(quad, NULL, &numPoints, &points, &weights);CHKERRQ(ierr);
4183113607cSMatthew G. Knepley           ierr = PetscFEGetTabulation(fe, numPoints, points, &basisField[f], &basisFieldDer[f], NULL);CHKERRQ(ierr);
4190f2d7e86SMatthew G. Knepley           for (q = 0; q < numPoints; ++q) {
4200f2d7e86SMatthew G. Knepley             CoordinatesRefToReal(dim, dim, v0, J, &points[q*dim], x);
421326413afSMatthew G. Knepley             ierr = EvaluateFieldJets(prob,    PETSC_FALSE, q, invJ, coefficients,    NULL, u, u_x, NULL);CHKERRQ(ierr);
422326413afSMatthew G. Knepley             ierr = EvaluateFieldJets(probAux, PETSC_FALSE, q, invJ, coefficientsAux, NULL, a, a_x, NULL);CHKERRQ(ierr);
4233bc3b0a0SMatthew G. Knepley             (*funcs[f])(u, NULL, u_x, a, NULL, a_x, x, &values[v]);
4240f2d7e86SMatthew G. Knepley           }
4253113607cSMatthew G. Knepley           ierr = PetscFERestoreTabulation(fe, numPoints, points, &basisField[f], &basisFieldDer[f], NULL);CHKERRQ(ierr);
4260f2d7e86SMatthew G. Knepley         } else {
4270f2d7e86SMatthew G. Knepley           for (comp = 0; comp < Ncf; ++comp) values[v+comp] = 0.0;
4280f2d7e86SMatthew G. Knepley         }
4290f2d7e86SMatthew G. Knepley         v += Ncf;
4300f2d7e86SMatthew G. Knepley       }
4310f2d7e86SMatthew G. Knepley     }
432326413afSMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, localU, c, NULL, &coefficients);CHKERRQ(ierr);
433326413afSMatthew G. Knepley     if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &coefficientsAux);CHKERRQ(ierr);}
4340f2d7e86SMatthew G. Knepley     ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr);
4350f2d7e86SMatthew G. Knepley   }
4360f2d7e86SMatthew G. Knepley   ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr);
4370f2d7e86SMatthew G. Knepley   ierr = PetscFree3(v0,J,invJ);CHKERRQ(ierr);
4380f2d7e86SMatthew G. Knepley   PetscFunctionReturn(0);
4390f2d7e86SMatthew G. Knepley }
4400f2d7e86SMatthew G. Knepley 
4410f2d7e86SMatthew G. Knepley #undef __FUNCT__
4420f2d7e86SMatthew G. Knepley #define __FUNCT__ "DMPlexProjectField"
4430f2d7e86SMatthew G. Knepley /*@C
4440f2d7e86SMatthew G. Knepley   DMPlexProjectField - This projects the given function of the fields into the function space provided.
4450f2d7e86SMatthew G. Knepley 
4460f2d7e86SMatthew G. Knepley   Input Parameters:
4470f2d7e86SMatthew G. Knepley + dm      - The DM
4480f2d7e86SMatthew G. Knepley . U       - The input field vector
4490f2d7e86SMatthew G. Knepley . funcs   - The functions to evaluate, one per field
4500f2d7e86SMatthew G. Knepley - mode    - The insertion mode for values
4510f2d7e86SMatthew G. Knepley 
4520f2d7e86SMatthew G. Knepley   Output Parameter:
4530f2d7e86SMatthew G. Knepley . X       - The output vector
4540f2d7e86SMatthew G. Knepley 
4550f2d7e86SMatthew G. Knepley   Level: developer
4560f2d7e86SMatthew G. Knepley 
4570f2d7e86SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff()
4580f2d7e86SMatthew G. Knepley @*/
4593bc3b0a0SMatthew G. Knepley PetscErrorCode DMPlexProjectField(DM dm, Vec U, void (**funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal [], PetscScalar []), InsertMode mode, Vec X)
4600f2d7e86SMatthew G. Knepley {
4610f2d7e86SMatthew G. Knepley   Vec            localX, localU;
4620f2d7e86SMatthew G. Knepley   PetscErrorCode ierr;
4630f2d7e86SMatthew G. Knepley 
4640f2d7e86SMatthew G. Knepley   PetscFunctionBegin;
4650f2d7e86SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4660f2d7e86SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
4670f2d7e86SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localU);CHKERRQ(ierr);
4680f2d7e86SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, U, INSERT_VALUES, localU);CHKERRQ(ierr);
4690f2d7e86SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, U, INSERT_VALUES, localU);CHKERRQ(ierr);
4703113607cSMatthew G. Knepley   ierr = DMPlexProjectFieldLocal(dm, localU, funcs, mode, localX);CHKERRQ(ierr);
4710f2d7e86SMatthew G. Knepley   ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr);
4720f2d7e86SMatthew G. Knepley   ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr);
4730f2d7e86SMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
4740f2d7e86SMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &localU);CHKERRQ(ierr);
4750f2d7e86SMatthew G. Knepley   PetscFunctionReturn(0);
4760f2d7e86SMatthew G. Knepley }
4770f2d7e86SMatthew G. Knepley 
4780f2d7e86SMatthew G. Knepley #undef __FUNCT__
479*d7ddef95SMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValues_FEM_Internal"
480*d7ddef95SMatthew G. Knepley static PetscErrorCode DMPlexInsertBoundaryValues_FEM_Internal(DM dm, PetscInt field, DMLabel label, PetscInt numids, const PetscInt ids[], void (*func)(const PetscReal[], PetscScalar *, void *), void *ctx, Vec locX)
48155f2e967SMatthew G. Knepley {
48255f2e967SMatthew G. Knepley   void        (**funcs)(const PetscReal x[], PetscScalar *u, void *ctx);
48355f2e967SMatthew G. Knepley   void         **ctxs;
484*d7ddef95SMatthew G. Knepley   PetscInt       numFields;
485*d7ddef95SMatthew G. Knepley   PetscErrorCode ierr;
486*d7ddef95SMatthew G. Knepley 
487*d7ddef95SMatthew G. Knepley   PetscFunctionBegin;
488*d7ddef95SMatthew G. Knepley   ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr);
489*d7ddef95SMatthew G. Knepley   ierr = PetscCalloc2(numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr);
490*d7ddef95SMatthew G. Knepley   funcs[field] = func;
491*d7ddef95SMatthew G. Knepley   ctxs[field]  = ctx;
492*d7ddef95SMatthew G. Knepley   ierr = DMPlexProjectFunctionLabelLocal(dm, label, numids, ids, funcs, ctxs, INSERT_BC_VALUES, locX);CHKERRQ(ierr);
493*d7ddef95SMatthew G. Knepley   ierr = PetscFree2(funcs,ctxs);CHKERRQ(ierr);
494*d7ddef95SMatthew G. Knepley   PetscFunctionReturn(0);
495*d7ddef95SMatthew G. Knepley }
496*d7ddef95SMatthew G. Knepley 
497*d7ddef95SMatthew G. Knepley #undef __FUNCT__
498*d7ddef95SMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValues_FVM_Internal"
499*d7ddef95SMatthew G. Knepley static PetscErrorCode DMPlexInsertBoundaryValues_FVM_Internal(DM dm, PetscReal time, Vec faceGeometry, Vec cellGeometry, Vec Grad,
500*d7ddef95SMatthew G. Knepley                                                               PetscInt field, DMLabel label, PetscInt numids, const PetscInt ids[], PetscErrorCode (*func)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*), void *ctx, Vec locX)
501*d7ddef95SMatthew G. Knepley {
502*d7ddef95SMatthew G. Knepley   PetscFV            fv;
503*d7ddef95SMatthew G. Knepley   DM                 dmFace, dmCell, dmGrad;
504*d7ddef95SMatthew G. Knepley   const PetscScalar *facegeom, *cellgeom, *grad;
505*d7ddef95SMatthew G. Knepley   PetscScalar       *x, *fx;
506*d7ddef95SMatthew G. Knepley   PetscInt           dim, fStart, fEnd, pdim, i;
507*d7ddef95SMatthew G. Knepley   PetscErrorCode     ierr;
508*d7ddef95SMatthew G. Knepley 
509*d7ddef95SMatthew G. Knepley   PetscFunctionBegin;
510*d7ddef95SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
511*d7ddef95SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
512*d7ddef95SMatthew G. Knepley   ierr = DMGetField(dm, field, (PetscObject *) &fv);CHKERRQ(ierr);
513*d7ddef95SMatthew G. Knepley   ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr);
514*d7ddef95SMatthew G. Knepley   ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
515*d7ddef95SMatthew G. Knepley   ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr);
516*d7ddef95SMatthew G. Knepley   ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
517*d7ddef95SMatthew G. Knepley   if (Grad) {
518*d7ddef95SMatthew G. Knepley     ierr = VecGetDM(Grad, &dmGrad);CHKERRQ(ierr);
519*d7ddef95SMatthew G. Knepley     ierr = VecGetArrayRead(Grad, &grad);CHKERRQ(ierr);
520*d7ddef95SMatthew G. Knepley     ierr = PetscFVGetNumComponents(fv, &pdim);CHKERRQ(ierr);
521*d7ddef95SMatthew G. Knepley     ierr = DMGetWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr);
522*d7ddef95SMatthew G. Knepley   }
523*d7ddef95SMatthew G. Knepley   ierr = VecGetArray(locX, &x);CHKERRQ(ierr);
524*d7ddef95SMatthew G. Knepley   for (i = 0; i < numids; ++i) {
525*d7ddef95SMatthew G. Knepley     IS              faceIS;
526*d7ddef95SMatthew G. Knepley     const PetscInt *faces;
527*d7ddef95SMatthew G. Knepley     PetscInt        numFaces, f;
528*d7ddef95SMatthew G. Knepley 
529*d7ddef95SMatthew G. Knepley     ierr = DMLabelGetStratumIS(label, ids[i], &faceIS);CHKERRQ(ierr);
530*d7ddef95SMatthew G. Knepley     if (!faceIS) continue; /* No points with that id on this process */
531*d7ddef95SMatthew G. Knepley     ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr);
532*d7ddef95SMatthew G. Knepley     ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr);
533*d7ddef95SMatthew G. Knepley     for (f = 0; f < numFaces; ++f) {
534*d7ddef95SMatthew G. Knepley       const PetscInt         face = faces[f], *cells;
535*d7ddef95SMatthew G. Knepley       const PetscFVFaceGeom *fg;
536*d7ddef95SMatthew G. Knepley 
537*d7ddef95SMatthew G. Knepley       if ((face < fStart) || (face >= fEnd)) continue; /* Refinement adds non-faces to labels */
538*d7ddef95SMatthew G. Knepley       ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr);
539*d7ddef95SMatthew G. Knepley       ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr);
540*d7ddef95SMatthew G. Knepley       if (Grad) {
541*d7ddef95SMatthew G. Knepley         const PetscFVCellGeom *cg;
542*d7ddef95SMatthew G. Knepley         const PetscScalar     *cx, *cgrad;
543*d7ddef95SMatthew G. Knepley         PetscScalar           *xG;
544*d7ddef95SMatthew G. Knepley         PetscReal              dx[3];
545*d7ddef95SMatthew G. Knepley         PetscInt               d;
546*d7ddef95SMatthew G. Knepley 
547*d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cg);CHKERRQ(ierr);
548*d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dm, cells[0], x, &cx);CHKERRQ(ierr);
549*d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dmGrad, cells[0], grad, &cgrad);CHKERRQ(ierr);
550*d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRef(dm, cells[1], x, &xG);CHKERRQ(ierr);
551*d7ddef95SMatthew G. Knepley         DMPlex_WaxpyD_Internal(dim, -1, cg->centroid, fg->centroid, dx);
552*d7ddef95SMatthew G. Knepley         for (d = 0; d < pdim; ++d) fx[d] = cx[d] + DMPlex_DotD_Internal(dim, &cgrad[d*dim], dx);
553*d7ddef95SMatthew G. Knepley         ierr = (*func)(time, fg->centroid, fg->normal, fx, xG, ctx);CHKERRQ(ierr);
554*d7ddef95SMatthew G. Knepley       } else {
555*d7ddef95SMatthew G. Knepley         const PetscScalar *xI;
556*d7ddef95SMatthew G. Knepley         PetscScalar       *xG;
557*d7ddef95SMatthew G. Knepley 
558*d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRead(dm, cells[0], x, &xI);CHKERRQ(ierr);
559*d7ddef95SMatthew G. Knepley         ierr = DMPlexPointLocalRef(dm, cells[1], x, &xG);CHKERRQ(ierr);
560*d7ddef95SMatthew G. Knepley         ierr = (*func)(time, fg->centroid, fg->normal, xI, xG, ctx);CHKERRQ(ierr);
561*d7ddef95SMatthew G. Knepley       }
562*d7ddef95SMatthew G. Knepley     }
563*d7ddef95SMatthew G. Knepley     ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr);
564*d7ddef95SMatthew G. Knepley     ierr = ISDestroy(&faceIS);CHKERRQ(ierr);
565*d7ddef95SMatthew G. Knepley   }
566*d7ddef95SMatthew G. Knepley   ierr = VecRestoreArray(locX, &x);CHKERRQ(ierr);
567*d7ddef95SMatthew G. Knepley   if (Grad) {
568*d7ddef95SMatthew G. Knepley     ierr = DMRestoreWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr);
569*d7ddef95SMatthew G. Knepley     ierr = VecRestoreArrayRead(Grad, &grad);CHKERRQ(ierr);
570*d7ddef95SMatthew G. Knepley   }
571*d7ddef95SMatthew G. Knepley   ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);
572*d7ddef95SMatthew G. Knepley   ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr);
573*d7ddef95SMatthew G. Knepley   PetscFunctionReturn(0);
574*d7ddef95SMatthew G. Knepley }
575*d7ddef95SMatthew G. Knepley 
576*d7ddef95SMatthew G. Knepley #undef __FUNCT__
577*d7ddef95SMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValues"
578*d7ddef95SMatthew G. Knepley PetscErrorCode DMPlexInsertBoundaryValues(DM dm, Vec locX, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM)
579*d7ddef95SMatthew G. Knepley {
580*d7ddef95SMatthew G. Knepley   PetscInt       numBd, b;
58155f2e967SMatthew G. Knepley   PetscErrorCode ierr;
58255f2e967SMatthew G. Knepley 
58355f2e967SMatthew G. Knepley   PetscFunctionBegin;
58455f2e967SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
585*d7ddef95SMatthew G. Knepley   PetscValidHeaderSpecific(locX, VEC_CLASSID, 2);
586*d7ddef95SMatthew G. Knepley   if (faceGeomFVM) {PetscValidHeaderSpecific(faceGeomFVM, VEC_CLASSID, 4);}
587*d7ddef95SMatthew G. Knepley   if (cellGeomFVM) {PetscValidHeaderSpecific(cellGeomFVM, VEC_CLASSID, 5);}
588*d7ddef95SMatthew G. Knepley   if (gradFVM)     {PetscValidHeaderSpecific(gradFVM, VEC_CLASSID, 6);}
58955f2e967SMatthew G. Knepley   ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr);
59055f2e967SMatthew G. Knepley   for (b = 0; b < numBd; ++b) {
59155f2e967SMatthew G. Knepley     PetscBool       isEssential;
592*d7ddef95SMatthew G. Knepley     const char     *labelname;
593*d7ddef95SMatthew G. Knepley     DMLabel         label;
594*d7ddef95SMatthew G. Knepley     PetscInt        field;
595*d7ddef95SMatthew G. Knepley     PetscObject     obj;
596*d7ddef95SMatthew G. Knepley     PetscClassId    id;
59755f2e967SMatthew G. Knepley     void          (*func)();
598*d7ddef95SMatthew G. Knepley     PetscInt        numids;
599*d7ddef95SMatthew G. Knepley     const PetscInt *ids;
60055f2e967SMatthew G. Knepley     void           *ctx;
60155f2e967SMatthew G. Knepley 
60263d5297fSMatthew G. Knepley     ierr = DMPlexGetBoundary(dm, b, &isEssential, NULL, &labelname, &field, &func, &numids, &ids, &ctx);CHKERRQ(ierr);
603*d7ddef95SMatthew G. Knepley     if (!isEssential) continue;
60463d5297fSMatthew G. Knepley     ierr = DMPlexGetLabel(dm, labelname, &label);CHKERRQ(ierr);
605*d7ddef95SMatthew G. Knepley     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
606*d7ddef95SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
607*d7ddef95SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
608*d7ddef95SMatthew G. Knepley       ierr = DMPlexInsertBoundaryValues_FEM_Internal(dm, field, label, numids, ids, (void (*)(const PetscReal[], PetscScalar *, void *)) func, ctx, locX);CHKERRQ(ierr);
609*d7ddef95SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
610*d7ddef95SMatthew G. Knepley       ierr = DMPlexInsertBoundaryValues_FVM_Internal(dm, time, faceGeomFVM, cellGeomFVM, gradFVM,
611*d7ddef95SMatthew G. Knepley                                                      field, label, numids, ids, (PetscErrorCode (*)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*)) func, ctx, locX);CHKERRQ(ierr);
612*d7ddef95SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
61355f2e967SMatthew G. Knepley   }
61455f2e967SMatthew G. Knepley   PetscFunctionReturn(0);
61555f2e967SMatthew G. Knepley }
61655f2e967SMatthew G. Knepley 
617cb1e1211SMatthew G Knepley #undef __FUNCT__
618cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeL2Diff"
619cb1e1211SMatthew G Knepley /*@C
620cb1e1211SMatthew G Knepley   DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h.
621cb1e1211SMatthew G Knepley 
622cb1e1211SMatthew G Knepley   Input Parameters:
623cb1e1211SMatthew G Knepley + dm    - The DM
624cb1e1211SMatthew G Knepley . funcs - The functions to evaluate for each field component
62551259fa3SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
626cb1e1211SMatthew G Knepley - X     - The coefficient vector u_h
627cb1e1211SMatthew G Knepley 
628cb1e1211SMatthew G Knepley   Output Parameter:
629cb1e1211SMatthew G Knepley . diff - The diff ||u - u_h||_2
630cb1e1211SMatthew G Knepley 
631cb1e1211SMatthew G Knepley   Level: developer
632cb1e1211SMatthew G Knepley 
63315496722SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2FieldDiff(), DMPlexComputeL2GradientDiff()
634878cb397SSatish Balay @*/
6350f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeL2Diff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff)
636cb1e1211SMatthew G Knepley {
637cb1e1211SMatthew G Knepley   const PetscInt   debug = 0;
638cb1e1211SMatthew G Knepley   PetscSection     section;
639c5bbbd5bSMatthew G. Knepley   PetscQuadrature  quad;
640cb1e1211SMatthew G Knepley   Vec              localX;
64115496722SMatthew G. Knepley   PetscScalar     *funcVal, *interpolant;
642cb1e1211SMatthew G Knepley   PetscReal       *coords, *v0, *J, *invJ, detJ;
643cb1e1211SMatthew G Knepley   PetscReal        localDiff = 0.0;
64415496722SMatthew G. Knepley   const PetscReal *quadPoints, *quadWeights;
64515496722SMatthew G. Knepley   PetscInt         dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, c, field, fieldOffset;
646cb1e1211SMatthew G Knepley   PetscErrorCode   ierr;
647cb1e1211SMatthew G Knepley 
648cb1e1211SMatthew G Knepley   PetscFunctionBegin;
649c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
650cb1e1211SMatthew G Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
651cb1e1211SMatthew G Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
652cb1e1211SMatthew G Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
653cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
654cb1e1211SMatthew G Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
655cb1e1211SMatthew G Knepley   for (field = 0; field < numFields; ++field) {
65615496722SMatthew G. Knepley     PetscObject  obj;
65715496722SMatthew G. Knepley     PetscClassId id;
658c5bbbd5bSMatthew G. Knepley     PetscInt     Nc;
659c5bbbd5bSMatthew G. Knepley 
66015496722SMatthew G. Knepley     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
66115496722SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
66215496722SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
66315496722SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
66415496722SMatthew G. Knepley 
6650f2d7e86SMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
6660f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
66715496722SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
66815496722SMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
66915496722SMatthew G. Knepley 
67015496722SMatthew G. Knepley       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
67115496722SMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
67215496722SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
673c5bbbd5bSMatthew G. Knepley     numComponents += Nc;
674cb1e1211SMatthew G Knepley   }
67515496722SMatthew G. Knepley   ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
6760f2d7e86SMatthew G. Knepley   ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
67715496722SMatthew G. Knepley   ierr = PetscMalloc6(numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
678cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
679cb1e1211SMatthew G Knepley   for (c = cStart; c < cEnd; ++c) {
680a1e44745SMatthew G. Knepley     PetscScalar *x = NULL;
681cb1e1211SMatthew G Knepley     PetscReal    elemDiff = 0.0;
682cb1e1211SMatthew G Knepley 
6838e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
684cb1e1211SMatthew G Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
685cb1e1211SMatthew G Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
686cb1e1211SMatthew G Knepley 
68715496722SMatthew G. Knepley     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
68815496722SMatthew G. Knepley       PetscObject  obj;
68915496722SMatthew G. Knepley       PetscClassId id;
690c110b1eeSGeoffrey Irving       void * const ctx = ctxs ? ctxs[field] : NULL;
69115496722SMatthew G. Knepley       PetscInt     Nb, Nc, q, fc;
692cb1e1211SMatthew G Knepley 
69315496722SMatthew G. Knepley       ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
69415496722SMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
69515496722SMatthew G. Knepley       if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
69615496722SMatthew G. Knepley       else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
69715496722SMatthew G. Knepley       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
698cb1e1211SMatthew G Knepley       if (debug) {
699cb1e1211SMatthew G Knepley         char title[1024];
700cb1e1211SMatthew G Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
70115496722SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr);
702cb1e1211SMatthew G Knepley       }
703cb1e1211SMatthew G Knepley       for (q = 0; q < numQuadPoints; ++q) {
70415496722SMatthew G. Knepley         CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords);
705c110b1eeSGeoffrey Irving         (*funcs[field])(coords, funcVal, ctx);
70615496722SMatthew G. Knepley         if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
70715496722SMatthew G. Knepley         else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
70815496722SMatthew G. Knepley         else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
70915496722SMatthew G. Knepley         for (fc = 0; fc < Nc; ++fc) {
71015496722SMatthew G. Knepley           if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "    elem %d field %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ);CHKERRQ(ierr);}
71115496722SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ;
712cb1e1211SMatthew G Knepley         }
713cb1e1211SMatthew G Knepley       }
71415496722SMatthew G. Knepley       fieldOffset += Nb*Nc;
715cb1e1211SMatthew G Knepley     }
716cb1e1211SMatthew G Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
717cb1e1211SMatthew G Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
718cb1e1211SMatthew G Knepley     localDiff += elemDiff;
719cb1e1211SMatthew G Knepley   }
72015496722SMatthew G. Knepley   ierr  = PetscFree6(funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr);
721cb1e1211SMatthew G Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
72286a74ee0SMatthew G. Knepley   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
723cb1e1211SMatthew G Knepley   *diff = PetscSqrtReal(*diff);
724cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
725cb1e1211SMatthew G Knepley }
726cb1e1211SMatthew G Knepley 
727cb1e1211SMatthew G Knepley #undef __FUNCT__
72840e14135SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeL2GradientDiff"
72940e14135SMatthew G. Knepley /*@C
73040e14135SMatthew G. Knepley   DMPlexComputeL2GradientDiff - This function computes the L_2 difference between the gradient of a function u and an FEM interpolant solution grad u_h.
73140e14135SMatthew G. Knepley 
73240e14135SMatthew G. Knepley   Input Parameters:
73340e14135SMatthew G. Knepley + dm    - The DM
73440e14135SMatthew G. Knepley . funcs - The gradient functions to evaluate for each field component
73551259fa3SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
73640e14135SMatthew G. Knepley . X     - The coefficient vector u_h
73740e14135SMatthew G. Knepley - n     - The vector to project along
73840e14135SMatthew G. Knepley 
73940e14135SMatthew G. Knepley   Output Parameter:
74040e14135SMatthew G. Knepley . diff - The diff ||(grad u - grad u_h) . n||_2
74140e14135SMatthew G. Knepley 
74240e14135SMatthew G. Knepley   Level: developer
74340e14135SMatthew G. Knepley 
74440e14135SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff()
74540e14135SMatthew G. Knepley @*/
7460f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeL2GradientDiff(DM dm, void (**funcs)(const PetscReal [], const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff)
747cb1e1211SMatthew G Knepley {
74840e14135SMatthew G. Knepley   const PetscInt  debug = 0;
749cb1e1211SMatthew G Knepley   PetscSection    section;
75040e14135SMatthew G. Knepley   PetscQuadrature quad;
75140e14135SMatthew G. Knepley   Vec             localX;
75240e14135SMatthew G. Knepley   PetscScalar    *funcVal, *interpolantVec;
75340e14135SMatthew G. Knepley   PetscReal      *coords, *realSpaceDer, *v0, *J, *invJ, detJ;
75440e14135SMatthew G. Knepley   PetscReal       localDiff = 0.0;
75540e14135SMatthew G. Knepley   PetscInt        dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp;
756cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
757cb1e1211SMatthew G Knepley 
758cb1e1211SMatthew G Knepley   PetscFunctionBegin;
759c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
76040e14135SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
76140e14135SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
76240e14135SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
76340e14135SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
76440e14135SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
765652b88e8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) {
7660f2d7e86SMatthew G. Knepley     PetscFE  fe;
76740e14135SMatthew G. Knepley     PetscInt Nc;
768652b88e8SMatthew G. Knepley 
7690f2d7e86SMatthew G. Knepley     ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
7700f2d7e86SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
7710f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
77240e14135SMatthew G. Knepley     numComponents += Nc;
773652b88e8SMatthew G. Knepley   }
77440e14135SMatthew G. Knepley   /* ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */
77540e14135SMatthew G. Knepley   ierr = PetscMalloc7(numComponents,&funcVal,dim,&coords,dim,&realSpaceDer,dim,&v0,dim*dim,&J,dim*dim,&invJ,dim,&interpolantVec);CHKERRQ(ierr);
77640e14135SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
77740e14135SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
77840e14135SMatthew G. Knepley     PetscScalar *x = NULL;
77940e14135SMatthew G. Knepley     PetscReal    elemDiff = 0.0;
780652b88e8SMatthew G. Knepley 
7818e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
78240e14135SMatthew G. Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
78340e14135SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
78440e14135SMatthew G. Knepley 
78540e14135SMatthew G. Knepley     for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) {
7860f2d7e86SMatthew G. Knepley       PetscFE          fe;
78751259fa3SMatthew G. Knepley       void * const     ctx = ctxs ? ctxs[field] : NULL;
78821454ff5SMatthew G. Knepley       const PetscReal *quadPoints, *quadWeights;
78940e14135SMatthew G. Knepley       PetscReal       *basisDer;
79021454ff5SMatthew G. Knepley       PetscInt         numQuadPoints, Nb, Ncomp, q, d, e, fc, f, g;
79140e14135SMatthew G. Knepley 
7920f2d7e86SMatthew G. Knepley       ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr);
79321454ff5SMatthew G. Knepley       ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
7940f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
7950f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Ncomp);CHKERRQ(ierr);
7960f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr);
79740e14135SMatthew G. Knepley       if (debug) {
79840e14135SMatthew G. Knepley         char title[1024];
79940e14135SMatthew G. Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
80040e14135SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Ncomp, &x[fieldOffset]);CHKERRQ(ierr);
801652b88e8SMatthew G. Knepley       }
80240e14135SMatthew G. Knepley       for (q = 0; q < numQuadPoints; ++q) {
80340e14135SMatthew G. Knepley         for (d = 0; d < dim; d++) {
80440e14135SMatthew G. Knepley           coords[d] = v0[d];
80540e14135SMatthew G. Knepley           for (e = 0; e < dim; e++) {
80640e14135SMatthew G. Knepley             coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0);
807652b88e8SMatthew G. Knepley           }
80840e14135SMatthew G. Knepley         }
80951259fa3SMatthew G. Knepley         (*funcs[field])(coords, n, funcVal, ctx);
81040e14135SMatthew G. Knepley         for (fc = 0; fc < Ncomp; ++fc) {
81140e14135SMatthew G. Knepley           PetscScalar interpolant = 0.0;
81240e14135SMatthew G. Knepley 
81340e14135SMatthew G. Knepley           for (d = 0; d < dim; ++d) interpolantVec[d] = 0.0;
81440e14135SMatthew G. Knepley           for (f = 0; f < Nb; ++f) {
81540e14135SMatthew G. Knepley             const PetscInt fidx = f*Ncomp+fc;
81640e14135SMatthew G. Knepley 
81740e14135SMatthew G. Knepley             for (d = 0; d < dim; ++d) {
81840e14135SMatthew G. Knepley               realSpaceDer[d] = 0.0;
81940e14135SMatthew G. Knepley               for (g = 0; g < dim; ++g) {
82040e14135SMatthew G. Knepley                 realSpaceDer[d] += invJ[g*dim+d]*basisDer[(q*Nb*Ncomp+fidx)*dim+g];
82140e14135SMatthew G. Knepley               }
82240e14135SMatthew G. Knepley               interpolantVec[d] += x[fieldOffset+fidx]*realSpaceDer[d];
82340e14135SMatthew G. Knepley             }
82440e14135SMatthew G. Knepley           }
82540e14135SMatthew G. Knepley           for (d = 0; d < dim; ++d) interpolant += interpolantVec[d]*n[d];
82640e14135SMatthew G. Knepley           if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "    elem %d fieldDer %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ);CHKERRQ(ierr);}
82740e14135SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ;
82840e14135SMatthew G. Knepley         }
82940e14135SMatthew G. Knepley       }
83040e14135SMatthew G. Knepley       comp        += Ncomp;
83140e14135SMatthew G. Knepley       fieldOffset += Nb*Ncomp;
83240e14135SMatthew G. Knepley     }
83340e14135SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
83440e14135SMatthew G. Knepley     if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "  elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);}
83540e14135SMatthew G. Knepley     localDiff += elemDiff;
83640e14135SMatthew G. Knepley   }
83740e14135SMatthew G. Knepley   ierr  = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr);
83840e14135SMatthew G. Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
83940e14135SMatthew G. Knepley   ierr  = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
84040e14135SMatthew G. Knepley   *diff = PetscSqrtReal(*diff);
841cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
842cb1e1211SMatthew G Knepley }
843cb1e1211SMatthew G Knepley 
844a0845e3aSMatthew G. Knepley #undef __FUNCT__
84573d901b8SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeL2FieldDiff"
84615496722SMatthew G. Knepley /*@C
84715496722SMatthew G. Knepley   DMPlexComputeL2FieldDiff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h, separated into field components.
84815496722SMatthew G. Knepley 
84915496722SMatthew G. Knepley   Input Parameters:
85015496722SMatthew G. Knepley + dm    - The DM
85115496722SMatthew G. Knepley . funcs - The functions to evaluate for each field component
85215496722SMatthew G. Knepley . ctxs  - Optional array of contexts to pass to each function, or NULL.
85315496722SMatthew G. Knepley - X     - The coefficient vector u_h
85415496722SMatthew G. Knepley 
85515496722SMatthew G. Knepley   Output Parameter:
85615496722SMatthew G. Knepley . diff - The array of differences, ||u^f - u^f_h||_2
85715496722SMatthew G. Knepley 
85815496722SMatthew G. Knepley   Level: developer
85915496722SMatthew G. Knepley 
86015496722SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff(), DMPlexComputeL2GradientDiff()
86115496722SMatthew G. Knepley @*/
8620f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeL2FieldDiff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal diff[])
86373d901b8SMatthew G. Knepley {
86473d901b8SMatthew G. Knepley   const PetscInt   debug = 0;
86573d901b8SMatthew G. Knepley   PetscSection     section;
86673d901b8SMatthew G. Knepley   PetscQuadrature  quad;
86773d901b8SMatthew G. Knepley   Vec              localX;
86815496722SMatthew G. Knepley   PetscScalar     *funcVal, *interpolant;
86973d901b8SMatthew G. Knepley   PetscReal       *coords, *v0, *J, *invJ, detJ;
87073d901b8SMatthew G. Knepley   PetscReal       *localDiff;
87115496722SMatthew G. Knepley   const PetscReal *quadPoints, *quadWeights;
87215496722SMatthew G. Knepley   PetscInt         dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, c, field, fieldOffset;
87373d901b8SMatthew G. Knepley   PetscErrorCode   ierr;
87473d901b8SMatthew G. Knepley 
87573d901b8SMatthew G. Knepley   PetscFunctionBegin;
876c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
87773d901b8SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
87873d901b8SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr);
87973d901b8SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
88073d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
88173d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
88273d901b8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) {
88315496722SMatthew G. Knepley     PetscObject  obj;
88415496722SMatthew G. Knepley     PetscClassId id;
88573d901b8SMatthew G. Knepley     PetscInt     Nc;
88673d901b8SMatthew G. Knepley 
88715496722SMatthew G. Knepley     ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
88815496722SMatthew G. Knepley     ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
88915496722SMatthew G. Knepley     if (id == PETSCFE_CLASSID) {
89015496722SMatthew G. Knepley       PetscFE fe = (PetscFE) obj;
89115496722SMatthew G. Knepley 
8920f2d7e86SMatthew G. Knepley       ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr);
8930f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
89415496722SMatthew G. Knepley     } else if (id == PETSCFV_CLASSID) {
89515496722SMatthew G. Knepley       PetscFV fv = (PetscFV) obj;
89615496722SMatthew G. Knepley 
89715496722SMatthew G. Knepley       ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr);
89815496722SMatthew G. Knepley       ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr);
89915496722SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
90073d901b8SMatthew G. Knepley     numComponents += Nc;
90173d901b8SMatthew G. Knepley   }
90215496722SMatthew G. Knepley   ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr);
9030f2d7e86SMatthew G. Knepley   ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr);
90415496722SMatthew G. Knepley   ierr = PetscCalloc7(numFields,&localDiff,numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr);
90573d901b8SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
90673d901b8SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
90773d901b8SMatthew G. Knepley     PetscScalar *x = NULL;
90873d901b8SMatthew G. Knepley 
9098e0841e0SMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr);
91073d901b8SMatthew G. Knepley     if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c);
91173d901b8SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
91273d901b8SMatthew G. Knepley 
91315496722SMatthew G. Knepley     for (field = 0, fieldOffset = 0; field < numFields; ++field) {
91415496722SMatthew G. Knepley       PetscObject  obj;
91515496722SMatthew G. Knepley       PetscClassId id;
91673d901b8SMatthew G. Knepley       void * const ctx = ctxs ? ctxs[field] : NULL;
91715496722SMatthew G. Knepley       PetscInt     Nb, Nc, q, fc;
91873d901b8SMatthew G. Knepley 
91915496722SMatthew G. Knepley       PetscReal       elemDiff = 0.0;
92015496722SMatthew G. Knepley 
92115496722SMatthew G. Knepley       ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr);
92215496722SMatthew G. Knepley       ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr);
92315496722SMatthew G. Knepley       if (id == PETSCFE_CLASSID)      {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);}
92415496722SMatthew G. Knepley       else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;}
92515496722SMatthew G. Knepley       else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
92673d901b8SMatthew G. Knepley       if (debug) {
92773d901b8SMatthew G. Knepley         char title[1024];
92873d901b8SMatthew G. Knepley         ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr);
92915496722SMatthew G. Knepley         ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr);
93073d901b8SMatthew G. Knepley       }
93173d901b8SMatthew G. Knepley       for (q = 0; q < numQuadPoints; ++q) {
93215496722SMatthew G. Knepley         CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords);
93373d901b8SMatthew G. Knepley         (*funcs[field])(coords, funcVal, ctx);
93415496722SMatthew G. Knepley         if (id == PETSCFE_CLASSID)      {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
93515496722SMatthew G. Knepley         else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);}
93615496722SMatthew G. Knepley         else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field);
93715496722SMatthew G. Knepley         for (fc = 0; fc < Nc; ++fc) {
93815496722SMatthew G. Knepley           if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, "    elem %d field %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ);CHKERRQ(ierr);}
93915496722SMatthew G. Knepley           elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ;
94073d901b8SMatthew G. Knepley         }
94173d901b8SMatthew G. Knepley       }
94215496722SMatthew G. Knepley       fieldOffset += Nb*Nc;
94373d901b8SMatthew G. Knepley       localDiff[field] += elemDiff;
94473d901b8SMatthew G. Knepley     }
94573d901b8SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr);
94673d901b8SMatthew G. Knepley   }
94773d901b8SMatthew G. Knepley   ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
94873d901b8SMatthew G. Knepley   ierr = MPI_Allreduce(localDiff, diff, numFields, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
94973d901b8SMatthew G. Knepley   for (field = 0; field < numFields; ++field) diff[field] = PetscSqrtReal(diff[field]);
95015496722SMatthew G. Knepley   ierr = PetscFree7(localDiff,funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr);
95173d901b8SMatthew G. Knepley   PetscFunctionReturn(0);
95273d901b8SMatthew G. Knepley }
95373d901b8SMatthew G. Knepley 
95473d901b8SMatthew G. Knepley #undef __FUNCT__
95573d901b8SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeIntegralFEM"
95673d901b8SMatthew G. Knepley /*@
95773d901b8SMatthew G. Knepley   DMPlexComputeIntegralFEM - Form the local integral F from the local input X using pointwise functions specified by the user
95873d901b8SMatthew G. Knepley 
95973d901b8SMatthew G. Knepley   Input Parameters:
96073d901b8SMatthew G. Knepley + dm - The mesh
96173d901b8SMatthew G. Knepley . X  - Local input vector
96273d901b8SMatthew G. Knepley - user - The user context
96373d901b8SMatthew G. Knepley 
96473d901b8SMatthew G. Knepley   Output Parameter:
96573d901b8SMatthew G. Knepley . integral - Local integral for each field
96673d901b8SMatthew G. Knepley 
96773d901b8SMatthew G. Knepley   Level: developer
96873d901b8SMatthew G. Knepley 
96973d901b8SMatthew G. Knepley .seealso: DMPlexComputeResidualFEM()
97073d901b8SMatthew G. Knepley @*/
9710f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeIntegralFEM(DM dm, Vec X, PetscReal *integral, void *user)
97273d901b8SMatthew G. Knepley {
97373d901b8SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dm->data;
97473d901b8SMatthew G. Knepley   DM                dmAux;
97573d901b8SMatthew G. Knepley   Vec               localX, A;
9762764a2aaSMatthew G. Knepley   PetscDS           prob, probAux = NULL;
97773d901b8SMatthew G. Knepley   PetscQuadrature   q;
97873d901b8SMatthew G. Knepley   PetscSection      section, sectionAux;
979bbce034cSMatthew G. Knepley   PetscFECellGeom  *cgeom;
98073d901b8SMatthew G. Knepley   PetscScalar      *u, *a = NULL;
9810f2d7e86SMatthew G. Knepley   PetscInt          dim, Nf, f, numCells, cStart, cEnd, c;
9820f2d7e86SMatthew G. Knepley   PetscInt          totDim, totDimAux;
98373d901b8SMatthew G. Knepley   PetscErrorCode    ierr;
98473d901b8SMatthew G. Knepley 
98573d901b8SMatthew G. Knepley   PetscFunctionBegin;
98673d901b8SMatthew G. Knepley   /*ierr = PetscLogEventBegin(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);*/
98773d901b8SMatthew G. Knepley   ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr);
98873d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
98973d901b8SMatthew G. Knepley   ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr);
990c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
99173d901b8SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
9922764a2aaSMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
9932764a2aaSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr);
99473d901b8SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
99573d901b8SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
99673d901b8SMatthew G. Knepley   numCells = cEnd - cStart;
9970f2d7e86SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {integral[f]    = 0.0;}
99873d901b8SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr);
99973d901b8SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr);
100073d901b8SMatthew G. Knepley   if (dmAux) {
100173d901b8SMatthew G. Knepley     ierr = DMGetDefaultSection(dmAux, &sectionAux);CHKERRQ(ierr);
10022764a2aaSMatthew G. Knepley     ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr);
10032764a2aaSMatthew G. Knepley     ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr);
100473d901b8SMatthew G. Knepley   }
1005*d7ddef95SMatthew G. Knepley   ierr = DMPlexInsertBoundaryValues(dm, localX, 0.0, NULL, NULL, NULL);CHKERRQ(ierr);
1006bbce034cSMatthew G. Knepley   ierr = PetscMalloc2(numCells*totDim,&u,numCells,&cgeom);CHKERRQ(ierr);
10070f2d7e86SMatthew G. Knepley   if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);}
100873d901b8SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
100973d901b8SMatthew G. Knepley     PetscScalar *x = NULL;
101073d901b8SMatthew G. Knepley     PetscInt     i;
101173d901b8SMatthew G. Knepley 
1012bbce034cSMatthew G. Knepley     ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, cgeom[c].v0, cgeom[c].J, cgeom[c].invJ, &cgeom[c].detJ);CHKERRQ(ierr);
1013bbce034cSMatthew G. Knepley     if (cgeom[c].detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", cgeom[c].detJ, c);
101473d901b8SMatthew G. Knepley     ierr = DMPlexVecGetClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
10150f2d7e86SMatthew G. Knepley     for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i];
101673d901b8SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr);
101773d901b8SMatthew G. Knepley     if (dmAux) {
101873d901b8SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
10190f2d7e86SMatthew G. Knepley       for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i];
102073d901b8SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr);
102173d901b8SMatthew G. Knepley     }
102273d901b8SMatthew G. Knepley   }
102373d901b8SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
10240f2d7e86SMatthew G. Knepley     PetscFE  fe;
102573d901b8SMatthew G. Knepley     PetscInt numQuadPoints, Nb;
102673d901b8SMatthew G. Knepley     /* Conforming batches */
102773d901b8SMatthew G. Knepley     PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize;
102873d901b8SMatthew G. Knepley     /* Remainder */
102973d901b8SMatthew G. Knepley     PetscInt Nr, offset;
103073d901b8SMatthew G. Knepley 
10312764a2aaSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
10320f2d7e86SMatthew G. Knepley     ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr);
10330f2d7e86SMatthew G. Knepley     ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr);
10340f2d7e86SMatthew G. Knepley     ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr);
103573d901b8SMatthew G. Knepley     ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr);
103673d901b8SMatthew G. Knepley     blockSize = Nb*numQuadPoints;
103773d901b8SMatthew G. Knepley     batchSize = numBlocks * blockSize;
10380f2d7e86SMatthew G. Knepley     ierr =  PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr);
103973d901b8SMatthew G. Knepley     numChunks = numCells / (numBatches*batchSize);
104073d901b8SMatthew G. Knepley     Ne        = numChunks*numBatches*batchSize;
104173d901b8SMatthew G. Knepley     Nr        = numCells % (numBatches*batchSize);
104273d901b8SMatthew G. Knepley     offset    = numCells - Nr;
1043bbce034cSMatthew G. Knepley     ierr = PetscFEIntegrate(fe, prob, f, Ne, cgeom, u, probAux, a, integral);CHKERRQ(ierr);
1044bbce034cSMatthew G. Knepley     ierr = PetscFEIntegrate(fe, prob, f, Nr, &cgeom[offset], &u[offset*totDim], probAux, &a[offset*totDimAux], integral);CHKERRQ(ierr);
104573d901b8SMatthew G. Knepley   }
1046bbce034cSMatthew G. Knepley   ierr = PetscFree2(u,cgeom);CHKERRQ(ierr);
104773d901b8SMatthew G. Knepley   if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);}
104873d901b8SMatthew G. Knepley   if (mesh->printFEM) {
104973d901b8SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Local integral:");CHKERRQ(ierr);
105073d901b8SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), " %g", integral[f]);CHKERRQ(ierr);}
105173d901b8SMatthew G. Knepley     ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "\n");CHKERRQ(ierr);
105273d901b8SMatthew G. Knepley   }
105373d901b8SMatthew G. Knepley   ierr  = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr);
105473d901b8SMatthew G. Knepley   /* TODO: Allreduce for integral */
105573d901b8SMatthew G. Knepley   /*ierr = PetscLogEventEnd(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);*/
105673d901b8SMatthew G. Knepley   PetscFunctionReturn(0);
105773d901b8SMatthew G. Knepley }
105873d901b8SMatthew G. Knepley 
105973d901b8SMatthew G. Knepley #undef __FUNCT__
1060d69c5d34SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeInterpolatorFEM"
1061d69c5d34SMatthew G. Knepley /*@
1062d69c5d34SMatthew G. Knepley   DMPlexComputeInterpolatorFEM - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM.
1063d69c5d34SMatthew G. Knepley 
1064d69c5d34SMatthew G. Knepley   Input Parameters:
1065d69c5d34SMatthew G. Knepley + dmf  - The fine mesh
1066d69c5d34SMatthew G. Knepley . dmc  - The coarse mesh
1067d69c5d34SMatthew G. Knepley - user - The user context
1068d69c5d34SMatthew G. Knepley 
1069d69c5d34SMatthew G. Knepley   Output Parameter:
1070934789fcSMatthew G. Knepley . In  - The interpolation matrix
1071d69c5d34SMatthew G. Knepley 
1072d69c5d34SMatthew G. Knepley   Note:
1073d69c5d34SMatthew G. Knepley   The first member of the user context must be an FEMContext.
1074d69c5d34SMatthew G. Knepley 
1075d69c5d34SMatthew G. Knepley   We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator,
1076d69c5d34SMatthew G. Knepley   like a GPU, or vectorize on a multicore machine.
1077d69c5d34SMatthew G. Knepley 
1078d69c5d34SMatthew G. Knepley   Level: developer
1079d69c5d34SMatthew G. Knepley 
1080d69c5d34SMatthew G. Knepley .seealso: DMPlexComputeJacobianFEM()
1081d69c5d34SMatthew G. Knepley @*/
1082934789fcSMatthew G. Knepley PetscErrorCode DMPlexComputeInterpolatorFEM(DM dmc, DM dmf, Mat In, void *user)
1083d69c5d34SMatthew G. Knepley {
1084d69c5d34SMatthew G. Knepley   DM_Plex          *mesh  = (DM_Plex *) dmc->data;
1085d69c5d34SMatthew G. Knepley   const char       *name  = "Interpolator";
10862764a2aaSMatthew G. Knepley   PetscDS           prob;
1087d69c5d34SMatthew G. Knepley   PetscFE          *feRef;
1088d69c5d34SMatthew G. Knepley   PetscSection      fsection, fglobalSection;
1089d69c5d34SMatthew G. Knepley   PetscSection      csection, cglobalSection;
1090d69c5d34SMatthew G. Knepley   PetscScalar      *elemMat;
1091942a7a06SMatthew G. Knepley   PetscInt          dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, c;
10920f2d7e86SMatthew G. Knepley   PetscInt          cTotDim, rTotDim = 0;
1093d69c5d34SMatthew G. Knepley   PetscErrorCode    ierr;
1094d69c5d34SMatthew G. Knepley 
1095d69c5d34SMatthew G. Knepley   PetscFunctionBegin;
1096d69c5d34SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1097c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr);
1098d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
1099d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
1100d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
1101d69c5d34SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
1102d69c5d34SMatthew G. Knepley   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
1103d69c5d34SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
11042764a2aaSMatthew G. Knepley   ierr = DMGetDS(dmf, &prob);CHKERRQ(ierr);
1105d69c5d34SMatthew G. Knepley   ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr);
1106d69c5d34SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
11070f2d7e86SMatthew G. Knepley     PetscFE  fe;
11080f2d7e86SMatthew G. Knepley     PetscInt rNb, Nc;
1109d69c5d34SMatthew G. Knepley 
11102764a2aaSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
11110f2d7e86SMatthew G. Knepley     ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr);
1112d69c5d34SMatthew G. Knepley     ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr);
11130f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
11140f2d7e86SMatthew G. Knepley     rTotDim += rNb*Nc;
1115d69c5d34SMatthew G. Knepley   }
11162764a2aaSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr);
11170f2d7e86SMatthew G. Knepley   ierr = PetscMalloc1(rTotDim*cTotDim,&elemMat);CHKERRQ(ierr);
11180f2d7e86SMatthew G. Knepley   ierr = PetscMemzero(elemMat, rTotDim*cTotDim * sizeof(PetscScalar));CHKERRQ(ierr);
1119d69c5d34SMatthew G. Knepley   for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) {
1120d69c5d34SMatthew G. Knepley     PetscDualSpace   Qref;
1121d69c5d34SMatthew G. Knepley     PetscQuadrature  f;
1122d69c5d34SMatthew G. Knepley     const PetscReal *qpoints, *qweights;
1123d69c5d34SMatthew G. Knepley     PetscReal       *points;
1124d69c5d34SMatthew G. Knepley     PetscInt         npoints = 0, Nc, Np, fpdim, i, k, p, d;
1125d69c5d34SMatthew G. Knepley 
1126d69c5d34SMatthew G. Knepley     /* Compose points from all dual basis functionals */
1127d69c5d34SMatthew G. Knepley     ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr);
11280f2d7e86SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feRef[fieldI], &Nc);CHKERRQ(ierr);
1129d69c5d34SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr);
1130d69c5d34SMatthew G. Knepley     for (i = 0; i < fpdim; ++i) {
1131d69c5d34SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1132d69c5d34SMatthew G. Knepley       ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, NULL);CHKERRQ(ierr);
1133d69c5d34SMatthew G. Knepley       npoints += Np;
1134d69c5d34SMatthew G. Knepley     }
1135d69c5d34SMatthew G. Knepley     ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr);
1136d69c5d34SMatthew G. Knepley     for (i = 0, k = 0; i < fpdim; ++i) {
1137d69c5d34SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1138d69c5d34SMatthew G. Knepley       ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr);
1139d69c5d34SMatthew G. Knepley       for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d];
1140d69c5d34SMatthew G. Knepley     }
1141d69c5d34SMatthew G. Knepley 
1142d69c5d34SMatthew G. Knepley     for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) {
11430f2d7e86SMatthew G. Knepley       PetscFE    fe;
1144d69c5d34SMatthew G. Knepley       PetscReal *B;
114536a6d9c0SMatthew G. Knepley       PetscInt   NcJ, cpdim, j;
1146d69c5d34SMatthew G. Knepley 
1147d69c5d34SMatthew G. Knepley       /* Evaluate basis at points */
11482764a2aaSMatthew G. Knepley       ierr = PetscDSGetDiscretization(prob, fieldJ, (PetscObject *) &fe);CHKERRQ(ierr);
11490f2d7e86SMatthew G. Knepley       ierr = PetscFEGetNumComponents(fe, &NcJ);CHKERRQ(ierr);
11500f2d7e86SMatthew G. Knepley       ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr);
1151ffe73a53SMatthew G. Knepley       /* For now, fields only interpolate themselves */
1152ffe73a53SMatthew G. Knepley       if (fieldI == fieldJ) {
11537c927364SMatthew G. Knepley         if (Nc != NcJ) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in fine space field %d does not match coarse field %d", Nc, NcJ);
11540f2d7e86SMatthew G. Knepley         ierr = PetscFEGetTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);
1155d69c5d34SMatthew G. Knepley         for (i = 0, k = 0; i < fpdim; ++i) {
1156d69c5d34SMatthew G. Knepley           ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr);
1157d69c5d34SMatthew G. Knepley           ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr);
1158d69c5d34SMatthew G. Knepley           for (p = 0; p < Np; ++p, ++k) {
115936a6d9c0SMatthew G. Knepley             for (j = 0; j < cpdim; ++j) {
11600f2d7e86SMatthew G. Knepley               for (c = 0; c < Nc; ++c) elemMat[(offsetI + i*Nc + c)*cTotDim + offsetJ + j*NcJ + c] += B[k*cpdim*NcJ+j*Nc+c]*qweights[p];
116136a6d9c0SMatthew G. Knepley             }
1162d69c5d34SMatthew G. Knepley           }
1163d69c5d34SMatthew G. Knepley         }
11640f2d7e86SMatthew G. Knepley         ierr = PetscFERestoreTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr);
1165ffe73a53SMatthew G. Knepley       }
116636a6d9c0SMatthew G. Knepley       offsetJ += cpdim*NcJ;
1167d69c5d34SMatthew G. Knepley     }
1168d69c5d34SMatthew G. Knepley     offsetI += fpdim*Nc;
1169549a8adaSMatthew G. Knepley     ierr = PetscFree(points);CHKERRQ(ierr);
1170d69c5d34SMatthew G. Knepley   }
11710f2d7e86SMatthew G. Knepley   if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rTotDim, cTotDim, elemMat);CHKERRQ(ierr);}
11727f5b169aSMatthew G. Knepley   /* Preallocate matrix */
11737f5b169aSMatthew G. Knepley   {
11747f5b169aSMatthew G. Knepley     PetscHashJK ht;
11757f5b169aSMatthew G. Knepley     PetscLayout rLayout;
11767f5b169aSMatthew G. Knepley     PetscInt   *dnz, *onz, *cellCIndices, *cellFIndices;
11777f5b169aSMatthew G. Knepley     PetscInt    locRows, rStart, rEnd, cell, r;
11787f5b169aSMatthew G. Knepley 
11797f5b169aSMatthew G. Knepley     ierr = MatGetLocalSize(In, &locRows, NULL);CHKERRQ(ierr);
11807f5b169aSMatthew G. Knepley     ierr = PetscLayoutCreate(PetscObjectComm((PetscObject) In), &rLayout);CHKERRQ(ierr);
11817f5b169aSMatthew G. Knepley     ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr);
11827f5b169aSMatthew G. Knepley     ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr);
11837f5b169aSMatthew G. Knepley     ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr);
11847f5b169aSMatthew G. Knepley     ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr);
11857f5b169aSMatthew G. Knepley     ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr);
11867f5b169aSMatthew G. Knepley     ierr = PetscCalloc4(locRows,&dnz,locRows,&onz,cTotDim,&cellCIndices,rTotDim,&cellFIndices);CHKERRQ(ierr);
11877f5b169aSMatthew G. Knepley     ierr = PetscHashJKCreate(&ht);CHKERRQ(ierr);
11887f5b169aSMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
11897f5b169aSMatthew G. Knepley       ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, cell, cellCIndices, cellFIndices);CHKERRQ(ierr);
11907f5b169aSMatthew G. Knepley       for (r = 0; r < rTotDim; ++r) {
11917f5b169aSMatthew G. Knepley         PetscHashJKKey  key;
11927f5b169aSMatthew G. Knepley         PetscHashJKIter missing, iter;
11937f5b169aSMatthew G. Knepley 
11947f5b169aSMatthew G. Knepley         key.j = cellFIndices[r];
11957f5b169aSMatthew G. Knepley         if (key.j < 0) continue;
11967f5b169aSMatthew G. Knepley         for (c = 0; c < cTotDim; ++c) {
11977f5b169aSMatthew G. Knepley           key.k = cellCIndices[c];
11987f5b169aSMatthew G. Knepley           if (key.k < 0) continue;
11997f5b169aSMatthew G. Knepley           ierr = PetscHashJKPut(ht, key, &missing, &iter);CHKERRQ(ierr);
12007f5b169aSMatthew G. Knepley           if (missing) {
12017f5b169aSMatthew G. Knepley             ierr = PetscHashJKSet(ht, iter, 1);CHKERRQ(ierr);
12027f5b169aSMatthew G. Knepley             if ((key.k >= rStart) && (key.k < rEnd)) ++dnz[key.j-rStart];
12037f5b169aSMatthew G. Knepley             else                                     ++onz[key.j-rStart];
12047f5b169aSMatthew G. Knepley           }
12057f5b169aSMatthew G. Knepley         }
12067f5b169aSMatthew G. Knepley       }
12077f5b169aSMatthew G. Knepley     }
12087f5b169aSMatthew G. Knepley     ierr = PetscHashJKDestroy(&ht);CHKERRQ(ierr);
12097f5b169aSMatthew G. Knepley     ierr = MatXAIJSetPreallocation(In, 1, dnz, onz, NULL, NULL);CHKERRQ(ierr);
12107f5b169aSMatthew G. Knepley     ierr = MatSetOption(In, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
12117f5b169aSMatthew G. Knepley     ierr = PetscFree4(dnz,onz,cellCIndices,cellFIndices);CHKERRQ(ierr);
12127f5b169aSMatthew G. Knepley   }
12137f5b169aSMatthew G. Knepley   /* Fill matrix */
12147f5b169aSMatthew G. Knepley   ierr = MatZeroEntries(In);CHKERRQ(ierr);
1215d69c5d34SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
1216934789fcSMatthew G. Knepley     ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr);
1217d69c5d34SMatthew G. Knepley   }
1218549a8adaSMatthew G. Knepley   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);}
1219d69c5d34SMatthew G. Knepley   ierr = PetscFree(feRef);CHKERRQ(ierr);
1220549a8adaSMatthew G. Knepley   ierr = PetscFree(elemMat);CHKERRQ(ierr);
1221934789fcSMatthew G. Knepley   ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1222934789fcSMatthew G. Knepley   ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1223d69c5d34SMatthew G. Knepley   if (mesh->printFEM) {
1224d69c5d34SMatthew G. Knepley     ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr);
1225934789fcSMatthew G. Knepley     ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr);
1226934789fcSMatthew G. Knepley     ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1227d69c5d34SMatthew G. Knepley   }
1228d69c5d34SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1229d69c5d34SMatthew G. Knepley   PetscFunctionReturn(0);
1230d69c5d34SMatthew G. Knepley }
12316c73c22cSMatthew G. Knepley 
12326c73c22cSMatthew G. Knepley #undef __FUNCT__
12337c927364SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeInjectorFEM"
12347c927364SMatthew G. Knepley PetscErrorCode DMPlexComputeInjectorFEM(DM dmc, DM dmf, VecScatter *sc, void *user)
12357c927364SMatthew G. Knepley {
1236e9d4ef1bSMatthew G. Knepley   PetscDS        prob;
12377c927364SMatthew G. Knepley   PetscFE       *feRef;
12387c927364SMatthew G. Knepley   Vec            fv, cv;
12397c927364SMatthew G. Knepley   IS             fis, cis;
12407c927364SMatthew G. Knepley   PetscSection   fsection, fglobalSection, csection, cglobalSection;
12417c927364SMatthew G. Knepley   PetscInt      *cmap, *cellCIndices, *cellFIndices, *cindices, *findices;
12427c927364SMatthew G. Knepley   PetscInt       cTotDim, fTotDim = 0, Nf, f, field, cStart, cEnd, c, dim, d, startC, offsetC, offsetF, m;
12437c927364SMatthew G. Knepley   PetscErrorCode ierr;
12447c927364SMatthew G. Knepley 
12457c927364SMatthew G. Knepley   PetscFunctionBegin;
124675a69067SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1247c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr);
12487c927364SMatthew G. Knepley   ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr);
12497c927364SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr);
12507c927364SMatthew G. Knepley   ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr);
12517c927364SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr);
12527c927364SMatthew G. Knepley   ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr);
12537c927364SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr);
1254e9d4ef1bSMatthew G. Knepley   ierr = DMGetDS(dmc, &prob);CHKERRQ(ierr);
12557c927364SMatthew G. Knepley   ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr);
12567c927364SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
12577c927364SMatthew G. Knepley     PetscFE  fe;
12587c927364SMatthew G. Knepley     PetscInt fNb, Nc;
12597c927364SMatthew G. Knepley 
1260e9d4ef1bSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr);
12617c927364SMatthew G. Knepley     ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr);
12627c927364SMatthew G. Knepley     ierr = PetscFEGetDimension(feRef[f], &fNb);CHKERRQ(ierr);
12637c927364SMatthew G. Knepley     ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr);
12647c927364SMatthew G. Knepley     fTotDim += fNb*Nc;
12657c927364SMatthew G. Knepley   }
1266e9d4ef1bSMatthew G. Knepley   ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr);
12677c927364SMatthew G. Knepley   ierr = PetscMalloc1(cTotDim,&cmap);CHKERRQ(ierr);
12687c927364SMatthew G. Knepley   for (field = 0, offsetC = 0, offsetF = 0; field < Nf; ++field) {
12697c927364SMatthew G. Knepley     PetscFE        feC;
12707c927364SMatthew G. Knepley     PetscDualSpace QF, QC;
12717c927364SMatthew G. Knepley     PetscInt       NcF, NcC, fpdim, cpdim;
12727c927364SMatthew G. Knepley 
1273e9d4ef1bSMatthew G. Knepley     ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &feC);CHKERRQ(ierr);
12747c927364SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feC, &NcC);CHKERRQ(ierr);
12757c927364SMatthew G. Knepley     ierr = PetscFEGetNumComponents(feRef[field], &NcF);CHKERRQ(ierr);
12767c927364SMatthew G. Knepley     if (NcF != NcC) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in fine space field %d does not match coarse field %d", NcF, NcC);
12777c927364SMatthew G. Knepley     ierr = PetscFEGetDualSpace(feRef[field], &QF);CHKERRQ(ierr);
12787c927364SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(QF, &fpdim);CHKERRQ(ierr);
12797c927364SMatthew G. Knepley     ierr = PetscFEGetDualSpace(feC, &QC);CHKERRQ(ierr);
12807c927364SMatthew G. Knepley     ierr = PetscDualSpaceGetDimension(QC, &cpdim);CHKERRQ(ierr);
12817c927364SMatthew G. Knepley     for (c = 0; c < cpdim; ++c) {
12827c927364SMatthew G. Knepley       PetscQuadrature  cfunc;
12837c927364SMatthew G. Knepley       const PetscReal *cqpoints;
12847c927364SMatthew G. Knepley       PetscInt         NpC;
12857c927364SMatthew G. Knepley 
12867c927364SMatthew G. Knepley       ierr = PetscDualSpaceGetFunctional(QC, c, &cfunc);CHKERRQ(ierr);
12877c927364SMatthew G. Knepley       ierr = PetscQuadratureGetData(cfunc, NULL, &NpC, &cqpoints, NULL);CHKERRQ(ierr);
12887c927364SMatthew G. Knepley       if (NpC != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not know how to do injection for moments");
12897c927364SMatthew G. Knepley       for (f = 0; f < fpdim; ++f) {
12907c927364SMatthew G. Knepley         PetscQuadrature  ffunc;
12917c927364SMatthew G. Knepley         const PetscReal *fqpoints;
12927c927364SMatthew G. Knepley         PetscReal        sum = 0.0;
12937c927364SMatthew G. Knepley         PetscInt         NpF, comp;
12947c927364SMatthew G. Knepley 
12957c927364SMatthew G. Knepley         ierr = PetscDualSpaceGetFunctional(QF, f, &ffunc);CHKERRQ(ierr);
12967c927364SMatthew G. Knepley         ierr = PetscQuadratureGetData(ffunc, NULL, &NpF, &fqpoints, NULL);CHKERRQ(ierr);
12977c927364SMatthew G. Knepley         if (NpC != NpF) continue;
12987c927364SMatthew G. Knepley         for (d = 0; d < dim; ++d) sum += PetscAbsReal(cqpoints[d] - fqpoints[d]);
12997c927364SMatthew G. Knepley         if (sum > 1.0e-9) continue;
13007c927364SMatthew G. Knepley         for (comp = 0; comp < NcC; ++comp) {
13017c927364SMatthew G. Knepley           cmap[(offsetC+c)*NcC+comp] = (offsetF+f)*NcF+comp;
13027c927364SMatthew G. Knepley         }
13037c927364SMatthew G. Knepley         break;
13047c927364SMatthew G. Knepley       }
13057c927364SMatthew G. Knepley     }
13067c927364SMatthew G. Knepley     offsetC += cpdim*NcC;
13077c927364SMatthew G. Knepley     offsetF += fpdim*NcF;
13087c927364SMatthew G. Knepley   }
13097c927364SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);}
13107c927364SMatthew G. Knepley   ierr = PetscFree(feRef);CHKERRQ(ierr);
13117c927364SMatthew G. Knepley 
13127c927364SMatthew G. Knepley   ierr = DMGetGlobalVector(dmf, &fv);CHKERRQ(ierr);
13137c927364SMatthew G. Knepley   ierr = DMGetGlobalVector(dmc, &cv);CHKERRQ(ierr);
13147c927364SMatthew G. Knepley   ierr = VecGetOwnershipRange(cv, &startC, NULL);CHKERRQ(ierr);
13157c927364SMatthew G. Knepley   ierr = PetscSectionGetConstrainedStorageSize(cglobalSection, &m);CHKERRQ(ierr);
13167c927364SMatthew G. Knepley   ierr = PetscMalloc2(cTotDim,&cellCIndices,fTotDim,&cellFIndices);CHKERRQ(ierr);
1317aa7da3c4SMatthew G. Knepley   ierr = PetscMalloc1(m,&cindices);CHKERRQ(ierr);
1318aa7da3c4SMatthew G. Knepley   ierr = PetscMalloc1(m,&findices);CHKERRQ(ierr);
13197c927364SMatthew G. Knepley   for (d = 0; d < m; ++d) cindices[d] = findices[d] = -1;
13207c927364SMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
13217c927364SMatthew G. Knepley     ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, c, cellCIndices, cellFIndices);CHKERRQ(ierr);
13227c927364SMatthew G. Knepley     for (d = 0; d < cTotDim; ++d) {
13237c927364SMatthew G. Knepley       if (cellCIndices[d] < 0) continue;
13247c927364SMatthew G. Knepley       if ((findices[cellCIndices[d]-startC] >= 0) && (findices[cellCIndices[d]-startC] != cellFIndices[cmap[d]])) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Coarse dof %d maps to both %d and %d", cindices[cellCIndices[d]-startC], findices[cellCIndices[d]-startC], cellFIndices[cmap[d]]);
13257c927364SMatthew G. Knepley       cindices[cellCIndices[d]-startC] = cellCIndices[d];
13267c927364SMatthew G. Knepley       findices[cellCIndices[d]-startC] = cellFIndices[cmap[d]];
13277c927364SMatthew G. Knepley     }
13287c927364SMatthew G. Knepley   }
13297c927364SMatthew G. Knepley   ierr = PetscFree(cmap);CHKERRQ(ierr);
13307c927364SMatthew G. Knepley   ierr = PetscFree2(cellCIndices,cellFIndices);CHKERRQ(ierr);
13317c927364SMatthew G. Knepley 
13327c927364SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, m, cindices, PETSC_OWN_POINTER, &cis);CHKERRQ(ierr);
13337c927364SMatthew G. Knepley   ierr = ISCreateGeneral(PETSC_COMM_SELF, m, findices, PETSC_OWN_POINTER, &fis);CHKERRQ(ierr);
13347c927364SMatthew G. Knepley   ierr = VecScatterCreate(cv, cis, fv, fis, sc);CHKERRQ(ierr);
13357c927364SMatthew G. Knepley   ierr = ISDestroy(&cis);CHKERRQ(ierr);
13367c927364SMatthew G. Knepley   ierr = ISDestroy(&fis);CHKERRQ(ierr);
13377c927364SMatthew G. Knepley   ierr = DMRestoreGlobalVector(dmf, &fv);CHKERRQ(ierr);
13387c927364SMatthew G. Knepley   ierr = DMRestoreGlobalVector(dmc, &cv);CHKERRQ(ierr);
133975a69067SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr);
1340cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
1341cb1e1211SMatthew G Knepley }
1342