#include "dtypes.h"
#include "petscpc.h"
#include "petscmg.h"
#include "petscmat.h"
#include "private/pcimpl.h"
#include "omg.h"
#include "oda.h"
#include "odaUtils.h"
#include "parUtils.h"
#include <iostream>
#include "dendro.h"
Go to the source code of this file.
Namespaces | |
namespace | ot |
Private Functions | |
PetscErrorCode | createInterpolationType1 (DAMG damgc, DAMG damgf, Mat *R) |
PetscErrorCode | createInterpolationType2 (DAMG damgc, DAMG damgf, Mat *R) |
PetscErrorCode | DAMGSetUp (DAMG *damg) |
PetscErrorCode | DAMGSetUpLevel (DAMG *damg, KSP ksp, int nlevels) |
PetscErrorCode | DAMGSolveKSP (DAMG *damg, int level) |
PetscErrorCode | rpMatDestroy (Mat R) |
int | scatterValues (Vec in, Vec out, PetscInt inSz, PetscInt outSz, int *&sendSz, int *&sendOff, int *&recvSz, int *&recvOff, MPI_Comm comm) |
Redistributes the vector, preserving the relative ordering of the elements. | |
Global Function Handle for using KSP_Shell | |
This will be used at the coarsest grid if not all processors are active on the coarsest grid | |
void(* | getPrivateMatricesForKSP_Shell )(Mat mat, Mat *AmatPrivate, Mat *PmatPrivate, MatStructure *pFlag) = NULL |
Variables for storing the various stencils used in the oda and omg module | |
double ***** | RmatType1Stencil = NULL |
double **** | RmatType2Stencil = NULL |
unsigned short **** | VtxMap1 = NULL |
unsigned short ***** | VtxMap2 = NULL |
unsigned short ***** | VtxMap3 = NULL |
unsigned short ****** | VtxMap4 = NULL |
Defines | |
#define | iC(fun) {CHKERRQ(fun);} |
Functions | |
PetscErrorCode | DAMG_Finalize () |
Destroys the stencils used in R and P. | |
PetscErrorCode | DAMG_Initialize (MPI_Comm comm) |
Initializes the stencils used in R and P. | |
void | DAMG_InitPrivateType1 (MPI_Comm comm) |
Loads the stencils used in R and P. Processor 0 reads the files and broadcasts to others. | |
void | DAMG_InitPrivateType2 (MPI_Comm comm) |
Loads the stencils used in R and P. The comm is split into many groups (each group containing a maximum of 1000 processors) and the first processor in each group reads the file and broadcasts the stencils to the other processors in the group. | |
void | DAMG_InitPrivateType3 (MPI_Comm comm) |
Loads the stencils used in R and P. No communication. Each processor opens a unique file. | |
PetscErrorCode | DAMGCreateAndSetDA (MPI_Comm comm, int &nlevels, void *user, DAMG **damg, std::vector< ot::TreeNode > &finestOctree, unsigned int dof, double loadFac, bool compressLut, bool incCorner) |
Constructs the Multigrid object. | |
PetscErrorCode | DAMGCreateJMatrix (DAMG damg, PetscErrorCode(*crjac)(DAMG, Mat *)) |
int | DAMGCreateSuppressedDOFs (DAMG *damg) |
Call this function to allocate memory for vectors used to mark the dirichlet nodes. | |
PetscErrorCode | DAMGDestroy (DAMG *damg) |
Destroy the DAMG object. | |
PetscErrorCode | DAMGInitialGuessCurrent (DAMG damg, Vec vec) |
Use the current vector as the initial guess. | |
PetscErrorCode | DAMGSetInitialGuess (DAMG *damg, PetscErrorCode(*guess)(DAMG, Vec)) |
Set a function handle that will be used to generate the initial guess. | |
PetscErrorCode | DAMGSetKSP (DAMG *damg, PetscErrorCode(*crjac)(DAMG, Mat *), PetscErrorCode(*compJac)(DAMG, Mat, Mat), PetscErrorCode(*rhs)(DAMG, Vec)) |
PetscErrorCode | DAMGSetNullSpace (DAMG *damg, PetscTruth has_cnst, PetscInt n, PetscErrorCode(*func)(DAMG, Vec[])) |
PetscErrorCode | DAMGSolve (DAMG *damg) |
Solves the problem. | |
PetscErrorCode | PC_KSP_Shell_Apply (void *ctx, Vec rhs, Vec sol) |
PetscErrorCode | PC_KSP_Shell_Destroy (void *ctx) |
PetscErrorCode | PC_KSP_Shell_SetUp (void *ctx) |
void | PrintDAMG (DAMG *damg) |
Prints detailed information about the meshes for each level. |
Definition in file omg.C.
|
|
|
|
Destroys the stencils used in R and P.
Definition at line 664 of file omg.C. References ot::DA_Finalize(), ot::destroyRmatType1Stencil(), ot::destroyRmatType2Stencil(), ot::destroyVtxMaps(), PROF_MG_FINAL_BEGIN, PROF_MG_FINAL_END, ot::RmatType1Stencil, ot::RmatType2Stencil, ot::VtxMap1, ot::VtxMap2, ot::VtxMap3, and ot::VtxMap4. Referenced by main(). |
|
Initializes the stencils used in R and P.
Definition at line 62 of file omg.C. References ot::DA_Initialize(), ot::DAMG_InitPrivateType1(), ot::DAMG_InitPrivateType2(), ot::DAMG_InitPrivateType3(), PROF_MG_INIT_BEGIN, and PROF_MG_INIT_END. Referenced by main(). |
|
Loads the stencils used in R and P. Processor 0 reads the files and broadcasts to others.
Definition at line 391 of file omg.C. References ot::readRmatType1Stencil(), ot::readRmatType2Stencil(), ot::readVtxMaps(), ot::RmatType1Stencil, ot::RmatType2Stencil, ot::VtxMap1, ot::VtxMap2, ot::VtxMap3, and ot::VtxMap4. Referenced by ot::DAMG_Initialize(). |
|
Loads the stencils used in R and P. The comm is split into many groups (each group containing a maximum of 1000 processors) and the first processor in each group reads the file and broadcasts the stencils to the other processors in the group.
Definition at line 91 of file omg.C. References ot::IreadRmatType1Stencil(), ot::IreadRmatType2Stencil(), ot::IreadVtxMaps(), ot::RmatType1Stencil, ot::RmatType2Stencil, par::splitComm2way(), ot::VtxMap1, ot::VtxMap2, ot::VtxMap3, and ot::VtxMap4. Referenced by ot::DAMG_Initialize(). |
|
Loads the stencils used in R and P. No communication. Each processor opens a unique file.
Definition at line 80 of file omg.C. References ot::IreadRmatType1Stencil(), ot::IreadRmatType2Stencil(), ot::IreadVtxMaps(), ot::RmatType1Stencil, ot::RmatType2Stencil, ot::VtxMap1, ot::VtxMap2, ot::VtxMap3, and ot::VtxMap4. Referenced by ot::DAMG_Initialize(). |
|
|
Definition at line 755 of file omg.C. References ot::DAMG. Referenced by main(). |
|
Call this function to allocate memory for vectors used to mark the dirichlet nodes.
Definition at line 43 of file omg.C. References ot::_p_DAMG::da, ot::_p_DAMG::da_aux, ot::DAMG, ot::_p_DAMG::dof, ot::DA::getLocalBufferSize(), ot::_p_DAMG::nlevels, ot::_p_DAMG::suppressedDOF, and ot::_p_DAMG::suppressedDOFaux. Referenced by main(). |
|
Destroy the DAMG object.
Definition at line 676 of file omg.C. References ot::_p_DAMG::B, ot::_p_DAMG::b, ot::_p_DAMG::da, ot::_p_DAMG::da_aux, ot::DAMG, ot::_p_DAMG::J, ot::_p_DAMG::ksp, ot::_p_DAMG::nlevels, ot::_p_DAMG::r, ot::_p_DAMG::R, ot::_p_DAMG::suppressedDOF, ot::_p_DAMG::suppressedDOFaux, and ot::_p_DAMG::x. Referenced by main(). |
|
Use the current vector as the initial guess.
|
|
Set a function handle that will be used to generate the initial guess.
Definition at line 1094 of file omg.C. References ot::_p_DAMG::initialguess, and ot::_p_DAMG::nlevels. Referenced by main(). |
|
Definition at line 770 of file omg.C. References ot::_p_DAMG::B, ot::PC_KSP_Shell::commActive, ot::_p_DAMG::da, ot::DAMG, ot::DAMGSetUpLevel(), ot::DA::getCommActive(), ot::DA::getNpesActive(), ot::DA::getNpesAll(), ot::DA::iAmActive(), ot::PC_KSP_Shell::iAmActive, ot::_p_DAMG::J, ot::_p_DAMG::ksp, ot::PC_KSP_Shell::ksp_private, ot::_p_DAMG::nlevels, ot::PC_KSP_Shell::pc, ot::PC_KSP_Shell_Apply(), ot::PC_KSP_Shell_Destroy(), ot::PC_KSP_Shell_SetUp(), PROF_MG_SET_KSP_BEGIN, PROF_MG_SET_KSP_END, ot::_p_DAMG::rhs, ot::PC_KSP_Shell::rhs_private, ot::PC_KSP_Shell::sol_private, and ot::_p_DAMG::solve. Referenced by main(), solve_neumann(), and solve_neumann_oct(). |
|
Definition at line 1038 of file omg.C. References ot::_p_DAMG::ksp, and ot::_p_DAMG::nlevels. |
|
Definition at line 2145 of file omg.C. References ot::createInterpolationType1(), ot::createInterpolationType2(), ot::DA::createVector(), ot::_p_DAMG::da, ot::_p_DAMG::da_aux, ot::_p_DAMG::nlevels, PROF_MG_SETUP_BEGIN, and PROF_MG_SETUP_END. Referenced by ot::DAMGCreateAndSetDA(). |
|
Definition at line 1156 of file omg.C. References ot::_p_DAMG::comm. Referenced by ot::DAMGSetKSP(). |
|
Solves the problem.
Definition at line 1106 of file omg.C. References DAMGGetFine, ot::_p_DAMG::initialguess, ot::_p_DAMG::nlevels, ot::_p_DAMG::solve, and ot::_p_DAMG::x. Referenced by main(), solve_neumann(), and solve_neumann_oct(). |
|
Definition at line 1141 of file omg.C. References ot::_p_DAMG::b, and ot::_p_DAMG::rhs. |
|
Definition at line 2710 of file omg.C. References PROF_PC_KSP_SHELL_APPLY_BEGIN, and PROF_PC_KSP_SHELL_APPLY_END. Referenced by ot::DAMGSetKSP(). |
|
Definition at line 2681 of file omg.C. References PROF_PC_KSP_SHELL_DESTROY_BEGIN, and PROF_PC_KSP_SHELL_DESTROY_END. Referenced by ot::DAMGSetKSP(). |
|
Definition at line 2562 of file omg.C. References ot::DA::iAmActive(), PROF_PC_KSP_SHELL_SETUP_BEGIN, and PROF_PC_KSP_SHELL_SETUP_END. Referenced by ot::DAMGSetKSP(). |
|
Prints detailed information about the meshes for each level.
Definition at line 2386 of file omg.C. References ot::_p_DAMG::da, ot::_p_DAMG::da_aux, DendroIntL, ot::DA::getCommActive(), ot::DA::getElementSize(), ot::DA::getGhostedElementSize(), ot::getGlobalMaxLevel(), ot::getGlobalMinLevel(), ot::DA::getNodeSize(), ot::DA::getNpesActive(), ot::DA::getNumBlocks(), ot::DA::iAmActive(), and ot::_p_DAMG::nlevels. Referenced by main(). |
|
Definition at line 2325 of file omg.C. References iC. Referenced by ot::createInterpolationType1(), and ot::createInterpolationType2(). |
|
Redistributes the vector, preserving the relative ordering of the elements.
Definition at line 2740 of file omg.C. References PROF_MG_SCATTER_BEGIN, and PROF_MG_SCATTER_END. Referenced by main(). |
|
Definition at line 190 of file externVars.h. |
|
Definition at line 168 of file externVars.h. Referenced by ot::DAMG_Finalize(), ot::DAMG_InitPrivateType1(), ot::DAMG_InitPrivateType2(), and ot::DAMG_InitPrivateType3(). |
|
Definition at line 167 of file externVars.h. Referenced by ot::DAMG_Finalize(), ot::DAMG_InitPrivateType1(), ot::DAMG_InitPrivateType2(), and ot::DAMG_InitPrivateType3(). |
|
Definition at line 169 of file externVars.h. Referenced by ot::DAMG_Finalize(), ot::DAMG_InitPrivateType1(), ot::DAMG_InitPrivateType2(), and ot::DAMG_InitPrivateType3(). |
|
Definition at line 170 of file externVars.h. Referenced by ot::DAMG_Finalize(), ot::DAMG_InitPrivateType1(), ot::DAMG_InitPrivateType2(), and ot::DAMG_InitPrivateType3(). |
|
Definition at line 171 of file externVars.h. Referenced by ot::DAMG_Finalize(), ot::DAMG_InitPrivateType1(), ot::DAMG_InitPrivateType2(), and ot::DAMG_InitPrivateType3(). |
|
Definition at line 172 of file externVars.h. Referenced by ot::DAMG_Finalize(), ot::DAMG_InitPrivateType1(), ot::DAMG_InitPrivateType2(), and ot::DAMG_InitPrivateType3(). |