Actual source code: quadratic.c

  1: #ifdef PETSC_RCS_HEADER
  2: static char vcid[] = "$Id: quadratic.c,v 1.7 2000/01/10 03:54:16 knepley Exp $";
  3: #endif

  5: /*
  6:    Defines piecewise quadratic function space on a two dimensional 
  7:    grid. Suitable for finite element type discretization of a PDE.
  8: */

 10: #include "src/grid/discretization/discimpl.h"         /*I "discretization.h" I*/
 11: #include "src/mesh/impls/triangular/triimpl.h"

 13: /* For precomputed integrals, the table is structured as follows:

 15:      precompInt[op,i,j] = int_{SE} <op phi^i(xi,eta), phi^j(xi,eta)> |J^{-1}|

 17:    The Jacobian in this case may not be constant over the element in question.
 18:    The numbering of the nodes in the standard element is

 20:                  3
 21:                  |
 22:                  | 
 23:                  5  4
 24:                  |   
 25:                  |    
 26:                  1--6--2
 27: */

 29: static int DiscDestroy_Triangular_2D_Quadratic(Discretization disc) {
 31:   return(0);
 32: }

 34: static int DiscView_Triangular_2D_Quadratic_File(Discretization disc, PetscViewer viewer) {
 36:   PetscViewerASCIIPrintf(viewer, "Quadratic discretizationn");
 37:   PetscViewerASCIIPrintf(viewer, "    %d shape functions per componentn", disc->funcs);
 38:   PetscViewerASCIIPrintf(viewer, "    %d registered operatorsn", disc->numOps);
 39:   return(0);
 40: }

 42: static int DiscView_Triangular_2D_Quadratic(Discretization disc, PetscViewer viewer) {
 43:   PetscTruth isascii;
 44:   int        ierr;

 47:   PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_ASCII, &isascii);
 48:   if (isascii == PETSC_TRUE) {
 49:     DiscView_Triangular_2D_Quadratic_File(disc, viewer);
 50:   }
 51:   return(0);
 52: }

 54: int DiscEvaluateShapeFunctions_Triangular_2D_Quadratic_Private(double xi, double eta, double *coords, double *x,
 55:                                                                double *y, double *dxxi, double *dxet, double *dyxi, double *dyet)
 56: {
 57:   /* ASSUMPTION: The coordinates passed in are corrected for periodicity */
 59:   *x    = 0.0; *y    = 0.0;
 60:   *dxxi = 0.0; *dxet = 0.0;
 61:   *dyxi = 0.0; *dyet = 0.0;
 62:   /* phi^0: 1 - 3 (xi + eta) + 2 (xi + eta)^2 */
 63:   *x    += coords[0*2+0]*(1.0 - (xi + eta))*(1.0 - 2.0*(xi + eta));
 64:   *dxxi += coords[0*2+0]*(-3.0 + 4.0*(xi + eta));
 65:   *dxet += coords[0*2+0]*(-3.0 + 4.0*(xi + eta));
 66:   *y    += coords[0*2+1]*(1.0 - (xi + eta))*(1.0 - 2.0*(xi + eta));
 67:   *dyxi += coords[0*2+1]*(-3.0 + 4.0*(xi + eta));
 68:   *dyet += coords[0*2+1]*(-3.0 + 4.0*(xi + eta));
 69:   /* phi^1: xi (2xi - 1) */
 70:   *x    += coords[1*2+0]*(xi*(2.0*xi - 1.0));
 71:   *dxxi += coords[1*2+0]*(4.0*xi  - 1.0);
 72:   *dxet += 0.0;
 73:   *y    += coords[1*2+1]*(xi*(2.0*xi - 1.0));
 74:   *dyxi += coords[1*2+1]*(4.0*xi  - 1.0);
 75:   *dyet += 0.0;
 76:   /* phi^2: eta (2eta - 1) */
 77:   *x    += coords[2*2+0]*(eta*(2.0*eta - 1.0));
 78:   *dxxi += 0.0;
 79:   *dxet += coords[2*2+0]*(4.0*eta - 1.0);
 80:   *y    += coords[2*2+1]*(eta*(2.0*eta - 1.0));
 81:   *dyxi += 0.0;
 82:   *dyet += coords[2*2+1]*(4.0*eta - 1.0);
 83:   /* phi^3: 4 xi eta */
 84:   *x    += coords[3*2+0]*(4.0*xi*eta);
 85:   *dxxi += coords[3*2+0]*(4.0*eta);
 86:   *dxet += coords[3*2+0]*(4.0*xi);
 87:   *y    += coords[3*2+1]*(4.0*xi*eta);
 88:   *dyxi += coords[3*2+1]*(4.0*eta);
 89:   *dyet += coords[3*2+1]*(4.0*xi);
 90:   /* phi^4: 4 eta (1 - xi - eta) */
 91:   *x    += coords[4*2+0]*(4.0*eta*(1.0 - (xi + eta)));
 92:   *dxxi += coords[4*2+0]*(-4.0*eta);
 93:   *dxet += coords[4*2+0]*(-8.0*eta + 4.0*(1.0 - xi));
 94:   *y    += coords[4*2+1]*(4.0*eta*(1.0 - (xi + eta)));
 95:   *dyxi += coords[4*2+1]*(-4.0*eta);
 96:   *dyet += coords[4*2+1]*(-8.0*eta + 4.0*(1.0 - xi));
 97:   /* phi^5: 4 xi (1 - xi - eta) */
 98:   *x    += coords[5*2+0]*(4.0*xi*(1.0 - (xi + eta)));
 99:   *dxxi += coords[5*2+0]*(-8.0*xi  + 4.0*(1.0 - eta));
100:   *dxet += coords[5*2+0]*(-4.0*xi);
101:   *y    += coords[5*2+1]*(4.0*xi*(1.0 - (xi + eta)));
102:   *dyxi += coords[5*2+1]*(-8.0*xi  + 4.0*(1.0 - eta));
103:   *dyet += coords[5*2+1]*(-4.0*xi);
104:   PetscLogFlops(36+20+20+20+30+20);
105:   return(0);
106: }

108: int DiscTransformCoords_Triangular_2D_Quadratic(double x, double y, double *coords, double *newXi, double *newEta)
109: {
110:   /* ASSUMPTION: The coordinates passed in are corrected for periodicity */
111:   double xi, eta;     /* Canonical coordinates of the point */
112:   double dxix;        /* PartDer{xi}{x}  */
113:   double detx;        /* PartDer{eta}{x} */
114:   double dxiy;        /* PartDer{xi}{y}  */
115:   double dety;        /* PartDer{eta}{y} */
116:   double dxxi;        /* PartDer{x}{xi}  */
117:   double dxet;        /* PartDer{x}{eta} */
118:   double dyxi;        /* PartDer{y}{xi}  */
119:   double dyet;        /* PartDer{y}{eta} */
120:   double new_x;       /* x(xi,eta) */
121:   double new_y;       /* x(xi,eta) */
122:   double residual_x;  /* x(xi,eta) - x */
123:   double residual_y;  /* x(xi,eta) - y */
124:   double jac, invjac; /* The Jacobian determinant and its inverse */
125:   int    maxIters = 100;
126:   int    iter;
127:   int    ierr;

130:   /* We have to solve

132:        sum_f x(f) phi^f(xi,eta) = x
133:        sum_f y(f) phi^f(xi,eta) = y

135:      where f runs over nodes (each one has coordinates and a shape function). We
136:      will use Newton's method

138:      / sum_f x(f) PartDer{phi^f}{xi} sum_f x(f) PartDer{phi^f}{eta}  / dxi   = / sum_f x(f) phi^f(xi,eta) - x
139:       sum_f y(f) PartDer{phi^f}{xi} sum_f y(f) PartDer{phi^f}{eta} /  deta /    sum_f y(f) phi^f(xi,eta) - y/

141:      which can be rewritten more compactly as

143:      / PartDer{x}{xi} PartDer{x}{eta}   / dxi   = / x(xi,eta) - x 
144:       PartDer{y}{xi} PartDer{y}{eta} /   deta /    y(xi,eta) - y /

146:      The initial guess will be the linear solution.

148:      ASSUMPTION: The coordinates passed in are all on the same sheet as x,y
149:   */
150:   /* Form linear solution */
151:   dxxi = coords[1*2+0] - coords[0*2+0];
152:   dxet = coords[2*2+0] - coords[0*2+0];
153:   dyxi = coords[1*2+1] - coords[0*2+1];
154:   dyet = coords[2*2+1] - coords[0*2+1];
155:   jac  = PetscAbsReal(dxxi*dyet - dxet*dyxi);
156:   if (jac < 1.0e-14) SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
157:   invjac = 1/jac;
158:   dxix   =  dyet*invjac;
159:   dxiy   = -dxet*invjac;
160:   detx   = -dyxi*invjac;
161:   dety   =  dxxi*invjac;
162:   xi     = dxix*(x - coords[0*2+0]) + dxiy*(y - coords[0*2+1]);
163:   eta    = detx*(x - coords[0*2+0]) + dety*(y - coords[0*2+1]);
164:   for(iter = 0; iter < maxIters; iter++) {
165:     /* This is clumsy, but I can't think of anything better right now */
166:     DiscEvaluateShapeFunctions_Triangular_2D_Quadratic_Private(xi, eta, coords, &new_x, &new_y, &dxxi, &dxet, &dyxi, &dyet);
167: 

169:     /* Check for convergence -- I should maybe make the tolerance variable */
170:     residual_x = new_x - x;
171:     residual_y = new_y - y;
172:     if (PetscAbsReal(residual_x) + PetscAbsReal(residual_y) < 1.0e-6) break;

174:     /* Solve the system */
175:     jac = PetscAbsReal(dxxi*dyet - dxet*dyxi);
176:     if (jac < 1.0e-14) {
177:       iter = maxIters;
178:       break;
179:     }

181:     /* These are the elements of the inverse matrix */
182:     invjac = 1/jac;
183:     dxix   =  dyet*invjac;
184:     dxiy   = -dxet*invjac;
185:     detx   = -dyxi*invjac;
186:     dety   =  dxxi*invjac;
187:     xi    -= dxix*residual_x + dxiy*residual_y;
188:     eta   -= detx*residual_x + dety*residual_y;
189:   }
190:   if (iter == maxIters) {
191:     PetscLogInfo(PETSC_NULL, "DiscTransformCoords_Triangular_2D_Quadratic: Newton iteration did not convergen");
192:     PetscLogInfo(PETSC_NULL, "x: %g y: %g maxIters: %dn", x, y, maxIters);
193:     for(iter = 0; iter < 6; iter++) {
194:       PetscLogInfo(PETSC_NULL, "  x%d: %g y%d: %gn", iter, coords[iter*2+0], iter, coords[iter*2+1]);
195:     }
196:     /* Use linear interpolation */
197:     xi  = dxix*(x - coords[0*2+0]) + dxiy*(y - coords[0*2+1]);
198:     eta = detx*(x - coords[0*2+0]) + dety*(y - coords[0*2+1]);
199:   }

201:   *newXi  = xi;
202:   *newEta = eta;
203:   PetscLogFlops(7+15+19*iter);
204:   return(0);
205: }

