Actual source code: mathematica.c
1: #ifdef PETSC_RCS_HEADER
2: static char vcid[] = "$Id: mathematica.c,v 1.9 2000/01/26 15:46:22 baggag Exp $";
3: #endif
5: /*
6: Written by Matt Knepley, knepley@cs.purdue.edu 7/23/97
7: Major overhall for interactivity 11/14/97
8: Reorganized 11/8/98
9: */
10: #include src/sys/src/viewer/viewerimpl.h
11: #include src/gvec/gvecimpl.h
12: #include "src/mesh/impls/triangular/2d/2dimpl.h"
13: #include src/sles/pc/pcimpl.h
14: #if 0
15: #include src/sles/pc/impls/ml/ml.h
16: #endif
17: #include src/mat/impls/aij/seq/aij.h
18: #include mathematica.h
20: PetscViewer VIEWER_MATHEMATICA_WORLD_PRIVATE = PETSC_NULL;
21: #ifdef PETSC_HAVE_MATHEMATICA
22: static void *mathematicaEnv = PETSC_NULL;
23: #endif
25: /*@C
26: PetscViewerMathematicaInitializePackage - This function initializes everything in the Petsc interface to Mathematica. It is
27: called from PetscDLLibraryRegister() when using dynamic libraries, and on the call to PetscInitialize()
28: when using static libraries.
30: Input Parameter:
31: path - The dynamic library path, or PETSC_NULL
33: Level: developer
35: .keywords: Petsc, initialize, package, PLAPACK
36: .seealso: PetscInitializePackage(), PetscInitialize()
37: @*/
38: int PetscViewerMathematicaInitializePackage(char *path) {
39: static PetscTruth initialized = PETSC_FALSE;
42: if (initialized == PETSC_TRUE) return(0);
43: initialized = PETSC_TRUE;
44: #ifdef PETSC_HAVE_MATHEMATICA
45: mathematicaEnv = (void *) MLInitialize(0);
46: #endif
47: return(0);
48: }
50: /*@C
51: PetscViewerMathematicaDestroyPackage - This function destroys everything in the Petsc interface to Mathematica. It is
52: called from PetscFinalize().
54: Level: developer
56: .keywords: Petsc, destroy, package, mathematica
57: .seealso: PetscFinalize()
58: @*/
59: int PetscViewerMathematicaFinalizePackage(void) {
61: #ifdef PETSC_HAVE_MATHEMATICA
62: if (mathematicaEnv != PETSC_NULL) MLDeinitialize((MLEnvironment) mathematicaEnv);
63: #endif
64: return(0);
65: }
67: int PetscViewerInitializeMathematicaWorld_Private()
68: {
72: if (VIEWER_MATHEMATICA_WORLD_PRIVATE) return(0);
73: PetscViewerMathematicaOpen(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_NULL, PETSC_NULL, &VIEWER_MATHEMATICA_WORLD_PRIVATE);
74:
75: return(0);
76: }
78: static int PetscViewerDestroy_Mathematica(PetscViewer viewer)
79: {
80: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
81: int ierr;
84: #ifdef PETSC_HAVE_MATHEMATICA
85: MLClose(vmath->link);
86: #endif
87: PetscFree(vmath);
88: return(0);
89: }
91: int PetscViewerDestroyMathematica_Private(void)
92: {
96: if (VIEWER_MATHEMATICA_WORLD_PRIVATE) {
97: PetscViewerDestroy(VIEWER_MATHEMATICA_WORLD_PRIVATE);
98: }
99: return(0);
100: }
102: EXTERN_C_BEGIN
103: int PetscViewerCreate_Mathematica(PetscViewer v)
104: {
105: PetscViewer_Mathematica *vmath;
106: char type[256];
107: PetscTruth isMotif, isPS, isPSFile, isTri, isVecTri, isVec, isSurface;
108: PetscTruth opt;
109: int ierr;
113: PetscNew(PetscViewer_Mathematica, &vmath);
114: v->data = (void *) vmath;
115: v->ops->destroy = PetscViewerDestroy_Mathematica;
116: v->ops->flush = 0;
117: PetscStrallocpy(PETSC_VIEWER_MATHEMATICA, &v->type_name);
119: vmath->graphicsType = GRAPHICS_MOTIF;
120: PetscOptionsGetString("viewer_", "-math_graphics", type, 255, &opt);
121: if (opt == PETSC_TRUE)
122: {
123: PetscStrcasecmp(type, "Motif", &isMotif);
124: PetscStrcasecmp(type, "PS", &isPS);
125: PetscStrcasecmp(type, "PSFile", &isPSFile);
126: if (isMotif == PETSC_TRUE)
127: vmath->graphicsType = GRAPHICS_MOTIF;
128: else if (isPS == PETSC_TRUE)
129: vmath->graphicsType = GRAPHICS_PS_STDOUT;
130: else if (isPSFile == PETSC_TRUE)
131: vmath->graphicsType = GRAPHICS_PS_FILE;
132: }
133: vmath->plotType = MATHEMATICA_TRIANGULATION_PLOT;
134: PetscOptionsGetString("viewer_", "-math_type", type, 255, &opt);
135: if (opt == PETSC_TRUE)
136: {
137: PetscStrcasecmp(type, "Triangulation", &isTri);
138: PetscStrcasecmp(type, "VectorTriangulation", &isVecTri);
139: PetscStrcasecmp(type, "Vector", &isVec);
140: PetscStrcasecmp(type, "Surface", &isSurface);
141: if (isTri == PETSC_TRUE)
142: vmath->plotType = MATHEMATICA_TRIANGULATION_PLOT;
143: else if (isVecTri == PETSC_TRUE)
144: vmath->plotType = MATHEMATICA_VECTOR_TRIANGULATION_PLOT;
145: else if (isVec == PETSC_TRUE)
146: vmath->plotType = MATHEMATICA_VECTOR_PLOT;
147: else if (isSurface == PETSC_TRUE)
148: vmath->plotType = MATHEMATICA_SURFACE_PLOT;
149: }
150: vmath->objName = PETSC_NULL;
152: return(0);
153: }
154: EXTERN_C_END
156: int PetscViewerMathematicaSetConnection(PetscViewer v, int port, const char machine[], const char mode[])
157: {
158: #ifdef PETSC_HAVE_MATHEMATICA
159: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) v->data;
160: int argc = 6;
161: #endif
162: char *argv[6];
163: int numPorts, numHosts, numModes;
164: int *ports;
165: char **hosts;
166: char **modes;
167: int size, rank;
168: int i;
169: PetscTruth opt;
170: int ierr;
173: MPI_Comm_size(v->comm, &size);
174: MPI_Comm_rank(v->comm, &rank);
176: /* Link name */
177: argv[0] = "-linkname";
178: ierr = PetscMalloc(15 * sizeof(char), &argv[1]);
179: if (port == PETSC_DECIDE) {
180: ierr = PetscMalloc(size * sizeof(int), &ports);
181: numPorts = size;
182: PetscOptionsGetIntArray("viewer_", "-math_port", ports, &numPorts, &opt);
183: if (opt == PETSC_TRUE) {
184: if (numPorts > rank)
185: sprintf(argv[1], "%6d", ports[rank]);
186: else
187: sprintf(argv[1], "%6d", ports[0]);
188: } else {
189: sprintf(argv[1], "%s", "math -mathlink");
190: }
191: PetscFree(ports);
192: } else {
193: sprintf(argv[1], "%6d", port);
194: }
196: /* Link host */
197: argv[2] = "-linkhost";
198: ierr = PetscMalloc(256 * sizeof(char), &argv[3]);
199: if (machine == PETSC_NULL) {
200: ierr = PetscMalloc(size * sizeof(char *), &hosts);
201: numHosts = size;
202: PetscOptionsGetStringArray("viewer_", "-math_host", hosts, &numHosts, &opt);
203: if (opt == PETSC_FALSE) {
204: PetscGetHostName(argv[3], 255);
205: } else {
206: if (numHosts > rank)
207: PetscStrncpy(argv[3], hosts[rank], 255);
208: else
209: PetscStrncpy(argv[3], hosts[0], 255);
210: for(i = 0; i < numHosts; i++) {
211: PetscFree(hosts[i]);
212: }
213: }
214: PetscFree(hosts);
215: } else {
216: PetscStrncpy(argv[3], machine, 255);
217: }
218: argv[3][255] = '0';
220: /* Link mode */
221: argv[4] = "-linkmode";
222: ierr = PetscMalloc(14 * sizeof(char), &argv[5]);
223: if (mode == PETSC_NULL) {
224: ierr = PetscMalloc(size * sizeof(char *), &modes);
225: numModes = size;
226: PetscOptionsGetStringArray("viewer_", "-math_mode", modes, &numModes, &opt);
227: if (opt == PETSC_FALSE) {
228: PetscStrncpy(argv[5], "Launch", 13);
229: } else {
230: if (numModes > rank)
231: PetscStrncpy(argv[5], modes[rank], 13);
232: else
233: PetscStrncpy(argv[5], modes[0], 13);
234: for(i = 0; i < numModes; i++) {
235: PetscFree(modes[i]);
236: }
237: }
238: PetscFree(modes);
239: } else {
240: PetscStrncpy(argv[5], mode, 13);
241: }
242: argv[5][13] = '0';
243: #ifdef PETSC_HAVE_MATHEMATICA
244: vmath->link = MLOpen(argc, argv);
245: #endif
246: PetscFree(argv[1]);
247: PetscFree(argv[3]);
248: PetscFree(argv[5]);
250: return(0);
251: }
253: /*----------------------------------------- Public Functions --------------------------------------------------------*/
254: /*@C
255: PetscViewerMathematicaOpen - Communicates with Mathemtica using MathLink.
257: Collective on comm
259: Input Parameters:
260: + comm - The MPI communicator
261: . port - [optional] The port to connect on, or PETSC_DECIDE
262: . machine - [optional] The machine to run Mathematica on, or PETSC_NULL
263: - mode - [optional] The connection mode, or PETSC_NULL
265: Output Parameter:
266: . viewer - The Mathematica viewer
268: Level: intermediate
270: Notes:
271: Most users should employ the following commands to access the
272: Mathematica viewers
273: $
274: $ PetscViewerMathematicaOpen(MPI_Comm comm, int port, char *machine, char *mode, PetscViewer &viewer)
275: $ MatView(Mat matrix, PetscViewer viewer)
276: $
277: $ or
278: $
279: $ PetscViewerMathematicaOpen(MPI_Comm comm, int port, char *machine, char *mode, PetscViewer &viewer)
280: $ VecView(Vec vector, PetscViewer viewer)
282: Options Database Keys:
283: $ -viewer_math_machine <machine> - The host machine for the kernel
284: $ -viewer_math_port <port> - The port for connection
285: $ -viewer_math_mode <mode> - The mode, e.g. Launch, Connect
286: $ -viewer_math_type <type> - The plot type, e.g. Triangulation, Vector
287: $ -viewer_math_graphics <output> - The output type, e.g. Motif, PS, PSFile
289: .keywords: PetscViewer, Mathematica, open
291: .seealso: MatView(), VecView()
292: @*/
293: int PetscViewerMathematicaOpen(MPI_Comm comm, int port, const char machine[], const char mode[], PetscViewer *v)
294: {
298: PetscViewerCreate(comm, v);
299: PetscViewerSetType(*v, PETSC_VIEWER_MATHEMATICA);
300: PetscViewerMathematicaSetConnection(*v, port, machine, mode);
301: return(0);
302: }
304: #ifdef PETSC_HAVE_MATHEMATICA
305: /*@C
306: PetscViewerMathematicaGetLink - Returns the link to Mathematica
308: Input Parameters:
309: . viewer - The Mathematica viewer
310: . link - The link to Mathematica
312: Level: intermediate
314: .keywords PetscViewer, Mathematica, link
315: .seealso PetscViewerMathematicaOpen()
316: @*/
317: int PetscViewerMathematicaGetLink(PetscViewer viewer, MLINK *link)
318: {
319: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
323: *link = vmath->link;
324: return(0);
325: }
326: #endif
328: /*@C
329: PetscViewerMathematicaSkipPackets - Discard packets sent by Mathematica until a certain packet type is received
331: Input Parameters:
332: . viewer - The Mathematica viewer
333: . type - The packet type to search for, e.g RETURNPKT
335: Level: advanced
337: .keywords PetscViewer, Mathematica, packets
338: .seealso PetscViewerMathematicaSetName(), PetscViewerMathematicaGetVector()
339: @*/
340: int PetscViewerMathematicaSkipPackets(PetscViewer viewer, int type)
341: {
342: #ifdef PETSC_HAVE_MATHEMATICA
343: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
344: MLINK link = vmath->link; /* The link to Mathematica */
345: int pkt; /* The packet type */
348: while((pkt = MLNextPacket(link)) && (pkt != type))
349: MLNewPacket(link);
350: if (!pkt) {
351: MLClearError(link);
352: SETERRQ(PETSC_ERR_LIB, (char *) MLErrorMessage(link));
353: }
354: return(0);
355: #endif
357: return(0);
358: }
360: /*@C
361: PetscViewerMathematicaLoadGraphics - Change the type of graphics output for Mathematica
363: Input Parameters:
364: . viewer - The Mathematica viewer
365: . type - The type of graphics, e.g. GRAPHICS_MOTIF, GRAPHICS_PS_FILE, GRAPHICS_PS_STDOUT
367: Level: intermediate
369: .keywords PetscViewer, Mathematica, packets
370: .seealso PetscViewerMathematicaSetName(), PetscViewerMathematicaSkipPackets()
371: @*/
372: int PetscViewerMathematicaLoadGraphics(PetscViewer viewer, GraphicsType type)
373: {
374: #ifdef PETSC_HAVE_MATHEMATICA
375: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
376: MLINK link = vmath->link; /* The link to Mathematica */
377: char programName[256];
378: int ierr;
381: /* Load graphics package */
382: MLPutFunction(link, "EvaluatePacket", 1);
383: switch(type)
384: {
385: case GRAPHICS_MOTIF:
386: MLPutFunction(link, "Get", 1);
387: MLPutString(link, "Motif.m");
388: break;
389: case GRAPHICS_PS_FILE:
390: MLPutFunction(link, "CompoundExpression", 4);
391: MLPutFunction(link, "Set", 2);
392: MLPutSymbol(link, "PetscGraphicsCounter");
393: MLPutInteger(link, 0);
394: MLPutFunction(link, "SetDelayed", 2);
395: MLPutSymbol(link, "$Display");
396: MLPutSymbol(link, "$FileDisplay");
397: MLPutFunction(link, "Set", 2);
398: MLPutSymbol(link, "$FileDisplay");
399: MLPutFunction(link, "OpenWrite", 1);
400: MLPutFunction(link, "StringJoin", 3);
401: if (!PetscGetProgramName(programName, 255))
402: MLPutString(link, programName);
403: else
404: MLPutString(link, "GVec");
405: MLPutFunction(link, "ToString", 1);
406: MLPutSymbol(link, "PetscGraphicsCounter");
407: MLPutString(link, ".ps");
408: MLPutFunction(link, "Set", 2);
409: MLPutSymbol(link, "$DisplayFunction");
410: MLPutFunction(link, "Function", 1);
411: MLPutFunction(link, "CompoundExpression", 2);
412: MLPutFunction(link, "Display", 3);
413: MLPutSymbol(link, "$Display");
414: MLPutFunction(link, "Slot", 1);
415: MLPutInteger(link, 1);
416: MLPutString(link, "EPS");
417: MLPutFunction(link, "Increment", 1);
418: MLPutSymbol(link, "PetscGraphicsCounter");
419: break;
420: case GRAPHICS_PS_STDOUT:
421: MLPutFunction(link, "Get", 1);
422: MLPutString(link, "PSDirect.m");
423: break;
424: default:
425: SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid graphics type: %d", type);
426: }
427: MLEndPacket(link);
429: /* Skip packets until ReturnPacket */
430: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
431: /* Skip ReturnPacket */
432: MLNewPacket(link);
434: /* Load PlotField.m for vector plots */
435: MLPutFunction(link, "EvaluatePacket", 1);
436: MLPutFunction(link, "Get", 1);
437: MLPutString(link, "Graphics/PlotField.m");
438: MLEndPacket(link);
440: /* Skip packets until ReturnPacket */
441: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
442: /* Skip ReturnPacket */
443: MLNewPacket(link);
445: return(0);
446: #else
448: switch(type)
449: {
450: case GRAPHICS_MOTIF:
451: case GRAPHICS_PS_FILE:
452: case GRAPHICS_PS_STDOUT:
453: break;
454: default:
455: SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid graphics type: %d", type);
456: }
457: return(0);
458: #endif
459: }
461: /*@C
462: PetscViewerMathematicaGetName - Retrieve the default name for objects communicated to Mathematica
464: Input Parameter:
465: . viewer - The Mathematica viewer
467: Output Parameter:
468: . name - The name for new objects created in Mathematica
470: Level: intermediate
472: .keywords PetscViewer, Mathematica, name
473: .seealso PetscViewerMathematicaSetName(), PetscViewerMathematicaClearName()
474: @*/
475: int PetscViewerMathematicaGetName(PetscViewer viewer, const char **name)
476: {
477: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
482: *name = vmath->objName;
483: return(0);
484: }
486: /*@C
487: PetscViewerMathematicaSetName - Override the default name for objects communicated to Mathematica
489: Input Parameters:
490: . viewer - The Mathematica viewer
491: . name - The name for new objects created in Mathematica
493: Level: intermediate
495: .keywords PetscViewer, Mathematica, name
496: .seealso PetscViewerMathematicaSetName(), PetscViewerMathematicaClearName()
497: @*/
498: int PetscViewerMathematicaSetName(PetscViewer viewer, const char name[])
499: {
500: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
505: vmath->objName = name;
506: return(0);
507: }
509: /*@C
510: PetscViewerMathematicaClearName - Use the default name for objects communicated to Mathematica
512: Input Parameter:
513: . viewer - The Mathematica viewer
515: Level: intermediate
517: .keywords PetscViewer, Mathematica, name
518: .seealso PetscViewerMathematicaGetName(), PetscViewerMathematicaSetName()
519: @*/
520: int PetscViewerMathematicaClearName(PetscViewer viewer)
521: {
522: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
526: vmath->objName = PETSC_NULL;
527: return(0);
528: }
530: /*@C
531: PetscViewerMathematicaGetVector - Retrieve a vector from Mathematica
533: Input Parameter:
534: . viewer - The Mathematica viewer
536: Output Parameter:
537: . v - The vector
539: Level: intermediate
541: .keywords PetscViewer, Mathematica, vector
542: .seealso VecView()
543: @*/
544: int PetscViewerMathematicaGetVector(PetscViewer viewer, Vec v)
545: {
546: #ifdef PETSC_HAVE_MATHEMATICA
547: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
548: MLINK link; /* The link to Mathematica */
549: char *name;
550: Scalar *mArray;
551: Scalar *array;
552: long mSize;
553: int size;
554: int ierr;
560: /* Determine the object name */
561: if (vmath->objName == PETSC_NULL)
562: name = "vec";
563: else
564: name = vmath->objName;
566: link = vmath->link;
567: VecGetLocalSize(v, &size);
568: VecGetArray(v, &array);
569: MLPutFunction(link, "EvaluatePacket", 1);
570: MLPutSymbol(link, name);
571: MLEndPacket(link);
572: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
573: MLGetRealList(link, &mArray, &mSize);
574: if (size != mSize) SETERRQ(PETSC_ERR_ARG_WRONG, "Incompatible vector sizes");
575: PetscMemcpy(array, mArray, mSize * sizeof(double));
576: MLDisownRealList(link, mArray, mSize);
577: VecRestoreArray(v, &array);
579: return(0);
580: #endif
582: return(0);
583: }
585: /*----------------------------------------- Private Functions -------------------------------------------------------*/
586: int PetscViewerMathematicaPutArray_Private(PetscViewer viewer, int n, double *a)
587: {
588: #ifdef PETSC_HAVE_MATHEMATICA
589: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
590: MLINK link = vmath->link; /* The link to Mathematica */
591: char *name;
592: int ierr;
595: /* Determine the object name */
596: if (vmath->objName == PETSC_NULL)
597: name = "vec";
598: else
599: name = vmath->objName;
601: /* Send the CSRMatrix object */
602: MLPutFunction(link, "EvaluatePacket", 1);
603: MLPutFunction(link, "Set", 2);
604: MLPutSymbol(link, name);
605: MLPutRealList(link, a, n);
606: MLEndPacket(link);
607: /* Skip packets until ReturnPacket */
608: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
609: /* Skip ReturnPacket */
610: MLNewPacket(link);
612: return(0);
613: #endif
615: return(0);
616: }
618: int PetscViewerMathematicaPutMatrix_Private(PetscViewer viewer, int m, int n, double *a)
619: {
620: #ifdef PETSC_HAVE_MATHEMATICA
621: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
622: MLINK link = vmath->link; /* The link to Mathematica */
623: char *name;
624: int ierr;
627: /* Determine the object name */
628: if (vmath->objName == PETSC_NULL)
629: name = "mat";
630: else
631: name = vmath->objName;
633: /* Send the dense matrix object */
634: MLPutFunction(link, "EvaluatePacket", 1);
635: MLPutFunction(link, "Set", 2);
636: MLPutSymbol(link, name);
637: MLPutFunction(link, "Transpose", 1);
638: MLPutFunction(link, "Partition", 2);
639: MLPutRealList(link, a, m*n);
640: MLPutInteger(link, m);
641: MLEndPacket(link);
642: /* Skip packets until ReturnPacket */
643: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
644: /* Skip ReturnPacket */
645: MLNewPacket(link);
647: return(0);
648: #endif
650: return(0);
651: }
653: int PetscViewerMathematicaPutSparse_Private(PetscViewer viewer, int m, int n, int *i, int *j, double *a)
654: {
655: #ifdef PETSC_HAVE_MATHEMATICA
656: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
657: MLINK link = vmath->link; /* The link to Mathematica */
658: const char *symbol;
659: char *name;
660: PetscTruth match;
661: int ierr;
664: /* Determine the object name */
665: if (vmath->objName == PETSC_NULL)
666: name = "mat";
667: else
668: name = vmath->objName;
670: /* Make sure Mathematica recognizes sparse matrices */
671: MLPutFunction(link, "EvaluatePacket", 1);
672: MLPutFunction(link, "Needs", 1);
673: MLPutString(link, "LinearAlgebra`CSRMatrix`");
674: MLEndPacket(link);
675: /* Skip packets until ReturnPacket */
676: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
677: /* Skip ReturnPacket */
678: MLNewPacket(link);
680: /* Send the CSRMatrix object */
681: MLPutFunction(link, "EvaluatePacket", 1);
682: MLPutFunction(link, "Set", 2);
683: MLPutSymbol(link, name);
684: MLPutFunction(link, "CSRMatrix", 5);
685: MLPutInteger(link, m);
686: MLPutInteger(link, n);
687: MLPutFunction(link, "Plus", 2);
688: MLPutIntegerList(link, i, m+1);
689: MLPutInteger(link, 1);
690: MLPutFunction(link, "Plus", 2);
691: MLPutIntegerList(link, j, i[m]);
692: MLPutInteger(link, 1);
693: MLPutRealList(link, a, i[m]);
694: MLEndPacket(link);
695: /* Skip packets until ReturnPacket */
696: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
697: /* Skip ReturnPacket */
698: MLNewPacket(link);
700: /* Check that matrix is valid */
701: MLPutFunction(link, "EvaluatePacket", 1);
702: MLPutFunction(link, "ValidQ", 1);
703: MLPutSymbol(link, name);
704: MLEndPacket(link);
705: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
706: MLGetSymbol(link, &symbol);
707: PetscStrcmp("True", (char *) symbol, &match);
708: if (match == PETSC_FALSE) {
709: MLDisownSymbol(link, symbol);
710: SETERRQ(PETSC_ERR_PLIB, "Invalid CSR matrix in Mathematica");
711: }
712: MLDisownSymbol(link, symbol);
713: /* Skip ReturnPacket */
714: MLNewPacket(link);
716: return(0);
717: #endif
719: return(0);
720: }
722: /*------------------------------------------- ML Functions ----------------------------------------------------------*/
723: int PetscViewerMathematicaPartitionMesh(PetscViewer viewer, int numElements, int numVertices, double *vertices, int **mesh,
724: int *numParts, int *colPartition)
725: {
726: #ifdef PETSC_HAVE_MATHEMATICA
727: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
728: MLINK link; /* The link to Mathematica */
729: const char *symbol;
730: int numOptions;
731: int partSize;
732: int *part;
733: long numCols;
734: int col;
735: PetscTruth match, opt;
736: int ierr;
740: link = vmath->link;
742: /* Make sure that the reduce.m package is loaded */
743: MLPutFunction(link, "EvaluatePacket", 1);
744: MLPutFunction(link, "Get", 1);
745: MLPutFunction(link, "StringJoin", 2);
746: MLPutFunction(link, "Environment", 1);
747: MLPutString(link, "PETSC_DIR");
748: MLPutString(link, "/src/pc/impls/ml/reduce.m");
749: MLEndPacket(link);
750: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
751: MLGetSymbol(link, &symbol);
752: PetscStrcmp("$Failed", (char *) symbol, &match);
753: if (match == PETSC_TRUE) {
754: MLDisownSymbol(link, symbol);
755: SETERRQ(PETSC_ERR_FILE_OPEN, "Unable to load package reduce.m");
756: }
757: MLDisownSymbol(link, symbol);
759: /* Partition the mesh */
760: numOptions = 0;
761: partSize = 0;
762: PetscOptionsGetInt(PETSC_NULL, "-pc_ml_partition_size", &partSize, &opt);
763: if ((opt == PETSC_TRUE) && (partSize > 0))
764: numOptions++;
765: MLPutFunction(link, "EvaluatePacket", 1);
766: MLPutFunction(link, "PartitionMesh", 1 + numOptions);
767: MLPutFunction(link, "MeshData", 5);
768: MLPutInteger(link, numElements);
769: MLPutInteger(link, numVertices);
770: MLPutInteger(link, numVertices);
771: MLPutFunction(link, "Partition", 2);
772: MLPutRealList(link, vertices, numVertices*2);
773: MLPutInteger(link, 2);
774: MLPutFunction(link, "Partition", 2);
775: MLPutIntegerList(link, mesh[MESH_ELEM], numElements*3);
776: MLPutInteger(link, 3);
777: if (partSize > 0)
778: {
779: MLPutFunction(link, "Rule", 2);
780: MLPutSymbol(link, "PartitionSize");
781: MLPutInteger(link, partSize);
782: }
783: MLEndPacket(link);
784: /* Skip packets until ReturnPacket */
785: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
787: /* Get the vertex partiton */
788: MLGetIntegerList(link, &part, &numCols);
789: if (numCols != numVertices) SETERRQ(PETSC_ERR_PLIB, "Invalid partition");
790: for(col = 0, *numParts = 0; col < numCols; col++) {
791: colPartition[col] = part[col]-1;
792: *numParts = PetscMax(*numParts, part[col]);
793: }
794: MLDisownIntegerList(link, part, numCols);
796: return(0);
797: #endif
799: return(0);
800: }
802: int PetscViewerMathematicaReduce(PetscViewer viewer, PC pc, int thresh)
803: {
804: #ifdef PETSC_HAVE_MATHEMATICA
805: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
806: MLINK link; /* The link to Mathematica */
807: PC_Multilevel *ml;
808: int *range;
809: long numRange;
810: int *null;
811: long numNull;
812: const char *symbol;
813: int numOptions;
814: int partSize;
815: int row, col;
816: PetscTruth match, opt;
817: int ierr;
822: link = vmath->link;
823: ml = (PC_Multilevel *) pc->data;
825: /* Make sure that the reduce.m package is loaded */
826: MLPutFunction(link, "EvaluatePacket", 1);
827: MLPutFunction(link, "Get", 1);
828: MLPutFunction(link, "StringJoin", 2);
829: MLPutFunction(link, "Environment", 1);
830: MLPutString(link, "PETSC_DIR");
831: MLPutString(link, "/src/pc/impls/ml/reduce.m");
832: MLEndPacket(link);
833: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
834: MLGetSymbol(link, &symbol);
835: PetscStrcmp("$Failed", (char *) symbol, &match);
836: if (match == PETSC_TRUE) {
837: MLDisownSymbol(link, symbol);
838: SETERRQ(PETSC_ERR_FILE_OPEN, "Unable to load package reduce.m");
839: }
840: MLDisownSymbol(link, symbol);
842: /* mesh = MeshData[] */
843: MLPutFunction(link, "EvaluatePacket", 1);
844: MLPutFunction(link, "Set", 2);
845: MLPutSymbol(link, "mesh");
846: MLPutFunction(link, "MeshData", 5);
847: MLPutInteger(link, ml->numElements[0]);
848: MLPutInteger(link, ml->numVertices[0]);
849: MLPutInteger(link, ml->numVertices[0]);
850: MLPutFunction(link, "Partition", 2);
851: MLPutRealList(link, ml->vertices, ml->numVertices[0]*2);
852: MLPutInteger(link, 2);
853: MLPutFunction(link, "Partition", 2);
854: MLPutIntegerList(link, ml->meshes[0][MESH_ELEM], ml->numElements[0]*3);
855: MLPutInteger(link, 3);
856: MLEndPacket(link);
857: /* Skip packets until ReturnPacket */
858: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
859: /* Skip ReturnPacket */
860: MLNewPacket(link);
861: /* Check that mesh is valid */
862: MLPutFunction(link, "EvaluatePacket", 1);
863: MLPutFunction(link, "ValidQ", 1);
864: MLPutSymbol(link, "mesh");
865: MLEndPacket(link);
866: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
867: MLGetSymbol(link, &symbol);
868: PetscStrcmp("True", (char *) symbol, &match);
869: if (match == PETSC_FALSE) {
870: MLDisownSymbol(link, symbol);
871: SETERRQ(PETSC_ERR_PLIB, "Invalid mesh in Mathematica");
872: }
873: MLDisownSymbol(link, symbol);
875: /* mat = gradient matrix */
876: MatView(ml->locB, viewer);
878: /* mattML = ReduceMatrix[mat,ml->minNodes] */
879: numOptions = 0;
880: partSize = 0;
881: PetscOptionsGetInt(PETSC_NULL, "-pc_ml_partition_size", &partSize, &opt);
882: if ((opt == PETSC_TRUE) && (partSize > 0))
883: numOptions++;
884: MLPutFunction(link, "EvaluatePacket", 1);
885: MLPutFunction(link, "Set", 2);
886: MLPutSymbol(link, "mattML");
887: MLPutFunction(link, "ReduceMatrix", 3 + numOptions);
888: MLPutSymbol(link, "mesh");
889: MLPutSymbol(link, "mat");
890: MLPutInteger(link, thresh);
891: if (partSize > 0) {
892: MLPutFunction(link, "Rule", 2);
893: MLPutSymbol(link, "PartitionSize");
894: MLPutInteger(link, partSize);
895: }
896: MLEndPacket(link);
897: /* Skip packets until ReturnPacket */
898: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
899: /* Skip ReturnPacket */
900: MLNewPacket(link);
901: /* Check that mattML is valid */
902: MLPutFunction(link, "EvaluatePacket", 1);
903: MLPutFunction(link, "ValidQ", 1);
904: MLPutSymbol(link, "mattML");
905: MLEndPacket(link);
906: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
907: MLGetSymbol(link, &symbol);
908: PetscStrcmp("True", (char *) symbol, &match);
909: if (match == PETSC_FALSE) {
910: MLDisownSymbol(link, symbol);
911: SETERRQ(PETSC_ERR_PLIB, "Invalid MLData in Mathematica");
912: }
913: MLDisownSymbol(link, symbol);
915: /* Copy information to the preconditioner */
916: MLPutFunction(link, "EvaluatePacket", 1);
917: MLPutFunction(link, "Part", 2);
918: MLPutSymbol(link, "mattML");
919: MLPutInteger(link, 3);
920: MLEndPacket(link);
921: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
922: MLGetInteger(link, &ml->numLevels);
924: /* Create lists of the range and nullspace columns */
925: MLPutFunction(link, "EvaluatePacket", 1);
926: MLPutFunction(link, "GetRange", 1);
927: MLPutSymbol(link, "mattML");
928: MLEndPacket(link);
929: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
930: MLGetIntegerList(link, &range, &numRange);
931: if (numRange > ml->sOrder->numLocVars) SETERRQ(PETSC_ERR_PLIB, "Invalid size for range space");
932: ml->rank = numRange;
933: ml->globalRank = ml->rank;
934: if (ml->rank > 0) {
935: PetscMalloc(numRange * sizeof(int), &ml->range);
936: for(row = 0; row < numRange; row++)
937: ml->range[row] = range[row]-1;
938: }
939: MLDisownIntegerList(link, range, numRange);
941: MLPutFunction(link, "EvaluatePacket", 1);
942: MLPutFunction(link, "GetNullColumns", 1);
943: MLPutSymbol(link, "mattML");
944: MLEndPacket(link);
945: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
946: MLGetIntegerList(link, &null, &numNull);
947: if (numRange + numNull != ml->sOrder->numLocVars) SETERRQ(PETSC_ERR_PLIB, "Invalid size for range and null spaces");
948: ml->numLocNullCols = numNull;
949: if (numNull > 0)
950: {
951: PetscMalloc(numNull * sizeof(int), &ml->nullCols);
952: for(col = 0; col < numNull; col++)
953: ml->nullCols[col] = null[col] - 1;
954: }
955: MLDisownIntegerList(link, null, numNull);
957: return(0);
958: #endif
960: return(0);
961: }
963: int PetscViewerMathematicaMultiLevelConvert(PetscViewer viewer, PC pc)
964: {
965: #ifdef PETSC_HAVE_MATHEMATICA
966: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
967: MLINK link; /* The link to Mathematica */
968: PC_Multilevel *ml;
969: Mat_SeqAIJ *grad;
970: int *numPartitions;
971: int *numPartitionCols, *cols;
972: int *numPartitionRows, *rows;
973: double *U, *singVals, *V;
974: long *Udims, *Vdims;
975: char **Uheads, **Vheads;
976: int *nnz;
977: int *offsets;
978: double *vals;
979: long numLevels, numParts, numCols, numRows, Udepth, numSingVals, Vdepth, len;
980: int numBdRows, numBdCols;
981: int mesh, level, part, col, row;
982: int ierr;
987: link = vmath->link;
988: ml = (PC_Multilevel *) pc->data;
990: /* ml->numLevels = ml[[3]] */
991: MLPutFunction(link, "EvaluatePacket", 1);
992: MLPutFunction(link, "Part", 2);
993: MLPutSymbol(link, "mattML");
994: MLPutInteger(link, 3);
995: MLEndPacket(link);
996: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
997: MLGetInteger(link, &ml->numLevels);
999: /* ml->numMeshes = Length[ml[[4]]] */
1000: MLPutFunction(link, "EvaluatePacket", 1);
1001: MLPutFunction(link, "Length", 1);
1002: MLPutFunction(link, "Part", 2);
1003: MLPutSymbol(link, "mattML");
1004: MLPutInteger(link, 4);
1005: MLEndPacket(link);
1006: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1007: MLGetInteger(link, &ml->numMeshes);
1009: /* ml->numElements[level] = ml[[4,level,1]] */
1010: PetscMalloc(ml->numMeshes * sizeof(int), &ml->numElements);
1012: /* ml->numVertices[level] = ml[[4,level,2]] */
1013: PetscMalloc(ml->numMeshes * sizeof(int), &ml->numVertices);
1015: /* ml->meshes = ml[[4]] */
1016: PetscMalloc(ml->numMeshes * sizeof(int **), &ml->meshes);
1017: for(mesh = 0; mesh < ml->numMeshes; mesh++) {
1018: PetscMalloc(NUM_MESH_DIV * sizeof(int *), &ml->meshes[mesh]);
1019: /* Here we should get meshes */
1020: PetscMalloc(1 * sizeof(int), &ml->meshes[mesh][MESH_OFFSETS]);
1021: PetscMalloc(1 * sizeof(int), &ml->meshes[mesh][MESH_ADJ]);
1022: PetscMalloc(1 * sizeof(int), &ml->meshes[mesh][MESH_ELEM]);
1023: }
1025: /* ml->numPartitions = Map[Length,Drop[ml[[5]],-1]] */
1026: MLPutFunction(link, "EvaluatePacket", 1);
1027: MLPutFunction(link, "Map", 2);
1028: MLPutSymbol(link, "Length");
1029: MLPutFunction(link, "Drop", 2);
1030: MLPutFunction(link, "Part", 2);
1031: MLPutSymbol(link, "mattML");
1032: MLPutInteger(link, 5);
1033: MLPutInteger(link, -1);
1034: MLEndPacket(link);
1035: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1036: MLGetIntegerList(link, &numPartitions, &numLevels);
1037: if (numLevels != ml->numLevels) SETERRQ(PETSC_ERR_PLIB, "Invalid node partition in MLData object");
1038: if (numLevels > 0) {
1039: PetscMalloc(numLevels * sizeof(int), &ml->numPartitions);
1040: PetscMemcpy(ml->numPartitions, numPartitions, numLevels * sizeof(int));
1041: }
1042: MLDisownIntegerList(link, numPartitions, numLevels);
1044: if (ml->numLevels > 0) {
1045: /* ml->numPartitionCols = Map[Length,ml[[5,level]]] */
1046: PetscMalloc(ml->numLevels * sizeof(int *), &ml->numPartitionCols);
1047: for(level = 0; level < ml->numLevels; level++) {
1048: MLPutFunction(link, "EvaluatePacket", 1);
1049: MLPutFunction(link, "Map", 2);
1050: MLPutSymbol(link, "Length");
1051: MLPutFunction(link, "Part", 3);
1052: MLPutSymbol(link, "mattML");
1053: MLPutInteger(link, 5);
1054: MLPutInteger(link, level+1);
1055: MLEndPacket(link);
1056: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1057: MLGetIntegerList(link, &numPartitionCols, &numParts);
1058: if (numParts != ml->numPartitions[level]) SETERRQ(PETSC_ERR_PLIB, "Invalid node partition in MLData object");
1059: if (numParts > 0) {
1060: PetscMalloc(numParts * sizeof(int), &ml->numPartitionCols[level]);
1061: PetscMemcpy(ml->numPartitionCols[level], numPartitionCols, numParts * sizeof(int));
1062: }
1063: MLDisownIntegerList(link, numPartitionCols, numParts);
1064: }
1066: /* ml->colPartition[level][part] = ml[[5,level,part]] */
1067: PetscMalloc(ml->numLevels * sizeof(int **), &ml->colPartition);
1068: for(level = 0; level < ml->numLevels; level++) {
1069: if (ml->numPartitions[level] == 0) continue;
1070: PetscMalloc(ml->numPartitions[level] * sizeof(int *), &ml->colPartition[level]);
1071: for(part = 0; part < ml->numPartitions[level]; part++) {
1072: MLPutFunction(link, "EvaluatePacket", 1);
1073: MLPutFunction(link, "Part", 4);
1074: MLPutSymbol(link, "mattML");
1075: MLPutInteger(link, 5);
1076: MLPutInteger(link, level+1);
1077: MLPutInteger(link, part+1);
1078: MLEndPacket(link);
1079: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1080: MLGetIntegerList(link, &cols, &numCols);
1081: if (numCols != ml->numPartitionCols[level][part]) SETERRQ(PETSC_ERR_PLIB, "Invalid node partition in MLData object");
1082: if (numCols > 0) {
1083: PetscMalloc(numCols * sizeof(int), &ml->colPartition[level][part]);
1084: /* Convert to zero-based indices */
1085: for(col = 0; col < numCols; col++) ml->colPartition[level][part][col] = cols[col] - 1;
1086: }
1087: MLDisownIntegerList(link, cols, numCols);
1088: }
1089: }
1091: /* ml->numPartitionRows = Map[Length,FlattenAt[ml[[6,level]],1]] */
1092: PetscMalloc(ml->numLevels * sizeof(int *), &ml->numPartitionRows);
1093: for(level = 0; level < ml->numLevels; level++) {
1094: MLPutFunction(link, "EvaluatePacket", 1);
1095: MLPutFunction(link, "Map", 2);
1096: MLPutSymbol(link, "Length");
1097: MLPutFunction(link, "FlattenAt", 2);
1098: MLPutFunction(link, "Part", 3);
1099: MLPutSymbol(link, "mattML");
1100: MLPutInteger(link, 6);
1101: MLPutInteger(link, level+1);
1102: MLPutInteger(link, 1);
1103: MLEndPacket(link);
1104: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1105: MLGetIntegerList(link, &numPartitionRows, &numParts);
1106: if (numParts != ml->numPartitions[level] + NUM_PART_ROW_DIV - 1) {
1107: SETERRQ(PETSC_ERR_PLIB, "Invalid edge partition in MLData object");
1108: }
1109: PetscMalloc(numParts * sizeof(int), &ml->numPartitionRows[level]);
1110: PetscMemcpy(ml->numPartitionRows[level], numPartitionRows, numParts * sizeof(int));
1111: MLDisownIntegerList(link, numPartitionRows, numParts);
1112: }
1114: /* ml->rowPartition[level][PART_ROW_INT][part] = ml[[6,level,1,part]]
1115: ml->rowPartition[level][PART_ROW_BD] = ml[[6,level,2]]
1116: ml->rowPartition[level][PART_ROW_RES] = ml[[6,level,3]] */
1117: PetscMalloc(ml->numLevels * sizeof(int ***), &ml->rowPartition);
1118: for(level = 0; level < ml->numLevels; level++) {
1119: PetscMalloc(NUM_PART_ROW_DIV * sizeof(int **), &ml->rowPartition[level]);
1120: /* Interior rows */
1121: if (ml->numPartitions[level] > 0) {
1122: PetscMalloc(ml->numPartitions[level] * sizeof(int *), &ml->rowPartition[level][PART_ROW_INT]);
1123: for(part = 0; part < ml->numPartitions[level]; part++) {
1124: MLPutFunction(link, "EvaluatePacket", 1);
1125: MLPutFunction(link, "Part", 5);
1126: MLPutSymbol(link, "mattML");
1127: MLPutInteger(link, 6);
1128: MLPutInteger(link, level+1);
1129: MLPutInteger(link, 1);
1130: MLPutInteger(link, part+1);
1131: MLEndPacket(link);
1132: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1133: MLGetIntegerList(link, &rows, &numRows);
1134: if (numRows != ml->numPartitionRows[level][part]) {
1135: SETERRQ(PETSC_ERR_PLIB, "Invalid edge partition in MLData object");
1136: }
1137: if (numRows > 0) {
1138: PetscMalloc(numRows * sizeof(int), &ml->rowPartition[level][PART_ROW_INT][part]);
1139: /* Convert to zero-based indices */
1140: for(row = 0; row < numRows; row++) {
1141: ml->rowPartition[level][PART_ROW_INT][part][row] = rows[row] - 1;
1142: }
1143: }
1144: MLDisownIntegerList(link, rows, numRows);
1145: }
1146: }
1147: /* Boundary rows */
1148: PetscMalloc(1 * sizeof(int *), &ml->rowPartition[level][PART_ROW_BD]);
1149: MLPutFunction(link, "EvaluatePacket", 1);
1150: MLPutFunction(link, "Part", 4);
1151: MLPutSymbol(link, "mattML");
1152: MLPutInteger(link, 6);
1153: MLPutInteger(link, level+1);
1154: MLPutInteger(link, 2);
1155: MLEndPacket(link);
1156: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1157: MLGetIntegerList(link, &rows, &numRows);
1158: if (numRows != ml->numPartitionRows[level][ml->numPartitions[level]]) {
1159: SETERRQ(PETSC_ERR_PLIB, "Invalid edge partition in MLData object");
1160: }
1161: if (numRows > 0) {
1162: PetscMalloc(numRows * sizeof(int), &ml->rowPartition[level][PART_ROW_BD][0]);
1163: /* Convert to zero-based indices */
1164: for(row = 0; row < numRows; row++) {
1165: ml->rowPartition[level][PART_ROW_BD][0][row] = rows[row] - 1;
1166: }
1167: }
1168: MLDisownIntegerList(link, rows, numRows);
1169: /* Residual rows*/
1170: PetscMalloc(1 * sizeof(int *), &ml->rowPartition[level][PART_ROW_RES]);
1171: MLPutFunction(link, "EvaluatePacket", 1);
1172: MLPutFunction(link, "Part", 4);
1173: MLPutSymbol(link, "mattML");
1174: MLPutInteger(link, 6);
1175: MLPutInteger(link, level+1);
1176: MLPutInteger(link, 3);
1177: MLEndPacket(link);
1178: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1179: MLGetIntegerList(link, &rows, &numRows);
1180: if (numRows != ml->numPartitionRows[level][ml->numPartitions[level]+1]) {
1181: SETERRQ(PETSC_ERR_PLIB, "Invalid edge partition in MLData object");
1182: }
1183: if (numRows > 0) {
1184: PetscMalloc(numRows * sizeof(int), &ml->rowPartition[level][PART_ROW_RES][0]);
1185: /* Convert to zero-based indices */
1186: for(row = 0; row < numRows; row++) {
1187: ml->rowPartition[level][PART_ROW_RES][0][row] = rows[row] - 1;
1188: }
1189: }
1190: MLDisownIntegerList(link, rows, numRows);
1191: }
1192: } else {
1193: PetscMalloc(1 * sizeof(int), &ml->numPartitions);
1194: PetscMalloc(1 * sizeof(int *), &ml->numPartitionCols);
1195: PetscMalloc(1 * sizeof(int **), &ml->colPartition);
1196: PetscMalloc(1 * sizeof(int *), &ml->numPartitionRows);
1197: PetscMalloc(1 * sizeof(int ***), &ml->rowPartition);
1198: }
1200: /* ml->numRows = ml[[1]] */
1201: MLPutFunction(link, "EvaluatePacket", 1);
1202: MLPutFunction(link, "Part", 2);
1203: MLPutSymbol(link, "mattML");
1204: MLPutInteger(link, 1);
1205: MLEndPacket(link);
1206: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1207: MLGetInteger(link, &ml->numRows);
1209: /* ml->numCols = ml[[2]] */
1210: MLPutFunction(link, "EvaluatePacket", 1);
1211: MLPutFunction(link, "Part", 2);
1212: MLPutSymbol(link, "mattML");
1213: MLPutInteger(link, 2);
1214: MLEndPacket(link);
1215: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1216: MLGetInteger(link, &ml->numCols);
1218: /* ml->zeroTol = ml[[9]] */
1219: MLPutFunction(link, "EvaluatePacket", 1);
1220: MLPutFunction(link, "N", 1);
1221: MLPutFunction(link, "Part", 2);
1222: MLPutSymbol(link, "mattML");
1223: MLPutInteger(link, 9);
1224: MLEndPacket(link);
1225: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1226: MLGetReal(link, &ml->zeroTol);
1228: if (ml->numLevels > 0) {
1229: /* ml->factors[level][part][FACT_U] = ml[[7,level,part,1]]
1230: ml->factors[level][part][FACT_DINV] = Divide[1,Select[ml[[7,level,part,2]],(#>ml[[9]])&]]
1231: ml->factors[level][part][FACT_V] = ml[[7,level,part,3]] */
1232: PetscMalloc(ml->numLevels * sizeof(double ***), &ml->factors);
1233: for(level = 0; level < ml->numLevels; level++) {
1234: if (ml->numPartitions[level] == 0) continue;
1235: PetscMalloc(ml->numPartitions[level] * sizeof(double **), &ml->factors[level]);
1236: for(part = 0; part < ml->numPartitions[level]; part++) {
1237: PetscMalloc(NUM_FACT_DIV * sizeof(double *), &ml->factors[level][part]);
1238: /* U, or U^T in LAPACK terms */
1239: MLPutFunction(link, "EvaluatePacket", 1);
1240: MLPutFunction(link, "Part", 5);
1241: MLPutSymbol(link, "mattML");
1242: MLPutInteger(link, 7);
1243: MLPutInteger(link, level+1);
1244: MLPutInteger(link, part+1);
1245: MLPutInteger(link, 1);
1246: MLEndPacket(link);
1247: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1248: MLGetRealArray(link, &U, &Udims, &Uheads, &Udepth);
1249: if (Udepth > 1) {
1250: if (Udepth != 2) SETERRQ(PETSC_ERR_PLIB, "Invalid U matrix");
1251: if ((Udims[0] != ml->numPartitionRows[level][part]) || (Udims[0] != Udims[1])) {
1252: SETERRQ(PETSC_ERR_PLIB, "Incompatible dimensions for U matrix");
1253: }
1254: PetscMalloc(Udims[0]*Udims[0] * sizeof(double), &ml->factors[level][part][FACT_U]);
1255: /* Notice that LAPACK will think that this is U^T, or U in LAPACK terms */
1256: PetscMemcpy(ml->factors[level][part][FACT_U], U, Udims[0]*Udims[0] * sizeof(double));
1257: } else if (ml->numPartitionRows[level][part] != 0) {
1258: SETERRQ(PETSC_ERR_PLIB, "Missing U matrix");
1259: }
1260: MLDisownRealArray(link, U, Udims, Uheads, Udepth);
1261: /* D^{-1} */
1262: MLPutFunction(link, "EvaluatePacket", 1);
1263: MLPutFunction(link, "Divide", 2);
1264: MLPutReal(link, 1.0);
1265: MLPutFunction(link, "Select", 2);
1266: MLPutFunction(link, "Part", 5);
1267: MLPutSymbol(link, "mattML");
1268: MLPutInteger(link, 7);
1269: MLPutInteger(link, level+1);
1270: MLPutInteger(link, part+1);
1271: MLPutInteger(link, 2);
1272: MLPutFunction(link, "Function", 2);
1273: MLPutSymbol(link, "x");
1274: MLPutFunction(link, "Greater", 2);
1275: MLPutSymbol(link, "x");
1276: MLPutFunction(link, "Part", 2);
1277: MLPutSymbol(link, "mattML");
1278: MLPutInteger(link, 9);
1279: MLEndPacket(link);
1280: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1281: MLGetRealList(link, &singVals, &numSingVals);
1282: if (numSingVals > ml->numPartitionCols[level][part]) {
1283: SETERRQ(PETSC_ERR_PLIB, "Invalid factor list in MLData object");
1284: }
1285: if (ml->numPartitionCols[level][part] > 0) {
1286: PetscMalloc(ml->numPartitionCols[level][part] * sizeof(double), &ml->factors[level][part][FACT_DINV]);
1287: PetscMemzero(ml->factors[level][part][FACT_DINV], ml->numPartitionCols[level][part] * sizeof(double));
1288: PetscMemcpy(ml->factors[level][part][FACT_DINV], singVals, numSingVals * sizeof(double));
1289: }
1290: MLDisownRealList(link, singVals, numSingVals);
1291: /* V^T, but V in LAPACK terms */
1292: MLPutFunction(link, "EvaluatePacket", 1);
1293: MLPutFunction(link, "Transpose", 1);
1294: MLPutFunction(link, "Part", 5);
1295: MLPutSymbol(link, "mattML");
1296: MLPutInteger(link, 7);
1297: MLPutInteger(link, level+1);
1298: MLPutInteger(link, part+1);
1299: MLPutInteger(link, 3);
1300: MLEndPacket(link);
1301: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1302: MLGetRealArray(link, &V, &Vdims, &Vheads, &Vdepth);
1303: if (Vdepth > 1) {
1304: if (Vdepth != 2) SETERRQ(PETSC_ERR_PLIB, "Invalid V matrix");
1305: if ((Vdims[0] != ml->numPartitionCols[level][part]) || (Vdims[0] != Vdims[1])) {
1306: SETERRQ(PETSC_ERR_PLIB, "Incompatible dimensions for U matrix");
1307: }
1308: PetscMalloc(Vdims[0]*Vdims[0] * sizeof(double), &ml->factors[level][part][FACT_V]);
1309: /* Notice that LAPACK will think that this is V, or V^T in LAPACK terms */
1310: PetscMemcpy(ml->factors[level][part][FACT_V], V, Vdims[0]*Vdims[0] * sizeof(double));
1311: } else if (ml->numPartitionCols[level][part] != 0) {
1312: SETERRQ(PETSC_ERR_PLIB, "Missing V matrix");
1313: }
1314: MLDisownRealArray(link, V, Vdims, Vheads, Vdepth);
1315: }
1316: }
1318: /* ml->grads = ml[[8]] */
1319: PetscMalloc(ml->numLevels * sizeof(Mat), &ml->grads);
1320: PetscMalloc(ml->numLevels * sizeof(Vec), &ml->bdReduceVecs);
1321: PetscMalloc(ml->numLevels * sizeof(Vec), &ml->colReduceVecs);
1322: PetscMalloc(ml->numLevels * sizeof(Vec), &ml->colReduceVecs2);
1323: for(level = 0; level < ml->numLevels; level++) {
1324: if (ml->numPartitions[level] == 0) continue;
1325: /* numRows = ml[[8,level,1]] */
1326: MLPutFunction(link, "EvaluatePacket", 1);
1327: MLPutFunction(link, "Part", 4);
1328: MLPutSymbol(link, "mattML");
1329: MLPutInteger(link, 8);
1330: MLPutInteger(link, level+1);
1331: MLPutInteger(link, 1);
1332: MLEndPacket(link);
1333: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1334: MLGetInteger(link, &numBdRows);
1335: VecCreateSeq(PETSC_COMM_SELF, numBdRows, &ml->bdReduceVecs[level]);
1336: /* numCols = ml[[8,level,2]] */
1337: MLPutFunction(link, "EvaluatePacket", 1);
1338: MLPutFunction(link, "Part", 4);
1339: MLPutSymbol(link, "mattML");
1340: MLPutInteger(link, 8);
1341: MLPutInteger(link, level+1);
1342: MLPutInteger(link, 2);
1343: MLEndPacket(link);
1344: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1345: MLGetInteger(link, &numBdCols);
1346: VecCreateSeq(PETSC_COMM_SELF, numBdCols, &ml->colReduceVecs[level]);
1347: VecDuplicate(ml->colReduceVecs[level], &ml->colReduceVecs2[level]);
1348: /* nnz = Map[Apply[Subtract,Sort[#,Greater]]&, Partition[ml[[8,level,3]],2,1]] */
1349: MLPutFunction(link, "EvaluatePacket", 1);
1350: MLPutFunction(link, "Map", 2);
1351: MLPutFunction(link, "Function", 2);
1352: MLPutSymbol(link, "x");
1353: MLPutFunction(link, "Apply", 2);
1354: MLPutSymbol(link, "Subtract");
1355: MLPutFunction(link, "Sort", 2);
1356: MLPutSymbol(link, "x");
1357: MLPutSymbol(link, "Greater");
1358: MLPutFunction(link, "Partition", 3);
1359: MLPutFunction(link, "Part", 4);
1360: MLPutSymbol(link, "mattML");
1361: MLPutInteger(link, 8);
1362: MLPutInteger(link, level+1);
1363: MLPutInteger(link, 3);
1364: MLPutInteger(link, 2);
1365: MLPutInteger(link, 1);
1366: MLEndPacket(link);
1367: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1368: MLGetIntegerList(link, &nnz, &len);
1369: if (len != numBdRows) SETERRQ(PETSC_ERR_PLIB, "Invalid boundary gradient matrix");
1370: MatCreateSeqAIJ(PETSC_COMM_SELF, numBdRows, numBdCols, PETSC_DEFAULT, nnz, &ml->grads[level]);
1371: grad = (Mat_SeqAIJ *) ml->grads[level]->data;
1372: PetscMemcpy(grad->ilen, nnz, numBdRows * sizeof(int));
1373: MLDisownIntegerList(link, nnz, len);
1374: /* ml->grads[level]->i = ml[[8,level,3]] */
1375: MLPutFunction(link, "EvaluatePacket", 1);
1376: MLPutFunction(link, "Part", 4);
1377: MLPutSymbol(link, "mattML");
1378: MLPutInteger(link, 8);
1379: MLPutInteger(link, level+1);
1380: MLPutInteger(link, 3);
1381: MLEndPacket(link);
1382: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1383: MLGetIntegerList(link, &offsets, &len);
1384: if (len != numBdRows+1) SETERRQ(PETSC_ERR_PLIB, "Invalid boundary gradient matrix");
1385: for(row = 0; row <= numBdRows; row++)
1386: grad->i[row] = offsets[row]-1;
1387: MLDisownIntegerList(link, offsets, len);
1388: /* ml->grads[level]->j = ml[[8,level,4]] */
1389: MLPutFunction(link, "EvaluatePacket", 1);
1390: MLPutFunction(link, "Part", 4);
1391: MLPutSymbol(link, "mattML");
1392: MLPutInteger(link, 8);
1393: MLPutInteger(link, level+1);
1394: MLPutInteger(link, 4);
1395: MLEndPacket(link);
1396: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1397: MLGetIntegerList(link, &cols, &len);
1398: if (len != grad->i[numBdRows]) SETERRQ(PETSC_ERR_PLIB, "Invalid boundary gradient matrix");
1399: for(col = 0; col < len; col++)
1400: grad->j[col] = cols[col]-1;
1401: MLDisownIntegerList(link, cols, len);
1402: /* ml->grads[level]->i = ml[[8,level,5]] */
1403: MLPutFunction(link, "EvaluatePacket", 1);
1404: MLPutFunction(link, "Part", 4);
1405: MLPutSymbol(link, "mattML");
1406: MLPutInteger(link, 8);
1407: MLPutInteger(link, level+1);
1408: MLPutInteger(link, 5);
1409: MLEndPacket(link);
1410: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1411: MLGetRealList(link, &vals, &len);
1412: if (len != grad->i[numBdRows]) SETERRQ(PETSC_ERR_PLIB, "Invalid boundary gradient matrix");
1413: PetscMemcpy(grad->a, vals, len * sizeof(double));
1414: MLDisownRealList(link, vals, len);
1415: /* Fix up all the members */
1416: grad->nz += len;
1417: MatAssemblyBegin(ml->grads[level], MAT_FINAL_ASSEMBLY);
1418: MatAssemblyEnd(ml->grads[level], MAT_FINAL_ASSEMBLY);
1419: }
1420: } else {
1421: PetscMalloc(1 * sizeof(double ***), &ml->factors);
1422: PetscMalloc(1 * sizeof(Mat), &ml->grads);
1423: PetscMalloc(1 * sizeof(Vec), &ml->bdReduceVecs);
1424: PetscMalloc(1 * sizeof(Vec), &ml->colReduceVecs);
1425: PetscMalloc(1 * sizeof(Vec), &ml->colReduceVecs2);
1426: }
1428: ml->interiorWorkLen = 1;
1429: for(level = 0; level < ml->numLevels; level++) {
1430: for(part = 0; part < ml->numPartitions[level]; part++)
1431: ml->interiorWorkLen = PetscMax(ml->interiorWorkLen, ml->numPartitionRows[level][part]);
1432: }
1433: PetscMalloc(ml->interiorWorkLen * sizeof(double), &ml->interiorWork);
1434: PetscMalloc(ml->interiorWorkLen * sizeof(double), &ml->interiorWork2);
1436: return(0);
1437: #endif
1439: return(0);
1440: }
1442: /*------------------------------ Functions for Triangular 2d Meshes -------------------------------------------------*/
1443: int PetscViewerMathematicaCreateSamplePoints_Triangular_2D(PetscViewer viewer, GVec v)
1444: {
1445: #ifdef PETSC_HAVE_MATHEMATICA
1446: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
1447: MLINK link = vmath->link; /* The link to Mathematica */
1448: Grid grid;
1449: Mesh_Triangular *tri;
1450: int numNodes;
1451: int *classes;
1452: int *offsets;
1453: double *array;
1454: int *localStart;
1455: int localOffset, comp;
1456: int node, nclass;
1457: int ierr;
1460: ierr = GVecGetGrid(v, &grid);
1461: tri = (Mesh_Triangular *) grid->mesh->data;
1462: numNodes = tri->numNodes;
1463: comp = grid->viewComp;
1464: offsets = grid->order->offsets;
1465: localStart = grid->order->localStart[grid->viewField];
1466: classes = grid->cm->classes;
1468: /* Attach a value to each point in the mesh -- Extra values are put in for fields not
1469: defined on some nodes, but these values are never used */
1470: VecGetArray(v, &array);
1471: MLPutFunction(link, "ReplaceAll", 2);
1472: MLPutFunction(link, "Thread", 1);
1473: MLPutFunction(link, "f", 2);
1474: MLPutSymbol(link, "nodes");
1475: MLPutFunction(link, "List", numNodes);
1476: for(node = 0; node < numNodes; node++)
1477: {
1478: nclass = classes[node];
1479: localOffset = localStart[nclass] + comp;
1480: MLPutReal(link, array[offsets[node] + localOffset]);
1481: }
1482: MLPutFunction(link, "Rule", 2);
1483: MLPutSymbol(link, "f");
1484: MLPutSymbol(link, "Append");
1485: VecRestoreArray(v, &array);
1487: return(0);
1488: #endif
1490: return(0);
1491: }
1493: int PetscViewerMathematicaCreateVectorSamplePoints_Triangular_2D(PetscViewer viewer, GVec v)
1494: {
1495: #ifdef PETSC_HAVE_MATHEMATICA
1496: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
1497: MLINK link = vmath->link; /* The link to Mathematica */
1498: Grid grid;
1499: Mesh_Triangular *tri;
1500: int numNodes;
1501: int *classes;
1502: int *offsets;
1503: int *fieldClasses;
1504: double *array;
1505: int *localStart;
1506: int localOffset;
1507: int node, nclass;
1508: int ierr;
1511: ierr = GVecGetGrid(v, &grid);
1512: tri = (Mesh_Triangular *) grid->mesh->data;
1513: numNodes = tri->numNodes;
1514: fieldClasses = grid->cm->fieldClasses[grid->viewField];
1515: offsets = grid->order->offsets;
1516: localStart = grid->order->localStart[grid->viewField];
1517: classes = grid->cm->classes;
1519: /* Make a list {{{x_0,y_0},{f^0_x,f^0_y}},...} */
1520: VecGetArray(v, &array);
1521: MLPutFunction(link, "ReplaceAll", 2);
1522: MLPutFunction(link, "Thread", 1);
1523: MLPutFunction(link, "f", 2);
1524: MLPutSymbol(link, "nodes");
1525: MLPutFunction(link, "List", numNodes);
1526: for(node = 0; node < numNodes; node++)
1527: {
1528: nclass = classes[node];
1529: if (fieldClasses[nclass])
1530: {
1531: localOffset = localStart[nclass];
1532: MLPutFunction(link, "List", 2);
1533: MLPutReal(link, array[offsets[node] + localOffset]);
1534: MLPutReal(link, array[offsets[node] + localOffset + 1]);
1535: }
1536: else
1537: {
1538: /* Vectors of length zero are ignored */
1539: MLPutFunction(link, "List", 2);
1540: MLPutReal(link, 0.0);
1541: MLPutReal(link, 0.0);
1542: }
1543: }
1544: MLPutFunction(link, "Rule", 2);
1545: MLPutSymbol(link, "f");
1546: MLPutSymbol(link, "List");
1547: VecRestoreArray(v, &array);
1549: return(0);
1550: #endif
1552: return(0);
1553: }
1555: int PetscViewerMathematicaCreateInterpolatedSamplePoints_Triangular_2D(PetscViewer viewer, GVec v, int vComp)
1556: {
1557: #ifdef PETSC_HAVE_MATHEMATICA
1558: #if 0
1559: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
1560: MLINK link = vmath->link; /* The link to Mathematica */
1561: InterpolationContext iCtx;
1562: Grid grid;
1563: Mesh mesh;
1564: double *x, *y, *values;
1565: long m, n;
1566: double startX, endX, incX;
1567: double startY, endY, incY;
1568: int comp;
1569: int proc, row, col;
1570: PetscTruth opt;
1571: int ierr;
1572: #endif
1575: #if 0
1576: ierr = GVecGetGrid(v, &grid);
1577: ierr = GridGetMesh(grid, &mesh);
1578: comp = grid->comp[grid->viewField];
1580: /* This sucks, but I will fix it later (It is for GridReduceInterpolationElementVec_Triangular_2D) */
1581: grid->classesOld = grid->classes;
1583: /* Setup InterpolationContext */
1584: iCtx.vec = v;
1585: iCtx.mesh = mesh;
1586: iCtx.order = grid->order;
1587: iCtx.ghostVec = grid->ghostVec;
1588: iCtx.field = grid->viewField;
1589: iCtx.numProcs = mesh->part->numProcs;
1590: iCtx.rank = mesh->part->rank;
1591: PetscMalloc(iCtx.numProcs * sizeof(int), &iCtx.activeProcs);
1592: PetscMalloc(iCtx.numProcs * sizeof(int), &iCtx.calcProcs);
1593: PetscMalloc(iCtx.numProcs*3 * sizeof(Scalar), &iCtx.coords);
1594: PetscMalloc(iCtx.numProcs * sizeof(Scalar *), &iCtx.values);
1595: for(proc = 0; proc < iCtx.numProcs; proc++) {
1596: PetscMalloc(comp * sizeof(Scalar), &iCtx.values[proc]);
1597: }
1599: /* Setup domain */
1600: startX = 0.0;
1601: startY = 0.0;
1602: endX = 1.0;
1603: endY = 1.0;
1604: ierr = PetscOptionsGetDouble("viewer_", "-math_start_x", &startX, &opt);
1605: ierr = PetscOptionsGetDouble("viewer_", "-math_start_y", &startY, &opt);
1606: ierr = PetscOptionsGetDouble("viewer_", "-math_end_x", &endX, &opt);
1607: ierr = PetscOptionsGetDouble("viewer_", "-math_end_y", &endY, &opt);
1608: ierr = PetscOptionsGetInt("viewer_", "-math_div_x", (int *) &n, &opt);
1609: ierr = PetscOptionsGetInt("viewer_", "-math_div_y", (int *) &m, &opt);
1610: ierr = PetscMalloc((n+1) * sizeof(double), &x);
1611: ierr = PetscMalloc((n+1) * sizeof(double), &y);
1612: ierr = PetscMalloc((n+1)*comp * sizeof(double), &values);
1613: incX = (endX - startX)/n;
1614: incY = (endY - startY)/m;
1616: x[0] = startX;
1617: for(col = 1; col <= n; col++)
1618: x[col] = x[col-1] + incX;
1620: MLPutFunction(link, "List", m+1);
1621: for(row = 0; row <= m; row++)
1622: {
1623: PetscMemzero(values, (n+1)*comp * sizeof(double));
1624: for(col = 0; col <= n; col++)
1625: y[col] = startY + incY*row;
1626: PointFunctionInterpolateField(n+1, comp, x, y, x, values, &iCtx);
1627: if (vComp >= 0)
1628: {
1629: for(col = 0; col <= n; col++)
1630: values[col] = values[col*comp+vComp];
1631: MLPutRealList(link, values, n+1);
1632: }
1633: else
1634: {
1635: MLPutFunction(link, "Transpose", 1);
1636: MLPutFunction(link, "List", 2);
1637: MLPutFunction(link, "Transpose", 1);
1638: MLPutFunction(link, "List", 2);
1639: MLPutRealList(link, x, n+1);
1640: MLPutRealList(link, y, n+1);
1641: MLPutFunction(link, "Partition", 2);
1642: MLPutRealList(link, values, (n+1)*comp);
1643: MLPutInteger(link, comp);
1644: }
1645: }
1647: /* This sucks, but I will fix it later (It is for GridReduceInterpolationElementVec_Triangular_2D) */
1648: grid->classesOld = PETSC_NULL;
1650: /* Cleanup */
1651: PetscFree(x);
1652: PetscFree(y);
1653: PetscFree(values);
1654: PetscFree(iCtx.activeProcs);
1655: PetscFree(iCtx.calcProcs);
1656: PetscFree(iCtx.coords);
1657: for(proc = 0; proc < iCtx.numProcs; proc++) {
1658: PetscFree(iCtx.values[proc]);
1659: }
1660: PetscFree(iCtx.values);
1661: #endif
1663: return(0);
1664: #endif
1666: return(0);
1667: }
1669: int PetscViewerMathematica_Mesh_Triangular_2D(PetscViewer viewer, Mesh mesh)
1670: {
1671: #ifdef PETSC_HAVE_MATHEMATICA
1672: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
1673: MLINK link = vmath->link;
1674: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
1675: int numCorners = tri->numCorners;
1676: int numFaces = tri->numFaces;
1677: int *faces = tri->faces;
1678: int numNodes = tri->numNodes;
1679: double *nodes = tri->nodes;
1680: int node, face, corner;
1681: int ierr;
1684: /* Load package to view graphics in a window */
1685: PetscViewerMathematicaLoadGraphics(viewer, vmath->graphicsType);
1687: /* Send the node coordinates */
1688: MLPutFunction(link, "EvaluatePacket", 1);
1689: MLPutFunction(link, "Set", 2);
1690: MLPutSymbol(link, "nodes");
1691: MLPutFunction(link, "List", numNodes);
1692: for(node = 0; node < numNodes; node++) {
1693: MLPutFunction(link, "List", 2);
1694: MLPutReal(link, nodes[node*2]);
1695: MLPutReal(link, nodes[node*2+1]);
1696: }
1697: MLEndPacket(link);
1698: /* Skip packets until ReturnPacket */
1699: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1700: /* Skip ReturnPacket */
1701: MLNewPacket(link);
1703: /* Send the faces */
1704: MLPutFunction(link, "EvaluatePacket", 1);
1705: MLPutFunction(link, "Set", 2);
1706: MLPutSymbol(link, "faces");
1707: MLPutFunction(link, "List", numFaces);
1708: for(face = 0; face < numFaces; face++) {
1709: MLPutFunction(link, "List", numCorners);
1710: for(corner = 0; corner < numCorners; corner++) {
1711: MLPutReal(link, faces[face*numCorners+corner]);
1712: }
1713: }
1714: MLEndPacket(link);
1715: /* Skip packets until ReturnPacket */
1716: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1717: /* Skip ReturnPacket */
1718: MLNewPacket(link);
1720: return(0);
1721: #endif
1723: return(0);
1724: }
1726: int PetscViewerMathematicaCheckMesh_Triangular_2D(PetscViewer viewer, Mesh mesh)
1727: {
1728: #ifdef PETSC_HAVE_MATHEMATICA
1729: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
1730: MLINK link = vmath->link;
1731: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
1732: int numCorners = tri->numCorners;
1733: int numFaces = tri->numFaces;
1734: int numNodes = tri->numNodes;
1735: const char *symbol;
1736: long args;
1737: int dim, type;
1738: PetscTruth match;
1739: int ierr;
1742: /* Check that nodes are defined */
1743: MLPutFunction(link, "EvaluatePacket", 1);
1744: MLPutFunction(link, "ValueQ", 1);
1745: MLPutSymbol(link, "nodes");
1746: MLEndPacket(link);
1747: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1748: MLGetSymbol(link, &symbol);
1749: PetscStrcmp("True", (char *) symbol, &match);
1750: if (match == PETSC_FALSE) {
1751: MLDisownSymbol(link, symbol);
1752: goto redefineMesh;
1753: }
1754: MLDisownSymbol(link, symbol);
1756: /* Check the dimensions */
1757: MLPutFunction(link, "EvaluatePacket", 1);
1758: MLPutFunction(link, "MatrixQ", 2);
1759: MLPutSymbol(link, "nodes");
1760: MLPutSymbol(link, "NumberQ");
1761: MLEndPacket(link);
1762: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1763: MLGetSymbol(link, &symbol);
1764: PetscStrcmp("True", (char *) symbol, &match);
1765: if (match == PETSC_FALSE) {
1766: MLDisownSymbol(link, symbol);
1767: goto redefineMesh;
1768: }
1769: MLDisownSymbol(link, symbol);
1770: MLPutFunction(link, "EvaluatePacket", 1);
1771: MLPutFunction(link, "Dimensions", 1);
1772: MLPutSymbol(link, "nodes");
1773: MLEndPacket(link);
1774: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1775: args = 0;
1776: type = MLGetNext(link);
1777: MLGetArgCount(link, &args);
1778: if (args != 2) {
1779: MLNewPacket(link);
1780: goto redefineMesh;
1781: }
1782: MLGetSymbol(link, &symbol);
1783: MLDisownSymbol(link, symbol);
1784: MLGetInteger(link, &dim);
1785: if (dim != numNodes) {
1786: MLNewPacket(link);
1787: goto redefineMesh;
1788: }
1789: MLGetInteger(link, &dim);
1790: if (dim != mesh->dim) {
1791: MLNewPacket(link);
1792: goto redefineMesh;
1793: }
1795: /* Check that faces are defined */
1796: MLPutFunction(link, "EvaluatePacket", 1);
1797: MLPutFunction(link, "ValueQ", 1);
1798: MLPutSymbol(link, "faces");
1799: MLEndPacket(link);
1800: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1801: MLGetSymbol(link, &symbol);
1802: PetscStrcmp("True", (char *) symbol, &match);
1803: if (match == PETSC_FALSE) {
1804: MLDisownSymbol(link, symbol);
1805: goto redefineMesh;
1806: }
1807: MLDisownSymbol(link, symbol);
1809: /* Check the dimensions */
1810: MLPutFunction(link, "EvaluatePacket", 1);
1811: MLPutFunction(link, "MatrixQ", 2);
1812: MLPutSymbol(link, "faces");
1813: MLPutSymbol(link, "NumberQ");
1814: MLEndPacket(link);
1815: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1816: MLGetSymbol(link, &symbol);
1817: PetscStrcmp("True", (char *) symbol, &match);
1818: if (match == PETSC_FALSE) {
1819: MLDisownSymbol(link, symbol);
1820: goto redefineMesh;
1821: }
1822: MLDisownSymbol(link, symbol);
1823: MLPutFunction(link, "EvaluatePacket", 1);
1824: MLPutFunction(link, "Dimensions", 1);
1825: MLPutSymbol(link, "faces");
1826: MLEndPacket(link);
1827: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1828: args = 0;
1829: type = MLGetNext(link);
1830: MLGetArgCount(link, &args);
1831: if (args != 2) {
1832: MLNewPacket(link);
1833: goto redefineMesh;
1834: }
1835: MLGetSymbol(link, &symbol);
1836: MLDisownSymbol(link, symbol);
1837: MLGetInteger(link, &dim);
1838: if (dim != numFaces) {
1839: MLNewPacket(link);
1840: goto redefineMesh;
1841: }
1842: MLGetInteger(link, &dim);
1843: if (dim != numCorners) {
1844: MLNewPacket(link);
1845: goto redefineMesh;
1846: }
1848: return(0);
1850: redefineMesh:
1851: MeshView(mesh, viewer);
1852: return(0);
1853: #endif
1855: return(0);
1856: }
1858: int PetscViewerMathematicaTriangulationPlot_Triangular_2D(PetscViewer viewer, GVec v)
1859: {
1860: #ifdef PETSC_HAVE_MATHEMATICA
1861: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
1862: MLINK link = vmath->link; /* The link to Mathematica */
1863: Mesh mesh;
1864: Mesh_Triangular *tri;
1865: Grid grid;
1866: int numCorners;
1867: int ierr;
1870: ierr = GVecGetGrid(v, &grid);
1871: mesh = grid->mesh;
1872: tri = (Mesh_Triangular *) mesh->data;
1873: numCorners = tri->numCorners;
1875: MLPutFunction(link, "Show", 2);
1876: MLPutFunction(link, "Graphics3D", 1);
1877: MLPutFunction(link, "Flatten", 1);
1878: MLPutFunction(link, "Map", 2);
1879: MLPutFunction(link, "Function", 1);
1880: if ((numCorners == 6) && ((grid->cm->numClasses == 1) || (grid->cm->fieldClasses[grid->viewField][1])))
1881: {
1882: MLPutFunction(link, "List", 4);
1883: /* Triangle 0--5--4 */
1884: MLPutFunction(link, "Polygon", 1);
1885: MLPutFunction(link, "Part", 2);
1886: MLPutSymbol(link, "points");
1887: MLPutFunction(link, "Plus", 2);
1888: MLPutFunction(link, "Part", 2);
1889: MLPutFunction(link, "Slot", 1);
1890: MLPutInteger(link, 1);
1891: MLPutFunction(link, "List", 3);
1892: MLPutInteger(link, 1);
1893: MLPutInteger(link, 6);
1894: MLPutInteger(link, 5);
1895: MLPutInteger(link, 1);
1896: /* Triangle 1--3--5 */
1897: MLPutFunction(link, "Polygon", 1);
1898: MLPutFunction(link, "Part", 2);
1899: MLPutSymbol(link, "points");
1900: MLPutFunction(link, "Plus", 2);
1901: MLPutFunction(link, "Part", 2);
1902: MLPutFunction(link, "Slot", 1);
1903: MLPutInteger(link, 1);
1904: MLPutFunction(link, "List", 3);
1905: MLPutInteger(link, 2);
1906: MLPutInteger(link, 4);
1907: MLPutInteger(link, 6);
1908: MLPutInteger(link, 1);
1909: /* Triangle 2--4--3 */
1910: MLPutFunction(link, "Polygon", 1);
1911: MLPutFunction(link, "Part", 2);
1912: MLPutSymbol(link, "points");
1913: MLPutFunction(link, "Plus", 2);
1914: MLPutFunction(link, "Part", 2);
1915: MLPutFunction(link, "Slot", 1);
1916: MLPutInteger(link, 1);
1917: MLPutFunction(link, "List", 3);
1918: MLPutInteger(link, 3);
1919: MLPutInteger(link, 5);
1920: MLPutInteger(link, 4);
1921: MLPutInteger(link, 1);
1922: /* Triangle 3--4--5 */
1923: MLPutFunction(link, "Polygon", 1);
1924: MLPutFunction(link, "Part", 2);
1925: MLPutSymbol(link, "points");
1926: MLPutFunction(link, "Plus", 2);
1927: MLPutFunction(link, "Part", 2);
1928: MLPutFunction(link, "Slot", 1);
1929: MLPutInteger(link, 1);
1930: MLPutFunction(link, "List", 3);
1931: MLPutInteger(link, 4);
1932: MLPutInteger(link, 5);
1933: MLPutInteger(link, 6);
1934: MLPutInteger(link, 1);
1935: }
1936: else if ((numCorners == 3) || (!grid->cm->fieldClasses[grid->viewField][1]))
1937: {
1938: /* Triangle 0--1--2 */
1939: MLPutFunction(link, "Polygon", 1);
1940: MLPutFunction(link, "Part", 2);
1941: MLPutSymbol(link, "points");
1942: MLPutFunction(link, "Plus", 2);
1943: MLPutFunction(link, "Part", 2);
1944: MLPutFunction(link, "Slot", 1);
1945: MLPutInteger(link, 1);
1946: MLPutFunction(link, "List", 3);
1947: MLPutInteger(link, 1);
1948: MLPutInteger(link, 2);
1949: MLPutInteger(link, 3);
1950: MLPutInteger(link, 1);
1951: } else {
1952: SETERRQ(PETSC_ERR_ARG_WRONG, "Invalid number of local nodes");
1953: }
1954: MLPutSymbol(link, "faces");
1955: MLPutFunction(link, "Rule", 2);
1956: MLPutSymbol(link, "AspectRatio");
1957: MLPutReal(link, 1.0);
1958: return(0);
1959: #endif
1961: return(0);
1962: }
1964: int PetscViewerMathematicaVectorPlot_Triangular_2D(PetscViewer viewer, GVec v)
1965: {
1966: #ifdef PETSC_HAVE_MATHEMATICA
1967: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
1968: MLINK link = vmath->link; /* The link to Mathematica */
1969: Grid grid;
1970: Mesh mesh;
1971: double scale;
1972: PetscTruth opt;
1973: int ierr;
1976: GVecGetGrid(v, &grid);
1977: GridGetMesh(grid, &mesh);
1979: MLPutFunction(link, "ListPlotVectorField", 3);
1980: MLPutSymbol(link, "points");
1981: MLPutFunction(link, "Rule", 2);
1982: MLPutSymbol(link, "AspectRatio");
1983: MLPutReal(link, mesh->sizeY/mesh->sizeX);
1984: MLPutFunction(link, "Rule", 2);
1985: MLPutSymbol(link, "ScaleFactor");
1986: PetscOptionsGetDouble("viewer_", "-math_scale", &scale, &opt);
1987: if (opt == PETSC_TRUE) {
1988: MLPutReal(link, scale);
1989: } else {
1990: MLPutSymbol(link, "None");
1991: }
1992: return(0);
1993: #endif
1995: return(0);
1996: }
1998: int PetscViewerMathematicaSurfacePlot_Triangular_2D(PetscViewer viewer, GVec v)
1999: {
2000: #ifdef PETSC_HAVE_MATHEMATICA
2001: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
2002: MLINK link = vmath->link; /* The link to Mathematica */
2003: double startX, endX;
2004: double startY, endY;
2005: PetscTruth opt;
2006: int ierr;
2009: /* Setup domain */
2010: startX = 0.0;
2011: startY = 0.0;
2012: endX = 1.0;
2013: endY = 1.0;
2014: ierr = PetscOptionsGetDouble("viewer_", "-math_start_x", &startX, &opt);
2015: ierr = PetscOptionsGetDouble("viewer_", "-math_start_y", &startY, &opt);
2016: ierr = PetscOptionsGetDouble("viewer_", "-math_end_x", &endX, &opt);
2017: ierr = PetscOptionsGetDouble("viewer_", "-math_end_y", &endY, &opt);
2019: MLPutFunction(link, "Show", 1);
2020: MLPutFunction(link, "SurfaceGraphics", 6);
2021: MLPutSymbol(link, "points");
2022: MLPutFunction(link, "Rule", 2);
2023: MLPutSymbol(link, "ClipFill");
2024: MLPutSymbol(link, "None");
2025: MLPutFunction(link, "Rule", 2);
2026: MLPutSymbol(link, "Axes");
2027: MLPutSymbol(link, "True");
2028: MLPutFunction(link, "Rule", 2);
2029: MLPutSymbol(link, "PlotLabel");
2030: MLPutString(link, vmath->objName);
2031: MLPutFunction(link, "Rule", 2);
2032: MLPutSymbol(link, "MeshRange");
2033: MLPutFunction(link, "List", 2);
2034: MLPutFunction(link, "List", 2);
2035: MLPutReal(link, startX);
2036: MLPutReal(link, endX);
2037: MLPutFunction(link, "List", 2);
2038: MLPutReal(link, startY);
2039: MLPutReal(link, endY);
2040: MLPutFunction(link, "Rule", 2);
2041: MLPutSymbol(link, "AspectRatio");
2042: MLPutReal(link, (endY - startY)/(endX - startX));
2043: return(0);
2044: #endif
2046: return(0);
2047: }
2049: int PetscViewerMathematica_GVec_Triangular_2D(PetscViewer viewer, GVec v)
2050: {
2051: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
2052: #ifdef PETSC_HAVE_MATHEMATICA
2053: MLINK link = vmath->link; /* The link to Mathematica */
2054: Mesh mesh;
2055: Mesh_Triangular *tri;
2056: Grid grid;
2057: Mat P;
2058: GVec w;
2059: int numCorners;
2060: int ierr;
2063: ierr = GVecGetGrid(v, &grid);
2064: mesh = grid->mesh;
2065: tri = (Mesh_Triangular *) mesh->data;
2066: numCorners = tri->numCorners;
2068: /* Check that a field component has been specified */
2069: if ((grid->viewField < 0) || (grid->viewField >= grid->fields))
2070: return(0);
2072: if (grid->isConstrained) {
2073: GVecCreate(grid, &w);
2074: GridGetConstraintMatrix(grid, &P);
2075: MatMult(P, v, w);
2076: } else {
2077: w = v;
2078: }
2080: /* Check that the mesh is defined correctly */
2081: PetscViewerMathematicaCheckMesh_Triangular_2D(viewer, mesh);
2083: /* Send the first component of the first field */
2084: MLPutFunction(link, "EvaluatePacket", 1);
2085: switch(vmath->plotType)
2086: {
2087: case MATHEMATICA_TRIANGULATION_PLOT:
2088: MLPutFunction(link, "Module", 2);
2089: /* Make temporary points with each value of the field component in the vector */
2090: MLPutFunction(link, "List", 2);
2091: MLPutSymbol(link, "f");
2092: MLPutFunction(link, "Set", 2);
2093: MLPutSymbol(link, "points");
2094: PetscViewerMathematicaCreateSamplePoints_Triangular_2D(viewer, w);
2095: /* Display the picture */
2096: PetscViewerMathematicaTriangulationPlot_Triangular_2D(viewer, w);
2097: break;
2098: case MATHEMATICA_VECTOR_TRIANGULATION_PLOT:
2099: if (grid->comp[grid->viewField] != 2) {
2100: SETERRQ(PETSC_ERR_ARG_WRONG, "Field must be a 2D vector field for this plot type");
2101: }
2102: MLPutFunction(link, "Module", 2);
2103: /* Make temporary list with points and 2D vector field values */
2104: MLPutFunction(link, "List", 2);
2105: MLPutSymbol(link, "f");
2106: MLPutFunction(link, "Set", 2);
2107: MLPutSymbol(link, "points");
2108: PetscViewerMathematicaCreateVectorSamplePoints_Triangular_2D(viewer, w);
2109: /* Display the picture */
2110: PetscViewerMathematicaVectorPlot_Triangular_2D(viewer, w);
2111: break;
2112: case MATHEMATICA_VECTOR_PLOT:
2113: if (grid->comp[grid->viewField] != 2) {
2114: SETERRQ(PETSC_ERR_ARG_WRONG, "Field must be a 2D vector field for this plot type");
2115: }
2116: MLPutFunction(link, "Module", 2);
2117: /* Make temporary list with points and 2D vector field values */
2118: MLPutFunction(link, "List", 2);
2119: MLPutSymbol(link, "f");
2120: MLPutFunction(link, "Set", 2);
2121: MLPutSymbol(link, "points");
2122: PetscViewerMathematicaCreateInterpolatedSamplePoints_Triangular_2D(viewer, w, -1);
2123: /* Display the picture */
2124: PetscViewerMathematicaVectorPlot_Triangular_2D(viewer, w);
2125: break;
2126: case MATHEMATICA_SURFACE_PLOT:
2127: if (grid->comp[grid->viewField] != 2) {
2128: SETERRQ(PETSC_ERR_ARG_WRONG, "Field must be a 2D vector field for this plot type");
2129: }
2130: MLPutFunction(link, "Module", 2);
2131: /* Make temporary list with interpolated field values on a square mesh */
2132: MLPutFunction(link, "List", 1);
2133: MLPutFunction(link, "Set", 2);
2134: MLPutSymbol(link, "points");
2135: PetscViewerMathematicaCreateInterpolatedSamplePoints_Triangular_2D(viewer, w, grid->viewComp);
2136: /* Display the picture */
2137: PetscViewerMathematicaSurfacePlot_Triangular_2D(viewer, w);
2138: break;
2139: default:
2140: SETERRQ(PETSC_ERR_ARG_WRONG, "Invalid plot type");
2141: }
2142: MLEndPacket(link);
2143: /* Skip packets until ReturnPacket */
2144: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
2145: /* Skip ReturnPacket */
2146: MLNewPacket(link);
2148: if (grid->isConstrained) {
2149: VecDestroy(w);
2150: }
2151: #else
2153: switch(vmath->plotType)
2154: {
2155: case MATHEMATICA_TRIANGULATION_PLOT:
2156: break;
2157: case MATHEMATICA_VECTOR_TRIANGULATION_PLOT:
2158: break;
2159: case MATHEMATICA_VECTOR_PLOT:
2160: break;
2161: case MATHEMATICA_SURFACE_PLOT:
2162: break;
2163: default:
2164: SETERRQ(PETSC_ERR_ARG_WRONG, "Invalid plot type");
2165: }
2166: #endif
2167: return(0);
2168: }