void showQMatrix( QMatrix *M ) { int i, j; printf( "Matrix: %s\n", Q_GetName( M ) ); for ( i = 0; i < Q_GetDim( M ); i++ ) { printf( "Row %d: ", i ); for ( j = 0; j < Q_GetDim( M ); j++ ) if ( Equal( Q_GetEl( M, i+1, j+1 ), 0.0, 1e-5 )) continue; else printf( "(%d) %.3lf ", j, Q_GetEl( M, i+1, j+1 ) ); printf( "\n" ); } /* for i */ } /* showQMatrix */
double GetMaxEigenval(QMatrix *A, PrecondProcType PrecondProc, double OmegaPrecond) /* returns estimate for maximum eigenvalue of the matrix A */ { double MaxEigenval; EigenvalInfoType *EigenvalInfo; Q_Lock(A); if (LASResult() == LASOK) { EigenvalInfo = (EigenvalInfoType *)*(Q_EigenvalInfo(A)); /* if eigenvalues not estimated yet, ... */ if (EigenvalInfo == NULL) { EigenvalInfo = (EigenvalInfoType *)malloc(sizeof(EigenvalInfoType)); if (EigenvalInfo != NULL) { *(Q_EigenvalInfo(A)) = (void *)EigenvalInfo; EstimEigenvals(A, PrecondProc, OmegaPrecond); } else { LASError(LASMemAllocErr, "GetMaxEigenval", Q_GetName(A), NULL, NULL); } } /* if eigenvalues estimated with an other preconditioner, ... */ if (EigenvalInfo->PrecondProcUsed != PrecondProc || EigenvalInfo->OmegaPrecondUsed != OmegaPrecond) { EstimEigenvals(A, PrecondProc, OmegaPrecond); } if (LASResult() == LASOK) MaxEigenval = EigenvalInfo->MaxEigenval; else MaxEigenval = 1.0; } else { MaxEigenval = 1.0; } return(MaxEigenval); }
Vector *Test4_QV(QMatrix *Q, Vector *V) /* VRes = Q * V ... implementation using local variables and pointers, descended counting of matrix elements */ { Vector *VRes; char *VResName; size_t Dim, Row, Clm, Len, ElCount; size_t *QLen; ElType **QEl, *PtrEl; Real Sum, Cmp; Real *VCmp, *VResCmp; if (LASResult() == LASOK) { if (Q->Dim == V->Dim) { Dim = Q->Dim; VRes = (Vector *)malloc(sizeof(Vector)); VResName = (char *)malloc((strlen(Q_GetName(Q)) + strlen(V_GetName(V)) + 10) * sizeof(char)); if (VRes != NULL && VResName != NULL) { sprintf(VResName, "(%s) * (%s)", Q_GetName(Q), V_GetName(V)); V_Constr(VRes, VResName, Dim, Tempor, True); if (LASResult() == LASOK) { /* initialisation of vector VRes */ if (Q->Symmetry || Q->ElOrder == Clmws) V_SetAllCmp(VRes, 0.0); /* analysis of multipliers of matrix Q and vector V is not implemented yet */ /* multiplication of matrix elements by vector components */ VCmp = V->Cmp; VResCmp = VRes->Cmp; QLen = Q->Len; QEl = Q->El; if (!Q->Symmetry) { if (Q->ElOrder == Rowws) { for (Row = Dim; Row > 0; Row--) { Len = QLen[Row]; PtrEl = QEl[Row] + Len - 1; Sum = 0.0; for (ElCount = Len; ElCount > 0; ElCount--) { Sum += (*PtrEl).Val * VCmp[(*PtrEl).Pos]; PtrEl--; } VResCmp[Row] = Sum; } } if (Q->ElOrder == Clmws) { for (Clm = Dim; Clm > 0; Clm--) { Len = QLen[Clm]; PtrEl = QEl[Clm] + Len - 1; Cmp = VCmp[Clm]; for (ElCount = Len; ElCount > 0; ElCount--) { VResCmp[(*PtrEl).Pos] += (*PtrEl).Val * Cmp; PtrEl--;; } } } } else { /* multiplication by symmetric matrix is not implemented yet */ V_SetAllCmp(VRes, 0.0); } } } else { LASError(LASMemAllocErr, "Mul_QV", Q_GetName(Q), V_GetName(V), NULL); if (VRes != NULL) free(VRes); if (VResName != NULL) free(VResName); } } else { LASError(LASDimErr, "Mul_QV", Q_GetName(Q), V_GetName(V), NULL); VRes = NULL; } } else { VRes = NULL; } if (Q != NULL) { if (Q->Instance == Tempor) { Q_Destr(Q); free(Q); } } if (V != NULL) { if (V->Instance == Tempor) { V_Destr(V); free(V); } } return(VRes); }
Vector *Test1_QV(QMatrix *Q, Vector *V) /* VRes = Q * V ... very simple implementation */ { Vector *VRes; char *VResName; size_t Dim, Row, Clm, ElCount; ElType **QEl; if (LASResult() == LASOK) { if (Q->Dim == V->Dim) { Dim = Q->Dim; QEl = Q->El; VRes = (Vector *)malloc(sizeof(Vector)); VResName = (char *)malloc((strlen(Q_GetName(Q)) + strlen(V_GetName(V)) + 10) * sizeof(char)); if (VRes != NULL && VResName != NULL) { sprintf(VResName, "(%s) * (%s)", Q_GetName(Q), V_GetName(V)); V_Constr(VRes, VResName, Dim, Tempor, True); if (LASResult() == LASOK) { /* initialisation of vector VRes */ V_SetAllCmp(VRes, 0.0); /* analysis of multipliers of matrix Q and vector V is not implemented yet */ /* multiplication of matrix elements by vector components */ if (!Q->Symmetry) { if (Q->ElOrder == Rowws) { for (Row = 1; Row <= Dim; Row++) { for (ElCount = 0; ElCount < Q->Len[Row]; ElCount++) { Clm = QEl[Row][ElCount].Pos; VRes->Cmp[Row] += QEl[Row][ElCount].Val * V->Cmp[Clm]; } } } if (Q->ElOrder == Clmws) { for (Clm = 1; Clm <= Dim; Clm++) { for (ElCount = 0; ElCount < Q->Len[Clm]; ElCount++) { Row = QEl[Clm][ElCount].Pos; VRes->Cmp[Row] += QEl[Clm][ElCount].Val * V->Cmp[Clm]; } } } } else { /* multiplication by symmetric matrix is not implemented yet */ V_SetAllCmp(VRes, 0.0); } } } else { LASError(LASMemAllocErr, "Mul_QV", Q_GetName(Q), V_GetName(V), NULL); if (VRes != NULL) free(VRes); if (VResName != NULL) free(VResName); } } else { LASError(LASDimErr, "Mul_QV", Q_GetName(Q), V_GetName(V), NULL); VRes = NULL; } } else { VRes = NULL; } if (Q != NULL) { if (Q->Instance == Tempor) { Q_Destr(Q); free(Q); } } if (V != NULL) { if (V->Instance == Tempor) { V_Destr(V); free(V); } } return(VRes); }
static void EstimEigenvals(QMatrix *A, PrecondProcType PrecondProc, double OmegaPrecond) /* estimates extremal eigenvalues of the matrix A by means of the Lanczos method */ { /* * for details to the Lanczos algorithm see * * G. H. Golub, Ch. F. van Loan: * Matrix Computations; * North Oxford Academic, Oxford, 1986 * * (for modification for preconditioned matrices compare with sec. 10.3) * */ double LambdaMin = 0.0, LambdaMax = 0.0; double LambdaMinOld, LambdaMaxOld; double GershBoundMin = 0.0, GershBoundMax = 0.0; double *Alpha, *Beta; size_t Dim, j; Boolean Found; Vector q, qOld, h, p; Q_Lock(A); Dim = Q_GetDim(A); V_Constr(&q, "q", Dim, Normal, True); V_Constr(&qOld, "qOld", Dim, Normal, True); V_Constr(&h, "h", Dim, Normal, True); if (PrecondProc != NULL) V_Constr(&p, "p", Dim, Normal, True); if (LASResult() == LASOK) { Alpha = (double *)malloc((Dim + 1) * sizeof(double)); Beta = (double *)malloc((Dim + 1) * sizeof(double)); if (Alpha != NULL && Beta != NULL) { j = 0; V_SetAllCmp(&qOld, 0.0); V_SetRndCmp(&q); if (Q_KerDefined(A)) OrthoRightKer_VQ(&q, A); if (Q_GetSymmetry(A) && PrecondProc != NULL) { (*PrecondProc)(A, &p, &q, OmegaPrecond); MulAsgn_VS(&q, 1.0 / sqrt(Mul_VV(&q, &p))); } else { MulAsgn_VS(&q, 1.0 / l2Norm_V(&q)); } Beta[0] = 1.0; do { j++; if (Q_GetSymmetry(A) && PrecondProc != NULL) { /* p = M^(-1) q */ (*PrecondProc)(A, &p, &q, OmegaPrecond); /* h = A p */ Asgn_VV(&h, Mul_QV(A, &p)); if (Q_KerDefined(A)) OrthoRightKer_VQ(&h, A); /* Alpha = p . h */ Alpha[j] = Mul_VV(&p, &h); /* r = h - Alpha q - Beta qOld */ SubAsgn_VV(&h, Add_VV(Mul_SV(Alpha[j], &q), Mul_SV(Beta[j-1], &qOld))); /* z = M^(-1) r */ (*PrecondProc)(A, &p, &h, OmegaPrecond); /* Beta = sqrt(r . z) */ Beta[j] = sqrt(Mul_VV(&h, &p)); Asgn_VV(&qOld, &q); /* q = r / Beta */ Asgn_VV(&q, Mul_SV(1.0 / Beta[j], &h)); } else { /* h = A p */ if (Q_GetSymmetry(A)) { Asgn_VV(&h, Mul_QV(A, &q)); } else { if (PrecondProc != NULL) { (*PrecondProc)(A, &h, Mul_QV(A, &q), OmegaPrecond); (*PrecondProc)(Transp_Q(A), &h, &h, OmegaPrecond); Asgn_VV(&h, Mul_QV(Transp_Q(A), &h)); } else { Asgn_VV(&h, Mul_QV(Transp_Q(A), Mul_QV(A, &q))); } } if (Q_KerDefined(A)) OrthoRightKer_VQ(&h, A); /* Alpha = q . h */ Alpha[j] = Mul_VV(&q, &h); /* r = h - Alpha q - Beta qOld */ SubAsgn_VV(&h, Add_VV(Mul_SV(Alpha[j], &q), Mul_SV(Beta[j-1], &qOld))); /* Beta = || r || */ Beta[j] = l2Norm_V(&h); Asgn_VV(&qOld, &q); /* q = r / Beta */ Asgn_VV(&q, Mul_SV(1.0 / Beta[j], &h)); } LambdaMaxOld = LambdaMax; LambdaMinOld = LambdaMin; /* determination of extremal eigenvalues of the tridiagonal matrix (Beta[i-1] Alpha[i] Beta[i]) (where 1 <= i <= j) by means of the method of bisection; bounds for eigenvalues are determined after Gershgorin circle theorem */ if (j == 1) { GershBoundMin = Alpha[1] - fabs(Beta[1]); GershBoundMax = Alpha[1] + fabs(Beta[1]); LambdaMin = Alpha[1]; LambdaMax = Alpha[1]; } else { GershBoundMin = min(Alpha[j] - fabs(Beta[j]) - fabs(Beta[j - 1]), GershBoundMin); GershBoundMax = max(Alpha[j] + fabs(Beta[j]) + fabs(Beta[j - 1]), GershBoundMax); SearchEigenval(j, Alpha, Beta, 1, GershBoundMin, LambdaMin, &Found, &LambdaMin); if (!Found) SearchEigenval(j, Alpha, Beta, 1, GershBoundMin, GershBoundMax, &Found, &LambdaMin); SearchEigenval(j, Alpha, Beta, j, LambdaMax, GershBoundMax, &Found, &LambdaMax); if (!Found) SearchEigenval(j, Alpha, Beta, j, GershBoundMin, GershBoundMax, &Found, &LambdaMax); } } while (!IsZero(Beta[j]) && j < Dim && (fabs(LambdaMin - LambdaMinOld) > EigenvalEps * LambdaMin || fabs(LambdaMax - LambdaMaxOld) > EigenvalEps * LambdaMax) && LASResult() == LASOK); if (Q_GetSymmetry(A)) { LambdaMin = (1.0 - j * EigenvalEps) * LambdaMin; } else { LambdaMin = (1.0 - sqrt(j) * EigenvalEps) * sqrt(LambdaMin); } if (Alpha != NULL) free(Alpha); if (Beta != NULL) free(Beta); } else { LASError(LASMemAllocErr, "EstimEigenvals", Q_GetName(A), NULL, NULL); } } V_Destr(&q); V_Destr(&qOld); V_Destr(&h); if (PrecondProc != NULL) V_Destr(&p); if (LASResult() == LASOK) { ((EigenvalInfoType *)*(Q_EigenvalInfo(A)))->MinEigenval = LambdaMin; ((EigenvalInfoType *)*(Q_EigenvalInfo(A)))->MaxEigenval = LambdaMax; ((EigenvalInfoType *)*(Q_EigenvalInfo(A)))->PrecondProcUsed = PrecondProc; ((EigenvalInfoType *)*(Q_EigenvalInfo(A)))->OmegaPrecondUsed = OmegaPrecond; } else { ((EigenvalInfoType *)*(Q_EigenvalInfo(A)))->MinEigenval = 1.0; ((EigenvalInfoType *)*(Q_EigenvalInfo(A)))->MaxEigenval = 1.0; ((EigenvalInfoType *)*(Q_EigenvalInfo(A)))->PrecondProcUsed = NULL; ((EigenvalInfoType *)*(Q_EigenvalInfo(A)))->OmegaPrecondUsed = 1.0; } Q_Unlock(A); }