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

First order upwind Scheme (Component-wise application to vectors). More...

#include <stdio.h>
#include <stdlib.h>
#include <basic.h>
#include <arrayfunctions.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 Interp1PrimFirstOrderUpwind (double *fI, double *fC, double *u, double *x, int upw, int dir, void *s, void *m, int uflag)
 1st order upwind reconstruction (component-wise) on a uniform grid More...
 

Detailed Description

First order upwind Scheme (Component-wise application to vectors).

Author
Debojyoti Ghosh

Definition in file Interp1PrimFirstOrderUpwind.c.

Macro Definition Documentation

◆ _MINIMUM_GHOSTS_

#define _MINIMUM_GHOSTS_   1

Minimum number of ghost points required for this interpolation method.

Definition at line 23 of file Interp1PrimFirstOrderUpwind.c.

Function Documentation

◆ Interp1PrimFirstOrderUpwind()

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

1st order upwind reconstruction (component-wise) 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{\bf f}_{j+1/2} = \left\{\begin{array}{cc} {\bf f}_{j} & {\rm upw} > 0 \\ {\bf f}_{j+1} & {\rm upw} \le 0 \end{array}\right.. \end{equation}

Implementation Notes:

  • The scalar interpolation method is applied to the vector function in a component-wise manner.
  • 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 61 of file Interp1PrimFirstOrderUpwind.c.

72 {
73  HyPar *solver = (HyPar*) s;
74 
75  int ghosts = solver->ghosts;
76  int ndims = solver->ndims;
77  int nvars = solver->nvars;
78  int *dim = solver->dim_local;
79 
80  /* create index and bounds for the outer loop, i.e., to loop over all 1D lines along
81  dimension "dir" */
82  int indexC[ndims], indexI[ndims], index_outer[ndims], bounds_outer[ndims], bounds_inter[ndims];
83  _ArrayCopy1D_(dim,bounds_outer,ndims); bounds_outer[dir] = 1;
84  _ArrayCopy1D_(dim,bounds_inter,ndims); bounds_inter[dir] += 1;
85  int N_outer; _ArrayProduct1D_(bounds_outer,ndims,N_outer);
86 
87  int i;
88 #pragma omp parallel for schedule(auto) default(shared) private(i,index_outer,indexC,indexI)
89  for (i=0; i<N_outer; i++) {
90  _ArrayIndexnD_(ndims,i,bounds_outer,index_outer,0);
91  _ArrayCopy1D_(index_outer,indexC,ndims);
92  _ArrayCopy1D_(index_outer,indexI,ndims);
93  for (indexI[dir] = 0; indexI[dir] < dim[dir]+1; indexI[dir]++) {
94  indexC[dir] = (upw > 0 ? indexI[dir]-1 : indexI[dir]);
95  int p; _ArrayIndex1D_(ndims,bounds_inter,indexI,0 ,p);
96  int q; _ArrayIndex1D_(ndims,dim ,indexC,ghosts,q);
97  int v; for (v=0; v<nvars; v++) fI[p*nvars+v] = fC[q*nvars+v];
98  }
99  }
100 
101  return(0);
102 }
int nvars
Definition: hypar.h:29
#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
int ghosts
Definition: hypar.h:52
#define _ArrayCopy1D_(x, y, size)
#define _ArrayProduct1D_(x, size, p)