HyPar  1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
timeintegration_cpp.h File Reference

Contains C++ function declarations for time integration. More...

Go to the source code of this file.

Functions

int TimeInitialize (void *, int, int, int, void *)
 
int TimeCleanup (void *)
 
int TimePreStep (void *)
 
int TimeStep (void *)
 
int TimePostStep (void *)
 
int TimePrintStep (void *)
 
int TimeError (void *, void *, double *)
 
int TimeGetAuxSolutions (int *, double **, void *, int, int)
 

Detailed Description

Contains C++ function declarations for time integration.

Author
Debojyoti Ghosh

Definition in file timeintegration_cpp.h.

Function Documentation

◆ TimeInitialize()

int TimeInitialize ( void *  s,
int  nsims,
int  rank,
int  nproc,
void *  ts 
)

Initialize the time integration

Initialize time integration: This function initializes all that is required for time integration.

  • It sets the number of iterations, time step size, simulation time, etc.
  • It allocates solution, right-hand-side, and stage solution arrays needed by specific time integration methods.
  • It calls the method-specific initialization functions.
Parameters
sArray of simulation objects of type SimulationObject
nsimsnumber of simulation objects
rankMPI rank of this process
nprocnumber of MPI processes
tsTime integration object of type TimeIntegration

Definition at line 28 of file TimeInitialize.c.