207: static int DiscEvaluateFunctionGalerkin_Triangular_2D_Quadratic(Discretization disc, Mesh mesh, PointFunction f,
208:                                                                 PetscScalar alpha, int elem, PetscScalar *array, void *ctx)
209: {
210:   Mesh_Triangular *tri               = (Mesh_Triangular *) mesh->data;
211:   double          *nodes             = tri->nodes;
212:   int             *elements          = tri->faces;
213:   int              numCorners        = mesh->numCorners;
214:   int              comp              = disc->comp;              /* The number of components in this field */
215:   int              funcs             = disc->funcs;             /* The number of shape functions per component */
216:   PetscScalar     *funcVal           = disc->funcVal;           /* Function value at a quadrature point */
217:   int              numQuadPoints     = disc->numQuadPoints;     /* Number of points used for Gaussian quadrature */
218:   double          *quadWeights       = disc->quadWeights;       /* Weights in the standard element for Gaussian quadrature */
219:   double          *quadShapeFuncs    = disc->quadShapeFuncs;    /* Shape functions evaluated at quadrature points */
220:   double          *quadShapeFuncDers = disc->quadShapeFuncDers; /* Shape function derivatives at quadrature points */
221:   double           jac;                                         /* |J| for map to standard element */
222:   double           x, y;                                        /* The integration point */
223:   double           dxxi;                                        /* PartDer{x}{xi}  */
224:   double           dxet;                                        /* PartDer{x}{xi}  */
225:   double           dyxi;                                        /* PartDer{y}{eta} */
226:   double           dyet;                                        /* PartDer{y}{eta} */
227:   double           coords[MAX_CORNERS*2];
228:   int              rank = -1;
229:   int              i, j, k, func, p;
230: #ifdef PETSC_USE_BOPT_g
231:   PetscTruth       opt;
232: #endif
233:   int              ierr;

236:   MPI_Comm_rank(disc->comm, &rank);

238:   /* For dummy collective calls */
239:   if (array == PETSC_NULL) {
240:     for(p = 0; p < numQuadPoints; p++) {
241:       (*f)(0, 0, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, ctx);
242:     }
243:     return(0);
244:   }

246: #ifdef PETSC_USE_BOPT_g
247:   if ((elem < 0) || (elem >= mesh->part->numOverlapElements)) {
248:     SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE, "Invalid element %d should be in [0,%d)", elem, mesh->part->numOverlapElements);
249:   }
250: #endif
251:   /* Calculate the determinant of the inverse Jacobian of the map to the standard element */
252:   for(i = 0; i < numCorners; i++) {
253:     coords[i*2]   = nodes[elements[elem*numCorners+i]*2];
254:     coords[i*2+1] = nodes[elements[elem*numCorners+i]*2+1];
255:   }

257:   /* Check for constant jacobian here */
258:   if (PETSC_FALSE) {
259:     jac = PetscAbsReal((coords[2] - coords[0])*(coords[5] - coords[1]) - (coords[4] - coords[0])*(coords[3] - coords[1]));
260:     if (jac < 1.0e-14) {
261:       PetscPrintf(PETSC_COMM_SELF, "[%d], elem: %d x1: %g y1: %g x2: %g y2: %g x3: %g y3: %gn",
262:                   rank, elem, coords[0], coords[1], coords[2], coords[3], coords[4], coords[5]);
263:       SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
264:     }
265:   }
266: #ifdef PETSC_USE_BOPT_g
267:   PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
268:   if (opt == PETSC_TRUE) {
269:     PetscPrintf(PETSC_COMM_SELF, "[%d]elem: %d x1: %g y1: %g x2: %g y2: %g x3: %g y3: %gn",
270:                 rank, elem, coords[0], coords[1], coords[2], coords[3], coords[4], coords[5]);
271:     PetscPrintf(PETSC_COMM_SELF, "  x4: %g y4: %g x5: %g y5: %g x6: %g y6: %gn",
272:                 coords[6], coords[7], coords[8], coords[9], coords[10], coords[11]);
273:   }
274: #endif

276:   /* Calculate element vector entries by Gaussian quadrature */
277:   for(p = 0; p < numQuadPoints; p++)
278:   {
279:     /* x                    = sum^{funcs}_{f=1} x_f phi^f(p)
280:        PartDer{x}{xi}(p)  = sum^{funcs}_{f=1} x_f PartDer{phi^f(p)}{xi}
281:        PartDer{x}{eta}(p) = sum^{funcs}_{f=1} x_f PartDer{phi^f(p)}{eta}
282:        y                    = sum^{funcs}_{f=1} y_f phi^f(p)
283:        PartDer{y}{xi}(p)  = sum^{funcs}_{f=1} y_f PartDer{phi^f(p)}{xi}
284:        PartDer{y}{eta}(p) = sum^{funcs}_{f=1} y_f PartDer{phi^f(p)}{eta} */
285:     x    = 0.0; y    = 0.0;
286:     dxxi = 0.0; dxet = 0.0;
287:     dyxi = 0.0; dyet = 0.0;
288:     if (mesh->isPeriodic == PETSC_TRUE) {
289:       for(func = 0; func < funcs; func++)
290:       {
291:         x    += MeshPeriodicRelativeX(mesh, coords[func*2],   coords[0])*quadShapeFuncs[p*funcs+func];
292:         dxxi += MeshPeriodicRelativeX(mesh, coords[func*2],   coords[0])*quadShapeFuncDers[p*funcs*2+func*2];
293:         dxet += MeshPeriodicRelativeX(mesh, coords[func*2],   coords[0])*quadShapeFuncDers[p*funcs*2+func*2+1];
294:         y    += MeshPeriodicRelativeY(mesh, coords[func*2+1], coords[1])*quadShapeFuncs[p*funcs+func];
295:         dyxi += MeshPeriodicRelativeY(mesh, coords[func*2+1], coords[1])*quadShapeFuncDers[p*funcs*2+func*2];
296:         dyet += MeshPeriodicRelativeY(mesh, coords[func*2+1], coords[1])*quadShapeFuncDers[p*funcs*2+func*2+1];
297:       }
298:       x = MeshPeriodicX(mesh, x);
299:       y = MeshPeriodicY(mesh, y);
300:     } else {
301:       for(func = 0; func < funcs; func++)
302:       {
303:         x    += coords[func*2]  *quadShapeFuncs[p*funcs+func];
304:         dxxi += coords[func*2]  *quadShapeFuncDers[p*funcs*2+func*2];
305:         dxet += coords[func*2]  *quadShapeFuncDers[p*funcs*2+func*2+1];
306:         y    += coords[func*2+1]*quadShapeFuncs[p*funcs+func];
307:         dyxi += coords[func*2+1]*quadShapeFuncDers[p*funcs*2+func*2];
308:         dyet += coords[func*2+1]*quadShapeFuncDers[p*funcs*2+func*2+1];
309:       }
310:     }
311:     jac = PetscAbsReal(dxxi*dyet - dxet*dyxi);
312:     if (jac < 1.0e-14) {
313:       PetscPrintf(PETSC_COMM_SELF, "[%d]p: %d x1: %g y1: %g x2: %g y2: %g x3: %g y3: %gn",
314:                   rank, p, coords[0], coords[1], coords[2], coords[3], coords[4], coords[5]);
315:       SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
316:     }
317:     (*f)(1, comp, &x, &y, PETSC_NULL, funcVal, ctx);
318: #ifdef PETSC_USE_BOPT_g
319:     PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
320:     if (opt == PETSC_TRUE) {
321:       PetscPrintf(PETSC_COMM_SELF, "[%d]p: %d jac: %g", rank, p, jac);
322:      for(j = 0; j < comp; j++)
323:         PetscPrintf(PETSC_COMM_SELF, " func[%d]: %g", j, PetscRealPart(funcVal[j]));
324:       PetscPrintf(PETSC_COMM_SELF, "n");
325:     }
326: #endif

328:     for(i = 0, k = 0; i < funcs; i++) {
329:       for(j = 0; j < comp; j++, k++) {
330:         array[k] += alpha*funcVal[j]*quadShapeFuncs[p*funcs+i]*jac*quadWeights[p];
331: #ifdef PETSC_USE_BOPT_g
332:         PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
333:         if (opt == PETSC_TRUE) {
334:           PetscPrintf(PETSC_COMM_SELF, "[%d]  array[%d]: %gn", rank, k, PetscRealPart(array[k]));
335:         }
336: #endif
337:       }
338:     }
339:   }
340:   PetscLogFlops((3 + 12*funcs + 5*funcs*comp) * numQuadPoints);
341:   return(0);
342: }

344: static int DiscEvaluateOperatorGalerkin_Triangular_2D_Quadratic(Discretization disc, Mesh mesh, int elemSize,
345:                                                                 int rowStart, int colStart, int op, PetscScalar alpha,
346:                                                                 int elem, PetscScalar *field, PetscScalar *array, void *ctx)
347: {
348:   Mesh_Triangular *tri        = (Mesh_Triangular *) mesh->data;
349:   double          *nodes      = tri->nodes;          /* The node coordinates */
350:   int             *elements   = tri->faces;          /* The element corners */
351:   int              numCorners = mesh->numCorners;    /* The number of corners per element */
352:   Operator         oper       = disc->operators[op]; /* The operator to discretize */
353:   Discretization   test       = oper->test;          /* The space of test functions */
354:   OperatorFunction opFunc     = oper->opFunc;        /* Integrals of operators which depend on J */
355:   PetscScalar     *precompInt = oper->precompInt;    /* Precomputed integrals of the operator on shape functions */
356:   int              rowSize    = test->size;          /* The number of shape functions per element */
357:   int              colSize    = disc->size;          /* The number of test  functions per element */
358:   double           x21, x31, y21, y31;               /* Coordinates of the element, with point 1 at the origin */
359:   double           jac;                              /* |J| for map to standard element */
360:   double           coords[MAX_CORNERS*2];            /* Coordinates of the element */
361:   int              rank;
362:   int              i, j, f;
363:   int              ierr;

366:   MPI_Comm_rank(disc->comm, &rank);
367: #ifdef PETSC_USE_BOPT_g
368:   /* Check for valid operator */
369:   if ((op < 0) || (op >= disc->numOps) || (!disc->operators[op])) SETERRQ(PETSC_ERR_ARG_WRONG, "Invalid operator");
370: #endif

372:   if (precompInt != PETSC_NULL)
373:   {
374:     /* Calculate the determinant of the inverse Jacobian of the map to the standard element
375:        which has been specified as constant here - 1/|x_{21} y_{31} - x_{31} y_{21}| */
376:     if (mesh->isPeriodic == PETSC_TRUE) {
377:       x21 = MeshPeriodicDiffX(mesh, nodes[elements[elem*numCorners+1]*2]   - nodes[elements[elem*numCorners]*2]);
378:       x31 = MeshPeriodicDiffX(mesh, nodes[elements[elem*numCorners+2]*2]   - nodes[elements[elem*numCorners]*2]);
379:       y21 = MeshPeriodicDiffY(mesh, nodes[elements[elem*numCorners+1]*2+1] - nodes[elements[elem*numCorners]*2+1]);
380:       y31 = MeshPeriodicDiffY(mesh, nodes[elements[elem*numCorners+2]*2+1] - nodes[elements[elem*numCorners]*2+1]);
381:     } else {
382:       x21 = nodes[elements[elem*numCorners+1]*2]   - nodes[elements[elem*numCorners]*2];
383:       x31 = nodes[elements[elem*numCorners+2]*2]   - nodes[elements[elem*numCorners]*2];
384:       y21 = nodes[elements[elem*numCorners+1]*2+1] - nodes[elements[elem*numCorners]*2+1];
385:       y31 = nodes[elements[elem*numCorners+2]*2+1] - nodes[elements[elem*numCorners]*2+1];
386:     }
387:     jac = PetscAbsReal(x21*y31 - x31*y21);
388:     if (jac < 1.0e-14) {
389:       PetscPrintf(PETSC_COMM_SELF, "[%d]x21: %g y21: %g x31: %g y31: %g jac: %gn", rank, x21, y21, x31, y31, jac);
390:       SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
391:     }
392:     /* PetscPrintf(PETSC_COMM_SELF, "x21: %g y21: %g x31: %g y31: %gn", x21, y21, x31, y31, jac); */

394:     /* Calculate element matrix entries which are all precomputed */
395:     for(i = 0; i < rowSize; i++)
396:       for(j = 0; j < colSize; j++)
397:         array[(i+rowStart)*elemSize + j+colStart] += alpha*precompInt[i*colSize + j]*jac;
398:     PetscLogFlops(7 + 2*rowSize*colSize);
399:   }
400:   else
401:   {
402:     if (opFunc == PETSC_NULL) SETERRQ(PETSC_ERR_ARG_CORRUPT, "Invalid function");
403:     if (mesh->isPeriodic == PETSC_TRUE) {
404:       coords[0*2+0] = nodes[elements[elem*numCorners+0]*2+0];
405:       coords[0*2+1] = nodes[elements[elem*numCorners+0]*2+1];
406:       for(f = 1; f < PetscMax(disc->funcs, test->funcs); f++) {
407:         coords[f*2+0] = MeshPeriodicRelativeX(mesh, nodes[elements[elem*numCorners+f]*2+0], coords[0*2+0]);
408:         coords[f*2+1] = MeshPeriodicRelativeY(mesh, nodes[elements[elem*numCorners+f]*2+1], coords[0*2+1]);
409:       }
410:     } else {
411:       for(f = 0; f < PetscMax(disc->funcs, test->funcs); f++) {
412:         coords[f*2+0] = nodes[elements[elem*numCorners+f]*2+0];
413:         coords[f*2+1] = nodes[elements[elem*numCorners+f]*2+1];
414:       }
415:     }

417:     (*opFunc)(disc, test, rowSize, colSize, rowStart, colStart, elemSize, coords, alpha, field, array, ctx);
418: 
419:   }
420:   return(0);
421: }

