Actual source code: constant1d.c

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

  5: /*
  6:    Defines piecewise constant 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" /* 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 node to the standard element is linear.
 19:    The numbering of the nodes in the standard element is

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

 24: static int DiscDestroy_Triangular_1D_Constant(Discretization disc) {
 26:   return(0);
 27: }

 29: static int DiscView_Triangular_1D_Constant_File(Discretization disc, PetscViewer viewer) {
 31:   PetscViewerASCIIPrintf(viewer, "Constant discretizationn");
 32:   PetscViewerASCIIPrintf(viewer, "    %d shape functions per componentn", disc->funcs);
 33:   PetscViewerASCIIPrintf(viewer, "    %d registered operatorsn", disc->numOps);
 34:   return(0);
 35: }

 37: static int DiscView_Triangular_1D_Constant(Discretization disc, PetscViewer viewer) {
 38:   PetscTruth isascii;
 39:   int        ierr;

 42:   PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_ASCII, &isascii);
 43:   if (isascii == PETSC_TRUE) {
 44:     DiscView_Triangular_1D_Constant_File(disc, viewer);
 45:   }
 46:   return(0);
 47: }

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

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

 75:   /* For dummy collective calls */
 76:   if (array == PETSC_NULL) {
 77:     for(p = 0; p < numQuadPoints; p++) {
 78:       (*f)(0, 0, PETSC_NULL, PETSC_NULL, PETSC_NULL, PETSC_NULL, ctx);
 79:     }
 80:     return(0);
 81:   }

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

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

120:     for(i = 0, k = 0; i < funcs; i++) {
121:       for(j = 0; j < comp; j++, k++) {
122:         array[k] += alpha*funcVal[j]*quadShapeFuncs[p*funcs+i]*jac*quadWeights[p];
123: #ifdef PETSC_USE_BOPT_g
124:         PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
125:         if (opt == PETSC_TRUE) {
126:           PetscPrintf(PETSC_COMM_SELF, "[%d]  array[%d]: %gn", rank, k, PetscRealPart(array[k]));
127:         }
128: #endif
129:       }
130:     }
131:   }
132:   PetscLogFlops(1 + (1 + 5*funcs*comp) * numQuadPoints);
133:   return(0);
134: }

136: static int DiscEvaluateOperatorGalerkin_Triangular_1D_Constant(Discretization disc, Mesh mesh, int elemSize,
137:                                                                int rowStart, int colStart, int op, PetscScalar alpha,
138:                                                                int elem, PetscScalar *field, PetscScalar *mat, void *ctx)
139: {
140:   int              dim        = disc->dim;
141:   Operator         oper       = disc->operators[op]; /* The operator to discretize */
142:   Discretization   test       = oper->test;          /* The space of test functions */
143:   OperatorFunction opFunc     = oper->opFunc;        /* Integrals of operators which depend on J */
144:   PetscScalar     *precompInt = oper->precompInt;    /* Precomputed integrals of the operator on shape functions */
145:   int              rowSize    = test->size;          /* The number of shape functions per element */
146:   int              colSize    = disc->size;          /* The number of test  functions per element */
147:   double           x21;                              /* Coordinates of the element, with point 1 at the origin */
148:   double           jac;                              /* |J| for map to standard element */
149:   double           coords[MAX_CORNERS*2];            /* Coordinates of the element */
150:   double           x;
151:   int              rank, node0, node1;
152:   int              i, j, f;
153:   int              ierr;

156:   MPI_Comm_rank(disc->comm, &rank);
157: #ifdef PETSC_USE_BOPT_g
158:   /* Check for valid operator */
159:   if ((op < 0) || (op >= disc->numOps) || (!disc->operators[op])) SETERRQ(PETSC_ERR_ARG_WRONG, "Invalid operator");
160: #endif

162:   if (precompInt != PETSC_NULL) {
163:     /* Calculate the determinant of the inverse Jacobian of the map to the standard element
164:        which must be a constant for linear elements - 1/|x_{21} y_{31} - x_{31} y_{21}| */
165:     MeshGetNodeFromElement(mesh, elem, 0, &node0);
166:     MeshGetNodeFromElement(mesh, elem, 1, &node1);
167:     MeshGetNodeCoords(mesh, node0, &coords[0*dim+0], PETSC_NULL, PETSC_NULL);
168:     MeshGetNodeCoords(mesh, node1, &coords[1*dim+0], PETSC_NULL, PETSC_NULL);
169:     x21  = MeshPeriodicDiffX(mesh, coords[1*dim+0] - coords[0*dim+0]);
170:     jac  = PetscAbsReal(x21);
171:     if (jac < 1.0e-14) {
172:       PetscPrintf(PETSC_COMM_SELF, "[%d]x21: %g jac: %gn", rank, x21, jac);
173:       SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
174:     }

176:     /* Calculate element matrix entries which are all precomputed */
177:     for(i = 0; i < rowSize; i++) {
178:       for(j = 0; j < colSize; j++) {
179:         mat[(i+rowStart)*elemSize + j+colStart] += alpha*precompInt[i*colSize + j]*jac;
180:       }
181:     }
182:     PetscLogFlops(1 + 2*rowSize*colSize);
183:   } else {
184:     if (opFunc == PETSC_NULL) SETERRQ(PETSC_ERR_ARG_CORRUPT, "Invalid function");
185:     MeshGetNodeFromElement(mesh, elem, 0, &node0);
186:     MeshGetNodeCoords(mesh, node0, &coords[0*dim+0], PETSC_NULL, PETSC_NULL);
187:     for(f = 1; f < PetscMax(disc->funcs, test->funcs); f++) {
188:       MeshGetNodeFromElement(mesh, elem, f, &node1);
189:       MeshGetNodeCoords(mesh, node1, &x, PETSC_NULL, PETSC_NULL);
190:       coords[f*dim+0] = MeshPeriodicRelativeX(mesh, x, coords[0*dim+0]);
191:     }

193:     (*opFunc)(disc, test, rowSize, colSize, rowStart, colStart, elemSize, coords, alpha, field, mat, ctx);
194: 
195:   }
196:   return(0);
197: }

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

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