34 {
35  TimeIntegration* TS = (TimeIntegration*) ts;
37  int ns, d, i;
38 
39  if (!sim) return(1);
40 
41 
42  TS->simulation = sim;
43  TS->nsims = nsims;
44 
45  TS->rank = rank;
46  TS->nproc = nproc;
47 
48  TS->n_iter = sim[0].solver.n_iter;
49  TS->restart_iter = sim[0].solver.restart_iter;
50  TS->dt = sim[0].solver.dt;
51 
52  TS->waqt = (double) TS->restart_iter * TS->dt;
53  TS->max_cfl = 0.0;
54  TS->norm = 0.0;
55  TS->TimeIntegrate = sim[0].solver.TimeIntegrate;
56  TS->iter_wctime_total = 0.0;
57 
58  TS->u_offsets = (long*) calloc (nsims, sizeof(long));
59  TS->u_sizes = (long*) calloc (nsims, sizeof(long));
60 
61  TS->u_size_total = 0;
62  for (ns = 0; ns < nsims; ns++) {
63  TS->u_offsets[ns] = TS->u_size_total;
64  TS->u_sizes[ns] = sim[ns].solver.npoints_local_wghosts * sim[ns].solver.nvars ;
65  TS->u_size_total += TS->u_sizes[ns];
66  }
67 
68  TS->u = (double*) calloc (TS->u_size_total,sizeof(double));
69  TS->rhs = (double*) calloc (TS->u_size_total,sizeof(double));
70  _ArraySetValue_(TS->u ,TS->u_size_total,0.0);
71  _ArraySetValue_(TS->rhs,TS->u_size_total,0.0);
72 
73  /* initialize arrays to NULL, then allocate as necessary */
74  TS->U = NULL;
75  TS->Udot = NULL;
76  TS->BoundaryFlux = NULL;
77 
78  TS->bf_offsets = (int*) calloc (nsims, sizeof(int));
79  TS->bf_sizes = (int*) calloc (nsims, sizeof(int));
80  TS->bf_size_total = 0;
81  for (ns = 0; ns < nsims; ns++) {
82  TS->bf_offsets[ns] = TS->bf_size_total;
83  TS->bf_sizes[ns] = 2 * sim[ns].solver.ndims * sim[ns].solver.nvars;
84  TS->bf_size_total += TS->bf_sizes[ns];
85  }
86 
87 #if defined(HAVE_CUDA)
88  if (sim[0].solver.use_gpu) {
89 
90  if (!strcmp(sim[0].solver.time_scheme,_RK_)) {
91 
92  /* explicit Runge-Kutta methods */
93  ExplicitRKParameters *params = (ExplicitRKParameters*) sim[0].solver.msti;
94  int nstages = params->nstages;
95 
96  gpuMalloc((void**)&TS->gpu_U, TS->u_size_total*nstages*sizeof(double));
97  gpuMalloc((void**)&TS->gpu_Udot, TS->u_size_total*nstages*sizeof(double));
98  gpuMalloc((void**)&TS->gpu_BoundaryFlux, TS->bf_size_total*nstages*sizeof(double));
99  gpuMemset(TS->gpu_U, 0, TS->u_size_total*nstages*sizeof(double));
100  gpuMemset(TS->gpu_Udot, 0, TS->u_size_total*nstages*sizeof(double));
101  gpuMemset(TS->gpu_BoundaryFlux, 0, TS->bf_size_total*nstages*sizeof(double));
102 
103  } else {
104 
105  fprintf(stderr,"ERROR in TimeInitialize(): %s is not yet implemented on GPUs.\n",
106  sim[0].solver.time_scheme );
107  return 1;
108 
109  }
110 
111  } else {
112 #endif
113 
114  if (!strcmp(sim[0].solver.time_scheme,_RK_)) {
115 
116  /* explicit Runge-Kutta methods */
117  ExplicitRKParameters *params = (ExplicitRKParameters*) sim[0].solver.msti;
118  int nstages = params->nstages;
119  TS->U = (double**) calloc (nstages,sizeof(double*));
120  TS->Udot = (double**) calloc (nstages,sizeof(double*));
121  for (i = 0; i < nstages; i++) {
122  TS->U[i] = (double*) calloc (TS->u_size_total,sizeof(double));
123  TS->Udot[i] = (double*) calloc (TS->u_size_total,sizeof(double));
124  }
125 
126  TS->BoundaryFlux = (double**) calloc (nstages,sizeof(double*));
127  for (i=0; i<nstages; i++) {
128  TS->BoundaryFlux[i] = (double*) calloc (TS->bf_size_total,sizeof(double));
129  }
130 
131  } else if (!strcmp(sim[0].solver.time_scheme,_FORWARD_EULER_)) {
132 
133  int nstages = 1;
134  TS->BoundaryFlux = (double**) calloc (nstages,sizeof(double*));
135  for (i=0; i<nstages; i++) {
136  TS->BoundaryFlux[i] = (double*) calloc (TS->bf_size_total,sizeof(double));
137  }
138 
139  } else if (!strcmp(sim[0].solver.time_scheme,_GLM_GEE_)) {
140 
141  /* General Linear Methods with Global Error Estimate */
142  GLMGEEParameters *params = (GLMGEEParameters*) sim[0].solver.msti;
143  int nstages = params->nstages;
144  int r = params->r;
145  TS->U = (double**) calloc (2*r-1 ,sizeof(double*));
146  TS->Udot = (double**) calloc (nstages,sizeof(double*));
147  for (i=0; i<2*r-1; i++) TS->U[i] = (double*) calloc (TS->u_size_total,sizeof(double));
148  for (i=0; i<nstages; i++) TS->Udot[i] = (double*) calloc (TS->u_size_total,sizeof(double));
149 
150  TS->BoundaryFlux = (double**) calloc (nstages,sizeof(double*));
151  for (i=0; i<nstages; i++) {
152  TS->BoundaryFlux[i] = (double*) calloc (TS->bf_size_total,sizeof(double));
153  }
154 
155 
156  if (!strcmp(params->ee_mode,_GLM_GEE_YYT_)) {
157  for (ns = 0; ns < nsims; ns++) {
158  for (i=0; i<r-1; i++) {
159  _ArrayCopy1D_( (sim[ns].solver.u),
160  (TS->U[r+i] + TS->u_offsets[ns]),
161  (TS->u_sizes[ns]) );
162  }
163  }
164  } else {
165  for (i=0; i<r-1; i++) {
166  _ArraySetValue_(TS->U[r+i],TS->u_size_total,0.0);
167  }
168  }
169  }
170 #if defined(HAVE_CUDA)
171  }
172 #endif
173 
174  /* set right-hand side function pointer */
176 
177  /* open files for writing */
178  if (!rank) {
179  if (sim[0].solver.write_residual) TS->ResidualFile = (void*) fopen("residual.out","w");
180  else TS->ResidualFile = NULL;
181  } else TS->ResidualFile = NULL;
182 
183  for (ns = 0; ns < nsims; ns++) {
184  sim[ns].solver.time_integrator = TS;
185  }
186 
187  return 0;
188 }
Structure of variables/parameters and function pointers for time integration.
int nvars
Definition: hypar.h:29
void gpuMemset(void *, int, size_t)
int TimeRHSFunctionExplicit(double *, double *, void *, void *, double)
Structure defining a simulation.
void * time_integrator
Definition: hypar.h:165
int n_iter
Definition: hypar.h:55
#define _FORWARD_EULER_
int restart_iter
Definition: hypar.h:58
int(* RHSFunction)(double *, double *, void *, void *, double)
int ndims
Definition: hypar.h:26
int(* TimeIntegrate)(void *)
void gpuMalloc(void **, size_t)
#define _GLM_GEE_YYT_
char ee_mode[_MAX_STRING_SIZE_]
#define _ArraySetValue_(x, size, value)
#define _GLM_GEE_
double dt
Definition: hypar.h:67
#define _RK_
Structure containing the parameters for an explicit Runge-Kutta method.
Structure containing the parameters for the General Linear Methods with Global Error Estimators (GLM-...
int npoints_local_wghosts
Definition: hypar.h:42
#define _ArrayCopy1D_(x, y, size)

◆ TimeCleanup()

int TimeCleanup ( void *  ts)

Clean up variables related to time integration

Clean up all allocations related to time integration

Parameters
tsObject of type TimeIntegration

Definition at line 18 of file TimeCleanup.c.

19 {
20  TimeIntegration* TS = (TimeIntegration*) ts;
22  int ns, nsims = TS->nsims;
23 
24  /* close files opened for writing */
25  if (!TS->rank) if (sim[0].solver.write_residual) fclose((FILE*)TS->ResidualFile);
26 
27 #if defined(HAVE_CUDA)
28  if (sim[0].solver.use_gpu) {
29  if (!strcmp(sim[0].solver.time_scheme,_RK_)) {
30  gpuFree(TS->gpu_U);
31  gpuFree(TS->gpu_Udot);
33  }
34  } else {
35 #endif
36  if (!strcmp(sim[0].solver.time_scheme,_RK_)) {
37  int i;
38  ExplicitRKParameters *params = (ExplicitRKParameters*) sim[0].solver.msti;
39  for (i=0; i<params->nstages; i++) free(TS->U[i]); free(TS->U);
40  for (i=0; i<params->nstages; i++) free(TS->Udot[i]); free(TS->Udot);
41  for (i=0; i<params->nstages; i++) free(TS->BoundaryFlux[i]); free(TS->BoundaryFlux);
42  } else if (!strcmp(sim[0].solver.time_scheme,_FORWARD_EULER_)) {
43  int nstages = 1, i;
44  for (i=0; i<nstages; i++) free(TS->BoundaryFlux[i]); free(TS->BoundaryFlux);
45  } else if (!strcmp(sim[0].solver.time_scheme,_GLM_GEE_)) {
46  int i;
47  GLMGEEParameters *params = (GLMGEEParameters*) sim[0].solver.msti;
48  for (i=0; i<2*params->r-1 ; i++) free(TS->U[i]); free(TS->U);
49  for (i=0; i<params->nstages; i++) free(TS->Udot[i]); free(TS->Udot);
50  for (i=0; i<params->nstages; i++) free(TS->BoundaryFlux[i]); free(TS->BoundaryFlux);
51  }
52 #if defined(HAVE_CUDA)
53  }
54 #endif
55 
56  /* deallocate arrays */
57  free(TS->u_offsets);
58  free(TS->u_sizes);
59  free(TS->bf_offsets);
60  free(TS->bf_sizes);
61  free(TS->u );
62  free(TS->rhs);
63  for (ns = 0; ns < nsims; ns++) {
64  sim[ns].solver.time_integrator = NULL;
65  }
66  return(0);
67 }
Structure of variables/parameters and function pointers for time integration.
Structure defining a simulation.
void * time_integrator
Definition: hypar.h:165
#define _FORWARD_EULER_
#define _GLM_GEE_
#define _RK_
Structure containing the parameters for an explicit Runge-Kutta method.
Structure containing the parameters for the General Linear Methods with Global Error Estimators (GLM-...
void gpuFree(void *)

◆ TimePreStep()

int TimePreStep ( void *  ts)

Function called at the beginning of a time step

Pre-time-step function: This function is called before each time step. Some notable things this does are:

  • Computes CFL and diffusion numbers.
  • Call the physics-specific pre-time-step function, if defined.
Parameters
tsObject of type TimeIntegration

Definition at line 22 of file TimePreStep.c.

23 {
24  TimeIntegration* TS = (TimeIntegration*) ts;
26 
28  int ns, nsims = TS->nsims;
29 
30  gettimeofday(&TS->iter_start_time,NULL);
31 
32  for (ns = 0; ns < nsims; ns++) {
33 
34  HyPar* solver = &(sim[ns].solver);
35  MPIVariables* mpi = &(sim[ns].mpi);
36 
37  double *u = NULL;
38 #if defined(HAVE_CUDA)
39  if (solver->use_gpu) {
40  u = solver->gpu_u;
41  } else{
42 #endif
43  u = solver->u;
44 #if defined(HAVE_CUDA)
45  }
46 #endif
47 
48  /* apply boundary conditions and exchange data over MPI interfaces */
49 
50  solver->ApplyBoundaryConditions( solver,
51  mpi,
52  u,
53  NULL,
54  TS->waqt );
55 
56  solver->ApplyIBConditions( solver,
57  mpi,
58  u,
59  TS->waqt );
60 
61 #if defined(HAVE_CUDA)
62  if (solver->use_gpu) {
64  solver->nvars,
65  solver->gpu_dim_local,
66  solver->ghosts,
67  mpi,
68  u );
69  } else {
70 #endif
72  solver->nvars,
73  solver->dim_local,
74  solver->ghosts,
75  mpi,
76  u );
77 #if defined(HAVE_CUDA)
78  }
79 #endif
80 
81  if ((TS->iter+1)%solver->screen_op_iter == 0) {
82 
83 #if defined(HAVE_CUDA)
84  if (!solver->use_gpu) {
85 #endif
86  _ArrayCopy1D_( solver->u,
87  (TS->u + TS->u_offsets[ns]),
88  (solver->npoints_local_wghosts*solver->nvars) );
89 #if defined(HAVE_CUDA)
90  }
91 #endif
92 
93  /* compute max CFL and diffusion number over the domain */
94  if (solver->ComputeCFL) {
95  double local_max_cfl = -1.0;
96  local_max_cfl = solver->ComputeCFL (solver,mpi,TS->dt,TS->waqt);
97  MPIMax_double(&TS->max_cfl ,&local_max_cfl ,1,&mpi->world);
98  } else {
99  TS->max_cfl = -1;
100  }
101  if (solver->ComputeDiffNumber) {
102  double local_max_diff = -1.0;
103  local_max_diff = solver->ComputeDiffNumber (solver,mpi,TS->dt,TS->waqt);
104  MPIMax_double(&TS->max_diff,&local_max_diff,1,&mpi->world);
105  } else {
106  TS->max_diff = -1;
107  }
108 
109  }
110 
111  /* set the step boundary flux integral value to zero */
112 #if defined(HAVE_CUDA)
113  if (solver->use_gpu) {
114  gpuArraySetValue(solver->StepBoundaryIntegral,2*solver->ndims*solver->nvars,0.0);
115  } else {
116 #endif
117  _ArraySetValue_(solver->StepBoundaryIntegral,2*solver->ndims*solver->nvars,0.0);
118 #if defined(HAVE_CUDA)
119  }
120 #endif
121 
122  if (solver->PreStep) {
123  solver->PreStep(u,solver,mpi,TS->waqt);
124  }
125 
126  }
127 
128  return 0;
129 }
Structure of variables/parameters and function pointers for time integration.
int nvars
Definition: hypar.h:29
Structure defining a simulation.
int MPIExchangeBoundariesnD(int, int, int *, int, void *, double *)
int(* ApplyBoundaryConditions)(void *, void *, double *, double *, double)
Definition: hypar.h:214
double * u
Definition: hypar.h:116
int(* ApplyIBConditions)(void *, void *, double *, double)
Definition: hypar.h:217
void gpuArraySetValue(double *, int, double)
int ndims
Definition: hypar.h:26
double(* ComputeCFL)(void *, void *, double, double)
Definition: hypar.h:269
int * gpu_dim_local
Definition: hypar.h:455
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
int MPIMax_double(double *, double *, int, void *)
Definition: MPIMax.c:38
double * gpu_u
Definition: hypar.h:459
MPI_Comm world
int(* PreStep)(double *, void *, void *, double)
Definition: hypar.h:339
double(* ComputeDiffNumber)(void *, void *, double, double)
Definition: hypar.h:272
#define _ArraySetValue_(x, size, value)
int * dim_local
Definition: hypar.h:37
int gpuMPIExchangeBoundariesnD(int, int, const int *, int, void *, double *)
int ghosts
Definition: hypar.h:52
double * StepBoundaryIntegral
Definition: hypar.h:384
Structure of MPI-related variables.
int use_gpu
Definition: hypar.h:449
int npoints_local_wghosts
Definition: hypar.h:42
#define _DECLARE_IERR_
Definition: basic.h:17
#define _ArrayCopy1D_(x, y, size)
struct timeval iter_start_time
int screen_op_iter
Definition: hypar.h:168

