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:  #include src/sles/pc/impls/ml/ml.h
 15:  #include src/mat/impls/aij/seq/aij.h
 16:  #include mathematica.h
 17: #include "petscfix.h"

 19: #if defined (PETSC_HAVE__SNPRINTF)
 20: #define snprintf _snprintf
 21: #endif

 23: PetscViewer  VIEWER_MATHEMATICA_WORLD_PRIVATE = PETSC_NULL;
 24: #ifdef PETSC_HAVE_MATHEMATICA
 25: static void *mathematicaEnv                   = PETSC_NULL;
 26: #endif

 28: /*@C
 29:   PetscViewerMathematicaInitializePackage - This function initializes everything in the Petsc interface to Mathematica. It is
 30:   called from PetscDLLibraryRegister() when using dynamic libraries, and on the call to PetscInitialize()
 31:   when using static libraries.

 33:   Input Parameter:
 34:   path - The dynamic library path, or PETSC_NULL

 36:   Level: developer

 38: .keywords: Petsc, initialize, package, PLAPACK
 39: .seealso: PetscInitializePackage(), PetscInitialize()
 40: @*/
 41: int PetscViewerMathematicaInitializePackage(char *path) {
 42:   static PetscTruth initialized = PETSC_FALSE;

 45:   if (initialized == PETSC_TRUE) return(0);
 46:   initialized = PETSC_TRUE;
 47: #ifdef PETSC_HAVE_MATHEMATICA
 48:   mathematicaEnv = (void *) MLInitialize(0);
 49: #endif
 50:   return(0);
 51: }

 53: /*@C
 54:   PetscViewerMathematicaDestroyPackage - This function destroys everything in the Petsc interface to Mathematica. It is
 55:   called from PetscFinalize().

 57:   Level: developer

 59: .keywords: Petsc, destroy, package, mathematica
 60: .seealso: PetscFinalize()
 61: @*/
 62: int PetscViewerMathematicaFinalizePackage(void) {
 64: #ifdef PETSC_HAVE_MATHEMATICA
 65:   if (mathematicaEnv != PETSC_NULL) MLDeinitialize((MLEnvironment) mathematicaEnv);
 66: #endif
 67:   return(0);
 68: }

 70: int PetscViewerInitializeMathematicaWorld_Private()
 71: {

 75:   if (VIEWER_MATHEMATICA_WORLD_PRIVATE) return(0);
 76:   PetscViewerMathematicaOpen(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_NULL, PETSC_NULL, &VIEWER_MATHEMATICA_WORLD_PRIVATE);
 77: 
 78:   return(0);
 79: }

 81: static int PetscViewerDestroy_Mathematica(PetscViewer viewer)
 82: {
 83:   PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
 84:   int                 ierr;

 87: #ifdef PETSC_HAVE_MATHEMATICA
 88:   MLClose(vmath->link);
 89: #endif
 90:   if (vmath->linkname != PETSC_NULL) {
 91:     PetscFree(vmath->linkname);
 92:   }
 93:   if (vmath->linkhost != PETSC_NULL) {
 94:     PetscFree(vmath->linkhost);
 95:   }
 96:   PetscFree(vmath);
 97:   return(0);
 98: }

100: int PetscViewerDestroyMathematica_Private(void)
101: {

105:   if (VIEWER_MATHEMATICA_WORLD_PRIVATE) {
106:     PetscViewerDestroy(VIEWER_MATHEMATICA_WORLD_PRIVATE);
107:   }
108:   return(0);
109: }

111: int PetscViewerMathematicaSetupConnection_Private(PetscViewer v) {
112: #ifdef PETSC_HAVE_MATHEMATICA
113:   PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) v->data;
114: #ifdef MATHEMATICA_3_0
115:   int                      argc = 6;
116:   char                    *argv[6];
117: #else
118:   int                      argc = 5;
119:   char                    *argv[5];
120: #endif
121:   char                     hostname[256];
122:   long                     lerr;
123:   int                      ierr;
124: #endif

127: #ifdef PETSC_HAVE_MATHEMATICA
128:   /* Link name */
129:   argv[0] = "-linkname";
130:   if (vmath->linkname == PETSC_NULL) {
131:     argv[1] = "math -mathlink";
132:   } else {
133:     argv[1] = vmath->linkname;
134:   }

136:   /* Link host */
137:   argv[2] = "-linkhost";
138:   if (vmath->linkhost == PETSC_NULL) {
139:     PetscGetHostName(hostname, 255);
140:     argv[3] = hostname;
141:   } else {
142:     argv[3] = vmath->linkhost;
143:   }

145:   /* Link mode */
146: #ifdef MATHEMATICA_3_0
147:   argv[4] = "-linkmode";
148:   switch(vmath->linkmode) {
149:   case MATHEMATICA_LINK_CREATE:
150:     argv[5] = "Create";
151:     break;
152:   case MATHEMATICA_LINK_CONNECT:
153:     argv[5] = "Connect";
154:     break;
155:   case MATHEMATICA_LINK_LAUNCH:
156:     argv[5] = "Launch";
157:     break;
158:   }
159: #else
160:   switch(vmath->linkmode) {
161:   case MATHEMATICA_LINK_CREATE:
162:     argv[4] = "-linkcreate";
163:     break;
164:   case MATHEMATICA_LINK_CONNECT:
165:     argv[4] = "-linkconnect";
166:     break;
167:   case MATHEMATICA_LINK_LAUNCH:
168:     argv[4] = "-linklaunch";
169:     break;
170:   }
171: #endif

173:   vmath->link = MLOpenInEnv(mathematicaEnv, argc, argv, &lerr);
174: #endif

176:   return(0);
177: }

179: EXTERN_C_BEGIN
180: int PetscViewerCreate_Mathematica(PetscViewer v) {
181:   PetscViewer_Mathematica *vmath;
182:   int                      ierr;


186:   PetscNew(PetscViewer_Mathematica, &vmath);
187:   v->data         = (void *) vmath;
188:   v->ops->destroy = PetscViewerDestroy_Mathematica;
189:   v->ops->flush   = 0;
190:   PetscStrallocpy(PETSC_VIEWER_MATHEMATICA, &v->type_name);

192:   vmath->linkname         = PETSC_NULL;
193:   vmath->linkhost         = PETSC_NULL;
194:   vmath->linkmode         = MATHEMATICA_LINK_CONNECT;
195:   vmath->graphicsType     = GRAPHICS_MOTIF;
196:   vmath->plotType         = MATHEMATICA_TRIANGULATION_PLOT;
197:   vmath->objName          = PETSC_NULL;

199:   PetscViewerMathematicaSetFromOptions(v);
200:   PetscViewerMathematicaSetupConnection_Private(v);
201:   return(0);
202: }
203: EXTERN_C_END

205: int PetscViewerMathematicaParseLinkMode_Private(char *modename, LinkMode *mode) {
206:   PetscTruth isCreate, isConnect, isLaunch;
207:   int        ierr;

210:   PetscStrcasecmp(modename, "Create",  &isCreate);
211:   PetscStrcasecmp(modename, "Connect", &isConnect);
212:   PetscStrcasecmp(modename, "Launch",  &isLaunch);
213:   if (isCreate == PETSC_TRUE) {
214:     *mode = MATHEMATICA_LINK_CREATE;
215:   } else if (isConnect == PETSC_TRUE) {
216:     *mode = MATHEMATICA_LINK_CONNECT;
217:   } else if (isLaunch == PETSC_TRUE) {
218:     *mode = MATHEMATICA_LINK_LAUNCH;
219:   } else {
220:     SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid Mathematica link mode: %s", modename);
221:   }
222:   return(0);
223: }

225: int PetscViewerMathematicaSetFromOptions(PetscViewer v) {
226:   PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) v->data;
227:   char                     linkname[256];
228:   char                     modename[256];
229:   char                     hostname[256];
230:   char                     type[256];
231:   int                      numPorts;
232:   int                     *ports;
233:   int                      numHosts;
234:   char                   **hosts;
235:   int                      size, rank;
236:   int                      h;
237:   PetscTruth               opt;
238:   int                      ierr;

241:   MPI_Comm_size(v->comm, &size);
242:   MPI_Comm_rank(v->comm, &rank);

244:   /* Get link name */
245:   PetscOptionsGetString("viewer_", "-math_linkname", linkname, 255, &opt);
246:   if (opt == PETSC_TRUE) {
247:     PetscViewerMathematicaSetLinkName(v, linkname);
248:   }
249:   /* Get link port */
250:   numPorts = size;
251:   PetscMalloc(size * sizeof(int), &ports);
252:   PetscOptionsGetIntArray("viewer_", "-math_linkport", ports, &numPorts, &opt);
253:   if (opt == PETSC_TRUE) {
254:     if (numPorts > rank) {
255:       snprintf(linkname, 255, "%6d", ports[rank]);
256:     } else {
257:       snprintf(linkname, 255, "%6d", ports[0]);
258:     }
259:     PetscViewerMathematicaSetLinkName(v, linkname);
260:   }
261:   PetscFree(ports);
262:   /* Get link host */
263:   numHosts = size;
264:   PetscMalloc(size * sizeof(char *), &hosts);
265:   PetscOptionsGetStringArray("viewer_", "-math_linkhost", hosts, &numHosts, &opt);
266:   if (opt == PETSC_TRUE) {
267:     if (numHosts > rank) {
268:       PetscStrncpy(hostname, hosts[rank], 255);
269:     } else {
270:       PetscStrncpy(hostname, hosts[0], 255);
271:     }
272:     PetscViewerMathematicaSetLinkHost(v, hostname);
273:   }
274:   for(h = 0; h < numHosts; h++) {
275:     PetscFree(hosts[h]);
276:   }
277:   PetscFree(hosts);
278:   /* Get link mode */
279:   PetscOptionsGetString("viewer_", "-math_linkmode", modename, 255, &opt);
280:   if (opt == PETSC_TRUE) {
281:     LinkMode mode;

283:     PetscViewerMathematicaParseLinkMode_Private(modename, &mode);
284:     PetscViewerMathematicaSetLinkMode(v, mode);
285:   }
286:   /* Get graphics type */
287:   PetscOptionsGetString("viewer_", "-math_graphics", type, 255, &opt);
288:   if (opt == PETSC_TRUE) {
289:     PetscTruth isMotif, isPS, isPSFile;

291:     PetscStrcasecmp(type, "Motif",  &isMotif);
292:     PetscStrcasecmp(type, "PS",     &isPS);
293:     PetscStrcasecmp(type, "PSFile", &isPSFile);
294:     if (isMotif == PETSC_TRUE) {
295:       vmath->graphicsType = GRAPHICS_MOTIF;
296:     } else if (isPS == PETSC_TRUE) {
297:       vmath->graphicsType = GRAPHICS_PS_STDOUT;
298:     } else if (isPSFile == PETSC_TRUE) {
299:       vmath->graphicsType = GRAPHICS_PS_FILE;
300:     }
301:   }
302:   /* Get plot type */
303:   PetscOptionsGetString("viewer_", "-math_type", type, 255, &opt);
304:   if (opt == PETSC_TRUE) {
305:     PetscTruth isTri, isVecTri, isVec, isSurface;

307:     PetscStrcasecmp(type, "Triangulation",       &isTri);
308:     PetscStrcasecmp(type, "VectorTriangulation", &isVecTri);
309:     PetscStrcasecmp(type, "Vector",              &isVec);
310:     PetscStrcasecmp(type, "Surface",             &isSurface);
311:     if (isTri == PETSC_TRUE) {
312:       vmath->plotType     = MATHEMATICA_TRIANGULATION_PLOT;
313:     } else if (isVecTri == PETSC_TRUE) {
314:       vmath->plotType     = MATHEMATICA_VECTOR_TRIANGULATION_PLOT;
315:     } else if (isVec == PETSC_TRUE) {
316:       vmath->plotType     = MATHEMATICA_VECTOR_PLOT;
317:     } else if (isSurface == PETSC_TRUE) {
318:       vmath->plotType     = MATHEMATICA_SURFACE_PLOT;
319:     }
320:   }
321:   return(0);
322: }

324: int PetscViewerMathematicaSetLinkName(PetscViewer v, const char *name) {
325:   PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) v->data;
326:   int                      ierr;

330:   PetscStrallocpy(name, &vmath->linkname);
331:   return(0);
332: }

