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

Characteristic-based first order upwind scheme. More...

#include <stdio.h>
#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_   1
 

Functions

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

Detailed Description

Characteristic-based first order upwind scheme.

Author
Debojyoti Ghosh

Definition in file Interp1PrimFirstOrderUpwindChar.c.

Macro Definition Documentation

◆ _MINIMUM_GHOSTS_

#define _MINIMUM_GHOSTS_   1

Minimum number of ghost points required for this interpolation method.

Definition at line 25 of file Interp1PrimFirstOrderUpwindChar.c.

Function Documentation

◆ Interp1PrimFirstOrderUpwindChar()

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

1st 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 1st 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 1st order upwind numerical approximation \(\hat{\bf f}_{j+1/2} \approx {\bf h}_{j+1/2}\) as:

\begin{equation} \hat{\alpha}^k_{j+1/2} = \left\{\begin{array}{cc} {\alpha}^k_{j} & {\rm upw} > 0 \\ {\alpha}^k_{j+1} & {\rm upw} \le 0 \end{array}\right.. \end{equation}

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:

  • 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}\)).
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 77 of file Interp1PrimFirstOrderUpwindChar.c.

88 {
89  HyPar *solver = (HyPar*) s;
90  int i, k, v;
92 
93  int ghosts = solver->ghosts;
94  int ndims = solver->ndims;
95  int nvars = solver->nvars;
96  int *dim = solver->dim_local;
97 
98  /* create index and bounds for the outer loop, i.e., to loop over all 1D lines along
99  dimension "dir" */
100  int indexC[ndims], indexI[ndims], index_outer[ndims], bounds_outer[ndims], bounds_inter[ndims];
101  _ArrayCopy1D_(dim,bounds_outer,ndims); bounds_outer[dir] = 1;
102  _ArrayCopy1D_(dim,bounds_inter,ndims); bounds_inter[dir] += 1;
103  int N_outer; _ArrayProduct1D_(bounds_outer,ndims,N_outer);
104 
105  /* allocate arrays for the averaged state, eigenvectors and characteristic interpolated f */
106  double R[nvars*nvars], L[nvars*nvars], uavg[nvars], fchar[nvars];
107 
108 #pragma omp parallel for schedule(auto) default(shared) private(i,k,v,R,L,uavg,fchar,index_outer,indexC,indexI)
109  for (i=0; i<N_outer; i++) {
110  _ArrayIndexnD_(ndims,i,bounds_outer,index_outer,0);
111  _ArrayCopy1D_(index_outer,indexC,ndims);
112  _ArrayCopy1D_(index_outer,indexI,ndims);
113 
114  for (indexI[dir] = 0; indexI[dir] < dim[dir]+1; indexI[dir]++) {
115 
116  int pL, pR; /* 1D index of the left and right cells */
117  indexC[dir] = indexI[dir]-1; _ArrayIndex1D_(ndims,dim,indexC,ghosts,pL);
118  indexC[dir] = indexI[dir] ; _ArrayIndex1D_(ndims,dim,indexC,ghosts,pR);
119  int p; /* 1D index of the interface */
120  _ArrayIndex1D_(ndims,bounds_inter,indexI,0,p);
121 
122  /* find averaged state at this interface */
123  IERR solver->AveragingFunction(uavg,&u[nvars*pL],&u[nvars*pR],solver->physics); CHECKERR(ierr);
124 
125  /* Get the left and right eigenvectors */
126  IERR solver->GetLeftEigenvectors (uavg,L,solver->physics,dir); CHECKERR(ierr);
127  IERR solver->GetRightEigenvectors (uavg,R,solver->physics,dir); CHECKERR(ierr);
128 
129  /* For each characteristic field */
130  for (v = 0; v < nvars; v++) {
131 
132  /* calculate the characteristic flux components along this characteristic */
133  double fcL = 0, fcR = 0;
134  for (k = 0; k < nvars; k++) {
135  fcL += L[v*nvars+k] * fC[pL*nvars+k];
136  fcR += L[v*nvars+k] * fC[pR*nvars+k];
137  }
138 
139  /* first order upwind approximation of the characteristic flux */
140  fchar[v] = (upw > 0 ? fcL : fcR);
141 
142  }
143 
144  /* calculate the interface u from the characteristic u */
145  IERR MatVecMult(nvars,(fI+nvars*p),R,fchar); CHECKERR(ierr);
146 
147  }
148  }
149 
150  return(0);
151 }
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