◆ TimeStep()

int TimeStep ( void *  ts)

Take one step in time

Advance one time step.

Parameters
tsObject of type TimeIntegration

Definition at line 13 of file TimeStep.c.

14 {
15  TimeIntegration *TS = (TimeIntegration*) ts;
17 
18  if (TS->TimeIntegrate) { TS->TimeIntegrate(TS); }
19 
20  return(0);
21 }
Structure of variables/parameters and function pointers for time integration.
Structure defining a simulation.
int(* TimeIntegrate)(void *)

◆ TimePostStep()

int TimePostStep ( void *  ts)

Function called at the end of a time step

Post-time-step function: this function is called at the end of each time step.

  • It updates the current simulation time.
  • It calls functions to print information and to write transient solution to file.
  • It will also call any physics-specific post-time-step function, if defined.
Parameters
tsObject of type TimeIntegration

Definition at line 28 of file TimePostStep.c.

29 {
30  TimeIntegration* TS = (TimeIntegration*) ts;
32  int ns, nsims = TS->nsims;
33 
34  /* update current time */
35  TS->waqt += TS->dt;
36 
37  if ((TS->iter+1)%sim[0].solver.screen_op_iter == 0) {
38 
39 #if defined(HAVE_CUDA)
40  if (sim[0].solver.use_gpu) {
41  TS->norm = -1;
42  } else {
43 #endif
44  /* Calculate norm for this time step */
45  double sum = 0.0;
46  double npts = 0;
47  for (ns = 0; ns < nsims; ns++) {
48  _ArrayAXPY_(sim[ns].solver.u,-1.0,(TS->u+TS->u_offsets[ns]),TS->u_sizes[ns]);
49  sum += ArraySumSquarenD( sim[ns].solver.nvars,
50  sim[ns].solver.ndims,
51  sim[ns].solver.dim_local,
52  sim[ns].solver.ghosts,
53  sim[ns].solver.index,
54  (TS->u+TS->u_offsets[ns]) );
55  npts += (double)sim[ns].solver.npoints_global;
56  }
57 
58  double global_sum = 0;
59  MPISum_double( &global_sum,
60  &sum,1,
61  &(sim[0].mpi.world) );
62 
63  TS->norm = sqrt(global_sum/npts);
64 #if defined(HAVE_CUDA)
65  }
66 #endif
67 
68  /* write to file */
69  if (TS->ResidualFile) {
70  fprintf((FILE*)TS->ResidualFile,"%10d\t%E\t%E\n",TS->iter+1,TS->waqt,TS->norm);
71  }
72 
73  }
74 
75 
76 #if defined(HAVE_CUDA)
77  if (!sim[0].solver.use_gpu) {
78 #endif
79  for (ns = 0; ns < nsims; ns++) {
80 
81  if (!strcmp(sim[ns].solver.ConservationCheck,"yes")) {
82  /* calculate volume integral of the solution at this time step */
83  IERR sim[ns].solver.VolumeIntegralFunction( sim[ns].solver.VolumeIntegral,
84  sim[ns].solver.u,
85  &(sim[ns].solver),
86  &(sim[ns].mpi) ); CHECKERR(ierr);
87  /* calculate surface integral of the flux at this time step */
88  IERR sim[ns].solver.BoundaryIntegralFunction( &(sim[ns].solver),
89  &(sim[ns].mpi)); CHECKERR(ierr);
90  /* calculate the conservation error at this time step */
91  IERR sim[ns].solver.CalculateConservationError( &(sim[ns].solver),
92  &(sim[ns].mpi)); CHECKERR(ierr);
93  }
94 
95  if (sim[ns].solver.PostStep) {
96  sim[ns].solver.PostStep(sim[ns].solver.u,&(sim[ns].solver),&(sim[ns].mpi),TS->waqt,TS->iter);
97  }
98 
99  }
100 #if defined(HAVE_CUDA)
101  }
102 #endif
103 
104  gettimeofday(&TS->iter_end_time,NULL);
105  long long walltime;
106  walltime = ( (TS->iter_end_time.tv_sec * 1000000 + TS->iter_end_time.tv_usec)
107  - (TS->iter_start_time.tv_sec * 1000000 + TS->iter_start_time.tv_usec));
108  TS->iter_wctime = (double) walltime / 1000000.0;
109  TS->iter_wctime_total += TS->iter_wctime;
110 
111  double global_total = 0, global_wctime = 0, global_mpi_total = 0, global_mpi_wctime = 0;
112 
113  MPIMax_double(&global_wctime, &TS->iter_wctime, 1, &(sim[0].mpi.world));
114  MPIMax_double(&global_total, &TS->iter_wctime_total, 1, &(sim[0].mpi.world));
115 
116 #if defined(HAVE_CUDA)
117  MPIMax_double(&global_mpi_wctime, &sim[0].mpi.wctime, 1, &(sim[0].mpi.world));
118  MPIMax_double(&global_mpi_total, &sim[0].mpi.wctime_total, 1, &(sim[0].mpi.world));
119 #endif
120 
121  return(0);
122 }
Structure of variables/parameters and function pointers for time integration.
#define IERR
Definition: basic.h:16
#define CHECKERR(ierr)
Definition: basic.h:18
Structure defining a simulation.
int MPISum_double(double *, double *, int, void *)
Definition: MPISum.c:39
long sum(const std::vector< int > &a_iv)
Definition: std_vec_ops.h:18
int(* CalculateConservationError)(void *, void *)
Definition: hypar.h:392
int(* VolumeIntegralFunction)(double *, double *, void *, void *)
Definition: hypar.h:388
int MPIMax_double(double *, double *, int, void *)
Definition: MPIMax.c:38
#define _ArrayAXPY_(x, a, y, size)
int(* BoundaryIntegralFunction)(void *, void *)
Definition: hypar.h:390
int(* PostStep)(double *, void *, void *, double, int)
Definition: hypar.h:341
INLINE double ArraySumSquarenD(int, int, int *, int, int *, double *)
struct timeval iter_start_time
struct timeval iter_end_time
int screen_op_iter
Definition: hypar.h:168

