#include /*I "petscdmplex.h" I*/ #include #include /* TSMonitorLGCtxDestroy - Destroys line graph contexts that where created with TSMonitorLGCtxNetworkCreate(). Collective on TSMonitorLGCtx_Network Input Parameter: . ctx - the monitor context */ PetscErrorCode TSMonitorLGCtxNetworkDestroy(TSMonitorLGCtxNetwork *ctx) { PetscInt i; PetscFunctionBegin; for (i=0; i<(*ctx)->nlg; i++) { PetscCall(PetscDrawLGDestroy(&(*ctx)->lg[i])); } PetscCall(PetscFree((*ctx)->lg)); PetscCall(PetscFree(*ctx)); PetscFunctionReturn(0); } PetscErrorCode TSMonitorLGCtxNetworkCreate(TS ts,const char host[],const char label[],int x,int y,int m,int n,PetscInt howoften,TSMonitorLGCtxNetwork *ctx) { PetscDraw draw; MPI_Comm comm; DM dm; PetscInt i,Start,End,e,nvar; PetscFunctionBegin; PetscCall(TSGetDM(ts,&dm)); PetscCall(PetscObjectGetComm((PetscObject)ts,&comm)); PetscCall(PetscNew(ctx)); i = 0; /* loop over edges counting number of line graphs needed */ PetscCall(DMNetworkGetEdgeRange(dm,&Start,&End)); for (e=Start; enlg = i; PetscCall(PetscMalloc1(i,&(*ctx)->lg)); i = 0; /* loop over edges creating all needed line graphs*/ PetscCall(DMNetworkGetEdgeRange(dm,&Start,&End)); for (e=Start; elg[i])); PetscCall(PetscDrawLGSetFromOptions((*ctx)->lg[i])); PetscCall(PetscDrawDestroy(&draw)); i++; } /* loop over vertices */ PetscCall(DMNetworkGetVertexRange(dm,&Start,&End)); for (e=Start; elg[i])); PetscCall(PetscDrawLGSetFromOptions((*ctx)->lg[i])); PetscCall(PetscDrawDestroy(&draw)); i++; } PetscCall(PetscDrawDestroy(&draw)); (*ctx)->howoften = howoften; PetscFunctionReturn(0); } /* TSMonitorLGCtxNetworkSolution - Monitors progress of the TS solvers for a DMNetwork solution with one window for each vertex and each edge Collective on TS Input Parameters: + ts - the TS context . step - current time-step . ptime - current time . u - current solution - dctx - the TSMonitorLGCtxNetwork object that contains all the options for the monitoring, this is created with TSMonitorLGCtxCreateNetwork() Options Database: . -ts_monitor_lg_solution_variables Level: intermediate Notes: Each process in a parallel run displays its component solutions in a separate window */ PetscErrorCode TSMonitorLGCtxNetworkSolution(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dctx) { TSMonitorLGCtxNetwork ctx = (TSMonitorLGCtxNetwork)dctx; const PetscScalar *xv; PetscScalar *yv; PetscInt i,v,Start,End,offset,nvar,e; TSConvergedReason reason; DM dm; Vec uv; PetscFunctionBegin; if (step < 0) PetscFunctionReturn(0); /* -1 indicates interpolated solution */ if (!step) { PetscDrawAxis axis; for (i=0; inlg; i++) { PetscCall(PetscDrawLGGetAxis(ctx->lg[i],&axis)); PetscCall(PetscDrawAxisSetLabels(axis,"Solution as function of time","Time","Solution")); PetscCall(PetscDrawLGReset(ctx->lg[i])); } } if (ctx->semilogy) { PetscInt n,j; PetscCall(VecDuplicate(u,&uv)); PetscCall(VecCopy(u,uv)); PetscCall(VecGetArray(uv,&yv)); PetscCall(VecGetLocalSize(uv,&n)); for (j=0; jlg[i],ptime,(const PetscReal*)(xv+offset))); i++; } /* iterate over vertices */ PetscCall(DMNetworkGetVertexRange(dm,&Start,&End)); for (v=Start; vlg[i],ptime,(const PetscReal*)(xv+offset))); i++; } if (ctx->semilogy) { PetscCall(VecRestoreArray(uv,&yv)); PetscCall(VecDestroy(&uv)); } else { PetscCall(VecRestoreArrayRead(u,&xv)); } PetscCall(TSGetConvergedReason(ts,&reason)); if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && reason)) { for (i=0; inlg; i++) { PetscCall(PetscDrawLGDraw(ctx->lg[i])); PetscCall(PetscDrawLGSave(ctx->lg[i])); } } PetscFunctionReturn(0); }