423: static int DiscEvaluateNonlinearOperatorGalerkin_Triangular_2D_Quadratic(Discretization disc, Mesh mesh, NonlinearOperator f,
424:                                                                          PetscScalar alpha, int elem, int numArgs, PetscScalar **field,
425:                                                                          PetscScalar *vec, void *ctx)
426: {
427:   Mesh_Triangular *tri        = (Mesh_Triangular *) mesh->data;
428:   double          *nodes      = tri->nodes;
429:   int             *elements   = tri->faces;
430:   int              numCorners = mesh->numCorners;
431:   int              comp       = disc->comp;     /* The number of components in this field */
432:   int              funcs      = disc->funcs;    /* The number of shape functions per component */
433:   PetscScalar     *funcVal    = disc->funcVal;  /* Function value at a quadrature point */
434:   PetscScalar    **fieldVal   = disc->fieldVal; /* Field value and derivatives at a quadrature point */
435:   double           jac;                         /* |J| for map to standard element */
436:   double           invjac;                      /* |J^{-1}| for map from standard element */
437:   int              numQuadPoints;               /* Number of points used for Gaussian quadrature */
438:   double          *quadWeights;                 /* Weights in the standard element for Gaussian quadrature */
439:   double          *quadShapeFuncs;              /* Shape functions evaluated at quadrature points */
440:   double          *quadShapeFuncDers;           /* Shape function derivatives evaluated at quadrature points */
441:   double           x, y;                        /* The integration point */
442:   double           dxxi;                        /* PartDer{x}{xi}  */
443:   double           dxet;                        /* PartDer{x}{eta}  */
444:   double           dyxi;                        /* PartDer{y}{xi} */
445:   double           dyet;                        /* PartDer{y}{eta} */
446:   double           dxix;                        /* PartDer{xi}{x}  */
447:   double           detx;                        /* PartDer{eta}{x} */
448:   double           dxiy;                        /* PartDer{xi}{y}  */
449:   double           dety;                        /* PartDer{eta}{y} */
450:   PetscScalar      dfxi;                        /* PartDer{field}{xi}  */
451:   PetscScalar      dfet;                        /* PartDer{field}{eta} */
452:   double           coords[12];                  /* Coordinates of the element */
453:   int              rank = -1;
454:   int              i, j, k, func, p, arg;
455: #ifdef PETSC_USE_BOPT_g
456:   PetscTruth       opt;
457: #endif
458:   int              ierr;

461:   MPI_Comm_rank(disc->comm, &rank);
462:   numQuadPoints     = disc->numQuadPoints;
463:   quadWeights       = disc->quadWeights;
464:   quadShapeFuncs    = disc->quadShapeFuncs;
465:   quadShapeFuncDers = disc->quadShapeFuncDers;

467:   /* Calculate the determinant of the inverse Jacobian of the map to the standard element */
468:   if (mesh->isPeriodic == PETSC_TRUE) {
469:     coords[0*2+0] = nodes[elements[elem*numCorners+0]*2+0];
470:     coords[0*2+1] = nodes[elements[elem*numCorners+0]*2+1];
471:     for(func = 1; func < funcs; func++) {
472:       coords[func*2+0] = MeshPeriodicRelativeX(mesh, nodes[elements[elem*numCorners+func]*2+0], coords[0*2+0]);
473:       coords[func*2+1] = MeshPeriodicRelativeY(mesh, nodes[elements[elem*numCorners+func]*2+1], coords[0*2+1]);
474:     }
475:   } else {
476:     for(func = 0; func < funcs; func++) {
477:       coords[func*2+0] = nodes[elements[elem*numCorners+func]*2+0];
478:       coords[func*2+1] = nodes[elements[elem*numCorners+func]*2+1];
479:     }
480:   }
481:   /* Check for constant jacobian here */
482:   if (PETSC_FALSE) {
483:     jac = PetscAbsReal((coords[2] - coords[0])*(coords[5] - coords[1]) - (coords[4] - coords[0])*(coords[3] - coords[1]));
484:     if (jac < 1.0e-14) {
485:       PetscPrintf(PETSC_COMM_SELF, "[%d], elem: %d x1: %g y1: %g x2: %g y2: %g x3: %g y3: %gn",
486:                   rank, elem, coords[0], coords[1], coords[2], coords[3], coords[4], coords[5]);
487:       SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
488:     }
489:   }
490: #ifdef PETSC_USE_BOPT_g
491:   PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
492:   if (opt == PETSC_TRUE) {
493:     PetscPrintf(PETSC_COMM_SELF, "[%d]elem: %d x1: %g y1: %g x2: %g y2: %g x3: %g y3: %gn",
494:                 rank, elem, coords[0], coords[1], coords[2], coords[3], coords[4], coords[5]);
495:     PetscPrintf(PETSC_COMM_SELF, "  x4: %g y4: %g x5: %g y5: %g x6: %g y6: %gn",
496:                 coords[6], coords[7], coords[8], coords[9], coords[10], coords[11]);
497:   }
498: #endif

500:   /* Calculate element vector entries by Gaussian quadrature */
501:   for(p = 0; p < numQuadPoints; p++) {
502:     /* x                      = sum^{funcs}_{f=1} x_f   phi^f(p)
503:        PartDer{x}{xi}(p)    = sum^{funcs}_{f=1} x_f   PartDer{phi^f(p)}{xi}
504:        PartDer{x}{eta}(p)   = sum^{funcs}_{f=1} x_f   PartDer{phi^f(p)}{eta}
505:        y                      = sum^{funcs}_{f=1} y_f   phi^f(p)
506:        PartDer{y}{xi}(p)    = sum^{funcs}_{f=1} y_f   PartDer{phi^f(p)}{xi}
507:        PartDer{y}{eta}(p)   = sum^{funcs}_{f=1} y_f   PartDer{phi^f(p)}{eta}
508:        u^i                    = sum^{funcs}_{f=1} u^i_f phi^f(p)
509:        PartDer{u^i}{xi}(p)  = sum^{funcs}_{f=1} u^i_f PartDer{phi^f(p)}{xi}
510:        PartDer{u^i}{eta}(p) = sum^{funcs}_{f=1} u^i_f PartDer{phi^f(p)}{eta} */
511:     x    = 0.0; y    = 0.0;
512:     dxxi = 0.0; dyxi = 0.0;
513:     dxet = 0.0; dyet = 0.0;
514:     for(arg = 0; arg < numArgs; arg++)
515:       for(j = 0; j < comp*3; j++)
516:         fieldVal[arg][j] = 0.0;
517:     for(func = 0; func < funcs; func++) {
518:       x    += coords[func*2]  *quadShapeFuncs[p*funcs+func];
519:       dxxi += coords[func*2]  *quadShapeFuncDers[p*funcs*2+func*2];
520:       dxet += coords[func*2]  *quadShapeFuncDers[p*funcs*2+func*2+1];
521:       y    += coords[func*2+1]*quadShapeFuncs[p*funcs+func];
522:       dyxi += coords[func*2+1]*quadShapeFuncDers[p*funcs*2+func*2];
523:       dyet += coords[func*2+1]*quadShapeFuncDers[p*funcs*2+func*2+1];
524:       for(arg = 0; arg < numArgs; arg++) {
525:         for(j = 0; j < comp; j++) {
526:           fieldVal[arg][j*3]   += field[arg][func*comp+j]*quadShapeFuncs[p*funcs+func];
527:           fieldVal[arg][j*3+1] += field[arg][func*comp+j]*quadShapeFuncDers[p*funcs*2+func*2];
528:           fieldVal[arg][j*3+2] += field[arg][func*comp+j]*quadShapeFuncDers[p*funcs*2+func*2+1];
529:         }
530:       }
531:     }
532:     if (mesh->isPeriodic == PETSC_TRUE) {
533:       x = MeshPeriodicX(mesh, x);
534:       y = MeshPeriodicY(mesh, y);
535:     }
536:     jac = PetscAbsReal(dxxi*dyet - dxet*dyxi);
537:     if (jac < 1.0e-14) {
538:       PetscPrintf(PETSC_COMM_SELF, "[%d]p: %d x1: %g y1: %g x2: %g y2: %g x3: %g y3: %gn",
539:                   rank, p, coords[0], coords[1], coords[2], coords[3], coords[4], coords[5]);
540:       SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
541:     }
542:     /* These are the elements of the inverse matrix */
543:     invjac = 1/jac;
544:     dxix   =  dyet*invjac;
545:     dxiy   = -dxet*invjac;
546:     detx   = -dyxi*invjac;
547:     dety   =  dxxi*invjac;

549:     /* Convert the field derivatives to old coordinates */
550:     for(arg = 0; arg < numArgs; arg++)
551:       for(j = 0; j < comp; j++) {
552:         dfxi                 = fieldVal[arg][j*3+1];
553:         dfet                 = fieldVal[arg][j*3+2];
554:         fieldVal[arg][j*3+1] = dfxi*dxix + dfet*detx;
555:         fieldVal[arg][j*3+2] = dfxi*dxiy + dfet*dety;
556:       }

558:     (*f)(1, comp, &x, &y, PETSC_NULL, numArgs, fieldVal, funcVal, ctx);
559: #ifdef PETSC_USE_BOPT_g
560:     PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
561:     if (opt == PETSC_TRUE) {
562:       PetscPrintf(PETSC_COMM_SELF, "[%d]p: %d jac: %g", rank, p, jac);
563:       for(j = 0; j < comp; j++)
564:         PetscPrintf(PETSC_COMM_SELF, " func[%d]: %g", j, PetscRealPart(funcVal[j]));
565:       PetscPrintf(PETSC_COMM_SELF, "n");
566:     }
567: #endif

569:     for(i = 0, k = 0; i < funcs; i++) {
570:       for(j = 0; j < comp; j++, k++) {
571:         vec[k] += alpha*funcVal[j]*quadShapeFuncs[p*funcs+i]*jac*quadWeights[p];
572: #ifdef PETSC_USE_BOPT_g
573:         PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
574:         if (opt == PETSC_TRUE) {
575:           PetscPrintf(PETSC_COMM_SELF, "[%d]  vec[%d]: %gn", rank, k, PetscRealPart(vec[k]));
576:         }
577: #endif
578:       }
579:     }
580:   }
581:   PetscLogFlops(((12 + (6*numArgs + 5)*comp)*funcs + 8 + 6*numArgs*comp) * numQuadPoints);
582:   return(0);
583: }

