HyPar  1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
ParabolicFunctionNC1.5Stage.c File Reference

Evaluate the parabolic term using a "1.5"-stage discretization. More...

#include <stdlib.h>
#include <basic.h>
#include <arrayfunctions.h>
#include <mpivars.h>
#include <hypar.h>

Go to the source code of this file.

Functions

int ParabolicFunctionNC1_5Stage (double *par, double *u, void *s, void *m, double t)
 

Detailed Description

Evaluate the parabolic term using a "1.5"-stage discretization.

Author
Debojyoti Ghosh

Definition in file ParabolicFunctionNC1.5Stage.c.

Function Documentation

int ParabolicFunctionNC1_5Stage ( double *  par,
double *  u,
void *  s,
void *  m,
double  t 
)

Evaluate the parabolic term using a "1.5"-stage finite-difference spatial discretization: The parabolic term is assumed to be of the form:

\begin{equation} {\bf P}\left({\bf u}\right) = \sum_{d1=0}^{D-1}\sum_{d2=0}^{D-1} \frac {\partial^2 h_{d1,d2}\left(\bf u\right)} {\partial x_{d1} \partial x_{d2}}, \end{equation}

where \(d1\) and \(d2\) are spatial dimension indices, and \(D\) is the total number of spatial dimensions (HyPar::ndims). This term is discretized at a grid point as:

\begin{equation} \left.{\bf P}\left({\bf u}\right)\right|_j = \sum_{d1=0}^{D-1} \sum_{d2=0}^{D-1} \frac { \mathcal{D}_{d1}\mathcal{D}_{d2} \left[ {\bf h}_{d1,d2} \right] } {\Delta x_{d1} \Delta x_{d2}}, \end{equation}

where \(\mathcal{D}\) denotes the finite-difference approximation to the first derivative.

  • When \(d = d1 = d2\), the operator \(\mathcal{D}^2_d\left[\cdot\right] = \mathcal{D}_{d}\mathcal{D}_{d}\left[\cdot\right]\) is computed in one step as the finite-difference approximation to the second derivative (Laplacian) \(\mathcal{L}_d\left[\cdot\right]\) using HyPar::SecondDerivativePar.
  • When \(d1 \ne d2 \), each of the first derivative approximations are \(\mathcal{D}_{d1}\) and \(\mathcal{D}_{d2}\) are computed separately, and thus the cross-derivative is evaluated in two steps using HyPar::FirstDerivativePar.

Note: this form of the parabolic term does allow for cross-derivatives ( \( d1 \ne d2 \)).

To use this form of the parabolic term:

  • specify "par_space_type" in solver.inp as "nonconservative-1.5stage" (HyPar::spatial_type_par).
  • the physical model must specify \({\bf h}_{d1,d2}\left({\bf u}\right)\) through HyPar::HFunction.
See Also
ParabolicFunctionNC2Stage()
Parameters
pararray to hold the computed parabolic term
usolution
sSolver object of type HyPar
mMPI object of type MPIVariables
tCurrent simulation time

Definition at line 35 of file ParabolicFunctionNC1.5Stage.c.

42 {
43  HyPar *solver = (HyPar*) s;
44  MPIVariables *mpi = (MPIVariables*) m;
45  double *Func = solver->fluxC;
46  double *Deriv1 = solver->Deriv1;
47  double *Deriv2 = solver->Deriv2;
48  int d, d1, d2, v, p, done;
49  double dxinv1, dxinv2;
51 
52  int ndims = solver->ndims;
53  int nvars = solver->nvars;
54  int ghosts = solver->ghosts;
55  int *dim = solver->dim_local;
56  double *dxinv = solver->dxinv;
57  int size = solver->npoints_local_wghosts;
58 
59  if (!solver->HFunction) return(0); /* zero parabolic terms */
60  solver->count_par++;
61 
62  int index[ndims];
63  _ArraySetValue_(par,size*nvars,0.0);
64 
65  for (d1 = 0; d1 < ndims; d1++) {
66  for (d2 = 0; d2 < ndims; d2++) {
67 
68  /* calculate the diffusion function */
69  IERR solver->HFunction(Func,u,d1,d2,solver,t); CHECKERR(ierr);
70  if (d1 == d2) {
71  /* second derivative with respect to the same independent variable */
72  IERR solver->SecondDerivativePar(Deriv2,Func,d1,solver,mpi); CHECKERR(ierr);
73  } else {
74  /* cross derivative */
75  IERR solver->FirstDerivativePar(Deriv1,Func ,d1, 1,solver,mpi); CHECKERR(ierr);
76  IERR MPIExchangeBoundariesnD(ndims,nvars,dim,ghosts,mpi,Deriv1); CHECKERR(ierr);
77  IERR solver->FirstDerivativePar(Deriv2,Deriv1,d2,-1,solver,mpi); CHECKERR(ierr);
78  }
79 
80  /* calculate the final term - second derivative of the diffusion function */
81  done = 0; _ArraySetValue_(index,ndims,0);
82  while (!done) {
83  _ArrayIndex1D_(ndims,dim,index,ghosts,p);
84  _GetCoordinate_(d1,index[d1],dim,ghosts,dxinv,dxinv1);
85  _GetCoordinate_(d2,index[d2],dim,ghosts,dxinv,dxinv2);
86  for (v=0; v<nvars; v++) par[nvars*p+v] += (dxinv1*dxinv2 * Deriv2[nvars*p+v]);
87  _ArrayIncrementIndex_(ndims,dim,index,done);
88  }
89 
90  }
91  }
92 
93  if (solver->flag_ib) _ArrayBlockMultiply_(par,solver->iblank,size,nvars);
94  return(0);
95 }
int npoints_local_wghosts
Definition: hypar.h:42
#define _ArraySetValue_(x, size, value)
double * iblank
Definition: hypar.h:436
#define _ArrayIncrementIndex_(N, imax, i, done)
double * Deriv1
Definition: hypar.h:151
int(* HFunction)(double *, double *, int, int, void *, double)
Definition: hypar.h:313
double * Deriv2
Definition: hypar.h:151
int flag_ib
Definition: hypar.h:441
int * dim_local
Definition: hypar.h:37
#define _GetCoordinate_(dir, i, dim, ghosts, x, coord)
Definition: basic.h:31
int ghosts
Definition: hypar.h:52
#define _ArrayIndex1D_(N, imax, i, ghost, index)
int count_par
Definition: hypar.h:418
#define _ArrayBlockMultiply_(x, a, n, bs)
int(* FirstDerivativePar)(double *, double *, int, int, void *, void *)
Definition: hypar.h:243
double * fluxC
Definition: hypar.h:128
int nvars
Definition: hypar.h:29
#define CHECKERR(ierr)
Definition: basic.h:18
double * dxinv
Definition: hypar.h:110
int ndims
Definition: hypar.h:26
int(* SecondDerivativePar)(double *, double *, int, void *, void *)
Definition: hypar.h:247
#define IERR
Definition: basic.h:16
int MPIExchangeBoundariesnD(int, int, int *, int, void *, double *)
Structure of MPI-related variables.
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
#define _DECLARE_IERR_
Definition: basic.h:17