中文
注册
我要评分
文档获取效率
文档正确性
内容完整性
文档易理解
在线提单
论坛求助
鲲鹏小智

KmlScaissGmresGet?I?

获得迭代求解相关参数。

接口定义

C Interface:

int KmlScaissGmresGetSII(KmlScasolverTask **handle, enum KmlSolverParam param, int *data, int nd);

int KmlScaissGmresGetSIS(KmlScasolverTask **handle, enum KmlSolverParam param, float *data, int nd);

int KmlScaissGmresGetDII(KmlScasolverTask **handle, enum KmlSolverParam param, int *data, int nd);

int KmlScaissGmresGetDID(KmlScasolverTask **handle, enum KmlSolverParam param, double *data, int nd);

参数

参数名

类型

描述

输入/输出

handle

KmlScasolverTask **

求解器句柄,传入之前步骤的变量。

输入/输出

param

enum KmlSolverParam

  • KMLSS_ITERATION_COUNT表示迭代次数。
  • KMLSS_MAX_ITERATION_COUNT表示设置的最大迭代次数。
  • KMLSS_TOLERANCE表示最终的迭代相对残差。
  • KMLSS_THRESHOLD表示设置的相对残差。

输入

data

  • 在KmlScaissGmresGetSII/KmlScaissGmresGetDII中为int *。
  • 在KmlScaissGmresGetSIS中为float *。
  • 在KmlScaissGmresGetDID中为double *。
  • 在KmlScaissGmresGetSII/KmlScaissGmresGetDII中为迭代步数。
  • 在KmlScaissGmresGetSIS中为残差。
  • 在KmlScaissGmresGetDID中为残差。

输出

nd

int

data数组元素个数。

输入

返回值

返回值

类型

描述

KMLSS_NO_ERROR

int

正常执行。

KMLSS_DATA_SIZE

int

参数nd不等于1。

KMLSS_NULL_ARGUMENT

int

handle,data中存在空参数。

KMLSS_BAD_SELECTOR

int

param为无效参数。

依赖

#include "kml_scaiss.h"

示例

C Interface:
MPI_Init(NULL, NULL);
int size, rank;
MPI_Comm_size(MPI_COMM_WORLD, &size);//获取总进程数
MPI_Comm_rank(MPI_COMM_WORLD, &rank);//获取当前进程标识

int mat_size = 8;
int n = mat_size / size;
int n_beg = n * rank;
if (n * size != mat_size && rank == (size - 1)) {
    n = mat_size - n * rank;
}
int ja[26] = { 0, 3, 4, 1, 2, 3, 5, 1, 2, 7, 0, 1, 3, 6, 0, 4, 5, 1, 4, 5, 7, 3, 6, 2, 5, 7 };
double a[26] = { 1.0, 1.0, 2.0, 9.0, 2.0, 1.0, -3.0, 2.0, 3.0, 2.0, 1.0, 1.0, 9.0, -5.0, 2.0, 6.0, 1.0, -3.0, 1.0, 4.0, 1.0, -5.0, 7.0, 2.0, 1.0, 2.0 };
int ia[9] = { 0, 3, 7, 10, 14, 17, 21, 23, 26 };
int a_beg = ia[n_beg];
for (int i = n_beg; i < (n_beg + n + 1); i++) {
    ia[i] -= a_beg;
}

/* Internal KML_SCAISS structure */
KmlScasolverTask *handle;

/* KML_SCAISS control parameters */
int error;            /* Output error handle */

/* Create data structures */
const double *a_holder = &a[a_beg];
const int *ja_holder = &ja[a_beg];
const int *ia_holder = &ia[n_beg];
error = KmlScaissGmresInitStripesDI(&handle, mat_size, 1, &n, &n_beg, &a_holder, &ja_holder, &ia_holder, MPI_COMM_WORLD);

double eps = 1e-4;
error = KmlScaissGmresSetDID(&handle, KMLSS_THRESHOLD, &eps, 1); 
if (error != 0) {
    printf("\nERROR in KmlScaissGmresSetDID with KMLSS_THRESHOLD: %d", error);
    return 1;
}
int max_iters = 2000;
error = KmlScaissGmresSetDII(&handle, KMLSS_MAX_ITERATION_COUNT, &max_iters, 1);
if (error != 0) {
    printf("\nERROR in KmlScaissGmresSetDII with KMLSS_MAX_ITERATION_COUNT: %d", error);
    return 1;
}
int precond = KMLSS_BJACOBI;
error = KmlScaissGmresSetDII(&handle, KMLSS_PRECONDITIONER_TYPE, &precond, 1);
if (error != 0) {
    printf("\nERROR in KmlScaissGmresSetDII with KMLSS_PRECONDITIONER_TYPE: %d", error);
    return 1;
}
int num_blocks = 1;
error = KmlScaissGmresPcSetDII(&handle, KMLSS_NUM_BLOCKS, &num_blocks, 1);
if (error != 0) {
	printf("ERROR in KmlScaissGmresPcSetDII: %d\n", error);
	return 1;
}
int subprecond = KMLSS_ICC;
error = KmlScaissGmresPcSetDII(&handle, KMLSS_BLOCK_METHOD, &subprecond, 1);
if (error != 0) {
    printf("ERROR in KmlScaissGmresPcSetDII: %d\n", error);
    return 1;
}
error = KmlScaissGmresAnalyzeDI(&handle);
if (error != 0) {
    printf("ERROR in KmlScaissGmresAnalyzeDI: %d\n", error);
    return 1;
}
error = KmlScaissGmresFactorizeDI(&handle);
if (error != 0) {
    printf("ERROR in KmlScaissGmresFactorizeDI: %d\n", error);
    return 1;
}

int nrhs = 1;
int ldx = n, ldb = n;
double b[8] = { 4.0, 9.0, 7.0, 6.0, 9.0, 3.0, 2.0, 5.0 };
double x[8] = { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0 };
error = KmlScaissGmresSolveDxDI(&handle, nrhs, &x[n_beg], ldx, &b[n_beg], ldb);
if (error != 0) {
    printf("ERROR in KmlScaissGmresSolveDI: %d\n", error);
    return 1;
}

int idata[] = {0};
error = KmlScaissGmresGetDII(&handle, KMLSS_ITERATION_COUNT, idata, 1); 
if (error != 0) {
    printf("\nERROR in KmlScaissGmresGetDII: %d", error);
    return 1;
    }
double ddata[] = {0.0};
error = KmlScaissGmresGetDID(&handle, KMLSS_TOLERANCE, ddata, 1);
if (error != 0) {
    printf("\nERROR in KmlScaissGmresGetDID: %d", error);
    return 1;
    }
if (rank == 0) printf("[Solve info] Res= \t%.12g\t, iteration = \t%d\t\n", sqrt(ddata[0]), idata[0]);

运行结果(2个进程):

[Solve info] Res=      2.04475505752e-15       , iteration =   8