HyPar  1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
Interp1PrimFifthOrderUpwindChar.c File Reference

Characteristic-based Fifth-order Upwind Scheme. More...

#include <stdlib.h>
#include <basic.h>
#include <arrayfunctions.h>
#include <matmult_native.h>
#include <mathfunctions.h>
#include <interpolation.h>
#include <mpivars.h>
#include <hypar.h>

Go to the source code of this file.

Macros

#define _MINIMUM_GHOSTS_   3
 

Functions

int Interp1PrimFifthOrderUpwindChar (double *fI, double *fC, double *u, double *x, int upw, int dir, void *s, void *m, int uflag)
 5th order upwind reconstruction (characteristic-based) on a uniform grid More...
 

Detailed Description

Characteristic-based Fifth-order Upwind Scheme.

Author
Debojyoti Ghosh

Definition in file Interp1PrimFifthOrderUpwindChar.c.

Macro Definition Documentation

◆ _MINIMUM_GHOSTS_

#define _MINIMUM_GHOSTS_   3

Minimum number of ghost points required for this interpolation method.

Definition at line 24 of file Interp1PrimFifthOrderUpwindChar.c.

Function Documentation

◆ Interp1PrimFifthOrderUpwindChar()

int Interp1PrimFifthOrderUpwindChar ( double *  fI,
double *  fC,
double *  u,
double *  x,
int  upw,
int  dir,
void *  s,
void *  m,
int  uflag 
)

5th order upwind reconstruction (characteristic-based) on a uniform grid

Computes the interpolated values of the first primitive of a function \({\bf f}\left({\bf u}\right)\) at the interfaces from the cell-centered values of the function using the fifth order upwind scheme on a uniform grid. The first primitive is defined as a function \({\bf h}\left({\bf u}\right)\) that satisfies:

\begin{equation} {\bf f}\left({\bf u}\left(x\right)\right) = \frac{1}{\Delta x} \int_{x-\Delta x/2}^{x+\Delta x/2} {\bf h}\left({\bf u}\left(\zeta\right)\right)d\zeta, \end{equation}

where \(x\) is the spatial coordinate along the dimension of the interpolation. This function computes the 5th order upwind numerical approximation \(\hat{\bf f}_{j+1/2} \approx {\bf h}_{j+1/2}\) as:

\begin{align} \hat{\alpha}^k_{j+1/2} = \frac{1}{30} {\alpha}^k_{j-2} - \frac{13}{60}{\alpha}^k_{j-1} + \frac{47}{60}{\alpha}^k_j + \frac{27}{60}{\alpha}^k_{j+1} - \frac{1}{20}{\alpha}^k_{j+2}, \end{align}

where

\begin{equation} \alpha^k = {\bf l}_k \cdot {\bf f},\ k=1,\cdots,n \end{equation}

is the \(k\)-th characteristic quantity, and \({\bf l}_k\) is the \(k\)-th left eigenvector, \({\bf r}_k\) is the \(k\)-th right eigenvector, and \(n\) is HyPar::nvars. The final interpolated function is computed from the interpolated characteristic quantities as:

\begin{equation} \hat{\bf f}_{j+1/2} = \sum_{k=1}^n \alpha^k_{j+1/2} {\bf r}_k \end{equation}

Implementation Notes:

  • This method assumes a uniform grid in the spatial dimension corresponding to the interpolation.
  • The method described above corresponds to a left-biased interpolation. The corresponding right-biased interpolation can be obtained by reflecting the equations about interface j+1/2.
  • The left and right eigenvectors are computed at an averaged quantity at j+1/2. Thus, this function requires functions to compute the average state, and the left and right eigenvectors. These are provided by the physical model through

    If these functions are not provided by the physical model, then a characteristic-based interpolation cannot be used.

  • The function computes the interpolant for the entire grid in one call. It loops over all the grid lines along the interpolation direction and carries out the 1D interpolation along these grid lines.
  • Location of cell-centers and cell interfaces along the spatial dimension of the interpolation is shown in the following figure:
    chap1_1Ddomain.png

Function arguments:

