VERB_code_2.2  2
 All Classes Namespaces Files Functions Variables Typedefs Macros Pages
MatrixSolver.h File Reference
#include "../Grid/Grid.h"
#include "../DiffusionCoefficient/DiffusionCoefficient.h"
#include "../GMRES/itlin.h"

Go to the source code of this file.

Functions

bool MakeMatrix (double *A, double *B1, double *C, double *R, double *x, double tau, double n, double f_lower, double f_upper, int nx, double dt, double *Dxx, double taulc, double alc, int g_flag, string lower_border_condition_type, string upper_border_condition_type, string approximationMethod="AM_Split_C")
 Make matrix for 1d diffusion. More...
 
bool SolveMatrix (double *f, double *A, double *B1, double *C, double *R, double f_lower, double f_upper, int nx, double dt)
 Solver for 1d general equation. More...
 
void over_relaxation_diag (DiagMatrix &A, Matrix1D< double > &B, Matrix1D< double > &X, int max_steps=1e5, double EE=1e-3)
 
void Lapack (DiagMatrix &A, Matrix1D< double > &B, Matrix1D< double > &X)
 
void gauss_solve (double *a, double *b, int n)
 Gauss inversion. More...
 
bool tridag (double a[], double b[], double c[], double r[], double u[], long n)
 Solve the AU=R system of equations, where A - tridiagonal matrix nxn with diagonals a[], b[], c[]. More...
 
void gmres_wrapout (CalculationMatrix &A, Matrix1D< double > &B, Matrix1D< double > &X, Parameters_structure::PSD::GMRES_parameters_structure GMRES_parameters)
 
void gmres2_wrapout (CalculationMatrix &A, Matrix1D< double > &B, Matrix1D< double > &X, Parameters_structure::PSD::GMRES_parameters_structure GMRES_parameters)
 
bool MakeModelMatrix_3D (CalculationMatrix &matr_A, CalculationMatrix &matr_B, CalculationMatrix &matr_C, Matrix3D< double > &L, Matrix3D< double > &pc, Matrix3D< double > &alpha, int L_size, int pc_size, int alpha_size, Matrix2D< double > &L_lowerBoundaryCondition, Matrix2D< double > &L_upperBoundaryCondition, Matrix2D< double > &pc_lowerBoundaryCondition, Matrix2D< double > &pc_upperBoundaryCondition, Matrix2D< double > &alpha_lowerBoundaryCondition, Matrix2D< double > &alpha_upperBoundaryCondition, string L_lowerBoundaryCondition_calculationType, string L_upperBoundaryCondition_calculationType, string pc_lowerBoundaryCondition_calculationType, string pc_upperBoundaryCondition_calculationType, string alpha_lowerBoundaryCondition_calculationType, string alpha_upperBoundaryCondition_calculationType, Matrix3D< double > &DLL, Matrix3D< double > &Dpcpc, Matrix3D< double > &DpcpcLpp, Matrix3D< double > &Daa, Matrix3D< double > &DaaLpp, Matrix3D< double > &Dpca, Matrix3D< double > &DpcaLpp, Matrix3D< double > &Jacobian, double dt, double Lpp, double tau=1e99, double tauLpp=1e99)
 
bool MakeModelMatrix_3D_KC (CalculationMatrix &matr_A, CalculationMatrix &matr_B, CalculationMatrix &matr_C, Matrix3D< double > &L, Matrix3D< double > &pc, Matrix3D< double > &alpha, int L_size, int pc_size, int alpha_size, Matrix2D< double > &L_lowerBoundaryCondition, Matrix2D< double > &L_upperBoundaryCondition, Matrix2D< double > &pc_lowerBoundaryCondition, Matrix2D< double > &pc_upperBoundaryCondition, Matrix2D< double > &alpha_lowerBoundaryCondition, Matrix2D< double > &alpha_upperBoundaryCondition, string L_lowerBoundaryCondition_calculationType, string L_upperBoundaryCondition_calculationType, string pc_lowerBoundaryCondition_calculationType, string pc_upperBoundaryCondition_calculationType, string alpha_lowerBoundaryCondition_calculationType, string alpha_upperBoundaryCondition_calculationType, Matrix3D< double > &DLL, Matrix3D< double > &Dpcpc, Matrix3D< double > &DpcpcLpp, Matrix3D< double > &Daa, Matrix3D< double > &DaaLpp, Matrix3D< double > &Dpca, Matrix3D< double > &DpcaLpp, Matrix3D< double > &Jacobian, double dt, double Lpp, double tau=1e99, double tauLpp=1e99)
 
void SecondDerivativeApproximation_3D (CalculationMatrix &matr_A, int il, int im, int ia, string FirstDerivative, string SecondDerivative, Matrix3D< double > &L, Matrix3D< double > &pc, Matrix3D< double > &alpha, Matrix3D< double > &Dxx, Matrix3D< double > &Jacobian, double multiplicator)
 
