OOFEM 3.0
Loading...
Searching...
No Matches
subspaceit.C
Go to the documentation of this file.
1/*
2 *
3 * ##### ##### ###### ###### ### ###
4 * ## ## ## ## ## ## ## ### ##
5 * ## ## ## ## #### #### ## # ##
6 * ## ## ## ## ## ## ## ##
7 * ## ## ## ## ## ## ## ##
8 * ##### ##### ## ###### ## ##
9 *
10 *
11 * OOFEM : Object Oriented Finite Element Code
12 *
13 * Copyright (C) 1993 - 2025 Borek Patzak
14 *
15 *
16 *
17 * Czech Technical University, Faculty of Civil Engineering,
18 * Department of Structural Mechanics, 166 29 Prague, Czech Republic
19 *
20 * This library is free software; you can redistribute it and/or
21 * modify it under the terms of the GNU Lesser General Public
22 * License as published by the Free Software Foundation; either
23 * version 2.1 of the License, or (at your option) any later version.
24 *
25 * This program is distributed in the hope that it will be useful,
26 * but WITHOUT ANY WARRANTY; without even the implied warranty of
27 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
28 * Lesser General Public License for more details.
29 *
30 * You should have received a copy of the GNU Lesser General Public
31 * License along with this library; if not, write to the Free Software
32 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
33 */
34
35//#define DETAILED_REPORT
36
37#include "subspaceit.h"
38#include "engngm.h"
39#include "domain.h"
40#include "floatmatrix.h"
41#include "floatarray.h"
42#include "mathfem.h"
43#include "gjacobi.h"
44#include "sparselinsystemnm.h"
45#include "classfactory.h"
46
47namespace oofem {
49
50SubspaceIteration :: SubspaceIteration(Domain *d, EngngModel *m) :
52 nitem(40)
53{
54}
55
56
58SubspaceIteration :: solve(SparseMtrx &a, SparseMtrx &b, FloatArray &_eigv, FloatMatrix &_r, double rtol, int nroot)
59//
60// this function solve the generalized eigenproblem using the Generalized
61// jacobi iteration
62//
63{
65 OOFEM_ERROR("matrices size mismatch");
66 }
67
68 FloatArray temp, w, d, tt, f, rtolv, eigv;
70 int nc1, ij = 0;
71 FloatMatrix ar, br, vec;
72 std :: unique_ptr< SparseLinearSystemNM > solver( GiveClassFactory().createSparseLinSolver(ST_Direct, domain, engngModel) );
74
76 int nc = min(2 * nroot, nroot + 8);
77 int nn = a.giveNumberOfColumns();
78 if ( nc > nn ) {
79 nc = nn;
80 }
81
82 ar.resize(nc, nc);
83 ar.zero();
84 br.resize(nc, nc);
85 br.zero();
86
87 //
88 // creation of initial iteration vectors
89 //
90 nc1 = nc - 1;
91
92 w.resize(nn);
93 w.zero();
94 d.resize(nc);
95 d.zero();
96 tt.resize(nn);
97 tt.zero();
98 rtolv.resize(nc);
99 rtolv.zero();
100 vec.resize(nc, nc);
101 vec.zero(); // eigen vectors of reduced problem
102
103 //
104 // create work arrays
105 //
106 r.resize(nn, nc);
107 r.zero();
108 eigv.resize(nc);
109 eigv.zero();
110
111 FloatArray h(nn);
112 for ( int i = 1; i <= nn; i++ ) {
113 h.at(i) = 1.0;
114 w.at(i) = b.at(i, i) / a.at(i, i);
115 }
116
117 b.times(h, tt);
118 r.setColumn(tt, 1);
119
120 for ( int j = 2; j <= nc; j++ ) {
121 double rt = 0.0;
122 for ( int i = 1; i <= nn; i++ ) {
123 if ( fabs( w.at(i) ) >= rt ) {
124 rt = fabs( w.at(i) );
125 ij = i;
126 }
127 }
128
129 tt.at(j) = ij;
130 w.at(ij) = 0.;
131 for ( int i = 1; i <= nn; i++ ) {
132 if ( i == ij ) {
133 h.at(i) = 1.0;
134 } else {
135 h.at(i) = 0.0;
136 }
137 }
138
139 b.times(h, tt);
140 r.setColumn(tt, j);
141 } // (r = z)
142
143# ifdef DETAILED_REPORT
144 OOFEM_LOG_INFO("SubspaceIteration :: solveYourselfAt: Degrees of freedom invoked by initial vectors :\n");
145 tt.printYourself();
146 OOFEM_LOG_INFO("SubspaceIteration :: solveYourselfAt: initial vectors for iteration:\n");
147 r.printYourself();
148# endif
149
150 //ish = 0;
151 a.factorized();
152 //
153 // start of iteration loop
154 //
155 for ( int nite = 0; ; ++nite ) { // label 100
156# ifdef DETAILED_REPORT
157 printf("SubspaceIteration :: solveYourselfAt: Iteration loop no. %d\n", nite);
158# endif
159 //
160 // compute projection ar and br of matrices a , b
161 //
162 for ( int j = 1; j <= nc; j++ ) {
163 f.beColumnOf(r, j);
164
165 solver->solve(a, f, tt);
166
167 for ( int i = j; i <= nc; i++ ) {
168 double art = 0.;
169 for ( int k = 1; k <= nn; k++ ) {
170 art += r.at(k, i) * tt.at(k);
171 }
172
173 ar.at(j, i) = art;
174 }
175
176 r.setColumn(tt, j); // (r = xbar)
177 }
178
179 ar.symmetrized(); // label 110
180#ifdef DETAILED_REPORT
181 OOFEM_LOG_INFO("SubspaceIteration :: solveYourselfAt: Printing projection matrix ar\n");
182 ar.printYourself();
183#endif
184 //
185 for ( int j = 1; j <= nc; j++ ) {
186 tt.beColumnOf(r, j);
187
188 b.times(tt, temp);
189 for ( int i = j; i <= nc; i++ ) {
190 double brt = 0.;
191 for ( int k = 1; k <= nn; k++ ) {
192 brt += r.at(k, i) * temp.at(k);
193 }
194
195 br.at(j, i) = brt;
196 } // label 180
197
198 r.setColumn(temp, j); // (r=zbar)
199 } // label 160
200
201 br.symmetrized();
202#ifdef DETAILED_REPORT
203 OOFEM_LOG_INFO("SubspaceIteration :: solveYourselfAt: Printing projection matrix br\n");
204 br.printYourself();
205#endif
206
207 //
208 // solution of reduced eigenvalue problem
209 //
210 mtd.solve(ar, br, eigv, vec);
211
212 // START EXPERIMENTAL
213#if 0
214 // solve the reduced problem by Inverse iteration
215 {
216 FloatMatrix x(nc,nc), z(nc,nc), zz(nc,nc), arinv;
217 FloatArray w(nc), ww(nc), tt(nc), t(nc);
218 double c;
219
220 // initial setting
221 for ( int i = 1;i <= nc; i++ ) {
222 ww.at(i)=1.0;
223 }
224
225
226 for ( int i = 1;i <= nc; i++ )
227 for ( int j = 1; j <= nc;j++ )
228 z.at(i,j)=1.0;
229
230 arinv.beInverseOf (ar);
231
232 for ( int i = 0;i < nitem; i++ ) {
233 // copy zz=z
234 zz = z;
235
236 // solve matrix equation K.X = M.X
237 x.beProductOf(arinv, z);
238 // evaluation of Rayleigh quotients
239 for ( int j = 1;j <= nc; j++ ) {
240 w.at(j) = 0.0;
241 for (k = 1; k<= nc; k++) w.at(j) += zz.at(k,j) * x.at(k,j);
242 }
243
244 z.beProductOf (br, x);
245
246 for ( int j = 1;j <= nc; j++ ) {
247 c = 0;
248 for ( int k = 1; k<= nc; k++ ) c += z.at(k,j) * x.at(k,j);
249 w.at(j) /= c;
250 }
251
252 // check convergence
253 int ac = 0;
254 for ( int j = 1;j <= nc; j++ ) {
255 if (fabs((ww.at(j)-w.at(j))/w.at(j))< rtol) ac++;
256 ww.at(j) = w.at(j);
257 }
258
259 //printf ("\n iterace cislo %d %d",i,ac);
260 //w.printYourself();
261
262 // Gramm-Schmidt ortogonalization
263 for ( int j = 1;j <= nc;j++ ) {
264 for ( int k = 1; k<= nc; k++ ) tt.at(k) = x.at(k,j);
265 t.beProductOf(br,tt) ;
266 for ( int ii = 1;ii < j; ii++ ) {
267 c = 0.0;
268 for ( int k = 1; k<= nc; k++ ) c += x.at(k,ii) * t.at(k);
269 for ( int k = 1; k<= nc; k++ ) x.at(k,j) -= x.at(k,ii) * c;
270 }
271 for ( int k = 1; k<= nc; k++) tt.at(k) = x.at(k,j);
272 t.beProductOf(br, tt);
273 c = 0.0;
274 for ( int k = 1; k<= nc; k++) c += x.at(k,j)*t.at(k);
275 for ( int k = 1; k<= nc; k++) x.at(k,j) /= sqrt(c);
276 }
277
278 if ( ac > nroot ) {
279 break;
280 }
281
282 // compute new approximation of Z
283 z.beProductOf(br,x);
284 }
285
286 eigv = w;
287 vec = x;
288 }
289#endif
290
291
292 //
293 // sorting eigenvalues according to their values
294 //
295 int is;
296 do {
297 is = 0; // label 350
298 for ( int i = 1; i <= nc1; i++ ) {
299 if ( fabs( eigv.at(i + 1) ) < fabs( eigv.at(i) ) ) {
300 is++;
301 std::swap(eigv.at(i), eigv.at(i + 1));
302 for ( int k = 1; k <= nc; k++ ) {
303 std::swap(vec.at(k, i), vec.at(k, i + 1));
304 }
305 }
306 } // label 360
307 } while ( is != 0 );
308
309# ifdef DETAILED_REPORT
310 OOFEM_LOG_INFO("SubspaceIteration :: solveYourselfAt: current eigen values of reduced problem \n");
311 eigv.printYourself();
312 OOFEM_LOG_INFO("SubspaceIteration :: solveYourselfAt: current eigen vectors of reduced problem \n");
313 vec.printYourself();
314# endif
315 //
316 // compute eigenvectors
317 //
318 for ( int i = 1; i <= nn; i++ ) { // label 375
319 for ( int j = 1; j <= nc; j++ ) {
320 tt.at(j) = r.at(i, j);
321 }
322
323 for ( int k = 1; k <= nc; k++ ) {
324 double rt = 0.;
325 for ( int j = 1; j <= nc; j++ ) {
326 rt += tt.at(j) * vec.at(j, k);
327 }
328
329 r.at(i, k) = rt;
330 }
331 } // label 420 (r = z)
332
333 //
334 // convergency check
335 //
336 for ( int i = 1; i <= nc; i++ ) {
337 double dif = ( eigv.at(i) - d.at(i) );
338 rtolv.at(i) = fabs( dif / eigv.at(i) );
339 }
340
341# ifdef DETAILED_REPORT
342 OOFEM_LOG_INFO("SubspaceIteration :: solveYourselfAt: Reached precision of eigenvalues:\n");
343 rtolv.printYourself();
344# endif
345 for ( int i = 1; i <= nroot; i++ ) {
346 if ( rtolv.at(i) > rtol ) {
347 goto label400;
348 }
349 }
350
351 OOFEM_LOG_INFO("SubspaceIteration :: solveYourselfAt: Convergence reached for RTOL=%20.15f\n", rtol);
352 status = CR_CONVERGED;
353 break;
354label400:
355 if ( nite >= nitem ) {
356 OOFEM_WARNING("SubspaceIteration :: solveYourselfAt: Convergence not reached in %d iteration - using current values", nitem);
357 status = CR_DIVERGED_ITS;
358 break;
359 }
360
361 d = eigv; // label 410 and 440
362
363 continue;
364 }
365
366
367 // compute eigenvectors
368 for ( int j = 1; j <= nc; j++ ) {
369 tt.beColumnOf(r, j);
370
372 r.setColumn(tt, j); // r = xbar
373 }
374
375 // one cad add a normalization of eigen-vectors here
376
377 // initialize original index locations
378 _r.resize(nn, nroot);
379 _eigv.resize(nroot);
380 for ( int i = 1; i <= nroot; i++ ) {
381 _eigv.at(i) = eigv.at(i);
382 for ( int j = 1; j <= nn; j++ ) {
383 _r.at(j, i) = r.at(j, i);
384 }
385 }
386
387 return status;
388}
389} // end namespace oofem
#define REGISTER_GeneralizedEigenValueSolver(class, type)
void resize(Index s)
Definition floatarray.C:94
double & at(Index i)
Definition floatarray.h:202
void beColumnOf(const FloatMatrix &mat, int col)
virtual void printYourself() const
Definition floatarray.C:762
void zero()
Zeroes all coefficients of receiver.
Definition floatarray.C:683
void beProductOf(const FloatMatrix &aMatrix, const FloatArray &anArray)
Definition floatarray.C:689
void resize(Index rows, Index cols)
Definition floatmatrix.C:79
void setColumn(const FloatArray &src, int c)
void zero()
Zeroes all coefficient of receiver.
bool beInverseOf(const FloatMatrix &src)
*Prints matrix to stdout Useful for debugging void printYourself() const
double at(std::size_t i, std::size_t j) const
ConvergedReason solve(FloatMatrix &K, FloatMatrix &M, FloatArray &w, FloatMatrix &x)
Definition gjacobi.C:54
Domain * domain
Pointer to domain.
Definition nummet.h:84
EngngModel * engngModel
Pointer to engineering model.
Definition nummet.h:86
SparseGeneralEigenValueSystemNM(Domain *d, EngngModel *m)
Constructor.
virtual double & at(int i, int j)=0
Returns coefficient at position (i,j).
virtual SparseMtrx * factorized()
Definition sparsemtrx.h:244
virtual FloatArray * backSubstitutionWith(FloatArray &y) const
Definition sparsemtrx.h:252
int giveNumberOfColumns() const
Returns number of columns of receiver.
Definition sparsemtrx.h:116
virtual void times(const FloatArray &x, FloatArray &answer) const
Definition sparsemtrx.h:137
int nitem
Max number of iterations.
Definition subspaceit.h:92
#define OOFEM_WARNING(...)
Definition error.h:80
#define OOFEM_ERROR(...)
Definition error.h:79
#define OOFEM_LOG_INFO(...)
Definition logger.h:143
FloatArrayF< N > min(const FloatArrayF< N > &a, const FloatArrayF< N > &b)
ClassFactory & GiveClassFactory()
@ CR_DIVERGED_ITS

This page is part of the OOFEM-3.0 documentation. Copyright Copyright (C) 1994-2025 Borek Patzak Bořek Patzák
Project e-mail: oofem@fsv.cvut.cz
Generated at for OOFEM by doxygen 1.15.0 written by Dimitri van Heesch, © 1997-2011