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