585: static int DiscEvaluateALEOperatorGalerkin_Triangular_2D_Quadratic(Discretization disc, Mesh mesh, int elemSize,
586:                                                                    int rowStart, int colStart, int op, PetscScalar alpha,
587:                                                                    int elem, PetscScalar *field, PetscScalar *ALEfield, PetscScalar *array,
588:                                                                    void *ctx)
589: {
590:   Mesh_Triangular    *tri        = (Mesh_Triangular *) mesh->data;
591:   double             *nodes      = tri->nodes;
592:   int                *elements   = tri->faces;
593:   int                 numCorners = mesh->numCorners;
594:   Discretization      test;                  /* The space of test functions */
595:   Operator            oper;                  /* The operator to discretize */
596:   int                 rowSize;               /* The number of shape functions per element */
597:   int                 colSize;               /* The number of test  functions per element */
598:   ALEOperatorFunction opFunc;                /* Integrals of operators which depend on J */
599:   double              coords[MAX_CORNERS*2]; /* Coordinates of the element */
600:   int                 f;
601:   int                 ierr;

604: #ifdef PETSC_USE_BOPT_g
605:   /* Check for valid operator */
606:   if ((op < 0) || (op >= disc->numOps) || (!disc->operators[op])) SETERRQ(PETSC_ERR_ARG_WRONG, "Invalid operator");
607: #endif
608:   /* Get discretization info */
609:   oper    = disc->operators[op];
610:   opFunc  = oper->ALEOpFunc;
611:   test    = oper->test;
612:   rowSize = test->size;
613:   colSize = disc->size;

615:   if (opFunc == PETSC_NULL) SETERRQ(PETSC_ERR_ARG_CORRUPT, "Invalid function");
616:   if (mesh->isPeriodic == PETSC_TRUE) {
617:     coords[0*2+0] = nodes[elements[elem*numCorners+0]*2+0];
618:     coords[0*2+1] = nodes[elements[elem*numCorners+0]*2+1];
619:     for(f = 1; f < PetscMax(disc->funcs, test->funcs); f++) {
620:       coords[f*2+0] = MeshPeriodicRelativeX(mesh, nodes[elements[elem*numCorners+f]*2+0], coords[0*2+0]);
621:       coords[f*2+1] = MeshPeriodicRelativeY(mesh, nodes[elements[elem*numCorners+f]*2+1], coords[0*2+1]);
622:     }
623:   } else {
624:     for(f = 0; f < PetscMax(disc->funcs, test->funcs); f++) {
625:       coords[f*2+0] = nodes[elements[elem*numCorners+f]*2+0];
626:       coords[f*2+1] = nodes[elements[elem*numCorners+f]*2+1];
627:     }
628:   }

630:   (*opFunc)(disc, test, rowSize, colSize, rowStart, colStart, elemSize, coords, alpha, field, ALEfield, array, ctx);
631: 
632:   return(0);
633: }

635: static int DiscEvaluateNonlinearALEOperatorGalerkin_Triangular_2D_Quadratic(Discretization disc, Mesh mesh, NonlinearOperator f,
636:                                                                             PetscScalar alpha, int elem, int numArgs, PetscScalar **field,
637:                                                                             PetscScalar *ALEfield, PetscScalar *vec, void *ctx)
638: {
639:   Mesh_Triangular *tri        = (Mesh_Triangular *) mesh->data;
640:   double          *nodes      = tri->nodes;
641:   int             *elements   = tri->faces;
642:   int              numCorners = mesh->numCorners;
643:   int              comp       = disc->comp;     /* The number of components in this field */
644:   int              funcs      = disc->funcs;    /* The number of shape functions per component */
645:   PetscScalar     *funcVal    = disc->funcVal;  /* Function value at a quadrature point */
646:   PetscScalar    **fieldVal   = disc->fieldVal; /* Field value and derivatives at a quadrature point */
647:   double           jac;                         /* |J| for map to standard element */
648:   double           invjac;                      /* |J^{-1}| for map from standard element */
649:   int              numQuadPoints;               /* Number of points used for Gaussian quadrature */
650:   double          *quadWeights;                 /* Weights in the standard element for Gaussian quadrature */
651:   double          *quadShapeFuncs;              /* Shape functions evaluated at quadrature points */
652:   double          *quadShapeFuncDers;           /* Shape function derivatives evaluated at quadrature points */
653:   double           x, y;                        /* The integration point */
654:   double           dxxi;                        /* PartDer{x}{xi}  */
655:   double           dxet;                        /* PartDer{x}{eta}  */
656:   double           dyxi;                        /* PartDer{y}{xi} */
657:   double           dyet;                        /* PartDer{y}{eta} */
658:   double           dxix;                        /* PartDer{xi}{x}  */
659:   double           detx;                        /* PartDer{eta}{x} */
660:   double           dxiy;                        /* PartDer{xi}{y}  */
661:   double           dety;                        /* PartDer{eta}{y} */
662:   PetscScalar      dfxi;                        /* PartDer{field}{xi}  */
663:   PetscScalar      dfet;                        /* PartDer{field}{eta} */
664:   double           coords[12];                  /* Coordinates of the element */
665:   int              rank;
666:   int              i, j, k, func, p, arg;
667: #ifdef PETSC_USE_BOPT_g
668:   PetscTruth       opt;
669: #endif
670:   int              ierr;

673:   MPI_Comm_rank(disc->comm, &rank);
674:   numQuadPoints     = disc->numQuadPoints;
675:   quadWeights       = disc->quadWeights;
676:   quadShapeFuncs    = disc->quadShapeFuncs;
677:   quadShapeFuncDers = disc->quadShapeFuncDers;

679:   /* Calculate the determinant of the inverse Jacobian of the map to the standard element */
680:   if (mesh->isPeriodic == PETSC_TRUE) {
681:     coords[0*2+0] = nodes[elements[elem*numCorners+0]*2+0];
682:     coords[0*2+1] = nodes[elements[elem*numCorners+0]*2+1];
683:     for(func = 1; func < funcs; func++) {
684:       coords[func*2+0] = MeshPeriodicRelativeX(mesh, nodes[elements[elem*numCorners+func]*2+0], coords[0*2+0]);
685:       coords[func*2+1] = MeshPeriodicRelativeY(mesh, nodes[elements[elem*numCorners+func]*2+1], coords[0*2+1]);
686:     }
687:   } else {
688:     for(func = 0; func < funcs; func++) {
689:       coords[func*2+0] = nodes[elements[elem*numCorners+func]*2+0];
690:       coords[func*2+1] = nodes[elements[elem*numCorners+func]*2+1];
691:     }
692:   }
693:   /* Check for constant jacobian here */
694:   if (PETSC_FALSE) {
695:     jac = PetscAbsReal((coords[2] - coords[0])*(coords[5] - coords[1]) - (coords[4] - coords[0])*(coords[3] - coords[1]));
696:     if (jac < 1.0e-14) {
697:       PetscPrintf(PETSC_COMM_SELF, "[%d], elem: %d x1: %g y1: %g x2: %g y2: %g x3: %g y3: %gn",
698:                   rank, elem, coords[0], coords[1], coords[2], coords[3], coords[4], coords[5]);
699:       SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
700:     }
701:   }
702: #ifdef PETSC_USE_BOPT_g
703:   PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
704:   if (opt == PETSC_TRUE) {
705:     PetscPrintf(PETSC_COMM_SELF, "[%d]elem: %d x1: %g y1: %g x2: %g y2: %g x3: %g y3: %gn",
706:                 rank, elem, coords[0], coords[1], coords[2], coords[3], coords[4], coords[5]);
707:     PetscPrintf(PETSC_COMM_SELF, "  x4: %g y4: %g x5: %g y5: %g x6: %g y6: %gn",
708:                 coords[6], coords[7], coords[8], coords[9], coords[10], coords[11]);
709:   }
710: #endif

712:   /* Calculate element vector entries by Gaussian quadrature */
713:   for(p = 0; p < numQuadPoints; p++) {
714:     /* x                      = sum^{funcs}_{f=1} x_f   phi^f(p)
715:        PartDer{x}{xi}(p)    = sum^{funcs}_{f=1} x_f   PartDer{phi^f(p)}{xi}
716:        PartDer{x}{eta}(p)   = sum^{funcs}_{f=1} x_f   PartDer{phi^f(p)}{eta}
717:        y                      = sum^{funcs}_{f=1} y_f   phi^f(p)
718:        PartDer{y}{xi}(p)    = sum^{funcs}_{f=1} y_f   PartDer{phi^f(p)}{xi}
719:        PartDer{y}{eta}(p)   = sum^{funcs}_{f=1} y_f   PartDer{phi^f(p)}{eta}
720:        u^i                    = sum^{funcs}_{f=1} u^i_f phi^f(p)
721:        PartDer{u^i}{xi}(p)  = sum^{funcs}_{f=1} u^i_f PartDer{phi^f(p)}{xi}
722:        PartDer{u^i}{eta}(p) = sum^{funcs}_{f=1} u^i_f PartDer{phi^f(p)}{eta} */
723:     x    = 0.0; y    = 0.0;
724:     dxxi = 0.0; dyxi = 0.0;
725:     dxet = 0.0; dyet = 0.0;
726:     for(arg = 0; arg < numArgs; arg++)
727:       for(j = 0; j < comp*3; j++)
728:         fieldVal[arg][j] = 0.0;
729:     for(func = 0; func < funcs; func++)
730:     {
731:       x    += coords[func*2]  *quadShapeFuncs[p*funcs+func];
732:       dxxi += coords[func*2]  *quadShapeFuncDers[p*funcs*2+func*2];
733:       dxet += coords[func*2]  *quadShapeFuncDers[p*funcs*2+func*2+1];
734:       y    += coords[func*2+1]*quadShapeFuncs[p*funcs+func];
735:       dyxi += coords[func*2+1]*quadShapeFuncDers[p*funcs*2+func*2];
736:       dyet += coords[func*2+1]*quadShapeFuncDers[p*funcs*2+func*2+1];
737:       for(arg = 0; arg < numArgs; arg++) {
738:         for(j = 0; j < comp; j++) {
739:           fieldVal[arg][j*3]   += (field[arg][func*comp+j] - ALEfield[func*comp+j])*quadShapeFuncs[p*funcs+func];
740:           fieldVal[arg][j*3+1] += field[arg][func*comp+j]*quadShapeFuncDers[p*funcs*2+func*2];
741:           fieldVal[arg][j*3+2] += field[arg][func*comp+j]*quadShapeFuncDers[p*funcs*2+func*2+1];
742:         }
743:       }
744:     }
745:     if (mesh->isPeriodic == PETSC_TRUE) {
746:       x = MeshPeriodicX(mesh, x);
747:       y = MeshPeriodicY(mesh, y);
748:     }
749:     jac = PetscAbsReal(dxxi*dyet - dxet*dyxi);
750:     if (jac < 1.0e-14) {
751:       PetscPrintf(PETSC_COMM_SELF, "[%d]p: %d x1: %g y1: %g x2: %g y2: %g x3: %g y3: %gn",
752:                   rank, p, coords[0], coords[1], coords[2], coords[3], coords[4], coords[5]);
753:       SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
754:     }
755:     /* These are the elements of the inverse matrix */
756:     invjac = 1/jac;
757:     dxix   =  dyet*invjac;
758:     dxiy   = -dxet*invjac;
759:     detx   = -dyxi*invjac;
760:     dety   =  dxxi*invjac;

762:     /* Convert the field derivatives to old coordinates */
763:     for(arg = 0; arg < numArgs; arg++) {
764:       for(j = 0; j < comp; j++) {
765:         dfxi                 = fieldVal[arg][j*3+1];
766:         dfet                 = fieldVal[arg][j*3+2];
767:         fieldVal[arg][j*3+1] = dfxi*dxix + dfet*detx;
768:         fieldVal[arg][j*3+2] = dfxi*dxiy + dfet*dety;
769:       }
770:     }

772:     (*f)(1, comp, &x, &y, PETSC_NULL, numArgs, fieldVal, funcVal, ctx);
773: #ifdef PETSC_USE_BOPT_g
774:     PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
775:     if (opt == PETSC_TRUE) {
776:       PetscPrintf(PETSC_COMM_SELF, "[%d]p: %d jac: %g", rank, p, jac);
777:       for(j = 0; j < comp; j++)
778:         PetscPrintf(PETSC_COMM_SELF, " func[%d]: %g", j, PetscRealPart(funcVal[j]));
779:       PetscPrintf(PETSC_COMM_SELF, "n");
780:     }
781: #endif

783:     for(i = 0, k = 0; i < funcs; i++) {
784:       for(j = 0; j < comp; j++, k++) {
785:         vec[k] += alpha*funcVal[j]*quadShapeFuncs[p*funcs+i]*jac*quadWeights[p];
786: #ifdef PETSC_USE_BOPT_g
787:         PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
788:         if (opt == PETSC_TRUE) {
789:           PetscPrintf(PETSC_COMM_SELF, "[%d]  vec[%d]: %gn", rank, k, PetscRealPart(vec[k]));
790:         }
791: #endif
792:       }
793:     }
794:   }
795:   PetscLogFlops(((12 + (7*numArgs + 5)*comp)*funcs + 8 + 6*numArgs*comp) * numQuadPoints);
796:   return(0);
797: }

