HyPar  1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
TransferToPETSc.cpp File Reference

Copy from a HyPar array to a PETSc vector. More...

#include <stdlib.h>
#include <vector>
#include <basic.h>
#include <arrayfunctions.h>
#include <simulation_object.h>
#include <bandedmatrix.h>
#include <petscinterface_struct.h>

Go to the source code of this file.

Macros

#define __FUNCT__   "TransferVecToPETSc"
 

Functions

int TransferVecToPETSc (const double *const u, Vec Y, void *ctxt, const int sim_idx, const int offset)
 
int TransferMatToPETSc (void *J, Mat A, void *ctxt)
 

Detailed Description

Copy from a HyPar array to a PETSc vector.

Author
Debojyoti Ghosh

Definition in file TransferToPETSc.cpp.

Macro Definition Documentation

◆ __FUNCT__

#define __FUNCT__   "TransferVecToPETSc"

Definition at line 17 of file TransferToPETSc.cpp.

Function Documentation

◆ TransferVecToPETSc()

int TransferVecToPETSc ( const double *const  u,
Vec  Y,
void *  ctxt,
const int  sim_idx,
const int  offset 
)

Copy data to a PETSc vector (used by PETSc time integrators, and with no ghost points) from a HyPar::u array (with ghost points).

See also
TransferVecFromPETSc()
Parameters
uHyPar::u type array (with ghost points)
YPETSc vector
ctxtObject of type PETScContext
sim_idxSimulation object index
offsetOffset

Definition at line 24 of file TransferToPETSc.cpp.

29 {
30  PETScContext* context = (PETScContext*) ctxt;
31  SimulationObject* sim = (SimulationObject*) context->simobj;
32  double* Yarr;
33 
34  PetscFunctionBegin;
35  VecGetArray(Y,&Yarr);
36  std::vector<int> index(sim[sim_idx].solver.ndims,0);
37  ArrayCopynD( sim[sim_idx].solver.ndims,
38  u,
39  (Yarr+offset),
40  sim[sim_idx].solver.dim_local,
41  sim[sim_idx].solver.ghosts,
42  0,
43  index.data(),
44  sim[sim_idx].solver.nvars );
45  VecRestoreArray(Y,&Yarr);
46 
47  PetscFunctionReturn(0);
48 }
int nvars
Definition: hypar.h:29
Structure defining a simulation.
INLINE int ArrayCopynD(int, const double *, double *, int *, int, int, int *, int)
Structure containing the variables for time-integration with PETSc.

◆ TransferMatToPETSc()

int TransferMatToPETSc ( void *  J,
Mat  A,
void *  ctxt 
)

Copy a matrix of type BandedMatrix to a PETSc matrix.

Parameters
JMatrix of type BandedMatrix
APETSc matrix
ctxtObject of type PETScContext

Definition at line 53 of file TransferToPETSc.cpp.

56 {
57  BandedMatrix *M = (BandedMatrix*) J;
58  PetscErrorCode ierr = 0;
59  int bs = M->BlockSize, nbands = M->nbands, bs2 = bs*bs;
60 
61  for (int i=0; i<M->nrows_local; i++) {
62  int colind[nbands];
63  double val[bs][bs*nbands];
64  for (int n=0; n<nbands; n++) {
65  colind[n] = M->ncol[nbands*i+n];
66  for (int p=0; p<bs; p++) {
67  for (int q = 0; q<bs; q++) {
68  val[p][n*bs+q] = M->data[i*nbands*bs2+n*bs2+p*bs+q];
69  }
70  }
71  }
72  MatSetValuesBlocked(A,1,&M->nrow[i],M->nbands,&colind[0],&val[0][0],INSERT_VALUES);
73  }
74 
75  MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
76  MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY);
77 
78  return(0);
79 }
double * data
Definition: bandedmatrix.h:24
Structure for defining a banded block matrix.
Definition: bandedmatrix.h:18