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

Compute the explicitly-treated part of the right-hand side for IMEX time integration. More...

#include <stdlib.h>
#include <basic.h>
#include <arrayfunctions.h>
#include <mpivars_cpp.h>
#include <simulation_object.h>
#include <petscinterface.h>

Go to the source code of this file.

Macros

#define __FUNCT__   "PetscRHSFunctionIMEX"
 

Functions

PetscErrorCode PetscRHSFunctionIMEX (TS ts, PetscReal t, Vec Y, Vec F, void *ctxt)
 

Detailed Description

Compute the explicitly-treated part of the right-hand side for IMEX time integration.

Author
Debojyoti Ghosh

Definition in file PetscRHSFunctionIMEX.cpp.

Macro Definition Documentation

#define __FUNCT__   "PetscRHSFunctionIMEX"

Definition at line 16 of file PetscRHSFunctionIMEX.cpp.

Function Documentation

PetscErrorCode PetscRHSFunctionIMEX ( TS  ts,
PetscReal  t,
Vec  Y,
Vec  F,
void *  ctxt 
)

Compute the explicitly-treated right-hand-side (RHS) for the implicit-explicit (IMEX) time integration of the governing equations: The ODE, obtained after discretizing the governing PDE in space, is expressed as follows (for the purpose of IMEX time integration):

\begin{eqnarray} \frac {d{\bf U}}{dt} &=& {\bf F}\left({\bf U}\right) + {\bf G}\left({\bf U}\right), \\ \Rightarrow \dot{\bf U} - {\bf G}\left({\bf U}\right) &=& {\bf F}\left({\bf U}\right), \end{eqnarray}

where \({\bf F}\) is non-stiff and integrated in time explicitly, and \({\bf G}\) is stiff and integrated in time implicitly, and \({\bf U}\) represents the entire solution vector (state vector).

Note that \({\bf F}\left({\bf U}\right)\) represents all terms that the user has indicated to be integrated in time explicitly (PETScContext::flag_hyperbolic_f, PETScContext::flag_hyperbolic_df, PETScContext::flag_hyperbolic, PETScContext::flag_parabolic, and PETScContext::flag_source).

This function computes \({\bf F}\left({\bf U}\right)\), given \({\bf U}\).

