Actual source code: bdlinear.c

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

  5: /*
  6:    Defines piecewise linear 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----2
 22: */

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

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

 42: #undef  __FUNCT__
 44: static int BoundaryDiscView_Triangular_2D_Linear(Discretization disc, PetscViewer viewer) {
 45:   PetscTruth isascii;
 46:   int        ierr;

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

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

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

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

 94: #ifdef PETSC_USE_BOPT_g
 95:   if ((edge < 0) || (edge >= mesh->numEdges)) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Invalid edge");
 96: #endif
 97:   /* Calculate the determinant of the inverse Jacobian of the map to the standard element
 98:      which must be a constant for linear elements */
 99:   x11 = nodes[edges[edge*2]  *2];
100:   y11 = nodes[edges[edge*2]  *2+1];
101:   if (mesh->isPeriodic == PETSC_TRUE) {
102:     x21 = MeshPeriodicDiffX(mesh, nodes[edges[edge*2+1]*2]   - x11);
103:     y21 = MeshPeriodicDiffY(mesh, nodes[edges[edge*2+1]*2+1] - y11);
104:   } else {
105:     x21 = nodes[edges[edge*2+1]*2]   - x11;
106:     y21 = nodes[edges[edge*2+1]*2+1] - y11;
107:   }
108:   jac = sqrt(x21*x21 + y21*y21);
109:   if (jac < 1.0e-14) {
110:     PetscPrintf(PETSC_COMM_SELF, "[%d]edge: %d x21: %g y21: %gn", rank, edge, x21, y21);
111:     SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
112:   }
113: #ifdef PETSC_USE_BOPT_g
114:   PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
115:   if (opt == PETSC_TRUE) {
116:     PetscPrintf(PETSC_COMM_SELF, "[%d]edge: %d x21: %g y21: %g jac: %gn", rank, edge, x21, y21, jac);
117:   }
118: #endif

120:   /* Calculate element vector entries by Gaussian quadrature */
121:   for(p = 0; p < numQuadPoints; p++) {
122:     x = MeshPeriodicX(mesh, x21*quadPoints[p] + x11);
123:     y = MeshPeriodicY(mesh, y21*quadPoints[p] + y11);
124:     (*f)(1, comp, &x, &y, PETSC_NULL, funcVal, uCtx);
125: #ifdef PETSC_USE_BOPT_g
126:     PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
127:     if (opt == PETSC_TRUE) {
128:       PetscPrintf(PETSC_COMM_SELF, "[%d]p:%d jac: %g", rank, p, jac);
129:       for(j = 0; j < comp; j++)
130:         PetscPrintf(PETSC_COMM_SELF, " func[%d]: %g", j, PetscRealPart(funcVal[j]));
131:       PetscPrintf(PETSC_COMM_SELF, "n");
132:     }
133: #endif

135:     for(i = 0, k = 0; i < funcs; i++) {
136:       for(j = 0; j < comp; j++, k++) {
137:         array[k] += alpha*funcVal[j]*quadShapeFuncs[p*funcs+i]*jac*quadWeights[p];
138: #ifdef PETSC_USE_BOPT_g
139:         PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
140:         if (opt == PETSC_TRUE) {
141:           PetscPrintf(PETSC_COMM_SELF, "[%d]  array[%d]: %gn", rank, k, PetscRealPart(array[k]));
142:         }
143: #endif
144:       }
145:     }
146:   }
147:   PetscLogFlops(6 + (4 + 5*funcs*comp) * numQuadPoints);
148:   return(0);
149: }

