Actual source code: bdquadratic.c

  1: #ifdef PETSC_RCS_HEADER
  2: static char vcid[] = "$Id: bdquadratic.c,v 1.4 2000/01/31 17:56:20 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), phi^j(xi)> |J^{-1}|

 17:    where we recall that |J| is a constant for linear affine maps,
 18:    and the map of any triangle to the standard element is linear.
 19:    The numbering of the nodes in the standard element is

 21:                  1--3--2
 22: */

 24: #undef  __FUNCT__
 26: static int BoundaryDiscDestroy_Triangular_2D_Quadratic(Discretization disc) {
 28:   return(0);
 29: }

 31: #undef  __FUNCT__
 33: static int BoundaryDiscView_Triangular_2D_Quadratic_File(Discretization disc, PetscViewer viewer) {
 35:   PetscViewerASCIIPrintf(viewer, "Quadratic boundary discretizationn");
 36:   PetscViewerASCIIPrintf(viewer, "    %d shape functions per componentn", disc->funcs);
 37:   PetscViewerASCIIPrintf(viewer, "    %d registered operatorsn", disc->numOps);
 38:   return(0);
 39: }

 41: #undef  __FUNCT__
 43: static int BoundaryDiscView_Triangular_2D_Quadratic(Discretization disc, PetscViewer viewer) {
 44:   PetscTruth isascii;
 45:   int        ierr;

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

 55: #undef  __FUNCT__
 57: static int BoundaryDiscEvaluateFunctionGalerkin_Triangular_2D_Quadratic(Discretization disc, Mesh mesh, PointFunction f,
 58:                                                                         PetscScalar alpha, int edge, PetscScalar *array, void *ctx)
 59: {
 60:   Mesh_Triangular *tri            = (Mesh_Triangular *) mesh->data;
 61:   EdgeContext     *bdCtx          = (EdgeContext *)         ctx;
 62:   void            *uCtx           = bdCtx->ctx;
 63:   double          *nodes          = tri->nodes;
 64:   int             *edges          = tri->edges;
 65:   int              comp           = disc->comp;                 /* The number of components in this field */
 66:   int              funcs          = disc->funcs;                /* The number of shape functions per component */
 67:   PetscScalar     *funcVal        = disc->funcVal;              /* Function value at a quadrature point */
 68:   int              numQuadPoints  = disc->numQuadPoints;        /* Number of points used for Gaussian quadrature */
 69:   double          *quadWeights    = disc->quadWeights;          /* Weights in the standard element for Gaussian quadrature */
 70:   double          *quadShapeFuncs = disc->quadShapeFuncs;       /* Shape function evaluated at quadrature points */
 71:   double          *quadShapeFuncDers = disc->quadShapeFuncDers; /* Shape function derivatives at quadrature points */
 72:   double           jac;                                         /* |J| for map to standard element */
 73:   double           x, y;                                        /* The integration point */
 74:   double           dxxi;                                        /* PartDer{x}{xi}  */
 75:   double           dyxi;                                        /* PartDer{y}{eta} */
 76:   double           coords[3*2];
 77:   int              rank;
 78:   int              i, j, k, p, func;
 79: #ifdef PETSC_USE_BOPT_g
 80:   PetscTruth       opt;
 81: #endif
 82:   int              ierr;

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

 87:   /* For dummy collective calls */
 88:   if (array == PETSC_NULL) {
 89:     for(p = 0; p < numQuadPoints; p++) {
 90:       (*f)(0, 0, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, ctx);
 91:     }
 92:     return(0);
 93:   }

 95: #ifdef PETSC_USE_BOPT_g
 96:   if ((edge < 0) || (edge >= mesh->numEdges)) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Invalid edge");
 97: #endif
 98:   /* Calculate the determinant of the inverse Jacobian of the map to the standard element */
 99:   if (mesh->isPeriodic == PETSC_TRUE) {
100:     coords[0*2]   = nodes[edges[edge*2+0]*2];
101:     coords[0*2+1] = nodes[edges[edge*2+0]*2+1];
102:     coords[1*2+0] = MeshPeriodicRelativeX(mesh, nodes[edges[edge*2+1]*2],   coords[0*2+0]);
103:     coords[1*2+1] = MeshPeriodicRelativeY(mesh, nodes[edges[edge*2+1]*2+1], coords[0*2+1]);
104:     coords[2*2+0] = MeshPeriodicRelativeX(mesh, nodes[bdCtx->midnode *2],   coords[0*2+0]);
105:     coords[2*2+1] = MeshPeriodicRelativeY(mesh, nodes[bdCtx->midnode *2+1], coords[0*2+1]);
106:   } else {
107:     coords[0*2]   = nodes[edges[edge*2+0]*2];
108:     coords[0*2+1] = nodes[edges[edge*2+0]*2+1];
109:     coords[1*2]   = nodes[edges[edge*2+1]*2];
110:     coords[1*2+1] = nodes[edges[edge*2+1]*2+1];
111:     coords[2*2]   = nodes[bdCtx->midnode *2];
112:     coords[2*2+1] = nodes[bdCtx->midnode *2+1];
113:   }
114: #ifdef PETSC_USE_BOPT_g
115:   PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
116:   if (opt == PETSC_TRUE) {
117:     PetscPrintf(PETSC_COMM_SELF, "[%d]edge: %d x1: %g y1: %g x2: %g y2: %g x3: %g y3: %gn",
118:                 rank, edge, coords[0], coords[1], coords[2], coords[3], coords[4], coords[5]);
119:   }
120: #endif

122:   /* Calculate element vector entries by Gaussian quadrature */
123:   for(p = 0; p < numQuadPoints; p++)
124:   {
125:     /* x                    = sum^{funcs}_{f=1} x_f phi^f(p)
126:        PartDer{x}{xi}(p)  = sum^{funcs}_{f=1} x_f PartDer{phi^f(p)}{xi}
127:        y                    = sum^{funcs}_{f=1} y_f phi^f(p)
128:        PartDer{y}{xi}(p)  = sum^{funcs}_{f=1} y_f PartDer{phi^f(p)}{xi} */
129:     x    = 0.0; y    = 0.0;
130:     dxxi = 0.0; dyxi = 0.0;
131:     for(func = 0; func < funcs; func++)
132:     {
133:       x    += coords[func*2]  *quadShapeFuncs[p*funcs+func];
134:       dxxi += coords[func*2]  *quadShapeFuncDers[p*funcs+func];
135:       y    += coords[func*2+1]*quadShapeFuncs[p*funcs+func];
136:       dyxi += coords[func*2+1]*quadShapeFuncDers[p*funcs+func];
137:     }
138:     if (mesh->isPeriodic == PETSC_TRUE) {
139:       x = MeshPeriodicX(mesh, x);
140:       y = MeshPeriodicY(mesh, y);
141:     }
142:     jac = sqrt(dxxi*dxxi + dyxi*dyxi);
143:     if (jac < 1.0e-14) {
144:       PetscPrintf(PETSC_COMM_SELF, "[%d]p: %d x1: %g y1: %g x2: %g y2: %g x3: %g y3: %gn",
145:                   rank, p, coords[0], coords[1], coords[2], coords[3], coords[4], coords[5]);
146:       SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
147:     }
148:     (*f)(1, comp, &x, &y, PETSC_NULL, funcVal, uCtx);
149: #ifdef PETSC_USE_BOPT_g
150:     PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
151:     if (opt == PETSC_TRUE) {
152:       PetscPrintf(PETSC_COMM_SELF, "[%d]p:%d jac: %g", rank, p, jac);
153:       for(j = 0; j < comp; j++)
154:         PetscPrintf(PETSC_COMM_SELF, " func[%d]: %g", j, PetscRealPart(funcVal[j]));
155:       PetscPrintf(PETSC_COMM_SELF, "n");
156:     }
157: #endif

159:     for(i = 0, k = 0; i < funcs; i++) {
160:       for(j = 0; j < comp; j++, k++) {
161:         array[k] += alpha*funcVal[j]*quadShapeFuncs[p*funcs+i]*jac*quadWeights[p];
162: #ifdef PETSC_USE_BOPT_g
163:         PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
164:         if (opt == PETSC_TRUE) {
165:           PetscPrintf(PETSC_COMM_SELF, "[%d]  array[%d]: %gn", rank, k, PetscRealPart(array[k]));
166:         }
167: #endif
168:       }
169:     }
170:   }
171:   PetscLogFlops((4 + 8*funcs + 5*funcs*comp) * numQuadPoints);
172:   return(0);
173: }

175: #undef  __FUNCT__
177: static int BoundaryDiscEvaluateOperatorGalerkin_Triangular_2D_Quadratic(Discretization disc, Mesh mesh, int elemSize,
178:                                                                      int rowStart, int colStart, int op, PetscScalar alpha,
179:                                                                      int edge, PetscScalar *field, PetscScalar *mat, void *ctx)
180: {
181:   Mesh_Triangular *tri        = (Mesh_Triangular *) mesh->data;
182:   EdgeContext     *bdCtx      = (EdgeContext *)         ctx;
183:   void            *uCtx       = bdCtx->ctx;
184:   double          *nodes      = tri->nodes;
185:   int             *edges      = tri->edges;
186:   Operator         oper       = disc->operators[op]; /* The operator to discretize */
187:   Discretization   test       = oper->test;          /* The space of test functions */
188:   int              rowSize    = test->size;          /* The number of shape functions per element */
189:   int              colSize    = disc->size;          /* The number of test  functions per element */
190:   OperatorFunction opFunc     = oper->opFunc;        /* Integrals of operators which depend on J */
191:   double           coords[MAX_CORNERS*2];            /* Coordinates of the element */
192:   int              ierr;

195: #ifdef PETSC_USE_BOPT_g
196:   /* Check for valid operator */
197:   if ((op < 0) || (op >= disc->numOps) || (!disc->operators[op])) SETERRQ(PETSC_ERR_ARG_WRONG, "Invalid operator");
198: #endif

200:   if (opFunc == PETSC_NULL) SETERRQ(PETSC_ERR_ARG_CORRUPT, "Invalid function");
201:   if (mesh->isPeriodic == PETSC_TRUE) {
202:     coords[0*2]   = nodes[edges[edge*2+0]*2];
203:     coords[0*2+1] = nodes[edges[edge*2+0]*2+1];
204:     coords[1*2+0] = MeshPeriodicRelativeX(mesh, nodes[edges[edge*2+1]*2],   coords[0*2+0]);
205:     coords[1*2+1] = MeshPeriodicRelativeY(mesh, nodes[edges[edge*2+1]*2+1], coords[0*2+1]);
206:     coords[2*2+0] = MeshPeriodicRelativeX(mesh, nodes[bdCtx->midnode *2],   coords[0*2+0]);
207:     coords[2*2+1] = MeshPeriodicRelativeY(mesh, nodes[bdCtx->midnode *2+1], coords[0*2+1]);
208:   } else {
209:     coords[0*2]   = nodes[edges[edge*2+0]*2];
210:     coords[0*2+1] = nodes[edges[edge*2+0]*2+1];
211:     coords[1*2]   = nodes[edges[edge*2+1]*2];
212:     coords[1*2+1] = nodes[edges[edge*2+1]*2+1];
213:     coords[2*2]   = nodes[bdCtx->midnode *2];
214:     coords[2*2+1] = nodes[bdCtx->midnode *2+1];
215:   }

217:   (*opFunc)(disc, test, rowSize, colSize, rowStart, colStart, elemSize, coords, alpha, field, mat, uCtx);
218: 

220:   return(0);
221: }

223: #undef  __FUNCT__
225: static int BoundaryDiscEvaluateALEOperatorGalerkin_Triangular_2D_Quadratic(Discretization disc, Mesh mesh, int elemSize,
226:                                                                         int rowStart, int colStart, int op, PetscScalar alpha,
227:                                                                         int elem, PetscScalar *field, PetscScalar *ALEfield, PetscScalar *mat,
228:                                                                         void *ctx)
229: {
230:   SETERRQ(PETSC_ERR_SUP, " ");
231: }

233: #undef  __FUNCT__
235: static int BoundaryDiscEvaluateNonlinearOperatorGalerkin_Triangular_2D_Quadratic(Discretization disc, Mesh mesh,
236:                                                                                  NonlinearOperator f, PetscScalar alpha, int elem,
237:                                                                                  int numArgs, PetscScalar **field, PetscScalar *vec, void *ctx)
238: {
239:   SETERRQ(PETSC_ERR_SUP, " ");
240: }

242: #undef  __FUNCT__
244: static int BoundaryDiscEvaluateNonlinearALEOperatorGalerkin_Triangular_2D_Quadratic(Discretization disc, Mesh mesh,
245:                                                                                     NonlinearOperator f, PetscScalar alpha, int elem,
246:                                                                                     int numArgs, PetscScalar **field, PetscScalar *ALEfield,
247:                                                                                     PetscScalar *vec, void *ctx)
248: {
249:   SETERRQ(PETSC_ERR_SUP, " ");
250: }

252: #undef  __FUNCT__
254: int BoundaryIdentity_Triangular_2D_Quadratic(Discretization disc, Discretization test, int rowSize, int colSize,
255:                                              int globalRowStart, int globalColStart, int globalSize, double *coords,
256:                                              PetscScalar alpha, PetscScalar *field, PetscScalar *array, void *ctx)
257: {
258:   int     numQuadPoints     = disc->numQuadPoints;     /* Number of points used for Gaussian quadrature */
259:   double *quadWeights       = disc->quadWeights;       /* Weights in the standard element for Gaussian quadrature */
260:   double *quadShapeFuncs    = disc->quadShapeFuncs;    /* Shape functions evaluated at quadrature points */
261:   double *quadTestFuncs     = test->quadShapeFuncs;    /* Test  functions evaluated at quadrature points */
262:   double *quadShapeFuncDers = disc->quadShapeFuncDers; /* Shape function derivatives evaluated at quadrature points */
263:   int     comp              = disc->comp;              /* The number of field components */
264:   int     tcomp             = test->comp;              /* The number of field components */
265:   int     funcs             = disc->funcs;             /* The number of shape functions */
266:   int     tfuncs            = test->funcs;             /* The number of test functions */
267:   double  dxxi;                                        /* PartDer{x}{xi}  */
268:   double  dyxi;                                        /* PartDer{y}{xi}  */
269:   double  jac;                                         /* |J| for map to standard element */
270:   int     i, j, c, f, p;

273: #ifdef PETSC_USE_BOPT_g
274:   if (comp != tcomp) SETERRQ(PETSC_ERR_ARG_INCOMP, "Shape and test fields must have the same number of components");
275: #endif
276:   /* Calculate element matrix entries by Gaussian quadrature */
277:   for(p = 0; p < numQuadPoints; p++) {
278:     /* PartDer{x}{xi}(p)  = sum^{funcs}_{f=1} x_f PartDer{phi^f(p)}{xi}
279:        PartDer{y}{xi}(p)  = sum^{funcs}_{f=1} y_f PartDer{phi^f(p)}{xi} */
280:     dxxi = 0.0; dyxi = 0.0;
281:     for(f = 0; f < funcs; f++) {
282:       dxxi += coords[f*2]  *quadShapeFuncDers[p*funcs+f];
283:       dyxi += coords[f*2+1]*quadShapeFuncDers[p*funcs+f];
284:     }
285:     jac = sqrt(dxxi*dxxi + dyxi*dyxi);
286: #ifdef PETSC_USE_BOPT_g
287:     if (jac < 1.0e-14) {
288:       PetscPrintf(PETSC_COMM_SELF, "p: %d x1: %g y1: %g x2: %g y2: %g x3: %g y3: %gn",
289:                   p, coords[0], coords[1], coords[2], coords[3], coords[4], coords[5]);
290:       SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
291:     }
292: #endif

294:     for(i = 0; i < tfuncs; i++) {
295:       for(j = 0; j < funcs; j++) {
296:         for(c = 0; c < comp; c++) {
297:           array[(i*tcomp+c+globalRowStart)*globalSize + j*comp+c+globalColStart] +=
298:             alpha*quadTestFuncs[p*tfuncs+i]*quadShapeFuncs[p*funcs+j]*jac*quadWeights[p];
299:         }
300:       }
301:     }
302:   }
303:   PetscLogFlops((4*funcs + 4 + 5*funcs*funcs*comp) * numQuadPoints);

305:   return(0);
306: }

308: #undef  __FUNCT__
310: int BoundaryDiscInterpolateField_Triangular_2D_Quadratic(Discretization disc, Mesh oldMesh, int elem, double x, double y, double z,
311:                                                       PetscScalar *oldFieldVal, PetscScalar *newFieldVal, InterpolationType type)
312: {
313:   SETERRQ(PETSC_ERR_SUP, " ");
314: }

316: #undef  __FUNCT__
318: int BoundaryDiscInterpolateElementVec_Triangular_2D_Quadratic(Discretization disc, ElementVec vec,
319:                                                               Discretization newDisc, ElementVec newVec)
320: {
321:   SETERRQ(PETSC_ERR_SUP, " ");
322: }

324: #undef  __FUNCT__
326: /*
327:   BoundaryDiscSetupQuadrature_Triangular_2D_Linear - Setup Gaussian quadrature with a 5 point integration rule

329:   Input Parameter:
330: . disc - The Discretization
331: */
332: int BoundaryDiscSetupQuadrature_Triangular_2D_Quadratic(Discretization disc) {
333:   int dim   = disc->dim;
334:   int funcs = disc->funcs;
335:   int p;

339:   disc->numQuadPoints = 5;
340:   PetscMalloc(disc->numQuadPoints*dim       * sizeof(double), &disc->quadPoints);
341:   PetscMalloc(disc->numQuadPoints           * sizeof(double), &disc->quadWeights);
342:   PetscMalloc(disc->numQuadPoints*funcs     * sizeof(double), &disc->quadShapeFuncs);
343:   PetscMalloc(disc->numQuadPoints*funcs*dim * sizeof(double), &disc->quadShapeFuncDers);
344:   PetscLogObjectMemory(disc, (disc->numQuadPoints*(funcs*(dim+1) + dim+1)) * sizeof(double));
345:   disc->quadPoints[0]  = 0.0469101;
346:   disc->quadWeights[0] = 0.118463;
347:   disc->quadPoints[1]  = 0.230765;
348:   disc->quadWeights[1] = 0.239314;
349:   disc->quadPoints[2]  = 0.500000;
350:   disc->quadWeights[2] = 0.284444;
351:   disc->quadPoints[3]  = 0.769235;
352:   disc->quadWeights[3] = 0.239314;
353:   disc->quadPoints[4]  = 0.953090;
354:   disc->quadWeights[4] = 0.118463;
355:   for(p = 0; p < disc->numQuadPoints; p++) {
356:     /* phi^0: (1 - xi) (1 - 2 xi) */
357:     disc->quadShapeFuncs[p*funcs+0]    = (1.0 - disc->quadPoints[p])*(1.0 - 2.0*disc->quadPoints[p]);
358:     disc->quadShapeFuncDers[p*funcs+0] =  4.0*disc->quadPoints[p] - 3.0;
359:     /* phi^1: -xi (1 - 2 xi) */
360:     disc->quadShapeFuncs[p*funcs+1]    = -disc->quadPoints[p]*(1.0 - 2.0*disc->quadPoints[p]);
361:     disc->quadShapeFuncDers[p*funcs+1] =  4.0*disc->quadPoints[p] - 1.0;
362:     /* phi^2: 4 xi (1 - xi) */
363:     disc->quadShapeFuncs[p*funcs+2]    =  4.0*disc->quadPoints[p]*(1.0 - disc->quadPoints[p]);
364:     disc->quadShapeFuncDers[p*funcs+2] =  4.0 - 8.0*disc->quadPoints[p];
365:   }
366:   return(0);
367: }

369: #undef  __FUNCT__
371: /*
372:   BoundaryDiscSetupOperators_Triangular_2D_Quadratic - Setup the default operators

374:   Input Parameter:
375: . disc - The Discretization
376: */
377: int BoundaryDiscSetupOperators_Triangular_2D_Quadratic(Discretization disc) {
378:   int newOp;

382:   /* The Identity operator I -- the matrix is symmetric */
383:   DiscretizationRegisterOperator(disc, BoundaryIdentity_Triangular_2D_Quadratic, &newOp);
384:   if (newOp != IDENTITY) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "Default operator %d not setup correctly", IDENTITY);
385:   /* The Laplacian operator Delta -- the matrix is symmetric */
386:   DiscretizationRegisterOperator(disc, PETSC_NULL, &newOp);
387:   if (newOp != LAPLACIAN) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "Default operator %d not setup correctly", LAPLACIAN);
388:   /* The Gradient operator nabla -- the matrix is rectangular */
389:   DiscretizationRegisterOperator(disc, PETSC_NULL, &newOp);
390:   if (newOp != GRADIENT) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "Default operator %d not setup correctly", GRADIENT);
391:   /* The Divergence operator nablacdot -- the matrix is rectangular */
392:   DiscretizationRegisterOperator(disc, PETSC_NULL, &newOp);
393:   if (newOp != DIVERGENCE) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "Default operator %d not setup correctly", DIVERGENCE);
394:   /* The weighted Laplacian operator -- the matrix is symmetric */
395:   DiscretizationRegisterOperator(disc, PETSC_NULL, &newOp);
396:   if (newOp != WEIGHTED_LAP) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "Default operator %d not setup correctly", WEIGHTED_LAP);
397:   return(0);
398: }

