System Preprocessors
|
#include <stdlib.h>
#include <stdio.h>
#include "petsc.h"
#include "petscis.h"
#include "syspro.h"
#include "sysprotransform.h"
#include "sysprolinear.h"
#include "linear_impl.h"
Go to the source code of this file.
Defines | |
#define | PREPROCESSOR "approximation" |
Functions | |
static PetscErrorCode | MatSymmetricPart (NumericalProblem inproblem, NumericalProblem outproblem) |
static PetscErrorCode | MatGustafssonMod (NumericalProblem inproblem, NumericalProblem outproblem) |
static PetscErrorCode | setup_approximation_choices () |
static PetscErrorCode | specific_approximation_choices (NumericalProblem inproblem, SalsaTransform transform) |
static PetscErrorCode | approximate_system (const char *type, int nopt, PetscBool overwrite, NumericalProblem inproblem, NumericalProblem *outproblem, void *gctx, void **ctx, PetscBool *success) |
static PetscErrorCode | unapproximate_system (const char *type, PetscBool overwrite, void *gctx, void *ctx, NumericalProblem problem, NumericalProblem nextproblem, NumericalSolution before, NumericalSolution after) |
PetscErrorCode | DeclareApproximationPreprocessor (void) |
Definition in file approximating.c.
#define PREPROCESSOR "approximation" |
Definition at line 39 of file approximating.c.
Referenced by DeclareApproximationPreprocessor(), setup_approximation_choices(), and specific_approximation_choices().
static PetscErrorCode approximate_system | ( | const char * | type, |
int | nopt, | ||
PetscBool | overwrite, | ||
NumericalProblem | inproblem, | ||
NumericalProblem * | outproblem, | ||
void * | gctx, | ||
void ** | ctx, | ||
PetscBool * | success | ||
) | [static] |
Definition at line 220 of file approximating.c.
References CHKERRQ(), ierr, LinearSystemDuplicatePointers(), LinearSystemGetParts(), MatGustafssonMod(), and MatSymmetricPart().
Referenced by DeclareApproximationPreprocessor().
PetscErrorCode DeclareApproximationPreprocessor | ( | void | ) |
Definition at line 280 of file approximating.c.
References approximate_system(), CHKERRQ(), DeclarePreprocessor(), ierr, PREPROCESSOR, PreprocessorSetPreservedCategories(), setup_approximation_choices(), specific_approximation_choices(), and unapproximate_system().
static PetscErrorCode MatGustafssonMod | ( | NumericalProblem | inproblem, |
NumericalProblem | outproblem | ||
) | [static] |
Definition at line 92 of file approximating.c.
References CHKERRQ(), ierr, LinearSystemGetParts(), and LinearSystemSetParts().
Referenced by approximate_system().
static PetscErrorCode MatSymmetricPart | ( | NumericalProblem | inproblem, |
NumericalProblem | outproblem | ||
) | [static] |
Definition at line 44 of file approximating.c.
References CHKERRQ(), ierr, LinearSystemGetParts(), LinearSystemSetParts(), SysProComputeQuantity(), and SysProRetrieveQuantity().
Referenced by approximate_system().
static PetscErrorCode setup_approximation_choices | ( | ) | [static] |
Definition at line 146 of file approximating.c.
References CHKERRQ(), ierr, NewTransformObject(), PREPROCESSOR, TransformGetByName(), and TransformObjectSetExplanation().
Referenced by DeclareApproximationPreprocessor().
static PetscErrorCode specific_approximation_choices | ( | NumericalProblem | inproblem, |
SalsaTransform | transform | ||
) | [static] |
This is the 'specific setup' phase of the approximation preprocessor. See Usage modes for details.
This routine eliminates the Gustafsson approximation for diagonally dominant systems, and the symmetric for symmetric systems.
Definition at line 178 of file approximating.c.
References CHKERRQ(), ierr, LinearSystemGetParts(), PREPROCESSOR, SysProRetrieveQuantity(), TransformObjectGetByName(), and TransformObjectMark().
Referenced by DeclareApproximationPreprocessor().
static PetscErrorCode unapproximate_system | ( | const char * | type, |
PetscBool | overwrite, | ||
void * | gctx, | ||
void * | ctx, | ||
NumericalProblem | problem, | ||
NumericalProblem | nextproblem, | ||
NumericalSolution | before, | ||
NumericalSolution | after | ||
) | [static] |
Definition at line 265 of file approximating.c.
References CHKERRQ(), ierr, and LinearSolutionCopy().
Referenced by DeclareApproximationPreprocessor().