799: int Identity_Triangular_2D_Quadratic(Discretization disc, Discretization test, int rowSize, int colSize,
800:                                      int globalRowStart, int globalColStart, int globalSize, double *coords,
801:                                      PetscScalar alpha, PetscScalar *field, PetscScalar *array, void *ctx)
802: {
803:   int     numQuadPoints;     /* Number of points used for Gaussian quadrature */
804:   double *quadWeights;       /* Weights in the standard element for Gaussian quadrature */
805:   double *quadShapeFuncs;    /* Shape functions evaluated at quadrature points */
806:   double *quadTestFuncs;     /* Test  functions evaluated at quadrature points */
807:   double *quadShapeFuncDers; /* Shape function derivatives evaluated at quadrature points */
808:   double  dxxi;              /* PartDer{x}{xi}  */
809:   double  dxet;              /* PartDer{x}{eta} */
810:   double  dyxi;              /* PartDer{y}{xi}  */
811:   double  dyet;              /* PartDer{y}{eta} */
812:   double  jac;               /* |J| for map to standard element */
813:   int     comp;              /* The number of field components */
814:   int     funcs;             /* The number of shape functions */
815:   int     i, j, c, f, p;

818:   /* Calculate element matrix entries by Gaussian quadrature */
819:   comp              = disc->comp;
820:   funcs             = disc->funcs;
821:   numQuadPoints     = disc->numQuadPoints;
822:   quadWeights       = disc->quadWeights;
823:   quadShapeFuncs    = disc->quadShapeFuncs;
824:   quadTestFuncs     = test->quadShapeFuncs;
825:   quadShapeFuncDers = disc->quadShapeFuncDers;
826:   for(p = 0; p < numQuadPoints; p++)
827:   {
828:     /* PartDer{x}{xi}(p)  = sum^{funcs}_{f=1} x_f PartDer{phi^f(p)}{xi}
829:        PartDer{x}{eta}(p) = sum^{funcs}_{f=1} x_f PartDer{phi^f(p)}{eta}
830:        PartDer{y}{xi}(p)  = sum^{funcs}_{f=1} y_f PartDer{phi^f(p)}{xi}
831:        PartDer{y}{eta}(p) = sum^{funcs}_{f=1} y_f PartDer{phi^f(p)}{eta} */
832:     dxxi = 0.0; dxet = 0.0;
833:     dyxi = 0.0; dyet = 0.0;
834:     for(f = 0; f < funcs; f++)
835:     {
836:       dxxi += coords[f*2]  *quadShapeFuncDers[p*funcs*2+f*2];
837:       dxet += coords[f*2]  *quadShapeFuncDers[p*funcs*2+f*2+1];
838:       dyxi += coords[f*2+1]*quadShapeFuncDers[p*funcs*2+f*2];
839:       dyet += coords[f*2+1]*quadShapeFuncDers[p*funcs*2+f*2+1];
840:     }
841:     jac = PetscAbsReal(dxxi*dyet - dxet*dyxi);
842: #ifdef PETSC_USE_BOPT_g
843:     if (jac < 1.0e-14) {
844:       PetscPrintf(PETSC_COMM_SELF, "p: %d x1: %g y1: %g x2: %g y2: %g x3: %g y3: %gn",
845:                   p, coords[0], coords[1], coords[2], coords[3], coords[4], coords[5]);
846:       SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
847:     }
848: #endif

850:     for(i = 0; i < funcs; i++)
851:     {
852:       for(j = 0; j < funcs; j++)
853:       {
854:                                 for(c = 0; c < comp; c++)
855:                                 {
856:                                         array[(i*comp+c+globalRowStart)*globalSize + j*comp+c+globalColStart] +=
857:                                                 alpha*quadTestFuncs[p*funcs+i]*quadShapeFuncs[p*funcs+j]*jac*quadWeights[p];
858:                                         /* PetscPrintf(PETSC_COMM_SELF, "  array[%d,%d]: %gn", i*comp+c+globalRowStart, j*comp+c+globalColStart,
859:                                                                                         array[(i*comp+c+globalRowStart)*globalSize + j*comp+c+globalColStart]); */
860:                                 }
861:       }
862:     }
863:   }
864:   PetscLogFlops((8*funcs + 3 + 5*funcs*funcs*comp) * numQuadPoints);

866:   return(0);
867: }

869: int Laplacian_Triangular_2D_Quadratic(Discretization disc, Discretization test, int rowSize, int colSize,
870:                                       int globalRowStart, int globalColStart, int globalSize, double *coords,
871:                                       PetscScalar alpha, PetscScalar *field, PetscScalar *array, void *ctx)
872: {
873:   int     numQuadPoints;     /* Number of points used for Gaussian quadrature */
874:   double *quadWeights;       /* Weights in the standard element for Gaussian quadrature */
875:   double *quadShapeFuncDers; /* Shape function derivatives evaluated at quadrature points */
876:   double *quadTestFuncDers;  /* Test  function derivatives evaluated at quadrature points */
877:   double  dxxi;              /* PartDer{x}{xi}  */
878:   double  dxet;              /* PartDer{x}{eta} */
879:   double  dyxi;              /* PartDer{y}{xi}  */
880:   double  dyet;              /* PartDer{y}{eta} */
881:   double  dxix;              /* PartDer{xi}{x}  */
882:   double  detx;              /* PartDer{eta}{x} */
883:   double  dxiy;              /* PartDer{xi}{y}  */
884:   double  dety;              /* PartDer{eta}{y} */
885:   double  dphix;             /* PartDer{phi_i}{x} times PartDer{phi_j}{x} */
886:   double  dphiy;             /* PartDer{phi_i}{y} times PartDer{phi_j}{y} */
887:   double  jac;               /* |J| for map to standard element */
888:   double  invjac;            /* |J^{-1}| for map from standard element */
889:   int     comp;              /* The number of field components */
890:   int     funcs;             /* The number of shape functions */
891:   int     i, j, c, f, p;

894:   /* Calculate element matrix entries by Gaussian quadrature */
895:   comp              = disc->comp;
896:   funcs             = disc->funcs;
897:   numQuadPoints     = disc->numQuadPoints;
898:   quadWeights       = disc->quadWeights;
899:   quadShapeFuncDers = disc->quadShapeFuncDers;
900:   quadTestFuncDers  = test->quadShapeFuncDers;
901:   for(p = 0; p < numQuadPoints; p++) {
902:     /* PartDer{x}{xi}(p)  = sum^{funcs}_{f=1} x_f PartDer{phi^f(p)}{xi}
903:        PartDer{x}{eta}(p) = sum^{funcs}_{f=1} x_f PartDer{phi^f(p)}{eta}
904:        PartDer{y}{xi}(p)  = sum^{funcs}_{f=1} y_f PartDer{phi^f(p)}{xi}
905:        PartDer{y}{eta}(p) = sum^{funcs}_{f=1} y_f PartDer{phi^f(p)}{eta} */
906:     dxxi = 0.0; dxet = 0.0;
907:     dyxi = 0.0; dyet = 0.0;
908:     for(f = 0; f < funcs; f++) {
909:       dxxi += coords[f*2]  *quadShapeFuncDers[p*funcs*2+f*2];
910:       dxet += coords[f*2]  *quadShapeFuncDers[p*funcs*2+f*2+1];
911:       dyxi += coords[f*2+1]*quadShapeFuncDers[p*funcs*2+f*2];
912:       dyet += coords[f*2+1]*quadShapeFuncDers[p*funcs*2+f*2+1];
913:     }
914:     jac  = PetscAbsReal(dxxi*dyet - dxet*dyxi);
915: #ifdef PETSC_USE_BOPT_g
916:     if (jac < 1.0e-14) {
917:       PetscPrintf(PETSC_COMM_SELF, "p: %d x1: %g y1: %g x2: %g y2: %g x3: %g y3: %gn",
918:                   p, coords[0], coords[1], coords[2], coords[3], coords[4], coords[5]);
919:       SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
920:     }
921: #endif
922:     /* These are the elements of the inverse matrix */
923:     invjac =  1.0/jac;
924:     dxix   =  dyet*invjac;
925:     dxiy   = -dxet*invjac;
926:     detx   = -dyxi*invjac;
927:     dety   =  dxxi*invjac;

929:     for(i = 0; i < funcs; i++) {
930:       for(j = 0; j < funcs; j++) {
931:         dphix = (quadTestFuncDers[p*funcs*2+i*2]*dxix + quadTestFuncDers[p*funcs*2+i*2+1]*detx)*
932:           (quadShapeFuncDers[p*funcs*2+j*2]*dxix + quadShapeFuncDers[p*funcs*2+j*2+1]*detx);
933:         dphiy = (quadTestFuncDers[p*funcs*2+i*2]*dxiy + quadTestFuncDers[p*funcs*2+i*2+1]*dety)*
934:           (quadShapeFuncDers[p*funcs*2+j*2]*dxiy + quadShapeFuncDers[p*funcs*2+j*2+1]*dety);
935:         for(c = 0; c < comp; c++) {
936:           array[(i*comp+c+globalRowStart)*globalSize + j*comp+c+globalColStart] +=
937:             -alpha*(dphix + dphiy)*jac*quadWeights[p];
938:           /* PetscPrintf(PETSC_COMM_SELF, "  array[%d,%d]: %gn", i*comp+c+globalRowStart, j*comp+c+globalColStart,
939:              array[(i*comp+c+globalRowStart)*globalSize + j*comp+c+globalColStart]); */
940:         }
941:       }
942:     }
943:   }
944:   PetscLogFlops((8*funcs + 8 + 19*funcs*funcs*comp) * numQuadPoints);

946:   return(0);
947: }