Argument Type Explanation ------—
fI double* Array to hold the computed interpolant at the grid interfaces. This array must have the same layout as the solution, but with no ghost points. Its size should be the same as u in all dimensions, except dir (the dimension along which to interpolate) along which it should be larger by 1 (number of interfaces is 1 more than the number of interior cell centers).
fC double* Array with the cell-centered values of the flux function \({\bf f}\left({\bf u}\right)\). This array must have the same layout and size as the solution, with ghost points.
u double* The solution array \({\bf u}\) (with ghost points). If the interpolation is characteristic based, this is needed to compute the eigendecomposition. For a multidimensional problem, the layout is as follows: u is a contiguous 1D array of size (nvars*dim[0]*dim[1]*...*dim[D-1]) corresponding to the multi-dimensional solution, with the following ordering - nvars, dim[0], dim[1], ..., dim[D-1], where nvars is the number of solution components (HyPar::nvars), dim is the local size (HyPar::dim_local), D is the number of spatial dimensions.
x double* The grid array (with ghost points). This is used only by non-uniform-grid interpolation methods. For multidimensional problems, the layout is as follows: x is a contiguous 1D array of size (dim[0]+dim[1]+...+dim[D-1]), with the spatial coordinates along dim[0] stored from 0,...,dim[0]-1, the spatial coordinates along dim[1] stored along dim[0],...,dim[0]+dim[1]-1, and so forth.
upw int Upwinding direction: if positive, a left-biased interpolant will be computed; if negative, a right-biased interpolant will be computed. If the interpolation method is central, then this has no effect.
dir int Spatial dimension along which to interpolate (eg: 0 for 1D; 0 or 1 for 2D; 0,1 or 2 for 3D)
s void* Solver object of type HyPar: the following variables are needed - HyPar::ghosts, HyPar::ndims, HyPar::nvars, HyPar::dim_local.
m void* MPI object of type MPIVariables: this is needed only by compact interpolation method that need to solve a global implicit system across MPI ranks.
uflag int A flag indicating if the function being interpolated \({\bf f}\) is the solution itself \({\bf u}\) (if 1, \({\bf f}\left({\bf u}\right) \equiv {\bf u}\)).

Reference:

Parameters
fIArray of interpolated function values at the interfaces
fCArray of cell-centered values of the function \({\bf f}\left({\bf u}\right)\)
uArray of cell-centered values of the solution \({\bf u}\)
xGrid coordinates
upwUpwind direction (left or right biased)
dirSpatial dimension along which to interpolation
sObject of type HyPar containing solver-related variables
mObject of type MPIVariables containing MPI-related variables
uflagFlag to indicate if \(f(u) \equiv u\), i.e, if the solution is being reconstructed

Definition at line 84 of file Interp1PrimFifthOrderUpwindChar.c.