254:   /* These are the elements of the inverse matrix */
255:   invjac = 1/jac;
256:   dxix   = invjac;

258:   /* Calculate element vector entries by Gaussian quadrature */
259:   for(p = 0; p < numQuadPoints; p++) {
260:     x = MeshPeriodicX(mesh, x21*quadPoints[p*dim] + x11);
261:     /* Can this be simplified? */
262:     for(arg = 0; arg < numArgs; arg++) {
263:       for(j = 0; j < comp*3; j++) fieldVal[arg][j] = 0.0;
264:       for(func = 0; func < funcs; func++) {
265:         for(j = 0; j < comp; j++) {
266:           fieldVal[arg][j*(dim+1)]   += field[arg][func*comp+j]*quadShapeFuncs[p*funcs+func];
267:           fieldVal[arg][j*(dim+1)+1] += field[arg][func*comp+j]*quadShapeFuncDers[p*funcs*dim+func*dim];
268:         }
269:       }
270:     }

272:     /* Convert the field derivatives to old coordinates */
273:     for(arg = 0; arg < numArgs; arg++) {
274:       for(j = 0; j < comp; j++) {
275:         dfxi                       = fieldVal[arg][j*(dim+1)+1];
276:         fieldVal[arg][j*(dim+1)+1] = dfxi*dxix;
277:       }
278:     }

280:     (*f)(1, comp, &x, PETSC_NULL, PETSC_NULL, numArgs, fieldVal, funcVal, ctx);
281: #ifdef PETSC_USE_BOPT_g
282:     PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
283:     if (opt == PETSC_TRUE) {
284:       PetscPrintf(PETSC_COMM_SELF, "[%d]p:%d jac: %g", rank, p, jac);
285:       for(j = 0; j < comp; j++)
286:         PetscPrintf(PETSC_COMM_SELF, " func[%d]: %g", j, PetscRealPart(funcVal[j]));
287:       PetscPrintf(PETSC_COMM_SELF, "n");
288:     }
289: #endif

291:     for(i = 0, k = 0; i < funcs; i++) {
292:       for(j = 0; j < comp; j++, k++) {
293:         vec[k] += alpha*funcVal[j]*quadShapeFuncs[p*funcs+i]*jac*quadWeights[p];
294: #ifdef PETSC_USE_BOPT_g
295:         PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
296:         if (opt == PETSC_TRUE) {
297:           PetscPrintf(PETSC_COMM_SELF, "[%d]  vec[%d]: %gn", rank, k, PetscRealPart(vec[k]));
298:         }
299: #endif
300:       }
301:     }
302:   }
303:   PetscLogFlops(2 + (1 + (4*numArgs + 5)*funcs*comp + numArgs*comp) * numQuadPoints);
304:   return(0);
305: }

