HyPar
1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
|
Compute the implicitly-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__ "PetscIFunctionIMEX" |
Functions | |
PetscErrorCode | PetscIFunctionIMEX (TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *ctxt) |
Compute the implicitly-treated part of the right-hand-side for IMEX time integration.
Definition in file PetscIFunctionIMEX.cpp.
#define __FUNCT__ "PetscIFunctionIMEX" |
Definition at line 16 of file PetscIFunctionIMEX.cpp.
PetscErrorCode PetscIFunctionIMEX | ( | TS | ts, |
PetscReal | t, | ||
Vec | Y, | ||
Vec | Ydot, | ||
Vec | F, | ||
void * | ctxt | ||
) |
Compute the implicitly-treated part of the right-hand-side 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 G}\left({\bf U}\right)\) represents all terms that the user has indicated to be integrated in time implicitly (PETScContext::flag_hyperbolic_f, PETScContext::flag_hyperbolic_df, PETScContext::flag_hyperbolic, PETScContext::flag_parabolic, and PETScContext::flag_source).
This function computes the left-hand-side of the above equation:
\begin{equation} \mathcal{G}\left(\dot{\bf U},{\bf U},t\right) = \dot{\bf U} - {\bf G}\left({\bf U}\right) \end{equation}
given \(\dot{\bf U}\) and \({\bf U}\).
Notes:
ts | The time integration object |
t | Current solution time |
Y | State vector (input) |
Ydot | Time derivative of the state vector (input) |
F | The computed function vector |
ctxt | Object of type PETScContext |
Definition at line 53 of file PetscIFunctionIMEX.cpp.