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

Compute the right-hand-side for implicit 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__   "PetscIFunctionImpl"
 

Functions

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

Detailed Description

Compute the right-hand-side for implicit time integration.

Author
Debojyoti Ghosh

Definition in file PetscIFunctionImpl.cpp.

Macro Definition Documentation

#define __FUNCT__   "PetscIFunctionImpl"

Definition at line 16 of file PetscIFunctionImpl.cpp.

Function Documentation

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

Compute the left-hand-side for the implicit time integration of the governing equations: The spatially discretized ODE can be expressed as

\begin{equation} \frac {d{\bf U}} {dt} = {\bf F}\left({\bf U}\right). \end{equation}

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

See Also
TSSetIFunction (https://petsc.org/release/docs/manualpages/TS/TSSetIFunction.html)

Note:

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

Definition at line 37 of file PetscIFunctionImpl.cpp.

43 {
44  PETScContext* context = (PETScContext*) ctxt;
45  SimulationObject* sim = (SimulationObject*) context->simobj;
46  int nsims = context->nsims;
47 
48  PetscFunctionBegin;
49  for (int ns = 0; ns < nsims; ns++) {
50 
51  HyPar* solver = &(sim[ns].solver);
52  MPIVariables* mpi = &(sim[ns].mpi);
53 
54  solver->count_RHSFunction++;
55 
56  int size = solver->npoints_local_wghosts;
57  double* u = solver->u;
58  double* rhs = solver->rhs;
59 
60  /* copy solution from PETSc vector */
61  TransferVecFromPETSc(u,Y,context,ns,context->offsets[ns]);
62  /* apply boundary conditions and exchange data over MPI interfaces */
63  solver->ApplyBoundaryConditions(solver,mpi,u,NULL,t);
65  solver->nvars,
66  solver->dim_local,
67  solver->ghosts,
68  mpi,
69  u );
70 
71  /* Evaluate hyperbolic, parabolic and source terms and the RHS */
72  solver->HyperbolicFunction(solver->hyp,u,solver,mpi,t,1,solver->FFunction,solver->Upwind);
73  solver->ParabolicFunction (solver->par,u,solver,mpi,t);
74  solver->SourceFunction (solver->source,u,solver,mpi,t);
75 
76  _ArraySetValue_(rhs,size*solver->nvars,0.0);
77  _ArrayAXPY_(solver->hyp ,-1.0,rhs,size*solver->nvars);
78  _ArrayAXPY_(solver->par , 1.0,rhs,size*solver->nvars);
79  _ArrayAXPY_(solver->source, 1.0,rhs,size*solver->nvars);
80 
81  /* save a copy of the solution and RHS for use in IJacobian */
82  _ArrayCopy1D_(u ,solver->uref ,(size*solver->nvars));
83  _ArrayCopy1D_(rhs,solver->rhsref,(size*solver->nvars));
84 
85  /* Transfer RHS to PETSc vector */
86  TransferVecToPETSc(rhs,F,context,ns,context->offsets[ns]);
87 
88  }
89 
90  /* LHS = Ydot - F(u) */
91  VecAYPX(F,-1.0,Ydot);
92 
93  PetscFunctionReturn(0);
94 }
double * hyp
Definition: hypar.h:119
int npoints_local_wghosts
Definition: hypar.h:42
#define _ArraySetValue_(x, size, value)
double * u
Definition: hypar.h:116
int TransferVecFromPETSc(double *const, const Vec, void *, const int, const int)
int(* Upwind)(double *, double *, double *, double *, double *, double *, int, void *, double)
Definition: hypar.h:295
int * dim_local
Definition: hypar.h:37
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
double * rhsref
Definition: hypar.h:398
int TransferVecToPETSc(const double *const, Vec, void *, const int, const int)
double * uref
Definition: hypar.h:397
#define _ArrayCopy1D_(x, y, size)
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.
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