949: int Weighted_Laplacian_Triangular_2D_Quadratic(Discretization disc, Discretization test, int rowSize, int colSize,
950:                                                int globalRowStart, int globalColStart, int globalSize, double *coords,
951:                                                PetscScalar alpha, PetscScalar *field, PetscScalar *array, void *ctx)
952: {
953:   int     numQuadPoints;     /* Number of points used for Gaussian quadrature */
954:   double *quadWeights;       /* Weights in the standard element for Gaussian quadrature */
955:   double *quadShapeFuncDers; /* Shape function derivatives evaluated at quadrature points */
956:   double *quadTestFuncDers;  /* Test  function derivatives evaluated at quadrature points */
957:   double  dxxi;              /* PartDer{x}{xi}  */
958:   double  dxet;              /* PartDer{x}{eta} */
959:   double  dyxi;              /* PartDer{y}{xi}  */
960:   double  dyet;              /* PartDer{y}{eta} */
961:   double  dxix;              /* PartDer{xi}{x}  */
962:   double  detx;              /* PartDer{eta}{x} */
963:   double  dxiy;              /* PartDer{xi}{y}  */
964:   double  dety;              /* PartDer{eta}{y} */
965:   double  dphix;             /* PartDer{phi_i}{x} times PartDer{phi_j}{x} */
966:   double  dphiy;             /* PartDer{phi_i}{y} times PartDer{phi_j}{y} */
967:   double  jac;               /* |J| for map to standard element */
968:   double  invjac;            /* |J^{-1}| for map from standard element */
969:         int     comp;              /* The number of field components */
970:         int     funcs;             /* The number of shape functions */
971:   int     i, j, c, f, p;

973:   /* Each element is weighted by its Jacobian. This is supposed to make smaller elements "stiffer". */
975:   /* Calculate element matrix entries by Gaussian quadrature */
976:         comp              = disc->comp;
977:         funcs             = disc->funcs;
978:   numQuadPoints     = disc->numQuadPoints;
979:   quadWeights       = disc->quadWeights;
980:   quadShapeFuncDers = disc->quadShapeFuncDers;
981:   quadTestFuncDers  = test->quadShapeFuncDers;
982:   for(p = 0; p < numQuadPoints; p++)
983:   {
984:     /* PartDer{x}{xi}(p)  = sum^{funcs}_{f=1} x_f PartDer{phi^f(p)}{xi}
985:        PartDer{x}{eta}(p) = sum^{funcs}_{f=1} x_f PartDer{phi^f(p)}{eta}
986:        PartDer{y}{xi}(p)  = sum^{funcs}_{f=1} y_f PartDer{phi^f(p)}{xi}
987:        PartDer{y}{eta}(p) = sum^{funcs}_{f=1} y_f PartDer{phi^f(p)}{eta} */
988:     dxxi = 0.0; dxet = 0.0;
989:     dyxi = 0.0; dyet = 0.0;
990:     for(f = 0; f < funcs; f++)
991:     {
992:       dxxi += coords[f*2]  *quadShapeFuncDers[p*funcs*2+f*2];
993:       dxet += coords[f*2]  *quadShapeFuncDers[p*funcs*2+f*2+1];
994:       dyxi += coords[f*2+1]*quadShapeFuncDers[p*funcs*2+f*2];
995:       dyet += coords[f*2+1]*quadShapeFuncDers[p*funcs*2+f*2+1];
996:     }
997:     jac  = PetscAbsReal(dxxi*dyet - dxet*dyxi);
998: #ifdef PETSC_USE_BOPT_g
999:     if (jac < 1.0e-14) {
1000:       PetscPrintf(PETSC_COMM_SELF, "p: %d x1: %g y1: %g x2: %g y2: %g x3: %g y3: %gn",
1001:                   p, coords[0], coords[1], coords[2], coords[3], coords[4], coords[5]);
1002:       SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
1003:     }
1004: #endif
1005:                 /* These are the elements of the inverse matrix */
1006:     invjac =  1.0/jac;
1007:                 dxix   =  dyet*invjac;
1008:                 dxiy   = -dxet*invjac;
1009:                 detx   = -dyxi*invjac;
1010:                 dety   =  dxxi*invjac;

1012:     for(i = 0; i < funcs; i++)
1013:     {
1014:       for(j = 0; j < funcs; j++)
1015:       {
1016:         dphix = (quadTestFuncDers[p*funcs*2+i*2]*dxix + quadTestFuncDers[p*funcs*2+i*2+1]*detx)*
1017:           (quadShapeFuncDers[p*funcs*2+j*2]*dxix + quadShapeFuncDers[p*funcs*2+j*2+1]*detx);
1018:         dphiy = (quadTestFuncDers[p*funcs*2+i*2]*dxiy + quadTestFuncDers[p*funcs*2+i*2+1]*dety)*
1019:           (quadShapeFuncDers[p*funcs*2+j*2]*dxiy + quadShapeFuncDers[p*funcs*2+j*2+1]*dety);
1020:                                 for(c = 0; c < comp; c++)
1021:                                 {
1022:                                         array[(i*comp+c+globalRowStart)*globalSize + j*comp+c+globalColStart] +=
1023:                                                 -alpha*(dphix + dphiy)*quadWeights[p];
1024:                                 }
1025:       }
1026:     }
1027:   }
1028:   PetscLogFlops((8*funcs + 8 + 18*funcs*funcs*comp) * numQuadPoints);

1030:   return(0);
1031: }

1033: int Divergence_Triangular_2D_Quadratic(Discretization disc, Discretization test, int rowSize, int colSize,
1034:                                        int globalRowStart, int globalColStart, int globalSize, double *coords,
1035:                                        PetscScalar alpha, PetscScalar *field, PetscScalar *array, void *ctx)
1036: {
1037:         /* We are using the convention that

1039:                    nabla matrix{v_1 cr v_2 cr vdots cr v_n} =
1040:                            matrix{v^{(1)}_1 cr vdots cr v^{(d)}_1 cr v^{(1)}_2 cr vdots cr v^{(d)}_n}

1042:                  and

1044:                    nabla cdot matrix{v^{(1)}_1 cr vdots cr v^{(d)}_1 cr v^{(1)}_2 cr vdots cr v^{(d)}_n} =
1045:                            matrix{v_1 cr v_2 cr vdots cr v_n}

1047:                  where $d$ is the number of space dimensions. This agrees with the convention which allows
1048:                  $Delta matrix{u_1 cr u_2} = 0$ to denote a set of scalar equations        This also requires that
1049:      the dimension of a vector must be divisible by the space dimension in order to be acted upon by
1050:      the divergence operator */
1051:   int     numQuadPoints;     /* Number of points used for Gaussian quadrature */
1052:   double *quadWeights;       /* Weights in the standard element for Gaussian quadrature */
1053:   double *quadTestFuncs;     /* Test  functions evaluated at quadrature points */
1054:   double *quadShapeFuncDers; /* Shape function derivatives evaluated at quadrature points */
1055:   double  dxxi;              /* PartDer{x}{xi}  */
1056:   double  dxet;              /* PartDer{x}{eta} */
1057:   double  dyxi;              /* PartDer{y}{xi}  */
1058:   double  dyet;              /* PartDer{y}{eta} */
1059:   double  dxix;              /* PartDer{xi}{x}  */
1060:   double  detx;              /* PartDer{eta}{x} */
1061:   double  dxiy;              /* PartDer{xi}{y}  */
1062:   double  dety;              /* PartDer{eta}{y} */
1063:   double  dphix;             /* PartDer{phi_i}{x} times PartDer{phi_j}{x} */
1064:   double  dphiy;             /* PartDer{phi_i}{y} times PartDer{phi_j}{y} */
1065:   double  jac;               /* |J| for map to standard element */
1066:   double  invjac;            /* |J^{-1}| for map from standard element */
1067:   int     comp;              /* The number of field components */
1068:   int     tcomp;             /* The number of field components for the test field */
1069:   int     funcs;             /* The number of shape functions */
1070:   int     tfuncs;            /* The number of test functions */
1071:   int     i, j, c, tc, f, p;

1074:   /* Calculate element matrix entries by Gaussian quadrature */
1075:   comp              = disc->comp;
1076:   tcomp             = test->comp;
1077:   funcs             = disc->funcs;
1078:   tfuncs            = test->funcs;
1079:   numQuadPoints     = disc->numQuadPoints;
1080:   quadWeights       = disc->quadWeights;
1081:   quadTestFuncs     = test->quadShapeFuncs;
1082:   quadShapeFuncDers = disc->quadShapeFuncDers;
1083:   for(p = 0; p < numQuadPoints; p++) {
1084:     /* PartDer{x}{xi}(p)  = sum^{funcs}_{f=1} x_f PartDer{phi^f(p)}{xi}
1085:        PartDer{x}{eta}(p) = sum^{funcs}_{f=1} x_f PartDer{phi^f(p)}{eta}
1086:        PartDer{y}{xi}(p)  = sum^{funcs}_{f=1} y_f PartDer{phi^f(p)}{xi}
1087:        PartDer{y}{eta}(p) = sum^{funcs}_{f=1} y_f PartDer{phi^f(p)}{eta} */
1088:     dxxi = 0.0; dxet = 0.0;
1089:     dyxi = 0.0; dyet = 0.0;
1090:     for(f = 0; f < funcs; f++)
1091:     {
1092:       dxxi += coords[f*2]  *quadShapeFuncDers[p*funcs*2+f*2];
1093:       dxet += coords[f*2]  *quadShapeFuncDers[p*funcs*2+f*2+1];
1094:       dyxi += coords[f*2+1]*quadShapeFuncDers[p*funcs*2+f*2];
1095:       dyet += coords[f*2+1]*quadShapeFuncDers[p*funcs*2+f*2+1];
1096:     }
1097:     jac  = PetscAbsReal(dxxi*dyet - dxet*dyxi);
1098: #ifdef PETSC_USE_BOPT_g
1099:     if (jac < 1.0e-14) {
1100:       PetscPrintf(PETSC_COMM_SELF, "p: %d x1: %g y1: %g x2: %g y2: %g x3: %g y3: %gn",
1101:                                                                         p, coords[0], coords[1], coords[2], coords[3], coords[4], coords[5]);
1102:       SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
1103:     }
1104: #endif
1105:                 /* PetscPrintf(PETSC_COMM_SELF, "p: %d x1: %g y1: %g x2: %g y2: %g x3: %g y3: %gn",
1106:                                                                 p, coords[0], coords[1], coords[2], coords[3], coords[4], coords[5]); */
1107:                 /* These are the elements of the inverse matrix */
1108:     invjac =  1.0/jac;
1109:     dxix   =  dyet*invjac;
1110:     dxiy   = -dxet*invjac;
1111:     detx   = -dyxi*invjac;
1112:     dety   =  dxxi*invjac;

1114:                 /* The rows are test functions */
1115:     for(i = 0; i < tfuncs; i++)
1116:     {
1117:       for(tc = 0; tc < tcomp; tc++)
1118:       {
1119:         /* The columns are shape functions */
1120:         for(j = 0; j < funcs; j++)
1121:         {
1122:           dphix = quadShapeFuncDers[p*funcs*2+j*2]*dxix + quadShapeFuncDers[p*funcs*2+j*2+1]*detx;
1123:           dphiy = quadShapeFuncDers[p*funcs*2+j*2]*dxiy + quadShapeFuncDers[p*funcs*2+j*2+1]*dety;
1124:           /* We divide by the number of space dimensions */
1125:           for(c = 0; c < comp/2; c++)
1126:           {
1127:             array[(i*tcomp+tc+globalRowStart)*globalSize + j*comp+c*2+globalColStart] +=
1128:               alpha*dphix*quadTestFuncs[p*tfuncs+i]*jac*quadWeights[p];
1129:                                                         /* PetscPrintf(PETSC_COMM_SELF, "  array[%d,%d]: %gn", i*tcomp+tc+globalRowStart, j*comp+c*2+globalColStart,
1130:                                                                                                         array[(i*tcomp+tc+globalRowStart)*globalSize + j*comp+c*2+globalColStart]); */
1131:             array[(i*tcomp+tc+globalRowStart)*globalSize + j*comp+c*2+1+globalColStart] +=
1132:               alpha*dphiy*quadTestFuncs[p*tfuncs+i]*jac*quadWeights[p];
1133:                                                         /* PetscPrintf(PETSC_COMM_SELF, "  array[%d,%d]: %gn", i*tcomp+tc+globalRowStart, j*comp+c*2+1+globalColStart,
1134:                                                                                                         array[(i*tcomp+tc+globalRowStart)*globalSize + j*comp+c*2+1+globalColStart]); */
1135:                                         }
1136:                                 }
1137:                         }
1138:     }
1139:   }
1140:   PetscLogFlops((8*funcs + 8 + 8*tfuncs*tcomp*funcs*comp) * numQuadPoints);

1142:   return(0);
1143: }

