HyPar  1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
PetscRHSFunctionExpl.cpp
Go to the documentation of this file.
1 
6 #ifdef with_petsc
7 
8 #include <stdlib.h>
9 #include <basic.h>
10 #include <arrayfunctions.h>
11 #include <mpivars_cpp.h>
12 #include <simulation_object.h>
13 #include <petscinterface.h>
14 
15 #undef __FUNCT__
16 #define __FUNCT__ "PetscRHSFunctionExpl"
17 
36 PetscErrorCode PetscRHSFunctionExpl( TS ts,
37  PetscReal t,
38  Vec Y,
39  Vec F,
40  void *ctxt )
41 {
42  PETScContext* context = (PETScContext*) ctxt;
43  SimulationObject* sim = (SimulationObject*) context->simobj;
44  int nsims = context->nsims;
45 
46  PetscFunctionBegin;
47 
48  for (int ns = 0; ns < nsims; ns++) {
49 
50  HyPar* solver = &(sim[ns].solver);
51  MPIVariables* mpi = &(sim[ns].mpi);
52 
53  solver->count_RHSFunction++;
54 
55  int size = solver->npoints_local_wghosts;
56  double* u = solver->u;
57  double* rhs = solver->rhs;
58 
59  /* copy solution from PETSc vector */
60  TransferVecFromPETSc(u,Y,context,ns,context->offsets[ns]);
61  /* apply boundary conditions and exchange data over MPI interfaces */
62  solver->ApplyBoundaryConditions(solver,mpi,u,NULL,t);
64  solver->nvars,
65  solver->dim_local,
66  solver->ghosts,
67  mpi,
68  u );
69 
70  /* Evaluate hyperbolic, parabolic and source terms and the RHS */
71  solver->HyperbolicFunction(solver->hyp,u,solver,mpi,t,1,solver->FFunction,solver->Upwind);
72 
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  /* Transfer RHS to PETSc vector */
82  TransferVecToPETSc(rhs,F,context,ns,context->offsets[ns]);
83 
84  }
85 
86  PetscFunctionReturn(0);
87 }
88 
89 #endif
int count_RHSFunction
Definition: hypar.h:422
int TransferVecToPETSc(const double *const, Vec, void *, const int, const int)
int nvars
Definition: hypar.h:29
int(* ParabolicFunction)(double *, double *, void *, void *, double)
Definition: hypar.h:256
Structure defining a simulation.
C++ declarations for MPI-related functions.
double * par
Definition: hypar.h:122
int MPIExchangeBoundariesnD(int, int, int *, int, void *, double *)
int(* ApplyBoundaryConditions)(void *, void *, double *, double *, double)
Definition: hypar.h:214
double * u
Definition: hypar.h:116
Some basic definitions and macros.
Simulation object.
int(* Upwind)(double *, double *, double *, double *, double *, double *, int, void *, double)
Definition: hypar.h:295
int(* SourceFunction)(double *, double *, void *, void *, double)
Definition: hypar.h:259
int ndims
Definition: hypar.h:26
PetscErrorCode PetscRHSFunctionExpl(TS ts, PetscReal t, Vec Y, Vec F, void *ctxt)
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(* FFunction)(double *, double *, int, void *, double)
Definition: hypar.h:276
double * rhs
Definition: hypar.h:399
Structure containing all solver-specific variables and functions.
Definition: hypar.h:23
#define _ArrayAXPY_(x, a, y, size)
Structure containing the variables for time-integration with PETSc.
#define _ArraySetValue_(x, size, value)
int * dim_local
Definition: hypar.h:37
int ghosts
Definition: hypar.h:52
Structure of MPI-related variables.
int npoints_local_wghosts
Definition: hypar.h:42
Contains macros and function definitions for common array operations.
int TransferVecFromPETSc(double *const, const Vec, void *, const int, const int)
double * source
Definition: hypar.h:125
double * hyp
Definition: hypar.h:119