334: int PetscViewerMathematicaSetLinkPort(PetscViewer v, int port) {
335:   char name[16];
336:   int  ierr;

339:   snprintf(name, 16, "%6d", port);
340:   PetscViewerMathematicaSetLinkName(v, name);
341:   return(0);
342: }

344: int PetscViewerMathematicaSetLinkHost(PetscViewer v, const char *host) {
345:   PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) v->data;
346:   int                      ierr;

350:   PetscStrallocpy(host, &vmath->linkhost);
351:   return(0);
352: }

354: int PetscViewerMathematicaSetLinkMode(PetscViewer v, LinkMode mode) {
355:   PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) v->data;

358:   vmath->linkmode = mode;
359:   return(0);
360: }

362: /*----------------------------------------- Public Functions --------------------------------------------------------*/
363: /*@C
364:   PetscViewerMathematicaOpen - Communicates with Mathemtica using MathLink.

366:   Collective on comm

368:   Input Parameters:
369: + comm    - The MPI communicator
370: . port    - [optional] The port to connect on, or PETSC_DECIDE
371: . machine - [optional] The machine to run Mathematica on, or PETSC_NULL
372: - mode    - [optional] The connection mode, or PETSC_NULL

374:   Output Parameter:
375: . viewer  - The Mathematica viewer

377:   Level: intermediate

379:   Notes:
380:   Most users should employ the following commands to access the 
381:   Mathematica viewers
382: $
383: $    PetscViewerMathematicaOpen(MPI_Comm comm, int port, char *machine, char *mode, PetscViewer &viewer)
384: $    MatView(Mat matrix, PetscViewer viewer)
385: $
386: $                or
387: $
388: $    PetscViewerMathematicaOpen(MPI_Comm comm, int port, char *machine, char *mode, PetscViewer &viewer)
389: $    VecView(Vec vector, PetscViewer viewer)

391:    Options Database Keys:
392: $    -viewer_math_linkhost <machine> - The host machine for the kernel
393: $    -viewer_math_linkname <name>    - The full link name for the connection
394: $    -viewer_math_linkport <port>    - The port for the connection
395: $    -viewer_math_mode <mode>        - The mode, e.g. Launch, Connect
396: $    -viewer_math_type <type>        - The plot type, e.g. Triangulation, Vector
397: $    -viewer_math_graphics <output>  - The output type, e.g. Motif, PS, PSFile

399: .keywords: PetscViewer, Mathematica, open

401: .seealso: MatView(), VecView()
402: @*/
403: int PetscViewerMathematicaOpen(MPI_Comm comm, int port, const char machine[], const char mode[], PetscViewer *v) {
404:   int      ierr;

407:   PetscViewerCreate(comm, v);
408: #if 0
409:   LinkMode linkmode;
410:   PetscViewerMathematicaSetLinkPort(*v, port);
411:   PetscViewerMathematicaSetLinkHost(*v, machine);
412:   PetscViewerMathematicaParseLinkMode_Private(mode, &linkmode);
413:   PetscViewerMathematicaSetLinkMode(*v, linkmode);
414: #endif
415:   PetscViewerSetType(*v, PETSC_VIEWER_MATHEMATICA);
416:   return(0);
417: }

419: #ifdef PETSC_HAVE_MATHEMATICA
420: /*@C
421:   PetscViewerMathematicaGetLink - Returns the link to Mathematica

423:   Input Parameters:
424: . viewer - The Mathematica viewer
425: . link   - The link to Mathematica

427:   Level: intermediate

429: .keywords PetscViewer, Mathematica, link
430: .seealso PetscViewerMathematicaOpen()
431: @*/
432: int PetscViewerMathematicaGetLink(PetscViewer viewer, MLINK *link)
433: {
434:   PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;

438:   *link = vmath->link;
439:   return(0);
440: }
441: #endif

443: /*@C
444:   PetscViewerMathematicaSkipPackets - Discard packets sent by Mathematica until a certain packet type is received

446:   Input Parameters:
447: . viewer - The Mathematica viewer
448: . type   - The packet type to search for, e.g RETURNPKT

450:   Level: advanced

452: .keywords PetscViewer, Mathematica, packets
453: .seealso PetscViewerMathematicaSetName(), PetscViewerMathematicaGetVector()
454: @*/
455: int PetscViewerMathematicaSkipPackets(PetscViewer viewer, int type)
456: {
457: #ifdef PETSC_HAVE_MATHEMATICA
458:   PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
459:   MLINK               link  = vmath->link; /* The link to Mathematica */
460:   int                 pkt;                 /* The packet type */

463:   while((pkt = MLNextPacket(link)) && (pkt != type))
464:     MLNewPacket(link);
465:   if (!pkt) {
466:     MLClearError(link);
467:     SETERRQ(PETSC_ERR_LIB, (char *) MLErrorMessage(link));
468:   }
469:   return(0);
470: #endif
472:   return(0);
473: }

475: /*@C
476:   PetscViewerMathematicaLoadGraphics - Change the type of graphics output for Mathematica

478:   Input Parameters:
479: . viewer - The Mathematica viewer
480: . type   - The type of graphics, e.g. GRAPHICS_MOTIF, GRAPHICS_PS_FILE, GRAPHICS_PS_STDOUT

482:   Level: intermediate

484: .keywords PetscViewer, Mathematica, packets
485: .seealso PetscViewerMathematicaSetName(), PetscViewerMathematicaSkipPackets()
486: @*/
487: int PetscViewerMathematicaLoadGraphics(PetscViewer viewer, GraphicsType type)
488: {
489: #ifdef PETSC_HAVE_MATHEMATICA
490:   PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
491:   MLINK               link  = vmath->link; /* The link to Mathematica */
492:   char                programName[256];
493:   int                 ierr;

496:   /* Load graphics package */
497:   MLPutFunction(link, "EvaluatePacket", 1);
498:     switch(type)
499:     {
500:     case GRAPHICS_MOTIF:
501:         MLPutFunction(link, "Get", 1);
502:           MLPutString(link, "Motif.m");
503:       break;
504:     case GRAPHICS_PS_FILE:
505:     MLPutFunction(link, "CompoundExpression", 4);
506:       MLPutFunction(link, "Set", 2);
507:         MLPutSymbol(link, "PetscGraphicsCounter");
508:         MLPutInteger(link, 0);
509:       MLPutFunction(link, "SetDelayed", 2);
510:         MLPutSymbol(link, "$Display");
511:         MLPutSymbol(link, "$FileDisplay");
512:       MLPutFunction(link, "Set", 2);
513:         MLPutSymbol(link, "$FileDisplay");
514:         MLPutFunction(link, "OpenWrite", 1);
515:           MLPutFunction(link, "StringJoin", 3);
516:           if (!PetscGetProgramName(programName, 255))
517:             MLPutString(link, programName);
518:           else
519:             MLPutString(link, "GVec");
520:             MLPutFunction(link, "ToString", 1);
521:               MLPutSymbol(link, "PetscGraphicsCounter");
522:             MLPutString(link, ".ps");
523:       MLPutFunction(link, "Set", 2);
524:         MLPutSymbol(link, "$DisplayFunction");
525:         MLPutFunction(link, "Function", 1);
526:           MLPutFunction(link, "CompoundExpression", 2);
527:             MLPutFunction(link, "Display", 3);
528:               MLPutSymbol(link, "$Display");
529:               MLPutFunction(link, "Slot", 1);
530:                 MLPutInteger(link, 1);
531:               MLPutString(link, "EPS");
532:             MLPutFunction(link, "Increment", 1);
533:               MLPutSymbol(link, "PetscGraphicsCounter");
534:       break;
535:     case GRAPHICS_PS_STDOUT:
536:       MLPutFunction(link, "Get", 1);
537:         MLPutString(link, "PSDirect.m");
538:       break;
539:     default:
540:       SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid graphics type: %d", type);
541:     }
542:   MLEndPacket(link);

544:   /* Skip packets until ReturnPacket */
545:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
546:   /* Skip ReturnPacket */
547:   MLNewPacket(link);

549:   /* Load PlotField.m for vector plots */
550:   MLPutFunction(link, "EvaluatePacket", 1);
551:     MLPutFunction(link, "Get", 1);
552:       MLPutString(link, "Graphics/PlotField.m");
553:   MLEndPacket(link);

555:   /* Skip packets until ReturnPacket */
556:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
557:   /* Skip ReturnPacket */
558:   MLNewPacket(link);

560:   return(0);
561: #else
563:   switch(type)
564:   {
565:   case GRAPHICS_MOTIF:
566:   case GRAPHICS_PS_FILE:
567:   case GRAPHICS_PS_STDOUT:
568:     break;
569:   default:
570:     SETERRQ1(PETSC_ERR_ARG_WRONG, "Invalid graphics type: %d", type);
571:   }
572:   return(0);
573: #endif
574: }

576: /*@C
577:   PetscViewerMathematicaGetName - Retrieve the default name for objects communicated to Mathematica

579:   Input Parameter:
580: . viewer - The Mathematica viewer

582:   Output Parameter:
583: . name   - The name for new objects created in Mathematica

585:   Level: intermediate

587: .keywords PetscViewer, Mathematica, name
588: .seealso PetscViewerMathematicaSetName(), PetscViewerMathematicaClearName()
589: @*/
590: int PetscViewerMathematicaGetName(PetscViewer viewer, const char **name)
591: {
592:   PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;

597:   *name = vmath->objName;
598:   return(0);
599: }

601: /*@C
602:   PetscViewerMathematicaSetName - Override the default name for objects communicated to Mathematica

604:   Input Parameters:
605: . viewer - The Mathematica viewer
606: . name   - The name for new objects created in Mathematica

608:   Level: intermediate

610: .keywords PetscViewer, Mathematica, name
611: .seealso PetscViewerMathematicaSetName(), PetscViewerMathematicaClearName()
612: @*/
613: int PetscViewerMathematicaSetName(PetscViewer viewer, const char name[])
614: {
615:   PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;

620:   vmath->objName = name;
621:   return(0);
622: }

624: /*@C
625:   PetscViewerMathematicaClearName - Use the default name for objects communicated to Mathematica

627:   Input Parameter:
628: . viewer - The Mathematica viewer

630:   Level: intermediate

632: .keywords PetscViewer, Mathematica, name
633: .seealso PetscViewerMathematicaGetName(), PetscViewerMathematicaSetName()
634: @*/
635: int PetscViewerMathematicaClearName(PetscViewer viewer)
636: {
637:   PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;

641:   vmath->objName = PETSC_NULL;
642:   return(0);
643: }

645: /*@C
646:   PetscViewerMathematicaGetVector - Retrieve a vector from Mathematica

648:   Input Parameter:
649: . viewer - The Mathematica viewer

651:   Output Parameter:
652: . v      - The vector

654:   Level: intermediate

656: .keywords PetscViewer, Mathematica, vector
657: .seealso VecView(), PetscViewerMathematicaPutVector()
658: @*/
659: int PetscViewerMathematicaGetVector(PetscViewer viewer, Vec v) {
660: #ifdef PETSC_HAVE_MATHEMATICA
661:   PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
662:   MLINK               link;   /* The link to Mathematica */
663:   char               *name;
664:   PetscScalar        *mArray;
665:   PetscScalar        *array;
666:   long                mSize;
667:   int                 size;
668:   int                 ierr;


674:   /* Determine the object name */
675:   if (vmath->objName == PETSC_NULL) {
676:     name = "vec";
677:   } else {
678:     name = (char *) vmath->objName;
679:   }

681:   link = vmath->link;
682:   VecGetLocalSize(v, &size);
683:   VecGetArray(v, &array);
684:   MLPutFunction(link, "EvaluatePacket", 1);
685:     MLPutSymbol(link, name);
686:   MLEndPacket(link);
687:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
688:   MLGetRealList(link, &mArray, &mSize);
689:   if (size != mSize) SETERRQ(PETSC_ERR_ARG_WRONG, "Incompatible vector sizes");
690:   PetscMemcpy(array, mArray, mSize * sizeof(double));
691:   MLDisownRealList(link, mArray, mSize);
692:   VecRestoreArray(v, &array);

694:   return(0);
695: #endif
697:   return(0);
698: }

