1e669de00SBarry Smith #include <petsc/private/dmnetworkimpl.h> /*I "petscdmplex.h" I*/ 2e669de00SBarry Smith #include <petscts.h> 3e669de00SBarry Smith #include <petscdraw.h> 4e669de00SBarry Smith 5e669de00SBarry Smith /* 64500feddSGetnet Betrie TSMonitorLGCtxDestroy - Destroys line graph contexts that where created with TSMonitorLGCtxNetworkCreate(). 7e669de00SBarry Smith 8e669de00SBarry Smith Collective on TSMonitorLGCtx_Network 9e669de00SBarry Smith 10e669de00SBarry Smith Input Parameter: 11e669de00SBarry Smith . ctx - the monitor context 12e669de00SBarry Smith 13e669de00SBarry Smith */ 14e669de00SBarry Smith PetscErrorCode TSMonitorLGCtxNetworkDestroy(TSMonitorLGCtxNetwork *ctx) 15e669de00SBarry Smith { 16e669de00SBarry Smith PetscInt i; 17e669de00SBarry Smith 18e669de00SBarry Smith PetscFunctionBegin; 19e669de00SBarry Smith for (i=0; i<(*ctx)->nlg; i++) { 20*9566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDestroy(&(*ctx)->lg[i])); 21e669de00SBarry Smith } 22*9566063dSJacob Faibussowitsch PetscCall(PetscFree((*ctx)->lg)); 23*9566063dSJacob Faibussowitsch PetscCall(PetscFree(*ctx)); 24e669de00SBarry Smith PetscFunctionReturn(0); 25e669de00SBarry Smith } 26e669de00SBarry Smith 27e669de00SBarry Smith PetscErrorCode TSMonitorLGCtxNetworkCreate(TS ts,const char host[],const char label[],int x,int y,int m,int n,PetscInt howoften,TSMonitorLGCtxNetwork *ctx) 28e669de00SBarry Smith { 29e669de00SBarry Smith PetscDraw draw; 30e669de00SBarry Smith MPI_Comm comm; 31e669de00SBarry Smith DM dm; 32e669de00SBarry Smith PetscInt i,Start,End,e,nvar; 33e669de00SBarry Smith 34e669de00SBarry Smith PetscFunctionBegin; 35*9566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&dm)); 36*9566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)ts,&comm)); 37*9566063dSJacob Faibussowitsch PetscCall(PetscNew(ctx)); 38e669de00SBarry Smith i = 0; 39e669de00SBarry Smith /* loop over edges counting number of line graphs needed */ 40*9566063dSJacob Faibussowitsch PetscCall(DMNetworkGetEdgeRange(dm,&Start,&End)); 41e669de00SBarry Smith for (e=Start; e<End; e++) { 42*9566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(dm,e,ALL_COMPONENTS,NULL,NULL,&nvar)); 43e669de00SBarry Smith if (!nvar) continue; 44e669de00SBarry Smith i++; 45e669de00SBarry Smith } 46e669de00SBarry Smith /* loop over vertices */ 47*9566063dSJacob Faibussowitsch PetscCall(DMNetworkGetVertexRange(dm,&Start,&End)); 48e669de00SBarry Smith for (e=Start; e<End; e++) { 49*9566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(dm,e,ALL_COMPONENTS,NULL,NULL,&nvar)); 50e669de00SBarry Smith if (!nvar) continue; 51e669de00SBarry Smith i++; 52e669de00SBarry Smith } 53e669de00SBarry Smith (*ctx)->nlg = i; 54*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(i,&(*ctx)->lg)); 55e669de00SBarry Smith 56e669de00SBarry Smith i = 0; 57e669de00SBarry Smith /* loop over edges creating all needed line graphs*/ 58*9566063dSJacob Faibussowitsch PetscCall(DMNetworkGetEdgeRange(dm,&Start,&End)); 59e669de00SBarry Smith for (e=Start; e<End; e++) { 60*9566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(dm,e,ALL_COMPONENTS,NULL,NULL,&nvar)); 61e669de00SBarry Smith if (!nvar) continue; 62*9566063dSJacob Faibussowitsch PetscCall(PetscDrawCreate(comm,host,label,x,y,m,n,&draw)); 63*9566063dSJacob Faibussowitsch PetscCall(PetscDrawSetFromOptions(draw)); 64*9566063dSJacob Faibussowitsch PetscCall(PetscDrawLGCreate(draw,nvar,&(*ctx)->lg[i])); 65*9566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetFromOptions((*ctx)->lg[i])); 66*9566063dSJacob Faibussowitsch PetscCall(PetscDrawDestroy(&draw)); 67e669de00SBarry Smith i++; 68e669de00SBarry Smith } 69e669de00SBarry Smith /* loop over vertices */ 70*9566063dSJacob Faibussowitsch PetscCall(DMNetworkGetVertexRange(dm,&Start,&End)); 71e669de00SBarry Smith for (e=Start; e<End; e++) { 72*9566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(dm,e,ALL_COMPONENTS,NULL,NULL,&nvar)); 73e669de00SBarry Smith if (!nvar) continue; 74*9566063dSJacob Faibussowitsch PetscCall(PetscDrawCreate(comm,host,label,x,y,m,n,&draw)); 75*9566063dSJacob Faibussowitsch PetscCall(PetscDrawSetFromOptions(draw)); 76*9566063dSJacob Faibussowitsch PetscCall(PetscDrawLGCreate(draw,nvar,&(*ctx)->lg[i])); 77*9566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetFromOptions((*ctx)->lg[i])); 78*9566063dSJacob Faibussowitsch PetscCall(PetscDrawDestroy(&draw)); 79e669de00SBarry Smith i++; 80e669de00SBarry Smith } 81*9566063dSJacob Faibussowitsch PetscCall(PetscDrawDestroy(&draw)); 82e669de00SBarry Smith (*ctx)->howoften = howoften; 83e669de00SBarry Smith PetscFunctionReturn(0); 84e669de00SBarry Smith } 85e669de00SBarry Smith 86e669de00SBarry Smith /* 87e669de00SBarry Smith TSMonitorLGCtxNetworkSolution - Monitors progress of the TS solvers for a DMNetwork solution with one window for each vertex and each edge 88e669de00SBarry Smith 89e669de00SBarry Smith Collective on TS 90e669de00SBarry Smith 91e669de00SBarry Smith Input Parameters: 92e669de00SBarry Smith + ts - the TS context 93e669de00SBarry Smith . step - current time-step 94e669de00SBarry Smith . ptime - current time 95e669de00SBarry Smith . u - current solution 964500feddSGetnet Betrie - dctx - the TSMonitorLGCtxNetwork object that contains all the options for the monitoring, this is created with TSMonitorLGCtxCreateNetwork() 97e669de00SBarry Smith 98e669de00SBarry Smith Options Database: 99e669de00SBarry Smith . -ts_monitor_lg_solution_variables 100e669de00SBarry Smith 101e669de00SBarry Smith Level: intermediate 102e669de00SBarry Smith 103e669de00SBarry Smith Notes: 104e669de00SBarry Smith Each process in a parallel run displays its component solutions in a separate window 105e669de00SBarry Smith 106e669de00SBarry Smith */ 107e669de00SBarry Smith PetscErrorCode TSMonitorLGCtxNetworkSolution(TS ts,PetscInt step,PetscReal ptime,Vec u,void *dctx) 108e669de00SBarry Smith { 109e669de00SBarry Smith TSMonitorLGCtxNetwork ctx = (TSMonitorLGCtxNetwork)dctx; 110e669de00SBarry Smith const PetscScalar *xv; 111e669de00SBarry Smith PetscScalar *yv; 112e669de00SBarry Smith PetscInt i,v,Start,End,offset,nvar,e; 113e669de00SBarry Smith TSConvergedReason reason; 114e669de00SBarry Smith DM dm; 115e669de00SBarry Smith Vec uv; 116e669de00SBarry Smith 117e669de00SBarry Smith PetscFunctionBegin; 118e669de00SBarry Smith if (step < 0) PetscFunctionReturn(0); /* -1 indicates interpolated solution */ 119e669de00SBarry Smith if (!step) { 120e669de00SBarry Smith PetscDrawAxis axis; 121e669de00SBarry Smith 122e669de00SBarry Smith for (i=0; i<ctx->nlg; i++) { 123*9566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetAxis(ctx->lg[i],&axis)); 124*9566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis,"Solution as function of time","Time","Solution")); 125*9566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(ctx->lg[i])); 126e669de00SBarry Smith } 127e669de00SBarry Smith } 128e669de00SBarry Smith 129e669de00SBarry Smith if (ctx->semilogy) { 130e669de00SBarry Smith PetscInt n,j; 131e669de00SBarry Smith 132*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u,&uv)); 133*9566063dSJacob Faibussowitsch PetscCall(VecCopy(u,uv)); 134*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(uv,&yv)); 135*9566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(uv,&n)); 136e669de00SBarry Smith for (j=0; j<n; j++) { 1374500feddSGetnet Betrie if (PetscRealPart(yv[j]) <= 0.0) yv[j] = -12; 1384e936a74SSatish Balay else yv[j] = PetscLog10Real(PetscRealPart(yv[j])); 139e669de00SBarry Smith } 140e669de00SBarry Smith xv = yv; 141e669de00SBarry Smith } else { 142*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(u,&xv)); 143e669de00SBarry Smith } 144e669de00SBarry Smith /* iterate over edges */ 145*9566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&dm)); 146e669de00SBarry Smith i = 0; 147*9566063dSJacob Faibussowitsch PetscCall(DMNetworkGetEdgeRange(dm,&Start,&End)); 148e669de00SBarry Smith for (e=Start; e<End; e++) { 149*9566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(dm,e,ALL_COMPONENTS,NULL,NULL,&nvar)); 150e669de00SBarry Smith if (!nvar) continue; 151e669de00SBarry Smith 152*9566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(dm,e,ALL_COMPONENTS,&offset)); 153*9566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddCommonPoint(ctx->lg[i],ptime,(const PetscReal*)(xv+offset))); 154e669de00SBarry Smith i++; 155e669de00SBarry Smith } 156e669de00SBarry Smith 157e669de00SBarry Smith /* iterate over vertices */ 158*9566063dSJacob Faibussowitsch PetscCall(DMNetworkGetVertexRange(dm,&Start,&End)); 159e669de00SBarry Smith for (v=Start; v<End; v++) { 160*9566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(dm,v,ALL_COMPONENTS,NULL,NULL,&nvar)); 161e669de00SBarry Smith if (!nvar) continue; 162e669de00SBarry Smith 163*9566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(dm,v,ALL_COMPONENTS,&offset)); 164*9566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddCommonPoint(ctx->lg[i],ptime,(const PetscReal*)(xv+offset))); 165e669de00SBarry Smith i++; 166e669de00SBarry Smith } 167e669de00SBarry Smith if (ctx->semilogy) { 168*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(uv,&yv)); 169*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&uv)); 170e669de00SBarry Smith } else { 171*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(u,&xv)); 172e669de00SBarry Smith } 173e669de00SBarry Smith 174*9566063dSJacob Faibussowitsch PetscCall(TSGetConvergedReason(ts,&reason)); 175e669de00SBarry Smith if (((ctx->howoften > 0) && (!(step % ctx->howoften))) || ((ctx->howoften == -1) && reason)) { 176e669de00SBarry Smith for (i=0; i<ctx->nlg; i++) { 177*9566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(ctx->lg[i])); 178*9566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSave(ctx->lg[i])); 179e669de00SBarry Smith } 180e669de00SBarry Smith } 181e669de00SBarry Smith PetscFunctionReturn(0); 182e669de00SBarry Smith } 183