/*****************************************************************************1*2* Elmer, A Finite Element Software for Multiphysical Problems3*4* Copyright 1st April 1995 - , CSC - IT Center for Science Ltd., Finland5*6* This library is free software; you can redistribute it and/or7* modify it under the terms of the GNU Lesser General Public8* License as published by the Free Software Foundation; either9* version 2.1 of the License, or (at your option) any later version.10*11* This library is distributed in the hope that it will be useful,12* but WITHOUT ANY WARRANTY; without even the implied warranty of13* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU14* Lesser General Public License for more details.15*16* You should have received a copy of the GNU Lesser General Public17* License along with this library (in file ../LGPL-2.1); if not, write18* to the Free Software Foundation, Inc., 51 Franklin Street,19* Fifth Floor, Boston, MA 02110-1301 USA20*21*****************************************************************************/2223/*******************************************************************************24*25* Solve symmetric positive definite eigenvalue problem.26*27*******************************************************************************28*29* Author: Juha Ruokolainen30*31* Address: CSC - IT Center for Science Ltd.32* Keilaranta 14, P.O. BOX 40533* 02101 Espoo, Finland34* Tel. +358 0 457 272335* Telefax: +358 0 457 230236* EMail: [email protected]37*38* Date: 30 May 199639*40* Modified by:41*42* Date of modification:43*44******************************************************************************/4546/*47* $Id: jacobi.c,v 1.1.1.1 2005/04/14 13:29:14 vierinen Exp $48*49* $Log: jacobi.c,v $50* Revision 1.1.1.1 2005/04/14 13:29:14 vierinen51* initial matc automake package52*53* Revision 1.2 1998/08/01 12:34:43 jpr54*55* Added Id, started Log.56*57*58*/5960#include "elmer/matc.h"6162/************************************************************************63P R O G R A M64To solve the generalized eigenproblem using the65generalized Jacobi iteration6667INPUT6869a(n,n) = Stiffness matrix (assumed positive definite)70b(n,n) = Mass matrix (assumed positive definite)71x(n,n) = Matrix storing eigenvectors on solution exit72eigv(n) = Vector storing eigenvalues on solution exit73d(n) = Working vector74n = Order of the matrices a and b75rtol = Convergence tolerance (usually set to 10.**-12)76nsmax = Maximum number of sweeps allowed (usually set to 15)7778OUTPUT7980a(n,n) = Diagonalized stiffness matrix81b(n,n) = Diagonalized mass matrix82x(n,n) = Eigenvectors stored columnwise83eigv(n) = Eigenvalues8485*************************************************************************/8687int matc_jacobi(double a[], double b[], double x[], double eigv[], double d[],88int n, double rtol)89{90register int i,j,k,ii,jj;91int nsmax=50, /* Max number of sweeps */92nsweep, /* Current sweep number */93nr,94jp1,jm1,kp1,km1,95convergence;9697double eps,98eptola,eptolb,99akk,ajj,ab,100check,sqch,101d1,d2,den,102ca,cg,103aj,bj,ak,bk,104xj,xk,105tol,dif,106epsa,epsb,107bb; /* Scale mass matrix */108109/************************************************************************110Initialize eigenvalue and eigenvector matrices111*************************************************************************/112for( i=0 ; i<n ; i++)113{114ii = i*n+i; /* address of diagonal element */115if (a[ii]<=0.0 || b[ii]<=0.0) return 0;116eigv[i] = d[i] = a[ii]/b[ii];117x[ii] = 1.0; /* make an unit matrix */118}119if(n==1) return 1; /* Return if single degree of freedom system */120121/************************************************************************122Initialize sweep counter and begin iteration123*************************************************************************/124nr=n - 1;125for ( nsweep = 0; nsweep < nsmax; nsweep++)126{127/************************************************************************128Check if present off-diagonal element is large enough to require zeroing129*************************************************************************/130131eps = pow(0.01, 2.0 * (nsweep + 1));132for(j=0 ; j<nr ; j++)133{134jj=j+1;135for(k=jj ; k<n ; k++)136{137eptola=( a[j*n+k]*a[j*n+k] / ( a[j*n+j]*a[k*n+k] ) );138eptolb=( b[j*n+k]*b[j*n+k] / ( b[j*n+j]*b[k*n+k] ) );139if ( eptola >= eps || eptolb >=eps )140{141142/************************************************************************143if zeroing is required, calculate the rotation matrix elements144*************************************************************************/145146akk=a[k*n+k]*b[j*n+k] - b[k*n+k]*a[j*n+k];147ajj=a[j*n+j]*b[j*n+k] - b[j*n+j]*a[j*n+k];148ab =a[j*n+j]*b[k*n+k] - a[k*n+k]*b[j*n+j];149check=(ab*ab+4.0*akk*ajj)/4.0;150151if (check<=0.0)152{153printf("***Error solution stop in *jacobi*\n");154printf(" check = %20.14e\n", check);155return 1;156}157sqch= sqrt(check);158d1 = ab/2.0+sqch;159d2 = ab/2.0-sqch;160den = d1;161if ( abs(d2) > abs(d1) ) den=d2;162163if (den == 0.0)164{165ca=0.0;166cg= -a[j*n+k] / a[k*n+k];167}168else169{170ca= akk/den;171cg= -ajj/den;172}173/************************************************************************174Perform the generalized rotation to zero the present off-diagonal element175*************************************************************************/176if (n != 2)177{178jp1=j+1;179jm1=j-1;180kp1=k+1;181km1=k-1;182/**************************************/183if ( jm1 >= 0)184for (i=0 ; i<=jm1 ; i++)185{186aj = a[i*n+j];187bj = b[i*n+j];188ak = a[i*n+k];189bk = b[i*n+k];190a[i*n+j] = aj+cg*ak;191b[i*n+j] = bj+cg*bk;192a[i*n+k] = ak+ca*aj;193b[i*n+k] = bk+ca*bj;194};195/**************************************/196if ((kp1-n+1) <= 0)197for (i=kp1 ; i<n ; i++)198{199aj = a[j*n+i];200bj = b[j*n+i];201ak = a[k*n+i];202bk = b[k*n+i];203a[j*n+i] = aj+cg*ak;204b[j*n+i] = bj+cg*bk;205a[k*n+i] = ak+ca*aj;206b[k*n+i] = bk+ca*bj;207};208/**************************************/209if ((jp1-km1) <= 0.0)210for (i=jp1 ; i<=km1 ; i++)211{212aj = a[j*n+i];213bj = b[j*n+i];214ak = a[i*n+k];215bk = b[i*n+k];216a[j*n+i] = aj+cg*ak;217b[j*n+i] = bj+cg*bk;218a[i*n+k] = ak+ca*aj;219b[i*n+k] = bk+ca*bj;220};221};222223ak = a[k*n+k];224bk = b[k*n+k];225a[k*n+k] = ak+2.0*ca*a[j*n+k]+ca*ca*a[j*n+j];226b[k*n+k] = bk+2.0*ca*b[j*n+k]+ca*ca*b[j*n+j];227a[j*n+j] = a[j*n+j]+2.0*cg*a[j*n+k]+cg*cg*ak;228b[j*n+j] = b[j*n+j]+2.0*cg*b[j*n+k]+cg*cg*bk;229a[j*n+k] = 0.0;230b[j*n+k] = 0.0;231232/************************************************************************233Update the eigenvector matrix after each rotation234*************************************************************************/235for (i=0 ; i<n ; i++)236{237xj = x[i*n+j];238xk = x[i*n+k];239x[i*n+j] = xj+cg*xk;240x[i*n+k] = xk+ca*xj;241};242};243};244};245/************************************************************************246Update the eigenvalues after each sweep247*************************************************************************/248for (i=0 ; i<n ; i++)249{250ii=i*n+i;251if (a[ii]<=0.0 || b[ii]<=0.0)252{253error( "*** Error solution stop in *jacobi*\n Matrix not positive definite.");254};255eigv[i] = a[ii] / b[ii];256};257258/************************************************************************259check for convergence260*************************************************************************/261convergence = 1;262for(i=0 ; i<n ; i++)263{264tol = rtol*d[i];265dif = abs(eigv[i]-d[i]);266if (dif > tol) convergence = 0;267if (! convergence) break;268};269/************************************************************************270check if all off-diag elements are satisfactorily small271*************************************************************************/272if (convergence)273{274eps=rtol*rtol;275for (j=0 ; j<nr ; j++)276{277jj=j+1;278for (k=jj ; k<n ; k++)279{280epsa = ( a[j*n+k]*a[j*n+k] / ( a[j*n+j]*a[k*n+k] ) );281epsb = ( b[j*n+k]*b[j*n+k] / ( b[j*n+j]*b[k*n+k] ) );282if ( epsa >= eps || epsb >=eps ) convergence = 0;283if (! convergence) break;284};285if (! convergence) break;286};287};288if (! convergence)289{290for (i=0 ; i<n ; i++)291d[i] = eigv[i];292};293if (convergence) break;294};295/************************************************************************296Fill out bottom triangle of resultant matrices and scale eigenvectors297*************************************************************************/298for (i=0 ; i<n ; i++)299for (j=i ; j<n ; j++)300{301b[j*n+i] = b[i*n+j];302a[j*n+i] = a[i*n+j];303};304305for (j=0 ; j<n ; j++)306{307bb = sqrt( b[j*n+j] );308for (k=0 ; k<n ; k++)309x[k*n+j] /= bb;310};311PrintOut( "jacobi: nsweeps %d\n", nsweep );312return 1 ;313}314315VARIABLE *mtr_jacob(VARIABLE *a)316{317VARIABLE *x, *ev;318319double *b, *d, rtol;320int i, n;321322if (NROW(a) != NCOL(a))323error("Jacob: Matrix must be square.\n");324325b = MATR(NEXT(a));326n = NROW(a);327328if (NROW(NEXT(a)) != NCOL(NEXT(a)) || n != NROW(NEXT(a)))329error("Jacob: Matrix dimensions incompatible.\n");330331rtol = *MATR(NEXT(NEXT(a)));332333x = var_new("eigv", TYPE_DOUBLE, NROW(a), NCOL(a));334d = (double *)ALLOCMEM(n * sizeof(double));335ev = var_temp_new(TYPE_DOUBLE, 1, n);336337(void)matc_jacobi(MATR(a), b, MATR(x), MATR(ev), d, n, rtol);338FREEMEM((char *)d);339340return ev;341}342343344