700: /*@C
701:   PetscViewerMathematicaPutVector - Send a vector to Mathematica

703:   Input Parameters:
704: + viewer - The Mathematica viewer
705: - v      - The vector

707:   Level: intermediate

709: .keywords PetscViewer, Mathematica, vector
710: .seealso VecView(), PetscViewerMathematicaGetVector()
711: @*/
712: int PetscViewerMathematicaPutVector(PetscViewer viewer, Vec v) {
713: #ifdef PETSC_HAVE_MATHEMATICA
714:   PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
715:   MLINK               link  = vmath->link; /* The link to Mathematica */
716:   char               *name;
717:   PetscScalar        *array;
718:   int                 size;
719:   int                 ierr;

722:   /* Determine the object name */
723:   if (vmath->objName == PETSC_NULL) {
724:     name = "vec";
725:   } else {
726:     name = (char *) vmath->objName;
727:   }
728:   VecGetLocalSize(v, &size);
729:   VecGetArray(v, &array);

731:   /* Send the Vector object */
732:   MLPutFunction(link, "EvaluatePacket", 1);
733:     MLPutFunction(link, "Set", 2);
734:       MLPutSymbol(link, name);
735:       MLPutRealList(link, array, size);
736:   MLEndPacket(link);
737:   /* Skip packets until ReturnPacket */
738:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
739:   /* Skip ReturnPacket */
740:   MLNewPacket(link);

742:   return(0);
743: #endif
745:   return(0);
746: }

748: int PetscViewerMathematicaPutMatrix(PetscViewer viewer, int m, int n, PetscReal *a) {
749: #ifdef PETSC_HAVE_MATHEMATICA
750:   PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
751:   MLINK               link  = vmath->link; /* The link to Mathematica */
752:   char               *name;
753:   int                 ierr;

756:   /* Determine the object name */
757:   if (vmath->objName == PETSC_NULL) {
758:     name = "mat";
759:   } else {
760:     name = (char *) vmath->objName;
761:   }

763:   /* Send the dense matrix object */
764:   MLPutFunction(link, "EvaluatePacket", 1);
765:     MLPutFunction(link, "Set", 2);
766:       MLPutSymbol(link, name);
767:       MLPutFunction(link, "Transpose", 1);
768:         MLPutFunction(link, "Partition", 2);
769:           MLPutRealList(link, a, m*n);
770:           MLPutInteger(link, m);
771:   MLEndPacket(link);
772:   /* Skip packets until ReturnPacket */
773:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
774:   /* Skip ReturnPacket */
775:   MLNewPacket(link);

777:   return(0);
778: #endif
780:   return(0);
781: }

783: int PetscViewerMathematicaPutCSRMatrix(PetscViewer viewer, int m, int n, int *i, int *j, PetscReal *a) {
784: #ifdef PETSC_HAVE_MATHEMATICA
785:   PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
786:   MLINK               link  = vmath->link; /* The link to Mathematica */
787:   const char         *symbol;
788:   char               *name;
789:   PetscTruth          match;
790:   int                 ierr;

793:   /* Determine the object name */
794:   if (vmath->objName == PETSC_NULL) {
795:     name = "mat";
796:   } else {
797:     name = (char *) vmath->objName;
798:   }

800:   /* Make sure Mathematica recognizes sparse matrices */
801:   MLPutFunction(link, "EvaluatePacket", 1);
802:     MLPutFunction(link, "Needs", 1);
803:       MLPutString(link, "LinearAlgebra`CSRMatrix`");
804:   MLEndPacket(link);
805:   /* Skip packets until ReturnPacket */
806:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
807:   /* Skip ReturnPacket */
808:   MLNewPacket(link);

810:   /* Send the CSRMatrix object */
811:   MLPutFunction(link, "EvaluatePacket", 1);
812:     MLPutFunction(link, "Set", 2);
813:       MLPutSymbol(link, name);
814:       MLPutFunction(link, "CSRMatrix", 5);
815:         MLPutInteger(link, m);
816:         MLPutInteger(link, n);
817:         MLPutFunction(link, "Plus", 2);
818:           MLPutIntegerList(link, i, m+1);
819:           MLPutInteger(link, 1);
820:         MLPutFunction(link, "Plus", 2);
821:           MLPutIntegerList(link, j, i[m]);
822:           MLPutInteger(link, 1);
823:         MLPutRealList(link, a, i[m]);
824:   MLEndPacket(link);
825:   /* Skip packets until ReturnPacket */
826:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
827:   /* Skip ReturnPacket */
828:   MLNewPacket(link);

830:   /* Check that matrix is valid */
831:   MLPutFunction(link, "EvaluatePacket", 1);
832:     MLPutFunction(link, "ValidQ", 1);
833:       MLPutSymbol(link, name);
834:   MLEndPacket(link);
835:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
836:   MLGetSymbol(link, &symbol);
837:   PetscStrcmp("True", (char *) symbol, &match);
838:   if (match == PETSC_FALSE) {
839:     MLDisownSymbol(link, symbol);
840:     SETERRQ(PETSC_ERR_PLIB, "Invalid CSR matrix in Mathematica");
841:   }
842:   MLDisownSymbol(link, symbol);
843:   /* Skip ReturnPacket */
844:   MLNewPacket(link);

846:   return(0);
847: #endif
849:   return(0);
850: }

852: /*------------------------------------------- ML Functions ----------------------------------------------------------*/
853: int PetscViewerMathematicaPartitionMesh(PetscViewer viewer, int numElements, int numVertices, double *vertices, int **mesh,
854:                                    int *numParts, int *colPartition)
855: {
856: #ifdef PETSC_HAVE_MATHEMATICA
857:   PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
858:   MLINK               link;   /* The link to Mathematica */
859:   const char         *symbol;
860:   int                 numOptions;
861:   int                 partSize;
862:   int                *part;
863:   long                numCols;
864:   int                 col;
865:   PetscTruth          match, opt;
866:   int                 ierr;

870:   link = vmath->link;

872:   /* Make sure that the reduce.m package is loaded */
873:   MLPutFunction(link, "EvaluatePacket", 1);
874:     MLPutFunction(link, "Get", 1);
875:       MLPutFunction(link, "StringJoin", 2);
876:         MLPutFunction(link, "Environment", 1);
877:           MLPutString(link, "PETSC_DIR");
878:         MLPutString(link, "/src/pc/impls/ml/reduce.m");
879:   MLEndPacket(link);
880:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
881:   MLGetSymbol(link, &symbol);
882:   PetscStrcmp("$Failed", (char *) symbol, &match);
883:   if (match == PETSC_TRUE) {
884:     MLDisownSymbol(link, symbol);
885:     SETERRQ(PETSC_ERR_FILE_OPEN, "Unable to load package reduce.m");
886:   }
887:   MLDisownSymbol(link, symbol);

889:   /* Partition the mesh */
890:   numOptions = 0;
891:   partSize   = 0;
892:   PetscOptionsGetInt(PETSC_NULL, "-pc_ml_partition_size", &partSize, &opt);
893:   if ((opt == PETSC_TRUE) && (partSize > 0))
894:     numOptions++;
895:   MLPutFunction(link, "EvaluatePacket", 1);
896:     MLPutFunction(link, "PartitionMesh", 1 + numOptions);
897:       MLPutFunction(link, "MeshData", 5);
898:         MLPutInteger(link, numElements);
899:         MLPutInteger(link, numVertices);
900:         MLPutInteger(link, numVertices);
901:         MLPutFunction(link, "Partition", 2);
902:           MLPutRealList(link, vertices, numVertices*2);
903:           MLPutInteger(link, 2);
904:         MLPutFunction(link, "Partition", 2);
905:           MLPutIntegerList(link, mesh[MESH_ELEM], numElements*3);
906:           MLPutInteger(link, 3);
907:       if (partSize > 0)
908:       {
909:         MLPutFunction(link, "Rule", 2);
910:           MLPutSymbol(link, "PartitionSize");
911:           MLPutInteger(link, partSize);
912:       }
913:   MLEndPacket(link);
914:   /* Skip packets until ReturnPacket */
915:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);

917:   /* Get the vertex partiton */
918:   MLGetIntegerList(link, &part, &numCols);
919:   if (numCols != numVertices) SETERRQ(PETSC_ERR_PLIB, "Invalid partition");
920:   for(col = 0, *numParts = 0; col < numCols; col++) {
921:     colPartition[col] = part[col]-1;
922:     *numParts = PetscMax(*numParts, part[col]);
923:   }
924:   MLDisownIntegerList(link, part, numCols);

926:   return(0);
927: #endif
929:   return(0);
930: }

932: int PetscViewerMathematicaReduce(PetscViewer viewer, PC pc, int thresh)
933: {
934: #ifdef PETSC_HAVE_MATHEMATICA
935:   PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
936:   MLINK               link;   /* The link to Mathematica */
937:   PC_Multilevel      *ml;
938:   int                *range;
939:   long                numRange;
940:   int                *null;
941:   long                numNull;
942:   const char         *symbol;
943:   int                 numOptions;
944:   int                 partSize;
945:   int                 row, col;
946:   PetscTruth          match, opt;
947:   int                 ierr;

952:   link = vmath->link;
953:   ml   = (PC_Multilevel *) pc->data;

955:   /* Make sure that the reduce.m package is loaded */
956:   MLPutFunction(link, "EvaluatePacket", 1);
957:     MLPutFunction(link, "Get", 1);
958:       MLPutFunction(link, "StringJoin", 2);
959:         MLPutFunction(link, "Environment", 1);
960:           MLPutString(link, "PETSC_DIR");
961:         MLPutString(link, "/src/pc/impls/ml/reduce.m");
962:   MLEndPacket(link);
963:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
964:   MLGetSymbol(link, &symbol);
965:   PetscStrcmp("$Failed", (char *) symbol, &match);
966:   if (match == PETSC_TRUE) {
967:     MLDisownSymbol(link, symbol);
968:     SETERRQ(PETSC_ERR_FILE_OPEN, "Unable to load package reduce.m");
969:   }
970:   MLDisownSymbol(link, symbol);

972:   /* mesh = MeshData[] */
973:   MLPutFunction(link, "EvaluatePacket", 1);
974:     MLPutFunction(link, "Set", 2);
975:       MLPutSymbol(link, "mesh");
976:       MLPutFunction(link, "MeshData", 5);
977:         MLPutInteger(link, ml->numElements[0]);
978:         MLPutInteger(link, ml->numVertices[0]);
979:         MLPutInteger(link, ml->numVertices[0]);
980:         MLPutFunction(link, "Partition", 2);
981:           MLPutRealList(link, ml->vertices, ml->numVertices[0]*2);
982:           MLPutInteger(link, 2);
983:         MLPutFunction(link, "Partition", 2);
984:           MLPutIntegerList(link, ml->meshes[0][MESH_ELEM], ml->numElements[0]*3);
985:           MLPutInteger(link, 3);
986:   MLEndPacket(link);
987:   /* Skip packets until ReturnPacket */
988:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
989:   /* Skip ReturnPacket */
990:   MLNewPacket(link);
991:   /* Check that mesh is valid */
992:   MLPutFunction(link, "EvaluatePacket", 1);
993:     MLPutFunction(link, "ValidQ", 1);
994:       MLPutSymbol(link, "mesh");
995:   MLEndPacket(link);
996:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
997:   MLGetSymbol(link, &symbol);
998:   PetscStrcmp("True", (char *) symbol, &match);
999:   if (match == PETSC_FALSE) {
1000:     MLDisownSymbol(link, symbol);
1001:     SETERRQ(PETSC_ERR_PLIB, "Invalid mesh in Mathematica");
1002:   }
1003:   MLDisownSymbol(link, symbol);

1005:   /* mat = gradient matrix */
1006:   MatView(ml->locB, viewer);

1008:   /* mattML = ReduceMatrix[mat,ml->minNodes] */
1009:   numOptions = 0;
1010:   partSize   = 0;
1011:   PetscOptionsGetInt(PETSC_NULL, "-pc_ml_partition_size", &partSize, &opt);
1012:   if ((opt == PETSC_TRUE) && (partSize > 0))
1013:     numOptions++;
1014:   MLPutFunction(link, "EvaluatePacket", 1);
1015:     MLPutFunction(link, "Set", 2);
1016:       MLPutSymbol(link, "mattML");
1017:       MLPutFunction(link, "ReduceMatrix", 3 + numOptions);
1018:         MLPutSymbol(link, "mesh");
1019:         MLPutSymbol(link, "mat");
1020:         MLPutInteger(link, thresh);
1021:         if (partSize > 0) {
1022:           MLPutFunction(link, "Rule", 2);
1023:             MLPutSymbol(link, "PartitionSize");
1024:             MLPutInteger(link, partSize);
1025:         }
1026:   MLEndPacket(link);
1027:   /* Skip packets until ReturnPacket */
1028:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1029:   /* Skip ReturnPacket */
1030:   MLNewPacket(link);
1031:   /* Check that mattML is valid */
1032:   MLPutFunction(link, "EvaluatePacket", 1);
1033:     MLPutFunction(link, "ValidQ", 1);
1034:       MLPutSymbol(link, "mattML");
1035:   MLEndPacket(link);
1036:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1037:   MLGetSymbol(link, &symbol);
1038:   PetscStrcmp("True", (char *) symbol, &match);
1039:   if (match == PETSC_FALSE) {
1040:     MLDisownSymbol(link, symbol);
1041:     SETERRQ(PETSC_ERR_PLIB, "Invalid MLData in Mathematica");
1042:   }
1043:   MLDisownSymbol(link, symbol);

1045:   /* Copy information to the preconditioner */
1046:   MLPutFunction(link, "EvaluatePacket", 1);
1047:     MLPutFunction(link, "Part", 2);
1048:       MLPutSymbol(link, "mattML");
1049:       MLPutInteger(link, 3);
1050:   MLEndPacket(link);
1051:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1052:   MLGetInteger(link, &ml->numLevels);

1054:   /* Create lists of the range and nullspace columns */
1055:   MLPutFunction(link, "EvaluatePacket", 1);
1056:     MLPutFunction(link, "GetRange", 1);
1057:       MLPutSymbol(link, "mattML");
1058:   MLEndPacket(link);
1059:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1060:   MLGetIntegerList(link, &range, &numRange);
1061:   if (numRange > ml->sOrder->numLocVars) SETERRQ(PETSC_ERR_PLIB, "Invalid size for range space");
1062:   ml->rank       = numRange;
1063:   ml->globalRank = ml->rank;
1064:   if (ml->rank > 0) {
1065:     PetscMalloc(numRange * sizeof(int), &ml->range);
1066:     for(row = 0; row < numRange; row++)
1067:       ml->range[row] = range[row]-1;
1068:   }
1069:   MLDisownIntegerList(link, range, numRange);

1071:   MLPutFunction(link, "EvaluatePacket", 1);
1072:     MLPutFunction(link, "GetNullColumns", 1);
1073:       MLPutSymbol(link, "mattML");
1074:   MLEndPacket(link);
1075:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1076:   MLGetIntegerList(link, &null, &numNull);
1077:   if (numRange + numNull != ml->sOrder->numLocVars) SETERRQ(PETSC_ERR_PLIB, "Invalid size for range and null spaces");
1078:   ml->numLocNullCols = numNull;
1079:   if (numNull > 0)
1080:   {
1081:     PetscMalloc(numNull * sizeof(int), &ml->nullCols);
1082:     for(col = 0; col < numNull; col++)
1083:       ml->nullCols[col] = null[col] - 1;
1084:   }
1085:   MLDisownIntegerList(link, null, numNull);

1087:   return(0);
1088: #endif
1090:   return(0);
1091: }