400: static struct _DiscretizationOps DOps = {PETSC_NULL/* DiscretizationSetup */,
401:                                          BoundaryDiscSetupOperators_Triangular_2D_Quadratic,
402:                                          PETSC_NULL/* DiscretizationSetFromOptions */,
403:                                          BoundaryDiscView_Triangular_2D_Quadratic,
404:                                          BoundaryDiscDestroy_Triangular_2D_Quadratic,
405:                                          BoundaryDiscEvaluateFunctionGalerkin_Triangular_2D_Quadratic,
406:                                          BoundaryDiscEvaluateOperatorGalerkin_Triangular_2D_Quadratic,
407:                                          BoundaryDiscEvaluateALEOperatorGalerkin_Triangular_2D_Quadratic,
408:                                          BoundaryDiscEvaluateNonlinearOperatorGalerkin_Triangular_2D_Quadratic,
409:                                          BoundaryDiscEvaluateNonlinearALEOperatorGalerkin_Triangular_2D_Quadratic,
410:                                          BoundaryDiscInterpolateField_Triangular_2D_Quadratic,
411:                                          BoundaryDiscInterpolateElementVec_Triangular_2D_Quadratic};

413: EXTERN_C_BEGIN
414: #undef   __FUNCT__
416: int BoundaryDiscCreate_Triangular_2D_Quadratic(Discretization disc) {
417:   int sarg;

421:   if (disc->comp <= 0) {
422:     SETERRQ(PETSC_ERR_ARG_WRONG, "Discretization must have at least 1 component. Call DiscretizationSetNumComponents() to set this.");
423:   }
424:   PetscMemcpy(disc->ops, &DOps, sizeof(struct _DiscretizationOps));
425:   disc->dim   = 2;
426:   disc->funcs = 3;
427:   disc->size  = disc->funcs*disc->comp;

429:   DiscretizationSetupDefaultOperators(disc);
430:   BoundaryDiscSetupQuadrature_Triangular_2D_Quadratic(disc);

432:   /* Storage */
433:   PetscMalloc(disc->comp * sizeof(PetscScalar),   &disc->funcVal);
434:   PetscMalloc(2          * sizeof(PetscScalar *), &disc->fieldVal);
435:   for(sarg = 0; sarg < 2; sarg++) {
436:     PetscMalloc(disc->comp*(disc->dim+1) * sizeof(PetscScalar), &disc->fieldVal[sarg]);
437:   }
438:   return(0);
439: }
440: EXTERN_C_END