151: #undef  __FUNCT__
153: static int BoundaryDiscEvaluateOperatorGalerkin_Triangular_2D_Linear(Discretization disc, Mesh mesh, int elemSize,
154:                                                                      int rowStart, int colStart, int op, PetscScalar alpha,
155:                                                                      int edge, PetscScalar *field, PetscScalar *mat, void *ctx)
156: {
157:   Mesh_Triangular *tri        = (Mesh_Triangular *) mesh->data;
158:   EdgeContext     *bdCtx      = (EdgeContext *)     ctx;
159:   void            *uCtx       = bdCtx->ctx;
160:   double          *nodes      = tri->nodes;
161:   int             *edges      = tri->edges;
162:   Operator         oper       = disc->operators[op]; /* The operator to discretize */
163:   Discretization   test       = oper->test;          /* The space of test functions */
164:   int              rowSize    = test->size;          /* The number of shape functions per element */
165:   int              colSize    = disc->size;          /* The number of test  functions per element */
166:   PetscScalar     *precompInt = oper->precompInt;    /* Precomputed integrals of the operator on shape functions */
167:   OperatorFunction opFunc     = oper->opFunc;        /* Integrals of operators which depend on J */
168:   double           x21, y21;                         /* Coordinates of the element, with point 1 at the origin */
169:   double           jac;                              /* |J| for map to standard element */
170:   double           coords[MAX_CORNERS*2];            /* Coordinates of the element */
171:   int              i, j;
172:   int              rank;
173:   int              ierr;

176:   MPI_Comm_rank(disc->comm, &rank);
177: #ifdef PETSC_USE_BOPT_g
178:   /* Check for valid operator */
179:   if ((op < 0) || (op >= disc->numOps) || (!disc->operators[op])) SETERRQ(PETSC_ERR_ARG_WRONG, "Invalid operator");
180: #endif

182:   if (precompInt != PETSC_NULL) {
183:     /* Calculate the determinant of the inverse Jacobian of the map to the standard element
184:        which must be a constant for linear elements */
185:     if (mesh->isPeriodic == PETSC_TRUE) {
186:       x21 = MeshPeriodicDiffX(mesh, nodes[edges[edge*2+1]*2]   - nodes[edges[edge*2]*2]);
187:       y21 = MeshPeriodicDiffY(mesh, nodes[edges[edge*2+1]*2+1] - nodes[edges[edge*2]*2+1]);
188:     } else {
189:       x21 = nodes[edges[edge*2+1]*2]   - nodes[edges[edge*2]*2];
190:       y21 = nodes[edges[edge*2+1]*2+1] - nodes[edges[edge*2]*2+1];
191:     }
192:     jac = sqrt(x21*x21 + y21*y21);
193:     if (jac < 1.0e-14) {
194:       PetscPrintf(PETSC_COMM_SELF, "[%d]x21: %g y21: %g jac: %gn", rank, x21, y21, jac);
195:       SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
196:     }

198:     /* Calculate element matrix entries which are all precomputed */
199:     for(i = 0; i < rowSize; i++)
200:       for(j = 0; j < colSize; j++)
201:         mat[(i+rowStart)*elemSize + j+colStart] += alpha*precompInt[i*colSize + j]*jac;
202:     PetscLogFlops(6 + 2*rowSize*colSize);
203:   } else {
204:     if (opFunc == PETSC_NULL) SETERRQ(PETSC_ERR_ARG_CORRUPT, "Invalid function");
205:     if (mesh->isPeriodic == PETSC_TRUE) {
206:       coords[0*2+0] = nodes[edges[edge*2+0]*2];
207:       coords[0*2+1] = nodes[edges[edge*2+0]*2+1];
208:       coords[1*2+0] = MeshPeriodicRelativeX(mesh, nodes[edges[edge*2+1]*2],   coords[0*2+0]);
209:       coords[1*2+1] = MeshPeriodicRelativeY(mesh, nodes[edges[edge*2+1]*2+1], coords[0*2+1]);
210:     } else {
211:       coords[0*2+0] = nodes[edges[edge*2+0]*2];
212:       coords[0*2+1] = nodes[edges[edge*2+0]*2+1];
213:       coords[1*2+0] = nodes[edges[edge*2+1]*2];
214:       coords[1*2+1] = nodes[edges[edge*2+1]*2+1];
215:     }

217:     (*opFunc)(disc, test, rowSize, colSize, rowStart, colStart, elemSize, coords, alpha, field, mat, uCtx);
218: 
219:   }
220:   return(0);
221: }