1093: int PetscViewerMathematicaMultiLevelConvert(PetscViewer viewer, PC pc)
1094: {
1095: #ifdef PETSC_HAVE_MATHEMATICA
1096:   PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
1097:   MLINK               link;   /* The link to Mathematica */
1098:   PC_Multilevel      *ml;
1099:   Mat_SeqAIJ         *grad;
1100:   int                *numPartitions;
1101:   int                *numPartitionCols, *cols;
1102:   int                *numPartitionRows, *rows;
1103:   double             *U, *singVals, *V;
1104:   long               *Udims, *Vdims;
1105:   char              **Uheads, **Vheads;
1106:   int                *nnz;
1107:   int                *offsets;
1108:   double             *vals;
1109:   long                numLevels, numParts, numCols, numRows, Udepth, numSingVals, Vdepth, len;
1110:   int                 numBdRows, numBdCols;
1111:   int                 mesh, level, part, col, row;
1112:   int                 ierr;

1117:   link = vmath->link;
1118:   ml   = (PC_Multilevel *) pc->data;

1120:   /* ml->numLevels = ml[[3]] */
1121:   MLPutFunction(link, "EvaluatePacket", 1);
1122:     MLPutFunction(link, "Part", 2);
1123:       MLPutSymbol(link, "mattML");
1124:       MLPutInteger(link, 3);
1125:   MLEndPacket(link);
1126:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1127:   MLGetInteger(link, &ml->numLevels);

1129:   /* ml->numMeshes = Length[ml[[4]]] */
1130:   MLPutFunction(link, "EvaluatePacket", 1);
1131:     MLPutFunction(link, "Length", 1);
1132:       MLPutFunction(link, "Part", 2);
1133:         MLPutSymbol(link, "mattML");
1134:         MLPutInteger(link, 4);
1135:   MLEndPacket(link);
1136:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1137:   MLGetInteger(link, &ml->numMeshes);

1139:   /* ml->numElements[level] = ml[[4,level,1]] */
1140:   PetscMalloc(ml->numMeshes * sizeof(int), &ml->numElements);

1142:   /* ml->numVertices[level] = ml[[4,level,2]] */
1143:   PetscMalloc(ml->numMeshes * sizeof(int), &ml->numVertices);

1145:   /* ml->meshes = ml[[4]] */
1146:   PetscMalloc(ml->numMeshes * sizeof(int **), &ml->meshes);
1147:   for(mesh = 0; mesh < ml->numMeshes; mesh++) {
1148:     PetscMalloc(NUM_MESH_DIV * sizeof(int *), &ml->meshes[mesh]);
1149:     /* Here we should get meshes */
1150:     PetscMalloc(1            * sizeof(int),   &ml->meshes[mesh][MESH_OFFSETS]);
1151:     PetscMalloc(1            * sizeof(int),   &ml->meshes[mesh][MESH_ADJ]);
1152:     PetscMalloc(1            * sizeof(int),   &ml->meshes[mesh][MESH_ELEM]);
1153:   }

1155:   /* ml->numPartitions = Map[Length,Drop[ml[[5]],-1]] */
1156:   MLPutFunction(link, "EvaluatePacket", 1);
1157:     MLPutFunction(link, "Map", 2);
1158:       MLPutSymbol(link, "Length");
1159:       MLPutFunction(link, "Drop", 2);
1160:         MLPutFunction(link, "Part", 2);
1161:           MLPutSymbol(link, "mattML");
1162:           MLPutInteger(link, 5);
1163:         MLPutInteger(link, -1);
1164:   MLEndPacket(link);
1165:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1166:   MLGetIntegerList(link, &numPartitions, &numLevels);
1167:   if (numLevels != ml->numLevels) SETERRQ(PETSC_ERR_PLIB, "Invalid node partition in MLData object");
1168:   if (numLevels > 0) {
1169:     PetscMalloc(numLevels * sizeof(int), &ml->numPartitions);
1170:     PetscMemcpy(ml->numPartitions, numPartitions, numLevels * sizeof(int));
1171:   }
1172:   MLDisownIntegerList(link, numPartitions, numLevels);

1174:   if (ml->numLevels > 0) {
1175:     /* ml->numPartitionCols = Map[Length,ml[[5,level]]] */
1176:     PetscMalloc(ml->numLevels * sizeof(int *), &ml->numPartitionCols);
1177:     for(level = 0; level < ml->numLevels; level++) {
1178:       MLPutFunction(link, "EvaluatePacket", 1);
1179:         MLPutFunction(link, "Map", 2);
1180:           MLPutSymbol(link, "Length");
1181:           MLPutFunction(link, "Part", 3);
1182:             MLPutSymbol(link, "mattML");
1183:             MLPutInteger(link, 5);
1184:             MLPutInteger(link, level+1);
1185:       MLEndPacket(link);
1186:       PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1187:       MLGetIntegerList(link, &numPartitionCols, &numParts);
1188:       if (numParts != ml->numPartitions[level]) SETERRQ(PETSC_ERR_PLIB, "Invalid node partition in MLData object");
1189:       if (numParts > 0) {
1190:         PetscMalloc(numParts * sizeof(int), &ml->numPartitionCols[level]);
1191:         PetscMemcpy(ml->numPartitionCols[level], numPartitionCols, numParts * sizeof(int));
1192:       }
1193:       MLDisownIntegerList(link, numPartitionCols, numParts);
1194:     }

1196:     /* ml->colPartition[level][part] = ml[[5,level,part]] */
1197:     PetscMalloc(ml->numLevels * sizeof(int **), &ml->colPartition);
1198:     for(level = 0; level < ml->numLevels; level++) {
1199:       if (ml->numPartitions[level] == 0) continue;
1200:       PetscMalloc(ml->numPartitions[level] * sizeof(int *), &ml->colPartition[level]);
1201:       for(part = 0; part < ml->numPartitions[level]; part++) {
1202:         MLPutFunction(link, "EvaluatePacket", 1);
1203:           MLPutFunction(link, "Part", 4);
1204:             MLPutSymbol(link, "mattML");
1205:             MLPutInteger(link, 5);
1206:             MLPutInteger(link, level+1);
1207:             MLPutInteger(link, part+1);
1208:         MLEndPacket(link);
1209:         PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1210:         MLGetIntegerList(link, &cols, &numCols);
1211:         if (numCols != ml->numPartitionCols[level][part]) SETERRQ(PETSC_ERR_PLIB, "Invalid node partition in MLData object");
1212:         if (numCols > 0) {
1213:           PetscMalloc(numCols * sizeof(int), &ml->colPartition[level][part]);
1214:           /* Convert to zero-based indices */
1215:           for(col = 0; col < numCols; col++) ml->colPartition[level][part][col] = cols[col] - 1;
1216:         }
1217:         MLDisownIntegerList(link, cols, numCols);
1218:       }
1219:     }

1221:     /* ml->numPartitionRows = Map[Length,FlattenAt[ml[[6,level]],1]] */
1222:     PetscMalloc(ml->numLevels * sizeof(int *), &ml->numPartitionRows);
1223:     for(level = 0; level < ml->numLevels; level++) {
1224:       MLPutFunction(link, "EvaluatePacket", 1);
1225:         MLPutFunction(link, "Map", 2);
1226:           MLPutSymbol(link, "Length");
1227:           MLPutFunction(link, "FlattenAt", 2);
1228:             MLPutFunction(link, "Part", 3);
1229:               MLPutSymbol(link, "mattML");
1230:               MLPutInteger(link, 6);
1231:               MLPutInteger(link, level+1);
1232:             MLPutInteger(link, 1);
1233:       MLEndPacket(link);
1234:       PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1235:       MLGetIntegerList(link, &numPartitionRows, &numParts);
1236:       if (numParts != ml->numPartitions[level] + NUM_PART_ROW_DIV - 1) {
1237:         SETERRQ(PETSC_ERR_PLIB, "Invalid edge partition in MLData object");
1238:       }
1239:       PetscMalloc(numParts * sizeof(int), &ml->numPartitionRows[level]);
1240:       PetscMemcpy(ml->numPartitionRows[level], numPartitionRows, numParts * sizeof(int));
1241:       MLDisownIntegerList(link, numPartitionRows, numParts);
1242:     }

1244:     /* ml->rowPartition[level][PART_ROW_INT][part] = ml[[6,level,1,part]]
1245:        ml->rowPartition[level][PART_ROW_BD]        = ml[[6,level,2]]
1246:        ml->rowPartition[level][PART_ROW_RES]       = ml[[6,level,3]] */
1247:     PetscMalloc(ml->numLevels * sizeof(int ***), &ml->rowPartition);
1248:     for(level = 0; level < ml->numLevels; level++) {
1249:       PetscMalloc(NUM_PART_ROW_DIV * sizeof(int **), &ml->rowPartition[level]);
1250:       /* Interior rows */
1251:       if (ml->numPartitions[level] > 0) {
1252:         PetscMalloc(ml->numPartitions[level] * sizeof(int *), &ml->rowPartition[level][PART_ROW_INT]);
1253:         for(part = 0; part < ml->numPartitions[level]; part++) {
1254:           MLPutFunction(link, "EvaluatePacket", 1);
1255:             MLPutFunction(link, "Part", 5);
1256:               MLPutSymbol(link, "mattML");
1257:               MLPutInteger(link, 6);
1258:               MLPutInteger(link, level+1);
1259:               MLPutInteger(link, 1);
1260:               MLPutInteger(link, part+1);
1261:           MLEndPacket(link);
1262:           PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1263:           MLGetIntegerList(link, &rows, &numRows);
1264:           if (numRows != ml->numPartitionRows[level][part]) {
1265:             SETERRQ(PETSC_ERR_PLIB, "Invalid edge partition in MLData object");
1266:           }
1267:           if (numRows > 0) {
1268:             PetscMalloc(numRows * sizeof(int), &ml->rowPartition[level][PART_ROW_INT][part]);
1269:             /* Convert to zero-based indices */
1270:             for(row = 0; row < numRows; row++) {
1271:               ml->rowPartition[level][PART_ROW_INT][part][row] = rows[row] - 1;
1272:             }
1273:           }
1274:           MLDisownIntegerList(link, rows, numRows);
1275:         }
1276:       }
1277:       /* Boundary rows */
1278:       PetscMalloc(1 * sizeof(int *), &ml->rowPartition[level][PART_ROW_BD]);
1279:       MLPutFunction(link, "EvaluatePacket", 1);
1280:         MLPutFunction(link, "Part", 4);
1281:           MLPutSymbol(link, "mattML");
1282:           MLPutInteger(link, 6);
1283:           MLPutInteger(link, level+1);
1284:           MLPutInteger(link, 2);
1285:       MLEndPacket(link);
1286:       PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1287:       MLGetIntegerList(link, &rows, &numRows);
1288:       if (numRows != ml->numPartitionRows[level][ml->numPartitions[level]]) {
1289:         SETERRQ(PETSC_ERR_PLIB, "Invalid edge partition in MLData object");
1290:       }
1291:       if (numRows > 0) {
1292:         PetscMalloc(numRows * sizeof(int), &ml->rowPartition[level][PART_ROW_BD][0]);
1293:         /* Convert to zero-based indices */
1294:         for(row = 0; row < numRows; row++) {
1295:           ml->rowPartition[level][PART_ROW_BD][0][row] = rows[row] - 1;
1296:         }
1297:       }
1298:       MLDisownIntegerList(link, rows, numRows);
1299:       /* Residual rows*/
1300:       PetscMalloc(1 * sizeof(int *), &ml->rowPartition[level][PART_ROW_RES]);
1301:       MLPutFunction(link, "EvaluatePacket", 1);
1302:         MLPutFunction(link, "Part", 4);
1303:           MLPutSymbol(link, "mattML");
1304:           MLPutInteger(link, 6);
1305:           MLPutInteger(link, level+1);
1306:           MLPutInteger(link, 3);
1307:       MLEndPacket(link);
1308:       PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1309:       MLGetIntegerList(link, &rows, &numRows);
1310:       if (numRows != ml->numPartitionRows[level][ml->numPartitions[level]+1]) {
1311:         SETERRQ(PETSC_ERR_PLIB, "Invalid edge partition in MLData object");
1312:       }
1313:       if (numRows > 0) {
1314:         PetscMalloc(numRows * sizeof(int), &ml->rowPartition[level][PART_ROW_RES][0]);
1315:         /* Convert to zero-based indices */
1316:         for(row = 0; row < numRows; row++) {
1317:           ml->rowPartition[level][PART_ROW_RES][0][row] = rows[row] - 1;
1318:         }
1319:       }
1320:       MLDisownIntegerList(link, rows, numRows);
1321:     }
1322:   } else {
1323:     PetscMalloc(1 * sizeof(int),     &ml->numPartitions);
1324:     PetscMalloc(1 * sizeof(int *),   &ml->numPartitionCols);
1325:     PetscMalloc(1 * sizeof(int **),  &ml->colPartition);
1326:     PetscMalloc(1 * sizeof(int *),   &ml->numPartitionRows);
1327:     PetscMalloc(1 * sizeof(int ***), &ml->rowPartition);
1328:   }

1330:   /* ml->numRows = ml[[1]] */
1331:   MLPutFunction(link, "EvaluatePacket", 1);
1332:     MLPutFunction(link, "Part", 2);
1333:       MLPutSymbol(link, "mattML");
1334:       MLPutInteger(link, 1);
1335:   MLEndPacket(link);
1336:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1337:   MLGetInteger(link, &ml->numRows);

1339:   /* ml->numCols = ml[[2]] */
1340:   MLPutFunction(link, "EvaluatePacket", 1);
1341:     MLPutFunction(link, "Part", 2);
1342:       MLPutSymbol(link, "mattML");
1343:       MLPutInteger(link, 2);
1344:   MLEndPacket(link);
1345:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1346:   MLGetInteger(link, &ml->numCols);

1348:   /* ml->zeroTol = ml[[9]] */
1349:   MLPutFunction(link, "EvaluatePacket", 1);
1350:     MLPutFunction(link, "N", 1);
1351:       MLPutFunction(link, "Part", 2);
1352:         MLPutSymbol(link, "mattML");
1353:         MLPutInteger(link, 9);
1354:   MLEndPacket(link);
1355:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1356:   MLGetReal(link, &ml->zeroTol);

1358:   if (ml->numLevels > 0) {
1359:     /* ml->factors[level][part][FACT_U]    = ml[[7,level,part,1]]
1360:        ml->factors[level][part][FACT_DINV] = Divide[1,Select[ml[[7,level,part,2]],(#>ml[[9]])&]]
1361:        ml->factors[level][part][FACT_V]    = ml[[7,level,part,3]] */
1362:     PetscMalloc(ml->numLevels * sizeof(double ***), &ml->factors);
1363:     for(level = 0; level < ml->numLevels; level++) {
1364:       if (ml->numPartitions[level] == 0) continue;
1365:       PetscMalloc(ml->numPartitions[level] * sizeof(double **), &ml->factors[level]);
1366:       for(part = 0; part < ml->numPartitions[level]; part++) {
1367:         PetscMalloc(NUM_FACT_DIV * sizeof(double *), &ml->factors[level][part]);
1368:         /* U, or U^T in LAPACK terms */
1369:         MLPutFunction(link, "EvaluatePacket", 1);
1370:           MLPutFunction(link, "Part", 5);
1371:             MLPutSymbol(link, "mattML");
1372:             MLPutInteger(link, 7);
1373:             MLPutInteger(link, level+1);
1374:             MLPutInteger(link, part+1);
1375:             MLPutInteger(link, 1);
1376:         MLEndPacket(link);
1377:         PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1378:         MLGetRealArray(link, &U, &Udims, &Uheads, &Udepth);
1379:         if (Udepth > 1) {
1380:           if (Udepth != 2) SETERRQ(PETSC_ERR_PLIB, "Invalid U matrix");
1381:           if ((Udims[0] != ml->numPartitionRows[level][part]) || (Udims[0] != Udims[1])) {
1382:             SETERRQ(PETSC_ERR_PLIB, "Incompatible dimensions for U matrix");
1383:           }
1384:           PetscMalloc(Udims[0]*Udims[0] * sizeof(double), &ml->factors[level][part][FACT_U]);
1385:           /* Notice that LAPACK will think that this is U^T, or U in LAPACK terms */
1386:           PetscMemcpy(ml->factors[level][part][FACT_U], U, Udims[0]*Udims[0] * sizeof(double));
1387:         } else if (ml->numPartitionRows[level][part] != 0) {
1388:           SETERRQ(PETSC_ERR_PLIB, "Missing U matrix");
1389:         }
1390:         MLDisownRealArray(link, U, Udims, Uheads, Udepth);
1391:         /* D^{-1} */
1392:         MLPutFunction(link, "EvaluatePacket", 1);
1393:           MLPutFunction(link, "Divide", 2);
1394:             MLPutReal(link, 1.0);
1395:             MLPutFunction(link, "Select", 2);
1396:               MLPutFunction(link, "Part", 5);
1397:                 MLPutSymbol(link, "mattML");
1398:                 MLPutInteger(link, 7);
1399:                 MLPutInteger(link, level+1);
1400:                 MLPutInteger(link, part+1);
1401:                 MLPutInteger(link, 2);
1402:               MLPutFunction(link, "Function", 2);
1403:                 MLPutSymbol(link, "x");
1404:                 MLPutFunction(link, "Greater", 2);
1405:                   MLPutSymbol(link, "x");
1406:                   MLPutFunction(link, "Part", 2);
1407:                     MLPutSymbol(link, "mattML");
1408:                     MLPutInteger(link, 9);
1409:         MLEndPacket(link);
1410:         PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1411:         MLGetRealList(link, &singVals, &numSingVals);
1412:         if (numSingVals > ml->numPartitionCols[level][part]) {
1413:           SETERRQ(PETSC_ERR_PLIB, "Invalid factor list in MLData object");
1414:         }
1415:         if (ml->numPartitionCols[level][part] > 0) {
1416:           PetscMalloc(ml->numPartitionCols[level][part] * sizeof(double), &ml->factors[level][part][FACT_DINV]);
1417:           PetscMemzero(ml->factors[level][part][FACT_DINV], ml->numPartitionCols[level][part] * sizeof(double));
1418:           PetscMemcpy(ml->factors[level][part][FACT_DINV], singVals, numSingVals * sizeof(double));
1419:         }
1420:         MLDisownRealList(link, singVals, numSingVals);
1421:         /* V^T, but V in LAPACK terms */
1422:         MLPutFunction(link, "EvaluatePacket", 1);
1423:           MLPutFunction(link, "Transpose", 1);
1424:             MLPutFunction(link, "Part", 5);
1425:               MLPutSymbol(link, "mattML");
1426:               MLPutInteger(link, 7);
1427:               MLPutInteger(link, level+1);
1428:               MLPutInteger(link, part+1);
1429:               MLPutInteger(link, 3);
1430:         MLEndPacket(link);
1431:         PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1432:         MLGetRealArray(link, &V, &Vdims, &Vheads, &Vdepth);
1433:         if (Vdepth > 1) {
1434:           if (Vdepth != 2) SETERRQ(PETSC_ERR_PLIB, "Invalid V matrix");
1435:           if ((Vdims[0] != ml->numPartitionCols[level][part]) || (Vdims[0] != Vdims[1])) {
1436:             SETERRQ(PETSC_ERR_PLIB, "Incompatible dimensions for U matrix");
1437:           }
1438:           PetscMalloc(Vdims[0]*Vdims[0] * sizeof(double), &ml->factors[level][part][FACT_V]);
1439:           /* Notice that LAPACK will think that this is V, or V^T in LAPACK terms */
1440:           PetscMemcpy(ml->factors[level][part][FACT_V], V, Vdims[0]*Vdims[0] * sizeof(double));
1441:         } else if (ml->numPartitionCols[level][part] != 0) {
1442:           SETERRQ(PETSC_ERR_PLIB, "Missing V matrix");
1443:         }
1444:         MLDisownRealArray(link, V, Vdims, Vheads, Vdepth);
1445:       }
1446:     }

1448:     /* ml->grads = ml[[8]] */
1449:     PetscMalloc(ml->numLevels * sizeof(Mat), &ml->grads);
1450:     PetscMalloc(ml->numLevels * sizeof(Vec), &ml->bdReduceVecs);
1451:     PetscMalloc(ml->numLevels * sizeof(Vec), &ml->colReduceVecs);
1452:     PetscMalloc(ml->numLevels * sizeof(Vec), &ml->colReduceVecs2);
1453:     for(level = 0; level < ml->numLevels; level++) {
1454:       if (ml->numPartitions[level] == 0) continue;
1455:       /* numRows = ml[[8,level,1]] */
1456:       MLPutFunction(link, "EvaluatePacket", 1);
1457:         MLPutFunction(link, "Part", 4);
1458:           MLPutSymbol(link, "mattML");
1459:           MLPutInteger(link, 8);
1460:           MLPutInteger(link, level+1);
1461:           MLPutInteger(link, 1);
1462:       MLEndPacket(link);
1463:       PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1464:       MLGetInteger(link, &numBdRows);
1465:       VecCreateSeq(PETSC_COMM_SELF, numBdRows, &ml->bdReduceVecs[level]);
1466:       /* numCols = ml[[8,level,2]] */
1467:       MLPutFunction(link, "EvaluatePacket", 1);
1468:         MLPutFunction(link, "Part", 4);
1469:           MLPutSymbol(link, "mattML");
1470:           MLPutInteger(link, 8);
1471:           MLPutInteger(link, level+1);
1472:           MLPutInteger(link, 2);
1473:       MLEndPacket(link);
1474:       PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1475:       MLGetInteger(link, &numBdCols);
1476:       VecCreateSeq(PETSC_COMM_SELF, numBdCols, &ml->colReduceVecs[level]);
1477:       VecDuplicate(ml->colReduceVecs[level], &ml->colReduceVecs2[level]);
1478:       /* nnz = Map[Apply[Subtract,Sort[#,Greater]]&, Partition[ml[[8,level,3]],2,1]] */
1479:       MLPutFunction(link, "EvaluatePacket", 1);
1480:         MLPutFunction(link, "Map", 2);
1481:           MLPutFunction(link, "Function", 2);
1482:             MLPutSymbol(link, "x");
1483:             MLPutFunction(link, "Apply", 2);
1484:               MLPutSymbol(link, "Subtract");
1485:               MLPutFunction(link, "Sort", 2);
1486:                 MLPutSymbol(link, "x");
1487:                 MLPutSymbol(link, "Greater");
1488:           MLPutFunction(link, "Partition", 3);
1489:             MLPutFunction(link, "Part", 4);
1490:               MLPutSymbol(link, "mattML");
1491:               MLPutInteger(link, 8);
1492:               MLPutInteger(link, level+1);
1493:               MLPutInteger(link, 3);
1494:             MLPutInteger(link, 2);
1495:             MLPutInteger(link, 1);
1496:       MLEndPacket(link);
1497:       PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1498:       MLGetIntegerList(link, &nnz, &len);
1499:       if (len != numBdRows) SETERRQ(PETSC_ERR_PLIB, "Invalid boundary gradient matrix");
1500:       MatCreateSeqAIJ(PETSC_COMM_SELF, numBdRows, numBdCols, PETSC_DEFAULT, nnz, &ml->grads[level]);
1501:       grad = (Mat_SeqAIJ *) ml->grads[level]->data;
1502:       PetscMemcpy(grad->ilen, nnz, numBdRows * sizeof(int));
1503:       MLDisownIntegerList(link, nnz, len);
1504:       /* ml->grads[level]->i = ml[[8,level,3]] */
1505:       MLPutFunction(link, "EvaluatePacket", 1);
1506:         MLPutFunction(link, "Part", 4);
1507:           MLPutSymbol(link, "mattML");
1508:           MLPutInteger(link, 8);
1509:           MLPutInteger(link, level+1);
1510:           MLPutInteger(link, 3);
1511:       MLEndPacket(link);
1512:       PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1513:       MLGetIntegerList(link, &offsets, &len);
1514:       if (len != numBdRows+1) SETERRQ(PETSC_ERR_PLIB, "Invalid boundary gradient matrix");
1515:       for(row = 0; row <= numBdRows; row++)
1516:         grad->i[row] = offsets[row]-1;
1517:       MLDisownIntegerList(link, offsets, len);
1518:       /* ml->grads[level]->j = ml[[8,level,4]] */
1519:       MLPutFunction(link, "EvaluatePacket", 1);
1520:         MLPutFunction(link, "Part", 4);
1521:           MLPutSymbol(link, "mattML");
1522:           MLPutInteger(link, 8);
1523:           MLPutInteger(link, level+1);
1524:           MLPutInteger(link, 4);
1525:       MLEndPacket(link);
1526:       PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1527:       MLGetIntegerList(link, &cols, &len);
1528:       if (len != grad->i[numBdRows]) SETERRQ(PETSC_ERR_PLIB, "Invalid boundary gradient matrix");
1529:       for(col = 0; col < len; col++)
1530:         grad->j[col] = cols[col]-1;
1531:       MLDisownIntegerList(link, cols, len);
1532:       /* ml->grads[level]->i = ml[[8,level,5]] */
1533:       MLPutFunction(link, "EvaluatePacket", 1);
1534:         MLPutFunction(link, "Part", 4);
1535:           MLPutSymbol(link, "mattML");
1536:           MLPutInteger(link, 8);
1537:           MLPutInteger(link, level+1);
1538:           MLPutInteger(link, 5);
1539:       MLEndPacket(link);
1540:       PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1541:       MLGetRealList(link, &vals, &len);
1542:       if (len != grad->i[numBdRows]) SETERRQ(PETSC_ERR_PLIB, "Invalid boundary gradient matrix");
1543:       PetscMemcpy(grad->a, vals, len * sizeof(double));
1544:       MLDisownRealList(link, vals, len);
1545:       /* Fix up all the members */
1546:       grad->nz += len;
1547:       MatAssemblyBegin(ml->grads[level], MAT_FINAL_ASSEMBLY);
1548:       MatAssemblyEnd(ml->grads[level], MAT_FINAL_ASSEMBLY);
1549:     }
1550:   } else {
1551:     PetscMalloc(1 * sizeof(double ***), &ml->factors);
1552:     PetscMalloc(1 * sizeof(Mat), &ml->grads);
1553:     PetscMalloc(1 * sizeof(Vec), &ml->bdReduceVecs);
1554:     PetscMalloc(1 * sizeof(Vec), &ml->colReduceVecs);
1555:     PetscMalloc(1 * sizeof(Vec), &ml->colReduceVecs2);
1556:   }

1558:   ml->interiorWorkLen = 1;
1559:   for(level = 0; level < ml->numLevels; level++) {
1560:     for(part = 0; part < ml->numPartitions[level]; part++)
1561:       ml->interiorWorkLen = PetscMax(ml->interiorWorkLen, ml->numPartitionRows[level][part]);
1562:   }
1563:   PetscMalloc(ml->interiorWorkLen * sizeof(double), &ml->interiorWork);
1564:   PetscMalloc(ml->interiorWorkLen * sizeof(double), &ml->interiorWork2);

1566:   return(0);
1567: #endif
1569:   return(0);
1570: }

