Actual source code: linear1d.c

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

  5: /*
  6:    Defines piecewise linear function space on a one 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" /* Just for MAX_CORNERS */

 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:    where we recall that |J| is a constant for linear affine maps,
 18:    and the map of any segment to the standard element is linear.
 19:    The numbering of the nodes in the standard element is

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

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

 31: #undef  __FUNCT__
 33: static int DiscView_Triangular_1D_Linear_File(Discretization disc, PetscViewer viewer) {
 35:   PetscViewerASCIIPrintf(viewer, "Linear 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 DiscView_Triangular_1D_Linear(Discretization disc, PetscViewer viewer) {
 44:   PetscTruth isascii;
 45:   int        ierr;

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

 55: #undef  __FUNCT__
 57: static int DiscEvaluateFunctionGalerkin_Triangular_1D_Linear(Discretization disc, Mesh mesh, PointFunction f, PetscScalar alpha,
 58:                                                              int elem, PetscScalar *array, void *ctx)
 59: {
 60:   int          dim            = disc->dim;
 61:   int          comp           = disc->comp;           /* The number of components in this field */
 62:   int          funcs          = disc->funcs;          /* The number of shape functions per component */
 63:   PetscScalar *funcVal        = disc->funcVal;        /* Function value at a quadrature point */
 64:   int          numQuadPoints  = disc->numQuadPoints;  /* Number of points used for Gaussian quadrature */
 65:   double      *quadPoints     = disc->quadPoints;     /* Points in the standard element for Gaussian quadrature */
 66:   double      *quadWeights    = disc->quadWeights;    /* Weights in the standard element for Gaussian quadrature */
 67:   double      *quadShapeFuncs = disc->quadShapeFuncs; /* Shape function evaluated at quadrature points */
 68:   double       jac;                                   /* |J| for map to standard element */
 69:   double       x;                                     /* The integration point */
 70:   double       x11, x21;
 71:   int          rank, node0, node1;
 72:   int          i, j, k, p;
 73: #ifdef PETSC_USE_BOPT_g
 74:   Partition    part;
 75:   int          numOverlapElements;
 76:   PetscTruth   opt;
 77: #endif
 78:   int          ierr;

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

 83:   /* For dummy collective calls */
 84:   if (array == PETSC_NULL) {
 85:     for(p = 0; p < numQuadPoints; p++) {
 86:       (*f)(0, 0, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, ctx);
 87:     }
 88:     return(0);
 89:   }

 91: #ifdef PETSC_USE_BOPT_g
 92:   MeshGetPartition(mesh, &part);
 93:   PartitionGetNumOverlapElements(part, &numOverlapElements);
 94:   if ((elem < 0) || (elem >= numOverlapElements)) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Invalid element");
 95: #endif
 96:   /* Calculate the determinant of the inverse Jacobian of the map to the standard element
 97:      which must be a constant for linear elements */
 98:   MeshGetNodeFromElement(mesh, elem, 0, &node0);
 99:   MeshGetNodeFromElement(mesh, elem, 1, &node1);
100:   MeshGetNodeCoords(mesh, node0, &x11, PETSC_NULL, PETSC_NULL);
101:   MeshGetNodeCoords(mesh, node1, &x,   PETSC_NULL, PETSC_NULL);
102:   x21  = MeshPeriodicDiffX(mesh, x - x11);
103:   jac  = PetscAbsReal(x21);
104:   if (jac < 1.0e-14) {
105:     PetscPrintf(PETSC_COMM_SELF, "[%d]elem: %d x21: %gn", rank, elem, x21);
106:     SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
107:   }
108: #ifdef PETSC_USE_BOPT_g
109:   PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
110:   if (opt == PETSC_TRUE) {
111:     PetscPrintf(PETSC_COMM_SELF, "[%d]elem: %d x21: %g jac: %gn", rank, elem, x21, jac);
112:   }
113: #endif

115:   /* Calculate element vector entries by Gaussian quadrature */
116:   for(p = 0; p < numQuadPoints; p++) {
117:     x    = MeshPeriodicX(mesh, x21*quadPoints[p*dim] + x11);
118:     (*f)(1, comp, &x, PETSC_NULL, PETSC_NULL, funcVal, ctx);
119: #ifdef PETSC_USE_BOPT_g
120:     PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
121:     if (opt == PETSC_TRUE) {
122:       PetscPrintf(PETSC_COMM_SELF, "[%d]p:%d jac: %g", rank, p, jac);
123:       for(j = 0; j < comp; j++) PetscPrintf(PETSC_COMM_SELF, " func[%d]: %g", j, PetscRealPart(funcVal[j]));
124:       PetscPrintf(PETSC_COMM_SELF, "n");
125:   }
126: #endif

128:     for(i = 0, k = 0; i < funcs; i++) {
129:       for(j = 0; j < comp; j++, k++) {
130:         array[k] += alpha*funcVal[j]*quadShapeFuncs[p*funcs+i]*jac*quadWeights[p];
131: #ifdef PETSC_USE_BOPT_g
132:         PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
133:         if (opt == PETSC_TRUE) {
134:           PetscPrintf(PETSC_COMM_SELF, "[%d]  array[%d]: %gn", rank, k, PetscRealPart(array[k]));
135:         }
136: #endif
137:       }
138:     }
139:   }
140:   PetscLogFlops(1 + (1 + 5*funcs*comp) * numQuadPoints);
141:   return(0);
142: }

144: #undef  __FUNCT__
146: static int DiscEvaluateOperatorGalerkin_Triangular_1D_Linear(Discretization disc, Mesh mesh, int elemSize,
147:                                                              int rowStart, int colStart, int op, PetscScalar alpha,
148:                                                              int elem, PetscScalar *field, PetscScalar *mat, void *ctx)
149: {
150:   int              dim        = disc->dim;
151:   Operator         oper       = disc->operators[op]; /* The operator to discretize */
152:   Discretization   test       = oper->test;          /* The space of test functions */
153:   OperatorFunction opFunc     = oper->opFunc;        /* Integrals of operators which depend on J */
154:   PetscScalar     *precompInt = oper->precompInt;    /* Precomputed integrals of the operator on shape functions */
155:   int              rowSize    = test->size;          /* The number of shape functions per element */
156:   int              colSize    = disc->size;          /* The number of test  functions per element */
157:   double           x21;                              /* Coordinates of the element, with point 1 at the origin */
158:   double           jac;                              /* |J| for map to standard element */
159:   double           coords[MAX_CORNERS*2];            /* Coordinates of the element */
160:   double           x;
161:   int              rank, node0, node1;
162:   int              i, j, f;
163:   int              ierr;

166:   MPI_Comm_rank(disc->comm, &rank);
167: #ifdef PETSC_USE_BOPT_g
168:   /* Check for valid operator */
169:   if ((op < 0) || (op >= disc->numOps) || (!disc->operators[op])) SETERRQ(PETSC_ERR_ARG_WRONG, "Invalid operator");
170: #endif

172:   if (precompInt != PETSC_NULL) {
173:     /* Calculate the determinant of the inverse Jacobian of the map to the standard element
174:        which must be a constant for linear elements - 1/|x_{21} y_{31} - x_{31} y_{21}| */
175:     MeshGetNodeFromElement(mesh, elem, 0, &node0);
176:     MeshGetNodeFromElement(mesh, elem, 1, &node1);
177:     MeshGetNodeCoords(mesh, node0, &coords[0*dim+0], PETSC_NULL, PETSC_NULL);
178:     MeshGetNodeCoords(mesh, node1, &coords[1*dim+0], PETSC_NULL, PETSC_NULL);
179:     x21  = MeshPeriodicDiffX(mesh, coords[1*dim+0] - coords[0*dim+0]);
180:     jac  = PetscAbsReal(x21);
181:     if (jac < 1.0e-14) {
182:       PetscPrintf(PETSC_COMM_SELF, "[%d]x21: %g jac: %gn", rank, x21, jac);
183:       SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
184:     }

186:     /* Calculate element matrix entries which are all precomputed */
187:     for(i = 0; i < rowSize; i++) {
188:       for(j = 0; j < colSize; j++) {
189:         mat[(i+rowStart)*elemSize + j+colStart] += alpha*precompInt[i*colSize + j]*jac;
190:       }
191:     }
192:     PetscLogFlops(1 + 2*rowSize*colSize);
193:   } else {
194:     if (opFunc == PETSC_NULL) SETERRQ(PETSC_ERR_ARG_CORRUPT, "Invalid function");
195:     MeshGetNodeFromElement(mesh, elem, 0, &node0);
196:     MeshGetNodeCoords(mesh, node0, &coords[0*dim+0], PETSC_NULL, PETSC_NULL);
197:     for(f = 1; f < PetscMax(disc->funcs, test->funcs); f++) {
198:       MeshGetNodeFromElement(mesh, elem, f, &node1);
199:       MeshGetNodeCoords(mesh, node1, &x, PETSC_NULL, PETSC_NULL);
200:       coords[f*dim+0] = MeshPeriodicRelativeX(mesh, x, coords[0*dim+0]);
201:     }

203:     (*opFunc)(disc, test, rowSize, colSize, rowStart, colStart, elemSize, coords, alpha, field, mat, ctx);
204: 
205:   }
206:   return(0);
207: }

209: #undef  __FUNCT__
211: static int DiscEvaluateNonlinearOperatorGalerkin_Triangular_1D_Linear(Discretization disc, Mesh mesh, NonlinearOperator f,
212:                                                                       PetscScalar alpha, int elem, int numArgs, PetscScalar **field,
213:                                                                       PetscScalar *vec, void *ctx)
214: {
215:   int              dim        = disc->dim;
216:   int              comp       = disc->comp;      /* The number of components in this field */
217:   int              funcs      = disc->funcs;     /* The number of shape functions per component */
218:   PetscScalar     *funcVal    = disc->funcVal;   /* Function value at a quadrature point */
219:   PetscScalar    **fieldVal   = disc->fieldVal;  /* Field value and derivatives at a quadrature point */
220:   double           jac;                          /* |J| for map to standard element */
221:   double           invjac;                       /* |J^{-1}| for map from standard element */
222:   int              numQuadPoints;                /* Number of points used for Gaussian quadrature */
223:   double          *quadPoints;                   /* Points in the standard element for Gaussian quadrature */
224:   double          *quadWeights;                  /* Weights in the standard element for Gaussian quadrature */
225:   double          *quadShapeFuncs;               /* Shape function evaluated at quadrature points */
226:   double          *quadShapeFuncDers;            /* Shape function derivatives evaluated at quadrature points */
227:   double           x;                            /* The integration point */
228:   double           dxix;                         /* PartDer{xi}{x}  */
229:   PetscScalar      dfxi;                         /* PartDer{field}{xi}  */
230:   double           x11, x21;
231:   int              rank, node0, node1;
232:   int              i, j, k, p, func, arg;
233: #ifdef PETSC_USE_BOPT_g
234:   PetscTruth       opt;
235: #endif
236:   int              ierr;

239:   if (numArgs > 2) SETERRQ(PETSC_ERR_SUP, "Only configured to handle two nonlinear arguments");
240:   MPI_Comm_rank(disc->comm, &rank);
241:   numQuadPoints     = disc->numQuadPoints;
242:   quadPoints        = disc->quadPoints;
243:   quadWeights       = disc->quadWeights;
244:   quadShapeFuncs    = disc->quadShapeFuncs;
245:   quadShapeFuncDers = disc->quadShapeFuncDers;
246: 
247:   /* Calculate the determinant of the inverse Jacobian of the map to the standard element
248:      which must be a constant for linear elements */
249:   MeshGetNodeFromElement(mesh, elem, 0, &node0);
250:   MeshGetNodeFromElement(mesh, elem, 1, &node1);
251:   MeshGetNodeCoords(mesh, node0, &x11, PETSC_NULL, PETSC_NULL);
252:   MeshGetNodeCoords(mesh, node1, &x,   PETSC_NULL, PETSC_NULL);
253:   x21  = MeshPeriodicDiffX(mesh, x - x11);
254:   jac  = PetscAbsReal(x21);
255:   if (jac < 1.0e-14) {
256:     PetscPrintf(PETSC_COMM_SELF, "[%d]elem: %d x21: %gn", rank, elem, x21);
257:     SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
258:   }
259: #ifdef PETSC_USE_BOPT_g
260:   PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
261:   if (opt == PETSC_TRUE) {
262:     PetscPrintf(PETSC_COMM_SELF, "[%d]elem: %d x21: %g jac: %gn", rank, elem, x21, jac);
263:   }
264: #endif

266:   /* These are the elements of the inverse matrix */
267:   invjac = 1/jac;
268:   dxix   = invjac;

270:   /* Calculate element vector entries by Gaussian quadrature */
271:   for(p = 0; p < numQuadPoints; p++) {
272:     x = MeshPeriodicX(mesh, x21*quadPoints[p*dim] + x11);
273:     /* Can this be simplified? */
274:     for(arg = 0; arg < numArgs; arg++) {
275:       for(j = 0; j < comp*(dim+1); j++) fieldVal[arg][j] = 0.0;
276:       for(func = 0; func < funcs; func++) {
277:         for(j = 0; j < comp; j++) {
278:           fieldVal[arg][j*(dim+1)]   += field[arg][func*comp+j]*quadShapeFuncs[p*funcs+func];
279:           fieldVal[arg][j*(dim+1)+1] += field[arg][func*comp+j]*quadShapeFuncDers[p*funcs*dim+func*dim];
280:         }
281:       }
282:     }

284:     /* Convert the field derivatives to old coordinates */
285:     for(arg = 0; arg < numArgs; arg++) {
286:       for(j = 0; j < comp; j++) {
287:         dfxi                       = fieldVal[arg][j*(dim+1)+1];
288:         fieldVal[arg][j*(dim+1)+1] = dfxi*dxix;
289:       }
290:     }

292:     (*f)(1, comp, &x, PETSC_NULL, PETSC_NULL, numArgs, fieldVal, funcVal, ctx);
293: #ifdef PETSC_USE_BOPT_g
294:     PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
295:     if (opt == PETSC_TRUE) {
296:       PetscPrintf(PETSC_COMM_SELF, "[%d]p:%d jac: %g", rank, p, jac);
297:       for(j = 0; j < comp; j++)
298:         PetscPrintf(PETSC_COMM_SELF, " func[%d]: %g", j, PetscRealPart(funcVal[j]));
299:       PetscPrintf(PETSC_COMM_SELF, "n");
300:     }
301: #endif

303:     for(i = 0, k = 0; i < funcs; i++) {
304:       for(j = 0; j < comp; j++, k++) {
305:         vec[k] += alpha*funcVal[j]*quadShapeFuncs[p*funcs+i]*jac*quadWeights[p];
306: #ifdef PETSC_USE_BOPT_g
307:         PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
308:         if (opt == PETSC_TRUE) {
309:           PetscPrintf(PETSC_COMM_SELF, "[%d]  vec[%d]: %gn", rank, k, PetscRealPart(vec[k]));
310:         }
311: #endif
312:       }
313:     }
314:   }
315:   PetscLogFlops(2 + (1 + (4*numArgs + 5)*funcs*comp + numArgs*comp) * numQuadPoints);
316:   return(0);
317: }

319: #undef  __FUNCT__
321: static int DiscEvaluateALEOperatorGalerkin_Triangular_1D_Linear(Discretization disc, Mesh mesh, int elemSize,
322:                                                                 int rowStart, int colStart, int op, PetscScalar alpha,
323:                                                                 int elem, PetscScalar *field, PetscScalar *ALEfield, PetscScalar *mat,
324:                                                                 void *ctx)
325: {
326:   int                 dim     = disc->dim;
327:   Operator            oper    = disc->operators[op]; /* The operator to discretize */
328:   Discretization      test    = oper->test;          /* The space of test functions */
329:   ALEOperatorFunction opFunc  = oper->ALEOpFunc;     /* Integrals of operators which depend on J */
330:   int                 rowSize = test->size;          /* The number of shape functions per element */
331:   int                 colSize = disc->size;          /* The number of test  functions per element */
332:   double              coords[MAX_CORNERS*2];         /* Coordinates of the element */
333:   double              x;
334:   int                 node0, node1;
335:   int                 f;
336:   int                 ierr;

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

344:   if (opFunc == PETSC_NULL) SETERRQ(PETSC_ERR_ARG_CORRUPT, "Invalid function");
345:   MeshGetNodeFromElement(mesh, elem, 0, &node0);
346:   MeshGetNodeCoords(mesh, node0, &coords[0*dim+0], PETSC_NULL, PETSC_NULL);
347:   for(f = 1; f < PetscMax(disc->funcs, test->funcs); f++) {
348:     MeshGetNodeFromElement(mesh, elem, f, &node1);
349:     MeshGetNodeCoords(mesh, node1, &x, PETSC_NULL, PETSC_NULL);
350:     coords[f*dim+0] = MeshPeriodicRelativeX(mesh, x, coords[0*dim+0]);
351:   }

353:   (*opFunc)(disc, test, rowSize, colSize, rowStart, colStart, elemSize, coords, alpha, field, ALEfield, mat, ctx);
354: 
355:   return(0);
356: }

358: #undef  __FUNCT__
360: static int DiscEvaluateNonlinearALEOperatorGalerkin_Triangular_1D_Linear(Discretization disc, Mesh mesh, NonlinearOperator f,
361:                                                                          PetscScalar alpha, int elem, int numArgs, PetscScalar **field,
362:                                                                          PetscScalar *ALEfield, PetscScalar *vec, void *ctx)
363: {
364:   int              dim        = disc->dim;
365:   int              comp       = disc->comp;      /* The number of components in this field */
366:   int              funcs      = disc->funcs;     /* The number of shape functions per component */
367:   PetscScalar     *funcVal    = disc->funcVal;   /* Function value at a quadrature point */
368:   PetscScalar    **fieldVal   = disc->fieldVal;  /* Field value and derivatives at a quadrature point */
369:   double           jac;                          /* |J| for map to standard element */
370:   double           invjac;                       /* |J^{-1}| for map from standard element */
371:   int              numQuadPoints;                /* Number of points used for Gaussian quadrature */
372:   double          *quadPoints;                   /* Points in the standard element for Gaussian quadrature */
373:   double          *quadWeights;                  /* Weights in the standard element for Gaussian quadrature */
374:   double          *quadShapeFuncs;               /* Shape function evaluated at quadrature points */
375:   double          *quadShapeFuncDers;            /* Shape function derivatives evaluated at quadrature points */
376:   double           x;                            /* The integration point */
377:   double           dxix;                         /* PartDer{xi}{x}  */
378:   PetscScalar      dfxi;                         /* PartDer{field}{xi}  */
379:   double           x11, x21;
380:   int              rank, node0, node1;
381:   int              i, j, k, p, func, arg;
382: #ifdef PETSC_USE_BOPT_g
383:   PetscTruth       opt;
384: #endif
385:   int              ierr;

388:   if (numArgs > 2) SETERRQ(PETSC_ERR_SUP, "Only configured to handle two nonlinear arguments");
389:   MPI_Comm_rank(disc->comm, &rank);

391:   numQuadPoints     = disc->numQuadPoints;
392:   quadPoints        = disc->quadPoints;
393:   quadWeights       = disc->quadWeights;
394:   quadShapeFuncs    = disc->quadShapeFuncs;
395:   quadShapeFuncDers = disc->quadShapeFuncDers;
396: 
397:   /* Calculate the determinant of the inverse Jacobian of the map to the standard element
398:      which must be a constant for linear elements */
399:   MeshGetNodeFromElement(mesh, elem, 0, &node0);
400:   MeshGetNodeFromElement(mesh, elem, 1, &node1);
401:   MeshGetNodeCoords(mesh, node0, &x11, PETSC_NULL, PETSC_NULL);
402:   MeshGetNodeCoords(mesh, node1, &x,   PETSC_NULL, PETSC_NULL);
403:   x21  = MeshPeriodicDiffX(mesh, x - x11);
404:   jac  = PetscAbsReal(x21);
405:   if (jac < 1.0e-14) {
406:     PetscPrintf(PETSC_COMM_SELF, "[%d]elem: %d x21: %gn", rank, elem, x21);
407:     SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
408:   }
409: #ifdef PETSC_USE_BOPT_g
410:   PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
411:   if (opt == PETSC_TRUE) {
412:     PetscPrintf(PETSC_COMM_SELF, "[%d]elem: %d x21: %g jac: %gn", rank, elem, x21, jac);
413:   }
414: #endif

416:   /* These are the elements of the inverse matrix */
417:   invjac = 1/jac;
418:   dxix   = invjac;

420:   /* Calculate element vector entries by Gaussian quadrature */
421:   for(p = 0; p < numQuadPoints; p++) {
422:     x = MeshPeriodicX(mesh, x21*quadPoints[p*dim] + x11);
423:     /* Can this be simplified? */
424:     for(arg = 0; arg < numArgs; arg++) {
425:       for(j = 0; j < comp*(dim+1); j++) fieldVal[arg][j] = 0.0;
426:       for(func = 0; func < funcs; func++)
427:         for(j = 0; j < comp; j++) {
428:           fieldVal[arg][j*(dim+1)]   += (field[arg][func*comp+j] - ALEfield[func*comp+j])*quadShapeFuncs[p*funcs+func];
429:           fieldVal[arg][j*(dim+1)+1] += field[arg][func*comp+j]*quadShapeFuncDers[p*funcs*dim+func*dim];
430:         }
431:     }

433:     /* Convert the field derivatives to old coordinates */
434:     for(arg = 0; arg < numArgs; arg++) {
435:       for(j = 0; j < comp; j++) {
436:         dfxi                       = fieldVal[arg][j*(dim+1)+1];
437:         fieldVal[arg][j*(dim+1)+1] = dfxi*dxix;
438:       }
439:     }

441:     (*f)(1, comp, &x, PETSC_NULL, PETSC_NULL, numArgs, fieldVal, funcVal, ctx);
442: #ifdef PETSC_USE_BOPT_g
443:     PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
444:     if (opt == PETSC_TRUE) {
445:       PetscPrintf(PETSC_COMM_SELF, "[%d]p:%d jac: %g", rank, p, jac);
446:       for(j = 0; j < comp; j++)
447:         PetscPrintf(PETSC_COMM_SELF, " func[%d]: %g", j, PetscRealPart(funcVal[j]));
448:       PetscPrintf(PETSC_COMM_SELF, "n");
449:     }
450: #endif

452:     for(i = 0, k = 0; i < funcs; i++) {
453:       for(j = 0; j < comp; j++, k++) {
454:         vec[k] += alpha*funcVal[j]*quadShapeFuncs[p*funcs+i]*jac*quadWeights[p];
455: #ifdef PETSC_USE_BOPT_g
456:         PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
457:         if (opt == PETSC_TRUE) {
458:           PetscPrintf(PETSC_COMM_SELF, "[%d]  vec[%d]: %gn", rank, k, PetscRealPart(vec[k]));
459:         }
460: #endif
461:       }
462:     }
463:   }
464:   PetscLogFlops(2 + (1 + (5*numArgs + 5)*funcs*comp + numArgs*comp) * numQuadPoints);
465:   return(0);
466: }

468: #undef  __FUNCT__
470: int Laplacian_Triangular_1D_Linear(Discretization disc, Discretization test, int rowSize, int colSize,
471:                                    int globalRowStart, int globalColStart, int globalSize, double *coords,
472:                                    PetscScalar alpha, PetscScalar *field, PetscScalar *array, void *ctx)
473: {
474:   double      x21;          /* Coordinates of the element, with point 1 at the origin */
475:   double      jac, invjac;  /* |J| and |J^{-1}| for map to standard element */
476:   PetscScalar entry;
477:   int         comp;         /* Number of components */
478:   int         i;

481:   /* Calculate the determinant of the inverse Jacobian of the map to the standard element
482:      which must be a constant for linear elements - 1/|x_{21}| */
483:   x21 = coords[1] - coords[0];
484:   jac = PetscAbsReal(x21);
485: #ifdef PETSC_USE_BOPT_g
486:   if (jac < 1.0e-14) {
487:     PetscPrintf(PETSC_COMM_SELF, "x21: %g jac: %gn", x21, jac);
488:     SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
489:   }
490: #endif
491:   invjac = 1.0/jac;

493:   comp = rowSize/disc->funcs;
494:   /* alpha PartDer{phi}{x}^2 |J| = alpha PartDer{xi}{x}^2 |J| = alpha |J^{-1}|^2 |J| = alpha |J^{-1}| */
495:   entry = alpha*invjac;
496:   for(i = 0; i < comp; i++) {
497:     /* phi^1 phi^1 */
498:     array[(0*comp+i+globalRowStart)*globalSize + 0*comp+i+globalColStart] = -entry;
499:     /* phi^1 phi^2 */
500:     array[(0*comp+i+globalRowStart)*globalSize + 1*comp+i+globalColStart] =  entry;
501:     /* phi^2 phi^1 */
502:     array[(1*comp+i+globalRowStart)*globalSize + 0*comp+i+globalColStart] =  entry;
503:     /* phi^2 phi^2 */
504:     array[(1*comp+i+globalRowStart)*globalSize + 1*comp+i+globalColStart] = -entry;
505:   }
506:   PetscLogFlops(4);

508:   return(0);
509: }

511: #undef  __FUNCT__
513: int Weighted_Laplacian_Triangular_1D_Linear(Discretization disc, Discretization test, int rowSize, int colSize,
514:                                             int globalRowStart, int globalColStart, int globalSize, double *coords,
515:                                             PetscScalar alpha, PetscScalar *field, PetscScalar *array, void *ctx)
516: {
517:   double      x21;          /* Coordinates of the element, with point 1 at the origin */
518:   double      jac, invjac;  /* |J| and |J^{-1}| for map to standard element */
519:   PetscScalar entry;
520:   int         comp;         /* Number of components */
521:   int         i;

524:   /* Calculate the determinant of the inverse Jacobian of the map to the standard element
525:      which must be a constant for linear elements - 1/|x_{21}| */
526:   x21 = coords[1] - coords[0];
527:   jac = PetscAbsReal(x21);
528: #ifdef PETSC_USE_BOPT_g
529:   if (jac < 1.0e-14) {
530:     PetscPrintf(PETSC_COMM_SELF, "x21: %g jac: %gn", x21, jac);
531:     SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
532:   }
533: #endif
534:   invjac = 1.0/jac;

536:   comp = rowSize/disc->funcs;
537:   /* alpha PartDer{phi}{x}^2 = alpha PartDer{xi}{x}^2 = alpha |J^{-1}|^2 */
538:   entry = alpha*invjac*invjac;
539:   for(i = 0; i < comp; i++) {
540:     /* phi^1 phi^1 */
541:     array[(0*comp+i+globalRowStart)*globalSize + 0*comp+i+globalColStart] = -entry;
542:     /* phi^1 phi^2 */
543:     array[(0*comp+i+globalRowStart)*globalSize + 1*comp+i+globalColStart] =  entry;
544:     /* phi^2 phi^1 */
545:     array[(1*comp+i+globalRowStart)*globalSize + 0*comp+i+globalColStart] =  entry;
546:     /* phi^2 phi^2 */
547:     array[(1*comp+i+globalRowStart)*globalSize + 1*comp+i+globalColStart] = -entry;
548:   }
549:   PetscLogFlops(4);

551:   return(0);
552: }

554: #undef  __FUNCT__
556: int Gradient_Triangular_1D_Linear(Discretization disc, Discretization test, int rowSize, int colSize,
557:                                   int globalRowStart, int globalColStart, int globalSize, double *coords,
558:                                   PetscScalar alpha, PetscScalar *field, PetscScalar *array, void *ctx)
559: {
560:   /* We are using the convention that

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

565:      and

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

570:      where $d$ is the number of space dimensions. This agrees with the convention which allows
571:      $Delta matrix{u_1 cr u_2} = 0$ to denote a set of scalar equations. This also means that
572:      the dimension of the test function vector must be divisible by the number of space dimensions */
573:   int     numQuadPoints;     /* Number of points used for Gaussian quadrature */
574:   double *quadWeights;       /* Weights in the standard element for Gaussian quadrature */
575:   double *quadShapeFuncs;    /* Shape functions evaluated at quadrature points */
576:   double *quadShapeFuncDers; /* Shape function derivatives evaluated at quadrature points */
577:   double *quadTestFuncDers;  /* Test function derivatives evaluated at quadrature points */
578:   double  dxxi;              /* PartDer{x}{xi}  */
579:   double  dxix;              /* PartDer{xi}{x}  */
580:   double  dphix;             /* PartDer{phi_i}{x} times PartDer{phi_j}{x} */
581:   double  jac;               /* |J| for map to standard element */
582:   double  invjac;            /* |J^{-1}| for map from standard element */
583:   int     dim;               /* The problem dimension */
584:   int     comp;              /* The number of field components */
585:   int     tcomp;             /* The number of field components for the test field */
586:   int     funcs;             /* The number of shape functions */
587:   int     tfuncs;            /* The number of test functions */
588:   int     i, j, c, tc, f, p;

591:   /* Calculate element matrix entries by Gaussian quadrature --
592:            Since we integrate by parts here, the test and shape functions are switched */
593:   dim               = disc->dim;
594:   comp              = disc->comp;
595:   tcomp             = test->comp;
596:   funcs             = disc->funcs;
597:   tfuncs            = test->funcs;
598:   numQuadPoints     = disc->numQuadPoints;
599:   quadWeights       = disc->quadWeights;
600:   quadShapeFuncs    = disc->quadShapeFuncs;
601:   quadShapeFuncDers = disc->quadShapeFuncDers;
602:   quadTestFuncDers  = test->quadShapeFuncDers;
603:   for(p = 0; p < numQuadPoints; p++) {
604:     /* PartDer{x}{xi}(p)  = sum^{funcs}_{f=1} x_f PartDer{phi^f(p)}{xi} */
605:     dxxi = 0.0;
606:     if (tfuncs >= funcs) {
607:       for(f = 0; f < tfuncs; f++) {
608:         dxxi += coords[f*dim]*quadTestFuncDers[p*tfuncs*dim+f*dim];
609:       }
610:     } else {
611:       for(f = 0; f < funcs; f++) {
612:         dxxi += coords[f*dim]*quadShapeFuncDers[p*funcs*dim+f*dim];
613:       }
614:     }
615:     jac = PetscAbsReal(dxxi);
616: #ifdef PETSC_USE_BOPT_g
617:     if (jac < 1.0e-14) {
618:       PetscPrintf(PETSC_COMM_SELF, "p: %d x1: %g x2: %gn", p, coords[0], coords[1]);
619:       SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
620:     }
621: #endif
622:     /* These are the elements of the inverse matrix */
623:     invjac = 1.0/jac;
624:     dxix   = invjac;

626:     /* The rows are test functions */
627:     for(i = 0; i < tfuncs; i++) {
628:       /* We divide by the space dimension */
629:       for(tc = 0; tc < tcomp/dim; tc++) {
630:         /* The columns are shape functions */
631:         for(j = 0; j < funcs; j++) {
632:           dphix = quadTestFuncDers[p*tfuncs*dim+i*dim]*dxix;
633:           for(c = 0; c < comp; c++) {
634:             array[(i*tcomp+tc*dim+globalRowStart)*globalSize + j*comp+c+globalColStart] +=
635:               -alpha*dphix*quadShapeFuncs[p*funcs+j]*jac*quadWeights[p];
636:           }
637:         }
638:       }
639:     }
640:   }
641:   PetscLogFlops((2*tfuncs + 1 + 4*tfuncs*tcomp/dim*funcs*comp) * numQuadPoints);

643:   return(0);
644: }

646: #undef  __FUNCT__
648: int Divergence_Triangular_1D_Linear(Discretization disc, Discretization test, int rowSize, int colSize,
649:                                     int globalRowStart, int globalColStart, int globalSize, double *coords,
650:                                     PetscScalar alpha, PetscScalar *field, PetscScalar *array, void *ctx)
651: {
652:         /* We are using the convention that

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

657:                  and

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

662:                  where $d$ is the number of space dimensions. This agrees with the convention which allows
663:                  $Delta matrix{u_1 cr u_2} = 0$ to denote a set of scalar equations        This also requires that
664:      the dimension of a vector must be divisible by the space dimension in order to be acted upon by
665:      the divergence operator */
666:   int     numQuadPoints;     /* Number of points used for Gaussian quadrature */
667:   double *quadWeights;       /* Weights in the standard element for Gaussian quadrature */
668:   double *quadTestFuncs;     /* Test  functions evaluated at quadrature points */
669:   double *quadShapeFuncDers; /* Shape function derivatives evaluated at quadrature points */
670:   double  dxxi;              /* PartDer{x}{xi}  */
671:   double  dxix;              /* PartDer{xi}{x}  */
672:   double  dphix;             /* PartDer{phi_i}{x} times PartDer{phi_j}{x} */
673:   double  jac;               /* |J| for map to standard element */
674:   double  invjac;            /* |J^{-1}| for map from standard element */
675:   int     dim;               /* The problem dimension */
676:   int     comp;              /* The number of field components */
677:   int     tcomp;             /* The number of field components for the test field */
678:   int     funcs;             /* The number of shape functions */
679:   int     tfuncs;            /* The number of test functions */
680:   int     i, j, c, tc, f, p;

683:   /* Calculate element matrix entries by Gaussian quadrature */
684:   dim              = disc->dim;
685:   comp              = disc->comp;
686:   tcomp             = test->comp;
687:   funcs             = disc->funcs;
688:   tfuncs            = test->funcs;
689:   numQuadPoints     = disc->numQuadPoints;
690:   quadWeights       = disc->quadWeights;
691:   quadTestFuncs     = test->quadShapeFuncs;
692:   quadShapeFuncDers = disc->quadShapeFuncDers;
693:   for(p = 0; p < numQuadPoints; p++) {
694:     /* PartDer{x}{xi}(p)  = sum^{funcs}_{f=1} x_f PartDer{phi^f(p)}{xi} */
695:     dxxi = 0.0;
696:     for(f = 0; f < funcs; f++) {
697:       dxxi += coords[f*dim]*quadShapeFuncDers[p*funcs*dim+f*dim];
698:     }
699:     jac = PetscAbsReal(dxxi);
700: #ifdef PETSC_USE_BOPT_g
701:     if (jac < 1.0e-14) {
702:       PetscPrintf(PETSC_COMM_SELF, "p: %d x1: %g x2: %gn", p, coords[0], coords[1]);
703:       SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
704:     }
705: #endif
706:     /* These are the elements of the inverse matrix */
707:     invjac = 1.0/jac;
708:     dxix   = invjac;

710:     /* The rows are test functions */
711:     for(i = 0; i < tfuncs; i++) {
712:       for(tc = 0; tc < tcomp; tc++) {
713:         /* The columns are shape functions */
714:         for(j = 0; j < funcs; j++) {
715:           dphix = quadShapeFuncDers[p*funcs*dim+j*dim]*dxix;
716:           /* We divide by the number of space dimensions */
717:           for(c = 0; c < comp/dim; c++) {
718:             array[(i*tcomp+tc+globalRowStart)*globalSize + j*comp+c*dim+globalColStart] +=
719:               alpha*dphix*quadTestFuncs[p*tfuncs+i]*jac*quadWeights[p];
720:           }
721:         }
722:       }
723:     }
724:   }
725:   PetscLogFlops((2*funcs + 1 + 4*tfuncs*tcomp*funcs*comp/dim) * numQuadPoints);

727:   return(0);
728: }

730: #undef  __FUNCT__
732: int DiscInterpolateField_Triangular_1D_Linear(Discretization disc, Mesh oldMesh, int elem, double x, double y, double z,
733:                                               PetscScalar *oldFieldVal, PetscScalar *newFieldVal, InterpolationType type)
734: {
735:   double           x11, x22;    /* Coordinates of vertex 0 and 1 */
736:   double           xi;          /* Canonical coordinates of the interpolation point */
737:   double           dxix;        /* PartDer{xi}{x}  */
738:   double           dxxi;        /* PartDer{x}{xi}  */
739:   double           jac, invjac; /* The Jacobian determinant and its inverse */
740:   int              comp = disc->comp;
741:   int              rank, node0, node1;
742:   int              neighbor, corner, c;
743: #ifdef PETSC_USE_BOPT_g
744:   PetscTruth       opt;
745: #endif
746:   int              ierr;

749:   MPI_Comm_rank(disc->comm, &rank);
750:   /* No scheme in place for boundary elements */
751:   for(corner = 0; corner < 2; corner++) {
752:     MeshGetElementNeighbor(oldMesh, elem, corner, &neighbor);
753:     if (neighbor < 0) {
754:       type = INTERPOLATION_LOCAL;
755:       break;
756:     }
757:   }

759:   switch (type) {
760:   case INTERPOLATION_LOCAL:
761:     MeshGetNodeFromElement(oldMesh, elem, 0, &node0);
762:     MeshGetNodeFromElement(oldMesh, elem, 1, &node1);
763:     MeshGetNodeCoords(oldMesh, node0, &x11, PETSC_NULL, PETSC_NULL);
764:     MeshGetNodeCoords(oldMesh, node1, &x22, PETSC_NULL, PETSC_NULL);
765:     dxxi = MeshPeriodicDiffX(oldMesh, x22 - x11);
766:     jac  = PetscAbsReal(dxxi);
767:     if (jac < 1.0e-14) {
768:       PetscPrintf(PETSC_COMM_SELF, "[%d]elem: %d x21: %gn", rank, elem, dxxi);
769:       SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
770:     }
771: #ifdef PETSC_USE_BOPT_g
772:     PetscOptionsHasName(PETSC_NULL, "-trace_interpolation", &opt);
773:     if (opt == PETSC_TRUE) {
774:       PetscPrintf(PETSC_COMM_SELF, "[%d]elem: %d x21: %g jac: %gn", rank, elem, dxxi, jac);
775:     }
776: #endif

778:     /* These are the elements of the inverse matrix */
779:     invjac = 1/jac;
780:     dxix   = invjac;
781:     xi     = dxix*MeshPeriodicDiffX(oldMesh, x - x11);
782:     for(c = 0 ; c < comp; c++) {
783:       newFieldVal[c] = oldFieldVal[0*comp+c]*(1.0 - xi) + oldFieldVal[1*comp+c]*xi;
784:     }
785:     PetscLogFlops(4+3*comp);
786:     break;
787:   default:
788:     SETERRQ1(PETSC_ERR_ARG_WRONG, "Unknown interpolation type %d", type);
789:   }
790: 
791:   return(0);
792: }

794: #undef  __FUNCT__
796: int DiscInterpolateElementVec_Triangular_1D_Linear(Discretization disc, ElementVec vec, Discretization newDisc, ElementVec newVec) {
797:   int          funcs = disc->funcs;
798:   int          comp  = disc->comp;
799:   int          size  = disc->size;
800:   PetscScalar *array, *newArray;
801:   PetscTruth   islin, isquad;
802:   int          f, c;
803:   int          ierr;

806:   ElementVecGetArray(vec,    &array);
807:   ElementVecGetArray(newVec, &newArray);
808:   PetscTypeCompare((PetscObject) newDisc, DISCRETIZATION_TRIANGULAR_1D_LINEAR,    &islin);
809:   PetscTypeCompare((PetscObject) newDisc, DISCRETIZATION_TRIANGULAR_1D_QUADRATIC, &isquad);
810:   if (islin == PETSC_TRUE) {
811:     PetscMemcpy(newArray, array, size * sizeof(PetscScalar));
812:   } else if (isquad == PETSC_TRUE) {
813:     for(f = 0; f < newDisc->funcs; f++) {
814:       for(c = 0; c < comp; c++) {
815:         if (f < funcs) {
816:           newArray[f*comp+c] = array[f*comp+c];
817:         } else {
818:           newArray[f*comp+c] = 0.5*(array[((f+1)%funcs)*comp+c] + array[((f+2)%funcs)*comp+c]);
819:         }
820:       }
821:     }
822:   } else {
823:     SETERRQ(PETSC_ERR_SUP, "Discretization not supported");
824:   }
825:   ElementVecRestoreArray(vec,    &array);
826:   ElementVecRestoreArray(newVec, &newArray);
827:   return(0);
828: }

830: #undef  __FUNCT__
832: /*
833:   DiscSetupQuadrature_Triangular_1D_Linear - Setup Gaussian quadrature with a 7 point integration rule

835:   Input Parameter:
836: . disc - The Discretization
837: */
838: int DiscSetupQuadrature_Triangular_1D_Linear(Discretization disc) {
839:   int dim   = disc->dim;
840:   int funcs = disc->funcs;
841:   int p;

845:   disc->numQuadPoints = 7;
846:   PetscMalloc(disc->numQuadPoints*dim       * sizeof(double), &disc->quadPoints);
847:   PetscMalloc(disc->numQuadPoints           * sizeof(double), &disc->quadWeights);
848:   PetscMalloc(disc->numQuadPoints*funcs     * sizeof(double), &disc->quadShapeFuncs);
849:   PetscMalloc(disc->numQuadPoints*funcs*dim * sizeof(double), &disc->quadShapeFuncDers);
850:   PetscLogObjectMemory(disc, (disc->numQuadPoints*(funcs*(dim+1) + dim+1)) * sizeof(double));
851:   disc->quadPoints[0]  = 0.0254460438286207377369052;
852:   disc->quadWeights[0] = 0.0647424830844348466353057;
853:   disc->quadPoints[1]  = 0.1292344072003027800680676;
854:   disc->quadWeights[1] = 0.1398526957446383339507339;
855:   disc->quadPoints[2]  = 0.29707742431130141654669679;
856:   disc->quadWeights[2] = 0.1909150252525594724751849;
857:   disc->quadPoints[3]  = 0.5000000000000000000000000;
858:   disc->quadWeights[3] = 0.2089795918367346938775510;
859:   disc->quadPoints[4]  = 0.70292257568869858345330321;
860:   disc->quadWeights[4] = disc->quadWeights[2];
861:   disc->quadPoints[5]  = 0.8707655927996972199319324;
862:   disc->quadWeights[5] = disc->quadWeights[1];
863:   disc->quadPoints[6]  = 0.9745539561713792622630948;
864:   disc->quadWeights[6] = disc->quadWeights[0];
865:   for(p = 0; p < disc->numQuadPoints; p++) {
866:     /* phi^0: 1 - xi */
867:     disc->quadShapeFuncs[p*funcs]              =  1.0 - disc->quadPoints[p*dim];
868:     disc->quadShapeFuncDers[p*funcs*dim+0*dim] = -1.0;
869:     /* phi^1: xi */
870:     disc->quadShapeFuncs[p*funcs+1]            =  disc->quadPoints[p*dim];
871:     disc->quadShapeFuncDers[p*funcs*dim+1*dim] =  1.0;
872:   }
873:   return(0);
874: }

876: #undef  __FUNCT__
878: /*
879:   DiscSetupOperators_Triangular_1D_Linear - Setup the default operators

881:   Input Parameter:
882: . disc - The Discretization
883: */
884: int DiscSetupOperators_Triangular_1D_Linear(Discretization disc) {
885:   int          comp = disc->comp;
886:   int          size = disc->size;
887:   PetscScalar *precompInt;
888:   int          newOp;
889:   int          c, i, j;
890:   int          ierr;

893:   /* The Identity operator I -- the matrix is symmetric */
894:   PetscMalloc(size*size * sizeof(PetscScalar), &precompInt);
895:   PetscLogObjectMemory(disc, size*size * sizeof(PetscScalar));
896:   PetscMemzero(precompInt, size*size * sizeof(PetscScalar));
897:   for(c = 0; c < comp; c++) {
898:     precompInt[(0*comp+c)*size + 0*comp+c] = 1.0/3.0;
899:     precompInt[(0*comp+c)*size + 1*comp+c] = 1.0/6.0;
900:     precompInt[(1*comp+c)*size + 1*comp+c] = 1.0/3.0;
901:   }
902:   for(i = 0; i < size; i++) {
903:     for(j = 0; j < i; j++) {
904:       precompInt[i*size + j] = precompInt[j*size + i];
905:     }
906:   }
907:   DiscretizationRegisterPrecomputedOperator(disc, precompInt, &newOp);
908:   if (newOp != IDENTITY) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "Default operator %d not setup correctly", IDENTITY);
909:   /* The Laplacian operator Delta -- the matrix is symmetric */
910:   DiscretizationRegisterOperator(disc, Laplacian_Triangular_1D_Linear, &newOp);
911:   if (newOp != LAPLACIAN) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "Default operator %d not setup correctly", LAPLACIAN);
912:   /* The Gradient operator nabla -- the matrix is rectangular */
913:   DiscretizationRegisterOperator(disc, Gradient_Triangular_1D_Linear, &newOp);
914:   if (newOp != GRADIENT) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "Default operator %d not setup correctly", GRADIENT);
915:   /* The Divergence operator nablacdot -- the matrix is rectangular */
916:   DiscretizationRegisterOperator(disc, Divergence_Triangular_1D_Linear, &newOp);
917:   if (newOp != DIVERGENCE) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "Default operator %d not setup correctly", DIVERGENCE);
918:   /* The weighted Laplacian operator -- the matrix is symmetric */
919:   DiscretizationRegisterOperator(disc, Weighted_Laplacian_Triangular_1D_Linear, &newOp);
920:   if (newOp != WEIGHTED_LAP) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "Default operator %d not setup correctly", WEIGHTED_LAP);
921:   return(0);
922: }