void SecondDerivativeApproximation_3D_KC_Diagonal (CalculationMatrix &matr_A, int il, int im, int ia, string FirstDerivative, string SecondDerivative, Matrix3D< double > &L, Matrix3D< double > &pc, Matrix3D< double > &alpha, Matrix3D< double > &Dxx, Matrix3D< double > &Jacobian, double multiplicator)
 
void SecondDerivativeApproximation_3D_KC_Mixed (CalculationMatrix &matr_A, int il, int im, int ia, string FirstDerivative, string SecondDerivative, Matrix3D< double > &L, Matrix3D< double > &pc, Matrix3D< double > &alpha, Matrix3D< double > &Dxx, Matrix3D< double > &Jacobian, double multiplicator)
 

Detailed Description

Making model matrixes, solving model matrixes.

Matrix form of linear equations: A*X[t+1] = B*X[t], where A - model matrix, B - RHS, X[t] - known values of function (PSD), X[t+1] - unknown values of function.

In that file there are procedures for making model matrix for 1d-diffusion, 2d-diffusion, some ideas of 3d-diffusion and mixed terms. Solver for tridiagonal matrix, solver by gauss method and iteration method (upper relaxation).

This file is under development, and has a lot of commented code, unfinished etc code.

There is a checked version for 1d diffusion (or split method of 2d, 3d diffusions) - 1d_universal_solver.

Author
Developed under supervision of the PI Yuri Shprits

Definition in file MatrixSolver.h.

Function Documentation

bool MakeMatrix ( double *  A,
double *  B1,
double *  C,
double *  R,
double *  x,
double  tau,
double  n,
double  f_lower,
double  f_upper,
int  nx,
double  dt,
double *  Dxx,
double  taulc,
double  alc,
int  g_flag,
string  lower_border_condition_type,
string  upper_border_condition_type,
string  approximationMethod 
)

Make matrix for 1d diffusion.

Make matrix for 1D diffusion.

Definition at line 71 of file MatrixSolver.cpp.

References B(), VF::G(), and VC::pi.

bool SolveMatrix ( double *  f,
double *  A,
double *  B1,
double *  C,
double *  R,
double  f_lower,
double  f_upper,
int  nx,
double  dt 
)

Solver for 1d general equation.

Solver for 1d-diffusion matrix (tridiagonal).

Definition at line 241 of file MatrixSolver.cpp.

References tridag().

void over_relaxation_diag ( DiagMatrix A,
Matrix1D< double > &  B,
Matrix1D< double > &  X,
int  max_steps,
double  EE 
)

Over relaxation iteration method.

Definition at line 1046 of file MatrixSolver.cpp.

References Output::echo(), and VF::max().

void Lapack ( DiagMatrix A,
Matrix1D< double > &  B,
Matrix1D< double > &  X 
)

Lapack inversion.

A * X = B - equation

Definition at line 953 of file MatrixSolver.cpp.

References Output::echo(), and VF::max().

void gauss_solve ( double *  a,
double *  b,
int  n 
)

Gauss inversion.

Definition at line 1206 of file MatrixSolver.cpp.

References EE, and I.

bool tridag ( double  a[],
double  b[],
double  c[],
double  r[],
double  u[],
long  n 
)

Solve the AU=R system of equations, where A - tridiagonal matrix nxn with diagonals a[], b[], c[].

Solver for a system of equations with 3-diagonal matrix.

Au = r Where A = diag(a, b, c),

Parameters
a[]- array, diagonal '-1' of the matrix
b[]- array, diagonal '0' of the matrix
c[]- array, diagonal '+1' of the matrix
r[]- array, r-vector
u[]- array, result
n- size of the matrix

Definition at line 1256 of file MatrixSolver.cpp.

bool MakeModelMatrix_3D ( CalculationMatrix matr_A,
CalculationMatrix matr_B,
CalculationMatrix matr_C,
Matrix3D< double > &  L,
Matrix3D< double > &  pc,
Matrix3D< double > &  alpha,
int  L_size,
int  pc_size,
int  alpha_size,
Matrix2D< double > &  L_lowerBoundaryCondition,
Matrix2D< double > &  L_upperBoundaryCondition,
Matrix2D< double > &  pc_lowerBoundaryCondition,
Matrix2D< double > &  pc_upperBoundaryCondition,
Matrix2D< double > &  alpha_lowerBoundaryCondition,
Matrix2D< double > &  alpha_upperBoundaryCondition,
string  L_lowerBoundaryCondition_calculationType,
string  L_upperBoundaryCondition_calculationType,
string  pc_lowerBoundaryCondition_calculationType,
string  pc_upperBoundaryCondition_calculationType,
string  alpha_lowerBoundaryCondition_calculationType,
string  alpha_upperBoundaryCondition_calculationType,
Matrix3D< double > &  DLL,
Matrix3D< double > &  Dpcpc,
Matrix3D< double > &  DpcpcLpp,
Matrix3D< double > &  Daa,
Matrix3D< double > &  DaaLpp,
Matrix3D< double > &  Dpca,
Matrix3D< double > &  DpcaLpp,
Matrix3D< double > &  Jacobian,
double  dt,
double  Lpp,
double  tau,
double  tauLpp 
)