◆ TimePrintStep()

int TimePrintStep ( void *  ts)

Print time integration related information

Print information to screen (also calls any physics-specific printing function, if defined).

Parameters
tsObject of type TimeIntegration

Definition at line 16 of file TimePrintStep.c.

17 {
18  TimeIntegration* TS = (TimeIntegration*) ts;
20  int ns, nsims = TS->nsims;
21 
22  if ((!TS->rank) && ((TS->iter+1)%sim[0].solver.screen_op_iter == 0)) {
23  if (nsims > 1) {
24  printf("--\n");
25  printf("iter=%7d, t=%1.3e\n", TS->iter+1, TS->waqt);
26  if (TS->max_cfl >= 0) printf(" CFL=%1.3E\n", TS->max_cfl);
27  if (TS->norm >= 0) printf(" norm=%1.4E\n", TS->norm);
28  printf(" wctime=%1.1E (s)\n",TS->iter_wctime);
29  } else {
30  printf("iter=%7d ",TS->iter+1 );
31  printf("t=%1.3E ",TS->waqt );
32  if (TS->max_cfl >= 0) printf("CFL=%1.3E ",TS->max_cfl );
33  if (TS->norm >= 0) printf("norm=%1.4E ",TS->norm );
34  printf("wctime: %1.1E (s) ",TS->iter_wctime);
35  }
36 
37  /* calculate and print conservation error */
38  if (!strcmp(sim[0].solver.ConservationCheck,"yes")) {
39 
40  double error = 0;
41  for (ns = 0; ns < nsims; ns++) {
42  int v;
43  for (v=0; v<sim[ns].solver.nvars; v++) {
44  error += (sim[ns].solver.ConservationError[v] * sim[ns].solver.ConservationError[v]);
45  }
46  }
47  error = sqrt(error);
48  printf(" cons_err=%1.4E\n", error);
49 
50  } else {
51 
52  if (nsims == 1) printf("\n");
53 
54  }
55 
56  /* print physics-specific info, if available */
57  for (ns = 0; ns < nsims; ns++) {
58  if (sim[ns].solver.PrintStep) {
59  if (nsims > 1) printf("Physics-specific output for domain %d:\n", ns);
60  sim[ns].solver.PrintStep( &(sim[ns].solver),
61  &(sim[ns].mpi),
62  TS->waqt );
63  }
64  }
65  if (nsims > 1) {
66  printf("--\n");
67  printf("\n");
68  }
69  }
70 
71  return(0);
72 }
int(* PrintStep)(void *, void *, double)
Definition: hypar.h:344
Structure of variables/parameters and function pointers for time integration.
int nvars
Definition: hypar.h:29
Structure defining a simulation.
double * ConservationError
Definition: hypar.h:374

