HyPar  1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
PetscSetInitialGuessROM.cpp
Go to the documentation of this file.
1 
6 #ifdef with_petsc
7 
8 #include <simulation_object.h>
9 #include <petscinterface.h>
10 #ifdef with_librom
11 #include <librom_interface.h>
12 #endif
13 
14 #ifdef with_librom
15 
18 PetscErrorCode PetscSetInitialGuessROM( SNES snes,
19  Vec X,
20  void* ctxt )
21 {
22  PETScContext* context = (PETScContext*) ctxt;
23  SimulationObject* sim = (SimulationObject*) context->simobj;
24  int nsims = context->nsims;
25 
26  PetscFunctionBegin;
27 
28  double stage_time = context->t_start;
29  if (context->stage_times.size() > (context->stage_index+1)) {
30  stage_time += context->stage_times[context->stage_index+1] * context->dt;
31  } else {
32  stage_time += context->dt;
33  }
34 
35  ((libROMInterface*)context->rom_interface)->predict(sim, stage_time);
36  if (!context->rank) {
37  printf( "libROM: Predicted ROM intial guess (time %1.4e), wallclock time: %f.\n",
38  stage_time, ((libROMInterface*)context->rom_interface)->predictWallclockTime() );
39  }
40 
41  for (int ns = 0; ns < nsims; ns++) {
42  TransferVecToPETSc(sim[ns].solver.u,X,context,ns,context->offsets[ns]);
43  }
44 
45 
46  PetscFunctionReturn(0);
47 }
48 
49 #endif
50 
51 #endif
int TransferVecToPETSc(const double *const, Vec, void *, const int, const int)
Structure defining a simulation.
std::vector< double > stage_times
Simulation object.
Structure containing the variables for time-integration with PETSc.
libROM interface class
Class implementing interface with libROM.
PetscErrorCode PetscSetInitialGuessROM(SNES snes, Vec X, void *ctxt)