Komplex  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Functions
azk_extract_solution.c File Reference

Extraction routine for getting the solution of a Komplex system. More...

#include <stdio.h>
#include <stdlib.h>
#include <math.h>
#include "az_aztec.h"
#include "azk_komplex.h"
Include dependency graph for azk_extract_solution.c:

Functions

void AZK_extract_solution_k2c (int *options, double *params, int *proc_config, AZ_MATRIX *Amat_komplex, AZ_PRECOND *Prec, double *vk, double *vc)
 Extract a Complex vector from a Komplex vector. More...
 
void AZK_extract_solution_k2g (int *options, double *params, int *proc_config, AZ_MATRIX *Amat_komplex, AZ_PRECOND *Prec, double *vk, double *vr, double *vi)
 
void AZK_extract_solution_k2ri (int *options, double *params, int *proc_config, AZ_MATRIX *Amat_komplex, AZ_PRECOND *Prec, double *vk, double *vr, double *vi)
 Extract real/imaginary parts of a complex vector from a Komplex vector. More...
 

Detailed Description

Extraction routine for getting the solution of a Komplex system.

KOMPLEX is an add-on module to AZTEC that allows users to solve complex-valued linear systems.

KOMPLEX solves a complex-valued linear system Ax = b by solving an equivalent real-valued system of twice the dimension. Specifically, writing in terms of real and imaginary parts, we have

\[ (A_r + i*A_i)*(x_r + i*x_i) = (b_r + i*b_i) \]

or by separating into real and imaginary equations we have

\[ \left( \begin{array}{rr} A_r & -A_i\\ A_i & A_r \end{array} \right) \left( \begin{array}{r} x_r\\ x_i \end{array} \right) = \left( \begin{array}{r} b_r\\ b_i \end{array} \right) \]

which is a real-valued system of twice the size. If we find xr and xi, we can form the solution to the original system as x = xr +i*xi.

KOMPLEX accept user linear systems in three forms with either global or local index values.

1) The first form is true complex. The user passes in an MSR or VBR format matrix where the values are stored like Fortran complex numbers. Thus, the values array is of type double that is twice as long as the number of complex values. Each complex entry is stored with real part followed by imaginary part (as in Fortran).

2) The second form stores real and imaginary parts separately, but the pattern for each is identical. Thus only the values of the imaginary part are passed to the creation routines.

3) The third form accepts two real-valued matrices with no assumption about the structure of the matrices. Each matrix is multiplied by a user-supplied complex constant. This is the most general form.

Each of the above forms supports a global or local index set. By this we mean that the index values (stored in bindx) refer to the global problem indices, or the local indices (for example after calling AZ_transform).

Function Documentation

void AZK_extract_solution_k2c ( int *  options,
double *  params,
int *  proc_config,
AZ_MATRIX *  Amat_komplex,
AZ_PRECOND *  Prec,
double *  vk,
double *  vc 
)

Extract a Complex vector from a Komplex vector.

Transforms a komplex vector to a complex vector.

Parameters
options(In) Determines specific solution method and other parameters.
params(In) Drop tolerance and convergence tolerance info.
proc_config(In) Machine configuration. proc_config[AZ_node] is the node number. proc_config[AZ_N_procs] is the number of processors.
Amat_komplex(In) Komplex version of matrix stored as an AZ_MATRIX structure.
Prec(In) Preconditioner for Amat stored as an AZ_PRECOND structure.
vk(In) Komplex version of vector.
vc(Out) Contains a complex vector with the real/imag parts interleaved as in Fortran complex format. Note that the user must allocate sufficient storage for results.
void AZK_extract_solution_k2ri ( int *  options,
double *  params,
int *  proc_config,
AZ_MATRIX *  Amat_komplex,
AZ_PRECOND *  Prec,
double *  vk,
double *  vr,
double *  vi 
)

Extract real/imaginary parts of a complex vector from a Komplex vector.

Transforms a komplex vector to real and imaginary parts.

Parameters
options(In) Determines specific solution method and other parameters.
params(In) Drop tolerance and convergence tolerance info.
proc_config(In) Machine configuration. proc_config[AZ_node] is the node number. proc_config[AZ_N_procs] is the number of processors.
Amat_komplex(In) Komplex version of matrix stored as an AZ_MATRIX structure.
Prec(In) Preconditioner for Amat stored as an AZ_PRECOND structure.
vk(In) Komplex version of vector.
vc(Out) Contains a complex vector with the real/imag parts interleaved as in Fortran complex format. Note that the user must allocate sufficient storage for results.