223: #undef  __FUNCT__
225: static int BoundaryDiscEvaluateALEOperatorGalerkin_Triangular_2D_Linear(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_Linear(Discretization disc, Mesh mesh,
236:                                                                               NonlinearOperator f, PetscScalar alpha,
237:                                                                               int elem, int numArgs, PetscScalar **field, PetscScalar *vec, void *ctx)
238: {
239:   SETERRQ(PETSC_ERR_SUP, " ");
240: }

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

252: #undef  __FUNCT__
254: int BoundaryLaplacian_Triangular_2D_Linear(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:   double      x21, y21; /* Coordinates of the element, with point 1 at the origin */
259:   double      jac;      /* |J| for map to standard element */
260:   PetscScalar w;        /* 1/(2 jac) */
261:   int         comp;     /* Number of components */
262:   int         i;

265:   /* Calculate the determinant of the inverse Jacobian of the map to the standard element
266:      which must be a constant for linear elements - 1/|x_{21} y_{31} - x_{31} y_{21}| */
267:   x21 = coords[2] - coords[0];
268:   y21 = coords[3] - coords[1];
269:   jac = sqrt(x21*x21 + y21*y21);
270: #ifdef PETSC_USE_BOPT_g
271:   if (jac < 1.0e-14) {
272:     PetscPrintf(PETSC_COMM_SELF, "x21: %g y21: %g jac: %gn", x21, y21, jac);
273:     SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
274:   }
275: #endif

277:   comp = rowSize/2;
278:   w  = 1.0/(2.0*jac);
279:   w *= alpha;
280:   for(i = 0; i < comp; i++)
281:   {
282:     /* phi^1 phi^1 */
283:     array[(0*comp+i+globalRowStart)*globalSize + 0*comp+i+globalColStart] = 0.0;
284:     /* phi^1 phi^2 */
285:     array[(0*comp+i+globalRowStart)*globalSize + 1*comp+i+globalColStart] = 0.0;
286:     /* phi^2 phi^1 */
287:     array[(1*comp+i+globalRowStart)*globalSize + 0*comp+i+globalColStart] =
288:                         array[(0*comp+i+globalRowStart)*globalSize + 1*comp+i+globalColStart];
289:     /* phi^2 phi^2 */
290:     array[(1*comp+i+globalRowStart)*globalSize + 1*comp+i+globalColStart] = 0.0;
291:   }
292:   PetscLogFlops(47);

294:   return(0);
295: }

297: #undef  __FUNCT__
299: int BoundaryDiscInterpolateField_Triangular_2D_Linear(Discretization disc, Mesh oldMesh, int elem, double x, double y, double z,
300:                                                       PetscScalar *oldFieldVal, PetscScalar *newFieldVal, InterpolationType type)
301: {
302:   SETERRQ(PETSC_ERR_SUP, " ");
303: }

305: #undef  __FUNCT__
307: int BoundaryDiscInterpolateElementVec_Triangular_2D_Linear(Discretization disc, ElementVec vec,
308:                                                            Discretization newDisc, ElementVec newVec)
309: {
310:   SETERRQ(PETSC_ERR_SUP, " ");
311: }


314: #undef  __FUNCT__
316: /*
317:   BoundaryDiscSetupQuadrature_Triangular_2D_Linear - Setup Gaussian quadrature with a 5 point integration rule

319:   Input Parameter:
320: . disc - The Discretization
321: */
322: int BoundaryDiscSetupQuadrature_Triangular_2D_Linear(Discretization disc) {
323:   int dim   = disc->dim;
324:   int funcs = disc->funcs;
325:   int p;

329:   disc->numQuadPoints = 5;
330:   PetscMalloc(disc->numQuadPoints*dim       * sizeof(double), &disc->quadPoints);
331:   PetscMalloc(disc->numQuadPoints           * sizeof(double), &disc->quadWeights);
332:   PetscMalloc(disc->numQuadPoints*funcs     * sizeof(double), &disc->quadShapeFuncs);
333:   PetscMalloc(disc->numQuadPoints*funcs*dim * sizeof(double), &disc->quadShapeFuncDers);
334:   PetscLogObjectMemory(disc, (disc->numQuadPoints*(funcs*(dim+1) + dim+1)) * sizeof(double));
335:   disc->quadPoints[0]  = 0.0469101;
336:   disc->quadWeights[0] = 0.118463;
337:   disc->quadPoints[1]  = 0.230765;
338:   disc->quadWeights[1] = 0.239314;
339:   disc->quadPoints[2]  = 0.500000;
340:   disc->quadWeights[2] = 0.284444;
341:   disc->quadPoints[3]  = 0.769235;
342:   disc->quadWeights[3] = 0.239314;
343:   disc->quadPoints[4]  = 0.953090;
344:   disc->quadWeights[4] = 0.118463;
345:   for(p = 0; p < disc->numQuadPoints; p++) {
346:     /* phi^0: 1 - xi */
347:     disc->quadShapeFuncs[p*funcs+0]    =  1.0 - disc->quadPoints[p];
348:     disc->quadShapeFuncDers[p*funcs+0] = -1.0;
349:     /* phi^1: xi */
350:     disc->quadShapeFuncs[p*funcs+1]    =  disc->quadPoints[p];
351:     disc->quadShapeFuncDers[p*funcs+1] =  1.0;
352:   }
353:   return(0);
354: }

356: #undef  __FUNCT__
358: /*
359:   BoundaryDiscSetupOperators_Triangular_2D_Linear - Setup the default operators

361:   Input Parameter:
362: . disc - The Discretization
363: */
364: int BoundaryDiscSetupOperators_Triangular_2D_Linear(Discretization disc) {
365:   int          comp = disc->comp;
366:   int          size = disc->size;
367:   PetscScalar *precompInt;
368:   int          newOp;
369:   int          c, i, j;
370:   int          ierr;

373:   /* The Identity operator I -- the matrix is symmetric */
374:   PetscMalloc(size*size * sizeof(PetscScalar), &precompInt);
375:   PetscLogObjectMemory(disc, size*size * sizeof(PetscScalar));
376:   PetscMemzero(precompInt, size*size * sizeof(PetscScalar));
377:   for(c = 0; c < comp; c++) {
378:     precompInt[(0*comp+c)*size + 0*comp+c] = 1.0/3.0;
379:     precompInt[(0*comp+c)*size + 1*comp+c] = 1.0/6.0;
380:     precompInt[(1*comp+c)*size + 1*comp+c] = 1.0/3.0;
381:   }
382:   for(i = 0; i < size; i++) {
383:     for(j = 0; j < i; j++) {
384:       precompInt[i*size + j] = precompInt[j*size + i];
385:     }
386:   }
387:   DiscretizationRegisterPrecomputedOperator(disc, precompInt, &newOp);
388:   if (newOp != IDENTITY) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "Default operator %d not setup correctly", IDENTITY);
389:   /* The Laplacian operator Delta -- the matrix is symmetric */
390:   DiscretizationRegisterOperator(disc, BoundaryLaplacian_Triangular_2D_Linear, &newOp);
391:   if (newOp != LAPLACIAN) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "Default operator %d not setup correctly", LAPLACIAN);
392:   /* The Gradient operator nabla -- the matrix is rectangular */
393:   DiscretizationRegisterOperator(disc, PETSC_NULL, &newOp);
394:   if (newOp != GRADIENT) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "Default operator %d not setup correctly", GRADIENT);
395:   /* The Divergence operator nablacdot -- the matrix is rectangular */
396:   DiscretizationRegisterOperator(disc, PETSC_NULL, &newOp);
397:   if (newOp != DIVERGENCE) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "Default operator %d not setup correctly", DIVERGENCE);
398:   /* The weighted Laplacian operator -- the matrix is symmetric */
399:   DiscretizationRegisterOperator(disc, PETSC_NULL, &newOp);
400:   if (newOp != WEIGHTED_LAP) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "Default operator %d not setup correctly", WEIGHTED_LAP);
401:   return(0);
402: }