307: static int DiscEvaluateALEOperatorGalerkin_Triangular_1D_Constant(Discretization disc, Mesh mesh, int elemSize,
308:                                                                   int rowStart, int colStart, int op, PetscScalar alpha,
309:                                                                   int elem, PetscScalar *field, PetscScalar *ALEfield, PetscScalar *mat,
310:                                                                   void *ctx)
311: {
312:   int                 dim     = disc->dim;
313:   Operator            oper    = disc->operators[op]; /* The operator to discretize */
314:   Discretization      test    = oper->test;          /* The space of test functions */
315:   ALEOperatorFunction opFunc  = oper->ALEOpFunc;     /* Integrals of operators which depend on J */
316:   int                 rowSize = test->size;          /* The number of shape functions per element */
317:   int                 colSize = disc->size;          /* The number of test  functions per element */
318:   double              coords[MAX_CORNERS*2];         /* Coordinates of the element */
319:   double              x;
320:   int                 node0, node1;
321:   int                 f;
322:   int                 ierr;

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

330:   if (opFunc == PETSC_NULL) SETERRQ(PETSC_ERR_ARG_CORRUPT, "Invalid function");
331:   MeshGetNodeFromElement(mesh, elem, 0, &node0);
332:   MeshGetNodeCoords(mesh, node0, &coords[0*dim+0], PETSC_NULL, PETSC_NULL);
333:   for(f = 1; f < PetscMax(disc->funcs, test->funcs); f++) {
334:     MeshGetNodeFromElement(mesh, elem, f, &node1);
335:     MeshGetNodeCoords(mesh, node1, &x, PETSC_NULL, PETSC_NULL);
336:     coords[f*dim+0] = MeshPeriodicRelativeX(mesh, x, coords[0*dim+0]);
337:   }

339:   (*opFunc)(disc, test, rowSize, colSize, rowStart, colStart, elemSize, coords, alpha, field, ALEfield, mat, ctx);
340: 
341:   return(0);
342: }

344: static int DiscEvaluateNonlinearALEOperatorGalerkin_Triangular_1D_Constant(Discretization disc, Mesh mesh, NonlinearOperator f,
345:                                                                            PetscScalar alpha, int elem, int numArgs, PetscScalar **field,
346:                                                                            PetscScalar *ALEfield, PetscScalar *vec, void *ctx)
347: {
348:   int              dim        = disc->dim;
349:   int              comp       = disc->comp;      /* The number of components in this field */
350:   int              funcs      = disc->funcs;     /* The number of shape functions per component */
351:   PetscScalar     *funcVal    = disc->funcVal;   /* Function value at a quadrature point */
352:   PetscScalar    **fieldVal   = disc->fieldVal;  /* Field value and derivatives at a quadrature point */
353:   double           jac;                          /* |J| for map to standard element */
354:   double           invjac;                       /* |J^{-1}| for map from standard element */
355:   int              numQuadPoints;                /* Number of points used for Gaussian quadrature */
356:   double          *quadPoints;                   /* Points in the standard element for Gaussian quadrature */
357:   double          *quadWeights;                  /* Weights in the standard element for Gaussian quadrature */
358:   double          *quadShapeFuncs;               /* Shape function evaluated at quadrature points */
359:   double          *quadShapeFuncDers;            /* Shape function derivatives evaluated at quadrature points */
360:   double           x;                            /* The integration point */
361:   double           dxix;                         /* PartDer{xi}{x}  */
362:   PetscScalar      dfxi;                         /* PartDer{field}{xi}  */
363:   double           x11, x21;
364:   int              rank, node0, node1;
365:   int              i, j, k, p, func, arg;
366: #ifdef PETSC_USE_BOPT_g
367:   PetscTruth       opt;
368: #endif
369:   int              ierr;

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

375:   numQuadPoints     = disc->numQuadPoints;
376:   quadPoints        = disc->quadPoints;
377:   quadWeights       = disc->quadWeights;
378:   quadShapeFuncs    = disc->quadShapeFuncs;
379:   quadShapeFuncDers = disc->quadShapeFuncDers;
380: 
381:   /* Calculate the determinant of the inverse Jacobian of the map to the standard element
382:      which must be a constant for linear elements */
383:   MeshGetNodeFromElement(mesh, elem, 0, &node0);
384:   MeshGetNodeFromElement(mesh, elem, 1, &node1);
385:   MeshGetNodeCoords(mesh, node0, &x11, PETSC_NULL, PETSC_NULL);
386:   MeshGetNodeCoords(mesh, node1, &x,   PETSC_NULL, PETSC_NULL);
387:   x21  = MeshPeriodicDiffX(mesh, x - x11);
388:   jac  = PetscAbsReal(x21);
389:   if (jac < 1.0e-14) {
390:     PetscPrintf(PETSC_COMM_SELF, "[%d]elem: %d x21: %gn", rank, elem, x21);
391:     SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
392:   }
393: #ifdef PETSC_USE_BOPT_g
394:   PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
395:   if (opt == PETSC_TRUE) {
396:     PetscPrintf(PETSC_COMM_SELF, "[%d]elem: %d x21: %g jac: %gn", rank, elem, x21, jac);
397:   }
398: #endif

400:   /* These are the elements of the inverse matrix */
401:   invjac = 1/jac;
402:   dxix   = invjac;

404:   /* Calculate element vector entries by Gaussian quadrature */
405:   for(p = 0; p < numQuadPoints; p++) {
406:     x = MeshPeriodicX(mesh, x21*quadPoints[p*dim] + x11);
407:     /* Can this be simplified? */
408:     for(arg = 0; arg < numArgs; arg++) {
409:       for(j = 0; j < comp*(dim+1); j++) fieldVal[arg][j] = 0.0;
410:       for(func = 0; func < funcs; func++)
411:         for(j = 0; j < comp; j++) {
412:           fieldVal[arg][j*(dim+1)]   += (field[arg][func*comp+j] - ALEfield[func*comp+j])*quadShapeFuncs[p*funcs+func];
413:           fieldVal[arg][j*(dim+1)+1] += field[arg][func*comp+j]*quadShapeFuncDers[p*funcs*dim+func*dim];
414:         }
415:     }

417:     /* Convert the field derivatives to old coordinates */
418:     for(arg = 0; arg < numArgs; arg++) {
419:       for(j = 0; j < comp; j++) {
420:         dfxi                       = fieldVal[arg][j*(dim+1)+1];
421:         fieldVal[arg][j*(dim+1)+1] = dfxi*dxix;
422:       }
423:     }

425:     (*f)(1, comp, &x, PETSC_NULL, PETSC_NULL, numArgs, fieldVal, funcVal, ctx);
426: #ifdef PETSC_USE_BOPT_g
427:     PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
428:     if (opt == PETSC_TRUE) {
429:       PetscPrintf(PETSC_COMM_SELF, "[%d]p:%d jac: %g", rank, p, jac);
430:       for(j = 0; j < comp; j++)
431:         PetscPrintf(PETSC_COMM_SELF, " func[%d]: %g", j, PetscRealPart(funcVal[j]));
432:       PetscPrintf(PETSC_COMM_SELF, "n");
433:     }
434: #endif

436:     for(i = 0, k = 0; i < funcs; i++) {
437:       for(j = 0; j < comp; j++, k++) {
438:         vec[k] += alpha*funcVal[j]*quadShapeFuncs[p*funcs+i]*jac*quadWeights[p];
439: #ifdef PETSC_USE_BOPT_g
440:         PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
441:         if (opt == PETSC_TRUE) {
442:           PetscPrintf(PETSC_COMM_SELF, "[%d]  vec[%d]: %gn", rank, k, PetscRealPart(vec[k]));
443:         }
444: #endif
445:       }
446:     }
447:   }
448:   PetscLogFlops(2 + (1 + (5*numArgs + 5)*funcs*comp + numArgs*comp) * numQuadPoints);
449:   return(0);
450: }