1572: #if 0
1573: /*------------------------------ Functions for Triangular 2d Meshes -------------------------------------------------*/
1574: int PetscViewerMathematicaCreateSamplePoints_Triangular_2D(PetscViewer viewer, GVec v)
1575: {
1576: #ifdef PETSC_HAVE_MATHEMATICA
1577:   PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
1578:   MLINK               link  = vmath->link; /* The link to Mathematica */
1579:   Grid                grid;
1580:   int                 numNodes;
1581:   int                *classes;
1582:   int                *offsets;
1583:   double             *array;
1584:   int                *localStart;
1585:   int                 localOffset, comp;
1586:   int                 node, nclass;
1587:   int                 ierr;

1590:   ierr       = GVecGetGrid(v, &grid);
1591:   numNodes   = grid->mesh->numNodes;
1592:   comp       = grid->viewComp;
1593:   offsets    = grid->order->offsets;
1594:   localStart = grid->order->localStart[grid->viewField];
1595:   classes    = grid->cm->classes;

1597:   /* Attach a value to each point in the mesh -- Extra values are put in for fields not
1598:      defined on some nodes, but these values are never used */
1599:   VecGetArray(v, &array);
1600:   MLPutFunction(link, "ReplaceAll", 2);
1601:     MLPutFunction(link, "Thread", 1);
1602:       MLPutFunction(link, "f", 2);
1603:         MLPutSymbol(link, "nodes");
1604:         MLPutFunction(link, "List", numNodes);
1605:         for(node = 0; node < numNodes; node++)
1606:         {
1607:           nclass      = classes[node];
1608:           localOffset = localStart[nclass] + comp;
1609:           MLPutReal(link, array[offsets[node] + localOffset]);
1610:         }
1611:     MLPutFunction(link, "Rule", 2);
1612:       MLPutSymbol(link, "f");
1613:       MLPutSymbol(link, "Append");
1614:   VecRestoreArray(v, &array);

1616:   return(0);
1617: #endif
1619:   return(0);
1620: }