1145: int DiscInterpolateField_Triangular_2D_Quadratic(Discretization disc, Mesh oldMesh, int elem, double x, double y, double z,
1146:                                                  PetscScalar *oldFieldVal, PetscScalar *newFieldVal, InterpolationType type)
1147: {
1148:   Mesh_Triangular *tri        = (Mesh_Triangular *) oldMesh->data;
1149:   int              numCorners = oldMesh->numCorners;
1150:   int             *elements   = tri->faces;
1151:   int             *neighbors  = tri->neighbors;
1152:   double          *nodes      = tri->nodes;
1153:   double           coords[24];  /* Coordinates of our "big element" */
1154:   double           xi, eta;     /* Canonical coordinates of the point */
1155:   double           x21, x31, y21, y31, jac, invjac, dx, dy, dxix, dxiy, detx, dety, xiOld, etaOld;
1156:   int              comp = disc->comp;
1157:   int              neighbor, corner, nelem, node, c;
1158:   int              ierr;

1161:   /* No scheme in place for boundary elements */
1162:   for(neighbor = 0; neighbor < 3; neighbor++)
1163:     if (neighbors[elem*3+neighbor] < 0)
1164:       type = INTERPOLATION_LOCAL;

1166:   switch (type)
1167:   {
1168:   case INTERPOLATION_LOCAL:
1169:     if (oldMesh->isPeriodic == PETSC_TRUE) {
1170:       coords[0*2+0] = MeshPeriodicRelativeX(oldMesh, nodes[elements[elem*numCorners+0]*2+0], x);
1171:       coords[0*2+1] = MeshPeriodicRelativeY(oldMesh, nodes[elements[elem*numCorners+0]*2+1], y);
1172:       coords[1*2+0] = MeshPeriodicRelativeX(oldMesh, nodes[elements[elem*numCorners+1]*2+0], x);
1173:       coords[1*2+1] = MeshPeriodicRelativeY(oldMesh, nodes[elements[elem*numCorners+1]*2+1], y);
1174:       coords[2*2+0] = MeshPeriodicRelativeX(oldMesh, nodes[elements[elem*numCorners+2]*2+0], x);
1175:       coords[2*2+1] = MeshPeriodicRelativeY(oldMesh, nodes[elements[elem*numCorners+2]*2+1], y);
1176:       coords[3*2+0] = MeshPeriodicRelativeX(oldMesh, nodes[elements[elem*numCorners+3]*2+0], x);
1177:       coords[3*2+1] = MeshPeriodicRelativeY(oldMesh, nodes[elements[elem*numCorners+3]*2+1], y);
1178:       coords[4*2+0] = MeshPeriodicRelativeX(oldMesh, nodes[elements[elem*numCorners+4]*2+0], x);
1179:       coords[4*2+1] = MeshPeriodicRelativeY(oldMesh, nodes[elements[elem*numCorners+4]*2+1], y);
1180:       coords[5*2+0] = MeshPeriodicRelativeX(oldMesh, nodes[elements[elem*numCorners+5]*2+0], x);
1181:       coords[5*2+1] = MeshPeriodicRelativeY(oldMesh, nodes[elements[elem*numCorners+5]*2+1], y);
1182:     } else {
1183:       coords[0*2+0] = nodes[elements[elem*numCorners+0]*2+0];
1184:       coords[0*2+1] = nodes[elements[elem*numCorners+0]*2+1];
1185:       coords[1*2+0] = nodes[elements[elem*numCorners+1]*2+0];
1186:       coords[1*2+1] = nodes[elements[elem*numCorners+1]*2+1];
1187:       coords[2*2+0] = nodes[elements[elem*numCorners+2]*2+0];
1188:       coords[2*2+1] = nodes[elements[elem*numCorners+2]*2+1];
1189:       coords[3*2+0] = nodes[elements[elem*numCorners+3]*2+0];
1190:       coords[3*2+1] = nodes[elements[elem*numCorners+3]*2+1];
1191:       coords[4*2+0] = nodes[elements[elem*numCorners+4]*2+0];
1192:       coords[4*2+1] = nodes[elements[elem*numCorners+4]*2+1];
1193:       coords[5*2+0] = nodes[elements[elem*numCorners+5]*2+0];
1194:       coords[5*2+1] = nodes[elements[elem*numCorners+5]*2+1];
1195:     }
1196:     /* Get the (xi,eta) coordinates of the point */
1197:     DiscTransformCoords_Triangular_2D_Quadratic(x, y, coords, &xi, &eta);
1198:     if ((xi < -1.0e-02) || (eta < -1.0e-02) || (xi > 1.01) || (eta > 1.01)) {
1199:       xiOld  = xi;
1200:       etaOld = eta;
1201:       /* Use linear approximation */
1202:       x21    = coords[1*2+0] - coords[0*2+0];
1203:       x31    = coords[2*2+0] - coords[0*2+0];
1204:       y21    = coords[1*2+1] - coords[0*2+1];
1205:       y31    = coords[2*2+1] - coords[0*2+1];
1206:       dx     = x - coords[0*2+0];
1207:       dy     = y - coords[0*2+1];
1208:       jac    = PetscAbsReal(x21*y31 - x31*y21);
1209:       invjac = 1/jac;
1210:       dxix   =  y31*invjac;
1211:       dxiy   = -x31*invjac;
1212:       detx   = -y21*invjac;
1213:       dety   =  x21*invjac;

1215:       xi     = dxix*dx + dxiy*dy;
1216:       eta    = detx*dx + dety*dy;
1217:       PetscPrintf(PETSC_COMM_SELF, "elem: %d x: %g y: %g xiOld: %g etaOld: %g xi: %g eta: %gn", elem, x, y, xiOld, etaOld, xi, eta);
1218:     }
1219:     for(c = 0; c < comp; c++) {
1220:       newFieldVal[c] = oldFieldVal[0*comp+c]*(1.0 - xi - eta)*(1.0 - 2.0*xi - 2.0*eta) +
1221:         oldFieldVal[1*comp+c]*xi *(2.0*xi  - 1.0)      +
1222:         oldFieldVal[2*comp+c]*eta*(2.0*eta - 1.0)      +
1223:         oldFieldVal[3*comp+c]*4.0*xi*eta               +
1224:         oldFieldVal[4*comp+c]*4.0*eta*(1.0 - xi - eta) +
1225:         oldFieldVal[5*comp+c]*4.0*xi *(1.0 - xi - eta);
1226:     }
1227:     PetscLogFlops(34*comp);
1228:     break;
1229:   case INTERPOLATION_HALO:
1230:     /* Here is our "big element" where numbers in parantheses represent
1231:        the numbering on the old little element:

1233:            2
1234:            |
1235:            | 
1236:            |  
1237:            6   5
1238:            |    
1239:            |     
1240:            |      
1241:        (1) 7---*---4 (0)
1242:            |      |
1243:            |      | 
1244:            |      |  
1245:            8   *   *   3
1246:            |      |    
1247:            |      |     
1248:            |      |      
1249:            0---9--10--11---1
1250:                   (2)

1252:        We search for the neighbor node by looking for the vertex not a member of the original element.
1253:     */
1254:     for(neighbor = 0; neighbor < 3; neighbor++)
1255:     {
1256:       nelem = neighbors[elem*3+neighbor];
1257:       for(corner = 0; corner < 3; corner++)
1258:       {
1259:         node = elements[nelem*numCorners+corner];
1260:         if ((node != elements[elem*numCorners+((neighbor+1)%3)]) && (node != elements[elem*numCorners+((neighbor+2)%3)]))
1261:         {
1262:           /* The neighboring elements give the vertices */
1263:           coords[neighbor*2]   = nodes[node*2];
1264:           coords[neighbor*2+1] = nodes[node*2+1];
1265:           break;
1266:         }
1267:       }
1268:     }
1269:     /* Element vertices form midnodes */
1270:     coords[3*2]   = nodes[elements[elem*numCorners]*2];
1271:     coords[3*2+1] = nodes[elements[elem*numCorners]*2+1];
1272:     coords[4*2]   = nodes[elements[elem*numCorners+1]*2];
1273:     coords[4*2+1] = nodes[elements[elem*numCorners+1]*2+1];
1274:     coords[5*2]   = nodes[elements[elem*numCorners+2]*2];
1275:     coords[5*2+1] = nodes[elements[elem*numCorners+2]*2+1];
1276:     /* Treat 4 triangles as one big element with quadratic shape functions */
1277:     SETERRQ(PETSC_ERR_SUP, "Unsupported interpolation type");
1278:   default:
1279:     SETERRQ(PETSC_ERR_ARG_WRONG, "Unknown interpolation type");
1280:   }
1281: 
1282:   return(0);
1283: }