452: int Laplacian_Triangular_1D_Constant(Discretization disc, Discretization test, int rowSize, int colSize,
453:                                    int globalRowStart, int globalColStart, int globalSize, double *coords,
454:                                    PetscScalar alpha, PetscScalar *field, PetscScalar *array, void *ctx)
455: {
456:   double      x21;          /* Coordinates of the element, with point 1 at the origin */
457:   double      jac, invjac;  /* |J| and |J^{-1}| for map to standard element */
458:   PetscScalar entry;
459:   int         comp;         /* Number of components */
460:   int         i;

463:   /* Calculate the determinant of the inverse Jacobian of the map to the standard element
464:      which must be a constant for linear elements - 1/|x_{21}| */
465:   x21 = coords[1] - coords[0];
466:   jac = PetscAbsReal(x21);
467: #ifdef PETSC_USE_BOPT_g
468:   if (jac < 1.0e-14) {
469:     PetscPrintf(PETSC_COMM_SELF, "x21: %g jac: %gn", x21, jac);
470:     SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
471:   }
472: #endif
473:   invjac = 1.0/jac;

475:   comp = rowSize/disc->funcs;
476:   /* alpha PartDer{phi}{x}^2 |J| = alpha PartDer{xi}{x}^2 |J| = alpha |J^{-1}|^2 |J| = alpha |J^{-1}| */
477:   entry = alpha*invjac;
478:   for(i = 0; i < comp; i++) {
479:     /* phi^1 phi^1 */
480:     array[(0*comp+i+globalRowStart)*globalSize + 0*comp+i+globalColStart] = entry;
481:   }
482:   PetscLogFlops(4);

484:   return(0);
485: }