1622: int PetscViewerMathematicaCreateVectorSamplePoints_Triangular_2D(PetscViewer viewer, GVec v)
1623: {
1624: #ifdef PETSC_HAVE_MATHEMATICA
1625:   PetscViewer_Mathematica *vmath = (PetscViewer_Mathematica *) viewer->data;
1626:   MLINK               link  = vmath->link; /* The link to Mathematica */
1627:   Grid                grid;
1628:   int                 numNodes;
1629:   int                *classes;
1630:   int                *offsets;
1631:   int                *fieldClasses;
1632:   double             *array;
1633:   int                *localStart;
1634:   int                 localOffset;
1635:   int                 node, nclass;
1636:   int                 ierr;

1639:   ierr         = GVecGetGrid(v, &grid);
1640:   numNodes     = grid->mesh->numNodes;
1641:   fieldClasses = grid->cm->fieldClasses[grid->viewField];
1642:   offsets      = grid->order->offsets;
1643:   localStart   = grid->order->localStart[grid->viewField];
1644:   classes      = grid->cm->classes;

1646:   /* Make a list {{{x_0,y_0},{f^0_x,f^0_y}},...} */
1647:   VecGetArray(v, &array);
1648:   MLPutFunction(link, "ReplaceAll", 2);
1649:     MLPutFunction(link, "Thread", 1);
1650:       MLPutFunction(link, "f", 2);
1651:         MLPutSymbol(link, "nodes");
1652:         MLPutFunction(link, "List", numNodes);
1653:         for(node = 0; node < numNodes; node++)
1654:         {
1655:           nclass = classes[node];
1656:           if (fieldClasses[nclass])
1657:           {
1658:             localOffset = localStart[nclass];
1659:             MLPutFunction(link, "List", 2);
1660:               MLPutReal(link, array[offsets[node] + localOffset]);
1661:               MLPutReal(link, array[offsets[node] + localOffset + 1]);
1662:           }
1663:           else
1664:           {
1665:             /* Vectors of length zero are ignored */
1666:             MLPutFunction(link, "List", 2);
1667:               MLPutReal(link, 0.0);
1668:               MLPutReal(link, 0.0);
1669:           }
1670:         }
1671:     MLPutFunction(link, "Rule", 2);
1672:       MLPutSymbol(link, "f");
1673:       MLPutSymbol(link, "List");
1674:   VecRestoreArray(v, &array);

1676:   return(0);
1677: #endif
1679:   return(0);
1680: }

