PLaSK library
Loading...
Searching...
No Matches
iterative_matrix.hpp
Go to the documentation of this file.
1/*
2 * This file is part of PLaSK (https://plask.app) by Photonics Group at TUL
3 * Copyright (c) 2022 Lodz University of Technology
4 *
5 * This program is free software: you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation, version 3.
8 *
9 * This program is distributed in the hope that it will be useful,
10 * but WITHOUT ANY WARRANTY; without even the implied warranty of
11 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12 * GNU General Public License for more details.
13 */
14#ifndef PLASK_COMMON_FEM_ITERATIVE_MATRIX_H
15#define PLASK_COMMON_FEM_ITERATIVE_MATRIX_H
16
17#include <algorithm>
18
19#include <nspcg/nspcg.hpp>
20
21#include "matrix.hpp"
22
23namespace plask {
24
95
97 typedef void (*NspcgFunc)(...);
98
99 protected:
100 int* icords;
101
103
104 const int nstore, ndim, mdim;
105
106 int ifact = 1;
107 int nw = 0, inw = 0;
108 double* wksp = nullptr;
109 int* iwksp = nullptr;
110 int kblsz = -1, nbl2d = -1;
111
113
114 virtual int get_maxnz() const { return mdim; }
115
116 public:
117 template <typename SolverT>
118 SparseMatrix(SolverT* solver, size_t rank, size_t size, size_t isiz)
121 params(&solver->iter_params),
122 nstore(2),
123 ndim(rank),
124 mdim(isiz) {}
125
126 template <typename SolverT>
130 params(&solver->iter_params),
131 nstore(4),
132 ndim(size),
133 mdim(size) {}
134
140
145
146 iparm.nstore = nstore;
147 iparm.itmax = params->maxit;
148 iparm.ipropa = 0;
149 iparm.ifact = (--ifact) ? 0 : 1;
150 if (ifact <= 0) ifact = params->nfact;
151 rparm.zeta = params->maxerr;
152
153 iparm.ns1 = params->ns1;
154 iparm.ns2 = params->ns2;
155 iparm.lvfill = params->lvfill;
156 iparm.ltrunc = params->ltrunc;
157 iparm.ndeg = params->ndeg;
158 rparm.omega = params->omega;
159
160 iparm.ns3 =
162 ? 40
163 : 0;
164 iparm.kblsz = kblsz;
165 iparm.nbl2d = nbl2d;
166
167 solver->writelog(LOG_DETAIL, "Iterating linear system");
168
169#ifdef NDEBUG
170 iparm.level = -1;
171#else
172 iparm.level = 3;
173#endif
174
175 int maxnz = get_maxnz();
176 size_t default_nw = 3 * rank + 2 * params->maxit + rank * maxnz + std::max(kblsz, 1);
177 size_t default_inw = maxnz + std::max(2 * int(rank), maxnz * maxnz + maxnz);
178
179 if (nw < default_nw) {
180 nw = default_nw;
183 iparm.ifact = 1; // we need to do factorization with new workspace
184 }
185
186 if (inw < default_inw) {
190 iparm.ifact = 1; // we need to do factorization with new workspace
191 }
192
193 // for (size_t r = 0; r < rank; ++r) {
194 // for (size_t c = 0; c < rank; ++c) {
195 // if (std::find(icords, A.icords+(ld+1), std::abs(int(r)-int(c))) == icords+(ld+1) )
196 // std::cout << " . ";
197 // else
198 // std::cout << str((*this))(r, c), "{:8.3f}") << " ";
199 // }
200 // std::cout << " " << str(B[r], "{:8.3f}") << std::endl;
201 // }
202
203 assert(B.size() == rank);
204
206 if (X.data() == nullptr || X.data() == B.data())
207 U.reset(B.size(), 1.);
208 else
209 U = X;
210 assert(U.size() == B.size());
211
212 int ier;
213
214 // TODO add choice of algorithms and predonditioners
215
217
221 throw BadInput(solver->getId(), "SOR oraccelerator must be used with SOR or LSOR preconditioner");
222 }
226 throw BadInput(solver->getId(), "SRCG accelerator must be used with SSOR or LSSOR preconditioner");
227 }
231 throw BadInput(solver->getId(), "SRSI accelerator must be used with SSOR or LSSOR preconditioner");
232 }
233 // clang-format off
235
236 switch (params->accelerator) {
257 };
258 // clang-format on
259
260 while (true) {
261 nspcg(precond_func, accel_func, ndim, mdim, rank, maxnz, data, icords, nullptr, nullptr, U.data(), nullptr, B.data(),
262 wksp, iwksp, nw, inw, iparm, rparm, ier);
263
264 // Increase workspace if needed
265 if (ier == -2 && nw) {
268 iparm.ifact = 1; // we need to do factorization with new workspace
269 } else if (ier == -3 && inw) {
272 iparm.ifact = 1; // we need to do factorization with new workspace
273 } else
274 break;
275 }
276
277 if (ier != 0) {
278 switch (ier) {
279 case -1: throw ComputationError(solver->getId(), "nonpositive matrix rank {}", rank);
280 case -2: throw ComputationError(solver->getId(), "insufficient real workspace ({} required)", nw);
281 case -3: throw ComputationError(solver->getId(), "insufficient integer workspace ({} required)", inw);
282 case -4: throw ComputationError(solver->getId(), "nonpositive diagonal element in stiffness matrix");
283 case -5: throw ComputationError(solver->getId(), "nonexistent diagonal element in stiffness matrix");
284 case -6: throw ComputationError(solver->getId(), "stiffness matrix A is not positive definite");
285 case -7: throw ComputationError(solver->getId(), "preconditioned matrix Q is not positive definite");
286 case -8: throw ComputationError(solver->getId(), "cannot permute stiffness matrix as requested");
287 case -9:
289 "Number of non-zero diagonals is not large enough to allow expansion of matrix");
290 case -10: throw ComputationError(solver->getId(), "inadmissible parameter encountered");
291 case -11: throw ComputationError(solver->getId(), "incorrect storage mode for block method");
292 case -12: throw ComputationError(solver->getId(), "zero pivot encountered in factorization");
293 case -13: throw ComputationError(solver->getId(), "breakdown in direction vector calculation");
294 case -14: throw ComputationError(solver->getId(), "breakdown in attempt to perform rotation");
295 case -15: throw ComputationError(solver->getId(), "breakdown in iterate calculation");
296 case -16: throw ComputationError(solver->getId(), "unimplemented combination of parameters");
297 case -18: throw ComputationError(solver->getId(), "unable to perform eigenvalue estimation");
298 case 1:
299 params->converged = false;
302 throw ComputationError(solver->getId(), "failed to converge in {} iterations (error {})", iparm.itmax,
303 rparm.zeta);
305 solver->writelog(LOG_WARNING, "Failed to converge in {} iterations (error {})", iparm.itmax,
306 rparm.zeta);
307 break;
309 solver->writelog(LOG_DETAIL, "Did not converge yen in {} iterations (error {})", iparm.itmax,
310 rparm.zeta);
311 break;
312 }
313 break;
314 case 2: solver->writelog(LOG_WARNING, "`maxerr` was too small, reset to {}", 3.55e-12); break;
315 case 3:
317 "NSPGS: `zbrent` failed to converge in the maximum number of {} iterations (signifies "
318 "difficulty in eigenvalue estimation)",
319 std::max(params->maxit, 50));
320 break;
321 case 4:
323 "NSPGS: In `zbrent`, f (a) and f (b) have the same sign (signifies difficulty in "
324 "eigenvalue estimation)");
325 break;
326 case 5: solver->writelog(LOG_DEBUG, "NSPGS: Negative pivot encountered in factorization"); break;
327 }
328 }
329 if (ier != 1) {
330 solver->writelog(LOG_DETAIL, "Converged after {} iterations (error {})", iparm.itmax, rparm.zeta);
331 params->converged = true;
332 }
333
334 params->iters = iparm.itmax;
335 params->err = rparm.zeta;
336
337 if (X.data() != U.data()) X = U;
338 }
339
341 std::fill(result.begin(), result.end(), 0.);
342 addmult(vector, result);
343 }
344};
345
353 template <typename SolverT>
355 icords[0] = 0;
356 icords[1] = 1;
357 icords[2] = major - 1;
358 icords[3] = major;
359 icords[4] = major + 1;
360 kblsz = major - 1;
361 nbl2d = major - 1;
362 }
363
371 template <typename SolverT>
372 SparseBandMatrix(SolverT* solver, size_t rank, size_t major, size_t minor) : SparseMatrix(solver, rank, 14 * rank, 14) {
373 icords[0] = 0;
374 icords[1] = 1;
375 icords[2] = minor - 1;
376 icords[3] = minor;
377 icords[4] = minor + 1;
378 icords[5] = major - minor - 1;
379 icords[6] = major - minor;
380 icords[7] = major - minor + 1;
381 icords[8] = major - 1;
382 icords[9] = major;
383 icords[10] = major + 1;
384 icords[11] = major + minor - 1;
385 icords[12] = major + minor;
386 icords[13] = major + minor + 1;
387 kblsz = minor - 1;
388 nbl2d = major - minor - 1;
389 }
390
397 template <typename SolverT>
398 SparseBandMatrix(SolverT* solver, size_t rank, std::initializer_list<int> bands)
400 size_t i = 0;
401 for (int band : bands) icords[i++] = band;
402 assert(icords[0] == 0);
403 }
404
411 double& operator()(size_t r, size_t c) override {
412 if (r == c) return data[r];
413 if (r < c) std::swap(r, c);
414 size_t i = std::find(icords, icords + mdim, r - c) - icords;
415 assert(i != mdim);
416 return data[c + rank * i];
417 }
418
420 for (size_t r = 0; r < rank; ++r) {
421 result[r] += data[r] * vector[r];
422 }
423 for (size_t d = 1; d < mdim; ++d) {
424 size_t sd = rank * d;
425 for (size_t r = 0; r < rank; ++r) {
426 size_t c = r + icords[d];
427 if (c >= rank) break;
428 result[r] += data[r + sd] * vector[c];
429 result[c] += data[r + sd] * vector[r];
430 }
431 }
432 }
433
434 protected:
458
459 public:
460 void setBC(DataVector<double>& B, size_t r, double val) override {
461 data[r] = 1.;
462 B[r] = val;
463 // above diagonal
464 for (ptrdiff_t i = mdim - 1; i > 0; --i) {
465 ptrdiff_t c = r - icords[i];
466 if (c >= 0) {
467 ptrdiff_t ii = c + rank * i;
468 assert(ii < size);
469 B[c] -= data[ii] * val;
470 data[ii] = 0.;
471 }
472 }
473 // below diagonal
474 for (ptrdiff_t i = 1; i < mdim; ++i) {
475 ptrdiff_t c = r + icords[i];
476 if (c < rank) {
477 size_t ii = r + rank * i;
478 assert(ii < size);
479 B[c] -= data[ii] * val;
480 data[ii] = 0.;
481 }
482 }
483 }
484};
485
487 int inz;
488 int* const ir;
489 int* const ic;
490
497 template <typename SolverT>
500 assert(maxnz >= rank);
501 for (size_t i = 0; i < rank; ++i) ir[i] = i + 1;
502 for (size_t i = 0; i < rank; ++i) ic[i] = i + 1;
505 }
506
513 double& operator()(size_t r, size_t c) override {
514 if (r == c) return data[r];
515 if (r > c) std::swap(r, c);
516 assert(inz < size);
517 ir[inz] = r + 1;
518 ic[inz] = c + 1;
519 return data[inz++];
520 }
521
522 void clear() override {
523 std::fill_n(data, size, 0.);
524 inz = rank;
525 }
526
528 for (size_t i = 0; i < rank; ++i) {
529 result[i] += data[i] * vector[i];
530 }
531 for (size_t i = rank; i < inz; ++i) {
532 size_t r = ir[i] - 1, c = ic[i] - 1;
533 result[r] += data[i] * vector[c];
534 result[c] += data[i] * vector[r];
535 }
536 }
537
538 protected:
540 switch (params->preconditioner) {
545 default: throw NotImplemented("preconditioner not implemented for non-rectangular or masked mesh");
546 };
547 assert(NULL);
548 }
549
550 int get_maxnz() const override { return inz; }
551
552 public:
553 void setBC(DataVector<double>& B, size_t r, double val) override {
554 data[r] = 1e32;
555 B[r] = val * 1e32;
556 }
557};
558
559} // namespace plask
560
561#endif // PLASK_COMMON_FEM_ITERATIVE_MATRIX_H