#include "petsc.h"#include "petscvec.h"#include "petscmat.h"#include "odaJac.h"Go to the source code of this file.
Defines | |
| #define | JACOBIAN_TYPE1_MULT_BLOCK |
| #define | MASS_MULT_BLOCK |
Functions | |
| PetscErrorCode | ComputeJacobian1 (ot::DA *da, Mat J) |
| PetscErrorCode | CreateActiveJacobian1 (ot::DA *da, Mat *jac) |
| PetscErrorCode | CreateAndComputeFullActiveJacobian1 (ot::DA *da, Mat *J) |
| PetscErrorCode | CreateAndComputeFullJacobian1 (ot::DA *da, Mat *J) |
| PetscErrorCode | CreateAndComputeMassMatrix (ot::DA *da, Mat *jac) |
| PetscErrorCode | CreateJacobian1 (ot::DA *da, Mat *jac) |
| PetscErrorCode | Jacobian1MatDestroy (Mat J) |
| PetscErrorCode | Jacobian1MatGetDiagonal (Mat J, Vec diag) |
| PetscErrorCode | Jacobian1MatMult (Mat J, Vec in, Vec out) |
| PetscErrorCode | Jacobian1ShellMatMult (Mat J, Vec in, Vec out) |
| PetscErrorCode | MassMatDestroy (Mat J) |
| PetscErrorCode | MassMatGetDiagonal (Mat J, Vec diag) |
| PetscErrorCode | MassMatMult (Mat J, Vec in, Vec out) |
Variables | |
| double **** | LaplacianType2Stencil |
| double **** | MassType2Stencil |
Definition in file odaJac.C.
|
|
Value: {\
unsigned int lev = data->da->getLevel(data->da->curr());\
double h = hFac*(1u << (maxD - lev));\
double fac1 = lapFac*h/2.0;\
double fac2 = massFac*h*h*h/8.0;\
unsigned int indices[8];\
data->da->getNodeIndices(indices);\
unsigned char childNum = data->da->getChildNumber();\
unsigned char hnMask = data->da->getHangingNodeIndex(data->da->curr());\
unsigned char elemType = 0;\
GET_ETYPE_BLOCK(elemType,hnMask,childNum)\
for(int k = 0;k < 8;k++) {\
for(int j=0;j<8;j++) {\
outArr[indices[k]] += (((fac1*(LaplacianType2Stencil[childNum][elemType][k][j])) +\
(fac2*(MassType2Stencil[childNum][elemType][k][j])))*inArr[indices[j]]);\
} /*end for j*/\
} /*end for k*/\
}
Definition at line 271 of file odaJac.C. Referenced by Jacobian1MatMult(). |
|
|
Value: {\
unsigned int lev = data->da->getLevel(data->da->curr());\
double h = hFac*(1u << (maxD - lev));\
double fac2 = h*h*h/8.0;\
unsigned int indices[8];\
data->da->getNodeIndices(indices);\
unsigned char childNum = data->da->getChildNumber();\
unsigned char hnMask = data->da->getHangingNodeIndex(data->da->curr());\
unsigned char elemType = 0;\
GET_ETYPE_BLOCK(elemType,hnMask,childNum)\
for(int k = 0;k < 8;k++) {\
for(int j=0;j<8;j++) {\
outArr[indices[k]] += ((fac2*(MassType2Stencil[childNum][elemType][k][j]))*inArr[indices[j]]);\
} /*end for j*/\
} /*end for k*/\
}
Definition at line 378 of file odaJac.C. Referenced by MassMatMult(). |
|
||||||||||||
|
|
|
||||||||||||
|
Definition at line 74 of file odaJac.C. References Jac1MFreeData::da, ot::DA::getCommActive(), ot::DA::getNodeSize(), ot::DA::iAmActive(), iC, Jac1MFreeData::inTmp, Jac1MFreeData::isFinestLevel, Jacobian1MatDestroy(), Jacobian1MatGetDiagonal(), Jacobian1MatMult(), Jac1MFreeData::Jmat_private, and Jac1MFreeData::outTmp. |
|
||||||||||||
|
||||||||||||
|
||||||||||||
|
Definition at line 103 of file odaJac.C. References Jac1MFreeData::da, ot::DA::getComm(), ot::DA::getNodeSize(), iC, Jac1MFreeData::inTmp, Jac1MFreeData::Jmat_private, MassMatDestroy(), MassMatGetDiagonal(), MassMatMult(), and Jac1MFreeData::outTmp. |
|
||||||||||||
|
Definition at line 54 of file odaJac.C. References Jac1MFreeData::da, ot::DA::getComm(), ot::DA::getNodeSize(), iC, Jac1MFreeData::inTmp, Jac1MFreeData::isFinestLevel, Jacobian1MatDestroy(), Jacobian1MatGetDiagonal(), Jacobian1MatMult(), Jac1MFreeData::Jmat_private, and Jac1MFreeData::outTmp. |
|
|
Definition at line 120 of file odaJac.C. References iC, Jac1MFreeData::inTmp, Jac1MFreeData::Jmat_private, and Jac1MFreeData::outTmp. Referenced by CreateActiveJacobian1(), CreateJacobian1(), and main(). |
|
||||||||||||
|
Definition at line 160 of file odaJac.C. References Jac1MFreeData::da, GET_ETYPE_BLOCK, ot::DA::getChildNumber(), ot::DA::getHangingNodeIndex(), ot::DA::getLevel(), ot::DA::getNodeIndices(), iC, LaplacianType2Stencil, MassType2Stencil, ot::DA::vecGetBuffer(), and ot::DA::vecRestoreBuffer(). Referenced by CreateActiveJacobian1(), CreateJacobian1(), and main(). |
|
||||||||||||||||
|
Definition at line 290 of file odaJac.C. References ot::DA::getMaxDepth(), iC, JACOBIAN_TYPE1_MULT_BLOCK, ot::DA::ReadFromGhostsBegin(), ot::DA::ReadFromGhostsEnd(), ot::DA::vecGetBuffer(), and ot::DA::vecRestoreBuffer(). Referenced by CreateActiveJacobian1(), CreateJacobian1(), and main(). |
|
||||||||||||||||
|
Definition at line 22 of file odaJac.C. Referenced by CreateJacobian1(). |
|
|
Definition at line 145 of file odaJac.C. References iC, and Jac1MFreeData::Jmat_private. Referenced by CreateAndComputeMassMatrix(). |
|
||||||||||||
|
Definition at line 225 of file odaJac.C. References GET_ETYPE_BLOCK, ot::DA::getChildNumber(), ot::DA::getHangingNodeIndex(), ot::DA::getLevel(), ot::DA::getNodeIndices(), iC, MassType2Stencil, ot::DA::vecGetBuffer(), and ot::DA::vecRestoreBuffer(). Referenced by CreateAndComputeMassMatrix(). |
|
||||||||||||||||
|
Definition at line 395 of file odaJac.C. References ot::DA::getMaxDepth(), iC, MASS_MULT_BLOCK, ot::DA::ReadFromGhostsBegin(), ot::DA::ReadFromGhostsEnd(), ot::DA::vecGetBuffer(), and ot::DA::vecRestoreBuffer(). Referenced by CreateAndComputeMassMatrix(). |
|
|
Definition at line 43 of file checkError.C. |
|
1.3.9.1