1682: int PetscViewerMathematicaCreateInterpolatedSamplePoints_Triangular_2D(PetscViewer viewer, GVec v, int vComp)
1683: {
1684: #ifdef PETSC_HAVE_MATHEMATICA
1685: #if 0
1686:   PetscViewer_Mathematica  *vmath = (PetscViewer_Mathematica *) viewer->data;
1687:   MLINK                link  = vmath->link; /* The link to Mathematica */
1688:   InterpolationContext iCtx;
1689:   Grid                 grid;
1690:   Mesh                 mesh;
1691:   double              *x, *y, *values;
1692:   long                 m, n;
1693:   double               startX, endX, incX;
1694:   double               startY, endY, incY;
1695:   int                  comp;
1696:   int                  proc, row, col;
1697:   PetscTruth           opt;
1698:   int                  ierr;
1699: #endif

1702: #if 0
1703:   ierr  = GVecGetGrid(v, &grid);
1704:   ierr  = GridGetMesh(grid, &mesh);
1705:   comp  = grid->comp[grid->viewField];

1707:   /* This sucks, but I will fix it later (It is for GridReduceInterpolationElementVec_Triangular_2D) */
1708:   grid->classesOld = grid->classes;

1710:   /* Setup InterpolationContext */
1711:   iCtx.vec         = v;
1712:   iCtx.mesh        = mesh;
1713:   iCtx.order       = grid->order;
1714:   iCtx.ghostVec    = grid->ghostVec;
1715:   iCtx.field       = grid->viewField;
1716:   iCtx.numProcs    = mesh->part->numProcs;
1717:   iCtx.rank        = mesh->part->rank;
1718:   PetscMalloc(iCtx.numProcs   * sizeof(int),      &iCtx.activeProcs);
1719:   PetscMalloc(iCtx.numProcs   * sizeof(int),      &iCtx.calcProcs);
1720:   PetscMalloc(iCtx.numProcs*3 * sizeof(PetscScalar),   &iCtx.coords);
1721:   PetscMalloc(iCtx.numProcs   * sizeof(PetscScalar *), &iCtx.values);
1722:   for(proc = 0; proc < iCtx.numProcs; proc++) {
1723:     PetscMalloc(comp * sizeof(PetscScalar), &iCtx.values[proc]);
1724:   }

1726:   /* Setup domain */
1727:   startX = 0.0;
1728:   startY = 0.0;
1729:   endX   = 1.0;
1730:   endY   = 1.0;
1731:   ierr   = PetscOptionsGetDouble("viewer_", "-math_start_x", &startX, &opt);
1732:   ierr   = PetscOptionsGetDouble("viewer_", "-math_start_y", &startY, &opt);
1733:   ierr   = PetscOptionsGetDouble("viewer_", "-math_end_x",   &endX,   &opt);
1734:   ierr   = PetscOptionsGetDouble("viewer_", "-math_end_y",   &endY,   &opt);
1735:   ierr   = PetscOptionsGetInt("viewer_", "-math_div_x", (int *) &n, &opt);
1736:   ierr   = PetscOptionsGetInt("viewer_", "-math_div_y", (int *) &m, &opt);
1737:   ierr   = PetscMalloc((n+1)      * sizeof(double), &x);
1738:   ierr   = PetscMalloc((n+1)      * sizeof(double), &y);
1739:   ierr   = PetscMalloc((n+1)*comp * sizeof(double), &values);
1740:   incX   = (endX - startX)/n;
1741:   incY   = (endY - startY)/m;

1743:   x[0] = startX;
1744:   for(col = 1; col <= n; col++)
1745:     x[col] = x[col-1] + incX;

1747:   MLPutFunction(link, "List", m+1);
1748:     for(row = 0; row <= m; row++)
1749:     {
1750:       PetscMemzero(values, (n+1)*comp * sizeof(double));
1751:       for(col = 0; col <= n; col++)
1752:         y[col] = startY + incY*row;
1753:       PointFunctionInterpolateField(n+1, comp, x, y, x, values, &iCtx);
1754:       if (vComp >= 0)
1755:       {
1756:         for(col = 0; col <= n; col++)
1757:           values[col] = values[col*comp+vComp];
1758:         MLPutRealList(link, values, n+1);
1759:       }
1760:       else
1761:       {
1762:         MLPutFunction(link, "Transpose", 1);
1763:         MLPutFunction(link, "List", 2);
1764:           MLPutFunction(link, "Transpose", 1);
1765:             MLPutFunction(link, "List", 2);
1766:               MLPutRealList(link, x, n+1);
1767:               MLPutRealList(link, y, n+1);
1768:           MLPutFunction(link, "Partition", 2);
1769:             MLPutRealList(link, values, (n+1)*comp);
1770:             MLPutInteger(link, comp);
1771:       }
1772:     }

1774:   /* This sucks, but I will fix it later (It is for GridReduceInterpolationElementVec_Triangular_2D) */
1775:   grid->classesOld = PETSC_NULL;

1777:   /* Cleanup */
1778:   PetscFree(x);
1779:   PetscFree(y);
1780:   PetscFree(values);
1781:   PetscFree(iCtx.activeProcs);
1782:   PetscFree(iCtx.calcProcs);
1783:   PetscFree(iCtx.coords);
1784:   for(proc = 0; proc < iCtx.numProcs; proc++) {
1785:     PetscFree(iCtx.values[proc]);
1786:   }
1787:   PetscFree(iCtx.values);
1788: #endif

1790:   return(0);
1791: #endif
1793:   return(0);
1794: }

1796: int PetscViewerMathematica_Mesh_Triangular_2D(PetscViewer viewer, Mesh mesh)
1797: {
1798: #ifdef PETSC_HAVE_MATHEMATICA
1799:   PetscViewer_Mathematica  *vmath = (PetscViewer_Mathematica *) viewer->data;
1800:   MLINK                link  = vmath->link;
1801:   Mesh_Triangular     *tri   = (Mesh_Triangular *) mesh->data;
1802:   int                  numCorners = mesh->numCorners;
1803:   int                  numFaces   = mesh->numFaces;
1804:   int                 *faces      = tri->faces;
1805:   int                  numNodes   = mesh->numNodes;
1806:   double              *nodes      = tri->nodes;
1807:   int                  node, face, corner;
1808:   int                  ierr;

1811:   /* Load package to view graphics in a window */
1812:   PetscViewerMathematicaLoadGraphics(viewer, vmath->graphicsType);

1814:   /* Send the node coordinates */
1815:   MLPutFunction(link, "EvaluatePacket", 1);
1816:     MLPutFunction(link, "Set", 2);
1817:       MLPutSymbol(link, "nodes");
1818:       MLPutFunction(link, "List", numNodes);
1819:       for(node = 0; node < numNodes; node++) {
1820:         MLPutFunction(link, "List", 2);
1821:           MLPutReal(link, nodes[node*2]);
1822:           MLPutReal(link, nodes[node*2+1]);
1823:       }
1824:   MLEndPacket(link);
1825:   /* Skip packets until ReturnPacket */
1826:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1827:   /* Skip ReturnPacket */
1828:   MLNewPacket(link);

1830:   /* Send the faces */
1831:   MLPutFunction(link, "EvaluatePacket", 1);
1832:     MLPutFunction(link, "Set", 2);
1833:       MLPutSymbol(link, "faces");
1834:       MLPutFunction(link, "List", numFaces);
1835:       for(face = 0; face < numFaces; face++) {
1836:         MLPutFunction(link, "List", numCorners);
1837:         for(corner = 0; corner < numCorners; corner++) {
1838:           MLPutReal(link, faces[face*numCorners+corner]);
1839:         }
1840:       }
1841:   MLEndPacket(link);
1842:   /* Skip packets until ReturnPacket */
1843:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1844:   /* Skip ReturnPacket */
1845:   MLNewPacket(link);

1847:   return(0);
1848: #endif
1850:   return(0);
1851: }

1853: int PetscViewerMathematicaCheckMesh_Triangular_2D(PetscViewer viewer, Mesh mesh)
1854: {
1855: #ifdef PETSC_HAVE_MATHEMATICA
1856:   PetscViewer_Mathematica  *vmath = (PetscViewer_Mathematica *) viewer->data;
1857:   MLINK                link  = vmath->link;
1858:   int                  numCorners = mesh->numCorners;
1859:   int                  numFaces   = mesh->numFaces;
1860:   int                  numNodes   = mesh->numNodes;
1861:   const char          *symbol;
1862:   long                 args;
1863:   int                  dim, type;
1864:   PetscTruth           match;
1865:   int                  ierr;

1868:   /* Check that nodes are defined */
1869:   MLPutFunction(link, "EvaluatePacket", 1);
1870:     MLPutFunction(link, "ValueQ", 1);
1871:       MLPutSymbol(link, "nodes");
1872:   MLEndPacket(link);
1873:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1874:   MLGetSymbol(link, &symbol);
1875:   PetscStrcmp("True", (char *) symbol, &match);
1876:   if (match == PETSC_FALSE) {
1877:     MLDisownSymbol(link, symbol);
1878:     goto redefineMesh;
1879:   }
1880:   MLDisownSymbol(link, symbol);

1882:   /* Check the dimensions */
1883:   MLPutFunction(link, "EvaluatePacket", 1);
1884:     MLPutFunction(link, "MatrixQ", 2);
1885:       MLPutSymbol(link, "nodes");
1886:       MLPutSymbol(link, "NumberQ");
1887:   MLEndPacket(link);
1888:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1889:   MLGetSymbol(link, &symbol);
1890:   PetscStrcmp("True", (char *) symbol, &match);
1891:   if (match == PETSC_FALSE) {
1892:     MLDisownSymbol(link, symbol);
1893:     goto redefineMesh;
1894:   }
1895:   MLDisownSymbol(link, symbol);
1896:   MLPutFunction(link, "EvaluatePacket", 1);
1897:     MLPutFunction(link, "Dimensions", 1);
1898:       MLPutSymbol(link, "nodes");
1899:   MLEndPacket(link);
1900:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1901:   args = 0;
1902:   type = MLGetNext(link);
1903:   MLGetArgCount(link, &args);
1904:   if (args != 2) {
1905:     MLNewPacket(link);
1906:     goto redefineMesh;
1907:   }
1908:   MLGetSymbol(link, &symbol);
1909:   MLDisownSymbol(link, symbol);
1910:   MLGetInteger(link, &dim);
1911:   if (dim != numNodes) {
1912:     MLNewPacket(link);
1913:     goto redefineMesh;
1914:   }
1915:   MLGetInteger(link, &dim);
1916:   if (dim != mesh->dim) {
1917:     MLNewPacket(link);
1918:     goto redefineMesh;
1919:   }

1921:   /* Check that faces are defined */
1922:   MLPutFunction(link, "EvaluatePacket", 1);
1923:     MLPutFunction(link, "ValueQ", 1);
1924:       MLPutSymbol(link, "faces");
1925:   MLEndPacket(link);
1926:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1927:   MLGetSymbol(link, &symbol);
1928:   PetscStrcmp("True", (char *) symbol, &match);
1929:   if (match == PETSC_FALSE) {
1930:     MLDisownSymbol(link, symbol);
1931:     goto redefineMesh;
1932:   }
1933:   MLDisownSymbol(link, symbol);

1935:   /* Check the dimensions */
1936:   MLPutFunction(link, "EvaluatePacket", 1);
1937:     MLPutFunction(link, "MatrixQ", 2);
1938:       MLPutSymbol(link, "faces");
1939:       MLPutSymbol(link, "NumberQ");
1940:   MLEndPacket(link);
1941:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1942:   MLGetSymbol(link, &symbol);
1943:   PetscStrcmp("True", (char *) symbol, &match);
1944:   if (match == PETSC_FALSE) {
1945:     MLDisownSymbol(link, symbol);
1946:     goto redefineMesh;
1947:   }
1948:   MLDisownSymbol(link, symbol);
1949:   MLPutFunction(link, "EvaluatePacket", 1);
1950:     MLPutFunction(link, "Dimensions", 1);
1951:       MLPutSymbol(link, "faces");
1952:   MLEndPacket(link);
1953:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
1954:   args = 0;
1955:   type = MLGetNext(link);
1956:   MLGetArgCount(link, &args);
1957:   if (args != 2) {
1958:     MLNewPacket(link);
1959:     goto redefineMesh;
1960:   }
1961:   MLGetSymbol(link, &symbol);
1962:   MLDisownSymbol(link, symbol);
1963:   MLGetInteger(link, &dim);
1964:   if (dim != numFaces) {
1965:     MLNewPacket(link);
1966:     goto redefineMesh;
1967:   }
1968:   MLGetInteger(link, &dim);
1969:   if (dim != numCorners) {
1970:     MLNewPacket(link);
1971:     goto redefineMesh;
1972:   }

1974:   return(0);

1976:  redefineMesh:
1977:   MeshView(mesh, viewer);
1978:   return(0);
1979: #endif
1981:   return(0);
1982: }