487: int Weighted_Laplacian_Triangular_1D_Constant(Discretization disc, Discretization test, int rowSize, int colSize,
488:                                             int globalRowStart, int globalColStart, int globalSize, double *coords,
489:                                             PetscScalar alpha, PetscScalar *field, PetscScalar *array, void *ctx)
490: {
491:   double      x21;          /* Coordinates of the element, with point 1 at the origin */
492:   double      jac, invjac;  /* |J| and |J^{-1}| for map to standard element */
493:   PetscScalar entry;
494:   int         comp;         /* Number of components */
495:   int         i;

498:   /* Calculate the determinant of the inverse Jacobian of the map to the standard element
499:      which must be a constant for linear elements - 1/|x_{21}| */
500:   x21 = coords[1] - coords[0];
501:   jac = PetscAbsReal(x21);
502: #ifdef PETSC_USE_BOPT_g
503:   if (jac < 1.0e-14) {
504:     PetscPrintf(PETSC_COMM_SELF, "x21: %g jac: %gn", x21, jac);
505:     SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
506:   }
507: #endif
508:   invjac = 1.0/jac;

510:   comp = rowSize/disc->funcs;
511:   /* alpha PartDer{phi}{x}^2 = alpha PartDer{xi}{x}^2 = alpha |J^{-1}|^2 */
512:   entry = alpha*invjac*invjac;
513:   for(i = 0; i < comp; i++) {
514:     /* phi^1 phi^1 */
515:     array[(0*comp+i+globalRowStart)*globalSize + 0*comp+i+globalColStart] = entry;
516:   }
517:   PetscLogFlops(4);

519:   return(0);
520: }

522: int Gradient_Triangular_1D_Constant(Discretization disc, Discretization test, int rowSize, int colSize,
523:                                   int globalRowStart, int globalColStart, int globalSize, double *coords,
524:                                   PetscScalar alpha, PetscScalar *field, PetscScalar *array, void *ctx)
525: {
526:   /* We are using the convention that

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

531:      and

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

536:      where $d$ is the number of space dimensions. This agrees with the convention which allows
537:      $Delta matrix{u_1 cr u_2} = 0$ to denote a set of scalar equations. This also means that
538:      the dimension of the test function vector must be divisible by the number of space dimensions */
539:   int     numQuadPoints;     /* Number of points used for Gaussian quadrature */
540:   double *quadWeights;       /* Weights in the standard element for Gaussian quadrature */
541:   double *quadShapeFuncs;    /* Shape functions evaluated at quadrature points */
542:   double *quadTestFuncDers;  /* Test function derivatives evaluated at quadrature points */
543:   double  dxxi;              /* PartDer{x}{xi}  */
544:   double  dxix;              /* PartDer{xi}{x}  */
545:   double  dphix;             /* PartDer{phi_i}{x} times PartDer{phi_j}{x} */
546:   double  jac;               /* |J| for map to standard element */
547:   double  invjac;            /* |J^{-1}| for map from standard element */
548:   int     dim;               /* The problem dimension */
549:   int     comp;              /* The number of field components */
550:   int     tcomp;             /* The number of field components for the test field */
551:   int     funcs;             /* The number of shape functions */
552:   int     tfuncs;            /* The number of test functions */
553:   int     i, j, c, tc, f, p;

556:   /* Calculate element matrix entries by Gaussian quadrature --
557:            Since we integrate by parts here, the test and shape functions are switched */
558:   dim              = disc->dim;
559:   comp             = disc->comp;
560:   tcomp            = test->comp;
561:   funcs            = disc->funcs;
562:   tfuncs           = test->funcs;
563:   numQuadPoints    = disc->numQuadPoints;
564:   quadWeights      = disc->quadWeights;
565:   quadShapeFuncs   = disc->quadShapeFuncs;
566:   quadTestFuncDers = test->quadShapeFuncDers;
567:   for(p = 0; p < numQuadPoints; p++) {
568:     /* PartDer{x}{xi}(p)  = sum^{funcs}_{f=1} x_f PartDer{phi^f(p)}{xi} */
569:     dxxi = 0.0;
570:     for(f = 0; f < tfuncs; f++) {
571:       dxxi += coords[f*dim]*quadTestFuncDers[p*tfuncs*dim+f*dim];
572:     }
573:     jac = PetscAbsReal(dxxi);
574: #ifdef PETSC_USE_BOPT_g
575:     if (jac < 1.0e-14) {
576:       PetscPrintf(PETSC_COMM_SELF, "p: %d x1: %g x2: %gn", p, coords[0], coords[1]);
577:       SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
578:     }
579: #endif
580:     /* These are the elements of the inverse matrix */
581:     invjac =  1.0/jac;
582:     dxix   =  invjac;

584:     /* The rows are test functions */
585:     for(i = 0; i < tfuncs; i++) {
586:       /* We divide by the space dimension */
587:       for(tc = 0; tc < tcomp/dim; tc++) {
588:         /* The columns are shape functions */
589:         for(j = 0; j < funcs; j++) {
590:           for(c = 0; c < comp; c++) {
591:             dphix = quadTestFuncDers[p*tfuncs*dim+i*dim]*dxix;
592:             array[(i*tcomp+tc*dim+globalRowStart)*globalSize + j*comp+c+globalColStart] +=
593:               -alpha*dphix*quadShapeFuncs[p*funcs+j]*jac*quadWeights[p];
594:           }
595:         }
596:       }
597:     }
598:   }
599:   PetscLogFlops((2*tfuncs + 1 + 4*tfuncs*tcomp/dim*funcs*comp) * numQuadPoints);

601:   return(0);
602: }