Create matrix form of the Fokker-Planck equation: matr_A * PSD(t+1) = matr_B * PSD(t) + matr_C

Calculation matr_A, matr_B, and matr_C, i.e. numerical approximation of the derivatives Takes boundary conditions and diffusion coefficients as an input. L, pc, alpha should be orthogonal for 3D!!!

Definition at line 302 of file MatrixSolver.cpp.

References AddBoundary(), VF::alc(), VF::bounce_time_new(), CalculationMatrix::change_ind, CalculationMatrix::index1d(), VC::pi, and SecondDerivativeApproximation_3D().

bool MakeModelMatrix_3D_KC ( CalculationMatrix matr_A,
CalculationMatrix matr_B,
CalculationMatrix matr_C,
Matrix3D< double > &  L,
Matrix3D< double > &  pc,
Matrix3D< double > &  alpha,
int  L_size,
int  pc_size,
int  alpha_size,
Matrix2D< double > &  L_lowerBoundaryCondition,
Matrix2D< double > &  L_upperBoundaryCondition,
Matrix2D< double > &  pc_lowerBoundaryCondition,
Matrix2D< double > &  pc_upperBoundaryCondition,
Matrix2D< double > &  alpha_lowerBoundaryCondition,
Matrix2D< double > &  alpha_upperBoundaryCondition,
string  L_lowerBoundaryCondition_calculationType,
string  L_upperBoundaryCondition_calculationType,
string  pc_lowerBoundaryCondition_calculationType,
string  pc_upperBoundaryCondition_calculationType,
string  alpha_lowerBoundaryCondition_calculationType,
string  alpha_upperBoundaryCondition_calculationType,
Matrix3D< double > &  DLL,
Matrix3D< double > &  Dpcpc,
Matrix3D< double > &  DpcpcLpp,
Matrix3D< double > &  Daa,
Matrix3D< double > &  DaaLpp,
Matrix3D< double > &  Dpca,
Matrix3D< double > &  DpcaLpp,
Matrix3D< double > &  Jacobian,
double  dt,
double  Lpp,
double  tau,
double  tauLpp 
)

kckim test Create matrix form of the Fokker-Planck equation: matr_A * PSD(t+1) = matr_B * PSD(t) + matr_C

Calculation matr_A, matr_B, and matr_C, i.e. numerical approximation of the derivatives Takes boundary conditions and diffusion coefficients as an input. L, pc, alpha should be orthogonal for 3D!!!

Definition at line 1299 of file MatrixSolver.cpp.

References AddBoundary(), CalculationMatrix::change_ind, CalculationMatrix::index1d(), SecondDerivativeApproximation_3D_KC_Diagonal(), and SecondDerivativeApproximation_3D_KC_Mixed().

void SecondDerivativeApproximation_3D ( CalculationMatrix matr_A,
int  il,
int  im,
int  ia,
string  FirstDerivative,
string  SecondDerivative,
Matrix3D< double > &  L,
Matrix3D< double > &  pc,
Matrix3D< double > &  alpha,
Matrix3D< double > &  Dxx,
Matrix3D< double > &  Jacobian,
double  multiplicator 
)

Second derivative approximation, returns coefficients to be putted into the model matrix. $L_{\alpha \beta}(y) = (D_{\alpha \beta} \cdot y_{\bar{x}_\alpha})_{x_{\beta}}$

Samarskiy, page 261

Returns coefficients to be put into model matrix for an approximation of a second derivative.

Definition at line 1122 of file MatrixSolver.cpp.

References GetDerivativeVector(), and CalculationMatrix::index1d().

void SecondDerivativeApproximation_3D_KC_Diagonal ( CalculationMatrix matr_A,
int  il,
int  im,
int  ia,
string  FirstDerivative,
string  SecondDerivative,
Matrix3D< double > &  L,
Matrix3D< double > &  pc,
Matrix3D< double > &  alpha,
Matrix3D< double > &  Dxx,
Matrix3D< double > &  Jacobian,
double  multiplicator 
)

kckim test Second derivative approximation, returns coefficients to be putted into the model matrix. $L_{\alpha \beta}(y) = (D_{\alpha \beta} \cdot y_{\bar{x}_\alpha})_{x_{\beta}}$

Samarskiy, page 261

Returns coefficients to be put into model matrix for an approximation of a second derivative.

Definition at line 1657 of file MatrixSolver.cpp.

References GetDerivativeVector(), and CalculationMatrix::index1d().

void SecondDerivativeApproximation_3D_KC_Mixed ( CalculationMatrix matr_A,
int  il,
int  im,
int  ia,
string  FirstDerivative,
string  SecondDerivative,
Matrix3D< double > &  L,
Matrix3D< double > &  pc,
Matrix3D< double > &  alpha,
Matrix3D< double > &  Dxx,
Matrix3D< double > &  Jacobian,
double  multiplicator 
)

Definition at line 1720 of file MatrixSolver.cpp.

References GetDerivativeVector(), and CalculationMatrix::index1d().