Actual source code: linear.c
1: #ifdef PETSC_RCS_HEADER
2: static char vcid[] = "$Id: linear.c,v 1.7 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: extern int DiscTransformCoords_Triangular_2D_Quadratic(double, double, double *, double *, double *);
15: /* For precomputed integrals, the table is structured as follows:
17: precompInt[op,i,j] = int_{SE} <op phi^i(xi,eta), phi^j(xi,eta)> |J^{-1}|
19: where we recall that |J| is a constant for linear affine maps,
20: and the map of any triangle to the standard element is linear.
21: The numbering of the nodes in the standard element is
23: 3
24: |
25: |
26: |
27: |
28: 1----2
29: */
31: static int DiscDestroy_Triangular_2D_Linear(Discretization disc) {
33: return(0);
34: }
36: static int DiscView_Triangular_2D_Linear_File(Discretization disc, PetscViewer viewer) {
38: PetscViewerASCIIPrintf(viewer, "Linear discretizationn");
39: PetscViewerASCIIPrintf(viewer, " %d shape functions per componentn", disc->funcs);
40: PetscViewerASCIIPrintf(viewer, " %d registered operatorsn", disc->numOps);
41: return(0);
42: }
44: static int DiscView_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: DiscView_Triangular_2D_Linear_File(disc, viewer);
52: }
53: return(0);
54: }
56: static int DiscEvaluateFunctionGalerkin_Triangular_2D_Linear(Discretization disc, Mesh mesh, PointFunction f, PetscScalar alpha,
57: int elem, PetscScalar *array, void *ctx)
58: {
59: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
60: double *nodes = tri->nodes;
61: int *elements = tri->faces;
62: int corners = mesh->numCorners;
63: int comp = disc->comp; /* The number of components in this field */
64: int funcs = disc->funcs; /* The number of shape functions per component */
65: PetscScalar *funcVal = disc->funcVal; /* Function value at a quadrature point */
66: int numQuadPoints = disc->numQuadPoints; /* Number of points used for Gaussian quadrature */
67: double *quadPoints = disc->quadPoints; /* Points in the standard element for Gaussian quadrature */
68: double *quadWeights = disc->quadWeights; /* Weights in the standard element for Gaussian quadrature */
69: double *quadShapeFuncs = disc->quadShapeFuncs; /* Shape function evaluated at quadrature points */
70: double jac; /* |J| for map to standard element */
71: double x, y; /* The integration point */
72: double x11, y11, x21, y21, x31, y31;
73: int rank = -1;
74: int i, j, k, p;
75: #ifdef PETSC_USE_BOPT_g
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: if ((elem < 0) || (elem >= mesh->part->numOverlapElements)) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Invalid element");
93: #endif
94: /* Calculate the determinant of the inverse Jacobian of the map to the standard element
95: which must be a constant for linear elements */
96: x11 = nodes[elements[elem*corners]*2];
97: y11 = nodes[elements[elem*corners]*2+1];
98: if (mesh->isPeriodic == PETSC_TRUE) {
99: x21 = MeshPeriodicDiffX(mesh, nodes[elements[elem*corners+1]*2] - x11);
100: x31 = MeshPeriodicDiffX(mesh, nodes[elements[elem*corners+2]*2] - x11);
101: y21 = MeshPeriodicDiffY(mesh, nodes[elements[elem*corners+1]*2+1] - y11);
102: y31 = MeshPeriodicDiffY(mesh, nodes[elements[elem*corners+2]*2+1] - y11);
103: } else {
104: x21 = nodes[elements[elem*corners+1]*2] - x11;
105: x31 = nodes[elements[elem*corners+2]*2] - x11;
106: y21 = nodes[elements[elem*corners+1]*2+1] - y11;
107: y31 = nodes[elements[elem*corners+2]*2+1] - y11;
108: }
109: jac = PetscAbsReal(x21*y31 - x31*y21);
110: if (jac < 1.0e-14) {
111: PetscPrintf(PETSC_COMM_SELF, "[%d]elem: %d x21: %g y21: %g x31: %g y31: %gn", rank, elem, x21, y21, x31, y31);
112: SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
113: }
114: #ifdef PETSC_USE_BOPT_g
115: PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
116: if (opt == PETSC_TRUE) {
117: PetscPrintf(PETSC_COMM_SELF, "[%d]elem: %d x21: %g y21: %g x31: %g y31: %g jac: %gn",
118: rank, elem, x21, y21, x31, y31, jac);
119: }
120: #endif
122: /* Calculate element vector entries by Gaussian quadrature */
123: for(p = 0; p < numQuadPoints; p++) {
124: x = MeshPeriodicX(mesh, x21*quadPoints[p*2] + x31*quadPoints[p*2+1] + x11);
125: y = MeshPeriodicY(mesh, y21*quadPoints[p*2] + y31*quadPoints[p*2+1] + y11);
126: (*f)(1, comp, &x, &y, PETSC_NULL, funcVal, ctx);
127: #ifdef PETSC_USE_BOPT_g
128: PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
129: if (opt == PETSC_TRUE) {
130: PetscPrintf(PETSC_COMM_SELF, "[%d]p:%d jac: %g", rank, p, jac);
131: for(j = 0; j < comp; j++)
132: PetscPrintf(PETSC_COMM_SELF, " func[%d]: %g", j, PetscRealPart(funcVal[j]));
133: PetscPrintf(PETSC_COMM_SELF, "n");
134: }
135: #endif
137: for(i = 0, k = 0; i < funcs; i++) {
138: for(j = 0; j < comp; j++, k++) {
139: array[k] += alpha*funcVal[j]*quadShapeFuncs[p*funcs+i]*jac*quadWeights[p];
140: #ifdef PETSC_USE_BOPT_g
141: PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
142: if (opt == PETSC_TRUE) {
143: PetscPrintf(PETSC_COMM_SELF, "[%d] array[%d]: %gn", rank, k, PetscRealPart(array[k]));
144: }
145: #endif
146: }
147: }
148: }
149: PetscLogFlops(7 + (8 + 5*funcs*comp) * numQuadPoints);
150: return(0);
151: }
153: static int DiscEvaluateOperatorGalerkin_Triangular_2D_Linear(Discretization disc, Mesh mesh, int elemSize,
154: int rowStart, int colStart, int op, PetscScalar alpha,
155: int elem, PetscScalar *field, PetscScalar *mat, void *ctx)
156: {
157: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
158: double *nodes = tri->nodes; /* The node coordinates */
159: int *elements = tri->faces; /* The element corners */
160: int numCorners = mesh->numCorners; /* The number of corners per element */
161: Operator oper = disc->operators[op]; /* The operator to discretize */
162: Discretization test = oper->test; /* The space of test functions */
163: OperatorFunction opFunc = oper->opFunc; /* Integrals of operators which depend on J */
164: PetscScalar *precompInt = oper->precompInt; /* Precomputed integrals of the operator on shape functions */
165: int rowSize = test->size; /* The number of shape functions per element */
166: int colSize = disc->size; /* The number of test functions per element */
167: double x21, x31, y21, y31; /* Coordinates of the element, with point 1 at the origin */
168: double jac; /* |J| for map to standard element */
169: double coords[MAX_CORNERS*2]; /* Coordinates of the element */
170: int rank;
171: int i, j, f;
172: int ierr;
175: MPI_Comm_rank(disc->comm, &rank);
176: #ifdef PETSC_USE_BOPT_g
177: /* Check for valid operator */
178: if ((op < 0) || (op >= disc->numOps) || (!disc->operators[op])) SETERRQ(PETSC_ERR_ARG_WRONG, "Invalid operator");
179: #endif
181: if (precompInt != PETSC_NULL)
182: {
183: /* Calculate the determinant of the inverse Jacobian of the map to the standard element
184: which must be a constant for linear elements - 1/|x_{21} y_{31} - x_{31} y_{21}| */
185: if (mesh->isPeriodic == PETSC_TRUE) {
186: x21 = MeshPeriodicDiffX(mesh, nodes[elements[elem*numCorners+1]*2] - nodes[elements[elem*numCorners]*2]);
187: x31 = MeshPeriodicDiffX(mesh, nodes[elements[elem*numCorners+2]*2] - nodes[elements[elem*numCorners]*2]);
188: y21 = MeshPeriodicDiffY(mesh, nodes[elements[elem*numCorners+1]*2+1] - nodes[elements[elem*numCorners]*2+1]);
189: y31 = MeshPeriodicDiffY(mesh, nodes[elements[elem*numCorners+2]*2+1] - nodes[elements[elem*numCorners]*2+1]);
190: } else {
191: x21 = nodes[elements[elem*numCorners+1]*2] - nodes[elements[elem*numCorners]*2];
192: x31 = nodes[elements[elem*numCorners+2]*2] - nodes[elements[elem*numCorners]*2];
193: y21 = nodes[elements[elem*numCorners+1]*2+1] - nodes[elements[elem*numCorners]*2+1];
194: y31 = nodes[elements[elem*numCorners+2]*2+1] - nodes[elements[elem*numCorners]*2+1];
195: }
196: jac = PetscAbsReal(x21*y31 - x31*y21);
197: if (jac < 1.0e-14) {
198: PetscPrintf(PETSC_COMM_SELF, "[%d]x21: %g y21: %g x31: %g y31: %g jac: %gn", rank, x21, y21, x31, y31, jac);
199: SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
200: }
202: /* Calculate element matrix entries which are all precomputed */
203: for(i = 0; i < rowSize; i++)
204: for(j = 0; j < colSize; j++)
205: mat[(i+rowStart)*elemSize + j+colStart] += alpha*precompInt[i*colSize + j]*jac;
206: PetscLogFlops(7 + 2*rowSize*colSize);
207: }
208: else
209: {
210: if (opFunc == PETSC_NULL) SETERRQ(PETSC_ERR_ARG_CORRUPT, "Invalid function");
211: if (mesh->isPeriodic == PETSC_TRUE) {
212: coords[0*2+0] = nodes[elements[elem*numCorners+0]*2+0];
213: coords[0*2+1] = nodes[elements[elem*numCorners+0]*2+1];
214: for(f = 1; f < PetscMax(disc->funcs, test->funcs); f++) {
215: coords[f*2+0] = MeshPeriodicRelativeX(mesh, nodes[elements[elem*numCorners+f]*2+0], coords[0*2+0]);
216: coords[f*2+1] = MeshPeriodicRelativeY(mesh, nodes[elements[elem*numCorners+f]*2+1], coords[0*2+1]);
217: }
218: } else {
219: for(f = 0; f < PetscMax(disc->funcs, test->funcs); f++) {
220: coords[f*2+0] = nodes[elements[elem*numCorners+f]*2+0];
221: coords[f*2+1] = nodes[elements[elem*numCorners+f]*2+1];
222: }
223: }
225: (*opFunc)(disc, test, rowSize, colSize, rowStart, colStart, elemSize, coords, alpha, field, mat, ctx);
226:
227: }
228: return(0);
229: }
231: static int DiscEvaluateNonlinearOperatorGalerkin_Triangular_2D_Linear(Discretization disc, Mesh mesh, NonlinearOperator f,
232: PetscScalar alpha, int elem, int numArgs, PetscScalar **field,
233: PetscScalar *vec, void *ctx)
234: {
235: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
236: double *nodes = tri->nodes; /* The node coordinates */
237: int *elements = tri->faces; /* The element corners */
238: int numCorners = mesh->numCorners; /* The number of corners per element */
239: int comp = disc->comp; /* The number of components in this field */
240: int funcs = disc->funcs; /* The number of shape functions per component */
241: PetscScalar *funcVal = disc->funcVal; /* Function value at a quadrature point */
242: PetscScalar **fieldVal = disc->fieldVal; /* Field value and derivatives at a quadrature point */
243: double jac; /* |J| for map to standard element */
244: double invjac; /* |J^{-1}| for map from standard element */
245: int numQuadPoints; /* Number of points used for Gaussian quadrature */
246: double *quadPoints; /* Points in the standard element for Gaussian quadrature */
247: double *quadWeights; /* Weights in the standard element for Gaussian quadrature */
248: double *quadShapeFuncs; /* Shape function evaluated at quadrature points */
249: double *quadShapeFuncDers; /* Shape function derivatives evaluated at quadrature points */
250: double x, y; /* The integration point */
251: double dxix; /* PartDer{xi}{x} */
252: double detx; /* PartDer{eta}{x} */
253: double dxiy; /* PartDer{xi}{y} */
254: double dety; /* PartDer{eta}{y} */
255: PetscScalar dfxi; /* PartDer{field}{xi} */
256: PetscScalar dfet; /* PartDer{field}{eta} */
257: double x11, y11, x21, y21, x31, y31;
258: int rank = -1;
259: int i, j, k, p, func, arg;
260: #ifdef PETSC_USE_BOPT_g
261: PetscTruth opt;
262: #endif
263: int ierr;
266: if (numArgs > 2) SETERRQ(PETSC_ERR_SUP, "Only configured to handle two nonlinear arguments");
267: MPI_Comm_rank(disc->comm, &rank);
268: numQuadPoints = disc->numQuadPoints;
269: quadPoints = disc->quadPoints;
270: quadWeights = disc->quadWeights;
271: quadShapeFuncs = disc->quadShapeFuncs;
272: quadShapeFuncDers = disc->quadShapeFuncDers;
273:
274: /* Calculate the determinant of the inverse Jacobian of the map to the standard element
275: which must be a constant for linear elements */
276: x11 = nodes[elements[elem*numCorners]*2];
277: y11 = nodes[elements[elem*numCorners]*2+1];
278: if (mesh->isPeriodic == PETSC_TRUE) {
279: x21 = MeshPeriodicDiffX(mesh, nodes[elements[elem*numCorners+1]*2] - x11);
280: x31 = MeshPeriodicDiffX(mesh, nodes[elements[elem*numCorners+2]*2] - x11);
281: y21 = MeshPeriodicDiffY(mesh, nodes[elements[elem*numCorners+1]*2+1] - y11);
282: y31 = MeshPeriodicDiffY(mesh, nodes[elements[elem*numCorners+2]*2+1] - y11);
283: } else {
284: x21 = nodes[elements[elem*numCorners+1]*2] - x11;
285: x31 = nodes[elements[elem*numCorners+2]*2] - x11;
286: y21 = nodes[elements[elem*numCorners+1]*2+1] - y11;
287: y31 = nodes[elements[elem*numCorners+2]*2+1] - y11;
288: }
289: jac = PetscAbsReal(x21*y31 - x31*y21);
290: if (jac < 1.0e-14) {
291: PetscPrintf(PETSC_COMM_SELF, "[%d]elem: %d x21: %g y21: %g x31: %g y31: %gn", rank, elem, x21, y21, x31, y31);
292: SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
293: }
294: #ifdef PETSC_USE_BOPT_g
295: PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
296: if (opt == PETSC_TRUE) {
297: PetscPrintf(PETSC_COMM_SELF, "[%d]elem: %d x21: %g y21: %g x31: %g y31: %g jac: %gn",
298: rank, elem, x21, y21, x31, y31, jac);
299: }
300: #endif
302: /* These are the elements of the inverse matrix */
303: invjac = 1/jac;
304: dxix = y31*invjac;
305: dxiy = -x31*invjac;
306: detx = -y21*invjac;
307: dety = x21*invjac;
309: /* Calculate element vector entries by Gaussian quadrature */
310: for(p = 0; p < numQuadPoints; p++)
311: {
312: x = MeshPeriodicX(mesh, x21*quadPoints[p*2] + x31*quadPoints[p*2+1] + x11);
313: y = MeshPeriodicY(mesh, y21*quadPoints[p*2] + y31*quadPoints[p*2+1] + y11);
314: /* Can this be simplified? */
315: for(arg = 0; arg < numArgs; arg++) {
316: for(j = 0; j < comp*3; j++)
317: fieldVal[arg][j] = 0.0;
318: for(func = 0; func < funcs; func++)
319: for(j = 0; j < comp; j++) {
320: fieldVal[arg][j*3] += field[arg][func*comp+j]*quadShapeFuncs[p*funcs+func];
321: fieldVal[arg][j*3+1] += field[arg][func*comp+j]*quadShapeFuncDers[p*funcs*2+func*2];
322: fieldVal[arg][j*3+2] += field[arg][func*comp+j]*quadShapeFuncDers[p*funcs*2+func*2+1];
323: }
324: }
326: /* Convert the field derivatives to old coordinates */
327: for(arg = 0; arg < numArgs; arg++) {
328: for(j = 0; j < comp; j++) {
329: dfxi = fieldVal[arg][j*3+1];
330: dfet = fieldVal[arg][j*3+2];
331: fieldVal[arg][j*3+1] = dfxi*dxix + dfet*detx;
332: fieldVal[arg][j*3+2] = dfxi*dxiy + dfet*dety;
333: }
334: }
336: (*f)(1, comp, &x, &y, PETSC_NULL, numArgs, fieldVal, funcVal, ctx);
337: #ifdef PETSC_USE_BOPT_g
338: PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
339: if (opt == PETSC_TRUE) {
340: PetscPrintf(PETSC_COMM_SELF, "[%d]p:%d jac: %g", rank, p, jac);
341: for(j = 0; j < comp; j++)
342: PetscPrintf(PETSC_COMM_SELF, " func[%d]: %g", j, PetscRealPart(funcVal[j]));
343: PetscPrintf(PETSC_COMM_SELF, "n");
344: }
345: #endif
347: for(i = 0, k = 0; i < funcs; i++) {
348: for(j = 0; j < comp; j++, k++) {
349: vec[k] += alpha*funcVal[j]*quadShapeFuncs[p*funcs+i]*jac*quadWeights[p];
350: #ifdef PETSC_USE_BOPT_g
351: PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
352: if (opt == PETSC_TRUE) {
353: PetscPrintf(PETSC_COMM_SELF, "[%d] vec[%d]: %gn", rank, k, PetscRealPart(vec[k]));
354: }
355: #endif
356: }
357: }
358: }
359: PetscLogFlops(12 + (8 + (6*numArgs + 5)*funcs*comp + 6*numArgs*comp) * numQuadPoints);
360: return(0);
361: }
363: static int DiscEvaluateALEOperatorGalerkin_Triangular_2D_Linear(Discretization disc, Mesh mesh, int elemSize,
364: int rowStart, int colStart, int op, PetscScalar alpha,
365: int elem, PetscScalar *field, PetscScalar *ALEfield, PetscScalar *mat,
366: void *ctx)
367: {
368: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
369: double *nodes = tri->nodes; /* The node coordinates */
370: int *elements = tri->faces; /* The element corners */
371: int numCorners = mesh->numCorners; /* The number of corners per element */
372: Operator oper = disc->operators[op]; /* The operator to discretize */
373: Discretization test = oper->test; /* The space of test functions */
374: ALEOperatorFunction opFunc = oper->ALEOpFunc; /* Integrals of operators which depend on J */
375: int rowSize = test->size; /* The number of shape functions per element */
376: int colSize = disc->size; /* The number of test functions per element */
377: double coords[MAX_CORNERS*2]; /* Coordinates of the element */
378: int f;
379: int ierr;
382: #ifdef PETSC_USE_BOPT_g
383: /* Check for valid operator */
384: if ((op < 0) || (op >= disc->numOps) || (!disc->operators[op])) SETERRQ(PETSC_ERR_ARG_WRONG, "Invalid operator");
385: #endif
387: if (opFunc == PETSC_NULL) SETERRQ(PETSC_ERR_ARG_CORRUPT, "Invalid function");
388: if (mesh->isPeriodic == PETSC_TRUE) {
389: coords[0*2+0] = nodes[elements[elem*numCorners+0]*2+0];
390: coords[0*2+1] = nodes[elements[elem*numCorners+0]*2+1];
391: for(f = 1; f < PetscMax(disc->funcs, test->funcs); f++) {
392: coords[f*2+0] = MeshPeriodicRelativeX(mesh, nodes[elements[elem*numCorners+f]*2+0], coords[0*2+0]);
393: coords[f*2+1] = MeshPeriodicRelativeY(mesh, nodes[elements[elem*numCorners+f]*2+1], coords[0*2+1]);
394: }
395: } else {
396: for(f = 0; f < PetscMax(disc->funcs, test->funcs); f++) {
397: coords[f*2+0] = nodes[elements[elem*numCorners+f]*2+0];
398: coords[f*2+1] = nodes[elements[elem*numCorners+f]*2+1];
399: }
400: }
402: (*opFunc)(disc, test, rowSize, colSize, rowStart, colStart, elemSize, coords, alpha, field, ALEfield, mat, ctx);
403:
404: return(0);
405: }
407: static int DiscEvaluateNonlinearALEOperatorGalerkin_Triangular_2D_Linear(Discretization disc, Mesh mesh, NonlinearOperator f,
408: PetscScalar alpha, int elem, int numArgs, PetscScalar **field,
409: PetscScalar *ALEfield, PetscScalar *vec, void *ctx)
410: {
411: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
412: double *nodes = tri->nodes; /* The node coordinates */
413: int *elements = tri->faces; /* The element corners */
414: int numCorners = mesh->numCorners; /* The number of corners per element */
415: int comp = disc->comp; /* The number of components in this field */
416: int funcs = disc->funcs; /* The number of shape functions per component */
417: PetscScalar *funcVal = disc->funcVal; /* Function value at a quadrature point */
418: PetscScalar **fieldVal = disc->fieldVal; /* Field value and derivatives at a quadrature point */
419: double jac; /* |J| for map to standard element */
420: double invjac; /* |J^{-1}| for map from standard element */
421: int numQuadPoints; /* Number of points used for Gaussian quadrature */
422: double *quadPoints; /* Points in the standard element for Gaussian quadrature */
423: double *quadWeights; /* Weights in the standard element for Gaussian quadrature */
424: double *quadShapeFuncs; /* Shape function evaluated at quadrature points */
425: double *quadShapeFuncDers; /* Shape function derivatives evaluated at quadrature points */
426: double x, y; /* The integration point */
427: double dxix; /* PartDer{xi}{x} */
428: double detx; /* PartDer{eta}{x} */
429: double dxiy; /* PartDer{xi}{y} */
430: double dety; /* PartDer{eta}{y} */
431: PetscScalar dfxi; /* PartDer{field}{xi} */
432: PetscScalar dfet; /* PartDer{field}{eta} */
433: double x11, y11, x21, y21, x31, y31;
434: int rank;
435: int i, j, k, p, func, arg;
436: #ifdef PETSC_USE_BOPT_g
437: PetscTruth opt;
438: #endif
439: int ierr;
442: if (numArgs > 2) SETERRQ(PETSC_ERR_SUP, "Only configured to handle two nonlinear arguments");
443: MPI_Comm_rank(disc->comm, &rank);
445: numQuadPoints = disc->numQuadPoints;
446: quadPoints = disc->quadPoints;
447: quadWeights = disc->quadWeights;
448: quadShapeFuncs = disc->quadShapeFuncs;
449: quadShapeFuncDers = disc->quadShapeFuncDers;
450:
451: /* Calculate the determinant of the inverse Jacobian of the map to the standard element
452: which must be a constant for linear elements */
453: x11 = nodes[elements[elem*numCorners]*2];
454: y11 = nodes[elements[elem*numCorners]*2+1];
455: if (mesh->isPeriodic == PETSC_TRUE) {
456: x21 = MeshPeriodicDiffX(mesh, nodes[elements[elem*numCorners+1]*2] - x11);
457: x31 = MeshPeriodicDiffX(mesh, nodes[elements[elem*numCorners+2]*2] - x11);
458: y21 = MeshPeriodicDiffY(mesh, nodes[elements[elem*numCorners+1]*2+1] - y11);
459: y31 = MeshPeriodicDiffY(mesh, nodes[elements[elem*numCorners+2]*2+1] - y11);
460: } else {
461: x21 = nodes[elements[elem*numCorners+1]*2] - x11;
462: x31 = nodes[elements[elem*numCorners+2]*2] - x11;
463: y21 = nodes[elements[elem*numCorners+1]*2+1] - y11;
464: y31 = nodes[elements[elem*numCorners+2]*2+1] - y11;
465: }
466: jac = PetscAbsReal(x21*y31 - x31*y21);
467: if (jac < 1.0e-14) {
468: PetscPrintf(PETSC_COMM_SELF, "[%d]elem: %d x21: %g y21: %g x31: %g y31: %gn", rank, elem, x21, y21, x31, y31);
469: SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
470: }
471: #ifdef PETSC_USE_BOPT_g
472: PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
473: if (opt == PETSC_TRUE) {
474: PetscPrintf(PETSC_COMM_SELF, "[%d]elem: %d x21: %g y21: %g x31: %g y31: %g jac: %gn",
475: rank, elem, x21, y21, x31, y31, jac);
476: }
477: #endif
479: /* These are the elements of the inverse matrix */
480: invjac = 1/jac;
481: dxix = y31*invjac;
482: dxiy = -x31*invjac;
483: detx = -y21*invjac;
484: dety = x21*invjac;
486: /* Calculate element vector entries by Gaussian quadrature */
487: for(p = 0; p < numQuadPoints; p++) {
488: x = MeshPeriodicX(mesh, x21*quadPoints[p*2] + x31*quadPoints[p*2+1] + x11);
489: y = MeshPeriodicY(mesh, y21*quadPoints[p*2] + y31*quadPoints[p*2+1] + y11);
490: /* Can this be simplified? */
491: for(arg = 0; arg < numArgs; arg++) {
492: for(j = 0; j < comp*3; j++) fieldVal[arg][j] = 0.0;
493: for(func = 0; func < funcs; func++)
494: for(j = 0; j < comp; j++) {
495: fieldVal[arg][j*3] += (field[arg][func*comp+j] - ALEfield[func*comp+j])*quadShapeFuncs[p*funcs+func];
496: fieldVal[arg][j*3+1] += field[arg][func*comp+j]*quadShapeFuncDers[p*funcs*2+func*2];
497: fieldVal[arg][j*3+2] += field[arg][func*comp+j]*quadShapeFuncDers[p*funcs*2+func*2+1];
498: }
499: }
501: /* Convert the field derivatives to old coordinates */
502: for(arg = 0; arg < numArgs; arg++) {
503: for(j = 0; j < comp; j++) {
504: dfxi = fieldVal[arg][j*3+1];
505: dfet = fieldVal[arg][j*3+2];
506: fieldVal[arg][j*3+1] = dfxi*dxix + dfet*detx;
507: fieldVal[arg][j*3+2] = dfxi*dxiy + dfet*dety;
508: }
509: }
511: (*f)(1, comp, &x, &y, PETSC_NULL, numArgs, fieldVal, funcVal, ctx);
512: #ifdef PETSC_USE_BOPT_g
513: PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
514: if (opt == PETSC_TRUE) {
515: PetscPrintf(PETSC_COMM_SELF, "[%d]p:%d jac: %g", rank, p, jac);
516: for(j = 0; j < comp; j++)
517: PetscPrintf(PETSC_COMM_SELF, " func[%d]: %g", j, PetscRealPart(funcVal[j]));
518: PetscPrintf(PETSC_COMM_SELF, "n");
519: }
520: #endif
522: for(i = 0, k = 0; i < funcs; i++) {
523: for(j = 0; j < comp; j++, k++) {
524: vec[k] += alpha*funcVal[j]*quadShapeFuncs[p*funcs+i]*jac*quadWeights[p];
525: #ifdef PETSC_USE_BOPT_g
526: PetscOptionsHasName(PETSC_NULL, "-trace_assembly", &opt);
527: if (opt == PETSC_TRUE) {
528: PetscPrintf(PETSC_COMM_SELF, "[%d] vec[%d]: %gn", rank, k, PetscRealPart(vec[k]));
529: }
530: #endif
531: }
532: }
533: }
534: PetscLogFlops(12 + (8 + (7*numArgs + 5)*funcs*comp + 6*numArgs*comp) * numQuadPoints);
535: return(0);
536: }
538: int Laplacian_Triangular_2D_Linear(Discretization disc, Discretization test, int rowSize, int colSize,
539: int globalRowStart, int globalColStart, int globalSize, double *coords,
540: PetscScalar alpha, PetscScalar *field, PetscScalar *array, void *ctx)
541: {
542: double x21, x31, y21, y31; /* Coordinates of the element, with point 1 at the origin */
543: double jac; /* |J| for map to standard element */
544: PetscScalar w; /* 1/(2 jac) */
545: int comp; /* Number of components */
546: int i;
549: /* Calculate the determinant of the inverse Jacobian of the map to the standard element
550: which must be a constant for linear elements - 1/|x_{21} y_{31} - x_{31} y_{21}| */
551: x21 = coords[2] - coords[0];
552: x31 = coords[4] - coords[0];
553: y21 = coords[3] - coords[1];
554: y31 = coords[5] - coords[1];
555: jac = PetscAbsReal(x21*y31 - x31*y21);
556: #ifdef PETSC_USE_BOPT_g
557: if (jac < 1.0e-14) {
558: PetscPrintf(PETSC_COMM_SELF, "x21: %g y21: %g x31: %g y31: %g jac: %gn", x21, y21, x31, y31, jac);
559: SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
560: }
561: #endif
562: /* PetscPrintf(PETSC_COMM_SELF, "x21: %g y21: %g x31: %g y31: %gn", x21, y21, x31, y31, jac); */
564: comp = rowSize/3;
565: w = 1.0/(2.0*jac);
566: w *= alpha;
567: for(i = 0; i < comp; i++) {
568: /* phi^1 phi^1 */
569: array[(0*comp+i+globalRowStart)*globalSize + 0*comp+i+globalColStart] =
570: (-x21*x21 + 2.0*x21*x31 - x31*x31 - (y21 - y31)*(y21 - y31))*w;
571: /* phi^1 phi^2 */
572: array[(0*comp+i+globalRowStart)*globalSize + 1*comp+i+globalColStart] = (-x21*x31 + x31*x31 + y31*(y31 - y21))*w;
573: /* phi^1 phi^3 */
574: array[(0*comp+i+globalRowStart)*globalSize + 2*comp+i+globalColStart] = (x21*x21 - x21*x31 + y21*(y21 - y31))*w;
575: /* phi^2 phi^1 */
576: array[(1*comp+i+globalRowStart)*globalSize + 0*comp+i+globalColStart] =
577: array[(0*comp+i+globalRowStart)*globalSize + 1*comp+i+globalColStart];
578: /* phi^2 phi^2 */
579: array[(1*comp+i+globalRowStart)*globalSize + 1*comp+i+globalColStart] = (-x31*x31 - y31*y31)*w;
580: /* phi^2 phi^3 */
581: array[(1*comp+i+globalRowStart)*globalSize + 2*comp+i+globalColStart] = (x21*x31 + y21*y31)*w;
582: /* phi^3 phi^1 */
583: array[(2*comp+i+globalRowStart)*globalSize + 0*comp+i+globalColStart] =
584: array[(0*comp+i+globalRowStart)*globalSize + 2*comp+i+globalColStart];
585: /* phi^3 phi^2 */
586: array[(2*comp+i+globalRowStart)*globalSize + 1*comp+i+globalColStart] =
587: array[(1*comp+i+globalRowStart)*globalSize + 2*comp+i+globalColStart];
588: /* phi^3 phi^3 */
589: array[(2*comp+i+globalRowStart)*globalSize + 2*comp+i+globalColStart] = (-x21*x21 - y21*y21)*w;
590: }
591: PetscLogFlops(47);
593: return(0);
594: }
596: int Weighted_Laplacian_Triangular_2D_Linear(Discretization disc, Discretization test, int rowSize, int colSize,
597: int globalRowStart, int globalColStart, int globalSize, double *coords,
598: PetscScalar alpha, PetscScalar *field, PetscScalar *array, void *ctx)
599: {
600: double x21, x31, y21, y31; /* Coordinates of the element, with point 1 at the origin */
601: #ifdef PETSC_USE_BOPT_g
602: double jac; /* |J| for map to standard element */
603: #endif
604: PetscScalar w; /* 1/2 */
605: int comp; /* Number of components */
606: int i;
608: /* Each element is weighted by its Jacobian. This is supposed to make smaller elements "stiffer". */
610: /* Calculate the determinant of the inverse Jacobian of the map to the standard element
611: which must be a constant for linear elements - 1/|x_{21} y_{31} - x_{31} y_{21}| */
612: x21 = coords[2] - coords[0];
613: x31 = coords[4] - coords[0];
614: y21 = coords[3] - coords[1];
615: y31 = coords[5] - coords[1];
616: #ifdef PETSC_USE_BOPT_g
617: jac = PetscAbsReal(x21*y31 - x31*y21);
618: if (jac < 1.0e-14) {
619: PetscPrintf(PETSC_COMM_SELF, "x21: %g y21: %g x31: %g y31: %g jac: %gn", x21, y21, x31, y31, jac);
620: SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
621: }
622: #endif
624: comp = rowSize/3;
625: w = 1.0/(2.0);
626: w *= alpha;
627: for(i = 0; i < comp; i++)
628: {
629: /* phi^1 phi^1 */
630: array[(0*comp+i+globalRowStart)*globalSize + 0*comp+i+globalColStart] =
631: (-x21*x21 + 2.0*x21*x31 - x31*x31 - (y21 - y31)*(y21 - y31))*w;
632: /* phi^1 phi^2 */
633: array[(0*comp+i+globalRowStart)*globalSize + 1*comp+i+globalColStart] = (-x21*x31 + x31*x31 + y31*(y31 - y21))*w;
634: /* phi^1 phi^3 */
635: array[(0*comp+i+globalRowStart)*globalSize + 2*comp+i+globalColStart] = (x21*x21 - x21*x31 + y21*(y21 - y31))*w;
636: /* phi^2 phi^1 */
637: array[(1*comp+i+globalRowStart)*globalSize + 0*comp+i+globalColStart] =
638: array[(0*comp+i+globalRowStart)*globalSize + 1*comp+i+globalColStart];
639: /* phi^2 phi^2 */
640: array[(1*comp+i+globalRowStart)*globalSize + 1*comp+i+globalColStart] = (-x31*x31 - y31*y31)*w;
641: /* phi^2 phi^3 */
642: array[(1*comp+i+globalRowStart)*globalSize + 2*comp+i+globalColStart] = (x21*x31 + y21*y31)*w;
643: /* phi^3 phi^1 */
644: array[(2*comp+i+globalRowStart)*globalSize + 0*comp+i+globalColStart] =
645: array[(0*comp+i+globalRowStart)*globalSize + 2*comp+i+globalColStart];
646: /* phi^3 phi^2 */
647: array[(2*comp+i+globalRowStart)*globalSize + 1*comp+i+globalColStart] =
648: array[(1*comp+i+globalRowStart)*globalSize + 2*comp+i+globalColStart];
649: /* phi^3 phi^3 */
650: array[(2*comp+i+globalRowStart)*globalSize + 2*comp+i+globalColStart] = (-x21*x21 - y21*y21)*w;
651: }
652: PetscLogFlops(47);
654: return(0);
655: }
657: int Gradient_Triangular_2D_Linear(Discretization disc, Discretization test, int rowSize, int colSize,
658: int globalRowStart, int globalColStart, int globalSize, double *coords,
659: PetscScalar alpha, PetscScalar *field, PetscScalar *array, void *ctx)
660: {
661: /* We are using the convention that
663: nabla matrix{v_1 cr v_2 cr vdots cr v_n} =
664: matrix{v^{(1)}_1 cr vdots cr v^{(d)}_1 cr v^{(1)}_2 cr vdots cr v^{(d)}_n}
666: and
668: nabla cdot matrix{v^{(1)}_1 cr vdots cr v^{(d)}_1 cr v^{(1)}_2 cr vdots cr v^{(d)}_n} =
669: matrix{v_1 cr v_2 cr vdots cr v_n}
671: where $d$ is the number of space dimensions. This agrees with the convention which allows
672: $Delta matrix{u_1 cr u_2} = 0$ to denote a set of scalar equations. This also means that
673: the dimension of the test function vector must be divisible by the number of space dimensions */
674: int numQuadPoints; /* Number of points used for Gaussian quadrature */
675: double *quadWeights; /* Weights in the standard element for Gaussian quadrature */
676: double *quadShapeFuncs; /* Shape functions evaluated at quadrature points */
677: double *quadTestFuncDers; /* Test function derivatives evaluated at quadrature points */
678: double dxxi; /* PartDer{x}{xi} */
679: double dxet; /* PartDer{x}{eta} */
680: double dyxi; /* PartDer{y}{xi} */
681: double dyet; /* PartDer{y}{eta} */
682: double dxix; /* PartDer{xi}{x} */
683: double detx; /* PartDer{eta}{x} */
684: double dxiy; /* PartDer{xi}{y} */
685: double dety; /* PartDer{eta}{y} */
686: double dphix; /* PartDer{phi_i}{x} times PartDer{phi_j}{x} */
687: double dphiy; /* PartDer{phi_i}{y} times PartDer{phi_j}{y} */
688: double jac; /* |J| for map to standard element */
689: double invjac; /* |J^{-1}| for map from standard element */
690: int comp; /* The number of field components */
691: int tcomp; /* The number of field components for the test field */
692: int funcs; /* The number of shape functions */
693: int tfuncs; /* The number of test functions */
694: int i, j, c, tc, f, p;
697: /* Calculate element matrix entries by Gaussian quadrature --
698: Since we integrate by parts here, the test and shape functions are switched */
699: comp = disc->comp;
700: tcomp = test->comp;
701: funcs = disc->funcs;
702: tfuncs = test->funcs;
703: numQuadPoints = disc->numQuadPoints;
704: quadWeights = disc->quadWeights;
705: quadShapeFuncs = disc->quadShapeFuncs;
706: quadTestFuncDers = test->quadShapeFuncDers;
707: for(p = 0; p < numQuadPoints; p++) {
708: /* PartDer{x}{xi}(p) = sum^{funcs}_{f=1} x_f PartDer{phi^f(p)}{xi}
709: PartDer{x}{eta}(p) = sum^{funcs}_{f=1} x_f PartDer{phi^f(p)}{eta}
710: PartDer{y}{xi}(p) = sum^{funcs}_{f=1} y_f PartDer{phi^f(p)}{xi}
711: PartDer{y}{eta}(p) = sum^{funcs}_{f=1} y_f PartDer{phi^f(p)}{eta} */
712: dxxi = 0.0; dxet = 0.0;
713: dyxi = 0.0; dyet = 0.0;
714: for(f = 0; f < tfuncs; f++) {
715: dxxi += coords[f*2] *quadTestFuncDers[p*tfuncs*2+f*2];
716: dxet += coords[f*2] *quadTestFuncDers[p*tfuncs*2+f*2+1];
717: dyxi += coords[f*2+1]*quadTestFuncDers[p*tfuncs*2+f*2];
718: dyet += coords[f*2+1]*quadTestFuncDers[p*tfuncs*2+f*2+1];
719: }
720: jac = PetscAbsReal(dxxi*dyet - dxet*dyxi);
721: #ifdef PETSC_USE_BOPT_g
722: if (jac < 1.0e-14) {
723: PetscPrintf(PETSC_COMM_SELF, "p: %d x1: %g y1: %g x2: %g y2: %g x3: %g y3: %gn",
724: p, coords[0], coords[1], coords[2], coords[3], coords[4], coords[5]);
725: SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
726: }
727: #endif
728: /* PetscPrintf(PETSC_COMM_SELF, "p: %d x1: %g y1: %g x2: %g y2: %g x3: %g y3: %gn",
729: p, coords[0], coords[1], coords[2], coords[3], coords[4], coords[5]); */
730: /* These are the elements of the inverse matrix */
731: invjac = 1.0/jac;
732: dxix = dyet*invjac;
733: dxiy = -dxet*invjac;
734: detx = -dyxi*invjac;
735: dety = dxxi*invjac;
736: /* PetscPrintf(PETSC_COMM_SELF, "dxix: %g dxiy: %g detx: %g dety: %g jac: %gn", dxix, dxiy, detx, dety, jac); */
738: /* The rows are test functions */
739: for(i = 0; i < tfuncs; i++) {
740: /* We divide by the space dimension */
741: for(tc = 0; tc < tcomp/2; tc++) {
742: /* The columns are shape functions */
743: for(j = 0; j < funcs; j++) {
744: for(c = 0; c < comp; c++) {
745: dphix = quadTestFuncDers[p*tfuncs*2+i*2]*dxix + quadTestFuncDers[p*tfuncs*2+i*2+1]*detx;
746: dphiy = quadTestFuncDers[p*tfuncs*2+i*2]*dxiy + quadTestFuncDers[p*tfuncs*2+i*2+1]*dety;
747: array[(i*tcomp+tc*2+globalRowStart)*globalSize + j*comp+c+globalColStart] +=
748: -alpha*dphix*quadShapeFuncs[p*funcs+j]*jac*quadWeights[p];
749: array[(i*tcomp+tc*2+1+globalRowStart)*globalSize + j*comp+c+globalColStart] +=
750: -alpha*dphiy*quadShapeFuncs[p*funcs+j]*jac*quadWeights[p];
751: }
752: }
753: }
754: }
755: }
756: PetscLogFlops((8*tfuncs + 8 + 8*tfuncs*tcomp*funcs*comp) * numQuadPoints);
758: return(0);
759: }
761: int DiscInterpolateField_Triangular_2D_Linear(Discretization disc, Mesh oldMesh, int elem, double x, double y, double z,
762: PetscScalar *oldFieldVal, PetscScalar *newFieldVal, InterpolationType type)
763: {
764: Mesh_Triangular *tri = (Mesh_Triangular *) oldMesh->data;
765: Partition p = oldMesh->part;
766: int numCorners = oldMesh->numCorners;
767: int *elements = tri->faces;
768: int *neighbors = tri->neighbors;
769: double *nodes = tri->nodes;
770: double coords[12]; /* Coordinates of our "big element" */
771: double x11, y11; /* Coordinates of vertex 0 */
772: double xi, eta; /* Canonical coordinates of the interpolation point */
773: double dxix; /* PartDer{xi}{x} */
774: double detx; /* PartDer{eta}{x} */
775: double dxiy; /* PartDer{xi}{y} */
776: double dety; /* PartDer{eta}{y} */
777: double dxxi; /* PartDer{x}{xi} */
778: double dxet; /* PartDer{x}{eta} */
779: double dyxi; /* PartDer{y}{xi} */
780: double dyet; /* PartDer{y}{eta} */
781: double jac, invjac; /* The Jacobian determinant and its inverse */
782: int comp = disc->comp;
783: int neighbor, corner, nelem, node, c;
784: #ifdef PETSC_USE_BOPT_g
785: PetscTruth opt;
786: #endif
787: int ierr;
790: /* No scheme in place for boundary elements */
791: for(neighbor = 0; neighbor < 3; neighbor++)
792: if (neighbors[elem*3+neighbor] < 0)
793: type = INTERPOLATION_LOCAL;
795: switch (type)
796: {
797: case INTERPOLATION_LOCAL:
798: x11 = nodes[elements[elem*numCorners]*2];
799: y11 = nodes[elements[elem*numCorners]*2+1];
800: if (oldMesh->isPeriodic == PETSC_TRUE) {
801: dxxi = MeshPeriodicDiffX(oldMesh, nodes[elements[elem*numCorners+1]*2] - x11);
802: dxet = MeshPeriodicDiffX(oldMesh, nodes[elements[elem*numCorners+2]*2] - x11);
803: dyxi = MeshPeriodicDiffY(oldMesh, nodes[elements[elem*numCorners+1]*2+1] - y11);
804: dyet = MeshPeriodicDiffY(oldMesh, nodes[elements[elem*numCorners+2]*2+1] - y11);
805: } else {
806: dxxi = nodes[elements[elem*numCorners+1]*2] - x11;
807: dxet = nodes[elements[elem*numCorners+2]*2] - x11;
808: dyxi = nodes[elements[elem*numCorners+1]*2+1] - y11;
809: dyet = nodes[elements[elem*numCorners+2]*2+1] - y11;
810: }
811: jac = PetscAbsReal(dxxi*dyet - dxet*dyxi);
812: if (jac < 1.0e-14) {
813: PetscPrintf(PETSC_COMM_SELF, "[%d]elem: %d x21: %g y21: %g x31: %g y31: %gn", p->rank, elem, dxxi, dyxi, dxet, dyet);
814: SETERRQ(PETSC_ERR_DISC_SING_JAC, "Singular Jacobian");
815: }
816: #ifdef PETSC_USE_BOPT_g
817: PetscOptionsHasName(PETSC_NULL, "-trace_interpolation", &opt);
818: if (opt == PETSC_TRUE) {
819: PetscPrintf(PETSC_COMM_SELF, "[%d]elem: %d x21: %g y21: %g x31: %g y31: %g jac: %gn",
820: p->rank, elem, dxxi, dyxi, dxet, dyet, jac);
821: }
822: #endif
824: /* These are the elements of the inverse matrix */
825: invjac = 1/jac;
826: dxix = dyet*invjac;
827: dxiy = -dxet*invjac;
828: detx = -dyxi*invjac;
829: dety = dxxi*invjac;
830: if (oldMesh->isPeriodic == PETSC_TRUE) {
831: xi = dxix*MeshPeriodicDiffX(oldMesh, x - x11) + dxiy*MeshPeriodicDiffY(oldMesh, y - y11);
832: eta = detx*MeshPeriodicDiffX(oldMesh, x - x11) + dety*MeshPeriodicDiffY(oldMesh, y - y11);
833: } else {
834: xi = dxix*(x - x11) + dxiy*(y - y11);
835: eta = detx*(x - x11) + dety*(y - y11);
836: }
837: for(c = 0 ; c < comp; c++) {
838: newFieldVal[c] = oldFieldVal[0*comp+c]*(1.0 - xi - eta) + oldFieldVal[1*comp+c]*xi + oldFieldVal[2*comp+c]*eta;
839: }
840: PetscLogFlops(7+15+7*comp);
841: break;
842: case INTERPOLATION_HALO:
843: /* Here is our "big element" where numbers in parantheses represent
844: the numbering on the old little element:
846: 2
847: |
848: |
849: |
850: (1) 4---3 (0)
851: | |
852: | |
853: | |
854: 0---5---1
855: (2)
857: We search for the neighbor node by looking for the vertex not a member of the original element.
858: */
859: for(neighbor = 0; neighbor < 3; neighbor++)
860: {
861: nelem = neighbors[elem*3+neighbor];
862: for(corner = 0; corner < 3; corner++)
863: {
864: node = elements[nelem*numCorners+corner];
865: if ((node != elements[elem*numCorners+((neighbor+1)%3)]) && (node != elements[elem*numCorners+((neighbor+2)%3)]))
866: {
867: /* The neighboring elements give the vertices */
868: coords[neighbor*2] = MeshPeriodicRelativeX(oldMesh, nodes[node*2+0], x);
869: coords[neighbor*2+1] = MeshPeriodicRelativeY(oldMesh, nodes[node*2+1], y);
870: break;
871: }
872: }
873: }
874: /* Element vertices form midnodes */
875: coords[3*2] = MeshPeriodicRelativeX(oldMesh, nodes[elements[elem*numCorners+0]*2+0], x);
876: coords[3*2+1] = MeshPeriodicRelativeY(oldMesh, nodes[elements[elem*numCorners+0]*2+1], y);
877: coords[4*2] = MeshPeriodicRelativeX(oldMesh, nodes[elements[elem*numCorners+1]*2+0], x);
878: coords[4*2+1] = MeshPeriodicRelativeY(oldMesh, nodes[elements[elem*numCorners+1]*2+1], y);
879: coords[5*2] = MeshPeriodicRelativeX(oldMesh, nodes[elements[elem*numCorners+2]*2+0], x);
880: coords[5*2+1] = MeshPeriodicRelativeY(oldMesh, nodes[elements[elem*numCorners+2]*2+1], y);
881: /* Get the (xi,eta) coordinates of the point */
882: DiscTransformCoords_Triangular_2D_Quadratic(x, y, coords, &xi, &eta);
883: if ((xi < -1.0e-02) || (eta < -1.0e-02)) {
884: PetscPrintf(PETSC_COMM_SELF, "Linear: elem: %d x: %g y: %g xi: %g eta: %gn", elem, x, y, xi, eta);
885: SETERRQ(PETSC_ERR_PLIB, "Standard element coordinates were negative");
886: }
887: for(c = 0 ; c < comp; c++) {
888: newFieldVal[c] = oldFieldVal[0*comp+c]*(1.0 - xi - eta)*(1.0 - 2.0*xi - 2.0*eta) +
889: oldFieldVal[1*comp+c]*xi *(2.0*xi - 1.0) +
890: oldFieldVal[2*comp+c]*eta*(2.0*eta - 1.0) +
891: oldFieldVal[3*comp+c]*4.0*xi*eta +
892: oldFieldVal[4*comp+c]*4.0*eta*(1.0 - xi - eta) +
893: oldFieldVal[5*comp+c]*4.0*xi *(1.0 - xi - eta);
894: }
895: PetscLogFlops(34*comp);
896: break;
897: default:
898: SETERRQ1(PETSC_ERR_ARG_WRONG, "Unknown interpolation type %d", type);
899: }
900:
901: return(0);
902: }
904: int DiscInterpolateElementVec_Triangular_2D_Linear(Discretization disc, ElementVec vec, Discretization newDisc, ElementVec newVec)
905: {
906: int funcs = disc->funcs;
907: int comp = disc->comp;
908: int size = disc->size;
909: PetscScalar *array, *newArray;
910: PetscTruth islin, isquad;
911: int f, c;
912: int ierr;
915: ElementVecGetArray(vec, &array);
916: ElementVecGetArray(newVec, &newArray);
917: PetscTypeCompare((PetscObject) newDisc, DISCRETIZATION_TRIANGULAR_2D_LINEAR, &islin);
918: PetscTypeCompare((PetscObject) newDisc, DISCRETIZATION_TRIANGULAR_2D_QUADRATIC, &isquad);
919: if (islin == PETSC_TRUE) {
920: PetscMemcpy(newArray, array, size * sizeof(PetscScalar));
921: } else if (isquad == PETSC_TRUE) {
922: for(f = 0; f < newDisc->funcs; f++) {
923: for(c = 0; c < comp; c++) {
924: if (f < funcs) {
925: newArray[f*comp+c] = array[f*comp+c];
926: } else {
927: newArray[f*comp+c] = 0.5*(array[((f+1)%funcs)*comp+c] + array[((f+2)%funcs)*comp+c]);
928: }
929: }
930: }
931: } else {
932: SETERRQ(PETSC_ERR_SUP, "Discretization not supported");
933: }
934: ElementVecRestoreArray(vec, &array);
935: ElementVecRestoreArray(newVec, &newArray);
936: return(0);
937: }
939: /*
940: DiscSetupQuadrature_Triangular_2D_Linear - Setup Gaussian quadrature with a 7 point integration rule
942: Input Parameter:
943: . disc - The Discretization
944: */
945: int DiscSetupQuadrature_Triangular_2D_Linear(Discretization disc) {
946: int dim = disc->dim;
947: int funcs = disc->funcs;
948: int p;
952: disc->numQuadPoints = 7;
953: PetscMalloc(disc->numQuadPoints*dim * sizeof(double), &disc->quadPoints);
954: PetscMalloc(disc->numQuadPoints * sizeof(double), &disc->quadWeights);
955: PetscMalloc(disc->numQuadPoints*funcs * sizeof(double), &disc->quadShapeFuncs);
956: PetscMalloc(disc->numQuadPoints*funcs*dim * sizeof(double), &disc->quadShapeFuncDers);
957: PetscLogObjectMemory(disc, (disc->numQuadPoints*(funcs*(dim+1) + dim+1)) * sizeof(double));
958: disc->quadPoints[0] = 1.0/3.0;
959: disc->quadPoints[1] = disc->quadPoints[0];
960: disc->quadWeights[0] = 0.11250000000000;
961: disc->quadPoints[2] = 0.797426985353087;
962: disc->quadPoints[3] = 0.101286507323456;
963: disc->quadWeights[1] = 0.0629695902724135;
964: disc->quadPoints[4] = disc->quadPoints[3];
965: disc->quadPoints[5] = disc->quadPoints[2];
966: disc->quadWeights[2] = disc->quadWeights[1];
967: disc->quadPoints[6] = disc->quadPoints[4];
968: disc->quadPoints[7] = disc->quadPoints[3];
969: disc->quadWeights[3] = disc->quadWeights[1];
970: disc->quadPoints[8] = 0.470142064105115;
971: disc->quadPoints[9] = 0.059715871789770;
972: disc->quadWeights[4] = 0.066197076394253;
973: disc->quadPoints[10] = disc->quadPoints[8];
974: disc->quadPoints[11] = disc->quadPoints[8];
975: disc->quadWeights[5] = disc->quadWeights[4];
976: disc->quadPoints[12] = disc->quadPoints[9];
977: disc->quadPoints[13] = disc->quadPoints[8];
978: disc->quadWeights[6] = disc->quadWeights[4];
979: for(p = 0; p < disc->numQuadPoints; p++) {
980: /* phi^0: 1 - xi - eta */
981: disc->quadShapeFuncs[p*funcs] = 1.0 - disc->quadPoints[p*2] - disc->quadPoints[p*2+1];
982: disc->quadShapeFuncDers[p*funcs*2+0*2] = -1.0;
983: disc->quadShapeFuncDers[p*funcs*2+0*2+1] = -1.0;
984: /* phi^1: xi */
985: disc->quadShapeFuncs[p*funcs+1] = disc->quadPoints[p*2];
986: disc->quadShapeFuncDers[p*funcs*2+1*2] = 1.0;
987: disc->quadShapeFuncDers[p*funcs*2+1*2+1] = 0.0;
988: /* phi^2: eta */
989: disc->quadShapeFuncs[p*funcs+2] = disc->quadPoints[p*2+1];
990: disc->quadShapeFuncDers[p*funcs*2+2*2] = 0.0;
991: disc->quadShapeFuncDers[p*funcs*2+2*2+1] = 1.0;
992: }
993: return(0);
994: }
996: /*
997: DiscSetupOperators_Triangular_2D_Linear - Setup the default operators
999: Input Parameter:
1000: . disc - The Discretization
1001: */
1002: int DiscSetupOperators_Triangular_2D_Linear(Discretization disc) {
1003: int comp = disc->comp;
1004: int size = disc->size;
1005: PetscScalar *precompInt;
1006: int newOp;
1007: int c, i, j;
1008: int ierr;
1011: /* The Identity operator I -- the matrix is symmetric */
1012: PetscMalloc(size*size * sizeof(PetscScalar), &precompInt);
1013: PetscLogObjectMemory(disc, size*size * sizeof(PetscScalar));
1014: PetscMemzero(precompInt, size*size * sizeof(PetscScalar));
1015: for(c = 0; c < comp; c++) {
1016: precompInt[(0*comp+c)*size + 0*comp+c] = 1.0/12.0;
1017: precompInt[(0*comp+c)*size + 1*comp+c] = 1.0/24.0;
1018: precompInt[(0*comp+c)*size + 2*comp+c] = 1.0/24.0;
1019: precompInt[(1*comp+c)*size + 1*comp+c] = 1.0/12.0;
1020: precompInt[(1*comp+c)*size + 2*comp+c] = 1.0/24.0;
1021: precompInt[(2*comp+c)*size + 2*comp+c] = 1.0/12.0;
1022: }
1023: for(i = 0; i < size; i++) {
1024: for(j = 0; j < i; j++) {
1025: precompInt[i*size + j] = precompInt[j*size + i];
1026: }
1027: }
1028: DiscretizationRegisterPrecomputedOperator(disc, precompInt, &newOp);
1029: if (newOp != IDENTITY) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "Default operator %d not setup correctly", IDENTITY);
1030: /* The Laplacian operator Delta -- the matrix is symmetric */
1031: DiscretizationRegisterOperator(disc, Laplacian_Triangular_2D_Linear, &newOp);
1032: if (newOp != LAPLACIAN) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "Default operator %d not setup correctly", LAPLACIAN);
1033: /* The Gradient operator nabla -- the matrix is rectangular */
1034: DiscretizationRegisterOperator(disc, Gradient_Triangular_2D_Linear, &newOp);
1035: if (newOp != GRADIENT) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "Default operator %d not setup correctly", GRADIENT);
1036: /* The Divergence operator nablacdot -- the matrix is rectangular */
1037: DiscretizationRegisterOperator(disc, PETSC_NULL, &newOp);
1038: if (newOp != DIVERGENCE) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "Default operator %d not setup correctly", DIVERGENCE);
1039: /* The weighted Laplacian operator -- the matrix is symmetric */
1040: DiscretizationRegisterOperator(disc, Weighted_Laplacian_Triangular_2D_Linear, &newOp);
1041: if (newOp != WEIGHTED_LAP) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE, "Default operator %d not setup correctly", WEIGHTED_LAP);
1042: return(0);
1043: }
1045: static struct _DiscretizationOps DOps = {PETSC_NULL/* DiscretizationSetup */,
1046: DiscSetupOperators_Triangular_2D_Linear,
1047: PETSC_NULL/* DiscretizationSetFromOptions */,
1048: DiscView_Triangular_2D_Linear,
1049: DiscDestroy_Triangular_2D_Linear,
1050: DiscEvaluateFunctionGalerkin_Triangular_2D_Linear,
1051: DiscEvaluateOperatorGalerkin_Triangular_2D_Linear,
1052: DiscEvaluateALEOperatorGalerkin_Triangular_2D_Linear,
1053: DiscEvaluateNonlinearOperatorGalerkin_Triangular_2D_Linear,
1054: DiscEvaluateNonlinearALEOperatorGalerkin_Triangular_2D_Linear,
1055: DiscInterpolateField_Triangular_2D_Linear,
1056: DiscInterpolateElementVec_Triangular_2D_Linear};
1058: EXTERN_C_BEGIN
1059: int DiscCreate_Triangular_2D_Linear(Discretization disc) {
1060: int arg;
1064: if (disc->comp <= 0) {
1065: SETERRQ(PETSC_ERR_ARG_WRONG, "Discretization must have at least 1 component. Call DiscretizationSetNumComponents() to set this.");
1066: }
1067: PetscMemcpy(disc->ops, &DOps, sizeof(struct _DiscretizationOps));
1068: disc->dim = 2;
1069: disc->funcs = 3;
1070: disc->size = disc->funcs*disc->comp;
1072: DiscretizationSetupDefaultOperators(disc);
1073: DiscSetupQuadrature_Triangular_2D_Linear(disc);
1075: DiscretizationCreate(disc->comm, &disc->bdDisc);
1076: DiscretizationSetNumComponents(disc->bdDisc, disc->comp);
1077: DiscretizationSetType(disc->bdDisc, BD_DISCRETIZATION_TRIANGULAR_2D_LINEAR);
1079: /* Storage */
1080: PetscMalloc(disc->comp * sizeof(PetscScalar), &disc->funcVal);
1081: PetscMalloc(2 * sizeof(PetscScalar *), &disc->fieldVal);
1082: for(arg = 0; arg < 2; arg++) {
1083: PetscMalloc(disc->comp*(disc->dim+1) * sizeof(PetscScalar), &disc->fieldVal[arg]);
1084: }
1085: return(0);
1086: }
1087: EXTERN_C_END