Actual source code: ex1.c


  2: static char help[] = "Used to benchmark changes to the PETSc VecScatter routines\n\n";
  3: #include <petscksp.h>
  4: extern PetscErrorCode PetscLogView_VecScatter(PetscViewer);

  6: int main(int argc, char **args)
  7: {
  8:   KSP         ksp;
  9:   Mat         A;
 10:   Vec         x, b;
 11:   PetscViewer fd;
 12:   char        file[PETSC_MAX_PATH_LEN];
 13:   PetscBool   flg, preload = PETSC_TRUE;

 15:   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
 16:   PetscCall(PetscLogDefaultBegin());
 17:   PetscCall(PetscOptionsGetString(NULL, NULL, "-f", file, sizeof(file), &flg));

 19:   PetscPreLoadBegin(preload, "Load system");

 21:   /*
 22:      Load the matrix and vector; then destroy the viewer.
 23:   */
 24:   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
 25:   PetscCall(MatSetFromOptions(A));
 26:   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &fd));
 27:   PetscCall(MatLoad(A, fd));
 28:   PetscCall(PetscViewerDestroy(&fd));

 30:   PetscCall(MatCreateVecs(A, &x, &b));
 31:   PetscCall(VecSetFromOptions(b));
 32:   PetscCall(VecSet(b, 1.0));

 34:   PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
 35:   PetscCall(KSPSetFromOptions(ksp));
 36:   PetscCall(KSPSetOperators(ksp, A, A));
 37:   PetscCall(KSPSetUp(ksp));
 38:   PetscCall(KSPSetUpOnBlocks(ksp));

 40:   PetscPreLoadStage("KSPSolve");
 41:   PetscCall(KSPSolve(ksp, b, x));

 43:   PetscCall(MatDestroy(&A));
 44:   PetscCall(VecDestroy(&b));
 45:   PetscCall(VecDestroy(&x));
 46:   PetscCall(KSPDestroy(&ksp));
 47:   PetscPreLoadEnd();
 48:   PetscCall(PetscLogView_VecScatter(PETSC_VIEWER_STDOUT_WORLD));

 50:   PetscCall(PetscFinalize());
 51:   return 0;
 52: }

 54: #include <petsctime.h>
 55: #include <petsc/private/petscimpl.h>
 56: #include <petsc/private/vecimpl.h>
 57: #include <petsc/private/kspimpl.h>
 58: #include <petscmachineinfo.h>
 59: #include <petscconfiginfo.h>
 60: /*
 61:    This is a special log viewer that prints out detailed information only for the VecScatter routines
 62: */
 63: typedef enum {
 64:   COUNT,
 65:   TIME,
 66:   NUMMESS,
 67:   MESSLEN,
 68:   REDUCT,
 69:   FLOPS
 70: } Stats;
 71: PetscErrorCode PetscLogView_VecScatter(PetscViewer viewer)
 72: {
 73:   MPI_Comm            comm      = PetscObjectComm((PetscObject)viewer);
 74:   PetscEventPerfInfo *eventInfo = NULL;
 75:   PetscLogDouble      locTotalTime, stats[6], maxstats[6], minstats[6], sumstats[6], avetime, ksptime;
 76:   PetscStageLog       stageLog;
 77:   const int           stage = 2;
 78:   int                 event, events[] = {VEC_ScatterBegin, VEC_ScatterEnd};
 79:   PetscMPIInt         rank, size;
 80:   PetscInt            i;
 81:   char                arch[128], hostname[128], username[128], pname[PETSC_MAX_PATH_LEN], date[128], version[256];

 83:   PetscFunctionBegin;
 84:   PetscCall(PetscTime(&locTotalTime));
 85:   locTotalTime -= petsc_BaseTime;
 86:   PetscCallMPI(MPI_Comm_size(comm, &size));
 87:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
 88:   PetscCall(PetscLogGetStageLog(&stageLog));
 89:   PetscCall(PetscViewerASCIIPrintf(viewer, "numProcs   = %d\n", size));

 91:   PetscCall(PetscGetArchType(arch, sizeof(arch)));
 92:   PetscCall(PetscGetHostName(hostname, sizeof(hostname)));
 93:   PetscCall(PetscGetUserName(username, sizeof(username)));
 94:   PetscCall(PetscGetProgramName(pname, sizeof(pname)));
 95:   PetscCall(PetscGetDate(date, sizeof(date)));
 96:   PetscCall(PetscGetVersion(version, sizeof(version)));
 97:   PetscCall(PetscViewerASCIIPrintf(viewer, "%s on a %s named %s with %d processors, by %s %s\n", pname, arch, hostname, size, username, date));
 98:   PetscCall(PetscViewerASCIIPrintf(viewer, "Using %s\n", version));
 99:   PetscCall(PetscViewerASCIIPrintf(viewer, "Configure options: %s", petscconfigureoptions));
100:   PetscCall(PetscViewerASCIIPrintf(viewer, "%s", petscmachineinfo));
101:   PetscCall(PetscViewerASCIIPrintf(viewer, "%s", petsccompilerinfo));
102:   PetscCall(PetscViewerASCIIPrintf(viewer, "%s", petsccompilerflagsinfo));
103:   PetscCall(PetscViewerASCIIPrintf(viewer, "%s", petsclinkerinfo));
104:   PetscCall(PetscViewerASCIIPrintf(viewer, "%s\n", PETSC_MPICC_SHOW));
105:   PetscCall(PetscOptionsView(NULL, viewer));
106: #if defined(PETSC_HAVE_HWLOC)
107:   PetscCall(PetscProcessPlacementView(viewer));
108: #endif
109:   PetscCall(PetscViewerASCIIPrintf(viewer, "----------------------------------------------------\n"));

111:   PetscCall(PetscViewerASCIIPrintf(viewer, "                Time     Min to Max Range   Proportion of KSP\n"));

113:   eventInfo = stageLog->stageInfo[stage].eventLog->eventInfo;
114:   PetscCall(MPIU_Allreduce(&eventInfo[KSP_Solve].time, &ksptime, 1, MPIU_PETSCLOGDOUBLE, MPI_SUM, PETSC_COMM_WORLD));
115:   ksptime = ksptime / size;

117:   for (i = 0; i < (int)(sizeof(events) / sizeof(int)); i++) {
118:     event          = events[i];
119:     stats[COUNT]   = eventInfo[event].count;
120:     stats[TIME]    = eventInfo[event].time;
121:     stats[NUMMESS] = eventInfo[event].numMessages;
122:     stats[MESSLEN] = eventInfo[event].messageLength;
123:     stats[REDUCT]  = eventInfo[event].numReductions;
124:     stats[FLOPS]   = eventInfo[event].flops;
125:     PetscCall(MPIU_Allreduce(stats, maxstats, 6, MPIU_PETSCLOGDOUBLE, MPI_MAX, PETSC_COMM_WORLD));
126:     PetscCall(MPIU_Allreduce(stats, minstats, 6, MPIU_PETSCLOGDOUBLE, MPI_MIN, PETSC_COMM_WORLD));
127:     PetscCall(MPIU_Allreduce(stats, sumstats, 6, MPIU_PETSCLOGDOUBLE, MPI_SUM, PETSC_COMM_WORLD));

129:     avetime = sumstats[1] / size;
130:     PetscCall(PetscViewerASCIIPrintf(viewer, "%s %4.2e   -%5.1f %% %5.1f %%   %4.2e %%\n", stageLog->eventLog->eventInfo[event].name, avetime, 100. * (avetime - minstats[1]) / avetime, 100. * (maxstats[1] - avetime) / avetime, 100. * avetime / ksptime));
131:   }
132:   PetscCall(PetscViewerFlush(viewer));
133:   PetscFunctionReturn(PETSC_SUCCESS);
134: }

136: /*TEST

138:    build:
139:      requires: defined(PETSC_USE_LOG)

141:    test:
142:      TODO: need to implement

144: TEST*/