404: static struct _DiscretizationOps DOps = {PETSC_NULL/* DiscretizationSetup */,
405:                                          BoundaryDiscSetupOperators_Triangular_2D_Linear,
406:                                          PETSC_NULL/* DiscretizationSetFromOptions */,
407:                                          BoundaryDiscView_Triangular_2D_Linear,
408:                                          BoundaryDiscDestroy_Triangular_2D_Linear,
409:                                          BoundaryDiscEvaluateFunctionGalerkin_Triangular_2D_Linear,
410:                                          BoundaryDiscEvaluateOperatorGalerkin_Triangular_2D_Linear,
411:                                          BoundaryDiscEvaluateALEOperatorGalerkin_Triangular_2D_Linear,
412:                                          BoundaryDiscEvaluateNonlinearOperatorGalerkin_Triangular_2D_Linear,
413:                                          BoundaryDiscEvaluateNonlinearALEOperatorGalerkin_Triangular_2D_Linear,
414:                                          BoundaryDiscInterpolateField_Triangular_2D_Linear,
415:                                          BoundaryDiscInterpolateElementVec_Triangular_2D_Linear};

417: EXTERN_C_BEGIN
418: #undef  __FUNCT__
420: int BoundaryDiscCreate_Triangular_2D_Linear(Discretization disc) {
421:   int arg;

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

433:   DiscretizationSetupDefaultOperators(disc);
434:   BoundaryDiscSetupQuadrature_Triangular_2D_Linear(disc);

436:   /* Storage */
437:   PetscMalloc(disc->comp * sizeof(PetscScalar),   &disc->funcVal);
438:   PetscMalloc(2          * sizeof(PetscScalar *), &disc->fieldVal);
439:   for(arg = 0; arg < 2; arg++) {
440:     PetscMalloc(disc->comp*(disc->dim+1) * sizeof(PetscScalar), &disc->fieldVal[arg]);
441:   }
442:   return(0);
443: }
444: EXTERN_C_END