1984: int PetscViewerMathematicaTriangulationPlot_Triangular_2D(PetscViewer viewer, GVec v)
1985: {
1986: #ifdef PETSC_HAVE_MATHEMATICA
1987:   PetscViewer_Mathematica  *vmath = (PetscViewer_Mathematica *) viewer->data;
1988:   MLINK                link  = vmath->link; /* The link to Mathematica */
1989:   Mesh                 mesh;
1990:   Grid                 grid;
1991:   int                  numCorners;
1992:   int                  ierr;

1995:   ierr       = GVecGetGrid(v, &grid);
1996:   mesh       = grid->mesh;
1997:   numCorners = mesh->numCorners;

1999:   MLPutFunction(link, "Show", 2);
2000:     MLPutFunction(link, "Graphics3D", 1);
2001:       MLPutFunction(link, "Flatten", 1);
2002:         MLPutFunction(link, "Map", 2);
2003:           MLPutFunction(link, "Function", 1);
2004:           if ((numCorners == 6) && ((grid->cm->numClasses == 1) || (grid->cm->fieldClasses[grid->viewField][1])))
2005:           {
2006:             MLPutFunction(link, "List", 4);
2007:             /* Triangle 0--5--4 */
2008:             MLPutFunction(link, "Polygon", 1);
2009:               MLPutFunction(link, "Part", 2);
2010:                 MLPutSymbol(link, "points");
2011:                 MLPutFunction(link, "Plus", 2);
2012:                   MLPutFunction(link, "Part", 2);
2013:                     MLPutFunction(link, "Slot", 1);
2014:                       MLPutInteger(link, 1);
2015:                     MLPutFunction(link, "List", 3);
2016:                       MLPutInteger(link, 1);
2017:                       MLPutInteger(link, 6);
2018:                       MLPutInteger(link, 5);
2019:                   MLPutInteger(link, 1);
2020:             /* Triangle 1--3--5 */
2021:             MLPutFunction(link, "Polygon", 1);
2022:               MLPutFunction(link, "Part", 2);
2023:                 MLPutSymbol(link, "points");
2024:                 MLPutFunction(link, "Plus", 2);
2025:                   MLPutFunction(link, "Part", 2);
2026:                     MLPutFunction(link, "Slot", 1);
2027:                       MLPutInteger(link, 1);
2028:                     MLPutFunction(link, "List", 3);
2029:                       MLPutInteger(link, 2);
2030:                       MLPutInteger(link, 4);
2031:                       MLPutInteger(link, 6);
2032:                   MLPutInteger(link, 1);
2033:             /* Triangle 2--4--3 */
2034:             MLPutFunction(link, "Polygon", 1);
2035:               MLPutFunction(link, "Part", 2);
2036:                 MLPutSymbol(link, "points");
2037:                 MLPutFunction(link, "Plus", 2);
2038:                   MLPutFunction(link, "Part", 2);
2039:                     MLPutFunction(link, "Slot", 1);
2040:                       MLPutInteger(link, 1);
2041:                     MLPutFunction(link, "List", 3);
2042:                       MLPutInteger(link, 3);
2043:                       MLPutInteger(link, 5);
2044:                       MLPutInteger(link, 4);
2045:                   MLPutInteger(link, 1);
2046:             /* Triangle 3--4--5 */
2047:             MLPutFunction(link, "Polygon", 1);
2048:               MLPutFunction(link, "Part", 2);
2049:                 MLPutSymbol(link, "points");
2050:                 MLPutFunction(link, "Plus", 2);
2051:                   MLPutFunction(link, "Part", 2);
2052:                     MLPutFunction(link, "Slot", 1);
2053:                       MLPutInteger(link, 1);
2054:                     MLPutFunction(link, "List", 3);
2055:                       MLPutInteger(link, 4);
2056:                       MLPutInteger(link, 5);
2057:                       MLPutInteger(link, 6);
2058:                   MLPutInteger(link, 1);
2059:           }
2060:           else if ((numCorners == 3) || (!grid->cm->fieldClasses[grid->viewField][1]))
2061:           {
2062:             /* Triangle 0--1--2 */
2063:             MLPutFunction(link, "Polygon", 1);
2064:               MLPutFunction(link, "Part", 2);
2065:                 MLPutSymbol(link, "points");
2066:                 MLPutFunction(link, "Plus", 2);
2067:                   MLPutFunction(link, "Part", 2);
2068:                     MLPutFunction(link, "Slot", 1);
2069:                       MLPutInteger(link, 1);
2070:                     MLPutFunction(link, "List", 3);
2071:                       MLPutInteger(link, 1);
2072:                       MLPutInteger(link, 2);
2073:                       MLPutInteger(link, 3);
2074:                   MLPutInteger(link, 1);
2075:           } else {
2076:             SETERRQ(PETSC_ERR_ARG_WRONG, "Invalid number of local nodes");
2077:           }
2078:           MLPutSymbol(link, "faces");
2079:     MLPutFunction(link, "Rule", 2);
2080:       MLPutSymbol(link, "AspectRatio");
2081:       MLPutReal(link, 1.0);
2082:   return(0);
2083: #endif
2085:   return(0);
2086: }

2088: int PetscViewerMathematicaVectorPlot_Triangular_2D(PetscViewer viewer, GVec v)
2089: {
2090: #ifdef PETSC_HAVE_MATHEMATICA
2091:   PetscViewer_Mathematica  *vmath = (PetscViewer_Mathematica *) viewer->data;
2092:   MLINK                link  = vmath->link; /* The link to Mathematica */
2093:   Grid                 grid;
2094:   Mesh                 mesh;
2095:   PetscReal            scale;
2096:   PetscTruth           opt;
2097:   int                  ierr;

2100:   GVecGetGrid(v, &grid);
2101:   GridGetMesh(grid, &mesh);

2103:   MLPutFunction(link, "ListPlotVectorField", 3);
2104:     MLPutSymbol(link, "points");
2105:     MLPutFunction(link, "Rule", 2);
2106:       MLPutSymbol(link, "AspectRatio");
2107:       MLPutReal(link, mesh->sizeY/mesh->sizeX);
2108:     MLPutFunction(link, "Rule", 2);
2109:       MLPutSymbol(link, "ScaleFactor");
2110:       PetscOptionsGetReal("viewer_", "-math_scale", &scale, &opt);
2111:       if (opt == PETSC_TRUE) {
2112:         MLPutReal(link, scale);
2113:       } else {
2114:         MLPutSymbol(link, "None");
2115:       }
2116:   return(0);
2117: #endif
2119:   return(0);
2120: }

2122: int PetscViewerMathematicaSurfacePlot_Triangular_2D(PetscViewer viewer, GVec v)
2123: {
2124: #ifdef PETSC_HAVE_MATHEMATICA
2125:   PetscViewer_Mathematica  *vmath = (PetscViewer_Mathematica *) viewer->data;
2126:   MLINK                link  = vmath->link; /* The link to Mathematica */
2127:   PetscReal            startX, endX;
2128:   PetscReal            startY, endY;
2129:   PetscTruth           opt;
2130:   int                  ierr;

2133:   /* Setup domain */
2134:   startX = 0.0;
2135:   startY = 0.0;
2136:   endX   = 1.0;
2137:   endY   = 1.0;
2138:   ierr   = PetscOptionsGetReal("viewer_", "-math_start_x", &startX, &opt);
2139:   ierr   = PetscOptionsGetReal("viewer_", "-math_start_y", &startY, &opt);
2140:   ierr   = PetscOptionsGetReal("viewer_", "-math_end_x",   &endX,   &opt);
2141:   ierr   = PetscOptionsGetReal("viewer_", "-math_end_y",   &endY,   &opt);

2143:   MLPutFunction(link, "Show", 1);
2144:     MLPutFunction(link, "SurfaceGraphics", 6);
2145:       MLPutSymbol(link, "points");
2146:       MLPutFunction(link, "Rule", 2);
2147:         MLPutSymbol(link, "ClipFill");
2148:         MLPutSymbol(link, "None");
2149:       MLPutFunction(link, "Rule", 2);
2150:         MLPutSymbol(link, "Axes");
2151:         MLPutSymbol(link, "True");
2152:       MLPutFunction(link, "Rule", 2);
2153:         MLPutSymbol(link, "PlotLabel");
2154:         MLPutString(link, vmath->objName);
2155:       MLPutFunction(link, "Rule", 2);
2156:         MLPutSymbol(link, "MeshRange");
2157:         MLPutFunction(link, "List", 2);
2158:           MLPutFunction(link, "List", 2);
2159:             MLPutReal(link, startX);
2160:             MLPutReal(link, endX);
2161:           MLPutFunction(link, "List", 2);
2162:             MLPutReal(link, startY);
2163:             MLPutReal(link, endY);
2164:       MLPutFunction(link, "Rule", 2);
2165:         MLPutSymbol(link, "AspectRatio");
2166:         MLPutReal(link, (endY - startY)/(endX - startX));
2167:   return(0);
2168: #endif
2170:   return(0);
2171: }

2173: int PetscViewerMathematica_GVec_Triangular_2D(PetscViewer viewer, GVec v)
2174: {
2175:   PetscViewer_Mathematica  *vmath = (PetscViewer_Mathematica *) viewer->data;
2176: #ifdef PETSC_HAVE_MATHEMATICA
2177:   MLINK                link  = vmath->link; /* The link to Mathematica */
2178:   Mesh                 mesh;
2179:   Grid                 grid;
2180:   Mat                  P;
2181:   GVec                 w;
2182:   int                  numCorners;
2183:   int                  ierr;

2186:   ierr       = GVecGetGrid(v, &grid);
2187:   mesh       = grid->mesh;
2188:   numCorners = mesh->numCorners;

2190:   /* Check that a field component has been specified */
2191:   if ((grid->viewField < 0) || (grid->viewField >= grid->numFields)) return(0);

2193:   if (grid->isConstrained) {
2194:     GVecCreate(grid, &w);
2195:     GridGetConstraintMatrix(grid, &P);
2196:     MatMult(P, v, w);
2197:   } else {
2198:     w = v;
2199:   }

2201:   /* Check that the mesh is defined correctly */
2202:   PetscViewerMathematicaCheckMesh_Triangular_2D(viewer, mesh);

2204:   /* Send the first component of the first field */
2205:   MLPutFunction(link, "EvaluatePacket", 1);
2206:   switch(vmath->plotType)
2207:   {
2208:   case MATHEMATICA_TRIANGULATION_PLOT:
2209:     MLPutFunction(link, "Module", 2);
2210:       /* Make temporary points with each value of the field component in the vector */
2211:       MLPutFunction(link, "List", 2);
2212:         MLPutSymbol(link, "f");
2213:         MLPutFunction(link, "Set", 2);
2214:           MLPutSymbol(link, "points");
2215:           PetscViewerMathematicaCreateSamplePoints_Triangular_2D(viewer, w);
2216:       /* Display the picture */
2217:       PetscViewerMathematicaTriangulationPlot_Triangular_2D(viewer, w);
2218:       break;
2219:   case MATHEMATICA_VECTOR_TRIANGULATION_PLOT:
2220:     if (grid->fields[grid->viewField].numComp != 2) {
2221:       SETERRQ(PETSC_ERR_ARG_WRONG, "Field must be a 2D vector field for this plot type");
2222:     }
2223:     MLPutFunction(link, "Module", 2);
2224:       /* Make temporary list with points and 2D vector field values */
2225:       MLPutFunction(link, "List", 2);
2226:         MLPutSymbol(link, "f");
2227:         MLPutFunction(link, "Set", 2);
2228:           MLPutSymbol(link, "points");
2229:           PetscViewerMathematicaCreateVectorSamplePoints_Triangular_2D(viewer, w);
2230:       /* Display the picture */
2231:       PetscViewerMathematicaVectorPlot_Triangular_2D(viewer, w);
2232:       break;
2233:   case MATHEMATICA_VECTOR_PLOT:
2234:     if (grid->fields[grid->viewField].numComp != 2) {
2235:       SETERRQ(PETSC_ERR_ARG_WRONG, "Field must be a 2D vector field for this plot type");
2236:     }
2237:     MLPutFunction(link, "Module", 2);
2238:       /* Make temporary list with points and 2D vector field values */
2239:       MLPutFunction(link, "List", 2);
2240:         MLPutSymbol(link, "f");
2241:         MLPutFunction(link, "Set", 2);
2242:           MLPutSymbol(link, "points");
2243:           PetscViewerMathematicaCreateInterpolatedSamplePoints_Triangular_2D(viewer, w, -1);
2244:       /* Display the picture */
2245:       PetscViewerMathematicaVectorPlot_Triangular_2D(viewer, w);
2246:       break;
2247:   case MATHEMATICA_SURFACE_PLOT:
2248:     if (grid->fields[grid->viewField].numComp != 2) {
2249:       SETERRQ(PETSC_ERR_ARG_WRONG, "Field must be a 2D vector field for this plot type");
2250:     }
2251:     MLPutFunction(link, "Module", 2);
2252:       /* Make temporary list with interpolated field values on a square mesh */
2253:       MLPutFunction(link, "List", 1);
2254:         MLPutFunction(link, "Set", 2);
2255:           MLPutSymbol(link, "points");
2256:           PetscViewerMathematicaCreateInterpolatedSamplePoints_Triangular_2D(viewer, w, grid->viewComp);
2257:       /* Display the picture */
2258:       PetscViewerMathematicaSurfacePlot_Triangular_2D(viewer, w);
2259:       break;
2260:   default:
2261:     SETERRQ(PETSC_ERR_ARG_WRONG, "Invalid plot type");
2262:   }
2263:   MLEndPacket(link);
2264:   /* Skip packets until ReturnPacket */
2265:   PetscViewerMathematicaSkipPackets(viewer, RETURNPKT);
2266:   /* Skip ReturnPacket */
2267:   MLNewPacket(link);

2269:   if (grid->isConstrained) {
2270:     VecDestroy(w);
2271:   }
2272: #else
2274:   switch(vmath->plotType)
2275:   {
2276:   case MATHEMATICA_TRIANGULATION_PLOT:
2277:       break;
2278:   case MATHEMATICA_VECTOR_TRIANGULATION_PLOT:
2279:       break;
2280:   case MATHEMATICA_VECTOR_PLOT:
2281:       break;
2282:   case MATHEMATICA_SURFACE_PLOT:
2283:       break;
2284:   default:
2285:     SETERRQ(PETSC_ERR_ARG_WRONG, "Invalid plot type");
2286:   }
2287: #endif
2288:   return(0);
2289: }
2290: #endif