604: int Divergence_Triangular_1D_Constant(Discretization disc, Discretization test, int rowSize, int colSize,
605:                                       int globalRowStart, int globalColStart, int globalSize, double *coords,
606:                                       PetscScalar alpha, PetscScalar *field, PetscScalar *array, void *ctx)
607: {
608:         /* We are using the convention that

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

613:                  and

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

618:                  where $d$ is the number of space dimensions. This agrees with the convention which allows
619:                  $Delta matrix{u_1 cr u_2} = 0$ to denote a set of scalar equations        This also requires that
620:      the dimension of a vector must be divisible by the space dimension in order to be acted upon by
621:      the divergence operator */
622:   int     numQuadPoints;     /* Number of points used for Gaussian quadrature */
623:   double *quadWeights;       /* Weights in the standard element for Gaussian quadrature */
624:   double *quadTestFuncs;     /* Test  functions evaluated at quadrature points */
625:   double *quadShapeFuncDers; /* Shape function derivatives evaluated at quadrature points */
626:   double *quadTestFuncDers;  /* Test function derivatives evaluated at quadrature points */
627:   double  dxxi;              /* PartDer{x}{xi}  */
628:   double  dxix;              /* PartDer{xi}{x}  */
629:   double  dphix;             /* PartDer{phi_i}{x} times PartDer{phi_j}{x} */
630:   double  jac;               /* |J| for map to standard element */
631:   double  invjac;            /* |J^{-1}| for map from standard element */
632:   int     dim;               /* The problem dimension */
633:   int     comp;              /* The number of field components */
634:   int     tcomp;             /* The number of field components for the test field */
635:   int     funcs;             /* The number of shape functions */
636:   int     tfuncs;            /* The number of test functions */
637:   int     i, j, c, tc, f, p;

640:   /* Calculate element matrix entries by Gaussian quadrature */
641:   dim              = disc->dim;
642:   comp              = disc->comp;
643:   tcomp             = test->comp;
644:   funcs             = disc->funcs;
645:   tfuncs            = test->funcs;
646:   numQuadPoints     = disc->numQuadPoints;
647:   quadWeights       = disc->quadWeights;
648:   quadTestFuncs     = test->quadShapeFuncs;
649:   quadShapeFuncDers = disc->quadShapeFuncDers;
650:   quadTestFuncDers  = test->quadShapeFuncDers;
651:   for(p = 0; p < numQuadPoints; p++) {
652:     /* PartDer{x}{xi}(p)  = sum^{funcs}_{f=1} x_f PartDer{phi^f(p)}{xi} */
653:     dxxi = 0.0;
654:     if (tfuncs >= funcs) {
655:       for(f = 0; f < tfuncs; f++) {
656:         dxxi += coords[f*dim]*quadTestFuncDers[p*tfuncs*dim+f*dim];
657:       }
658:     } else {
659:       for(f = 0; f < funcs; f++) {
660:         dxxi += coords[f*dim]*quadShapeFuncDers[p*funcs*dim+f*dim];
661:       }
662:     }
663:     jac = PetscAbsReal(dxxi);
664: #ifdef PETSC_USE_BOPT_g
665:     if (jac < 1.0e-14) {
666:       PetscPrintf(PETSC_COMM_SELF, "p: %d x1: %g x2: %gn", p, coords[0], coords[1]);
667:       SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
668:     }
669: #endif
670:     /* These are the elements of the inverse matrix */
671:     invjac = 1.0/jac;
672:     dxix   = invjac;

674:     /* The rows are test functions */
675:     for(i = 0; i < tfuncs; i++) {
676:       for(tc = 0; tc < tcomp; tc++) {
677:         /* The columns are shape functions */
678:         for(j = 0; j < funcs; j++) {
679:           dphix = quadShapeFuncDers[p*funcs*dim+j*dim]*dxix;
680:           /* We divide by the number of space dimensions */
681:           for(c = 0; c < comp/dim; c++) {
682:             array[(i*tcomp+tc+globalRowStart)*globalSize + j*comp+c*dim+globalColStart] +=
683:               alpha*dphix*quadTestFuncs[p*tfuncs+i]*jac*quadWeights[p];
684:           }
685:         }
686:       }
687:     }
688:   }
689:   PetscLogFlops((2*funcs + 1 + 4*tfuncs*tcomp*funcs*comp/dim) * numQuadPoints);

691:   return(0);
692: }