See Also
PetscIFunctionIMEX(), TSSetRHSFunction (https://petsc.org/release/docs/manualpages/TS/TSSetRHSFunction.html)

Note:

Parameters
tsTime integration object
tCurrent simulation time
YState vector (input)
FThe computed right-hand-side vector
ctxtObject of type PETScContext

Definition at line 49 of file PetscRHSFunctionIMEX.cpp.

54 {
55  PETScContext* context = (PETScContext*) ctxt;
56  SimulationObject* sim = (SimulationObject*) context->simobj;
57  int nsims = context->nsims;
58 
59  context->waqt = t;
60 
61  PetscFunctionBegin;
62  for (int ns = 0; ns < nsims; ns++) {
63 
64  HyPar* solver = &(sim[ns].solver);
65  MPIVariables* mpi = &(sim[ns].mpi);
66 
67  solver->count_RHSFunction++;
68 
69  int size = solver->npoints_local_wghosts;
70  double *u = solver->u;
71  double *rhs = solver->rhs;
72 
73  /* copy solution from PETSc vector */
74  TransferVecFromPETSc(u,Y,context,ns,context->offsets[ns]);
75  /* apply boundary conditions and exchange data over MPI interfaces */
76  solver->ApplyBoundaryConditions(solver,mpi,u,NULL,t);
77  MPIExchangeBoundariesnD(solver->ndims,solver->nvars,solver->dim_local,
78  solver->ghosts,mpi,u);
79 
80  /* initialize right-hand side to zero */
81  _ArraySetValue_(rhs,size*solver->nvars,0.0);
82 
83  /* Evaluate hyperbolic, parabolic and source terms and the RHS */
84  if ((!strcmp(solver->SplitHyperbolicFlux,"yes")) && solver->flag_fdf_specified) {
85  if (context->flag_hyperbolic_f == _EXPLICIT_) {
86  solver->HyperbolicFunction(solver->hyp,u,solver,mpi,t,0,solver->FdFFunction,solver->UpwindFdF);
87  _ArrayAXPY_(solver->hyp,-1.0,rhs,size*solver->nvars);
88  }
89  if (context->flag_hyperbolic_df == _EXPLICIT_) {
90  solver->HyperbolicFunction(solver->hyp,u,solver,mpi,t,0,solver->dFFunction,solver->UpwinddF);
91  _ArrayAXPY_(solver->hyp,-1.0,rhs,size*solver->nvars);
92  }
93  } else if (!strcmp(solver->SplitHyperbolicFlux,"yes")) {
94  if (context->flag_hyperbolic_f == _EXPLICIT_) {
95  solver->HyperbolicFunction(solver->hyp,u,solver,mpi,t,0,solver->FFunction,solver->Upwind);
96  _ArrayAXPY_(solver->hyp,-1.0,rhs,size*solver->nvars);
97  solver->HyperbolicFunction(solver->hyp,u,solver,mpi,t,0,solver->dFFunction,solver->UpwinddF);
98  _ArrayAXPY_(solver->hyp, 1.0,rhs,size*solver->nvars);
99  }
100  if (context->flag_hyperbolic_df == _EXPLICIT_) {
101  solver->HyperbolicFunction(solver->hyp,u,solver,mpi,t,0,solver->dFFunction,solver->UpwinddF);
102  _ArrayAXPY_(solver->hyp,-1.0,rhs,size*solver->nvars);
103  }
104  } else {
105  if (context->flag_hyperbolic == _EXPLICIT_) {
106  solver->HyperbolicFunction(solver->hyp,u,solver,mpi,t,0,solver->FFunction,solver->Upwind);
107  _ArrayAXPY_(solver->hyp,-1.0,rhs,size*solver->nvars);
108  }
109  }
110  if (context->flag_parabolic == _EXPLICIT_) {
111  solver->ParabolicFunction (solver->par,u,solver,mpi,t);
112  _ArrayAXPY_(solver->par, 1.0,rhs,size*solver->nvars);
113  }
114  if (context->flag_source == _EXPLICIT_) {
115  solver->SourceFunction (solver->source,u,solver,mpi,t);
116  _ArrayAXPY_(solver->source, 1.0,rhs,size*solver->nvars);
117  }
118 
119  /* Transfer RHS to PETSc vector */
120  TransferVecToPETSc(rhs,F,context,ns,context->offsets[ns]);
121  }
122 
123  PetscFunctionReturn(0);
124 }
double * hyp
Definition: hypar.h:119
int npoints_local_wghosts
Definition: hypar.h:42
#define _ArraySetValue_(x, size, value)
int(* dFFunction)(double *, double *, int, void *, double)
Definition: hypar.h:280
double * u
Definition: hypar.h:116
int TransferVecFromPETSc(double *const, const Vec, void *, const int, const int)
#define _EXPLICIT_
int(* Upwind)(double *, double *, double *, double *, double *, double *, int, void *, double)
Definition: hypar.h:295
int flag_fdf_specified
Definition: hypar.h:290
int(* UpwinddF)(double *, double *, double *, double *, double *, double *, int, void *, double)
Definition: hypar.h:300
int * dim_local
Definition: hypar.h:37
int(* UpwindFdF)(double *, double *, double *, double *, double *, double *, int, void *, double)
Definition: hypar.h:307
int(* FdFFunction)(double *, double *, int, void *, double)
Definition: hypar.h:286
int(* FFunction)(double *, double *, int, void *, double)
Definition: hypar.h:276
int(* ParabolicFunction)(double *, double *, void *, void *, double)
Definition: hypar.h:256
int ghosts
Definition: hypar.h:52
int count_RHSFunction
Definition: hypar.h:422
double * par
Definition: hypar.h:122
double * source
Definition: hypar.h:125
int(* SourceFunction)(double *, double *, void *, void *, double)
Definition: hypar.h:259
int TransferVecToPETSc(const double *const, Vec, void *, const int, const int)
int(* HyperbolicFunction)(double *, double *, void *, void *, double, int, int(*)(double *, double *, int, void *, double), int(*)(double *, double *, double *, double *, double *, double *, int, void *, double))
Definition: hypar.h:250
int nvars
Definition: hypar.h:29
Structure containing the variables for time-integration with PETSc.
Structure defining a simulation.
char SplitHyperbolicFlux[_MAX_STRING_SIZE_]
Definition: hypar.h:92
int ndims
Definition: hypar.h:26
int(* ApplyBoundaryConditions)(void *, void *, double *, double *, double)
Definition: hypar.h:214
int MPIExchangeBoundariesnD(int, int, int *, int, void *, double *)
Structure of MPI-related variables.
#define _ArrayAXPY_(x, a, y, size)
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
double * rhs
Definition: hypar.h:399