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

Pre-step function for 2D Navier Stokes equations. More...

#include <arrayfunctions.h>
#include <mathfunctions.h>
#include <matmult_native.h>
#include <physicalmodels/navierstokes2d.h>
#include <hypar.h>

Go to the source code of this file.

Functions

int NavierStokes2DPreStep (double *u, void *s, void *m, double waqt)
 

Detailed Description

Pre-step function for 2D Navier Stokes equations.

Author
Debojyoti Ghosh

Definition in file NavierStokes2DPreStep.c.

Function Documentation

◆ NavierStokes2DPreStep()

int NavierStokes2DPreStep ( double *  u,
void *  s,
void *  m,
double  waqt 
)

Pre-step function for the 2D Navier Stokes equations: This function is called at the beginning of each time step.

  • The solution at the beginning of the step is copied into NavierStokes2D::solution for the linearized flux partitioning.
  • The linearized "fast" Jacobian (representing the acoustic modes) NavierStokes2D::fast_jac is computed as:

    \begin{equation} A_f\left({\bf u}\right) = R\left({\bf u}\right)\Lambda_f\left({\bf u}\right)L\left({\bf u}\right) \end{equation}

    where \(R\) and \(L\) are the matrices of right and left eigenvectors, and,

    \begin{equation} \Lambda_f = diag\left[0,0,u+c,u-c \right] \end{equation}

    with \(c\) being the speed of sound.

    Reference:
  • Ghosh, D., Constantinescu, E. M., Semi-Implicit Time Integration of Atmospheric Flows with Characteristic-Based Flux Partitioning, SIAM Journal on Scientific Computing, 38 (3), 2016, A1848-A1875, http://dx.doi.org/10.1137/15M1044369.
Parameters
uSolution vector
sSolver object of type HyPar
mMPI object of type MPIVariables
waqtCurrent simulation time

Definition at line 33 of file NavierStokes2DPreStep.c.

39 {
40  HyPar *solver = (HyPar*) s;
41  NavierStokes2D *param = (NavierStokes2D*) solver->physics;
42  int *dim = solver->dim_local;
43  int ghosts = solver->ghosts, dir, p;
44  double *A;
45  static const int ndims = _MODEL_NDIMS_;
46  static const int JacSize = _MODEL_NVARS_*_MODEL_NVARS_;
47  static int index[_MODEL_NDIMS_], bounds[_MODEL_NDIMS_], offset[_MODEL_NDIMS_];
50 
51  /* set bounds for array index to include ghost points */
52  _ArrayAddCopy1D_(dim,(2*ghosts),bounds,ndims);
53  /* set offset such that index is compatible with ghost point arrangement */
54  _ArraySetValue_(offset,ndims,-ghosts);
55  /* copy the solution to act as a reference for linearization */
57 
58  int done = 0; _ArraySetValue_(index,ndims,0);
59  while (!done) {
60  _ArrayIndex1DWO_(ndims,dim,index,offset,ghosts,p);
61 
62  dir = _XDIR_;
63  A = (param->fast_jac + 2*JacSize*p + dir*JacSize);
64  /* get the eigenvalues, and left & right eigenvectors */
65  _NavierStokes2DEigenvalues_ ((u+_MODEL_NVARS_*p),D,param->gamma,dir);
68  /* remove the entropy modes (corresponding to eigenvalues u) */
69  D[2*_MODEL_NVARS_+2] = D[3*_MODEL_NVARS_+3] = 0.0;
70  /* assemble the Jacobian */
71  MatMult4(_MODEL_NVARS_,DL,D,L );
72  MatMult4(_MODEL_NVARS_,A ,R,DL);
73 
74  dir = _YDIR_;
75  A = (param->fast_jac + 2*JacSize*p + dir*JacSize);
76  /* get the eigenvalues, and left & right eigenvectors */
80  /* remove the entropy modes (corresponding to eigenvalues v) */
81  D[1*_MODEL_NVARS_+1] = D[3*_MODEL_NVARS_+3] = 0.0;
82  /* assemble the Jacobian */
83  MatMult4(_MODEL_NVARS_,DL,D,L );
84  MatMult4(_MODEL_NVARS_,A ,R,DL);
85 
86  _ArrayIncrementIndex_(ndims,bounds,index,done);
87  }
88 
89  return(0);
90 }
#define _NavierStokes2DLeftEigenvectors_(u, L, ga, dir)
#define _NavierStokes2DEigenvalues_(u, D, gamma, dir)
#define _ArrayAddCopy1D_(x, a, y, size)
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
#define _NavierStokes2DRightEigenvectors_(u, R, ga, dir)
#define _YDIR_
Definition: euler2d.h:41
#define _MODEL_NDIMS_
Definition: euler1d.h:56
#define _ArrayIndex1DWO_(N, imax, i, offset, ghost, index)
Structure containing variables and parameters specific to the 2D Navier Stokes equations. This structure contains the physical parameters, variables, and function pointers specific to the 2D Navier-Stokes equations.
#define MatMult4(N, A, X, Y)
#define _ArraySetValue_(x, size, value)
int * dim_local
Definition: hypar.h:37
#define _ArrayIncrementIndex_(N, imax, i, done)
void * physics
Definition: hypar.h:266
#define _XDIR_
Definition: euler1d.h:75
int ghosts
Definition: hypar.h:52
#define _MODEL_NVARS_
Definition: euler1d.h:58
int npoints_local_wghosts
Definition: hypar.h:42
#define _ArrayCopy1D_(x, y, size)