694: int DiscInterpolateField_Triangular_1D_Constant(Discretization disc, Mesh oldMesh, int elem, double x, double y, double z,
695:                                               PetscScalar *oldFieldVal, PetscScalar *newFieldVal, InterpolationType type)
696: {
697:   int comp = disc->comp;
698:   int rank;
699:   int neighbor, corner, c;

703:   MPI_Comm_rank(disc->comm, &rank);
704:   /* No scheme in place for boundary elements */
705:   for(corner = 0; corner < 2; corner++) {
706:     MeshGetElementNeighbor(oldMesh, elem, corner, &neighbor);
707:     if (neighbor < 0) {
708:       type = INTERPOLATION_LOCAL;
709:       break;
710:     }
711:   }

713:   switch (type) {
714:   case INTERPOLATION_LOCAL:
715:     for(c = 0 ; c < comp; c++) {
716:       newFieldVal[c] = oldFieldVal[0*comp+c];
717:     }
718:     PetscLogFlops(4+3*comp);
719:     break;
720:   default:
721:     SETERRQ1(PETSC_ERR_ARG_WRONG, "Unknown interpolation type %d", type);
722:   }
723: 
724:   return(0);
725: }

727: int DiscInterpolateElementVec_Triangular_1D_Constant(Discretization disc, ElementVec vec, Discretization newDisc, ElementVec newVec) {
728:   int          comp  = disc->comp;
729:   PetscScalar *array, *newArray;
730:   int          f, c;
731:   int          ierr;

734:   ElementVecGetArray(vec,    &array);
735:   ElementVecGetArray(newVec, &newArray);
736:   for(f = 0; f < newDisc->funcs; f++) {
737:     for(c = 0; c < comp; c++) {
738:       newArray[f*comp+c] = array[c];
739:     }
740:   }
741:   ElementVecRestoreArray(vec,    &array);
742:   ElementVecRestoreArray(newVec, &newArray);
743:   return(0);
744: }

746: /*
747:   DiscSetupQuadrature_Triangular_1D_Constant - Setup Gaussian quadrature with a 7 point integration rule

749:   Input Parameter:
750: . disc - The Discretization
751: */
752: int DiscSetupQuadrature_Triangular_1D_Constant(Discretization disc) {
753:   int dim   = disc->dim;
754:   int funcs = disc->funcs;
755:   int p;

759:   disc->numQuadPoints = 7;
760:   PetscMalloc(disc->numQuadPoints*dim       * sizeof(double), &disc->quadPoints);
761:   PetscMalloc(disc->numQuadPoints           * sizeof(double), &disc->quadWeights);
762:   PetscMalloc(disc->numQuadPoints*funcs     * sizeof(double), &disc->quadShapeFuncs);
763:   PetscMalloc(disc->numQuadPoints*funcs*dim * sizeof(double), &disc->quadShapeFuncDers);
764:   PetscLogObjectMemory(disc, (disc->numQuadPoints*(funcs*(dim+1) + dim+1)) * sizeof(double));
765:   disc->quadPoints[0]  = 0.0254460438286207377369052;
766:   disc->quadWeights[0] = 0.0647424830844348466353057;
767:   disc->quadPoints[1]  = 0.1292344072003027800680676;
768:   disc->quadWeights[1] = 0.1398526957446383339507339;
769:   disc->quadPoints[2]  = 0.29707742431130141654669679;
770:   disc->quadWeights[2] = 0.1909150252525594724751849;
771:   disc->quadPoints[3]  = 0.5000000000000000000000000;
772:   disc->quadWeights[3] = 0.2089795918367346938775510;
773:   disc->quadPoints[4]  = 0.70292257568869858345330321;
774:   disc->quadWeights[4] = disc->quadWeights[2];
775:   disc->quadPoints[5]  = 0.8707655927996972199319324;
776:   disc->quadWeights[5] = disc->quadWeights[1];
777:   disc->quadPoints[6]  = 0.9745539561713792622630948;
778:   disc->quadWeights[6] = disc->quadWeights[0];
779:   for(p = 0; p < disc->numQuadPoints; p++) {
780:     /* phi^0: 1 */
781:     disc->quadShapeFuncs[p*funcs]              = 1.0;
782:     disc->quadShapeFuncDers[p*funcs*dim+0*dim] = 0.0;
783:   }
784:   return(0);
785: }

