HyPar  1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
PetscIJacobian.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__ "PetscIJacobian"
17 
52 PetscErrorCode PetscIJacobian( TS ts,
53  PetscReal t,
54  Vec Y,
55  Vec Ydot,
56  PetscReal a,
57  Mat A,
58  Mat B,
59  void *ctxt )
60 {
61  PETScContext* context = (PETScContext*) ctxt;
62 
63  PetscFunctionBegin;
64  for (int ns = 0; ns < context->nsims; ns++) {
65  ((SimulationObject*)context->simobj)[ns].solver.count_IJacobian++;
66  }
67  context->shift = a;
68  context->waqt = t;
69  /* Construct preconditioning matrix */
70  if (context->flag_use_precon) PetscComputePreconMatImpl(B,Y,context);
71 
72  PetscFunctionReturn(0);
73 }
74 
75 #undef __FUNCT__
76 #define __FUNCT__ "PetscJacobianFunction_JFNK"
77 
111 PetscErrorCode PetscJacobianFunction_JFNK( Mat Jacobian,
112  Vec Y,
113  Vec F )
115 {
116  PETScContext* context(nullptr);
117 
118  PetscFunctionBegin;
119 
120  MatShellGetContext(Jacobian,&context);
121  SimulationObject* sim = (SimulationObject*) context->simobj;
122  int nsims = context->nsims;
123 
124  double normY;
125  VecNorm(Y,NORM_2,&normY);
126 
127  if (normY < 1e-16) {
128 
129  /* F = 0 */
130  VecZeroEntries(F);
131  /* [J]Y = aY - F(Y) */
132  VecAXPBY(F,context->shift,0,Y);
133 
134  } else {
135 
136  double epsilon = context->jfnk_eps / normY;
137  double t = context->waqt; /* current stage/step time */
138 
139  for (int ns = 0; ns < nsims; ns++) {
140 
141  HyPar* solver = &(sim[ns].solver);
142  MPIVariables* mpi = &(sim[ns].mpi);
143  solver->count_IJacFunction++;
144 
145  int size = solver->npoints_local_wghosts;
146 
147  double *u = solver->u;
148  double *uref = solver->uref;
149  double *rhsref = solver->rhsref;
150  double *rhs = solver->rhs;
151 
152  /* copy solution from PETSc vector */
153  TransferVecFromPETSc(u,Y,context,ns,context->offsets[ns]);
154  _ArrayAYPX_(uref,epsilon,u,size*solver->nvars);
155  /* apply boundary conditions and exchange data over MPI interfaces */
156  solver->ApplyBoundaryConditions(solver,mpi,u,NULL,t);
158  solver->nvars,
159  solver->dim_local,
160  solver->ghosts,
161  mpi,
162  u );
163 
164  /* Evaluate hyperbolic, parabolic and source terms and the RHS for U+dU */
165  _ArraySetValue_(rhs,size*solver->nvars,0.0);
166  solver->HyperbolicFunction( solver->hyp,u,solver,mpi,t,0,
167  solver->FFunction,solver->Upwind);
168  _ArrayAXPY_(solver->hyp,-1.0,rhs,size*solver->nvars);
169  solver->ParabolicFunction (solver->par,u,solver,mpi,t);
170  _ArrayAXPY_(solver->par, 1.0,rhs,size*solver->nvars);
171  solver->SourceFunction (solver->source,u,solver,mpi,t);
172  _ArrayAXPY_(solver->source, 1.0,rhs,size*solver->nvars);
173 
174  _ArrayAXPY_(rhsref,-1.0,rhs,size*solver->nvars);
175  /* Transfer RHS to PETSc vector */
176  TransferVecToPETSc(rhs,F,context,ns,context->offsets[ns]);
177  }
178 
179  /* [J]Y = aY - F(Y) */
180  VecAXPBY(F,context->shift,(-1.0/epsilon),Y);
181 
182  }
183 
184  PetscFunctionReturn(0);
185 }
186 
187 #undef __FUNCT__
188 #define __FUNCT__ "PetscJacobianFunction_Linear"
189 
236 PetscErrorCode PetscJacobianFunction_Linear( Mat Jacobian,
237  Vec Y,
238  Vec F )
240 {
241  PETScContext* context(nullptr);
242 
243  PetscFunctionBegin;
244 
245  MatShellGetContext(Jacobian,&context);
246  SimulationObject* sim = (SimulationObject*) context->simobj;
247  int nsims = context->nsims;
248 
249  double normY;
250  VecNorm(Y,NORM_2,&normY);
251 
252  if (normY < 1e-16) {
253 
254  /* F = 0 */
255  VecZeroEntries(F);
256  /* [J]Y = aY - F(Y) */
257  VecAXPBY(F,context->shift,0,Y);
258 
259  } else {
260 
261  double t = context->waqt; /* current stage/step time */
262 
263  for (int ns = 0; ns < nsims; ns++) {
264 
265  HyPar* solver = &(sim[ns].solver);
266  MPIVariables* mpi = &(sim[ns].mpi);
267  solver->count_IJacFunction++;
268 
269  int size = solver->npoints_local_wghosts;
270 
271  double *u = solver->u;
272  double *uref = solver->uref;
273  double *rhsref = solver->rhsref;
274  double *rhs = solver->rhs;
275 
276  /* copy solution from PETSc vector */
277  TransferVecFromPETSc(u,Y,context,ns,context->offsets[ns]);
278  _ArrayAYPX_(uref,1.0,u,size*solver->nvars);
279  /* apply boundary conditions and exchange data over MPI interfaces */
280  solver->ApplyBoundaryConditions(solver,mpi,u,NULL,t);
282  solver->nvars,
283  solver->dim_local,
284  solver->ghosts,
285  mpi,
286  u );
287 
288  /* Evaluate hyperbolic, parabolic and source terms and the RHS for U+dU */
289  _ArraySetValue_(rhs,size*solver->nvars,0.0);
290  solver->HyperbolicFunction( solver->hyp,u,solver,mpi,t,0,
291  solver->FFunction,solver->Upwind);
292  _ArrayAXPY_(solver->hyp,-1.0,rhs,size*solver->nvars);
293  solver->ParabolicFunction (solver->par,u,solver,mpi,t);
294  _ArrayAXPY_(solver->par, 1.0,rhs,size*solver->nvars);
295  solver->SourceFunction (solver->source,u,solver,mpi,t);
296  _ArrayAXPY_(solver->source, 1.0,rhs,size*solver->nvars);
297 
298  _ArrayAXPY_(rhsref,-1.0,rhs,size*solver->nvars);
299  /* Transfer RHS to PETSc vector */
300  TransferVecToPETSc(rhs,F,context,ns,context->offsets[ns]);
301  }
302 
303  /* [J]Y = aY - F(Y) */
304  VecAXPBY(F,context->shift,-1.0,Y);
305 
306  }
307 
308  PetscFunctionReturn(0);
309 }
310 
311 #endif
int TransferVecToPETSc(const double *const, Vec, void *, const int, const int)
int nvars
Definition: hypar.h:29
int count_IJacFunction
Definition: hypar.h:422
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.
#define _ArrayAYPX_(x, a, y, size)
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
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 PetscComputePreconMatImpl(Mat, Vec, void *)
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)
PetscErrorCode PetscJacobianFunction_JFNK(Mat Jacobian, Vec Y, Vec F)
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
PetscErrorCode PetscIJacobian(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *ctxt)
int ghosts
Definition: hypar.h:52
Structure of MPI-related variables.
int npoints_local_wghosts
Definition: hypar.h:42
double * uref
Definition: hypar.h:397
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
PetscErrorCode PetscJacobianFunction_Linear(Mat Jacobian, Vec Y, Vec F)