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/ksp/pc/pcimpl.h
12: #include src/mat/impls/aij/seq/aij.h
13: #include mathematica.h
14: #include "petscfix.h"
16: #if defined (PETSC_HAVE__SNPRINTF) && !defined(PETSC_HAVE_SNPRINTF)
17: #define snprintf _snprintf
18: #endif
20: PetscViewer PETSC_VIEWER_MATHEMATICA_WORLD_PRIVATE = PETSC_NULL;
21: #ifdef PETSC_HAVE_MATHEMATICA
22: static void *mathematicaEnv = PETSC_NULL;
23: #endif
27: /*@C
28: PetscViewerMathematicaInitializePackage - This function initializes everything in the Petsc interface to Mathematica. It is
29: called from PetscDLLibraryRegister() when using dynamic libraries, and on the call to PetscInitialize()
30: when using static libraries.
32: Input Parameter:
33: path - The dynamic library path, or PETSC_NULL
35: Level: developer
37: .keywords: Petsc, initialize, package, PLAPACK
38: .seealso: PetscInitializePackage(), PetscInitialize()
39: @*/
40: int PetscViewerMathematicaInitializePackage(char *path) {
41: static PetscTruth initialized = PETSC_FALSE;
44: if (initialized == PETSC_TRUE) return(0);
45: initialized = PETSC_TRUE;
46: #ifdef PETSC_HAVE_MATHEMATICA
47: mathematicaEnv = (void *) MLInitialize(0);
48: #endif
49: return(0);
50: }
54: /*@C
55: PetscViewerMathematicaDestroyPackage - This function destroys everything in the Petsc interface to Mathematica. It is
56: called from PetscFinalize().
58: Level: developer
60: .keywords: Petsc, destroy, package, mathematica
61: .seealso: PetscFinalize()
62: @*/
63: int PetscViewerMathematicaFinalizePackage(void) {
65: #ifdef PETSC_HAVE_MATHEMATICA
66: if (mathematicaEnv != PETSC_NULL) MLDeinitialize((MLEnvironment) mathematicaEnv);
67: #endif
68: return(0);
69: }
73: int PetscViewerInitializeMathematicaWorld_Private()
74: {
78: if (PETSC_VIEWER_MATHEMATICA_WORLD_PRIVATE) return(0);
79: PetscViewerMathematicaOpen(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_NULL, PETSC_NULL, &PETSC_VIEWER_MATHEMATICA_WORLD_PRIVATE);
80:
81: return(0);
82: }
86: static int PetscViewerDestroy_Mathematica(PetscViewer viewer)
87: {
88: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
89: int ierr;
92: #ifdef PETSC_HAVE_MATHEMATICA
93: MLClose(vmath->link);
94: #endif
95: if (vmath->linkname != PETSC_NULL) {
96: PetscFree(vmath->linkname);
97: }
98: if (vmath->linkhost != PETSC_NULL) {
99: PetscFree(vmath->linkhost);
100: }
101: PetscFree(vmath);
102: return(0);
103: }
107: int PetscViewerDestroyMathematica_Private(void)
108: {
112: if (PETSC_VIEWER_MATHEMATICA_WORLD_PRIVATE) {
113: PetscViewerDestroy(PETSC_VIEWER_MATHEMATICA_WORLD_PRIVATE);
114: }
115: return(0);
116: }
120: int PetscViewerMathematicaSetupConnection_Private(PetscViewer v) {
121: #ifdef PETSC_HAVE_MATHEMATICA
122: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) v->data;
123: #ifdef MATHEMATICA_3_0
124: int argc = 6;
125: char *argv[6];
126: #else
127: int argc = 5;
128: char *argv[5];
129: #endif
130: char hostname[256];
131: long lerr;
132: int ierr;
133: #endif
136: #ifdef PETSC_HAVE_MATHEMATICA
137: /* Link name */
138: argv[0] = "-linkname";
139: if (vmath->linkname == PETSC_NULL) {
140: argv[1] = "math -mathlink";
141: } else {
142: argv[1] = vmath->linkname;
143: }
145: /* Link host */
146: argv[2] = "-linkhost";
147: if (vmath->linkhost == PETSC_NULL) {
148: PetscGetHostName(hostname, 255);
149: argv[3] = hostname;
150: } else {
151: argv[3] = vmath->linkhost;
152: }
154: /* Link mode */
155: #ifdef MATHEMATICA_3_0
156: argv[4] = "-linkmode";
157: switch(vmath->linkmode) {
158: case MATHEMATICA_LINK_CREATE:
159: argv[5] = "Create";
160: break;
161: case MATHEMATICA_LINK_CONNECT:
162: argv[5] = "Connect";
163: break;
164: case MATHEMATICA_LINK_LAUNCH:
165: argv[5] = "Launch";
166: break;
167: }
168: #else
169: switch(vmath->linkmode) {
170: case MATHEMATICA_LINK_CREATE:
171: argv[4] = "-linkcreate";
172: break;
173: case MATHEMATICA_LINK_CONNECT:
174: argv[4] = "-linkconnect";
175: break;
176: case MATHEMATICA_LINK_LAUNCH:
177: argv[4] = "-linklaunch";
178: break;
179: }
180: #endif
182: vmath->link = MLOpenInEnv(mathematicaEnv, argc, argv, &lerr);
183: #endif
185: return(0);
186: }
188: EXTERN_C_BEGIN
191: int PetscViewerCreate_Mathematica(PetscViewer v) {
192: PetscViewer_Mathematica *vmath;
193: int ierr;
197: PetscNew(PetscViewer_Mathematica, &vmath);
198: v->data = (void *) vmath;
199: v->ops->destroy = PetscViewerDestroy_Mathematica;
200: v->ops->flush = 0;
201: PetscStrallocpy(PETSC_VIEWER_MATHEMATICA, &v->type_name);
203: vmath->linkname = PETSC_NULL;
204: vmath->linkhost = PETSC_NULL;
205: vmath->linkmode = MATHEMATICA_LINK_CONNECT;
206: vmath->graphicsType = GRAPHICS_MOTIF;
207: vmath->plotType = MATHEMATICA_TRIANGULATION_PLOT;
208: vmath->objName = PETSC_NULL;
210: PetscViewerMathematicaSetFromOptions(v);
211: PetscViewerMathematicaSetupConnection_Private(v);
212: return(0);
213: }
214: EXTERN_C_END
218: int PetscViewerMathematicaParseLinkMode_Private(char *modename, LinkMode *mode) {
219: PetscTruth isCreate, isConnect, isLaunch;
220: int ierr;
223: PetscStrcasecmp(modename, "Create", &isCreate);
224: PetscStrcasecmp(modename, "Connect", &isConnect);
225: PetscStrcasecmp(modename, "Launch", &isLaunch);
226: if (isCreate == PETSC_TRUE) {
227: *mode = MATHEMATICA_LINK_CREATE;
228: } else if (isConnect == PETSC_TRUE) {
229: *mode = MATHEMATICA_LINK_CONNECT;
230: } else if (isLaunch == PETSC_TRUE) {
231: *mode = MATHEMATICA_LINK_LAUNCH;
232: } else {
233: SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid Mathematica link mode: %s", modename);
234: }
235: return(0);
236: }
240: int PetscViewerMathematicaSetFromOptions(PetscViewer v) {
241: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) v->data;
242: char linkname[256];
243: char modename[256];
244: char hostname[256];
245: char type[256];
246: int numPorts;
247: int *ports;
248: int numHosts;
249: char **hosts;
250: int size, rank;
251: int h;
252: PetscTruth opt;
253: int ierr;
256: MPI_Comm_size(v->comm, &size);
257: MPI_Comm_rank(v->comm, &rank);
259: /* Get link name */
260: PetscOptionsGetString("viewer_", "-math_linkname", linkname, 255, &opt);
261: if (opt == PETSC_TRUE) {
262: PetscViewerMathematicaSetLinkName(v, linkname);
263: }
264: /* Get link port */
265: numPorts = size;
266: PetscMalloc(size * sizeof(int), &ports);
267: PetscOptionsGetIntArray("viewer_", "-math_linkport", ports, &numPorts, &opt);
268: if (opt == PETSC_TRUE) {
269: if (numPorts > rank) {
270: snprintf(linkname, 255, "%6d", ports[rank]);
271: } else {
272: snprintf(linkname, 255, "%6d", ports[0]);
273: }
274: PetscViewerMathematicaSetLinkName(v, linkname);
275: }
276: PetscFree(ports);
277: /* Get link host */
278: numHosts = size;
279: PetscMalloc(size * sizeof(char *), &hosts);
280: PetscOptionsGetStringArray("viewer_", "-math_linkhost", hosts, &numHosts, &opt);
281: if (opt == PETSC_TRUE) {
282: if (numHosts > rank) {
283: PetscStrncpy(hostname, hosts[rank], 255);
284: } else {
285: PetscStrncpy(hostname, hosts[0], 255);
286: }
287: PetscViewerMathematicaSetLinkHost(v, hostname);
288: }
289: for(h = 0; h < numHosts; h++) {
290: PetscFree(hosts[h]);
291: }
292: PetscFree(hosts);
293: /* Get link mode */
294: PetscOptionsGetString("viewer_", "-math_linkmode", modename, 255, &opt);
295: if (opt == PETSC_TRUE) {
296: LinkMode mode;
298: PetscViewerMathematicaParseLinkMode_Private(modename, &mode);
299: PetscViewerMathematicaSetLinkMode(v, mode);
300: }
301: /* Get graphics type */
302: PetscOptionsGetString("viewer_", "-math_graphics", type, 255, &opt);
303: if (opt == PETSC_TRUE) {
304: PetscTruth isMotif, isPS, isPSFile;
306: PetscStrcasecmp(type, "Motif", &isMotif);
307: PetscStrcasecmp(type, "PS", &isPS);
308: PetscStrcasecmp(type, "PSFile", &isPSFile);
309: if (isMotif == PETSC_TRUE) {
310: vmath->graphicsType = GRAPHICS_MOTIF;
311: } else if (isPS == PETSC_TRUE) {
312: vmath->graphicsType = GRAPHICS_PS_STDOUT;
313: } else if (isPSFile == PETSC_TRUE) {
314: vmath->graphicsType = GRAPHICS_PS_FILE;
315: }
316: }
317: /* Get plot type */
318: PetscOptionsGetString("viewer_", "-math_type", type, 255, &opt);
319: if (opt == PETSC_TRUE) {
320: PetscTruth isTri, isVecTri, isVec, isSurface;
322: PetscStrcasecmp(type, "Triangulation", &isTri);
323: PetscStrcasecmp(type, "VectorTriangulation", &isVecTri);
324: PetscStrcasecmp(type, "Vector", &isVec);
325: PetscStrcasecmp(type, "Surface", &isSurface);
326: if (isTri == PETSC_TRUE) {
327: vmath->plotType = MATHEMATICA_TRIANGULATION_PLOT;
328: } else if (isVecTri == PETSC_TRUE) {
329: vmath->plotType = MATHEMATICA_VECTOR_TRIANGULATION_PLOT;
330: } else if (isVec == PETSC_TRUE) {
331: vmath->plotType = MATHEMATICA_VECTOR_PLOT;
332: } else if (isSurface == PETSC_TRUE) {
333: vmath->plotType = MATHEMATICA_SURFACE_PLOT;
334: }
335: }
336: return(0);
337: }
341: int PetscViewerMathematicaSetLinkName(PetscViewer v, const char *name) {
342: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) v->data;
343: int ierr;
348: PetscStrallocpy(name, &vmath->linkname);
349: return(0);
350: }
354: int PetscViewerMathematicaSetLinkPort(PetscViewer v, int port) {
355: char name[16];
356: int ierr;
359: snprintf(name, 16, "%6d", port);
360: PetscViewerMathematicaSetLinkName(v, name);
361: return(0);
362: }
366: int PetscViewerMathematicaSetLinkHost(PetscViewer v, const char *host) {
367: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) v->data;
368: int ierr;
373: PetscStrallocpy(host, &vmath->linkhost);
374: return(0);
375: }
379: int PetscViewerMathematicaSetLinkMode(PetscViewer v, LinkMode mode) {
380: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) v->data;
383: vmath->linkmode = mode;
384: return(0);
385: }
387: /*----------------------------------------- Public Functions --------------------------------------------------------*/
390: /*@C
391: PetscViewerMathematicaOpen - Communicates with Mathemtica using MathLink.
393: Collective on comm
395: Input Parameters:
396: + comm - The MPI communicator
397: . port - [optional] The port to connect on, or PETSC_DECIDE
398: . machine - [optional] The machine to run Mathematica on, or PETSC_NULL
399: - mode - [optional] The connection mode, or PETSC_NULL
401: Output Parameter:
402: . viewer - The Mathematica viewer
404: Level: intermediate
406: Notes:
407: Most users should employ the following commands to access the
408: Mathematica viewers
409: $
410: $ PetscViewerMathematicaOpen(MPI_Comm comm, int port, char *machine, char *mode, PetscViewer &viewer)
411: $ MatView(Mat matrix, PetscViewer viewer)
412: $
413: $ or
414: $
415: $ PetscViewerMathematicaOpen(MPI_Comm comm, int port, char *machine, char *mode, PetscViewer &viewer)
416: $ VecView(Vec vector, PetscViewer viewer)
418: Options Database Keys:
419: $ -viewer_math_linkhost <machine> - The host machine for the kernel
420: $ -viewer_math_linkname <name> - The full link name for the connection
421: $ -viewer_math_linkport <port> - The port for the connection
422: $ -viewer_math_mode <mode> - The mode, e.g. Launch, Connect
423: $ -viewer_math_type <type> - The plot type, e.g. Triangulation, Vector
424: $ -viewer_math_graphics <output> - The output type, e.g. Motif, PS, PSFile
426: .keywords: PetscViewer, Mathematica, open
428: .seealso: MatView(), VecView()
429: @*/
430: int PetscViewerMathematicaOpen(MPI_Comm comm, int port, const char machine[], const char mode[], PetscViewer *v) {
431: int ierr;
434: PetscViewerCreate(comm, v);
435: #if 0
436: LinkMode linkmode;
437: PetscViewerMathematicaSetLinkPort(*v, port);
438: PetscViewerMathematicaSetLinkHost(*v, machine);
439: PetscViewerMathematicaParseLinkMode_Private(mode, &linkmode);
440: PetscViewerMathematicaSetLinkMode(*v, linkmode);
441: #endif
442: PetscViewerSetType(*v, PETSC_VIEWER_MATHEMATICA);
443: return(0);
444: }
446: #ifdef PETSC_HAVE_MATHEMATICA
449: /*@C
450: PetscViewerMathematicaGetLink - Returns the link to Mathematica
452: Input Parameters:
453: . viewer - The Mathematica viewer
454: . link - The link to Mathematica
456: Level: intermediate
458: .keywords PetscViewer, Mathematica, link
459: .seealso PetscViewerMathematicaOpen()
460: @*/
461: int PetscViewerMathematicaGetLink(PetscViewer viewer, MLINK *link)
462: {
463: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
467: *link = vmath->link;
468: return(0);
469: }
470: #endif
474: /*@C
475: PetscViewerMathematicaSkipPackets - Discard packets sent by Mathematica until a certain packet type is received
477: Input Parameters:
478: . viewer - The Mathematica viewer
479: . type - The packet type to search for, e.g RETURNPKT
481: Level: advanced
483: .keywords PetscViewer, Mathematica, packets
484: .seealso PetscViewerMathematicaSetName(), PetscViewerMathematicaGetVector()
485: @*/
486: int PetscViewerMathematicaSkipPackets(PetscViewer viewer, int type)
487: {
488: #ifdef PETSC_HAVE_MATHEMATICA
489: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
490: MLINK link = vmath->link; /* The link to Mathematica */
491: int pkt; /* The packet type */
494: while((pkt = MLNextPacket(link)) && (pkt != type))
495: MLNewPacket(link);
496: if (!pkt) {
497: MLClearError(link);
498: SETERRQ(PETSC_ERR_LIB, (char *) MLErrorMessage(link));
499: }
500: return(0);
501: #endif
503: return(0);
504: }
508: /*@C
509: PetscViewerMathematicaLoadGraphics - Change the type of graphics output for Mathematica
511: Input Parameters:
512: . viewer - The Mathematica viewer
513: . type - The type of graphics, e.g. GRAPHICS_MOTIF, GRAPHICS_PS_FILE, GRAPHICS_PS_STDOUT
515: Level: intermediate
517: .keywords PetscViewer, Mathematica, packets
518: .seealso PetscViewerMathematicaSetName(), PetscViewerMathematicaSkipPackets()
519: @*/
520: int PetscViewerMathematicaLoadGraphics(PetscViewer viewer, GraphicsType type)
521: {
522: #ifdef PETSC_HAVE_MATHEMATICA
523: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
524: MLINK link = vmath->link; /* The link to Mathematica */
525: char programName[256];
526: int ierr;
529: /* Load graphics package */
530: MLPutFunction(link, "EvaluatePacket", 1);
531: switch(type)
532: {
533: case GRAPHICS_MOTIF:
534: MLPutFunction(link, "Get", 1);
535: MLPutString(link, "Motif.m");
536: break;
537: case GRAPHICS_PS_FILE:
538: MLPutFunction(link, "CompoundExpression", 4);
539: MLPutFunction(link, "Set", 2);
540: MLPutSymbol(link, "PetscGraphicsCounter");
541: MLPutInteger(link, 0);
542: MLPutFunction(link, "SetDelayed", 2);
543: MLPutSymbol(link, "$Display");
544: MLPutSymbol(link, "$FileDisplay");
545: MLPutFunction(link, "Set", 2);
546: MLPutSymbol(link, "$FileDisplay");
547: MLPutFunction(link, "OpenWrite", 1);
548: MLPutFunction(link, "StringJoin", 3);
549: if (!PetscGetProgramName(programName, 255))
550: MLPutString(link, programName);
551: else
552: MLPutString(link, "GVec");
553: MLPutFunction(link, "ToString", 1);
554: MLPutSymbol(link, "PetscGraphicsCounter");
555: MLPutString(link, ".ps");
556: MLPutFunction(link, "Set", 2);
557: MLPutSymbol(link, "$DisplayFunction");
558: MLPutFunction(link, "Function", 1);
559: MLPutFunction(link, "CompoundExpression", 2);
560: MLPutFunction(link, "Display", 3);
561: MLPutSymbol(link, "$Display");
562: MLPutFunction(link, "Slot", 1);
563: MLPutInteger(link, 1);
564: MLPutString(link, "EPS");
565: MLPutFunction(link, "Increment", 1);
566: MLPutSymbol(link, "PetscGraphicsCounter");
567: break;
568: case GRAPHICS_PS_STDOUT:
569: MLPutFunction(link, "Get", 1);
570: MLPutString(link, "PSDirect.m");
571: break;
572: default:
573: SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid graphics type: %d", type);
574: }
575: MLEndPacket(link);
577: /* Skip packets until ReturnPacket */
578: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
579: /* Skip ReturnPacket */
580: MLNewPacket(link);
582: /* Load PlotField.m for vector plots */
583: MLPutFunction(link, "EvaluatePacket", 1);
584: MLPutFunction(link, "Get", 1);
585: MLPutString(link, "Graphics/PlotField.m");
586: MLEndPacket(link);
588: /* Skip packets until ReturnPacket */
589: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
590: /* Skip ReturnPacket */
591: MLNewPacket(link);
593: return(0);
594: #else
596: switch(type)
597: {
598: case GRAPHICS_MOTIF:
599: case GRAPHICS_PS_FILE:
600: case GRAPHICS_PS_STDOUT:
601: break;
602: default:
603: SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid graphics type: %d", type);
604: }
605: return(0);
606: #endif
607: }
611: /*@C
612: PetscViewerMathematicaGetName - Retrieve the default name for objects communicated to Mathematica
614: Input Parameter:
615: . viewer - The Mathematica viewer
617: Output Parameter:
618: . name - The name for new objects created in Mathematica
620: Level: intermediate
622: .keywords PetscViewer, Mathematica, name
623: .seealso PetscViewerMathematicaSetName(), PetscViewerMathematicaClearName()
624: @*/
625: int PetscViewerMathematicaGetName(PetscViewer viewer, const char **name)
626: {
627: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
632: *name = vmath->objName;
633: return(0);
634: }
638: /*@C
639: PetscViewerMathematicaSetName - Override the default name for objects communicated to Mathematica
641: Input Parameters:
642: . viewer - The Mathematica viewer
643: . name - The name for new objects created in Mathematica
645: Level: intermediate
647: .keywords PetscViewer, Mathematica, name
648: .seealso PetscViewerMathematicaSetName(), PetscViewerMathematicaClearName()
649: @*/
650: int PetscViewerMathematicaSetName(PetscViewer viewer, const char name[])
651: {
652: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
657: vmath->objName = name;
658: return(0);
659: }
663: /*@C
664: PetscViewerMathematicaClearName - Use the default name for objects communicated to Mathematica
666: Input Parameter:
667: . viewer - The Mathematica viewer
669: Level: intermediate
671: .keywords PetscViewer, Mathematica, name
672: .seealso PetscViewerMathematicaGetName(), PetscViewerMathematicaSetName()
673: @*/
674: int PetscViewerMathematicaClearName(PetscViewer viewer)
675: {
676: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
680: vmath->objName = PETSC_NULL;
681: return(0);
682: }
686: /*@C
687: PetscViewerMathematicaGetVector - Retrieve a vector from Mathematica
689: Input Parameter:
690: . viewer - The Mathematica viewer
692: Output Parameter:
693: . v - The vector
695: Level: intermediate
697: .keywords PetscViewer, Mathematica, vector
698: .seealso VecView(), PetscViewerMathematicaPutVector()
699: @*/
700: int PetscViewerMathematicaGetVector(PetscViewer viewer, Vec v) {
701: #ifdef PETSC_HAVE_MATHEMATICA
702: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
703: MLINK link; /* The link to Mathematica */
704: char *name;
705: PetscScalar *mArray;
706: PetscScalar *array;
707: long mSize;
708: int size;
709: int ierr;
715: /* Determine the object name */
716: if (vmath->objName == PETSC_NULL) {
717: name = "vec";
718: } else {
719: name = (char *) vmath->objName;
720: }
722: link = vmath->link;
723: VecGetLocalSize(v, &size);
724: VecGetArray(v, &array);
725: MLPutFunction(link, "EvaluatePacket", 1);
726: MLPutSymbol(link, name);
727: MLEndPacket(link);
728: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
729: MLGetRealList(link, &mArray, &mSize);
730: if (size != mSize) SETERRQ(PETSC_ERR_ARG_WRONG, "Incompatible vector sizes");
731: PetscMemcpy(array, mArray, mSize * sizeof(double));
732: MLDisownRealList(link, mArray, mSize);
733: VecRestoreArray(v, &array);
735: return(0);
736: #endif
738: return(0);
739: }
743: /*@C
744: PetscViewerMathematicaPutVector - Send a vector to Mathematica
746: Input Parameters:
747: + viewer - The Mathematica viewer
748: - v - The vector
750: Level: intermediate
752: .keywords PetscViewer, Mathematica, vector
753: .seealso VecView(), PetscViewerMathematicaGetVector()
754: @*/
755: int PetscViewerMathematicaPutVector(PetscViewer viewer, Vec v) {
756: #ifdef PETSC_HAVE_MATHEMATICA
757: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
758: MLINK link = vmath->link; /* The link to Mathematica */
759: char *name;
760: PetscScalar *array;
761: int size;
762: int ierr;
765: /* Determine the object name */
766: if (vmath->objName == PETSC_NULL) {
767: name = "vec";
768: } else {
769: name = (char *) vmath->objName;
770: }
771: VecGetLocalSize(v, &size);
772: VecGetArray(v, &array);
774: /* Send the Vector object */
775: MLPutFunction(link, "EvaluatePacket", 1);
776: MLPutFunction(link, "Set", 2);
777: MLPutSymbol(link, name);
778: MLPutRealList(link, array, size);
779: MLEndPacket(link);
780: /* Skip packets until ReturnPacket */
781: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
782: /* Skip ReturnPacket */
783: MLNewPacket(link);
785: VecRestoreArray(v, &array);
786: return(0);
787: #endif
789: return(0);
790: }
792: int PetscViewerMathematicaPutMatrix(PetscViewer viewer, int m, int n, PetscReal *a) {
793: #ifdef PETSC_HAVE_MATHEMATICA
794: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
795: MLINK link = vmath->link; /* The link to Mathematica */
796: char *name;
797: int ierr;
800: /* Determine the object name */
801: if (vmath->objName == PETSC_NULL) {
802: name = "mat";
803: } else {
804: name = (char *) vmath->objName;
805: }
807: /* Send the dense matrix object */
808: MLPutFunction(link, "EvaluatePacket", 1);
809: MLPutFunction(link, "Set", 2);
810: MLPutSymbol(link, name);
811: MLPutFunction(link, "Transpose", 1);
812: MLPutFunction(link, "Partition", 2);
813: MLPutRealList(link, a, m*n);
814: MLPutInteger(link, m);
815: MLEndPacket(link);
816: /* Skip packets until ReturnPacket */
817: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
818: /* Skip ReturnPacket */
819: MLNewPacket(link);
821: return(0);
822: #endif
824: return(0);
825: }
827: int PetscViewerMathematicaPutCSRMatrix(PetscViewer viewer, int m, int n, int *i, int *j, PetscReal *a) {
828: #ifdef PETSC_HAVE_MATHEMATICA
829: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
830: MLINK link = vmath->link; /* The link to Mathematica */
831: const char *symbol;
832: char *name;
833: PetscTruth match;
834: int ierr;
837: /* Determine the object name */
838: if (vmath->objName == PETSC_NULL) {
839: name = "mat";
840: } else {
841: name = (char *) vmath->objName;
842: }
844: /* Make sure Mathematica recognizes sparse matrices */
845: MLPutFunction(link, "EvaluatePacket", 1);
846: MLPutFunction(link, "Needs", 1);
847: MLPutString(link, "LinearAlgebra`CSRMatrix`");
848: MLEndPacket(link);
849: /* Skip packets until ReturnPacket */
850: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
851: /* Skip ReturnPacket */
852: MLNewPacket(link);
854: /* Send the CSRMatrix object */
855: MLPutFunction(link, "EvaluatePacket", 1);
856: MLPutFunction(link, "Set", 2);
857: MLPutSymbol(link, name);
858: MLPutFunction(link, "CSRMatrix", 5);
859: MLPutInteger(link, m);
860: MLPutInteger(link, n);
861: MLPutFunction(link, "Plus", 2);
862: MLPutIntegerList(link, i, m+1);
863: MLPutInteger(link, 1);
864: MLPutFunction(link, "Plus", 2);
865: MLPutIntegerList(link, j, i[m]);
866: MLPutInteger(link, 1);
867: MLPutRealList(link, a, i[m]);
868: MLEndPacket(link);
869: /* Skip packets until ReturnPacket */
870: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
871: /* Skip ReturnPacket */
872: MLNewPacket(link);
874: /* Check that matrix is valid */
875: MLPutFunction(link, "EvaluatePacket", 1);
876: MLPutFunction(link, "ValidQ", 1);
877: MLPutSymbol(link, name);
878: MLEndPacket(link);
879: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
880: MLGetSymbol(link, &symbol);
881: PetscStrcmp("True", (char *) symbol, &match);
882: if (match == PETSC_FALSE) {
883: MLDisownSymbol(link, symbol);
884: SETERRQ(PETSC_ERR_PLIB, "Invalid CSR matrix in Mathematica");
885: }
886: MLDisownSymbol(link, symbol);
887: /* Skip ReturnPacket */
888: MLNewPacket(link);
890: return(0);
891: #endif
893: return(0);
894: }
896: /*------------------------------------------- ML Functions ----------------------------------------------------------*/
899: int PetscViewerMathematicaPartitionMesh(PetscViewer viewer, int numElements, int numVertices, double *vertices, int **mesh,
900: int *numParts, int *colPartition)
901: {
902: #ifdef PETSC_HAVE_MATHEMATICA
903: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
904: MLINK link; /* The link to Mathematica */
905: const char *symbol;
906: int numOptions;
907: int partSize;
908: int *part;
909: long numCols;
910: int col;
911: PetscTruth match, opt;
912: int ierr;
916: link = vmath->link;
918: /* Make sure that the reduce.m package is loaded */
919: MLPutFunction(link, "EvaluatePacket", 1);
920: MLPutFunction(link, "Get", 1);
921: MLPutFunction(link, "StringJoin", 2);
922: MLPutFunction(link, "Environment", 1);
923: MLPutString(link, "PETSC_DIR");
924: MLPutString(link, "/src/pc/impls/ml/reduce.m");
925: MLEndPacket(link);
926: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
927: MLGetSymbol(link, &symbol);
928: PetscStrcmp("$Failed", (char *) symbol, &match);
929: if (match == PETSC_TRUE) {
930: MLDisownSymbol(link, symbol);
931: SETERRQ(PETSC_ERR_FILE_OPEN, "Unable to load package reduce.m");
932: }
933: MLDisownSymbol(link, symbol);
935: /* Partition the mesh */
936: numOptions = 0;
937: partSize = 0;
938: PetscOptionsGetInt(PETSC_NULL, "-pc_ml_partition_size", &partSize, &opt);
939: if ((opt == PETSC_TRUE) && (partSize > 0))
940: numOptions++;
941: MLPutFunction(link, "EvaluatePacket", 1);
942: MLPutFunction(link, "PartitionMesh", 1 + numOptions);
943: MLPutFunction(link, "MeshData", 5);
944: MLPutInteger(link, numElements);
945: MLPutInteger(link, numVertices);
946: MLPutInteger(link, numVertices);
947: MLPutFunction(link, "Partition", 2);
948: MLPutRealList(link, vertices, numVertices*2);
949: MLPutInteger(link, 2);
950: MLPutFunction(link, "Partition", 2);
951: MLPutIntegerList(link, mesh[MESH_ELEM], numElements*3);
952: MLPutInteger(link, 3);
953: if (partSize > 0)
954: {
955: MLPutFunction(link, "Rule", 2);
956: MLPutSymbol(link, "PartitionSize");
957: MLPutInteger(link, partSize);
958: }
959: MLEndPacket(link);
960: /* Skip packets until ReturnPacket */
961: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
963: /* Get the vertex partiton */
964: MLGetIntegerList(link, &part, &numCols);
965: if (numCols != numVertices) SETERRQ(PETSC_ERR_PLIB, "Invalid partition");
966: for(col = 0, *numParts = 0; col < numCols; col++) {
967: colPartition[col] = part[col]-1;
968: *numParts = PetscMax(*numParts, part[col]);
969: }
970: MLDisownIntegerList(link, part, numCols);
972: return(0);
973: #endif
975: return(0);
976: }
980: int PetscViewerMathematicaReduce(PetscViewer viewer, PC pc, int thresh)
981: {
982: #ifdef PETSC_HAVE_MATHEMATICA
983: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
984: MLINK link; /* The link to Mathematica */
985: PC_Multilevel *ml;
986: int *range;
987: long numRange;
988: int *null;
989: long numNull;
990: const char *symbol;
991: int numOptions;
992: int partSize;
993: int row, col;
994: PetscTruth match, opt;
995: int ierr;
1000: link = vmath->link;
1001: ml = (PC_Multilevel *) pc->data;
1003: /* Make sure that the reduce.m package is loaded */
1004: MLPutFunction(link, "EvaluatePacket", 1);
1005: MLPutFunction(link, "Get", 1);
1006: MLPutFunction(link, "StringJoin", 2);
1007: MLPutFunction(link, "Environment", 1);
1008: MLPutString(link, "PETSC_DIR");
1009: MLPutString(link, "/src/pc/impls/ml/reduce.m");
1010: MLEndPacket(link);
1011: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1012: MLGetSymbol(link, &symbol);
1013: PetscStrcmp("$Failed", (char *) symbol, &match);
1014: if (match == PETSC_TRUE) {
1015: MLDisownSymbol(link, symbol);
1016: SETERRQ(PETSC_ERR_FILE_OPEN, "Unable to load package reduce.m");
1017: }
1018: MLDisownSymbol(link, symbol);
1020: /* mesh = MeshData[] */
1021: MLPutFunction(link, "EvaluatePacket", 1);
1022: MLPutFunction(link, "Set", 2);
1023: MLPutSymbol(link, "mesh");
1024: MLPutFunction(link, "MeshData", 5);
1025: MLPutInteger(link, ml->numElements[0]);
1026: MLPutInteger(link, ml->numVertices[0]);
1027: MLPutInteger(link, ml->numVertices[0]);
1028: MLPutFunction(link, "Partition", 2);
1029: MLPutRealList(link, ml->vertices, ml->numVertices[0]*2);
1030: MLPutInteger(link, 2);
1031: MLPutFunction(link, "Partition", 2);
1032: MLPutIntegerList(link, ml->meshes[0][MESH_ELEM], ml->numElements[0]*3);
1033: MLPutInteger(link, 3);
1034: MLEndPacket(link);
1035: /* Skip packets until ReturnPacket */
1036: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1037: /* Skip ReturnPacket */
1038: MLNewPacket(link);
1039: /* Check that mesh is valid */
1040: MLPutFunction(link, "EvaluatePacket", 1);
1041: MLPutFunction(link, "ValidQ", 1);
1042: MLPutSymbol(link, "mesh");
1043: MLEndPacket(link);
1044: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1045: MLGetSymbol(link, &symbol);
1046: PetscStrcmp("True", (char *) symbol, &match);
1047: if (match == PETSC_FALSE) {
1048: MLDisownSymbol(link, symbol);
1049: SETERRQ(PETSC_ERR_PLIB, "Invalid mesh in Mathematica");
1050: }
1051: MLDisownSymbol(link, symbol);
1053: /* mat = gradient matrix */
1054: MatView(ml->locB, viewer);
1056: /* mattML = ReduceMatrix[mat,ml->minNodes] */
1057: numOptions = 0;
1058: partSize = 0;
1059: PetscOptionsGetInt(PETSC_NULL, "-pc_ml_partition_size", &partSize, &opt);
1060: if ((opt == PETSC_TRUE) && (partSize > 0))
1061: numOptions++;
1062: MLPutFunction(link, "EvaluatePacket", 1);
1063: MLPutFunction(link, "Set", 2);
1064: MLPutSymbol(link, "mattML");
1065: MLPutFunction(link, "ReduceMatrix", 3 + numOptions);
1066: MLPutSymbol(link, "mesh");
1067: MLPutSymbol(link, "mat");
1068: MLPutInteger(link, thresh);
1069: if (partSize > 0) {
1070: MLPutFunction(link, "Rule", 2);
1071: MLPutSymbol(link, "PartitionSize");
1072: MLPutInteger(link, partSize);
1073: }
1074: MLEndPacket(link);
1075: /* Skip packets until ReturnPacket */
1076: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1077: /* Skip ReturnPacket */
1078: MLNewPacket(link);
1079: /* Check that mattML is valid */
1080: MLPutFunction(link, "EvaluatePacket", 1);
1081: MLPutFunction(link, "ValidQ", 1);
1082: MLPutSymbol(link, "mattML");
1083: MLEndPacket(link);
1084: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1085: MLGetSymbol(link, &symbol);
1086: PetscStrcmp("True", (char *) symbol, &match);
1087: if (match == PETSC_FALSE) {
1088: MLDisownSymbol(link, symbol);
1089: SETERRQ(PETSC_ERR_PLIB, "Invalid MLData in Mathematica");
1090: }
1091: MLDisownSymbol(link, symbol);
1093: /* Copy information to the preconditioner */
1094: MLPutFunction(link, "EvaluatePacket", 1);
1095: MLPutFunction(link, "Part", 2);
1096: MLPutSymbol(link, "mattML");
1097: MLPutInteger(link, 3);
1098: MLEndPacket(link);
1099: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1100: MLGetInteger(link, &ml->numLevels);
1102: /* Create lists of the range and nullspace columns */
1103: MLPutFunction(link, "EvaluatePacket", 1);
1104: MLPutFunction(link, "GetRange", 1);
1105: MLPutSymbol(link, "mattML");
1106: MLEndPacket(link);
1107: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1108: MLGetIntegerList(link, &range, &numRange);
1109: if (numRange > ml->sOrder->numLocVars) SETERRQ(PETSC_ERR_PLIB, "Invalid size for range space");
1110: ml->rank = numRange;
1111: ml->globalRank = ml->rank;
1112: if (ml->rank > 0) {
1113: PetscMalloc(numRange * sizeof(int), &ml->range);
1114: for(row = 0; row < numRange; row++)
1115: ml->range[row] = range[row]-1;
1116: }
1117: MLDisownIntegerList(link, range, numRange);
1119: MLPutFunction(link, "EvaluatePacket", 1);
1120: MLPutFunction(link, "GetNullColumns", 1);
1121: MLPutSymbol(link, "mattML");
1122: MLEndPacket(link);
1123: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1124: MLGetIntegerList(link, &null, &numNull);
1125: if (numRange + numNull != ml->sOrder->numLocVars) SETERRQ(PETSC_ERR_PLIB, "Invalid size for range and null spaces");
1126: ml->numLocNullCols = numNull;
1127: if (numNull > 0)
1128: {
1129: PetscMalloc(numNull * sizeof(int), &ml->nullCols);
1130: for(col = 0; col < numNull; col++)
1131: ml->nullCols[col] = null[col] - 1;
1132: }
1133: MLDisownIntegerList(link, null, numNull);
1135: return(0);
1136: #endif
1138: return(0);
1139: }
1143: int PetscViewerMathematicaMultiLevelConvert(PetscViewer viewer, PC pc)
1144: {
1145: #ifdef PETSC_HAVE_MATHEMATICA
1146: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
1147: MLINK link; /* The link to Mathematica */
1148: PC_Multilevel *ml;
1149: Mat_SeqAIJ *grad;
1150: int *numPartitions;
1151: int *numPartitionCols, *cols;
1152: int *numPartitionRows, *rows;
1153: double *U, *singVals, *V;
1154: long *Udims, *Vdims;
1155: char **Uheads, **Vheads;
1156: int *nnz;
1157: int *offsets;
1158: double *vals;
1159: long numLevels, numParts, numCols, numRows, Udepth, numSingVals, Vdepth, len;
1160: int numBdRows, numBdCols;
1161: int mesh, level, part, col, row;
1162: int ierr;
1167: link = vmath->link;
1168: ml = (PC_Multilevel *) pc->data;
1170: /* ml->numLevels = ml[[3]] */
1171: MLPutFunction(link, "EvaluatePacket", 1);
1172: MLPutFunction(link, "Part", 2);
1173: MLPutSymbol(link, "mattML");
1174: MLPutInteger(link, 3);
1175: MLEndPacket(link);
1176: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1177: MLGetInteger(link, &ml->numLevels);
1179: /* ml->numMeshes = Length[ml[[4]]] */
1180: MLPutFunction(link, "EvaluatePacket", 1);
1181: MLPutFunction(link, "Length", 1);
1182: MLPutFunction(link, "Part", 2);
1183: MLPutSymbol(link, "mattML");
1184: MLPutInteger(link, 4);
1185: MLEndPacket(link);
1186: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1187: MLGetInteger(link, &ml->numMeshes);
1189: /* ml->numElements[level] = ml[[4,level,1]] */
1190: PetscMalloc(ml->numMeshes * sizeof(int), &ml->numElements);
1192: /* ml->numVertices[level] = ml[[4,level,2]] */
1193: PetscMalloc(ml->numMeshes * sizeof(int), &ml->numVertices);
1195: /* ml->meshes = ml[[4]] */
1196: PetscMalloc(ml->numMeshes * sizeof(int **), &ml->meshes);
1197: for(mesh = 0; mesh < ml->numMeshes; mesh++) {
1198: PetscMalloc(NUM_MESH_DIV * sizeof(int *), &ml->meshes[mesh]);
1199: /* Here we should get meshes */
1200: PetscMalloc(1 * sizeof(int), &ml->meshes[mesh][MESH_OFFSETS]);
1201: PetscMalloc(1 * sizeof(int), &ml->meshes[mesh][MESH_ADJ]);
1202: PetscMalloc(1 * sizeof(int), &ml->meshes[mesh][MESH_ELEM]);
1203: }
1205: /* ml->numPartitions = Map[Length,Drop[ml[[5]],-1]] */
1206: MLPutFunction(link, "EvaluatePacket", 1);
1207: MLPutFunction(link, "Map", 2);
1208: MLPutSymbol(link, "Length");
1209: MLPutFunction(link, "Drop", 2);
1210: MLPutFunction(link, "Part", 2);
1211: MLPutSymbol(link, "mattML");
1212: MLPutInteger(link, 5);
1213: MLPutInteger(link, -1);
1214: MLEndPacket(link);
1215: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1216: MLGetIntegerList(link, &numPartitions, &numLevels);
1217: if (numLevels != ml->numLevels) SETERRQ(PETSC_ERR_PLIB, "Invalid node partition in MLData object");
1218: if (numLevels > 0) {
1219: PetscMalloc(numLevels * sizeof(int), &ml->numPartitions);
1220: PetscMemcpy(ml->numPartitions, numPartitions, numLevels * sizeof(int));
1221: }
1222: MLDisownIntegerList(link, numPartitions, numLevels);
1224: if (ml->numLevels > 0) {
1225: /* ml->numPartitionCols = Map[Length,ml[[5,level]]] */
1226: PetscMalloc(ml->numLevels * sizeof(int *), &ml->numPartitionCols);
1227: for(level = 0; level < ml->numLevels; level++) {
1228: MLPutFunction(link, "EvaluatePacket", 1);
1229: MLPutFunction(link, "Map", 2);
1230: MLPutSymbol(link, "Length");
1231: MLPutFunction(link, "Part", 3);
1232: MLPutSymbol(link, "mattML");
1233: MLPutInteger(link, 5);
1234: MLPutInteger(link, level+1);
1235: MLEndPacket(link);
1236: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1237: MLGetIntegerList(link, &numPartitionCols, &numParts);
1238: if (numParts != ml->numPartitions[level]) SETERRQ(PETSC_ERR_PLIB, "Invalid node partition in MLData object");
1239: if (numParts > 0) {
1240: PetscMalloc(numParts * sizeof(int), &ml->numPartitionCols[level]);
1241: PetscMemcpy(ml->numPartitionCols[level], numPartitionCols, numParts * sizeof(int));
1242: }
1243: MLDisownIntegerList(link, numPartitionCols, numParts);
1244: }
1246: /* ml->colPartition[level][part] = ml[[5,level,part]] */
1247: PetscMalloc(ml->numLevels * sizeof(int **), &ml->colPartition);
1248: for(level = 0; level < ml->numLevels; level++) {
1249: if (ml->numPartitions[level] == 0) continue;
1250: PetscMalloc(ml->numPartitions[level] * sizeof(int *), &ml->colPartition[level]);
1251: for(part = 0; part < ml->numPartitions[level]; part++) {
1252: MLPutFunction(link, "EvaluatePacket", 1);
1253: MLPutFunction(link, "Part", 4);
1254: MLPutSymbol(link, "mattML");
1255: MLPutInteger(link, 5);
1256: MLPutInteger(link, level+1);
1257: MLPutInteger(link, part+1);
1258: MLEndPacket(link);
1259: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1260: MLGetIntegerList(link, &cols, &numCols);
1261: if (numCols != ml->numPartitionCols[level][part]) SETERRQ(PETSC_ERR_PLIB, "Invalid node partition in MLData object");
1262: if (numCols > 0) {
1263: PetscMalloc(numCols * sizeof(int), &ml->colPartition[level][part]);
1264: /* Convert to zero-based indices */
1265: for(col = 0; col < numCols; col++) ml->colPartition[level][part][col] = cols[col] - 1;
1266: }
1267: MLDisownIntegerList(link, cols, numCols);
1268: }
1269: }
1271: /* ml->numPartitionRows = Map[Length,FlattenAt[ml[[6,level]],1]] */
1272: PetscMalloc(ml->numLevels * sizeof(int *), &ml->numPartitionRows);
1273: for(level = 0; level < ml->numLevels; level++) {
1274: MLPutFunction(link, "EvaluatePacket", 1);
1275: MLPutFunction(link, "Map", 2);
1276: MLPutSymbol(link, "Length");
1277: MLPutFunction(link, "FlattenAt", 2);
1278: MLPutFunction(link, "Part", 3);
1279: MLPutSymbol(link, "mattML");
1280: MLPutInteger(link, 6);
1281: MLPutInteger(link, level+1);
1282: MLPutInteger(link, 1);
1283: MLEndPacket(link);
1284: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1285: MLGetIntegerList(link, &numPartitionRows, &numParts);
1286: if (numParts != ml->numPartitions[level] + NUM_PART_ROW_DIV - 1) {
1287: SETERRQ(PETSC_ERR_PLIB, "Invalid edge partition in MLData object");
1288: }
1289: PetscMalloc(numParts * sizeof(int), &ml->numPartitionRows[level]);
1290: PetscMemcpy(ml->numPartitionRows[level], numPartitionRows, numParts * sizeof(int));
1291: MLDisownIntegerList(link, numPartitionRows, numParts);
1292: }
1294: /* ml->rowPartition[level][PART_ROW_INT][part] = ml[[6,level,1,part]]
1295: ml->rowPartition[level][PART_ROW_BD] = ml[[6,level,2]]
1296: ml->rowPartition[level][PART_ROW_RES] = ml[[6,level,3]] */
1297: PetscMalloc(ml->numLevels * sizeof(int ***), &ml->rowPartition);
1298: for(level = 0; level < ml->numLevels; level++) {
1299: PetscMalloc(NUM_PART_ROW_DIV * sizeof(int **), &ml->rowPartition[level]);
1300: /* Interior rows */
1301: if (ml->numPartitions[level] > 0) {
1302: PetscMalloc(ml->numPartitions[level] * sizeof(int *), &ml->rowPartition[level][PART_ROW_INT]);
1303: for(part = 0; part < ml->numPartitions[level]; part++) {
1304: MLPutFunction(link, "EvaluatePacket", 1);
1305: MLPutFunction(link, "Part", 5);
1306: MLPutSymbol(link, "mattML");
1307: MLPutInteger(link, 6);
1308: MLPutInteger(link, level+1);
1309: MLPutInteger(link, 1);
1310: MLPutInteger(link, part+1);
1311: MLEndPacket(link);
1312: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1313: MLGetIntegerList(link, &rows, &numRows);
1314: if (numRows != ml->numPartitionRows[level][part]) {
1315: SETERRQ(PETSC_ERR_PLIB, "Invalid edge partition in MLData object");
1316: }
1317: if (numRows > 0) {
1318: PetscMalloc(numRows * sizeof(int), &ml->rowPartition[level][PART_ROW_INT][part]);
1319: /* Convert to zero-based indices */
1320: for(row = 0; row < numRows; row++) {
1321: ml->rowPartition[level][PART_ROW_INT][part][row] = rows[row] - 1;
1322: }
1323: }
1324: MLDisownIntegerList(link, rows, numRows);
1325: }
1326: }
1327: /* Boundary rows */
1328: PetscMalloc(1 * sizeof(int *), &ml->rowPartition[level][PART_ROW_BD]);
1329: MLPutFunction(link, "EvaluatePacket", 1);
1330: MLPutFunction(link, "Part", 4);
1331: MLPutSymbol(link, "mattML");
1332: MLPutInteger(link, 6);
1333: MLPutInteger(link, level+1);
1334: MLPutInteger(link, 2);
1335: MLEndPacket(link);
1336: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1337: MLGetIntegerList(link, &rows, &numRows);
1338: if (numRows != ml->numPartitionRows[level][ml->numPartitions[level]]) {
1339: SETERRQ(PETSC_ERR_PLIB, "Invalid edge partition in MLData object");
1340: }
1341: if (numRows > 0) {
1342: PetscMalloc(numRows * sizeof(int), &ml->rowPartition[level][PART_ROW_BD][0]);
1343: /* Convert to zero-based indices */
1344: for(row = 0; row < numRows; row++) {
1345: ml->rowPartition[level][PART_ROW_BD][0][row] = rows[row] - 1;
1346: }
1347: }
1348: MLDisownIntegerList(link, rows, numRows);
1349: /* Residual rows*/
1350: PetscMalloc(1 * sizeof(int *), &ml->rowPartition[level][PART_ROW_RES]);
1351: MLPutFunction(link, "EvaluatePacket", 1);
1352: MLPutFunction(link, "Part", 4);
1353: MLPutSymbol(link, "mattML");
1354: MLPutInteger(link, 6);
1355: MLPutInteger(link, level+1);
1356: MLPutInteger(link, 3);
1357: MLEndPacket(link);
1358: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1359: MLGetIntegerList(link, &rows, &numRows);
1360: if (numRows != ml->numPartitionRows[level][ml->numPartitions[level]+1]) {
1361: SETERRQ(PETSC_ERR_PLIB, "Invalid edge partition in MLData object");
1362: }
1363: if (numRows > 0) {
1364: PetscMalloc(numRows * sizeof(int), &ml->rowPartition[level][PART_ROW_RES][0]);
1365: /* Convert to zero-based indices */
1366: for(row = 0; row < numRows; row++) {
1367: ml->rowPartition[level][PART_ROW_RES][0][row] = rows[row] - 1;
1368: }
1369: }
1370: MLDisownIntegerList(link, rows, numRows);
1371: }
1372: } else {
1373: PetscMalloc(1 * sizeof(int), &ml->numPartitions);
1374: PetscMalloc(1 * sizeof(int *), &ml->numPartitionCols);
1375: PetscMalloc(1 * sizeof(int **), &ml->colPartition);
1376: PetscMalloc(1 * sizeof(int *), &ml->numPartitionRows);
1377: PetscMalloc(1 * sizeof(int ***), &ml->rowPartition);
1378: }
1380: /* ml->numRows = ml[[1]] */
1381: MLPutFunction(link, "EvaluatePacket", 1);
1382: MLPutFunction(link, "Part", 2);
1383: MLPutSymbol(link, "mattML");
1384: MLPutInteger(link, 1);
1385: MLEndPacket(link);
1386: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1387: MLGetInteger(link, &ml->numRows);
1389: /* ml->numCols = ml[[2]] */
1390: MLPutFunction(link, "EvaluatePacket", 1);
1391: MLPutFunction(link, "Part", 2);
1392: MLPutSymbol(link, "mattML");
1393: MLPutInteger(link, 2);
1394: MLEndPacket(link);
1395: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1396: MLGetInteger(link, &ml->numCols);
1398: /* ml->zeroTol = ml[[9]] */
1399: MLPutFunction(link, "EvaluatePacket", 1);
1400: MLPutFunction(link, "N", 1);
1401: MLPutFunction(link, "Part", 2);
1402: MLPutSymbol(link, "mattML");
1403: MLPutInteger(link, 9);
1404: MLEndPacket(link);
1405: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1406: MLGetReal(link, &ml->zeroTol);
1408: if (ml->numLevels > 0) {
1409: /* ml->factors[level][part][FACT_U] = ml[[7,level,part,1]]
1410: ml->factors[level][part][FACT_DINV] = Divide[1,Select[ml[[7,level,part,2]],(#>ml[[9]])&]]
1411: ml->factors[level][part][FACT_V] = ml[[7,level,part,3]] */
1412: PetscMalloc(ml->numLevels * sizeof(double ***), &ml->factors);
1413: for(level = 0; level < ml->numLevels; level++) {
1414: if (ml->numPartitions[level] == 0) continue;
1415: PetscMalloc(ml->numPartitions[level] * sizeof(double **), &ml->factors[level]);
1416: for(part = 0; part < ml->numPartitions[level]; part++) {
1417: PetscMalloc(NUM_FACT_DIV * sizeof(double *), &ml->factors[level][part]);
1418: /* U, or U^T in LAPACK terms */
1419: MLPutFunction(link, "EvaluatePacket", 1);
1420: MLPutFunction(link, "Part", 5);
1421: MLPutSymbol(link, "mattML");
1422: MLPutInteger(link, 7);
1423: MLPutInteger(link, level+1);
1424: MLPutInteger(link, part+1);
1425: MLPutInteger(link, 1);
1426: MLEndPacket(link);
1427: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1428: MLGetRealArray(link, &U, &Udims, &Uheads, &Udepth);
1429: if (Udepth > 1) {
1430: if (Udepth != 2) SETERRQ(PETSC_ERR_PLIB, "Invalid U matrix");
1431: if ((Udims[0] != ml->numPartitionRows[level][part]) || (Udims[0] != Udims[1])) {
1432: SETERRQ(PETSC_ERR_PLIB, "Incompatible dimensions for U matrix");
1433: }
1434: PetscMalloc(Udims[0]*Udims[0] * sizeof(double), &ml->factors[level][part][FACT_U]);
1435: /* Notice that LAPACK will think that this is U^T, or U in LAPACK terms */
1436: PetscMemcpy(ml->factors[level][part][FACT_U], U, Udims[0]*Udims[0] * sizeof(double));
1437: } else if (ml->numPartitionRows[level][part] != 0) {
1438: SETERRQ(PETSC_ERR_PLIB, "Missing U matrix");
1439: }
1440: MLDisownRealArray(link, U, Udims, Uheads, Udepth);
1441: /* D^{-1} */
1442: MLPutFunction(link, "EvaluatePacket", 1);
1443: MLPutFunction(link, "Divide", 2);
1444: MLPutReal(link, 1.0);
1445: MLPutFunction(link, "Select", 2);
1446: MLPutFunction(link, "Part", 5);
1447: MLPutSymbol(link, "mattML");
1448: MLPutInteger(link, 7);
1449: MLPutInteger(link, level+1);
1450: MLPutInteger(link, part+1);
1451: MLPutInteger(link, 2);
1452: MLPutFunction(link, "Function", 2);
1453: MLPutSymbol(link, "x");
1454: MLPutFunction(link, "Greater", 2);
1455: MLPutSymbol(link, "x");
1456: MLPutFunction(link, "Part", 2);
1457: MLPutSymbol(link, "mattML");
1458: MLPutInteger(link, 9);
1459: MLEndPacket(link);
1460: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1461: MLGetRealList(link, &singVals, &numSingVals);
1462: if (numSingVals > ml->numPartitionCols[level][part]) {
1463: SETERRQ(PETSC_ERR_PLIB, "Invalid factor list in MLData object");
1464: }
1465: if (ml->numPartitionCols[level][part] > 0) {
1466: PetscMalloc(ml->numPartitionCols[level][part] * sizeof(double), &ml->factors[level][part][FACT_DINV]);
1467: PetscMemzero(ml->factors[level][part][FACT_DINV], ml->numPartitionCols[level][part] * sizeof(double));
1468: PetscMemcpy(ml->factors[level][part][FACT_DINV], singVals, numSingVals * sizeof(double));
1469: }
1470: MLDisownRealList(link, singVals, numSingVals);
1471: /* V^T, but V in LAPACK terms */
1472: MLPutFunction(link, "EvaluatePacket", 1);
1473: MLPutFunction(link, "Transpose", 1);
1474: MLPutFunction(link, "Part", 5);
1475: MLPutSymbol(link, "mattML");
1476: MLPutInteger(link, 7);
1477: MLPutInteger(link, level+1);
1478: MLPutInteger(link, part+1);
1479: MLPutInteger(link, 3);
1480: MLEndPacket(link);
1481: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1482: MLGetRealArray(link, &V, &Vdims, &Vheads, &Vdepth);
1483: if (Vdepth > 1) {
1484: if (Vdepth != 2) SETERRQ(PETSC_ERR_PLIB, "Invalid V matrix");
1485: if ((Vdims[0] != ml->numPartitionCols[level][part]) || (Vdims[0] != Vdims[1])) {
1486: SETERRQ(PETSC_ERR_PLIB, "Incompatible dimensions for U matrix");
1487: }
1488: PetscMalloc(Vdims[0]*Vdims[0] * sizeof(double), &ml->factors[level][part][FACT_V]);
1489: /* Notice that LAPACK will think that this is V, or V^T in LAPACK terms */
1490: PetscMemcpy(ml->factors[level][part][FACT_V], V, Vdims[0]*Vdims[0] * sizeof(double));
1491: } else if (ml->numPartitionCols[level][part] != 0) {
1492: SETERRQ(PETSC_ERR_PLIB, "Missing V matrix");
1493: }
1494: MLDisownRealArray(link, V, Vdims, Vheads, Vdepth);
1495: }
1496: }
1498: /* ml->grads = ml[[8]] */
1499: PetscMalloc(ml->numLevels * sizeof(Mat), &ml->grads);
1500: PetscMalloc(ml->numLevels * sizeof(Vec), &ml->bdReduceVecs);
1501: PetscMalloc(ml->numLevels * sizeof(Vec), &ml->colReduceVecs);
1502: PetscMalloc(ml->numLevels * sizeof(Vec), &ml->colReduceVecs2);
1503: for(level = 0; level < ml->numLevels; level++) {
1504: if (ml->numPartitions[level] == 0) continue;
1505: /* numRows = ml[[8,level,1]] */
1506: MLPutFunction(link, "EvaluatePacket", 1);
1507: MLPutFunction(link, "Part", 4);
1508: MLPutSymbol(link, "mattML");
1509: MLPutInteger(link, 8);
1510: MLPutInteger(link, level+1);
1511: MLPutInteger(link, 1);
1512: MLEndPacket(link);
1513: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1514: MLGetInteger(link, &numBdRows);
1515: VecCreateSeq(PETSC_COMM_SELF, numBdRows, &ml->bdReduceVecs[level]);
1516: /* numCols = ml[[8,level,2]] */
1517: MLPutFunction(link, "EvaluatePacket", 1);
1518: MLPutFunction(link, "Part", 4);
1519: MLPutSymbol(link, "mattML");
1520: MLPutInteger(link, 8);
1521: MLPutInteger(link, level+1);
1522: MLPutInteger(link, 2);
1523: MLEndPacket(link);
1524: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1525: MLGetInteger(link, &numBdCols);
1526: VecCreateSeq(PETSC_COMM_SELF, numBdCols, &ml->colReduceVecs[level]);
1527: VecDuplicate(ml->colReduceVecs[level], &ml->colReduceVecs2[level]);
1528: /* nnz = Map[Apply[Subtract,Sort[#,Greater]]&, Partition[ml[[8,level,3]],2,1]] */
1529: MLPutFunction(link, "EvaluatePacket", 1);
1530: MLPutFunction(link, "Map", 2);
1531: MLPutFunction(link, "Function", 2);
1532: MLPutSymbol(link, "x");
1533: MLPutFunction(link, "Apply", 2);
1534: MLPutSymbol(link, "Subtract");
1535: MLPutFunction(link, "Sort", 2);
1536: MLPutSymbol(link, "x");
1537: MLPutSymbol(link, "Greater");
1538: MLPutFunction(link, "Partition", 3);
1539: MLPutFunction(link, "Part", 4);
1540: MLPutSymbol(link, "mattML");
1541: MLPutInteger(link, 8);
1542: MLPutInteger(link, level+1);
1543: MLPutInteger(link, 3);
1544: MLPutInteger(link, 2);
1545: MLPutInteger(link, 1);
1546: MLEndPacket(link);
1547: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1548: MLGetIntegerList(link, &nnz, &len);
1549: if (len != numBdRows) SETERRQ(PETSC_ERR_PLIB, "Invalid boundary gradient matrix");
1550: MatCreate(PETSC_COMM_SELF,numBdRows,numBdCols,numBdRows,numBdCols,&ml->grads[level]);
1551: MatSetType(ml->grads[level],MATSEQAIJ);
1552: MatSeqAIJSetPreallocation(ml->grads[level],PETSC_DEFAULT,nnz);
1553: grad = (Mat_SeqAIJ *) ml->grads[level]->data;
1554: PetscMemcpy(grad->ilen, nnz, numBdRows * sizeof(int));
1555: MLDisownIntegerList(link, nnz, len);
1556: /* ml->grads[level]->i = ml[[8,level,3]] */
1557: MLPutFunction(link, "EvaluatePacket", 1);
1558: MLPutFunction(link, "Part", 4);
1559: MLPutSymbol(link, "mattML");
1560: MLPutInteger(link, 8);
1561: MLPutInteger(link, level+1);
1562: MLPutInteger(link, 3);
1563: MLEndPacket(link);
1564: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1565: MLGetIntegerList(link, &offsets, &len);
1566: if (len != numBdRows+1) SETERRQ(PETSC_ERR_PLIB, "Invalid boundary gradient matrix");
1567: for(row = 0; row <= numBdRows; row++)
1568: grad->i[row] = offsets[row]-1;
1569: MLDisownIntegerList(link, offsets, len);
1570: /* ml->grads[level]->j = ml[[8,level,4]] */
1571: MLPutFunction(link, "EvaluatePacket", 1);
1572: MLPutFunction(link, "Part", 4);
1573: MLPutSymbol(link, "mattML");
1574: MLPutInteger(link, 8);
1575: MLPutInteger(link, level+1);
1576: MLPutInteger(link, 4);
1577: MLEndPacket(link);
1578: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1579: MLGetIntegerList(link, &cols, &len);
1580: if (len != grad->i[numBdRows]) SETERRQ(PETSC_ERR_PLIB, "Invalid boundary gradient matrix");
1581: for(col = 0; col < len; col++)
1582: grad->j[col] = cols[col]-1;
1583: MLDisownIntegerList(link, cols, len);
1584: /* ml->grads[level]->i = ml[[8,level,5]] */
1585: MLPutFunction(link, "EvaluatePacket", 1);
1586: MLPutFunction(link, "Part", 4);
1587: MLPutSymbol(link, "mattML");
1588: MLPutInteger(link, 8);
1589: MLPutInteger(link, level+1);
1590: MLPutInteger(link, 5);
1591: MLEndPacket(link);
1592: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1593: MLGetRealList(link, &vals, &len);
1594: if (len != grad->i[numBdRows]) SETERRQ(PETSC_ERR_PLIB, "Invalid boundary gradient matrix");
1595: PetscMemcpy(grad->a, vals, len * sizeof(double));
1596: MLDisownRealList(link, vals, len);
1597: /* Fix up all the members */
1598: grad->nz += len;
1599: MatAssemblyBegin(ml->grads[level], MAT_FINAL_ASSEMBLY);
1600: MatAssemblyEnd(ml->grads[level], MAT_FINAL_ASSEMBLY);
1601: }
1602: } else {
1603: PetscMalloc(1 * sizeof(double ***), &ml->factors);
1604: PetscMalloc(1 * sizeof(Mat), &ml->grads);
1605: PetscMalloc(1 * sizeof(Vec), &ml->bdReduceVecs);
1606: PetscMalloc(1 * sizeof(Vec), &ml->colReduceVecs);
1607: PetscMalloc(1 * sizeof(Vec), &ml->colReduceVecs2);
1608: }
1610: ml->interiorWorkLen = 1;
1611: for(level = 0; level < ml->numLevels; level++) {
1612: for(part = 0; part < ml->numPartitions[level]; part++)
1613: ml->interiorWorkLen = PetscMax(ml->interiorWorkLen, ml->numPartitionRows[level][part]);
1614: }
1615: PetscMalloc(ml->interiorWorkLen * sizeof(double), &ml->interiorWork);
1616: PetscMalloc(ml->interiorWorkLen * sizeof(double), &ml->interiorWork2);
1618: return(0);
1619: #endif
1621: return(0);
1622: }
1624: #if 0
1625: /*------------------------------ Functions for Triangular 2d Meshes -------------------------------------------------*/
1628: int PetscViewerMathematicaCreateSamplePoints_Triangular_2D(PetscViewer viewer, GVec v)
1629: {
1630: #ifdef PETSC_HAVE_MATHEMATICA
1631: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
1632: MLINK link = vmath->link; /* The link to Mathematica */
1633: Grid grid;
1634: int numNodes;
1635: int *classes;
1636: int *offsets;
1637: double *array;
1638: int *localStart;
1639: int localOffset, comp;
1640: int node, nclass;
1641: int ierr;
1644: GVecGetGrid(v, &grid);
1645: numNodes = grid->mesh->numNodes;
1646: comp = grid->viewComp;
1647: offsets = grid->order->offsets;
1648: localStart = grid->order->localStart[grid->viewField];
1649: classes = grid->cm->classes;
1651: /* Attach a value to each point in the mesh -- Extra values are put in for fields not
1652: defined on some nodes, but these values are never used */
1653: VecGetArray(v, &array);
1654: MLPutFunction(link, "ReplaceAll", 2);
1655: MLPutFunction(link, "Thread", 1);
1656: MLPutFunction(link, "f", 2);
1657: MLPutSymbol(link, "nodes");
1658: MLPutFunction(link, "List", numNodes);
1659: for(node = 0; node < numNodes; node++)
1660: {
1661: nclass = classes[node];
1662: localOffset = localStart[nclass] + comp;
1663: MLPutReal(link, array[offsets[node] + localOffset]);
1664: }
1665: MLPutFunction(link, "Rule", 2);
1666: MLPutSymbol(link, "f");
1667: MLPutSymbol(link, "Append");
1668: VecRestoreArray(v, &array);
1670: return(0);
1671: #endif
1673: return(0);
1674: }
1678: int PetscViewerMathematicaCreateVectorSamplePoints_Triangular_2D(PetscViewer viewer, GVec v)
1679: {
1680: #ifdef PETSC_HAVE_MATHEMATICA
1681: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
1682: MLINK link = vmath->link; /* The link to Mathematica */
1683: Grid grid;
1684: int numNodes;
1685: int *classes;
1686: int *offsets;
1687: int *fieldClasses;
1688: double *array;
1689: int *localStart;
1690: int localOffset;
1691: int node, nclass;
1692: int ierr;
1695: GVecGetGrid(v, &grid);
1696: numNodes = grid->mesh->numNodes;
1697: fieldClasses = grid->cm->fieldClasses[grid->viewField];
1698: offsets = grid->order->offsets;
1699: localStart = grid->order->localStart[grid->viewField];
1700: classes = grid->cm->classes;
1702: /* Make a list {{{x_0,y_0},{f^0_x,f^0_y}},...} */
1703: VecGetArray(v, &array);
1704: MLPutFunction(link, "ReplaceAll", 2);
1705: MLPutFunction(link, "Thread", 1);
1706: MLPutFunction(link, "f", 2);
1707: MLPutSymbol(link, "nodes");
1708: MLPutFunction(link, "List", numNodes);
1709: for(node = 0; node < numNodes; node++)
1710: {
1711: nclass = classes[node];
1712: if (fieldClasses[nclass])
1713: {
1714: localOffset = localStart[nclass];
1715: MLPutFunction(link, "List", 2);
1716: MLPutReal(link, array[offsets[node] + localOffset]);
1717: MLPutReal(link, array[offsets[node] + localOffset + 1]);
1718: }
1719: else
1720: {
1721: /* Vectors of length zero are ignored */
1722: MLPutFunction(link, "List", 2);
1723: MLPutReal(link, 0.0);
1724: MLPutReal(link, 0.0);
1725: }
1726: }
1727: MLPutFunction(link, "Rule", 2);
1728: MLPutSymbol(link, "f");
1729: MLPutSymbol(link, "List");
1730: VecRestoreArray(v, &array);
1732: return(0);
1733: #endif
1735: return(0);
1736: }
1740: int PetscViewerMathematicaCreateInterpolatedSamplePoints_Triangular_2D(PetscViewer viewer, GVec v, int vComp)
1741: {
1742: #ifdef PETSC_HAVE_MATHEMATICA
1743: #if 0
1744: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
1745: MLINK link = vmath->link; /* The link to Mathematica */
1746: InterpolationContext iCtx;
1747: Grid grid;
1748: Mesh mesh;
1749: double *x, *y, *values;
1750: long m, n;
1751: double startX, endX, incX;
1752: double startY, endY, incY;
1753: int comp;
1754: int proc, row, col;
1755: PetscTruth opt;
1756: int ierr;
1757: #endif
1760: #if 0
1761: GVecGetGrid(v, &grid);
1762: GridGetMesh(grid, &mesh);
1763: comp = grid->comp[grid->viewField];
1765: /* This sucks, but I will fix it later (It is for GridReduceInterpolationElementVec_Triangular_2D) */
1766: grid->classesOld = grid->classes;
1768: /* Setup InterpolationContext */
1769: iCtx.vec = v;
1770: iCtx.mesh = mesh;
1771: iCtx.order = grid->order;
1772: iCtx.ghostVec = grid->ghostVec;
1773: iCtx.field = grid->viewField;
1774: iCtx.numProcs = mesh->part->numProcs;
1775: iCtx.rank = mesh->part->rank;
1776: PetscMalloc(iCtx.numProcs * sizeof(int), &iCtx.activeProcs);
1777: PetscMalloc(iCtx.numProcs * sizeof(int), &iCtx.calcProcs);
1778: PetscMalloc(iCtx.numProcs*3 * sizeof(PetscScalar), &iCtx.coords);
1779: PetscMalloc(iCtx.numProcs * sizeof(PetscScalar *), &iCtx.values);
1780: for(proc = 0; proc < iCtx.numProcs; proc++) {
1781: PetscMalloc(comp * sizeof(PetscScalar), &iCtx.values[proc]);
1782: }
1784: /* Setup domain */
1785: startX = 0.0;
1786: startY = 0.0;
1787: endX = 1.0;
1788: endY = 1.0;
1789: PetscOptionsGetDouble("viewer_", "-math_start_x", &startX, &opt);
1790: PetscOptionsGetDouble("viewer_", "-math_start_y", &startY, &opt);
1791: PetscOptionsGetDouble("viewer_", "-math_end_x", &endX, &opt);
1792: PetscOptionsGetDouble("viewer_", "-math_end_y", &endY, &opt);
1793: PetscOptionsGetInt("viewer_", "-math_div_x", (int *) &n, &opt);
1794: PetscOptionsGetInt("viewer_", "-math_div_y", (int *) &m, &opt);
1795: PetscMalloc((n+1) * sizeof(double), &x);
1796: PetscMalloc((n+1) * sizeof(double), &y);
1797: PetscMalloc((n+1)*comp * sizeof(double), &values);
1798: incX = (endX - startX)/n;
1799: incY = (endY - startY)/m;
1801: x[0] = startX;
1802: for(col = 1; col <= n; col++)
1803: x[col] = x[col-1] + incX;
1805: MLPutFunction(link, "List", m+1);
1806: for(row = 0; row <= m; row++)
1807: {
1808: PetscMemzero(values, (n+1)*comp * sizeof(double));
1809: for(col = 0; col <= n; col++)
1810: y[col] = startY + incY*row;
1811: PointFunctionInterpolateField(n+1, comp, x, y, x, values, &iCtx);
1812: if (vComp >= 0)
1813: {
1814: for(col = 0; col <= n; col++)
1815: values[col] = values[col*comp+vComp];
1816: MLPutRealList(link, values, n+1);
1817: }
1818: else
1819: {
1820: MLPutFunction(link, "Transpose", 1);
1821: MLPutFunction(link, "List", 2);
1822: MLPutFunction(link, "Transpose", 1);
1823: MLPutFunction(link, "List", 2);
1824: MLPutRealList(link, x, n+1);
1825: MLPutRealList(link, y, n+1);
1826: MLPutFunction(link, "Partition", 2);
1827: MLPutRealList(link, values, (n+1)*comp);
1828: MLPutInteger(link, comp);
1829: }
1830: }
1832: /* This sucks, but I will fix it later (It is for GridReduceInterpolationElementVec_Triangular_2D) */
1833: grid->classesOld = PETSC_NULL;
1835: /* Cleanup */
1836: PetscFree(x);
1837: PetscFree(y);
1838: PetscFree(values);
1839: PetscFree(iCtx.activeProcs);
1840: PetscFree(iCtx.calcProcs);
1841: PetscFree(iCtx.coords);
1842: for(proc = 0; proc < iCtx.numProcs; proc++) {
1843: PetscFree(iCtx.values[proc]);
1844: }
1845: PetscFree(iCtx.values);
1846: #endif
1848: return(0);
1849: #endif
1851: return(0);
1852: }
1856: int PetscViewerMathematica_Mesh_Triangular_2D(PetscViewer viewer, Mesh mesh)
1857: {
1858: #ifdef PETSC_HAVE_MATHEMATICA
1859: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
1860: MLINK link = vmath->link;
1861: Mesh_Triangular *tri = (Mesh_Triangular *) mesh->data;
1862: int numCorners = mesh->numCorners;
1863: int numFaces = mesh->numFaces;
1864: int *faces = tri->faces;
1865: int numNodes = mesh->numNodes;
1866: double *nodes = tri->nodes;
1867: int node, face, corner;
1868: int ierr;
1871: /* Load package to view graphics in a window */
1872: PetscViewerMathematicaLoadGraphics(viewer, vmath->graphicsType);
1874: /* Send the node coordinates */
1875: MLPutFunction(link, "EvaluatePacket", 1);
1876: MLPutFunction(link, "Set", 2);
1877: MLPutSymbol(link, "nodes");
1878: MLPutFunction(link, "List", numNodes);
1879: for(node = 0; node < numNodes; node++) {
1880: MLPutFunction(link, "List", 2);
1881: MLPutReal(link, nodes[node*2]);
1882: MLPutReal(link, nodes[node*2+1]);
1883: }
1884: MLEndPacket(link);
1885: /* Skip packets until ReturnPacket */
1886: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1887: /* Skip ReturnPacket */
1888: MLNewPacket(link);
1890: /* Send the faces */
1891: MLPutFunction(link, "EvaluatePacket", 1);
1892: MLPutFunction(link, "Set", 2);
1893: MLPutSymbol(link, "faces");
1894: MLPutFunction(link, "List", numFaces);
1895: for(face = 0; face < numFaces; face++) {
1896: MLPutFunction(link, "List", numCorners);
1897: for(corner = 0; corner < numCorners; corner++) {
1898: MLPutReal(link, faces[face*numCorners+corner]);
1899: }
1900: }
1901: MLEndPacket(link);
1902: /* Skip packets until ReturnPacket */
1903: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1904: /* Skip ReturnPacket */
1905: MLNewPacket(link);
1907: return(0);
1908: #endif
1910: return(0);
1911: }
1915: int PetscViewerMathematicaCheckMesh_Triangular_2D(PetscViewer viewer, Mesh mesh)
1916: {
1917: #ifdef PETSC_HAVE_MATHEMATICA
1918: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
1919: MLINK link = vmath->link;
1920: int numCorners = mesh->numCorners;
1921: int numFaces = mesh->numFaces;
1922: int numNodes = mesh->numNodes;
1923: const char *symbol;
1924: long args;
1925: int dim, type;
1926: PetscTruth match;
1927: int ierr;
1930: /* Check that nodes are defined */
1931: MLPutFunction(link, "EvaluatePacket", 1);
1932: MLPutFunction(link, "ValueQ", 1);
1933: MLPutSymbol(link, "nodes");
1934: MLEndPacket(link);
1935: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1936: MLGetSymbol(link, &symbol);
1937: PetscStrcmp("True", (char *) symbol, &match);
1938: if (match == PETSC_FALSE) {
1939: MLDisownSymbol(link, symbol);
1940: goto redefineMesh;
1941: }
1942: MLDisownSymbol(link, symbol);
1944: /* Check the dimensions */
1945: MLPutFunction(link, "EvaluatePacket", 1);
1946: MLPutFunction(link, "MatrixQ", 2);
1947: MLPutSymbol(link, "nodes");
1948: MLPutSymbol(link, "NumberQ");
1949: MLEndPacket(link);
1950: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1951: MLGetSymbol(link, &symbol);
1952: PetscStrcmp("True", (char *) symbol, &match);
1953: if (match == PETSC_FALSE) {
1954: MLDisownSymbol(link, symbol);
1955: goto redefineMesh;
1956: }
1957: MLDisownSymbol(link, symbol);
1958: MLPutFunction(link, "EvaluatePacket", 1);
1959: MLPutFunction(link, "Dimensions", 1);
1960: MLPutSymbol(link, "nodes");
1961: MLEndPacket(link);
1962: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1963: args = 0;
1964: type = MLGetNext(link);
1965: MLGetArgCount(link, &args);
1966: if (args != 2) {
1967: MLNewPacket(link);
1968: goto redefineMesh;
1969: }
1970: MLGetSymbol(link, &symbol);
1971: MLDisownSymbol(link, symbol);
1972: MLGetInteger(link, &dim);
1973: if (dim != numNodes) {
1974: MLNewPacket(link);
1975: goto redefineMesh;
1976: }
1977: MLGetInteger(link, &dim);
1978: if (dim != mesh->dim) {
1979: MLNewPacket(link);
1980: goto redefineMesh;
1981: }
1983: /* Check that faces are defined */
1984: MLPutFunction(link, "EvaluatePacket", 1);
1985: MLPutFunction(link, "ValueQ", 1);
1986: MLPutSymbol(link, "faces");
1987: MLEndPacket(link);
1988: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1989: MLGetSymbol(link, &symbol);
1990: PetscStrcmp("True", (char *) symbol, &match);
1991: if (match == PETSC_FALSE) {
1992: MLDisownSymbol(link, symbol);
1993: goto redefineMesh;
1994: }
1995: MLDisownSymbol(link, symbol);
1997: /* Check the dimensions */
1998: MLPutFunction(link, "EvaluatePacket", 1);
1999: MLPutFunction(link, "MatrixQ", 2);
2000: MLPutSymbol(link, "faces");
2001: MLPutSymbol(link, "NumberQ");
2002: MLEndPacket(link);
2003: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
2004: MLGetSymbol(link, &symbol);
2005: PetscStrcmp("True", (char *) symbol, &match);
2006: if (match == PETSC_FALSE) {
2007: MLDisownSymbol(link, symbol);
2008: goto redefineMesh;
2009: }
2010: MLDisownSymbol(link, symbol);
2011: MLPutFunction(link, "EvaluatePacket", 1);
2012: MLPutFunction(link, "Dimensions", 1);
2013: MLPutSymbol(link, "faces");
2014: MLEndPacket(link);
2015: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
2016: args = 0;
2017: type = MLGetNext(link);
2018: MLGetArgCount(link, &args);
2019: if (args != 2) {
2020: MLNewPacket(link);
2021: goto redefineMesh;
2022: }
2023: MLGetSymbol(link, &symbol);
2024: MLDisownSymbol(link, symbol);
2025: MLGetInteger(link, &dim);
2026: if (dim != numFaces) {
2027: MLNewPacket(link);
2028: goto redefineMesh;
2029: }
2030: MLGetInteger(link, &dim);
2031: if (dim != numCorners) {
2032: MLNewPacket(link);
2033: goto redefineMesh;
2034: }
2036: return(0);
2038: redefineMesh:
2039: MeshView(mesh, viewer);
2040: return(0);
2041: #endif
2043: return(0);
2044: }
2048: int PetscViewerMathematicaTriangulationPlot_Triangular_2D(PetscViewer viewer, GVec v)
2049: {
2050: #ifdef PETSC_HAVE_MATHEMATICA
2051: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
2052: MLINK link = vmath->link; /* The link to Mathematica */
2053: Mesh mesh;
2054: Grid grid;
2055: int numCorners;
2056: int ierr;
2059: GVecGetGrid(v, &grid);
2060: mesh = grid->mesh;
2061: numCorners = mesh->numCorners;
2063: MLPutFunction(link, "Show", 2);
2064: MLPutFunction(link, "Graphics3D", 1);
2065: MLPutFunction(link, "Flatten", 1);
2066: MLPutFunction(link, "Map", 2);
2067: MLPutFunction(link, "Function", 1);
2068: if ((numCorners == 6) && ((grid->cm->numClasses == 1) || (grid->cm->fieldClasses[grid->viewField][1])))
2069: {
2070: MLPutFunction(link, "List", 4);
2071: /* Triangle 0--5--4 */
2072: MLPutFunction(link, "Polygon", 1);
2073: MLPutFunction(link, "Part", 2);
2074: MLPutSymbol(link, "points");
2075: MLPutFunction(link, "Plus", 2);
2076: MLPutFunction(link, "Part", 2);
2077: MLPutFunction(link, "Slot", 1);
2078: MLPutInteger(link, 1);
2079: MLPutFunction(link, "List", 3);
2080: MLPutInteger(link, 1);
2081: MLPutInteger(link, 6);
2082: MLPutInteger(link, 5);
2083: MLPutInteger(link, 1);
2084: /* Triangle 1--3--5 */
2085: MLPutFunction(link, "Polygon", 1);
2086: MLPutFunction(link, "Part", 2);
2087: MLPutSymbol(link, "points");
2088: MLPutFunction(link, "Plus", 2);
2089: MLPutFunction(link, "Part", 2);
2090: MLPutFunction(link, "Slot", 1);
2091: MLPutInteger(link, 1);
2092: MLPutFunction(link, "List", 3);
2093: MLPutInteger(link, 2);
2094: MLPutInteger(link, 4);
2095: MLPutInteger(link, 6);
2096: MLPutInteger(link, 1);
2097: /* Triangle 2--4--3 */
2098: MLPutFunction(link, "Polygon", 1);
2099: MLPutFunction(link, "Part", 2);
2100: MLPutSymbol(link, "points");
2101: MLPutFunction(link, "Plus", 2);
2102: MLPutFunction(link, "Part", 2);
2103: MLPutFunction(link, "Slot", 1);
2104: MLPutInteger(link, 1);
2105: MLPutFunction(link, "List", 3);
2106: MLPutInteger(link, 3);
2107: MLPutInteger(link, 5);
2108: MLPutInteger(link, 4);
2109: MLPutInteger(link, 1);
2110: /* Triangle 3--4--5 */
2111: MLPutFunction(link, "Polygon", 1);
2112: MLPutFunction(link, "Part", 2);
2113: MLPutSymbol(link, "points");
2114: MLPutFunction(link, "Plus", 2);
2115: MLPutFunction(link, "Part", 2);
2116: MLPutFunction(link, "Slot", 1);
2117: MLPutInteger(link, 1);
2118: MLPutFunction(link, "List", 3);
2119: MLPutInteger(link, 4);
2120: MLPutInteger(link, 5);
2121: MLPutInteger(link, 6);
2122: MLPutInteger(link, 1);
2123: }
2124: else if ((numCorners == 3) || (!grid->cm->fieldClasses[grid->viewField][1]))
2125: {
2126: /* Triangle 0--1--2 */
2127: MLPutFunction(link, "Polygon", 1);
2128: MLPutFunction(link, "Part", 2);
2129: MLPutSymbol(link, "points");
2130: MLPutFunction(link, "Plus", 2);
2131: MLPutFunction(link, "Part", 2);
2132: MLPutFunction(link, "Slot", 1);
2133: MLPutInteger(link, 1);
2134: MLPutFunction(link, "List", 3);
2135: MLPutInteger(link, 1);
2136: MLPutInteger(link, 2);
2137: MLPutInteger(link, 3);
2138: MLPutInteger(link, 1);
2139: } else {
2140: SETERRQ(PETSC_ERR_ARG_WRONG, "Invalid number of local nodes");
2141: }
2142: MLPutSymbol(link, "faces");
2143: MLPutFunction(link, "Rule", 2);
2144: MLPutSymbol(link, "AspectRatio");
2145: MLPutReal(link, 1.0);
2146: return(0);
2147: #endif
2149: return(0);
2150: }
2154: int PetscViewerMathematicaVectorPlot_Triangular_2D(PetscViewer viewer, GVec v)
2155: {
2156: #ifdef PETSC_HAVE_MATHEMATICA
2157: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
2158: MLINK link = vmath->link; /* The link to Mathematica */
2159: Grid grid;
2160: Mesh mesh;
2161: PetscReal scale;
2162: PetscTruth opt;
2163: int ierr;
2166: GVecGetGrid(v, &grid);
2167: GridGetMesh(grid, &mesh);
2169: MLPutFunction(link, "ListPlotVectorField", 3);
2170: MLPutSymbol(link, "points");
2171: MLPutFunction(link, "Rule", 2);
2172: MLPutSymbol(link, "AspectRatio");
2173: MLPutReal(link, mesh->sizeY/mesh->sizeX);
2174: MLPutFunction(link, "Rule", 2);
2175: MLPutSymbol(link, "ScaleFactor");
2176: PetscOptionsGetReal("viewer_", "-math_scale", &scale, &opt);
2177: if (opt == PETSC_TRUE) {
2178: MLPutReal(link, scale);
2179: } else {
2180: MLPutSymbol(link, "None");
2181: }
2182: return(0);
2183: #endif
2185: return(0);
2186: }
2190: int PetscViewerMathematicaSurfacePlot_Triangular_2D(PetscViewer viewer, GVec v)
2191: {
2192: #ifdef PETSC_HAVE_MATHEMATICA
2193: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
2194: MLINK link = vmath->link; /* The link to Mathematica */
2195: PetscReal startX, endX;
2196: PetscReal startY, endY;
2197: PetscTruth opt;
2198: int ierr;
2201: /* Setup domain */
2202: startX = 0.0;
2203: startY = 0.0;
2204: endX = 1.0;
2205: endY = 1.0;
2206: PetscOptionsGetReal("viewer_", "-math_start_x", &startX, &opt);
2207: PetscOptionsGetReal("viewer_", "-math_start_y", &startY, &opt);
2208: PetscOptionsGetReal("viewer_", "-math_end_x", &endX, &opt);
2209: PetscOptionsGetReal("viewer_", "-math_end_y", &endY, &opt);
2211: MLPutFunction(link, "Show", 1);
2212: MLPutFunction(link, "SurfaceGraphics", 6);
2213: MLPutSymbol(link, "points");
2214: MLPutFunction(link, "Rule", 2);
2215: MLPutSymbol(link, "ClipFill");
2216: MLPutSymbol(link, "None");
2217: MLPutFunction(link, "Rule", 2);
2218: MLPutSymbol(link, "Axes");
2219: MLPutSymbol(link, "True");
2220: MLPutFunction(link, "Rule", 2);
2221: MLPutSymbol(link, "PlotLabel");
2222: MLPutString(link, vmath->objName);
2223: MLPutFunction(link, "Rule", 2);
2224: MLPutSymbol(link, "MeshRange");
2225: MLPutFunction(link, "List", 2);
2226: MLPutFunction(link, "List", 2);
2227: MLPutReal(link, startX);
2228: MLPutReal(link, endX);
2229: MLPutFunction(link, "List", 2);
2230: MLPutReal(link, startY);
2231: MLPutReal(link, endY);
2232: MLPutFunction(link, "Rule", 2);
2233: MLPutSymbol(link, "AspectRatio");
2234: MLPutReal(link, (endY - startY)/(endX - startX));
2235: return(0);
2236: #endif
2238: return(0);
2239: }
2243: int PetscViewerMathematica_GVec_Triangular_2D(PetscViewer viewer, GVec v)
2244: {
2245: PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
2246: #ifdef PETSC_HAVE_MATHEMATICA
2247: MLINK link = vmath->link; /* The link to Mathematica */
2248: Mesh mesh;
2249: Grid grid;
2250: Mat P;
2251: GVec w;
2252: int numCorners;
2253: int ierr;
2256: GVecGetGrid(v, &grid);
2257: mesh = grid->mesh;
2258: numCorners = mesh->numCorners;
2260: /* Check that a field component has been specified */
2261: if ((grid->viewField < 0) || (grid->viewField >= grid->numFields)) return(0);
2263: if (grid->isConstrained) {
2264: GVecCreate(grid, &w);
2265: GridGetConstraintMatrix(grid, &P);
2266: MatMult(P, v, w);
2267: } else {
2268: w = v;
2269: }
2271: /* Check that the mesh is defined correctly */
2272: PetscViewerMathematicaCheckMesh_Triangular_2D(viewer, mesh);
2274: /* Send the first component of the first field */
2275: MLPutFunction(link, "EvaluatePacket", 1);
2276: switch(vmath->plotType)
2277: {
2278: case MATHEMATICA_TRIANGULATION_PLOT:
2279: MLPutFunction(link, "Module", 2);
2280: /* Make temporary points with each value of the field component in the vector */
2281: MLPutFunction(link, "List", 2);
2282: MLPutSymbol(link, "f");
2283: MLPutFunction(link, "Set", 2);
2284: MLPutSymbol(link, "points");
2285: PetscViewerMathematicaCreateSamplePoints_Triangular_2D(viewer, w);
2286: /* Display the picture */
2287: PetscViewerMathematicaTriangulationPlot_Triangular_2D(viewer, w);
2288: break;
2289: case MATHEMATICA_VECTOR_TRIANGULATION_PLOT:
2290: if (grid->fields[grid->viewField].numComp != 2) {
2291: SETERRQ(PETSC_ERR_ARG_WRONG, "Field must be a 2D vector field for this plot type");
2292: }
2293: MLPutFunction(link, "Module", 2);
2294: /* Make temporary list with points and 2D vector field values */
2295: MLPutFunction(link, "List", 2);
2296: MLPutSymbol(link, "f");
2297: MLPutFunction(link, "Set", 2);
2298: MLPutSymbol(link, "points");
2299: PetscViewerMathematicaCreateVectorSamplePoints_Triangular_2D(viewer, w);
2300: /* Display the picture */
2301: PetscViewerMathematicaVectorPlot_Triangular_2D(viewer, w);
2302: break;
2303: case MATHEMATICA_VECTOR_PLOT:
2304: if (grid->fields[grid->viewField].numComp != 2) {
2305: SETERRQ(PETSC_ERR_ARG_WRONG, "Field must be a 2D vector field for this plot type");
2306: }
2307: MLPutFunction(link, "Module", 2);
2308: /* Make temporary list with points and 2D vector field values */
2309: MLPutFunction(link, "List", 2);
2310: MLPutSymbol(link, "f");
2311: MLPutFunction(link, "Set", 2);
2312: MLPutSymbol(link, "points");
2313: PetscViewerMathematicaCreateInterpolatedSamplePoints_Triangular_2D(viewer, w, -1);
2314: /* Display the picture */
2315: PetscViewerMathematicaVectorPlot_Triangular_2D(viewer, w);
2316: break;
2317: case MATHEMATICA_SURFACE_PLOT:
2318: if (grid->fields[grid->viewField].numComp != 2) {
2319: SETERRQ(PETSC_ERR_ARG_WRONG, "Field must be a 2D vector field for this plot type");
2320: }
2321: MLPutFunction(link, "Module", 2);
2322: /* Make temporary list with interpolated field values on a square mesh */
2323: MLPutFunction(link, "List", 1);
2324: MLPutFunction(link, "Set", 2);
2325: MLPutSymbol(link, "points");
2326: PetscViewerMathematicaCreateInterpolatedSamplePoints_Triangular_2D(viewer, w, grid->viewComp);
2327: /* Display the picture */
2328: PetscViewerMathematicaSurfacePlot_Triangular_2D(viewer, w);
2329: break;
2330: default:
2331: SETERRQ(PETSC_ERR_ARG_WRONG, "Invalid plot type");
2332: }
2333: MLEndPacket(link);
2334: /* Skip packets until ReturnPacket */
2335: PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
2336: /* Skip ReturnPacket */
2337: MLNewPacket(link);
2339: if (grid->isConstrained) {
2340: VecDestroy(w);
2341: }
2342: #else
2344: switch(vmath->plotType)
2345: {
2346: case MATHEMATICA_TRIANGULATION_PLOT:
2347: break;
2348: case MATHEMATICA_VECTOR_TRIANGULATION_PLOT:
2349: break;
2350: case MATHEMATICA_VECTOR_PLOT:
2351: break;
2352: case MATHEMATICA_SURFACE_PLOT:
2353: break;
2354: default:
2355: SETERRQ(PETSC_ERR_ARG_WRONG, "Invalid plot type");
2356: }
2357: #endif
2358: return(0);
2359: }
2360: #endif