Actual source code: fwk.c
  2: #include <petscsys.h>
  3: #include <petscfwk.h>
  5: /* FIX: is it okay to include this directly? */
  6: #include <ctype.h>
  8:  PetscClassId PETSC_FWK_CLASSID;
  9: static char PETSC_FWK_CLASS_NAME[] = "PetscFwk";
 10: static PetscBool  PetscFwkPackageInitialized = PETSC_FALSE;
 13: typedef PetscErrorCode (*PetscFwkPythonLoadVTableFunction)(PetscFwk fwk, const char* path, const char* name, void **vtable);
 14: typedef PetscErrorCode (*PetscFwkPythonClearVTableFunction)(PetscFwk fwk, void **vtable);
 15: typedef PetscErrorCode (*PetscFwkPythonCallFunction)(PetscFwk fwk, const char* message, void *vtable);
 18: PetscFwkPythonLoadVTableFunction      PetscFwkPythonLoadVTable      = PETSC_NULL;
 19: PetscFwkPythonClearVTableFunction     PetscFwkPythonClearVTable     = PETSC_NULL;
 20: PetscFwkPythonCallFunction            PetscFwkPythonCall            = PETSC_NULL;
 23: #define PETSC_FWK_CHECKINIT_PYTHON()                                        \
 24:   if(PetscFwkPythonLoadVTable == PETSC_NULL) {                                \
 26:     PetscPythonInitialize(PETSC_NULL,PETSC_NULL);        \
 27:     if(PetscFwkPythonLoadVTable == PETSC_NULL) {                                \
 28:       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,                                \
 29:               "Couldn't initialize Python support for PetscFwk");        \
 30:     }                                                                        \
 31:   }                                                                        
 32: 
 33: #define PETSC_FWK_LOAD_VTABLE_PYTHON(fwk, path, name)                   \
 34:   PETSC_FWK_CHECKINIT_PYTHON();                                                \
 35:   {                                                                        \
 37:     PetscFwkPythonLoadVTable(fwk, path, name, &(fwk->vtable));   \
 38:     if (ierr) { PetscPythonPrintError(); SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB, "Python error"); } \
 39:   }
 41: #define PETSC_FWK_CLEAR_VTABLE_PYTHON(fwk)                              \
 42:   PETSC_FWK_CHECKINIT_PYTHON();                                                \
 43:   {                                                                        \
 45:     PetscFwkPythonClearVTable(fwk, &(fwk->vtable));              \
 46:     if (ierr) { PetscPythonPrintError(); SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB, "Python error"); } \
 47:   }
 48: 
 49: #define PETSC_FWK_CALL_PYTHON(fwk, message)                             \
 50:   PETSC_FWK_CHECKINIT_PYTHON();                                         \
 51:   {                                                                        \
 53:     PetscFwkPythonCall(fwk, message, fwk->vtable);                           \
 54:     if (ierr) { PetscPythonPrintError(); SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB, "Python error"); } \
 55:   }
 56: 
 58: /* ---------------------------------------------------------------------------------------------- */
 59: struct _n_PetscFwkGraph {
 60:   PetscInt vcount, vmax; /* actual and allocated number of vertices */
 61:   PetscInt *i, *j, *outdegree; /* (A)IJ structure for the underlying matrix: 
 62:                                   i[row]             the row offset into j, 
 63:                                   i[row+1] - i[row]  allocated number of entries for row,
 64:                                   outdegree[row]     actual number of entries in row
 65:                                */
 66:   PetscInt *indegree;
 67:   PetscInt nz, maxnz;
 68:   PetscInt rowreallocs, colreallocs;
 69: };
 71: typedef struct _n_PetscFwkGraph *PetscFwkGraph;
 73: #define CHUNKSIZE 5
 74: /*
 75:     Inserts the (row,col) entry, allocating larger arrays, if necessary. 
 76:     Does NOT check whether row and col are within range (< graph->vcount).
 77:     Does NOT check whether the entry is already present.
 78: */
 79: #undef  __FUNCT__
 81: PetscErrorCode PetscFwkGraphExpandRow_Private(PetscFwkGraph graph, PetscInt row)
 82: {
 84:   PetscInt       rowlen, rowmax, rowoffset;
 85:   PetscInt       ii;
 88:   rowlen = graph->outdegree[row];
 89:   rowmax = graph->i[row+1] - graph->i[row];
 90:   rowoffset = graph->i[row];
 91:   if (rowlen >= rowmax) {
 92:     /* there is no extra room in row, therefore enlarge */
 93:     PetscInt   new_nz = graph->i[graph->vcount] + CHUNKSIZE;
 94:     PetscInt   *new_i=0,*new_j=0;
 95: 
 96:     /* malloc new storage space */
 97:     PetscMalloc(new_nz*sizeof(PetscInt),&new_j);
 98:     PetscMalloc((graph->vmax+1)*sizeof(PetscInt),&new_i);
 99: 
100:     /* copy over old data into new slots */
101:     for (ii=0; ii<row+1; ii++) {new_i[ii] = graph->i[ii];}
102:     for (ii=row+1; ii<graph->vmax+1; ii++) {new_i[ii] = graph->i[ii]+CHUNKSIZE;}
103:     PetscMemcpy(new_j,graph->j,(rowoffset+rowlen)*sizeof(PetscInt));
104:     PetscMemcpy(new_j+rowoffset+rowlen+CHUNKSIZE,graph->j+rowoffset+rowlen,(new_nz - CHUNKSIZE - rowoffset - rowlen)*sizeof(PetscInt));
105:     /* free up old matrix storage */
106:     PetscFree(graph->j);
107:     PetscFree(graph->i);
108:     graph->i = new_i; graph->j = new_j;
109:     graph->maxnz     += CHUNKSIZE;
110:     graph->colreallocs++;
111:   }
112:   return(0);
113: }/* PetscFwkGraphExpandRow_Private() */
115: #undef  __FUNCT__
117: PetscErrorCode PetscFwkGraphAddVertex(PetscFwkGraph graph, PetscInt *v) {
118:   PetscInt ii;
121:   if(graph->vcount >= graph->vmax) {
122:     /* Add rows */
123:     PetscInt   *new_i=0, *new_outdegree=0, *new_indegree;
124: 
125:     /* malloc new storage space */
126:     PetscMalloc((graph->vmax+CHUNKSIZE+1)*sizeof(PetscInt),&new_i);
127:     PetscMalloc((graph->vmax+CHUNKSIZE)*sizeof(PetscInt),&new_outdegree);
128:     PetscMalloc((graph->vmax+CHUNKSIZE)*sizeof(PetscInt),&new_indegree);
129:     PetscMemzero(new_outdegree, (graph->vmax+CHUNKSIZE)*sizeof(PetscInt));
130:     PetscMemzero(new_indegree, (graph->vmax+CHUNKSIZE)*sizeof(PetscInt));
133:     /* copy over old data into new slots */
134:     PetscMemcpy(new_i,graph->i,(graph->vmax+1)*sizeof(PetscInt));
135:     PetscMemcpy(new_outdegree,graph->outdegree,(graph->vmax)*sizeof(PetscInt));
136:     PetscMemcpy(new_indegree,graph->indegree,(graph->vmax)*sizeof(PetscInt));
137:     for (ii=graph->vmax+1; ii<=graph->vmax+CHUNKSIZE; ++ii) {
138:       new_i[ii] = graph->i[graph->vmax];
139:     }
141:     /* free up old matrix storage */
142:     PetscFree(graph->i);
143:     PetscFree(graph->outdegree);
144:     PetscFree(graph->indegree);
145:     graph->i = new_i; graph->outdegree = new_outdegree; graph->indegree = new_indegree;
146: 
147:     graph->vmax += CHUNKSIZE;
148:     graph->rowreallocs++;
149:   }
150:   if(v) {*v = graph->vcount;}
151:   ++(graph->vcount);
152:   return(0);
153: }/* PetscFwkGraphAddVertex() */
155: #undef  __FUNCT__
157: PetscErrorCode PetscFwkGraphAddEdge(PetscFwkGraph graph, PetscInt row, PetscInt col) {
158:   PetscErrorCode        ierr;
159:   PetscInt              *rp,low,high,t,ii,i;
162:   if(row < 0 || row >= graph->vcount) {
163:     SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Source vertex %D out of range: min %D max %D",row, 0, graph->vcount);
164:   }
165:   if(col < 0 || col >= graph->vcount) {
166:     SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Target vertex %D out of range: min %D max %D",col, 0, graph->vcount);
167:   }
168:   rp   = graph->j + graph->i[row];
169:   low  = 0;
170:   high = graph->outdegree[row];
171:   while (high-low > 5) {
172:     t = (low+high)/2;
173:     if (rp[t] > col) high = t;
174:     else             low  = t;
175:   }
176:   for (i=low; i<high; ++i) {
177:     if (rp[i] > col) break;
178:     if (rp[i] == col) {
179:       goto we_are_done;
180:     }
181:   }
182:   PetscFwkGraphExpandRow_Private(graph, row);
183:   /* 
184:      If the graph was empty before, graph->j was NULL and rp was NULL as well.  
185:      Now that the row has been expanded, rp needs to be reset. 
186:   */
187:   rp = graph->j + graph->i[row];
188:   /* shift up all the later entries in this row */
189:   for (ii=graph->outdegree[row]; ii>=i; --ii) {
190:     rp[ii+1] = rp[ii];
191:   }
192:   rp[i] = col;
193:   ++(graph->outdegree[row]);
194:   ++(graph->indegree[col]);
195:   ++(graph->nz);
196: 
197:  we_are_done:
198:   return(0);
199: }/* PetscFwkGraphAddEdge() */
202: #undef  __FUNCT__
204: PetscErrorCode PetscFwkGraphTopologicalSort(PetscFwkGraph graph, PetscInt *n, PetscInt **queue) {
205:   PetscBool  *queued;
206:   PetscInt   *indegree;
207:   PetscInt ii, k, jj, Nqueued = 0;
208:   PetscBool  progress;
211:   if(!n || !queue) {
212:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG, "Invalid return argument pointers n or vertices");
213:   }
214:   *n = graph->vcount;
215:   PetscMalloc(sizeof(PetscInt)*graph->vcount, queue);
216:   PetscMalloc(sizeof(PetscBool)*graph->vcount, &queued);
217:   PetscMalloc(sizeof(PetscInt)*graph->vcount, &indegree);
218:   for(ii = 0; ii < graph->vcount; ++ii) {
219:     queued[ii]   = PETSC_FALSE;
220:     indegree[ii] = graph->indegree[ii];
221:   }
222:   while(Nqueued < graph->vcount) {
223:     progress = PETSC_FALSE;
224:     for(ii = 0; ii < graph->vcount; ++ii) {
225:       /* If ii is not queued yet, and the indegree is 0, queue it. */
226:       if(!queued[ii] && !indegree[ii]) {
227:         (*queue)[Nqueued] = ii;
228:         queued[ii] = PETSC_TRUE;
229:         ++Nqueued;
230:         progress = PETSC_TRUE;
231:         /* Reduce the indegree of all vertices in row ii */
232:         for(k = 0; k < graph->outdegree[ii]; ++k) {
233:           jj = graph->j[graph->i[ii]+k];
234:           --(indegree[jj]);
235:           /* 
236:              It probably would be more efficient to make a recursive call to the body of the for loop 
237:              with the jj in place of ii, but we use a simple-minded algorithm instead, since the graphs
238:              we anticipate encountering are tiny. 
239:           */
240:         }/* for(k) */
241:       }/* if(!queued) */
242:     }/* for(ii) */
243:     /* If no progress was made during this iteration, the graph must have a cycle */
244:     if(!progress) {
245:       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Cycle detected in the dependency graph");
246:     }
247:   }/* while(Nqueued) */
248:   PetscFree(queued);
249:   PetscFree(indegree);
250:   return(0);
251: }/* PetscFwkGraphTopologicalSort() */
253: #undef  __FUNCT__
255: PetscErrorCode PetscFwkGraphDestroy(PetscFwkGraph graph)
256: {
259:   PetscFree(graph->i);
260:   PetscFree(graph->j);
261:   PetscFree(graph->outdegree);
262:   PetscFree(graph->indegree);
263:   PetscFree(graph);
264:   return(0);
265: }/* PetscFwkGraphDestroy() */
267: #undef  __FUNCT__
269: PetscErrorCode PetscFwkGraphCreate(PetscFwkGraph *graph_p)
270: {
271:   PetscFwkGraph graph;
274:   PetscNew(struct _n_PetscFwkGraph, graph_p);
275:   graph = *graph_p;
276:   graph->vcount = graph->vmax = graph->nz = graph->maxnz = graph->rowreallocs = graph->colreallocs = 0;
277:   PetscMalloc(sizeof(PetscInt), &(graph->i));
278:   graph->j = graph->outdegree = graph->indegree = PETSC_NULL;
279:   return(0);
280: }/* PetscFwkGraphCreate() */
284: /* ------------------------------------------------------------------------------------------------------- */
286: typedef enum{PETSC_FWK_VTABLE_NONE, PETSC_FWK_VTABLE_SO, PETSC_FWK_VTABLE_PY} PetscFwkVTableType;
288: struct _n_PetscFwkVTable_SO {
289:   char           *path, *name;
290: };
293: struct _p_PetscFwk {
294:   PETSCHEADER(int);
295:   PetscFwkVTableType        vtable_type;
296:   void                      *vtable;
297:   char *                    url;
298:   PetscFwk                  parent;
299:   PetscInt                  N, maxN;
300:   PetscFwk                  *component;
301:   PetscFwkGraph             dep_graph;
302: };
304: static PetscFwk defaultFwk = PETSC_NULL;
308: /* ------------------------------------------------------------------------------------------------------- */
310: typedef PetscErrorCode (*PetscFwkCallFunction)(PetscFwk, const char*);
311: typedef PetscErrorCode (*PetscFwkMessageFunction)(PetscFwk);
312: typedef void (*QueryFunction)(void);
314: static PetscDLLibrary PetscFwkDLLibrariesLoaded = 0;
316: #undef  __FUNCT__
318: PetscErrorCode PetscFwkCall_SO(PetscFwk fwk, const char* path, const char* name, const char* message) {
319:   size_t    namelen, messagelen, msgfunclen, callfunclen;
320:   char *msgfunc  = PETSC_NULL, *callfunc = PETSC_NULL;
321:   PetscFwkCallFunction call = PETSC_NULL;
322:   PetscFwkMessageFunction msg = PETSC_NULL;
325:   PetscStrlen(name, &namelen);
326:   PetscStrlen(message, &messagelen);
327:   msgfunclen = namelen + messagelen;
328:   PetscMalloc(sizeof(char)*(msgfunclen+1), &msgfunc);
329:   msgfunc[0] = '\0';
330:   if(namelen){
331:     PetscStrcat(msgfunc, name);
332:   }
333:   PetscStrcat(msgfunc, message);
334:   if(namelen) {
335:     /* HACK: is 'toupper' part of the C standard? Looks like starting with C89. */
336:     msgfunc[namelen] = toupper(msgfunc[namelen]);
337:   }
338:   PetscDLLibrarySym(((PetscObject)fwk)->comm, &PetscFwkDLLibrariesLoaded, path, msgfunc, (void **)(&msg));
339:   PetscFree(msgfunc);
340:   if(msg) {
341:     (*msg)(fwk);
342:     return(0);
343:   }
344:   callfunclen        = namelen+4;
345:   PetscMalloc(sizeof(char)*(callfunclen+1), &callfunc);
346:   if(namelen) {
347:     PetscStrcpy(callfunc, name);
348:   }
349:   if(namelen){
350:     PetscStrcat(callfunc, "Call");
351:   }
352:   else {
353:     PetscStrcat(callfunc, "call");
354:   }
355:   PetscDLLibrarySym(((PetscObject)fwk)->comm, &PetscFwkDLLibrariesLoaded, path, callfunc, (void**)(&call));
356:   PetscFree(callfunc);
357:   if(call) {
358:     (*call)(fwk, message);
359:     return(0);
360:   }
361:   SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "PetscFwk '%s' cannot execute '%s'", ((PetscObject)fwk)->name, message);
362:   return(0);
363: }/* PetscFwkCall_SO() */
366: #undef  __FUNCT__
368: PetscErrorCode PetscFwkCall_NONE(PetscFwk fwk, const char* message) {
369:   PetscFwkCallFunction call = PETSC_NULL;
370:   PetscFwkMessageFunction msg = PETSC_NULL;
373:   PetscFListFind(((PetscObject)fwk)->qlist, ((PetscObject)fwk)->comm, message,PETSC_FALSE, (QueryFunction*)(&msg));
374:   if(msg) {
375:     (*msg)(fwk);
376:     return(0);
377:   }
378:   PetscFListFind(((PetscObject)fwk)->qlist, ((PetscObject)fwk)->comm, "call",PETSC_FALSE, (QueryFunction*)(&call));
379:   if(call) {
380:     (*call)(fwk, message);
381:     return(0);
382:   }
383:   SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "PetscFwk '%s' cannot execute '%s'", ((PetscObject)fwk)->name, message);
384:   return(0);
385: }/* PetscFwkCall_NONE() */
388: #undef  __FUNCT__
390: /*@C 
391:    PetscFwkCall -- send a string message to a PetscFwk object.  
393:    Logically collective on PetscFwk.
395:    Input paramters:
396: +  fwk     -- a PetscFwk object
397: -  message -- a character string
399:    Notes: In response to the message the object performs actions defined by its URL (see PetscFwkSetURL()). 
400:           Side effects may include, in particular, actions on the composed objects (see PetscObjectCompose()).
402:   Level: beginner.
404: .seealso: PetscFwkSetURL(), PetscFwkGetURL(), PetscObjectCompose(), PetscObjectQuery()
405: @*/
406: PetscErrorCode PetscFwkCall(PetscFwk fwk, const char* message) {
411:   if(!message || !message[0]) {
412:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Null or empty message string");
413:   }
414:   switch(fwk->vtable_type) {
415:   case PETSC_FWK_VTABLE_NONE:
416:     PetscFwkCall_NONE(fwk, message);
417:     break;
418:   case PETSC_FWK_VTABLE_SO:
419:     PetscFwkCall_SO(fwk,
420:                            ((struct _n_PetscFwkVTable_SO*)fwk->vtable)->path,
421:                            ((struct _n_PetscFwkVTable_SO*)fwk->vtable)->name,
422:                            message);
423: 
424:     break;
425:   case PETSC_FWK_VTABLE_PY:
426:     PETSC_FWK_CALL_PYTHON(fwk, message);
427:     break;
428:   default:
429:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown PetscFwkVTableType value");
430:     break;
431:   }
432:   return(0);
433: }/* PetscFwkCall() */
436: #define PETSC_FWK_MAX_URL_LENGTH 1024
437: /* 
438:    Normalize the url (by truncating to PETSC_FWK_MAX_URL_LENGTH) and parse it to find out the component type and location.
439:    Warning: the returned char pointers are borrowed and their contents must be copied elsewhere to be preserved.
440: */
441: #undef  __FUNCT__
443: PetscErrorCode  PetscFwkParseURL_Private(const char inurl[], char **outpath, char **outname, PetscFwkVTableType *outtype){
444:   char *n, *s;
445:   static PetscInt nlen = PETSC_FWK_MAX_URL_LENGTH;
446:   static char path[PETSC_FWK_MAX_URL_LENGTH+1], name[PETSC_FWK_MAX_URL_LENGTH+1];
447:   PetscFwkVTableType type = PETSC_FWK_VTABLE_NONE;
450:   /* FIX: this routine should replace the filesystem path by an abolute path for real normalization */
451:   /* Copy the inurl so we can manipulate it inplace and also truncate to the max allowable length */
452:   PetscStrncpy(path, inurl, nlen);
453:   /* Split url <path>:<name> into <path> and <name> */
454:   PetscStrrchr(path,':',&n);
455:   /* Make sure it's not the ":/" of the "://" separator */
456:   if(!n[0] || n[0] == '/') {
457:     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,
458:            "Could not locate component name within the URL.\n"
459:            "Must have url = [<path/><library>:]<name>.\n"
460:            "Instead got %s\n"
461:            "Remember that URL is always truncated to the max allowed length of %d",
462:            inurl, nlen);
463:   }
464:   /* Copy n to name */
465:   PetscStrcpy(name, n);
466:   /* If n isn't the whole path (i.e., there is a ':' separator), end 'path' right before the located ':' */
467:   if(n == path) {
468:     /* 
469:        No library is provided, so the component is assumed to be "local", that is
470:        defined in an already loaded executable. So we set type to .so, path to "",
471:        and look for the configure symbol among already loaded symbols 
472:        (or count on PetscDLXXX to do that.
473:     */
474:     type = PETSC_FWK_VTABLE_SO;
475:     path[0] = '\0';
476:   }
477:   else {
478:     n[-1] = '\0';
479:     /* Find the library suffix and determine the component library type: .so or .py */
480:     PetscStrrchr(path,'.',&s);
481:     /* FIX: we should really be using PETSc's internally defined suffices */
482:     if(s != path && s[-1] == '.') {
483:       if((s[0] == 'a' && s[1] == '\0') || (s[0] == 's' && s[1] == 'o' && s[2] == '\0')){
484:         type = PETSC_FWK_VTABLE_SO;
485:       }
486:       else if (s[0] == 'p' && s[1] == 'y' && s[2] == '\0'){
487:         type = PETSC_FWK_VTABLE_PY;
488:       }
489:       else {
490:         SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,
491:                  "Unknown library suffix within the URL.\n"
492:                  "Must have url = [<path/><library>:]<name>,\n"
493:                  "where library = <libname>.<suffix>, suffix = .a || .so || .py.\n"
494:                  "Instead got url %s\n"
495:                  "Remember that URL is always truncated to the max allowed length of %d",
496:                  inurl, s,nlen);
497:       }
498:     }
499:   }
500:   *outpath = path;
501:   *outname = name;
502:   *outtype = type;
503:   return(0);
504: }/* PetscFwkParseURL_Private() */
506: #undef  __FUNCT__
508: PetscErrorCode  PetscFwkClearURL_Private(PetscFwk fwk) {
511:   switch(fwk->vtable_type) {
512:   case PETSC_FWK_VTABLE_SO:
513:     {
514:       struct _n_PetscFwkVTable_SO *vt = (struct _n_PetscFwkVTable_SO*)(fwk->vtable);
515:       PetscFree(vt->path);
516:       PetscFree(vt->name);
517:       PetscFree(vt);
518:       fwk->vtable = PETSC_NULL;
519:       fwk->vtable_type = PETSC_FWK_VTABLE_NONE;
520:     }
521:     break;
522:   case PETSC_FWK_VTABLE_PY:
523:     PETSC_FWK_CLEAR_VTABLE_PYTHON(fwk);
524:     break;
525:   case PETSC_FWK_VTABLE_NONE:
526:     break;
527:   default:
528:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,
529:              "Unknown PetscFwk vtable type: %d", fwk->vtable_type);
530:   }
531:   PetscFree(fwk->url);
532:   fwk->url = PETSC_NULL;
533:   return(0);
534: }/* PetscFwkClearURL_Private() */
536: #undef  __FUNCT__
538: /*@C 
539:    PetscFwkSetURL -- set a backend implementing PetscFwk functionality from the URL string.
541:    Logically collective on PetscFwk.
543:    Input paramters:
544: +  fwk -- a PetscFwk object
545: -  url -- URL string
547:    Notes: URL can point to a backendn -- a .so file or a .py file.  
548:      The .so file must contain symbols for function 'void call(const char[])' or symbols 'void <msg>(void)'
549:    for any message <msg> that PetscFwk is expected to understand.  When PetscFwkCall() is invoked
550:    with <msg>, a symbol for 'void <msg>(void)' is sought and called, if found.  If not, a symbol for
551:    'void call(const char[])' is sought and called with <msg> as the argument. If neither symbol is found,
552:    an error occurs.
553:      The .py file must define a class that implements 'call(str)' or '<msg>()' methods, as above.
554:      If no URL has been set, fwk attempts to reponsd to the message using function '<msg>' (failing that,
555:    'call') retrieved from fwk using PetscObjectQueryFunction().  If neither '<msg>' nor 'call' have been
556:    previusly composed with fwk (see PetscObjectComposeFunction()), an error occurs.
558:    Level: intermediate.
560: .seealso: PetscFwkGetURL(), PetscObjectCompose(), PetscObjectQuery(), PetscObjectComposeFunction()
561: @*/
562: PetscErrorCode  PetscFwkSetURL(PetscFwk fwk, const char url[]) {
564:   char *path, *name;
568:   if(fwk->vtable) {
569:     PetscFwkClearURL_Private(fwk);
570:   }
571:   PetscStrallocpy(url,  &(fwk->url));
572:   PetscFwkParseURL_Private(url, &path, &name, &fwk->vtable_type);
573:   switch(fwk->vtable_type) {
574:   case PETSC_FWK_VTABLE_SO:
575:     {
576:       struct _n_PetscFwkVTable_SO *vt;
577:       PetscMalloc(sizeof(struct _n_PetscFwkVTable_SO), &(vt));
578:       fwk->vtable = (void*)vt;
579:       PetscStrallocpy(path, &vt->path);
580:       PetscStrallocpy(name, &vt->name);
581:     }
582:     break;
583:   case PETSC_FWK_VTABLE_PY:
584:     PETSC_FWK_LOAD_VTABLE_PYTHON(fwk, path, name);
585:     break;
586:   default:
587:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,
588:              "Unknown PetscFwk vtable type: %d", fwk->vtable_type);
589:   }
590:   return(0);
591: }/* PetscFwkSetURL() */
594: #undef  __FUNCT__
596: /*@C 
597:    PetscFwkGetURL -- retrieve the URL defining the backend that implements the PetscFwkCall() functionality.
599:    Not collective.
601:    Input paramters:
602: .  fwk -- a PetscFwk object
604:    Output parameters:
605: .  url -- the URL string
607:    Level: beginner.
609: .seealso: PetscFwkSetURL(), PetscFwkCall()
610: @*/
611: PetscErrorCode  PetscFwkGetURL(PetscFwk fwk, const char **url) {
615:   *url = fwk->url;
616:   return(0);
617: }/* PetscFwkGetURL() */
618: 
620: /* ------------------------------------------------------------------------------------------------------- */
621: #undef  __FUNCT__
623: PetscErrorCode PetscFwkView(PetscFwk fwk, PetscViewer viewer) {
624:   PetscInt *vertices, N;
625:   PetscInt i, id;
626:   PetscBool         iascii;
627:   PetscErrorCode    ierr;
630:   if (!viewer) {
631:     PetscViewerASCIIGetStdout(((PetscObject)fwk)->comm,&viewer);
632:   }
635:   PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
636:   if (iascii) {
637:     PetscViewerASCIIPrintf(viewer, "Framework url: %s\n", fwk->url);
638:     PetscViewerASCIIPrintf(viewer, "Components are visited in this order:\n");
639:     PetscFwkGraphTopologicalSort(fwk->dep_graph, &N, &vertices);
640:     for(i = 0; i < N; ++i) {
641:       if(i) {
642:         PetscViewerASCIIPrintf(viewer, "\n");
643:       }
644:       id = vertices[i];
645:       PetscViewerASCIIPrintf(viewer, "%d: name:%s, url:%s", id, ((PetscObject)(fwk->component[id]))->name, fwk->component[id]->url);
646:     }
647:     PetscViewerASCIIPrintf(viewer, "\n");
648:     PetscFree(vertices);
649:   }
650:   else {
651:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG, "PetscFwk can only be viewed with an ASCII viewer");
652:   }
653:   return(0);
654: }/* PetscFwkView() */
656: #undef  __FUNCT__
658: PetscErrorCode PetscFwkGetParent(PetscFwk fwk, PetscFwk *parent)
659: {
663:   *parent = fwk->parent;
664:   return(0);
665: }/* PetscFwkGetParent() */
667: #undef  __FUNCT__
669: PetscErrorCode PetscFwkVisit(PetscFwk fwk, const char* message){
670:   PetscInt i, id, N, *vertices;
671:   PetscFwk component;
676:   PetscFwkGraphTopologicalSort(fwk->dep_graph, &N, &vertices);
677:   for(i = 0; i < N; ++i) {
678:     id = vertices[i];
679:     component = fwk->component[id];
680:     /* Save the component's visitor */
681:     component->parent = fwk;
682:     /* Call "configure" */
683:     PetscFwkCall(component, message);
684:     /* Clear visitor */
685:     component->parent = PETSC_NULL;
686:   }
687:   PetscFree(vertices);
688:   return(0);
689: }/* PetscFwkVisit() */
691: #undef  __FUNCT__
693: PetscErrorCode  PetscFwkGetKeyID_Private(PetscFwk fwk, const char key[], PetscInt *_id, PetscBool  *_found){
694:   PetscInt i;
695:   PetscBool  eq;
698:   /* Check whether a component with the given key has already been registered. */
699:   if(_found){*_found = PETSC_FALSE;}
700:   for(i = 0; i < fwk->N; ++i) {
701:     PetscStrcmp(key, ((PetscObject)fwk->component[i])->name, &eq);
702:     if(eq) {
703:       if(_id) {*_id = i;}
704:       if(_found){*_found = PETSC_TRUE;}
705:     }
706:   }
707:   return(0);
708: }/* PetscFwkGetKeyID_Private() */
710: #undef  __FUNCT__
712: PetscErrorCode  PetscFwkRegisterKey_Private(PetscFwk fwk, const char key[], PetscInt *_id) {
713:   PetscInt v, id;
714:   PetscBool  found;
717:   /* Check whether a component with the given url has already been registered.  If so, return its id, if it has been requested. */
718:   PetscFwkGetKeyID_Private(fwk, key, &id, &found);
719:   if(found) {
720:     goto alldone;
721:   }
722:   /***/
723:   /* No such key found. */
724:   /* Create a new component for this key. */
725:   if(fwk->N >= fwk->maxN) {
726:     /* No more empty component slots, therefore, expand the component array */
727:     PetscFwk *new_components;
728:     PetscMalloc(sizeof(PetscFwk)*(fwk->maxN+CHUNKSIZE), &new_components);
729:     PetscMemcpy(new_components, fwk->component, sizeof(PetscFwk)*(fwk->maxN));
730:     PetscMemzero(new_components+fwk->maxN,sizeof(PetscFwk)*(CHUNKSIZE));
731:     PetscFree(fwk->component);
732:     fwk->component = new_components;
733:     fwk->maxN += CHUNKSIZE;
734:   }
735:   id = fwk->N;
736:   ++(fwk->N);
737:   /* Store key */
738:   /* Create the corresponding component. */
739:   PetscFwkCreate(((PetscObject)fwk)->comm, fwk->component+id);
740:   PetscObjectSetName((PetscObject)fwk->component[id], key);
741:   /* Add a new vertex to the dependence graph.  This vertex will correspond to the newly registered component. */
742:   PetscFwkGraphAddVertex(fwk->dep_graph, &v);
743:   /* v must equal id */
744:   if(v != id) {
745:     SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT, "New dependence graph vertex %d for key %s not the same as component id %d", v, key, id);
746:   }
747:  alldone:
748:   if(_id) {
749:     *_id = id;
750:   }
751:   return(0);
752: }/* PetscFwkRegisterKey_Private() */
755: #undef  __FUNCT__
757: PetscErrorCode  PetscFwkRegisterComponent(PetscFwk fwk, const char key[]){
762:   PetscFwkRegisterKey_Private(fwk, key, PETSC_NULL);
763:   return(0);
764: }/* PetscFwkRegisterComponent() */
766: #undef  __FUNCT__
768: PetscErrorCode  PetscFwkRegisterComponentURL(PetscFwk fwk, const char key[], const char url[]){
770:   PetscInt id;
775:   PetscFwkRegisterKey_Private(fwk, key, &id);
776:   PetscFwkSetURL(fwk->component[id], url);
777:   return(0);
778: }/* PetscFwkRegisterComponentURL() */
781: #undef  __FUNCT__
783: PetscErrorCode  PetscFwkRegisterDependence(PetscFwk fwk, const char clientkey[], const char serverkey[])
784: {
785:   PetscInt clientid, serverid;
791:   /* Register keys */
792:   PetscFwkRegisterKey_Private(fwk, clientkey, &clientid);
793:   PetscFwkRegisterKey_Private(fwk, serverkey, &serverid);
794:   /*
795:     Add the dependency edge to the dependence_graph as follows (serverurl, clienturl): 
796:      this means "server preceeds client", so server should be configured first.
797:   */
798:   PetscFwkGraphAddEdge(fwk->dep_graph, clientid, serverid);
799:   return(0);
800: }/* PetscFwkRegisterDependence() */
804: #undef  __FUNCT__
806: /*@C 
807:    PetscFwkDestroy -- destroy PetscFwk.
809:    Not collective.
811:    Input paramters:
812: .  fwk -- a PetscFwk object
814:    Level: beginner.
816: .seealso: PetscFwkCreate(), PetscFwkSetURL(), PetscFwkCall()
817: @*/
818: PetscErrorCode  PetscFwkDestroy(PetscFwk *fwk)
819: {
820:   PetscInt       i;
822: 
824:   if (!*fwk) return(0);
826:   if (--((PetscObject)(*fwk))->refct > 0) return(0);
827:   for(i = 0; i < (*fwk)->N; ++i){
828:     PetscObjectDestroy((PetscObject*)&(*fwk)->component[i]);
829:   }
830:   PetscFree((*fwk)->component);
831:   PetscFwkGraphDestroy((*fwk)->dep_graph);
832:   PetscHeaderDestroy(fwk);
833:   return(0);
834: }/* PetscFwkDestroy()*/
836: #undef  __FUNCT__
838: /*@C 
839:    PetscFwkCreate -- create an empty PetscFwk object.
841:    Logically collective on comm.
843:    Input paramters:
844: .  comm      -- the MPI_Comm to create the PetscFwk object on.
846:    Output parameters:
847: .  framework -- the created PetscFwk object.
849:    Level: beginner.
851: .seealso: PetscFwkDestroy(), PetscFwkSetURL(), PetscFwkCall()
852: @*/
853: PetscErrorCode  PetscFwkCreate(MPI_Comm comm, PetscFwk *framework){
854:   PetscFwk fwk;
857: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
858:   PetscFwkInitializePackage(PETSC_NULL);
859: #endif
861:   PetscHeaderCreate(fwk,_p_PetscFwk,PetscInt,PETSC_FWK_CLASSID,0,"PetscFwk","Object dependence Framework","fwk",comm,PetscFwkDestroy,PetscFwkView);
862:   fwk->parent      = PETSC_NULL;
863:   fwk->component   = PETSC_NULL;
864:   fwk->vtable_type = PETSC_FWK_VTABLE_NONE;
865:   fwk->vtable      = PETSC_NULL;
866:   fwk->N = fwk->maxN = 0;
867:   /* FIX: should only create a graph on demand */
868:   PetscFwkGraphCreate(&fwk->dep_graph);
869:   PetscObjectChangeTypeName((PetscObject)fwk,PETSCFWK);
870:   *framework = fwk;
871:   return(0);
872: }/* PetscFwkCreate() */
875: #undef  __FUNCT__
877: PetscErrorCode  PetscFwkGetComponent(PetscFwk fwk, const char key[], PetscFwk *_component, PetscBool  *_found) {
878:   PetscInt id;
879:   PetscBool  found;
886:   PetscFwkGetKeyID_Private(fwk, key, &id, &found);
887:   if(found && _component) {
888:     *_component = fwk->component[id];
889:   }
890:   if(_found) {*_found = found;}
891:   return(0);
892: }/* PetscFwkGetComponent() */
896: #undef  __FUNCT__
898: PetscErrorCode PetscFwkFinalizePackage(void)
899: {
902:   if (defaultFwk) {
903:     PetscFwkDestroy(&defaultFwk);
904:   }
905:   PetscFwkPackageInitialized = PETSC_FALSE;
906:   return(0);
907: }/* PetscFwkFinalizePackage() */
910: #undef  __FUNCT__
912: PetscErrorCode PetscFwkInitializePackage(const char path[]){
915:   if(PetscFwkPackageInitialized) return(0);
916:   PetscFwkPackageInitialized = PETSC_TRUE;
917:   /* Register classes */
918:   PetscClassIdRegister(PETSC_FWK_CLASS_NAME, &PETSC_FWK_CLASSID);
919:   /* Register finalization routine */
920:   PetscRegisterFinalize(PetscFwkFinalizePackage);
921:   return(0);
922: }/* PetscFwkInitializePackage() */
924: /* ---------------------------------------------------------------------*/
925: /*
926:     The variable Petsc_Fwk_default_keyval is used to indicate an MPI attribute that
927:   is attached to a communicator, in this case the attribute is a PetscFwk.
928: */
929: static PetscMPIInt Petsc_Fwk_default_keyval = MPI_KEYVAL_INVALID;
931: #undef  __FUNCT__
933: PetscFwk  PETSC_FWK_DEFAULT_(MPI_Comm comm) {
935:   PetscBool      flg;
936:   PetscFwk       fwk;
939:   if (Petsc_Fwk_default_keyval == MPI_KEYVAL_INVALID) {
940:     MPI_Keyval_create(MPI_NULL_COPY_FN,MPI_NULL_DELETE_FN,&Petsc_Fwk_default_keyval,0);
941:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_FWK_DEFAULT_",__FILE__,__SDIR__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," "); return(0);}
942:   }
943:   MPI_Attr_get(comm,Petsc_Fwk_default_keyval,(void **)(&fwk),(PetscMPIInt*)&flg);
944:   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_FWK_DEFAULT_",__FILE__,__SDIR__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," "); return(0);}
945:   if (!flg) { /* PetscFwk not yet created */
946:     PetscFwkCreate(comm, &fwk);
947:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_FWK_DEFAULT_",__FILE__,__SDIR__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," "); return(0);}
948:     PetscObjectRegisterDestroy((PetscObject)fwk);
949:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_FWK_DEFAULT_",__FILE__,__SDIR__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," "); return(0);}
950:     MPI_Attr_put(comm,Petsc_Fwk_default_keyval,(void*)fwk);
951:     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_FWK_DEFAULT_",__FILE__,__SDIR__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," "); return(0);}
952:   }
953:   PetscFunctionReturn(fwk);
954: }/* PETSC_FWK_DEFAULT_() */