1285: int DiscInterpolateElementVec_Triangular_2D_Quadratic(Discretization disc, ElementVec vec, Discretization newDisc, ElementVec newVec)
1286: {
1287:   int          comp  = disc->comp;
1288:   int          size  = disc->size;
1289:   PetscScalar *array, *newArray;
1290:   PetscTruth   islin, isquad;
1291:   int          f, c;
1292:   int          ierr;

1295:   ElementVecGetArray(vec,    &array);
1296:   ElementVecGetArray(newVec, &newArray);
1297:   PetscTypeCompare((PetscObject) newDisc, DISCRETIZATION_TRIANGULAR_2D_LINEAR,    &islin);
1298:   PetscTypeCompare((PetscObject) newDisc, DISCRETIZATION_TRIANGULAR_2D_QUADRATIC, &isquad);
1299:   if (isquad == PETSC_TRUE) {
1300:     PetscMemcpy(newArray, array, size * sizeof(PetscScalar));
1301:   } else if (islin == PETSC_TRUE) {
1302:     for(f = 0; f < newDisc->funcs; f++) {
1303:       for(c = 0; c < comp; c++) {
1304:         newArray[f*comp+c] = array[f*comp+c];
1305:       }
1306:     }
1307:   } else {
1308:     SETERRQ(PETSC_ERR_SUP, "Discretization not supported");
1309:   }
1310:   ElementVecRestoreArray(vec,    &array);
1311:   ElementVecRestoreArray(newVec, &newArray);
1312:   return(0);
1313: }

1315: /*
1316:   DiscSetupQuadrature_Triangular_2D_Quadratic - Setup Gaussian quadrature with a 7 point integration rule

1318:   Input Parameter:
1319: . disc - The Discretization
1320: */
1321: int DiscSetupQuadrature_Triangular_2D_Quadratic(Discretization disc) {
1322:   int    dim   = disc->dim;
1323:   int    funcs = disc->funcs;
1324:   double xi, eta;
1325:   int    p;
1326:   int    ierr;

1329:   disc->numQuadPoints = 7;
1330:   PetscMalloc(disc->numQuadPoints*dim       * sizeof(double), &disc->quadPoints);
1331:   PetscMalloc(disc->numQuadPoints           * sizeof(double), &disc->quadWeights);
1332:   PetscMalloc(disc->numQuadPoints*funcs     * sizeof(double), &disc->quadShapeFuncs);
1333:   PetscMalloc(disc->numQuadPoints*funcs*dim * sizeof(double), &disc->quadShapeFuncDers);
1334:   PetscLogObjectMemory(disc, (disc->numQuadPoints*(funcs*(dim+1) + dim+1)) * sizeof(double));
1335:   disc->quadPoints[0]  = 1.0/3.0;
1336:   disc->quadPoints[1]  = disc->quadPoints[0];
1337:   disc->quadWeights[0] = 0.11250000000000;
1338:   disc->quadPoints[2]  = 0.797426985353087;
1339:   disc->quadPoints[3]  = 0.101286507323456;
1340:   disc->quadWeights[1] = 0.0629695902724135;
1341:   disc->quadPoints[4]  = disc->quadPoints[3];
1342:   disc->quadPoints[5]  = disc->quadPoints[2];
1343:   disc->quadWeights[2] = disc->quadWeights[1];
1344:   disc->quadPoints[6]  = disc->quadPoints[4];
1345:   disc->quadPoints[7]  = disc->quadPoints[3];
1346:   disc->quadWeights[3] = disc->quadWeights[1];
1347:   disc->quadPoints[8]  = 0.470142064105115;
1348:   disc->quadPoints[9]  = 0.059715871789770;
1349:   disc->quadWeights[4] = 0.066197076394253;
1350:   disc->quadPoints[10] = disc->quadPoints[8];
1351:   disc->quadPoints[11] = disc->quadPoints[8];
1352:   disc->quadWeights[5] = disc->quadWeights[4];
1353:   disc->quadPoints[12] = disc->quadPoints[9];
1354:   disc->quadPoints[13] = disc->quadPoints[8];
1355:   disc->quadWeights[6] = disc->quadWeights[4];
1356:   for(p = 0; p < disc->numQuadPoints; p++) {
1357:     xi  = disc->quadPoints[p*2];
1358:     eta = disc->quadPoints[p*2+1];
1359:     /* phi^0: 1 - 3 (xi + eta) + 2 (xi + eta)^2 */
1360:     disc->quadShapeFuncs[p*funcs]            =  1.0 - 3.0*(xi + eta) + 2.0*(xi + eta)*(xi + eta);
1361:     disc->quadShapeFuncDers[p*funcs*2+0*2]   = -3.0 + 4.0*(xi + eta);
1362:     disc->quadShapeFuncDers[p*funcs*2+0*2+1] = -3.0 + 4.0*(xi + eta);
1363:     /* phi^1: xi (2xi - 1) */
1364:     disc->quadShapeFuncs[p*funcs+1]          =  xi*(2.0*xi - 1.0);
1365:     disc->quadShapeFuncDers[p*funcs*2+1*2]   =  4.0*xi - 1.0;
1366:     disc->quadShapeFuncDers[p*funcs*2+1*2+1] =  0.0;
1367:     /* phi^2: eta (2eta - 1) */
1368:     disc->quadShapeFuncs[p*funcs+2]          =  eta*(2.0*eta - 1.0);
1369:     disc->quadShapeFuncDers[p*funcs*2+2*2]   =  0.0;
1370:     disc->quadShapeFuncDers[p*funcs*2+2*2+1] =  4.0*eta - 1.0;
1371:     /* phi^3: 4 xi eta */
1372:     disc->quadShapeFuncs[p*funcs+3]          =  4.0*xi*eta;
1373:     disc->quadShapeFuncDers[p*funcs*2+3*2]   =  4.0*eta;
1374:     disc->quadShapeFuncDers[p*funcs*2+3*2+1] =  4.0*xi;
1375:     /* phi^4: 4 eta (1 - xi - eta) */
1376:     disc->quadShapeFuncs[p*funcs+4]          =  4.0*eta*(1.0 - xi - eta);
1377:     disc->quadShapeFuncDers[p*funcs*2+4*2]   = -4.0*eta;
1378:     disc->quadShapeFuncDers[p*funcs*2+4*2+1] = -8.0*eta + 4.0*(1.0 - xi);
1379:     /* phi^5: 4 xi (1 - xi - eta) */
1380:     disc->quadShapeFuncs[p*funcs+5]          =  4.0*xi*(1.0 - xi - eta);
1381:     disc->quadShapeFuncDers[p*funcs*2+5*2]   = -8.0*xi + 4.0*(1.0 - eta);
1382:     disc->quadShapeFuncDers[p*funcs*2+5*2+1] = -4.0*xi;
1383:   }
1384:   return(0);
1385: }

1387: /*
1388:   DiscSetupOperators_Triangular_2D_Quadratic - Setup the default operators

1390:   Input Parameter:
1391: . disc - The Discretization
1392: */
1393: int DiscSetupOperators_Triangular_2D_Quadratic(Discretization disc) {
1394:   int newOp;

1398:   /* The Identity operator I -- the matrix is symmetric */
1399:   DiscretizationRegisterOperator(disc, Identity_Triangular_2D_Quadratic, &newOp);
1400:   if (newOp != IDENTITY) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "Default operator %d not setup correctly", IDENTITY);
1401:   /* The Laplacian operator Delta -- the matrix is symmetric */
1402:   DiscretizationRegisterOperator(disc, Laplacian_Triangular_2D_Quadratic, &newOp);
1403:   if (newOp != LAPLACIAN) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "Default operator %d not setup correctly", LAPLACIAN);
1404:   /* The Gradient operator nabla -- the matrix is rectangular */
1405:   DiscretizationRegisterOperator(disc, PETSC_NULL, &newOp);
1406:   if (newOp != GRADIENT) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "Default operator %d not setup correctly", GRADIENT);
1407:   /* The Divergence operator nablacdot -- the matrix is rectangular */
1408:   DiscretizationRegisterOperator(disc, Divergence_Triangular_2D_Quadratic, &newOp);
1409:   if (newOp != DIVERGENCE) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "Default operator %d not setup correctly", DIVERGENCE);
1410:   /* The weighted Laplacian operator -- the matrix is symmetric */
1411:   DiscretizationRegisterOperator(disc, Weighted_Laplacian_Triangular_2D_Quadratic, &newOp);
1412:   if (newOp != WEIGHTED_LAP) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "Default operator %d not setup correctly", WEIGHTED_LAP);
1413:   return(0);
1414: }

1416: static struct _DiscretizationOps DOps = {PETSC_NULL/* DiscretizationSetup */,
1417:                                          DiscSetupOperators_Triangular_2D_Quadratic,
1418:                                          PETSC_NULL/* DiscretizationSetFromOptions */,
1419:                                          DiscView_Triangular_2D_Quadratic,
1420:                                          DiscDestroy_Triangular_2D_Quadratic,
1421:                                          DiscEvaluateFunctionGalerkin_Triangular_2D_Quadratic,
1422:                                          DiscEvaluateOperatorGalerkin_Triangular_2D_Quadratic,
1423:                                          DiscEvaluateALEOperatorGalerkin_Triangular_2D_Quadratic,
1424:                                          DiscEvaluateNonlinearOperatorGalerkin_Triangular_2D_Quadratic,
1425:                                          DiscEvaluateNonlinearALEOperatorGalerkin_Triangular_2D_Quadratic,
1426:                                          DiscInterpolateField_Triangular_2D_Quadratic,
1427:                                          DiscInterpolateElementVec_Triangular_2D_Quadratic};

1429: EXTERN_C_BEGIN
1430: int DiscCreate_Triangular_2D_Quadratic(Discretization disc) {
1431:   int    arg;
1432:   int    ierr;

1435:   if (disc->comp <= 0) {
1436:     SETERRQ(PETSC_ERR_ARG_WRONG, "Discretization must have at least 1 component. Call DiscretizationSetNumComponents() to set this.");
1437:   }
1438:   PetscMemcpy(disc->ops, &DOps, sizeof(struct _DiscretizationOps));
1439:   disc->dim   = 2;
1440:   disc->funcs = 6;
1441:   disc->size  = disc->funcs*disc->comp;

1443:   DiscretizationSetupDefaultOperators(disc);
1444:   DiscSetupQuadrature_Triangular_2D_Quadratic(disc);

1446:   DiscretizationCreate(disc->comm, &disc->bdDisc);
1447:   DiscretizationSetNumComponents(disc->bdDisc, disc->comp);
1448:   DiscretizationSetType(disc->bdDisc, BD_DISCRETIZATION_TRIANGULAR_2D_QUADRATIC);

1450:   /* Storage */
1451:   PetscMalloc(disc->comp * sizeof(PetscScalar),   &disc->funcVal);
1452:   PetscMalloc(2          * sizeof(PetscScalar *), &disc->fieldVal);
1453:   for(arg = 0; arg < 2; arg++) {
1454:     PetscMalloc(disc->comp*(disc->dim+1) * sizeof(PetscScalar), &disc->fieldVal[arg]);
1455:   }
1456:   return(0);
1457: }
1458: EXTERN_C_END