787: /*
788:   DiscSetupOperators_Triangular_1D_Constant - Setup the default operators

790:   Input Parameter:
791: . disc - The Discretization
792: */
793: int DiscSetupOperators_Triangular_1D_Constant(Discretization disc) {
794:   int          comp = disc->comp;
795:   int          size = disc->size;
796:   PetscScalar *precompInt;
797:   int          newOp;
798:   int          c, i, j;
799:   int          ierr;

802:   /* The Identity operator I -- the matrix is symmetric */
803:   PetscMalloc(size*size * sizeof(PetscScalar), &precompInt);
804:   PetscLogObjectMemory(disc, size*size * sizeof(PetscScalar));
805:   PetscMemzero(precompInt, size*size * sizeof(PetscScalar));
806:   for(c = 0; c < comp; c++) {
807:     precompInt[(0*comp+c)*size + 0*comp+c] = 1.0;
808:   }
809:   for(i = 0; i < size; i++) {
810:     for(j = 0; j < i; j++) {
811:       precompInt[i*size + j] = precompInt[j*size + i];
812:     }
813:   }
814:   DiscretizationRegisterPrecomputedOperator(disc, precompInt, &newOp);
815:   if (newOp != IDENTITY) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "Default operator %d not setup correctly", IDENTITY);
816:   /* The Laplacian operator Delta -- the matrix is symmetric */
817:   DiscretizationRegisterOperator(disc, Laplacian_Triangular_1D_Constant, &newOp);
818:   if (newOp != LAPLACIAN) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "Default operator %d not setup correctly", LAPLACIAN);
819:   /* The Gradient operator nabla -- the matrix is rectangular */
820:   DiscretizationRegisterOperator(disc, Gradient_Triangular_1D_Constant, &newOp);
821:   if (newOp != GRADIENT) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "Default operator %d not setup correctly", GRADIENT);
822:   /* The Divergence operator nablacdot -- the matrix is rectangular */
823:   DiscretizationRegisterOperator(disc, Divergence_Triangular_1D_Constant, &newOp);
824:   if (newOp != DIVERGENCE) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "Default operator %d not setup correctly", DIVERGENCE);
825:   /* The weighted Laplacian operator -- the matrix is symmetric */
826:   DiscretizationRegisterOperator(disc, Weighted_Laplacian_Triangular_1D_Constant, &newOp);
827:   if (newOp != WEIGHTED_LAP) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "Default operator %d not setup correctly", WEIGHTED_LAP);
828:   return(0);
829: }

831: static struct _DiscretizationOps DOps = {PETSC_NULL/* DiscretizationSetup */,
832:                                          DiscSetupOperators_Triangular_1D_Constant,
833:                                          PETSC_NULL/* DiscretizationSetFromOptions */,
834:                                          DiscView_Triangular_1D_Constant,
835:                                          DiscDestroy_Triangular_1D_Constant,
836:                                          DiscEvaluateFunctionGalerkin_Triangular_1D_Constant,
837:                                          DiscEvaluateOperatorGalerkin_Triangular_1D_Constant,
838:                                          DiscEvaluateALEOperatorGalerkin_Triangular_1D_Constant,
839:                                          DiscEvaluateNonlinearOperatorGalerkin_Triangular_1D_Constant,
840:                                          DiscEvaluateNonlinearALEOperatorGalerkin_Triangular_1D_Constant,
841:                                          DiscInterpolateField_Triangular_1D_Constant,
842:                                          DiscInterpolateElementVec_Triangular_1D_Constant};

844: EXTERN_C_BEGIN
845: int DiscCreate_Triangular_1D_Constant(Discretization disc) {
846:   int arg;

850:   if (disc->comp <= 0) {
851:     SETERRQ(PETSC_ERR_ARG_WRONG, "Discretization must have at least 1 component. Call DiscretizationSetNumComponents() to set this.");
852:   }
853:   PetscMemcpy(disc->ops, &DOps, sizeof(struct _DiscretizationOps));
854:   disc->dim   = 1;
855:   disc->funcs = 1;
856:   disc->size  = disc->funcs*disc->comp;

858:   DiscretizationSetupDefaultOperators(disc);
859:   DiscSetupQuadrature_Triangular_1D_Constant(disc);

861:   /* Storage */
862:   PetscMalloc(disc->comp * sizeof(PetscScalar),   &disc->funcVal);
863:   PetscMalloc(2          * sizeof(PetscScalar *), &disc->fieldVal);
864:   for(arg = 0; arg < 2; arg++) {
865:     PetscMalloc(disc->comp*(disc->dim+1) * sizeof(PetscScalar), &disc->fieldVal[arg]);
866:   }
867:   return(0);
868: }
869: EXTERN_C_END