accelerInt  v0.1
Functions | Variables
solver_main.cu File Reference

the generic main file for all GPU solvers More...

#include <stdlib.h>
#include <stdio.h>
#include <math.h>
#include <string.h>
#include <stdbool.h>
#include <complex.h>
#include <sys/types.h>
#include <sys/stat.h>
#include <cuda.h>
#include <cuda_runtime.h>
#include <helper_cuda.h>
#include <cuComplex.h>
#include "header.cuh"
#include "timer.h"
#include "solver.cuh"
#include "gpu_memory.cuh"
#include "read_initial_conditions.cuh"
#include "launch_bounds.cuh"
#include <assert.h>
Include dependency graph for solver_main.cu:

Go to the source code of this file.

Functions

void write_log (int NUM, double t, const double *y_host, FILE *pFile)
 Writes state vectors to file. More...
 
void memcpy2D_in (double *dst, const int pitch_dst, double const *src, const int pitch_src, const int offset, const size_t width, const int height)
 A convienience method to copy memory between host pointers of different pitches, widths and heights. Enables easier use of CUDA's cudaMemcpy2D functions. More...
 
void memcpy2D_out (double *dst, const int pitch_dst, double const *src, const int pitch_src, const int offset, const size_t width, const int height)
 A convienience method to copy memory between host pointers of different pitches, widths and heights. Enables easier use of CUDA's cudaMemcpy2D functions. More...
 
int main (int argc, char *argv[])
 

Variables

__device__ int integrator_steps [DIVERGENCE_TEST] = {0}
 If DIVERGENCE_TEST is defined, this creates a device array for tracking. More...
 

Detailed Description

the generic main file for all GPU solvers

Author
Nicholas Curtis
Date
03/09/2015

Contains main function, setup, initialization, logging, timing and driver functions

Definition in file solver_main.cu.

Function Documentation

◆ main()

int main ( int  argc,
char *  argv[] 
)

Main function

Parameters
[in]argccommand line argument count
[in]argvcommand line argument vector

This allows running the integrators from the command line. The syntax is as follows:
./solver-name [num_threads] [num_IVPs]

  • num_threads [Optional, Default:1]
    • The number OpenMP threads to utilize
    • The number of threads cannot be greater than recognized by OpenMP via omp_get_max_threads()
  • num_IVPs [Optional, Default:1]
    • The number of initial value problems to solve.
    • This must be less than the number of conditions in the data file if SAME_IC is not defined.
    • If SAME_IC is defined, then the initial conditions in the mechanism files will be used.

Number of independent systems

Definition at line 150 of file solver_main.cu.

◆ memcpy2D_in()

void memcpy2D_in ( double *  dst,
const int  pitch_dst,
double const *  src,
const int  pitch_src,
const int  offset,
const size_t  width,
const int  height 
)
inline

A convienience method to copy memory between host pointers of different pitches, widths and heights. Enables easier use of CUDA's cudaMemcpy2D functions.

Parameters
[out]dstThe destination array
[in]pitch_dstThe width (in bytes) of the destination array. This corresponds to the padded number of IVPs to be solved.
[in]srcThe source pointer
[in]pitch_srcThe width (in bytes) of the source array. This corresponds to the (non-padded) number of IVPs read by read_initial_conditions
[in]offsetThe offset within the source array (IVP index) to copy from. This is useful in the case (for large models) where the solver and state vector memory will not fit in device memory and the integration must be split into multiple kernel calls.
[in]widthThe size (in bytes) of memory to copy for each entry in the state vector
[in]heightThe number of entries in the state vector

Definition at line 96 of file solver_main.cu.

◆ memcpy2D_out()

void memcpy2D_out ( double *  dst,
const int  pitch_dst,
double const *  src,
const int  pitch_src,
const int  offset,
const size_t  width,
const int  height 
)
inline

A convienience method to copy memory between host pointers of different pitches, widths and heights. Enables easier use of CUDA's cudaMemcpy2D functions.

Parameters
[out]dstThe destination array
[in]pitch_dstThe width (in bytes) of the source array. This corresponds to the (non-padded) number of IVPs read by read_initial_conditions
[in]srcThe source pointer
[in]pitch_srcThe width (in bytes) of the destination array. This corresponds to the padded number of IVPs to be solved.
[in]offsetThe offset within the destination array (IVP index) to copy to. This is useful in the case (for large models) where the solver and state vector memory will not fit in device memory and the integration must be split into multiple kernel calls.
[in]widthThe size (in bytes) of memory to copy for each entry in the state vector
[in]heightThe number of entries in the state vector

Definition at line 122 of file solver_main.cu.

◆ write_log()

void write_log ( int  NUM,
double  t,
const double *  y_host,
FILE *  pFile 
)

Writes state vectors to file.

Parameters
[in]NUMthe number of state vectors to write
[in]tthe current system time
[in]y_hostthe current state vectors
[in]pFilethe opened binary file object

The resulting file is updated as: system time
temperature, mass fractions (State #1)
temperature, mass fractions (State #2)...

Definition at line 59 of file solver_main.cu.

Variable Documentation

◆ integrator_steps

__device__ int integrator_steps[DIVERGENCE_TEST] = {0}

If DIVERGENCE_TEST is defined, this creates a device array for tracking.

Definition at line 44 of file solver_main.cu.