◆ TimeError()

int TimeError ( void *  s,
void *  m,
double *  uex 
)

Compute/estimate error in solution

Computes the time integration error in the solution: If the time integration method chosen has a mechanism to compute the error in the numerical integration, it is computed in this function. In addition, if an exact solution is available (see function argument uex, the error of the computed solution (stored in HyPar::u) with respect to this exact solution is also computed. These errors are written to a text file whose name depends on the time integration method being used.

Time integration method with error estimation currently implemented:

Parameters
sSolver object of type HyPar
mMPI object of type MPIVariables
uexExact solution (stored with the same array layout as HyPar::u (may be NULL)

Definition at line 27 of file TimeError.c.

34 {
35  HyPar *solver = (HyPar*) s;
36  MPIVariables *mpi = (MPIVariables*) m;
38  int size = solver->npoints_local_wghosts * solver->nvars;
39  double sum = 0.0, global_sum = 0.0;
40  static const double tolerance = 1e-15;
41  if (!TS) return(0);
42 
43  if (!strcmp(solver->time_scheme,_GLM_GEE_)) {
44  /* For GLM-GEE methods, calculate the norm of the estimated global error */
45  GLMGEEParameters *params = (GLMGEEParameters*) solver->msti;
46  double error[6] = {0,0,0,0,0,0}, *Uerr;
47  if (!strcmp(params->ee_mode,_GLM_GEE_YEPS_)) Uerr = TS->U[params->r];
48  else {
49  Uerr = TS->U[0];
50  _ArraySubtract1D_(Uerr,solver->u,TS->U[params->r],size);
51  _ArrayScale1D_(Uerr,(1.0/(1.0-params->gamma)),size);
52  }
53 
54  /* calculate solution norm for relative errors */
55  double sol_norm[3] = {0.0,0.0,0.0};
56  /* L1 */
57  sum = ArraySumAbsnD (solver->nvars,solver->ndims,solver->dim_local,
58  solver->ghosts,solver->index,solver->u);
59  global_sum = 0; MPISum_double(&global_sum,&sum,1,&mpi->world);
60  sol_norm[0] = global_sum/((double)solver->npoints_global);
61  /* L2 */
62  sum = ArraySumSquarenD(solver->nvars,solver->ndims,solver->dim_local,
63  solver->ghosts,solver->index,solver->u);
64  global_sum = 0; MPISum_double(&global_sum,&sum,1,&mpi->world);
65  sol_norm[1] = sqrt(global_sum/((double)solver->npoints_global));
66  /* Linf */
67  sum = ArrayMaxnD (solver->nvars,solver->ndims,solver->dim_local,
68  solver->ghosts,solver->index,solver->u);
69  global_sum = 0; MPIMax_double(&global_sum,&sum,1,&mpi->world);
70  sol_norm[2] = global_sum;
71 
72  /* calculate L1 norm of error */
73  sum = ArraySumAbsnD (solver->nvars,solver->ndims,solver->dim_local,
74  solver->ghosts,solver->index,Uerr);
75  global_sum = 0; MPISum_double(&global_sum,&sum,1,&mpi->world);
76  error[0] = global_sum/((double)solver->npoints_global);
77  /* calculate L2 norm of error */
78  sum = ArraySumSquarenD(solver->nvars,solver->ndims,solver->dim_local,
79  solver->ghosts,solver->index,Uerr);
80  global_sum = 0; MPISum_double(&global_sum,&sum,1,&mpi->world);
81  error[1] = sqrt(global_sum/((double)solver->npoints_global));
82  /* calculate Linf norm of error */
83  sum = ArrayMaxnD (solver->nvars,solver->ndims,solver->dim_local,
84  solver->ghosts,solver->index,Uerr);
85  global_sum = 0; MPIMax_double(&global_sum,&sum,1,&mpi->world);
86  error[2] = global_sum;
87 
88  if (uex) {
89  _ArrayAXBY_(TS->Udot[0],1.0,solver->u,-1.0,uex,size);
90  _ArrayAXPY_(Uerr,-1.0,TS->Udot[0],size);
91  /* calculate L1 norm of error */
92  sum = ArraySumAbsnD (solver->nvars,solver->ndims,solver->dim_local,
93  solver->ghosts,solver->index,TS->Udot[0]);
94  global_sum = 0; MPISum_double(&global_sum,&sum,1,&mpi->world);
95  error[3] = global_sum/((double)solver->npoints_global);
96  /* calculate L2 norm of error */
97  sum = ArraySumSquarenD(solver->nvars,solver->ndims,solver->dim_local,
98  solver->ghosts,solver->index,TS->Udot[0]);
99  global_sum = 0; MPISum_double(&global_sum,&sum,1,&mpi->world);
100  error[4] = sqrt(global_sum/((double)solver->npoints_global));
101  /* calculate Linf norm of error */
102  sum = ArrayMaxnD (solver->nvars,solver->ndims,solver->dim_local,
103  solver->ghosts,solver->index,TS->Udot[0]);
104  global_sum = 0; MPIMax_double(&global_sum,&sum,1,&mpi->world);
105  error[5] = global_sum;
106  } else error[3] = error[4] = error[5] = -1;
107 
108  if ( (sol_norm[0] > tolerance)
109  && (sol_norm[1] > tolerance)
110  && (sol_norm[2] > tolerance) ) {
111  error[0] /= sol_norm[0];
112  error[1] /= sol_norm[1];
113  error[2] /= sol_norm[2];
114  if (uex) {
115  error[3] /= sol_norm[0];
116  error[4] /= sol_norm[1];
117  error[5] /= sol_norm[2];
118  }
119  }
120 
121  /* write to file */
122  if (!mpi->rank) {
123  FILE *out;
124  out = fopen("glm_err.dat","w");
125  fprintf(out,"%1.16E %1.16E %1.16E %1.16E ",TS->dt,error[0],error[1],error[2]);
126  fprintf(out,"%1.16E %1.16E %1.16E\n",error[3],error[4],error[5]);
127  fclose(out);
128  printf("Estimated time integration errors (GLM-GEE time-integration):-\n");
129  printf(" L1 Error : %1.16E\n",error[0]);
130  printf(" L2 Error : %1.16E\n",error[1]);
131  printf(" Linfinity Error : %1.16E\n",error[2]);
132  }
133  }
134 
135  return(0);
136 }
Structure of variables/parameters and function pointers for time integration.
int nvars
Definition: hypar.h:29
void * time_integrator
Definition: hypar.h:165
#define _ArrayAXBY_(z, a, x, b, y, size)
int MPISum_double(double *, double *, int, void *)
Definition: MPISum.c:39
long sum(const std::vector< int > &a_iv)
Definition: std_vec_ops.h:18
double * u
Definition: hypar.h:116
int ndims
Definition: hypar.h:26
#define _ArraySubtract1D_(x, a, b, size)
INLINE double ArraySumAbsnD(int, int, int *, int, int *, double *)
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
int MPIMax_double(double *, double *, int, void *)
Definition: MPIMax.c:38
#define _ArrayAXPY_(x, a, y, size)
int npoints_global
Definition: hypar.h:40
char time_scheme[_MAX_STRING_SIZE_]
Definition: hypar.h:78
MPI_Comm world
int * index
Definition: hypar.h:102
int * dim_local
Definition: hypar.h:37
#define _GLM_GEE_
INLINE double ArrayMaxnD(int, int, int *, int, int *, double *)
int ghosts
Definition: hypar.h:52
Structure of MPI-related variables.
#define _ArrayScale1D_(x, a, size)
void * msti
Definition: hypar.h:366
Structure containing the parameters for the General Linear Methods with Global Error Estimators (GLM-...
int npoints_local_wghosts
Definition: hypar.h:42
INLINE double ArraySumSquarenD(int, int, int *, int, int *, double *)
#define _GLM_GEE_YEPS_

◆ TimeGetAuxSolutions()

int TimeGetAuxSolutions ( int *  N,
double **  uaux,
void *  s,
int  nu,
int  ns 
)

Function to get auxiliary solutions if available (for example, in GLM-GEE methods)

Return auxiliary solution: Some time integrators may have the concept of auxiliary solutions that they evolve along with the main solution HyPar::u (these may be used for error estimation, for example). This function returns a pointer to such an auxiliary solution. Note that the auxiliary solution has the same dimensions and array layout as the main solution.

  • Call with the final argument n less than 0 to get the total number of auxiliary solutions (this value will be stored in the argument N at exit).
  • Call with the final argument n less than or equal to zero to get the n-th auxiliary solution (the argument uaux will point to this at exit).

Note that auxiliary solutions are numbered in the C convention: 0,1,...,N-1.

Time integration methods which use auxiliary solutions currently implemented:

Parameters
NNumber of auxiliary solutions
uauxPointer to the array holding the auxiliary solution
sSolver object of type HyPar
nuIndex of the auxiliary solution to return
nsIndex of the simulation domain of which the auxiliary solution to return

Definition at line 29 of file TimeGetAuxSolutions.c.

36 {
37  HyPar *solver = (HyPar*) s;
39 
40  if (nu >= 0) {
41  if (!strcmp(solver->time_scheme,_GLM_GEE_)) {
42  GLMGEEParameters *params = (GLMGEEParameters*) solver->msti;
43  *uaux = (TS->U[params->r+nu] + TS->u_offsets[ns]);
44  }
45  } else {
46  if (!TS) *N = 0;
47  else {
48  if (!strcmp(solver->time_scheme,_GLM_GEE_)) {
49  GLMGEEParameters *params = (GLMGEEParameters*) solver->msti;
50  *N = params->r - 1;
51  } else *N = 0;
52  }
53  }
54 
55  return(0);
56 }
Structure of variables/parameters and function pointers for time integration.
void * time_integrator
Definition: hypar.h:165
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
char time_scheme[_MAX_STRING_SIZE_]
Definition: hypar.h:78
#define _GLM_GEE_
void * msti
Definition: hypar.h:366
Structure containing the parameters for the General Linear Methods with Global Error Estimators (GLM-...