HyPar  1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
PetscPreTimeStep.cpp File Reference

Pre-time-step function. More...

#include <stdio.h>
#include <basic.h>
#include <arrayfunctions.h>
#include <mpivars_cpp.h>
#include <simulation_object.h>
#include <petscinterface.h>
#include <librom_interface.h>

Go to the source code of this file.

Macros

#define __FUNCT__   "PetscPreTimeStep"
 

Functions

int OutputSolution (void *, int, double)
 
PetscErrorCode PetscPreTimeStep (TS ts)
 

Detailed Description

Pre-time-step function.

Author
Debojyoti Ghosh

Definition in file PetscPreTimeStep.cpp.

Macro Definition Documentation

#define __FUNCT__   "PetscPreTimeStep"

Definition at line 21 of file PetscPreTimeStep.cpp.

Function Documentation

int OutputSolution ( void *  s,
int  nsims,
double  a_time 
)

Write solutions to file

Write out the solution to file

Parameters
sArray of simulation objects of type SimulationObject
nsimsNumber of simulation objects
a_timeCurrent simulation time

Definition at line 27 of file OutputSolution.cpp.

30 {
31  SimulationObject* simobj = (SimulationObject*) s;
32  int ns;
34 
35  for (ns = 0; ns < nsims; ns++) {
36 
37  HyPar* solver = &(simobj[ns].solver);
38  MPIVariables* mpi = &(simobj[ns].mpi);
39 
40  if ((!solver->WriteOutput) && (strcmp(solver->plot_solution,"yes"))) continue;
41 
42  /* time integration module may have auxiliary arrays to write out, so get them */
43  int NSolutions = 0;
44  IERR TimeGetAuxSolutions(&NSolutions,NULL,solver,-1,ns); CHECKERR(ierr);
45  if (NSolutions > 10) NSolutions = 10;
46 
47  int nu;
48  char fname_root[_MAX_STRING_SIZE_];
49  char aux_fname_root[_MAX_STRING_SIZE_];
50  strcpy(fname_root, solver->op_fname_root);
51  strcpy(aux_fname_root, solver->aux_op_fname_root);
52 
53  if (nsims > 1) {
54  char index[_MAX_STRING_SIZE_];
55  GetStringFromInteger(ns, index, (int)log10(nsims)+1);
56  strcat(fname_root, "_");
57  strcat(fname_root, index);
58  strcat(aux_fname_root, "_");
59  strcat(aux_fname_root, index);
60  }
61 
62  for (nu=0; nu<NSolutions; nu++) {
63 
64  double *uaux = NULL;
65  IERR TimeGetAuxSolutions(&NSolutions,&uaux,solver,nu,ns); CHECKERR(ierr);
66 
67  IERR WriteArray( solver->ndims,
68  solver->nvars,
69  solver->dim_global,
70  solver->dim_local,
71  solver->ghosts,
72  solver->x,
73  uaux,
74  solver,
75  mpi,
76  aux_fname_root ); CHECKERR(ierr);
77 
78  aux_fname_root[2]++;
79  }
80 
81 #if defined(HAVE_CUDA)
82  if (solver->use_gpu) {
83  /* Copy values from GPU memory to CPU memory for writing. */
84  gpuMemcpy(solver->x, solver->gpu_x, sizeof(double)*solver->size_x, gpuMemcpyDeviceToHost);
85 
86 #ifdef CUDA_VAR_ORDERDING_AOS
87  gpuMemcpy( solver->u,
88  solver->gpu_u,
89  sizeof(double)*solver->ndof_cells_wghosts,
91 #else
92  double *h_u = (double *) malloc(sizeof(double)*solver->ndof_cells_wghosts);
93  gpuMemcpy(h_u, solver->gpu_u, sizeof(double)*solver->ndof_cells_wghosts, gpuMemcpyDeviceToHost);
94  for (int i=0; i<solver->npoints_local_wghosts; i++) {
95  for (int v=0; v<solver->nvars; v++) {
96  solver->u[i*solver->nvars+v] = h_u[i+v*solver->npoints_local_wghosts];
97  }
98  }
99  free(h_u);
100 #endif
101  }
102 #endif
103 
104  WriteArray( solver->ndims,
105  solver->nvars,
106  solver->dim_global,
107  solver->dim_local,
108  solver->ghosts,
109  solver->x,
110  solver->u,
111  solver,
112  mpi,
113  fname_root );
114 
115  if (!strcmp(solver->plot_solution, "yes")) {
116  PlotArray( solver->ndims,
117  solver->nvars,
118  solver->dim_global,
119  solver->dim_local,
120  solver->ghosts,
121  solver->x,
122  solver->u,
123  a_time,
124  solver,
125  mpi,
126  fname_root );
127  }
128 
129  /* increment the index string, if required */
130  if ((!strcmp(solver->output_mode,"serial")) && (!strcmp(solver->op_overwrite,"no"))) {
132  }
133 
134  }
135 
136  return(0);
137 }
int npoints_local_wghosts
Definition: hypar.h:42
int * dim_global
Definition: hypar.h:33
double * u
Definition: hypar.h:116
char aux_op_fname_root[_MAX_STRING_SIZE_]
Definition: hypar.h:208
void GetStringFromInteger(int, char *, int)
char plot_solution[_MAX_STRING_SIZE_]
Definition: hypar.h:194
int * dim_local
Definition: hypar.h:37
void gpuMemcpy(void *, const void *, size_t, enum gpuMemcpyKind)
int TimeGetAuxSolutions(int *, double **, void *, int, int)
int ghosts
Definition: hypar.h:52
char * filename_index
Definition: hypar.h:197
char op_fname_root[_MAX_STRING_SIZE_]
Definition: hypar.h:206
int size_x
Definition: hypar.h:48
int WriteArray(int, int, int *, int *, int, double *, double *, void *, void *, char *)
Definition: WriteArray.c:27
#define _MAX_STRING_SIZE_
Definition: basic.h:14
double * gpu_u
Definition: hypar.h:459
char op_overwrite[_MAX_STRING_SIZE_]
Definition: hypar.h:191
int ndof_cells_wghosts
Definition: hypar.h:47
double * gpu_x
Definition: hypar.h:457
int nvars
Definition: hypar.h:29
#define CHECKERR(ierr)
Definition: basic.h:18
int PlotArray(int, int, int *, int *, int, double *, double *, double, void *, void *, char *)
Definition: PlotArray.c:31
int use_gpu
Definition: hypar.h:449
Structure defining a simulation.
int ndims
Definition: hypar.h:26
#define IERR
Definition: basic.h:16
int index_length
Definition: hypar.h:199
Structure of MPI-related variables.
double * x
Definition: hypar.h:107
char output_mode[_MAX_STRING_SIZE_]
Definition: hypar.h:183
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
void IncrementFilenameIndex(char *f, int len)
int(* WriteOutput)(int, int, int *, double *, double *, char *, int *)
Definition: hypar.h:211
#define _DECLARE_IERR_
Definition: basic.h:17
PetscErrorCode PetscPreTimeStep ( TS  ts)

Function called before a time step

Parameters
tsTime integration object

Definition at line 24 of file PetscPreTimeStep.cpp.

25 {
26  PETScContext* context(nullptr);
27 
28  PetscFunctionBegin;
29 
30  TSGetApplicationContext(ts,&context);
31  if (!context) {
32  fprintf(stderr,"Error in PetscPreTimeStep: Null context!\n");
33  return(1);
34  }
35  gettimeofday(&(context->iter_start_time),NULL);
36  SimulationObject* sim = (SimulationObject*) context->simobj;
37  int nsims = context->nsims;
38 
39  Vec Y;
40  TSGetSolution(ts,&Y);
41 
42  double waqt;
43  TSGetTime(ts,&waqt);
44  double dt;
45  TSGetTimeStep(ts,&dt);
46  int iter;
47  TSGetStepNumber(ts,&iter);
48 
49  context->dt = dt;
50  context->waqt = waqt;
51  context->t_start = waqt;
52 
53  TSType time_scheme;
54  TSGetType(ts,&time_scheme);
55 
56  for (int ns = 0; ns < nsims; ns++) {
57 
58  HyPar* solver = &(sim[ns].solver);
59  MPIVariables* mpi = &(sim[ns].mpi);
60 
61  /* get solution */
62  TransferVecFromPETSc(solver->u,Y,context,ns,context->offsets[ns]);
63 
64  /* save a copy of the solution to compute norm at end of time step */
65  _ArrayCopy1D_(solver->u,solver->u0,(solver->npoints_local_wghosts*solver->nvars));
66 
67  /* apply boundary conditions and exchange data over MPI interfaces */
68  solver->ApplyBoundaryConditions(solver,mpi,solver->u,NULL,waqt);
69  solver->ApplyIBConditions(solver,mpi,solver->u,waqt);
71  solver->nvars,
72  solver->dim_local,
73  solver->ghosts,
74  mpi,
75  solver->u );
76 
77  /* Call any physics-specific pre-step function */
78  if (solver->PreStep) solver->PreStep(solver->u,solver,mpi,waqt);
79 
80  /* If using a non-linear scheme with ARKIMEX methods,
81  compute the non-linear finite-difference operator */
82  if (!strcmp(time_scheme,TSARKIMEX)) {
83  solver->NonlinearInterp(solver->u,solver,mpi,waqt,solver->FFunction);
84  }
85 
86  /* set the step boundary flux integral value to zero */
87  _ArraySetValue_(solver->StepBoundaryIntegral,2*solver->ndims*solver->nvars,0.0);
88 
89  TransferVecToPETSc(solver->u,Y,context,ns,context->offsets[ns]);
90 
91  }
92 
93 
94  if (!iter) {
95  for (int ns = 0; ns < nsims; ns++) {
96  HyPar* solver = &(sim[ns].solver);
97  MPIVariables* mpi = &(sim[ns].mpi);
98  if (solver->PhysicsOutput) solver->PhysicsOutput(solver,mpi,waqt);
99  }
100  OutputSolution(sim, nsims,waqt);
101 #ifdef with_librom
102  context->op_times_arr.push_back(waqt);
103 #endif
104  }
105 
106 #ifdef with_librom
107  if ( (context->rom_mode == _ROM_MODE_TRAIN_)
108  && (iter%((libROMInterface*)context->rom_interface)->samplingFrequency() == 0) ) {
109  ((libROMInterface*)context->rom_interface)->takeSample( sim, waqt );
110  }
111 #endif
112 
113  PetscFunctionReturn(0);
114 }
int npoints_local_wghosts
Definition: hypar.h:42
#define _ArraySetValue_(x, size, value)
double * u
Definition: hypar.h:116
int TransferVecFromPETSc(double *const, const Vec, void *, const int, const int)
int * dim_local
Definition: hypar.h:37
double * u0
Definition: hypar.h:396
int(* PhysicsOutput)(void *, void *, double)
Definition: hypar.h:347
int(* FFunction)(double *, double *, int, void *, double)
Definition: hypar.h:276
int ghosts
Definition: hypar.h:52
int(* NonlinearInterp)(double *, void *, void *, double, int(*)(double *, double *, int, void *, double))
Definition: hypar.h:228
double * StepBoundaryIntegral
Definition: hypar.h:384
Class implementing interface with libROM.
int TransferVecToPETSc(const double *const, Vec, void *, const int, const int)
#define _ArrayCopy1D_(x, y, size)
int nvars
Definition: hypar.h:29
int OutputSolution(void *, int, double)
Structure containing the variables for time-integration with PETSc.
Structure defining a simulation.
int ndims
Definition: hypar.h:26
int(* ApplyBoundaryConditions)(void *, void *, double *, double *, double)
Definition: hypar.h:214
int MPIExchangeBoundariesnD(int, int, int *, int, void *, double *)
int(* PreStep)(double *, void *, void *, double)
Definition: hypar.h:339
Structure of MPI-related variables.
#define _ROM_MODE_TRAIN_
int(* ApplyIBConditions)(void *, void *, double *, double)
Definition: hypar.h:217
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23