924: static struct _DiscretizationOps DOps = {PETSC_NULL/* DiscretizationSetup */,
925:                                          DiscSetupOperators_Triangular_1D_Linear,
926:                                          PETSC_NULL/* DiscretizationSetFromOptions */,
927:                                          DiscView_Triangular_1D_Linear,
928:                                          DiscDestroy_Triangular_1D_Linear,
929:                                          DiscEvaluateFunctionGalerkin_Triangular_1D_Linear,
930:                                          DiscEvaluateOperatorGalerkin_Triangular_1D_Linear,
931:                                          DiscEvaluateALEOperatorGalerkin_Triangular_1D_Linear,
932:                                          DiscEvaluateNonlinearOperatorGalerkin_Triangular_1D_Linear,
933:                                          DiscEvaluateNonlinearALEOperatorGalerkin_Triangular_1D_Linear,
934:                                          DiscInterpolateField_Triangular_1D_Linear,
935:                                          DiscInterpolateElementVec_Triangular_1D_Linear};

937: EXTERN_C_BEGIN
938: #undef  __FUNCT__
940: int DiscCreate_Triangular_1D_Linear(Discretization disc) {
941:   int arg;

945:   if (disc->comp <= 0) {
946:     SETERRQ(PETSC_ERR_ARG_WRONG, "Discretization must have at least 1 component. Call DiscretizationSetNumComponents() to set this.");
947:   }
948:   PetscMemcpy(disc->ops, &DOps, sizeof(struct _DiscretizationOps));
949:   disc->dim   = 1;
950:   disc->funcs = 2;
951:   disc->size  = disc->funcs*disc->comp;

953:   DiscretizationSetupDefaultOperators(disc);
954:   DiscSetupQuadrature_Triangular_1D_Linear(disc);

956:   /* Storage */
957:   PetscMalloc(disc->comp * sizeof(PetscScalar),   &disc->funcVal);
958:   PetscMalloc(2          * sizeof(PetscScalar *), &disc->fieldVal);
959:   for(arg = 0; arg < 2; arg++) {
960:     PetscMalloc(disc->comp*(disc->dim+1) * sizeof(PetscScalar), &disc->fieldVal[arg]);
961:   }
962:   return(0);
963: }
964: EXTERN_C_END