HyPar  1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
PetscIFunctionIMEX.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__ "PetscIFunctionIMEX"
17 
53 PetscErrorCode PetscIFunctionIMEX( TS ts,
54  PetscReal t,
55  Vec Y,
56  Vec Ydot,
57  Vec F,
58  void *ctxt )
59 {
60  PETScContext* context = (PETScContext*) ctxt;
61  SimulationObject* sim = (SimulationObject*) context->simobj;
62  int nsims = context->nsims;
63 
64  PetscFunctionBegin;
65 
66  context->waqt = t;
67 
68  for (int ns = 0; ns < nsims; ns++) {
69 
70  HyPar* solver = &(sim[ns].solver);
71  MPIVariables* mpi = &(sim[ns].mpi);
72 
73  solver->count_IFunction++;
74 
75  int size = solver->npoints_local_wghosts;
76  double *u = solver->u;
77  double *rhs = solver->rhs;
78 
79  /* copy solution from PETSc vector */
80  TransferVecFromPETSc(u,Y,context,ns,context->offsets[ns]);
81  /* apply boundary conditions and exchange data over MPI interfaces */
82  solver->ApplyBoundaryConditions(solver,mpi,u,NULL,t);
84  solver->nvars,
85  solver->dim_local,
86  solver->ghosts,
87  mpi,
88  u );
89 
90  /* initialize right-hand side to zero */
91  _ArraySetValue_(rhs,size*solver->nvars,0.0);
92 
93  /* Evaluate hyperbolic, parabolic and source terms and the RHS */
94  if ((!strcmp(solver->SplitHyperbolicFlux,"yes")) && solver->flag_fdf_specified) {
95  if (context->flag_hyperbolic_f == _IMPLICIT_) {
96  solver->HyperbolicFunction(solver->hyp,u,solver,mpi,t,0,solver->FdFFunction,solver->UpwindFdF);
97  _ArrayAXPY_(solver->hyp,-1.0,rhs,size*solver->nvars);
98  }
99  if (context->flag_hyperbolic_df == _IMPLICIT_) {
100  solver->HyperbolicFunction(solver->hyp,u,solver,mpi,t,0,solver->dFFunction,solver->UpwinddF);
101  _ArrayAXPY_(solver->hyp,-1.0,rhs,size*solver->nvars);
102  }
103  } else if (!strcmp(solver->SplitHyperbolicFlux,"yes")) {
104  if (context->flag_hyperbolic_f == _IMPLICIT_) {
105  solver->HyperbolicFunction(solver->hyp,u,solver,mpi,t,0,solver->FFunction,solver->Upwind);
106  _ArrayAXPY_(solver->hyp,-1.0,rhs,size*solver->nvars);
107  solver->HyperbolicFunction(solver->hyp,u,solver,mpi,t,0,solver->dFFunction,solver->UpwinddF);
108  _ArrayAXPY_(solver->hyp, 1.0,rhs,size*solver->nvars);
109  }
110  if (context->flag_hyperbolic_df == _IMPLICIT_) {
111  solver->HyperbolicFunction(solver->hyp,u,solver,mpi,t,0,solver->dFFunction,solver->UpwinddF);
112  _ArrayAXPY_(solver->hyp,-1.0,rhs,size*solver->nvars);
113  }
114  } else {
115  if (context->flag_hyperbolic == _IMPLICIT_) {
116  solver->HyperbolicFunction(solver->hyp,u,solver,mpi,t,0,solver->FFunction,solver->Upwind);
117  _ArrayAXPY_(solver->hyp,-1.0,rhs,size*solver->nvars);
118  }
119  }
120  if (context->flag_parabolic == _IMPLICIT_) {
121  solver->ParabolicFunction (solver->par,u,solver,mpi,t);
122  _ArrayAXPY_(solver->par, 1.0,rhs,size*solver->nvars);
123  }
124  if (context->flag_source == _IMPLICIT_) {
125  solver->SourceFunction (solver->source,u,solver,mpi,t);
126  _ArrayAXPY_(solver->source, 1.0,rhs,size*solver->nvars);
127  }
128 
129  /* save a copy of the solution and RHS for use in IJacobian */
130  _ArrayCopy1D_(u ,solver->uref ,(size*solver->nvars));
131  _ArrayCopy1D_(rhs,solver->rhsref,(size*solver->nvars));
132 
133  /* Transfer RHS to PETSc vector */
134  TransferVecToPETSc(rhs,F,context,ns,context->offsets[ns]);
135 
136  }
137 
138  /* LHS = Ydot - F(u) */
139  VecAYPX(F,-1.0,Ydot);
140 
141  PetscFunctionReturn(0);
142 }
143 
144 #endif
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.
char SplitHyperbolicFlux[_MAX_STRING_SIZE_]
Definition: hypar.h:92
C++ declarations for MPI-related functions.
#define _IMPLICIT_
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(* FdFFunction)(double *, double *, int, void *, double)
Definition: hypar.h:286
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
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)
int(* UpwindFdF)(double *, double *, double *, double *, double *, double *, int, void *, double)
Definition: hypar.h:307
int flag_fdf_specified
Definition: hypar.h:290
Structure containing the variables for time-integration with PETSc.
#define _ArraySetValue_(x, size, value)
int * dim_local
Definition: hypar.h:37
double * rhsref
Definition: hypar.h:398
int(* dFFunction)(double *, double *, int, void *, double)
Definition: hypar.h:280
PetscErrorCode PetscIFunctionIMEX(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *ctxt)
int(* UpwinddF)(double *, double *, double *, double *, double *, double *, int, void *, double)
Definition: hypar.h:300
int ghosts
Definition: hypar.h:52
Structure of MPI-related variables.
int npoints_local_wghosts
Definition: hypar.h:42
int count_IFunction
Definition: hypar.h:422
double * uref
Definition: hypar.h:397
#define _ArrayCopy1D_(x, y, size)
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