95 {
96  HyPar *solver = (HyPar*) s;
97  int i, k, v;
99 
100  int ghosts = solver->ghosts;
101  int ndims = solver->ndims;
102  int nvars = solver->nvars;
103  int *dim = solver->dim_local;
104 
105  /* define some constants */
106  static const double one_by_thirty = 1.0/30.0,
107  thirteen_by_sixty = 13.0/60.0,
108  fortyseven_by_sixty = 47.0/60.0,
109  twentyseven_by_sixty = 27.0/60.0,
110  one_by_twenty = 1.0/20.0;
111 
112  /* create index and bounds for the outer loop, i.e., to loop over all 1D lines along
113  dimension "dir" */
114  int indexC[ndims], indexI[ndims], index_outer[ndims], bounds_outer[ndims], bounds_inter[ndims];
115  _ArrayCopy1D_(dim,bounds_outer,ndims); bounds_outer[dir] = 1;
116  _ArrayCopy1D_(dim,bounds_inter,ndims); bounds_inter[dir] += 1;
117  int N_outer; _ArrayProduct1D_(bounds_outer,ndims,N_outer);
118 
119  /* allocate arrays for the averaged state, eigenvectors and characteristic interpolated f */
120  double R[nvars*nvars], L[nvars*nvars], uavg[nvars], fchar[nvars];
121 
122 #pragma omp parallel for schedule(auto) default(shared) private(i,k,v,R,L,uavg,fchar,index_outer,indexC,indexI)
123  for (i=0; i<N_outer; i++) {
124  _ArrayIndexnD_(ndims,i,bounds_outer,index_outer,0);
125  _ArrayCopy1D_(index_outer,indexC,ndims);
126  _ArrayCopy1D_(index_outer,indexI,ndims);
127 
128  for (indexI[dir] = 0; indexI[dir] < dim[dir]+1; indexI[dir]++) {
129 
130  /* 1D indices of the stencil grid points */
131  int qm1,qm2,qm3,qp1,qp2;
132  if (upw > 0) {
133  indexC[dir] = indexI[dir]-3; _ArrayIndex1D_(ndims,dim,indexC,ghosts,qm3);
134  indexC[dir] = indexI[dir]-2; _ArrayIndex1D_(ndims,dim,indexC,ghosts,qm2);
135  indexC[dir] = indexI[dir]-1; _ArrayIndex1D_(ndims,dim,indexC,ghosts,qm1);
136  indexC[dir] = indexI[dir] ; _ArrayIndex1D_(ndims,dim,indexC,ghosts,qp1);
137  indexC[dir] = indexI[dir]+1; _ArrayIndex1D_(ndims,dim,indexC,ghosts,qp2);
138  } else {
139  indexC[dir] = indexI[dir]+2; _ArrayIndex1D_(ndims,dim,indexC,ghosts,qm3);
140  indexC[dir] = indexI[dir]+1; _ArrayIndex1D_(ndims,dim,indexC,ghosts,qm2);
141  indexC[dir] = indexI[dir] ; _ArrayIndex1D_(ndims,dim,indexC,ghosts,qm1);
142  indexC[dir] = indexI[dir]-1; _ArrayIndex1D_(ndims,dim,indexC,ghosts,qp1);
143  indexC[dir] = indexI[dir]-2; _ArrayIndex1D_(ndims,dim,indexC,ghosts,qp2);
144  }
145 
146  int p; /* 1D index of the interface */
147  _ArrayIndex1D_(ndims,bounds_inter,indexI,0,p);
148 
149  /* find averaged state at this interface */
150  IERR solver->AveragingFunction(uavg,&u[nvars*qm1],&u[nvars*qp1],solver->physics); CHECKERR(ierr);
151 
152  /* Get the left and right eigenvectors */
153  IERR solver->GetLeftEigenvectors (uavg,L,solver->physics,dir); CHECKERR(ierr);
154  IERR solver->GetRightEigenvectors (uavg,R,solver->physics,dir); CHECKERR(ierr);
155 
156  /* For each characteristic field */
157  for (v = 0; v < nvars; v++) {
158 
159  /* calculate the characteristic flux components along this characteristic */
160  double fm3, fm2, fm1, fp1, fp2;
161  fm3 = fm2 = fm1 = fp1 = fp2 = 0;
162  for (k = 0; k < nvars; k++) {
163  fm3 += L[v*nvars+k] * fC[qm3*nvars+k];
164  fm2 += L[v*nvars+k] * fC[qm2*nvars+k];
165  fm1 += L[v*nvars+k] * fC[qm1*nvars+k];
166  fp1 += L[v*nvars+k] * fC[qp1*nvars+k];
167  fp2 += L[v*nvars+k] * fC[qp2*nvars+k];
168  }
169 
170  /* fifth order upwind approximation of the characteristic flux */
171  fchar[v] = one_by_thirty * fm3
172  - thirteen_by_sixty * fm2
173  + fortyseven_by_sixty * fm1
174  + twentyseven_by_sixty * fp1
175  - one_by_twenty * fp2;
176  }
177 
178  /* calculate the interface u from the characteristic u */
179  IERR MatVecMult(nvars,(fI+nvars*p),R,fchar); CHECKERR(ierr);
180 
181  }
182  }
183 
184  return(0);
185 }
int nvars
Definition: hypar.h:29
#define IERR
Definition: basic.h:16
#define CHECKERR(ierr)
Definition: basic.h:18
int(* AveragingFunction)(double *, double *, double *, void *)
Definition: hypar.h:354
#define _ArrayIndexnD_(N, index, imax, i, ghost)
int ndims
Definition: hypar.h:26
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
#define _ArrayIndex1D_(N, imax, i, ghost, index)
int * dim_local
Definition: hypar.h:37
void * physics
Definition: hypar.h:266
int ghosts
Definition: hypar.h:52
#define MatVecMult(N, y, A, x)
#define _DECLARE_IERR_
Definition: basic.h:17
#define _ArrayCopy1D_(x, y, size)
int(* GetRightEigenvectors)(double *, double *, void *, int)
Definition: hypar.h:359
#define _ArrayProduct1D_(x, size, p)
int(* GetLeftEigenvectors)(double *, double *, void *, int)
Definition: hypar.h:357