OPAL (Object Oriented Parallel Accelerator Library) 2024.2
OPAL
AmrMultiGrid.cpp
Go to the documentation of this file.
1//
2// Class AmrMultiGrid
3// Main class of the AMR Poisson multigrid solver.
4// It implements the multigrid solver described in https://doi.org/10.1016/j.cpc.2019.106912
5//
6// Copyright (c) 2017 - 2020, Matthias Frey, Paul Scherrer Institut, Villigen PSI, Switzerland
7// All rights reserved
8//
9// Implemented as part of the PhD thesis
10// "Precise Simulations of Multibunches in High Intensity Cyclotrons"
11//
12// This file is part of OPAL.
13//
14// OPAL is free software: you can redistribute it and/or modify
15// it under the terms of the GNU General Public License as published by
16// the Free Software Foundation, either version 3 of the License, or
17// (at your option) any later version.
18//
19// You should have received a copy of the GNU General Public License
20// along with OPAL. If not, see <https://www.gnu.org/licenses/>.
21//
22#include "AmrMultiGrid.h"
23
24#include <algorithm>
25#include <filesystem>
26#include <functional>
27#include <map>
28#include <numeric>
29
31#include "OPALconfig.h"
32#include "Physics/Units.h"
34#include "Utilities/Timer.h"
35#include "Utilities/Util.h"
36
37#include <AMReX_ParallelDescriptor.H>
38
39#if AMR_MG_WRITE
40 #include <iomanip>
41#endif
42
44 const std::string& bsolver,
45 const std::string& prec,
46 const bool& rebalance,
47 const std::string& reuse,
48 const std::string& bcx,
49 const std::string& bcy,
50 const std::string& bcz,
51 const std::string& smoother,
52 const std::size_t& nSweeps,
53 const std::string& interp,
54 const std::string& norm)
55 : AmrPoissonSolver<AmrBoxLib>(itsAmrObject_p)
56 , comm_mp( new comm_t( amrex::ParallelDescriptor::Communicator() ) )
57 , nIter_m(0)
58 , bIter_m(0)
59 , maxiter_m(100)
60 , nSweeps_m(nSweeps)
61 , mglevel_m(0)
62 , lbase_m(0)
63 , lfine_m(0)
64 , nlevel_m(1)
65 , nBcPoints_m(0)
66 , snorm_m(norm)
67 , eps_m(1.0e-10)
68 , verbose_m(false)
69 , fname_m(OpalData::getInstance()->getInputBasename() + std::string(".solver"))
70 , flag_m(std::ios::out)
71{
72#if AMR_MG_TIMER
73 this->initTimer_m();
74#endif
75
76 const Boundary bcs[AMREX_SPACEDIM] = {
77 this->convertToEnumBoundary_m(bcx),
78 this->convertToEnumBoundary_m(bcy),
80 };
81
82 this->initPhysicalBoundary_m(&bcs[0]);
83
85
86 norm_m = this->convertToEnumNorm_m(norm);
87
88 const Interpolater interpolater = this->convertToEnumInterpolater_m(interp);
89 this->initInterpolater_m(interpolater);
90
91 // interpolater for crse-fine-interface
93
94 // preconditioner
95 const Preconditioner precond = this->convertToEnumPreconditioner_m(prec);
96 this->initPrec_m(precond, reuse);
97
98 // base level solver
99 const BaseSolver solver = this->convertToEnumBaseSolver_m(bsolver);
100 this->initBaseSolver_m(solver, rebalance, reuse);
101
102 if (std::filesystem::exists(fname_m)) {
103 flag_m = std::ios::app;
104 INFOMSG("Appending solver information to existing file: " << fname_m << endl);
105 } else {
106 INFOMSG("Creating new file for solver information: " << fname_m << endl);
107 }
108}
109
110
114 unsigned short baseLevel,
115 unsigned short finestLevel,
116 bool prevAsGuess)
117{
118 lbase_m = baseLevel;
119 lfine_m = finestLevel;
120 nlevel_m = lfine_m - lbase_m + 1;
121
122 /* we cannot use the previous solution
123 * if we have to regrid (AmrPoissonSolver::hasToRegrid())
124 *
125 * regrid_m is set in AmrBoxlib::regrid()
126 */
127 bool reset = !prevAsGuess;
128
129 if ( this->regrid_m )
130 reset = true;
131
132 this->initLevels_m(rho, itsAmrObject_mp->Geom(), this->regrid_m);
133
134 // build all necessary matrices and vectors
135 this->setup_m(rho, phi, this->regrid_m);
136
137 this->initGuess_m(reset);
138
139 // actual solve
140 scalar_t error = this->iterate_m();
141
142 for (int lev = nlevel_m - 1; lev > -1; --lev) {
143 averageDown_m(lev);
144 }
145
146 // write efield to AMReX
147 this->computeEfield_m(efield);
148
149 // copy solution back
150 for (int lev = 0; lev < nlevel_m; ++lev) {
151 int ilev = lbase_m + lev;
152
153 phi[ilev]->setVal(0.0, phi[ilev]->nGrow());
154
155 this->trilinos2amrex_m(lev, 0, *phi[ilev], mglevel_m[lev]->phi_p);
156 }
157
158 if ( verbose_m )
159 this->writeSDDSData_m(error);
160
161 // we can now reset
162 this->regrid_m = false;
163}
164
165
166void AmrMultiGrid::setNumberOfSweeps(const std::size_t& nSweeps) {
167 nSweeps_m = nSweeps;
168}
169
170
171void AmrMultiGrid::setMaxNumberOfIterations(const std::size_t& maxiter) {
172 if ( maxiter < 1 )
173 throw OpalException("AmrMultiGrid::setMaxNumberOfIterations()",
174 "The max. number of iterations needs to be positive!");
175
176 maxiter_m = maxiter;
177}
178
179
181 return nIter_m;
182}
183
184
188
189
190void AmrMultiGrid::setVerbose(bool verbose) {
191 verbose_m = verbose;
192}
193
194
196 eps_m = eps;
197}
198
199
201{
202 // make sure it's reset
203 nBcPoints_m = 0;
204
205 for (unsigned int i = 0; i < AMREX_SPACEDIM; ++i) {
206 switch ( bc[i] ) {
209 break;
210 case Boundary::OPEN:
212 break;
215 break;
216 default:
217 throw OpalException("AmrMultiGrid::initPhysicalBoundary_m()",
218 "This type of boundary is not supported");
219 }
220 // we use the maximum in order to build matrices
221 go_t tmp = bc_m[i]->getNumberOfPoints();
222 if ( nBcPoints_m < tmp )
223 nBcPoints_m = tmp;
224 }
225}
226
227
228void AmrMultiGrid::initLevels_m(const amrex::Vector<AmrField_u>& rho,
229 const amrex::Vector<AmrGeometry_t>& geom,
230 bool regrid)
231{
232 if ( !regrid )
233 return;
234
235 // although we do a resize afterwards, we do this to be safe
236 for (int lev = nlevel_m; lev < (int)mglevel_m.size(); ++lev) {
237 mglevel_m[lev].reset(nullptr);
238 }
239
240 mglevel_m.resize(nlevel_m);
241
242 amrex::Periodicity period(AmrIntVect_t(D_DECL(0, 0, 0)));
243
244 AmrIntVect_t rr = AmrIntVect_t(D_DECL(2, 2, 2));
245
246 for (int lev = 0; lev < nlevel_m; ++lev) {
247 int ilev = lbase_m + lev;
248
249 // do not initialize base level every time
250 if (mglevel_m[lev] == nullptr || lev > lbase_m) {
251 mglevel_m[lev].reset(new AmrMultiGridLevel_t(itsAmrObject_mp->getMeshScaling(),
252 rho[ilev]->boxArray(),
253 rho[ilev]->DistributionMap(),
254 geom[ilev],
255 rr,
256 bc_m,
257 comm_mp));
258 } else {
259 mglevel_m[lev]->buildLevelMask();
260 }
261
262 mglevel_m[lev]->refmask.reset(
264 mglevel_m[lev]->dmap, 1, 2)
265 );
266 mglevel_m[lev]->refmask->setVal(AmrMultiGridLevel_t::Refined::NO, 2);
267 mglevel_m[lev]->refmask->FillBoundary(period);
268
269 amrex::BoxArray ba = mglevel_m[lev]->grids;
270 ba.coarsen(rr);
271 mglevel_m[lev]->crsemask.reset(
273 mglevel_m[lev]->dmap, 1, 2)
274 );
275 mglevel_m[lev]->crsemask->setVal(AmrMultiGridLevel_t::Refined::NO, 2);
276 }
277
278 for (int lev = 1; lev < nlevel_m; ++lev) {
279 mglevel_m[lev]->crsemask->setVal(AmrMultiGridLevel_t::Refined::YES, 0);
280
281 // used for boundary interpolation --> replaces expensive calls to isBoundary
282 mglevel_m[lev]->crsemask->setDomainBndry(AmrMultiGridLevel_t::Mask::PHYSBNDRY,
283 mglevel_m[lev-1]->geom); //FIXME: geometry of lev - 1
284 // really needed ?
285 mglevel_m[lev]->crsemask->FillBoundary(period);
286 }
287
288 /* to complete initialization we need to fill
289 * the mask of refinement
290 */
291 for (int lev = 0; lev < nlevel_m-1; ++lev) {
292 // get boxarray with refined cells
293 amrex::BoxArray ba = mglevel_m[lev]->grids;
294 ba.refine(rr);
295 ba = amrex::intersect(mglevel_m[lev+1]->grids, ba);
296 ba.coarsen(rr);
297
298 // refined cells
299 amrex::DistributionMapping dmap(ba, comm_mp->getSize());
300 AmrMultiGridLevel_t::mask_t refined(ba, dmap, 1, 0);
301 refined.setVal(AmrMultiGridLevel_t::Refined::YES);
302// refined.setDomainBndry(AmrMultiGridLevel_t::Mask::PHYSBNDRY, mglevel_m[lev]->geom);
303
304 // fill intersection with YES
305 mglevel_m[lev]->refmask->copy(refined, 0, 0, 1, 0, 2);
306
307 /* physical boundary cells will never be refined cells
308 * since they are ghost cells
309 */
310 mglevel_m[lev]->refmask->setDomainBndry(AmrMultiGridLevel_t::Mask::PHYSBNDRY,
311 mglevel_m[lev]->geom);
312
313 mglevel_m[lev]->refmask->FillBoundary(period);
314 }
315}
316
317
319 for (int lev = 0; lev < nlevel_m; ++lev) {
320 mglevel_m[lev]->refmask.reset(nullptr);
321 mglevel_m[lev]->crsemask.reset(nullptr);
322 mglevel_m[lev]->mask.reset(nullptr);
323 }
324}
325
326
328 if ( !reset )
329 return;
330
331 // reset
332 for (int lev = 0; lev < nlevel_m; ++lev)
333 mglevel_m[lev]->phi_p->putScalar(0.0);
334}
335
336
338
339 // initial error
340 std::vector<scalar_t> rhsNorms;
341 std::vector<scalar_t> resNorms;
342
343 this->initResidual_m(rhsNorms, resNorms);
344
345 std::for_each(rhsNorms.begin(), rhsNorms.end(),
346 [this](double& val){ val *= eps_m; });
347
348 nIter_m = 0;
349 bIter_m = 0;
350
351 while ( !isConverged_m(rhsNorms, resNorms) && nIter_m < maxiter_m ) {
352
354
355// /* in contrast to algorithm, we average down now
356// * --> potential is valid also on coarse covered
357// * cells
358// * --> however, it may take 1-2 iterations longer
359// */
360// for (int lev = nlevel_m - 1; lev > -1; --lev) {
361// averageDown_m(lev);
362// }
363
364 // update residual
365 for (int lev = 0; lev < nlevel_m; ++lev) {
366
367 this->residual_m(lev,
368 mglevel_m[lev]->residual_p,
369 mglevel_m[lev]->rho_p,
370 mglevel_m[lev]->phi_p);
371 }
372
373
374 for (lo_t lev = 0; lev < nlevel_m; ++lev)
375 resNorms[lev] = getLevelResidualNorm(lev);
376
377 ++nIter_m;
378
379#if AMR_MG_WRITE
380 this->writeResidualNorm_m();
381#endif
382
383 bIter_m += solver_mp->getNumIters();
384 }
385
386 return std::accumulate(resNorms.begin(),
387 resNorms.end(), 0.0,
388 std::plus<scalar_t>());
389}
390
391
392bool AmrMultiGrid::isConverged_m(std::vector<scalar_t>& rhsNorms,
393 std::vector<scalar_t>& resNorms)
394{
395 return std::equal(resNorms.begin(), resNorms.end(),
396 rhsNorms.begin(),
397 std::less<scalar_t>());
398}
399
400
402 Teuchos::RCP<vector_t>& r,
403 const Teuchos::RCP<vector_t>& b,
404 const Teuchos::RCP<vector_t>& x)
405{
406 /*
407 * r = b - A*x
408 */
409 if ( level < lfine_m ) {
410
411 vector_t fine2crse(mglevel_m[level]->Awf_p->getDomainMap(), true);
412
413 // get boundary for
414 if ( mglevel_m[level]->Bfine_p != Teuchos::null ) {
415 mglevel_m[level]->Bfine_p->apply(*mglevel_m[level+1]->phi_p, fine2crse);
416 }
417
418 // operation: fine2crse += A * x
419 mglevel_m[level]->Awf_p->apply(*x, fine2crse,
420 Teuchos::NO_TRANS,
421 scalar_t(1.0),
422 scalar_t(1.0));
423
424 if ( mglevel_m[level]->Bcrse_p != Teuchos::null ) {
425 // operation: fine2crse += B * phi^(l-1)
426 mglevel_m[level]->Bcrse_p->apply(*mglevel_m[level-1]->phi_p,
427 fine2crse,
428 Teuchos::NO_TRANS,
429 scalar_t(1.0),
430 scalar_t(1.0));
431 }
432
433 Teuchos::RCP<vector_t> uncov_Ax = Teuchos::rcp( new vector_t(mglevel_m[level]->map_p, true) );
434 mglevel_m[level]->UnCovered_p->apply(fine2crse, *uncov_Ax);
435
436 Teuchos::RCP<vector_t> uncov_b = Teuchos::rcp( new vector_t(mglevel_m[level]->map_p, true) );
437
438 mglevel_m[level]->UnCovered_p->apply(*b, *uncov_b);
439
440 // ONLY subtract coarse rho
441// mglevel_m[level]->residual_p->putScalar(0.0);
442
443 r->update(1.0, *uncov_b, -1.0, *uncov_Ax, 0.0);
444
445 } else {
446 /* finest level: Awf_p == Anf_p
447 */
448 Teuchos::RCP<vector_t> Ax = Teuchos::rcp( new vector_t(mglevel_m[level]->map_p, true) );
449 mglevel_m[level]->Anf_p->apply(*x, *Ax);
450
451 if ( mglevel_m[level]->Bcrse_p != Teuchos::null ) {
452 // operationr: Ax += B * phi^(l-1)
453 mglevel_m[level]->Bcrse_p->apply(*mglevel_m[level-1]->phi_p,
454 *Ax,
455 Teuchos::NO_TRANS,
456 scalar_t(1.0),
457 scalar_t(1.0));
458 }
459 r->update(1.0, *b, -1.0, *Ax, 0.0);
460 }
461}
462
463
464void AmrMultiGrid::relax_m(const lo_t& level) {
465
466 if ( level == lfine_m ) {
467
468 if ( level == lbase_m ) {
469 /* Anf_p == Awf_p
470 */
471 Teuchos::RCP<vector_t> Ax = Teuchos::rcp( new vector_t(mglevel_m[level]->map_p, true) );
472 mglevel_m[level]->Anf_p->apply(*mglevel_m[level]->phi_p, *Ax);
473 mglevel_m[level]->residual_p->update(1.0, *mglevel_m[level]->rho_p, -1.0, *Ax, 0.0);
474
475 } else {
476 this->residual_no_fine_m(level,
477 mglevel_m[level]->residual_p,
478 mglevel_m[level]->phi_p,
479 mglevel_m[level-1]->phi_p,
480 mglevel_m[level]->rho_p);
481 }
482 }
483
484 if ( level > 0 ) {
485 // phi^(l, save) = phi^(l)
486 Teuchos::RCP<vector_t> phi_save = Teuchos::rcp( new vector_t(mglevel_m[level]->map_p) );
487 Tpetra::deep_copy(*phi_save, *mglevel_m[level]->phi_p);
488
489 mglevel_m[level-1]->error_p->putScalar(0.0);
490
491 // smoothing
492 this->smooth_m(level,
493 mglevel_m[level]->error_p,
494 mglevel_m[level]->residual_p);
495
496
497 // phi = phi + e
498 mglevel_m[level]->phi_p->update(1.0, *mglevel_m[level]->error_p, 1.0);
499
500 /*
501 * restrict
502 */
503 this->restrict_m(level);
504
505 this->relax_m(level - 1);
506
507 /*
508 * prolongate / interpolate
509 */
510 this->prolongate_m(level);
511
512 // residual update
513 Teuchos::RCP<vector_t> tmp = Teuchos::rcp( new vector_t(mglevel_m[level]->map_p) );
514 this->residual_no_fine_m(level, tmp,
515 mglevel_m[level]->error_p,
516 mglevel_m[level-1]->error_p,
517 mglevel_m[level]->residual_p);
518
519 Tpetra::deep_copy(*mglevel_m[level]->residual_p, *tmp);
520
521 // delta error
522 Teuchos::RCP<vector_t> derror = Teuchos::rcp( new vector_t(mglevel_m[level]->map_p, true) );
523
524 // smoothing
525 this->smooth_m(level, derror, mglevel_m[level]->residual_p);
526
527 // e^(l) += de^(l)
528 mglevel_m[level]->error_p->update(1.0, *derror, 1.0);
529
530 // phi^(l) = phi^(l, save) + e^(l)
531 mglevel_m[level]->phi_p->update(1.0, *phi_save, 1.0, *mglevel_m[level]->error_p, 0.0);
532
533 } else {
534 // e = A^(-1)r
535#if AMR_MG_TIMER
536 IpplTimings::startTimer(bottomTimer_m);
537#endif
538
539 solver_mp->solve(mglevel_m[level]->error_p,
540 mglevel_m[level]->residual_p);
541
542#if AMR_MG_TIMER
543 IpplTimings::stopTimer(bottomTimer_m);
544#endif
545 // phi = phi + e
546 mglevel_m[level]->phi_p->update(1.0, *mglevel_m[level]->error_p, 1.0);
547 }
548}
549
550
552 Teuchos::RCP<vector_t>& result,
553 const Teuchos::RCP<vector_t>& rhs,
554 const Teuchos::RCP<vector_t>& crs_rhs,
555 const Teuchos::RCP<vector_t>& b)
556{
557 vector_t crse2fine(mglevel_m[level]->Anf_p->getDomainMap(), true);
558
559 // get boundary for
560 if ( mglevel_m[level]->Bcrse_p != Teuchos::null ) {
561 mglevel_m[level]->Bcrse_p->apply(*crs_rhs, crse2fine);
562 }
563
564 // operation: crse2fine = 1.0 * crse2fine + 1.0 * A^(l) * rhs
565 mglevel_m[level]->Anf_p->apply(*rhs, crse2fine,
566 Teuchos::NO_TRANS,
567 scalar_t(1.0),
568 scalar_t(1.0));
569
570 result->update(1.0, *b, -1.0, crse2fine, 0.0);
571}
572
573
574#if AMR_MG_WRITE
575void AmrMultiGrid::writeResidualNorm_m() {
576 scalar_t err = 0.0;
577
578 std::ofstream out;
579 if ( comm_mp->getRank() == 0 )
580 out.open("residual.dat", std::ios::app);
581
582 for (int lev = 0; lev < nlevel_m; ++lev) {
583 scalar_t tmp = evalNorm_m(mglevel_m[lev]->residual_p);
584
585 if ( comm_mp->getRank() == 0 )
586 out << std::setw(15) << std::right << tmp;
587 }
588
589 if ( comm_mp->getRank() == 0 )
590 out.close();
591}
592#endif
593
594
596AmrMultiGrid::evalNorm_m(const Teuchos::RCP<const vector_t>& x)
597{
598 scalar_t norm = 0.0;
599
600 switch ( norm_m ) {
601 case Norm::L1:
602 {
603 norm = x->norm1();
604 break;
605 }
606 case Norm::L2:
607 {
608 norm = x->norm2();
609 break;
610 }
611 case Norm::LINF:
612 {
613 norm = x->normInf();
614 break;
615 }
616 default:
617 throw OpalException("AmrMultiGrid::evalNorm_m()",
618 "This type of norm not suppported.");
619 }
620 return norm;
621}
622
623
624void AmrMultiGrid::initResidual_m(std::vector<scalar_t>& rhsNorms,
625 std::vector<scalar_t>& resNorms)
626{
627 rhsNorms.clear();
628 resNorms.clear();
629
630#if AMR_MG_WRITE
631 std::ofstream out;
632
633 if ( comm_mp->getRank() == 0) {
634 out.open("residual.dat", std::ios::out);
635
636 for (int lev = 0; lev < nlevel_m; ++lev)
637 out << std::setw(14) << std::right << "level" << lev;
638 out << std::endl;
639 }
640#endif
641
642 for (int lev = 0; lev < nlevel_m; ++lev) {
643 this->residual_m(lev,
644 mglevel_m[lev]->residual_p,
645 mglevel_m[lev]->rho_p,
646 mglevel_m[lev]->phi_p);
647
648 resNorms.push_back(evalNorm_m(mglevel_m[lev]->residual_p));
649
650#if AMR_MG_WRITE
651 if ( comm_mp->getRank() == 0 )
652 out << std::setw(15) << std::right << resNorms.back();
653#endif
654
655 rhsNorms.push_back(evalNorm_m(mglevel_m[lev]->rho_p));
656 }
657
658#if AMR_MG_WRITE
659 if ( comm_mp->getRank() == 0 )
660 out.close();
661#endif
662}
663
664
666#if AMR_MG_TIMER
667 IpplTimings::startTimer(efieldTimer_m);
668#endif
669 Teuchos::RCP<vector_t> efield_p = Teuchos::null;
670 for (int lev = nlevel_m - 1; lev > -1; --lev) {
671 int ilev = lbase_m + lev;
672
673 efield_p = Teuchos::rcp( new vector_t(mglevel_m[lev]->map_p, false) );
674
675
676 for (int d = 0; d < AMREX_SPACEDIM; ++d) {
677 efield[ilev][d]->setVal(0.0, efield[ilev][d]->nGrow());
678 efield_p->putScalar(0.0);
679 mglevel_m[lev]->G_p[d]->apply(*mglevel_m[lev]->phi_p, *efield_p);
680 this->trilinos2amrex_m(lev, 0, *efield[ilev][d], efield_p);
681 }
682 }
683#if AMR_MG_TIMER
684 IpplTimings::stopTimer(efieldTimer_m);
685#endif
686}
687
688
689void AmrMultiGrid::setup_m(const amrex::Vector<AmrField_u>& rho,
690 const amrex::Vector<AmrField_u>& phi,
691 const bool& matrices)
692{
693#if AMR_MG_TIMER
694 IpplTimings::startTimer(buildTimer_m);
695#endif
696
697 if ( lbase_m == lfine_m )
698 this->buildSingleLevel_m(rho, phi, matrices);
699 else
700 this->buildMultiLevel_m(rho, phi, matrices);
701
702 mglevel_m[lfine_m]->error_p->putScalar(0.0);
703
704 if ( matrices ) {
705 this->clearMasks_m();
706 // set the bottom solve operator
707 if (!solver_mp->hasOperator()) {
708 solver_mp->setOperator(mglevel_m[lbase_m]->Anf_p, mglevel_m[0].get());
709 }
710 }
711
712#if AMR_MG_TIMER
713 IpplTimings::stopTimer(buildTimer_m);
714#endif
715}
716
717
718void AmrMultiGrid::buildSingleLevel_m(const amrex::Vector<AmrField_u>& rho,
719 const amrex::Vector<AmrField_u>& phi,
720 const bool& matrices)
721{
722 this->open_m(lbase_m, matrices);
723
724 const scalar_t* invdx = mglevel_m[lbase_m]->invCellSize();
725
726 const scalar_t invdx2[] = {
727 D_DECL( invdx[0] * invdx[0],
728 invdx[1] * invdx[1],
729 invdx[2] * invdx[2] )
730 };
731
732 if (matrices) {
733 for (amrex::MFIter mfi(*mglevel_m[lbase_m]->mask, true);
734 mfi.isValid(); ++mfi)
735 {
736 const box_t& tbx = mfi.tilebox();
737 const basefab_t& mfab = (*mglevel_m[lbase_m]->mask)[mfi];
738 const farraybox_t& rhofab = (*rho[lbase_m])[mfi];
739 const farraybox_t& pfab = (*phi[lbase_m])[mfi];
740
741 const int* lo = tbx.loVect();
742 const int* hi = tbx.hiVect();
743
744 for (int i = lo[0]; i <= hi[0]; ++i) {
745 for (int j = lo[1]; j <= hi[1]; ++j) {
746#if AMREX_SPACEDIM == 3
747 for (int k = lo[2]; k <= hi[2]; ++k) {
748#endif
749 AmrIntVect_t iv(D_DECL(i, j, k));
750 go_t gidx = mglevel_m[lbase_m]->serialize(iv);
751
752 if (!solver_mp->hasOperator()) {
753 this->buildNoFinePoissonMatrix_m(lbase_m, gidx, iv, mfab, invdx2);
754 this->buildGradientMatrix_m(lbase_m, gidx, iv, mfab, invdx);
755 }
756
757 mglevel_m[lbase_m]->rho_p->replaceGlobalValue(gidx, rhofab(iv, 0));
758 mglevel_m[lbase_m]->phi_p->replaceGlobalValue(gidx, pfab(iv, 0));
759
760#if AMREX_SPACEDIM == 3
761 }
762#endif
763 }
764 }
765 }
766 } else {
768 }
769
770 this->close_m(lbase_m, matrices);
771
772 if (matrices) {
773 mglevel_m[lbase_m]->Awf_p = Teuchos::null;
774 mglevel_m[lbase_m]->UnCovered_p = Teuchos::null;
775 }
776}
777
778
779void AmrMultiGrid::buildMultiLevel_m(const amrex::Vector<AmrField_u>& rho,
780 const amrex::Vector<AmrField_u>& phi,
781 const bool& matrices)
782{
783 // the base level has no smoother --> nlevel_m - 1
784 if ( matrices ) {
785 // although we do a resize afterwards, we do this to be safe
786 for (int lev = nlevel_m-1; lev < (int)smoother_m.size(); ++lev) {
787 smoother_m[lev].reset();
788 }
789 smoother_m.resize(nlevel_m-1);
790 }
791
792 for (int lev = 0; lev < nlevel_m; ++lev) {
793 this->open_m(lev, matrices);
794
795 int ilev = lbase_m + lev;
796
797 // find all coarse cells that are covered by fine cells
798// AmrIntVect_t rr = mglevel_m[lev]->refinement();
799
800 const scalar_t* invdx = mglevel_m[lev]->invCellSize();
801
802 const scalar_t invdx2[] = {
803 D_DECL( invdx[0] * invdx[0],
804 invdx[1] * invdx[1],
805 invdx[2] * invdx[2] )
806 };
807
808 if ( matrices ) {
809 for (amrex::MFIter mfi(*mglevel_m[lev]->mask, true);
810 mfi.isValid(); ++mfi)
811 {
812 const box_t& tbx = mfi.tilebox();
813 const basefab_t& mfab = (*mglevel_m[lev]->mask)[mfi];
814 const basefab_t& rfab = (*mglevel_m[lev]->refmask)[mfi];
815 const basefab_t& cfab = (*mglevel_m[lev]->crsemask)[mfi];
816 const farraybox_t& rhofab = (*rho[ilev])[mfi];
817 const farraybox_t& pfab = (*phi[ilev])[mfi];
818
819 const int* lo = tbx.loVect();
820 const int* hi = tbx.hiVect();
821
822 for (int i = lo[0]; i <= hi[0]; ++i) {
823 int ii = i << 1;
824 for (int j = lo[1]; j <= hi[1]; ++j) {
825 int jj = j << 1;
826#if AMREX_SPACEDIM == 3
827 for (int k = lo[2]; k <= hi[2]; ++k) {
828 int kk = k << 1;
829#endif
830 AmrIntVect_t iv(D_DECL(i, j, k));
831 go_t gidx = mglevel_m[lev]->serialize(iv);
832
833 this->buildRestrictionMatrix_m(lev, gidx, iv,
834 D_DECL(ii, jj, kk), rfab);
835
836 this->buildInterpolationMatrix_m(lev, gidx, iv, cfab);
837
838 this->buildCrseBoundaryMatrix_m(lev, gidx, iv, mfab,
839 cfab, invdx2);
840
841 this->buildFineBoundaryMatrix_m(lev, gidx, iv,
842 mfab, rfab);
843
844 this->buildCompositePoissonMatrix_m(lev, gidx, iv, mfab,
845 rfab, invdx2);
846
847 if (lev > lbase_m || (lev == lbase_m && !solver_mp->hasOperator())) {
848 this->buildNoFinePoissonMatrix_m(lev, gidx, iv, mfab, invdx2);
849 this->buildGradientMatrix_m(lev, gidx, iv, mfab, invdx);
850 }
851
852 mglevel_m[lev]->rho_p->replaceGlobalValue(gidx, rhofab(iv, 0));
853 mglevel_m[lev]->phi_p->replaceGlobalValue(gidx, pfab(iv, 0));
854#if AMREX_SPACEDIM == 3
855 }
856#endif
857 }
858 }
859 }
860 } else {
861 for (lo_t lev = 0; lev < nlevel_m; ++lev) {
862 int ilev = lbase_m + lev;
863 this->buildDensityVector_m(lev, *rho[ilev]);
864 }
865 }
866
867 this->close_m(lev, matrices);
868
869 if ( matrices && lev > lbase_m ) {
870 smoother_m[lev-1].reset( new AmrSmoother(mglevel_m[lev]->Anf_p,
872 }
873 }
874}
875
876
877void AmrMultiGrid::open_m(const lo_t& level,
878 const bool& matrices)
879{
880 if ( matrices ) {
881
882 if ( level > lbase_m ) {
883
884 /*
885 * interpolation matrix
886 */
887
888 int nNeighbours = (nBcPoints_m + 1) * interp_mp->getNumberOfPoints();
889
890 mglevel_m[level]->I_p = Teuchos::rcp( new matrix_t(mglevel_m[level]->map_p,
891 nNeighbours) );
892
893 /*
894 * coarse boundary matrix
895 */
896
897 nNeighbours = 2 * AMREX_SPACEDIM * nBcPoints_m *
898 2 * AMREX_SPACEDIM * interface_mp->getNumberOfPoints();
899
900 mglevel_m[level]->Bcrse_p = Teuchos::rcp(
901 new matrix_t(mglevel_m[level]->map_p, nNeighbours) );
902
903 }
904
905
906 if ( level < lfine_m ) {
907
908 /*
909 * restriction matrix
910 */
911
912 // refinement 2
913 int nNeighbours = AMREX_D_TERM(2, * 2, * 2);
914
915 mglevel_m[level]->R_p = Teuchos::rcp(
916 new matrix_t(mglevel_m[level]->map_p, nNeighbours) );
917
918 /*
919 * fine boundary matrix
920 */
921
922 // refinement 2
923 nNeighbours = 2 * AMREX_SPACEDIM * AMREX_D_TERM(2, * 2, * 2);
924
925 mglevel_m[level]->Bfine_p = Teuchos::rcp(
926 new matrix_t(mglevel_m[level]->map_p, nNeighbours) );
927
928 }
929
930 /*
931 * no-fine Poisson matrix
932 */
933
934 int nPhysBoundary = 2 * AMREX_SPACEDIM * nBcPoints_m;
935
936 // number of internal stencil points
937 int nIntBoundary = AMREX_SPACEDIM * interface_mp->getNumberOfPoints();
938
939 int nEntries = (AMREX_SPACEDIM << 1) + 2 /* plus boundaries */ + nPhysBoundary + nIntBoundary;
940
941 if (level > lbase_m || (level == lbase_m && !solver_mp->hasOperator())) {
942 mglevel_m[level]->Anf_p = Teuchos::rcp(
943 new matrix_t(mglevel_m[level]->map_p, nEntries) );
944 }
945
946 /*
947 * with-fine / composite Poisson matrix
948 */
949 if ( lbase_m != lfine_m ) {
950 nEntries = (AMREX_SPACEDIM << 1) + 5 /* plus boundaries */ + nPhysBoundary + nIntBoundary;
951
952 mglevel_m[level]->Awf_p = Teuchos::rcp(
953 new matrix_t(mglevel_m[level]->map_p, nEntries) );
954
955 /*
956 * uncovered cells matrix
957 */
958 mglevel_m[level]->UnCovered_p = Teuchos::rcp(
959 new matrix_t(mglevel_m[level]->map_p, 1) );
960 }
961
962 /*
963 * gradient matrices
964 */
965 nEntries = 11;
966
967 if (level > lbase_m || (level == lbase_m && !solver_mp->hasOperator())) {
968 for (int d = 0; d < AMREX_SPACEDIM; ++d) {
969 mglevel_m[level]->G_p[d] = Teuchos::rcp(
970 new matrix_t(mglevel_m[level]->map_p, nEntries) );
971 }
972 }
973 }
974
975 mglevel_m[level]->rho_p = Teuchos::rcp(
976 new vector_t(mglevel_m[level]->map_p, false) );
977
978 if ( matrices ) {
979 mglevel_m[level]->phi_p = Teuchos::rcp(
980 new vector_t(mglevel_m[level]->map_p, false) );
981 }
982}
983
984
985void AmrMultiGrid::close_m(const lo_t& level,
986 const bool& matrices)
987{
988 if ( matrices ) {
989 if ( level > lbase_m ) {
990
991 mglevel_m[level]->I_p->fillComplete(mglevel_m[level-1]->map_p, // col map (domain map)
992 mglevel_m[level]->map_p); // row map (range map)
993
994 mglevel_m[level]->Bcrse_p->fillComplete(mglevel_m[level-1]->map_p, // col map
995 mglevel_m[level]->map_p); // row map
996 }
997
998 if ( level < lfine_m ) {
999
1000 mglevel_m[level]->R_p->fillComplete(mglevel_m[level+1]->map_p,
1001 mglevel_m[level]->map_p);
1002
1003 mglevel_m[level]->Bfine_p->fillComplete(mglevel_m[level+1]->map_p,
1004 mglevel_m[level]->map_p);
1005 }
1006
1007 if (level > lbase_m || (level == lbase_m && !solver_mp->hasOperator())) {
1008 mglevel_m[level]->Anf_p->fillComplete();
1009
1010 for (int d = 0; d < AMREX_SPACEDIM; ++d)
1011 mglevel_m[level]->G_p[d]->fillComplete();
1012 }
1013
1014 if ( lbase_m != lfine_m ) {
1015 mglevel_m[level]->Awf_p->fillComplete();
1016
1017 mglevel_m[level]->UnCovered_p->fillComplete();
1018 }
1019 }
1020}
1021
1022
1024 const go_t& gidx,
1025 const AmrIntVect_t& iv,
1026 const basefab_t& mfab,
1027 const scalar_t* invdx2)
1028{
1029 /*
1030 * Laplacian of "no fine"
1031 */
1032
1033 /*
1034 * 1D not supported
1035 * 2D --> 5 elements per row
1036 * 3D --> 7 elements per row
1037 */
1038
1039 umap_t map;
1040 indices_t indices;
1041 coefficients_t values;
1042
1043 /*
1044 * check neighbours in all directions (Laplacian stencil --> cross)
1045 */
1046 for (int d = 0; d < AMREX_SPACEDIM; ++d) {
1047 for (int shift = -1; shift <= 1; shift += 2) {
1048 AmrIntVect_t biv = iv;
1049 biv[d] += shift;
1050
1051 switch ( mfab(biv) )
1052 {
1054 case AmrMultiGridLevel_t::Mask::COVERED: // covered --> interior cell
1055 {
1056 map[mglevel_m[level]->serialize(biv)] += invdx2[d];
1057 break;
1058 }
1060 {
1061 // boundary cell
1062 // only level > 0 have this kind of boundary
1063#ifndef NDEBUG
1064 if ( level == lbase_m )
1065 throw OpalException("AmrMultiGrid::buildNoFinePoissonMatrix_m()",
1066 "Error in mask for level "
1067 + std::to_string(level) + "!");
1068#endif
1069 /* Dirichlet boundary conditions from coarser level.
1070 */
1071 interface_mp->fine(biv, map, invdx2[d], d, -shift,
1072 mglevel_m[level].get());
1073 break;
1074 }
1076 {
1077 // physical boundary cell
1078 mglevel_m[level]->applyBoundary(biv, d, map,
1079 invdx2[d] /*matrix coefficient*/);
1080 break;
1081 }
1082 default:
1083 throw OpalException("AmrMultiGrid::buildNoFinePoissonMatrix_m()",
1084 "Error in mask for level "
1085 + std::to_string(level) + "!");
1086 }
1087 }
1088 }
1089
1090 // check center
1091 map[gidx] += AMREX_D_TERM(- 2.0 * invdx2[0],
1092 - 2.0 * invdx2[1],
1093 - 2.0 * invdx2[2]);
1094
1095 this->map2vector_m(map, indices, values);
1096
1097 mglevel_m[level]->Anf_p->insertGlobalValues(gidx,
1098 indices.size(),
1099 &values[0],
1100 &indices[0]);
1101}
1102
1103
1105 const go_t& gidx,
1106 const AmrIntVect_t& iv,
1107 const basefab_t& mfab,
1108 const basefab_t& rfab,
1109 const scalar_t* invdx2)
1110{
1111 /*
1112 * Laplacian of "with fine"
1113 *
1114 * For the finest level: Awf == Anf
1115 */
1116 if ( rfab(iv) == AmrMultiGridLevel_t::Refined::YES ) //|| lbase_m != lfine_m )
1117 return;
1118 /*
1119 * Only cells that are not refined
1120 */
1121
1122 /*
1123 * 1D not supported by AmrMultiGrid
1124 * 2D --> 5 elements per row
1125 * 3D --> 7 elements per row
1126 */
1127
1128 umap_t map;
1129 indices_t indices;
1130 coefficients_t values;
1131
1132 /*
1133 * Only cells that are not refined
1134 */
1135
1136 /*
1137 * check neighbours in all directions (Laplacian stencil --> cross)
1138 */
1139 for (int d = 0; d < AMREX_SPACEDIM; ++d) {
1140 for (int shift = -1; shift <= 1; shift += 2) {
1141 AmrIntVect_t biv = iv;
1142 biv[d] += shift;
1143
1144 if ( rfab(biv) != AmrMultiGridLevel_t::Refined::YES )
1145 {
1146 /*
1147 * It can't be a refined cell!
1148 */
1149 switch ( mfab(biv) )
1150 {
1152 // covered --> interior cell
1154 {
1155 map[mglevel_m[level]->serialize(biv)] += invdx2[d];
1156 map[gidx] -= invdx2[d]; // add center once
1157 break;
1158 }
1160 {
1161 // boundary cell
1162 // only level > 0 have this kind of boundary
1163#ifndef NDEBUG
1164 if ( level == lbase_m )
1165 throw OpalException("AmrMultiGrid::buildCompositePoissonMatrix_m()",
1166 "Error in mask for level "
1167 + std::to_string(level) + "!");
1168#endif
1169
1170 /* We are on the fine side of the crse-fine interface
1171 * --> normal stencil --> no flux matching required
1172 * --> interpolation of fine ghost cell required
1173 * (used together with Bcrse)
1174 */
1175
1176 /* Dirichlet boundary conditions from coarser level.
1177 */
1178 interface_mp->fine(biv, map, invdx2[d], d, -shift,
1179 mglevel_m[level].get());
1180
1181 // add center once
1182 map[gidx] -= invdx2[d];
1183 break;
1184 }
1186 {
1187 // physical boundary cell
1188 mglevel_m[level]->applyBoundary(biv, d, map,
1189 invdx2[d] /*matrix coefficient*/);
1190
1191 // add center once
1192 map[gidx] -= invdx2[d];
1193 break;
1194 }
1195 default:
1196 throw OpalException("AmrMultiGrid::buildCompositePoissonMatrix_m()",
1197 "Error in mask for level "
1198 + std::to_string(level) + "!");
1199 }
1200 } else {
1201 /*
1202 * If neighbour cell is refined, we are on the coarse
1203 * side of the crse-fine interface --> flux matching
1204 * required --> interpolation of fine ghost cell
1205 * (used together with Bfine)
1206 */
1207
1208
1209 // flux matching, coarse part
1210
1211 /* 2D --> 2 fine cells to compute flux per coarse-fine-interace --> avg = 2
1212 * 3D --> 4 fine cells to compute flux per coarse-fine-interace --> avg = 4
1213 *
1214 * @precondition: refinement of 2
1215 */
1216 // top and bottom for all directions
1217 const scalar_t* invcdx = mglevel_m[level]->invCellSize();
1218 const scalar_t* invfdx = mglevel_m[level+1]->invCellSize();
1219 scalar_t invavg = AMREX_D_PICK(1.0, 0.5, 0.25);
1220 scalar_t value = - invavg * invcdx[d] * invfdx[d];
1221
1222 for (int d1 = 0; d1 < 2; ++d1) {
1223#if AMREX_SPACEDIM == 3
1224 for (int d2 = 0; d2 < 2; ++d2) {
1225#endif
1226
1227 /* in order to get a top iv --> needs to be odd value in "d"
1228 * in order to get a bottom iv --> needs to be even value in "d"
1229 */
1230 AmrIntVect_t fake(D_DECL(0, 0, 0));
1231
1232 fake[(d+1)%AMREX_SPACEDIM] = d1;
1233#if AMREX_SPACEDIM == 3
1234 fake[(d+2)%AMREX_SPACEDIM] = d2;
1235#endif
1236 interface_mp->coarse(iv, map, value, d, shift, rfab,
1237 fake, mglevel_m[level].get());
1238
1239#if AMREX_SPACEDIM == 3
1240 }
1241#endif
1242 }
1243 }
1244 }
1245 }
1246
1247 this->map2vector_m(map, indices, values);
1248
1249 mglevel_m[level]->Awf_p->insertGlobalValues(gidx,
1250 indices.size(),
1251 &values[0],
1252 &indices[0]);
1253
1254 scalar_t vv = 1.0;
1255 mglevel_m[level]->UnCovered_p->insertGlobalValues(gidx,
1256 1,
1257 &vv,
1258 &gidx);
1259}
1260
1261
1263 const go_t& gidx,
1264 const AmrIntVect_t& iv,
1265 D_DECL(const go_t& ii,
1266 const go_t& jj,
1267 const go_t& kk),
1268 const basefab_t& rfab)
1269{
1270 /*
1271 * x^(l) = R * x^(l+1)
1272 */
1273
1274 // finest level does not need to have a restriction matrix
1275 if ( rfab(iv) == AmrMultiGridLevel_t::Refined::NO || level == lfine_m )
1276 return;
1277
1278 /* Difficulty: If a fine cell belongs to another processor than the underlying
1279 * coarse cell, we get an error when filling the matrix since the
1280 * cell (--> global index) does not belong to the same processor.
1281 * Solution: Find all coarse cells that are covered by fine cells, thus,
1282 * the distributionmap is correct.
1283 *
1284 *
1285 */
1286 indices_t indices;
1287 indices.reserve(2 << (AMREX_SPACEDIM - 1));
1288 coefficients_t values;
1289 values.reserve(2 << (AMREX_SPACEDIM -1));
1290
1291 // neighbours
1292 for (int iref = 0; iref < 2; ++iref) {
1293 for (int jref = 0; jref < 2; ++jref) {
1294#if AMREX_SPACEDIM == 3
1295 for (int kref = 0; kref < 2; ++kref) {
1296#endif
1297 AmrIntVect_t riv(D_DECL(ii + iref, jj + jref, kk + kref));
1298
1299 indices.push_back( mglevel_m[level+1]->serialize(riv) );
1300 values.push_back( AMREX_D_PICK(0.5, 0.25, 0.125) );
1301#if AMREX_SPACEDIM == 3
1302 }
1303#endif
1304 }
1305 }
1306
1307 mglevel_m[level]->R_p->insertGlobalValues(gidx,
1308 indices.size(),
1309 &values[0],
1310 &indices[0]);
1311}
1312
1313
1315 const go_t& gidx,
1316 const AmrIntVect_t& iv,
1317 const basefab_t& cfab)
1318{
1319 /* crse: level - 1
1320 * fine (this): level
1321 */
1322
1323 /*
1324 * This does not include ghost cells
1325 * --> no boundaries
1326 *
1327 * x^(l) = I * x^(l-1)
1328 */
1329
1330 if ( level == lbase_m )
1331 return;
1332
1333 umap_t map;
1334 indices_t indices;
1335 coefficients_t values;
1336
1337 /*
1338 * we need boundary + indices from coarser level
1339 */
1340 interp_mp->stencil(iv, cfab, map, 1.0, mglevel_m[level-1].get());
1341
1342 this->map2vector_m(map, indices, values);
1343
1344 mglevel_m[level]->I_p->insertGlobalValues(gidx,
1345 indices.size(),
1346 &values[0],
1347 &indices[0]);
1348}
1349
1350
1352 const go_t& gidx,
1353 const AmrIntVect_t& iv,
1354 const basefab_t& mfab,
1355 const basefab_t& cfab,
1356 const scalar_t* invdx2)
1357{
1358 /*
1359 * fine (this): level
1360 * coarse: level - 1
1361 */
1362
1363 // the base level has only physical boundaries
1364 if ( level == lbase_m )
1365 return;
1366
1367 // iv is a fine cell
1368
1369 umap_t map;
1370 indices_t indices;
1371 coefficients_t values;
1372
1373 // check its neighbours to see if at crse-fine interface
1374 for (int d = 0; d < AMREX_SPACEDIM; ++d) {
1375 for (int shift = -1; shift <= 1; shift += 2) {
1376 // neighbour
1377 AmrIntVect_t niv = iv;
1378 niv[d] += shift;
1379
1380 if ( mfab(niv) == AmrMultiGridLevel_t::Mask::BNDRY )
1381 {
1382 // neighbour does not belong to fine grids
1383 // includes no cells at physical boundary
1384
1385 // coarse cell that is not refined
1386 AmrIntVect_t civ = iv;
1387 civ[d] += shift;
1388 civ.coarsen(mglevel_m[level]->refinement());
1389
1390 // we need boundary + indices from coarser level
1391 // we need normalization by mesh size squared (of fine cell size)
1392 // --> Laplacian for fine cell
1393 // (Dirichlet boundary for coarse --> fine)
1394 interface_mp->coarse(civ, map, invdx2[d], d, shift, cfab,
1395 iv, mglevel_m[level-1].get());
1396 }
1397#ifndef NDEBUG
1398 else if ( mfab(niv) == AmrMultiGridLevel_t::Mask::PHYSBNDRY ) {
1399 throw OpalException("AmrMultiGrid::buildCrseBoundaryMatrix_m()",
1400 "Fine meshes shouldn't be connected "
1401 "to physical (i.e. mesh) boundary!");
1402 }
1403#endif
1404 }
1405 }
1406
1407 this->map2vector_m(map, indices, values);
1408
1409 mglevel_m[level]->Bcrse_p->insertGlobalValues(gidx,
1410 indices.size(),
1411 &values[0],
1412 &indices[0]);
1413}
1414
1415
1417 const go_t& gidx,
1418 const AmrIntVect_t& iv,
1419 const basefab_t& mfab,
1420 const basefab_t& rfab)
1421{
1422 /* fine: level + 1
1423 * coarse (this): level
1424 */
1425
1426 // the finest level does not need data from a finer level
1427 if ( rfab(iv) == AmrMultiGridLevel_t::Refined::YES || level == lfine_m )
1428 return;
1429
1430 const scalar_t* invcdx = mglevel_m[level]->invCellSize();
1431 const scalar_t* invfdx = mglevel_m[level+1]->invCellSize();
1432
1433 // inverse of number of fine cell gradients
1434 scalar_t invavg = AMREX_D_PICK(1, 0.5, 0.25);
1435
1436 umap_t map;
1437 indices_t indices;
1438 coefficients_t values;
1439
1440 auto fill = [&](umap_t& map,
1441 D_DECL(int ii, int jj, int kk),
1442 int* begin, int* end, int d,
1443 const AmrIntVect_t& iv, int shift)
1444 {
1445 for (int iref = ii - begin[0]; iref <= ii + end[0]; ++iref) {
1446 for (int jref = jj - begin[1]; jref <= jj + end[1]; ++jref) {
1447#if AMREX_SPACEDIM == 3
1448 for (int kref = kk - begin[2]; kref <= kk + end[2]; ++kref) {
1449#endif
1450 /* Since all fine cells on the not-refined cell are
1451 * outside of the "domain" --> we need to interpolate
1452 */
1453 AmrIntVect_t riv(D_DECL(iref, jref, kref));
1454
1455 if ( (riv[d] >> 1) /*refinement*/ == iv[d] ) {
1456 /* the fine cell is on the coarse side --> fine
1457 * ghost cell --> we need to interpolate
1458 */
1459 scalar_t value = - invavg * invcdx[d] * invfdx[d];
1460
1461 interface_mp->fine(riv, map, value, d, shift,
1462 mglevel_m[level+1].get());
1463 } else {
1464 scalar_t value = invavg * invcdx[d] * invfdx[d];
1465 map[mglevel_m[level+1]->serialize(riv)] += value;
1466 }
1467#if AMREX_SPACEDIM == 3
1468 }
1469#endif
1470 }
1471 }
1472 };
1473
1474
1475 /*
1476 * iv is a coarse cell that got not refined
1477 *
1478 * --> check all neighbours to see if at crse-fine
1479 * interface
1480 */
1481 for (int d = 0; d < AMREX_SPACEDIM; ++d) {
1482 for (int shift = -1; shift <= 1; shift += 2) {
1483 // neighbour
1484 AmrIntVect_t covered = iv;
1485 covered[d] += shift;
1486
1487 if ( rfab(covered) == AmrMultiGridLevel_t::Refined::YES &&
1488 mfab(covered) != AmrMultiGridLevel_t::PHYSBNDRY )
1489 {
1490 // neighbour is covered by fine cells
1491
1492 /*
1493 * "shift" is the amount to a coarse cell that got refined
1494 * "d" is the direction to shift
1495 *
1496 * --> check all covered neighbour cells
1497 */
1498
1499 /* we need to iterate over correct fine cells. It depends
1500 * on the orientation of the interface
1501 */
1502 int begin[AMREX_SPACEDIM] = { D_DECL( int(d == 0), int(d == 1), int(d == 2) ) };
1503 int end[AMREX_SPACEDIM] = { D_DECL( int(d != 0), int(d != 1), int(d != 2) ) };
1504
1505 /*
1506 * neighbour cell got refined but is not on physical boundary
1507 * --> we are at a crse-fine interface
1508 *
1509 * we need now to find out which fine cells
1510 * are required to satisfy the flux matching
1511 * condition
1512 */
1513
1514 switch ( shift ) {
1515 case -1:
1516 {
1517 // --> interface is on the lower face
1518 int ii = iv[0] << 1; // refinemet in x
1519 int jj = iv[1] << 1; // refinemet in y
1520#if AMREX_SPACEDIM == 3
1521 int kk = iv[2] << 1; // refinemet in z
1522#endif
1523 // iterate over all fine cells at the interface
1524 // start with lower cells --> cover coarse neighbour
1525 // cell
1526 fill(map, D_DECL(ii, jj, kk), &begin[0], &end[0], d, iv, shift);
1527 break;
1528 }
1529 case 1:
1530 default:
1531 {
1532 // --> interface is on the upper face
1533 int ii = covered[0] << 1; // refinemet in x
1534 int jj = covered[1] << 1; // refinemet in y
1535#if AMREX_SPACEDIM == 3
1536 int kk = covered[2] << 1; // refinemet in z
1537#endif
1538 fill(map, D_DECL(ii, jj, kk), &begin[0], &end[0], d, iv, shift);
1539 break;
1540 }
1541 }
1542 }
1543 }
1544 }
1545
1546 this->map2vector_m(map, indices, values);
1547
1548 // iv: not covered coarse cell at crse-fine interface
1549
1550 mglevel_m[level]->Bfine_p->insertGlobalValues(gidx,
1551 indices.size(),
1552 &values[0],
1553 &indices[0]);
1554}
1555
1556
1558 const AmrField_t& rho)
1559{
1560 this->amrex2trilinos_m(level, 0, rho, mglevel_m[level]->rho_p);
1561}
1562
1563
1565 const AmrField_t& phi)
1566{
1567 this->amrex2trilinos_m(level, 0, phi, mglevel_m[level]->phi_p);
1568}
1569
1570
1572 const go_t& gidx,
1573 const AmrIntVect_t& iv,
1574 const basefab_t& mfab,
1575 const scalar_t* invdx)
1576{
1577 umap_t map;
1578 indices_t indices;
1579 coefficients_t values;
1580
1581 auto check = [&](const AmrIntVect_t& iv,
1582 const basefab_t& mfab,
1583 int dir,
1584 scalar_t shift)
1585 {
1586 switch ( mfab(iv) )
1587 {
1589 // interior cells
1591 // covered --> interior cell
1592 map[mglevel_m[level]->serialize(iv)] -= shift * 0.5 * invdx[dir];
1593 break;
1595 {
1596 // interior boundary cells --> only level > 0
1597#ifndef NDEBUG
1598 if ( level == lbase_m )
1599 throw OpalException("AmrMultiGrid::buildGradientMatrix_m()",
1600 "Error in mask for level "
1601 + std::to_string(level) + "!");
1602#endif
1603
1604 scalar_t value = - shift * 0.5 * invdx[dir];
1605
1606 // use 1st order Lagrange --> only cells of this level required
1607 AmrIntVect_t tmp = iv;
1608 // first fine cell on refined coarse cell (closer to interface)
1609 tmp[dir] -= shift;
1610 map[mglevel_m[level]->serialize(tmp)] += 2.0 * value;
1611
1612 // second fine cell on refined coarse cell (further away from interface)
1613 tmp[dir] -= shift;
1614 map[mglevel_m[level]->serialize(tmp)] -= value;
1615 break;
1616 }
1618 {
1619 // physical boundary cells
1620
1621 scalar_t value = - shift * 0.5 * invdx[dir];
1622
1623 mglevel_m[level]->applyBoundary(iv, dir,
1624 map, value);
1625 break;
1626 }
1627 default:
1628 break;
1629 }
1630 };
1631
1632 for (int d = 0; d < AMREX_SPACEDIM; ++d) {
1633 for (int shift = -1; shift <= 1; shift += 2) {
1634 AmrIntVect_t niv = iv;
1635 niv[d] += shift;
1636 check(niv, mfab, d, shift);
1637 }
1638
1639 this->map2vector_m(map, indices, values);
1640
1641 mglevel_m[level]->G_p[d]->insertGlobalValues(gidx,
1642 indices.size(),
1643 &values[0],
1644 &indices[0]);
1645 }
1646}
1647
1648
1650 const lo_t& comp,
1651 const AmrField_t& mf,
1652 Teuchos::RCP<vector_t>& mv)
1653{
1654 if ( mv.is_null() )
1655 mv = Teuchos::rcp( new vector_t(mglevel_m[level]->map_p, false) );
1656
1657 for (amrex::MFIter mfi(mf, true); mfi.isValid(); ++mfi) {
1658 const amrex::Box& tbx = mfi.tilebox();
1659 const amrex::FArrayBox& fab = mf[mfi];
1660
1661 const int* lo = tbx.loVect();
1662 const int* hi = tbx.hiVect();
1663
1664 for (int i = lo[0]; i <= hi[0]; ++i) {
1665 for (int j = lo[1]; j <= hi[1]; ++j) {
1666#if AMREX_SPACEDIM == 3
1667 for (int k = lo[2]; k <= hi[2]; ++k) {
1668#endif
1669 AmrIntVect_t iv(D_DECL(i, j, k));
1670
1671 go_t gidx = mglevel_m[level]->serialize(iv);
1672
1673 mv->replaceGlobalValue(gidx, fab(iv, comp));
1674#if AMREX_SPACEDIM == 3
1675 }
1676#endif
1677 }
1678 }
1679 }
1680}
1681
1682
1684 const lo_t& comp,
1685 AmrField_t& mf,
1686 const Teuchos::RCP<vector_t>& mv)
1687{
1688 Teuchos::ArrayRCP<const amr::scalar_t> data = mv->get1dView();
1689
1690 for (amrex::MFIter mfi(mf, true); mfi.isValid(); ++mfi) {
1691 const amrex::Box& tbx = mfi.tilebox();
1692 amrex::FArrayBox& fab = mf[mfi];
1693
1694 const int* lo = tbx.loVect();
1695 const int* hi = tbx.hiVect();
1696
1697 for (int i = lo[0]; i <= hi[0]; ++i) {
1698 for (int j = lo[1]; j <= hi[1]; ++j) {
1699#if AMREX_SPACEDIM == 3
1700 for (int k = lo[2]; k <= hi[2]; ++k) {
1701#endif
1702 AmrIntVect_t iv(D_DECL(i, j, k));
1703
1704 go_t gidx = mglevel_m[level]->serialize(iv);
1705 lo_t lidx = mglevel_m[level]->map_p->getLocalElement(gidx);
1706
1707 fab(iv, comp) = data[lidx];
1708 }
1709#if AMREX_SPACEDIM == 3
1710 }
1711#endif
1712 }
1713 }
1714}
1715
1716
1717inline
1719 coefficients_t& values)
1720{
1721 indices.clear();
1722 values.clear();
1723
1724 indices.reserve(map.size());
1725 values.reserve(map.size());
1726
1727 std::for_each(map.begin(), map.end(),
1728 [&](const std::pair<const go_t, scalar_t>& entry)
1729 {
1730#ifndef NDEBUG
1731 if ( entry.first < 0 ) {
1732 throw OpalException("AmrMultiGrid::map2vector_m()",
1733 "Negative matrix index!");
1734 }
1735#endif
1736
1737 indices.push_back(entry.first);
1738 values.push_back(entry.second);
1739 }
1740 );
1741
1742 map.clear();
1743}
1744
1745
1747 Teuchos::RCP<vector_t>& e,
1748 Teuchos::RCP<vector_t>& r)
1749{
1750#if AMR_MG_TIMER
1751 IpplTimings::startTimer(smoothTimer_m);
1752#endif
1753
1754 // base level has no smoother --> l - 1
1755 smoother_m[level-1]->smooth(e, r);
1756
1757#if AMR_MG_TIMER
1758 IpplTimings::stopTimer(smoothTimer_m);
1759#endif
1760}
1761
1762
1763void AmrMultiGrid::restrict_m(const lo_t& level) {
1764
1765#if AMR_MG_TIMER
1766 IpplTimings::startTimer(restrictTimer_m);
1767#endif
1768
1769 Teuchos::RCP<vector_t> tmp = Teuchos::rcp( new vector_t(mglevel_m[level]->map_p) );
1770
1771 this->residual_no_fine_m(level, tmp,
1772 mglevel_m[level]->error_p,
1773 mglevel_m[level-1]->error_p,
1774 mglevel_m[level]->residual_p);
1775
1776 mglevel_m[level-1]->residual_p->putScalar(0.0);
1777
1778 // average down: residual^(l-1) = R^(l) * tmp
1779 mglevel_m[level-1]->R_p->apply(*tmp, *mglevel_m[level-1]->residual_p);
1780
1781
1782 // composite matrix, i.e. matrix without covered cells
1783 // r^(l-1) = rho^(l-1) - A * phi^(l-1)
1784
1785 vector_t fine2crse(mglevel_m[level-1]->Awf_p->getDomainMap(), true);
1786
1787 // get boundary for
1788 mglevel_m[level-1]->Bfine_p->apply(*mglevel_m[level]->phi_p, fine2crse);
1789
1790 // operation: fine2coarse += A * phi
1791 mglevel_m[level-1]->Awf_p->apply(*mglevel_m[level-1]->phi_p,
1792 fine2crse, Teuchos::NO_TRANS,
1793 scalar_t(1.0), scalar_t(1.0));
1794
1795 if ( mglevel_m[level-1]->Bcrse_p != Teuchos::null ) {
1796 // operation: fine2coarse += B * phi
1797 mglevel_m[level-1]->Bcrse_p->apply(*mglevel_m[level-2]->phi_p,
1798 fine2crse, Teuchos::NO_TRANS,
1799 scalar_t(1.0), scalar_t(1.0));
1800 }
1801
1802 Teuchos::RCP<vector_t> uncoveredRho = Teuchos::rcp( new vector_t(mglevel_m[level-1]->map_p, true) );
1803
1804 mglevel_m[level-1]->UnCovered_p->apply(*mglevel_m[level-1]->rho_p, *uncoveredRho);
1805
1806
1807 //FIXME tmp2 not needed
1808 Teuchos::RCP<vector_t> tmp2 = Teuchos::rcp( new vector_t(mglevel_m[level-1]->map_p, true) );
1809 mglevel_m[level-1]->UnCovered_p->apply(fine2crse, *tmp2);
1810
1811 // ONLY subtract coarse rho
1812 mglevel_m[level-1]->residual_p->update(1.0, *uncoveredRho, -1.0, *tmp2, 1.0);
1813
1814#if AMR_MG_TIMER
1815 IpplTimings::stopTimer(restrictTimer_m);
1816#endif
1817}
1818
1819
1821#if AMR_MG_TIMER
1822 IpplTimings::startTimer(interpTimer_m);
1823#endif
1824 // interpolate error from l-1 to l
1825 // operation: e^(l) = 1.0 * e^(l) + 1.0 * I^(l) * e^(l-1)
1826 mglevel_m[level]->I_p->apply(*mglevel_m[level-1]->error_p,
1827 *mglevel_m[level]->error_p,
1828 Teuchos::NO_TRANS,
1829 scalar_t(1.0),
1830 scalar_t(1.0));
1831#if AMR_MG_TIMER
1832 IpplTimings::stopTimer(interpTimer_m);
1833#endif
1834}
1835
1836
1838
1839 if (level == lfine_m )
1840 return;
1841#if AMR_MG_TIMER
1842 IpplTimings::startTimer(averageTimer_m);
1843#endif
1844 Teuchos::RCP<vector_t> phicrse = Teuchos::rcp( new vector_t(mglevel_m[level]->map_p, false) );
1845
1846 // operation: phicrse = 0.0 * phicrse + 1.0 * R^(l) * phi^(l+1)
1847 mglevel_m[level]->R_p->apply(*mglevel_m[level+1]->phi_p, *phicrse);
1848
1849 Teuchos::RCP<vector_t> uncov_phi = Teuchos::rcp( new vector_t(mglevel_m[level]->map_p, true) );
1850
1851 mglevel_m[level]->UnCovered_p->apply(*mglevel_m[level]->phi_p, *uncov_phi);
1852
1853 mglevel_m[level]->phi_p->update(1.0, *phicrse, 1.0, *uncov_phi, 0.0);
1854#if AMR_MG_TIMER
1855 IpplTimings::stopTimer(averageTimer_m);
1856#endif
1857}
1858
1859
1861 switch ( interp ) {
1864 break;
1866 throw OpalException("AmrMultiGrid::initInterpolater_m()",
1867 "Not yet implemented.");
1870 break;
1871 default:
1872 throw OpalException("AmrMultiGrid::initInterpolater_m()",
1873 "No such interpolater available.");
1874 }
1875}
1876
1877
1879 switch ( interface ) {
1882 break;
1886 break;
1889 break;
1890 default:
1891 throw OpalException("AmrMultiGrid::initCrseFineInterp_m()",
1892 "No such interpolater for the coarse-fine interface available.");
1893 }
1894}
1895
1896
1898 const bool& rebalance,
1899 const std::string& reuse)
1900{
1901 switch ( solver ) {
1902 // Belos solvers
1904 solver_mp.reset( new BelosSolver_t("BICGSTAB", prec_mp) );
1905 break;
1906 case BaseSolver::MINRES:
1907 solver_mp.reset( new BelosSolver_t("MINRES", prec_mp) );
1908 break;
1909 case BaseSolver::PCPG:
1910 solver_mp.reset( new BelosSolver_t("PCPG", prec_mp) );
1911 break;
1912 case BaseSolver::CG:
1913 solver_mp.reset( new BelosSolver_t("Pseudoblock CG", prec_mp) );
1914 break;
1915 case BaseSolver::GMRES:
1916 solver_mp.reset( new BelosSolver_t("Pseudoblock GMRES", prec_mp) );
1917 break;
1919 solver_mp.reset( new BelosSolver_t("Stochastic CG", prec_mp) );
1920 break;
1922 solver_mp.reset( new BelosSolver_t("RCG", prec_mp) );
1923 break;
1925 solver_mp.reset( new BelosSolver_t("GCRODR", prec_mp) );
1926 break;
1927 // Amesos2 solvers
1928#ifdef HAVE_AMESOS2_KLU2
1929 case BaseSolver::KLU2:
1930 solver_mp.reset( new Amesos2Solver_t("klu2") );
1931 break;
1932#endif
1933#if HAVE_AMESOS2_SUPERLU
1934 case BaseSolver::SUPERLU:
1935 solver_mp.reset( new Amesos2Solver_t("superlu") );
1936 break;
1937#endif
1938#ifdef HAVE_AMESOS2_UMFPACK
1939 case BaseSolver::UMFPACK:
1940 solver_mp.reset( new Amesos2Solver_t("umfpack") );
1941 break;
1942#endif
1943#ifdef HAVE_AMESOS2_PARDISO_MKL
1944 case BaseSolver::PARDISO_MKL:
1945 solver_mp.reset( new Amesos2Solver_t("pardiso_mkl") );
1946 break;
1947#endif
1948#ifdef HAVE_AMESOS2_MUMPS
1949 case BaseSolver::MUMPS:
1950 solver_mp.reset( new Amesos2Solver_t("mumps") );
1951 break;
1952#endif
1953#ifdef HAVE_AMESOS2_LAPACK
1954 case BaseSolver::LAPACK:
1955 solver_mp.reset( new Amesos2Solver_t("lapack") );
1956 break;
1957#endif
1958 case BaseSolver::SA:
1959 {
1960 std::string muelu = MueLuSolver_t::convertToMueLuReuseOption(reuse);
1961 solver_mp.reset( new MueLuSolver_t(rebalance, muelu) );
1962 break;
1963 }
1964 default:
1965 throw OpalException("AmrMultiGrid::initBaseSolver_m()",
1966 "No such bottom solver available.");
1967 }
1968}
1969
1970
1972 const std::string& reuse)
1973{
1974 switch ( prec ) {
1975 case Preconditioner::ILUT:
1976 case Preconditioner::CHEBYSHEV:
1977 case Preconditioner::RILUK:
1978 case Preconditioner::JACOBI:
1979 case Preconditioner::BLOCK_JACOBI:
1980 case Preconditioner::GS:
1981 case Preconditioner::BLOCK_GS:
1982 prec_mp.reset( new Ifpack2Preconditioner_t(prec) );
1983 break;
1984 case Preconditioner::SA:
1985 {
1986 std::string muelu = MueLuPreconditioner_t::convertToMueLuReuseOption(reuse);
1987 prec_mp.reset( new MueLuPreconditioner_t(muelu) );
1988 break;
1989 }
1990 case Preconditioner::NONE:
1991 prec_mp.reset();
1992 break;
1993 default:
1994 throw OpalException("AmrMultiGrid::initPrec_m()",
1995 "No such preconditioner available.");
1996 }
1997}
1998
1999
2002 std::map<std::string, Boundary> map;
2003
2004 map["DIRICHLET"] = Boundary::DIRICHLET;
2005 map["OPEN"] = Boundary::OPEN;
2006 map["PERIODIC"] = Boundary::PERIODIC;
2007
2008 auto boundary = map.find(bc);
2009
2010 if ( boundary == map.end() )
2011 throw OpalException("AmrMultiGrid::convertToEnumBoundary_m()",
2012 "No boundary type '" + bc + "'.");
2013 return boundary->second;
2014}
2015
2018 std::map<std::string, Interpolater> map;
2019
2020 map["TRILINEAR"] = Interpolater::TRILINEAR;
2021 map["LAGRANGE"] = Interpolater::LAGRANGE;
2023
2024 auto interpolater = map.find(interp);
2025
2026 if ( interpolater == map.end() )
2027 throw OpalException("AmrMultiGrid::convertToEnumInterpolater_m()",
2028 "No interpolater '" + interp + "'.");
2029 return interpolater->second;
2030}
2031
2032
2034AmrMultiGrid::convertToEnumBaseSolver_m(const std::string& bsolver) {
2035 std::map<std::string, BaseSolver> map;
2036
2037 map["BICGSTAB"] = BaseSolver::BICGSTAB;
2038 map["MINRES"] = BaseSolver::MINRES;
2039 map["PCPG"] = BaseSolver::PCPG;
2040 map["CG"] = BaseSolver::CG;
2041 map["GMRES"] = BaseSolver::GMRES;
2042 map["STOCHASTIC_CG"] = BaseSolver::STOCHASTIC_CG;
2043 map["RECYCLING_CG"] = BaseSolver::RECYCLING_GMRES;
2044 map["RECYCLING_GMRES"] = BaseSolver::RECYCLING_GMRES;
2045#ifdef HAVE_AMESOS2_KLU2
2046 map["KLU2"] = BaseSolver::KLU2;
2047#endif
2048#ifdef HAVE_AMESOS2_SUPERLU
2049 map["SUPERLU"] = BaseSolver::SUPERLU;
2050#endif
2051#ifdef HAVE_AMESOS2_UMFPACK
2052 map["UMFPACK"] = BaseSolver::UMFPACK;
2053#endif
2054#ifdef HAVE_AMESOS2_PARDISO_MKL
2055 map["PARDISO_MKL"] = BaseSolver::PARDISO_MKL;
2056#endif
2057#ifdef HAVE_AMESOS2_MUMPS
2058 map["MUMPS"] = BaseSolver::MUMPS;
2059#endif
2060#ifdef HAVE_AMESOS2_LAPACK
2061 map["LAPACK"] = BaseSolver::LAPACK;
2062#endif
2063 map["SA"] = BaseSolver::SA;
2064
2065 auto solver = map.find(bsolver);
2066
2067 if ( solver == map.end() )
2068 throw OpalException("AmrMultiGrid::convertToEnumBaseSolver_m()",
2069 "No bottom solver '" + bsolver + "'.");
2070 return solver->second;
2071}
2072
2073
2076 std::map<std::string, Preconditioner> map;
2077
2078 map["NONE"] = Preconditioner::NONE;
2079
2081
2083
2084 auto precond = map.find(prec);
2085
2086 if ( precond == map.end() )
2087 throw OpalException("AmrMultiGrid::convertToEnumPreconditioner_m()",
2088 "No preconditioner '" + prec + "'.");
2089 return precond->second;
2090}
2091
2092
2094AmrMultiGrid::convertToEnumSmoother_m(const std::string& smoother) {
2095 return AmrSmoother::convertToEnumSmoother(smoother);
2096}
2097
2098
2100AmrMultiGrid::convertToEnumNorm_m(const std::string& norm) {
2101 std::map<std::string, Norm> map;
2102
2103 map["L1_NORM"] = Norm::L1;
2104 map["L2_NORM"] = Norm::L2;
2105 map["LINF_NORM"] = Norm::LINF;
2106
2107 auto n = map.find(norm);
2108
2109 if ( n == map.end() )
2110 throw OpalException("AmrMultiGrid::convertToEnumNorm_m()",
2111 "No norm '" + norm + "'.");
2112 return n->second;
2113}
2114
2115
2116void AmrMultiGrid::writeSDDSHeader_m(std::ofstream& outfile) {
2117 OPALTimer::Timer simtimer;
2118
2119 std::string dateStr(simtimer.date());
2120 std::string timeStr(simtimer.time());
2121 std::string indent(" ");
2122
2123 outfile << "SDDS1" << std::endl;
2124 outfile << "&description\n"
2125 << indent << "text=\"Solver statistics '" << OpalData::getInstance()->getInputFn()
2126 << "' " << dateStr << "" << timeStr << "\",\n"
2127 << indent << "contents=\"solver info\"\n"
2128 << "&end\n";
2129 outfile << "&parameter\n"
2130 << indent << "name=processors,\n"
2131 << indent << "type=long,\n"
2132 << indent << "description=\"Number of Cores used\"\n"
2133 << "&end\n";
2134 outfile << "&parameter\n"
2135 << indent << "name=revision,\n"
2136 << indent << "type=string,\n"
2137 << indent << "description=\"git revision of opal\"\n"
2138 << "&end\n";
2139 outfile << "&parameter\n"
2140 << indent << "name=flavor,\n"
2141 << indent << "type=string,\n"
2142 << indent << "description=\"OPAL flavor that wrote file\"\n"
2143 << "&end\n";
2144 outfile << "&column\n"
2145 << indent << "name=t,\n"
2146 << indent << "type=double,\n"
2147 << indent << "units=ns,\n"
2148 << indent << "description=\"1 Time\"\n"
2149 << "&end\n";
2150 outfile << "&column\n"
2151 << indent << "name=mg_iter,\n"
2152 << indent << "type=long,\n"
2153 << indent << "units=1,\n"
2154 << indent << "description=\"2 Number of Multigrid Iterations\"\n"
2155 << "&end\n";
2156 outfile << "&column\n"
2157 << indent << "name=bottom_iter,\n"
2158 << indent << "type=long,\n"
2159 << indent << "units=1,\n"
2160 << indent << "description=\"3 Total Number of Bottom Solver Iterations\"\n"
2161 << "&end\n";
2162 outfile << "&column\n"
2163 << indent << "name=regrid,\n"
2164 << indent << "type=bool,\n"
2165 << indent << "units=1,\n"
2166 << indent << "description=\"4 Regrid Step\"\n"
2167 << "&end\n";
2168 outfile << "&column\n"
2169 << indent << "name=" + snorm_m + ",\n"
2170 << indent << "type=double,\n"
2171 << indent << "units=1,\n"
2172 << indent << "description=\"5 Error\"\n"
2173 << "&end\n"
2174 << "&data\n"
2175 << indent << "mode=ascii,\n"
2176 << indent << "no_row_counts=1\n"
2177 << "&end\n"
2178 << comm_mp->getSize() << '\n'
2179 << OPAL_PROJECT_NAME << " " << OPAL_PROJECT_VERSION << " git rev. #" << Util::getGitRevision() << '\n'
2180 << (OpalData::getInstance()->isInOPALTMode()? "opal-t":
2181 (OpalData::getInstance()->isInOPALCyclMode()? "opal-cycl": "opal-map")) << std::endl;
2182}
2183
2184
2186#if AMR_MG_TIMER
2187 IpplTimings::startTimer(dumpTimer_m);
2188#endif
2189 unsigned int pwi = 10;
2190
2191 std::ofstream outfile;
2192
2193 if ( comm_mp->getRank() == 0 ) {
2194 outfile.open(fname_m.c_str(), flag_m);
2195 outfile.precision(15);
2196 outfile.setf(std::ios::scientific, std::ios::floatfield);
2197
2198 if ( flag_m == std::ios::out ) {
2199 flag_m = std::ios::app;
2200 writeSDDSHeader_m(outfile);
2201 }
2202
2203 outfile << itsAmrObject_mp->getT() * Units::s2ns << std::setw(pwi) << '\t' // 1
2204 << this->nIter_m << std::setw(pwi) << '\t' // 2
2205 << this->bIter_m << std::setw(pwi) << '\t' // 3
2206 << this->regrid_m << std::setw(pwi) << '\t' // 4
2207 << error << '\n'; // 5
2208 }
2209#if AMR_MG_TIMER
2210 IpplTimings::stopTimer(dumpTimer_m);
2211#endif
2212}
2213
2214
2215#if AMR_MG_TIMER
2216void AmrMultiGrid::initTimer_m() {
2217 buildTimer_m = IpplTimings::getTimer("AMR MG matrix setup");
2218 restrictTimer_m = IpplTimings::getTimer("AMR MG restrict");
2219 smoothTimer_m = IpplTimings::getTimer("AMR MG smooth");
2220 interpTimer_m = IpplTimings::getTimer("AMR MG prolongate");
2221 efieldTimer_m = IpplTimings::getTimer("AMR MG e-field");
2222 averageTimer_m = IpplTimings::getTimer("AMR MG average down");
2223 bottomTimer_m = IpplTimings::getTimer("AMR MG bottom-solver");
2224 dumpTimer_m = IpplTimings::getTimer("AMR MG dump");
2225}
2226#endif
2227
2228
2229double AmrMultiGrid::getXRangeMin(unsigned short level) {
2230 return itsAmrObject_mp->Geom(level).ProbLo(0);
2231}
2232
2233
2234double AmrMultiGrid::getXRangeMax(unsigned short level) {
2235 return itsAmrObject_mp->Geom(level).ProbHi(0);
2236}
2237
2238
2239double AmrMultiGrid::getYRangeMin(unsigned short level) {
2240 return itsAmrObject_mp->Geom(level).ProbLo(1);
2241}
2242
2243
2244double AmrMultiGrid::getYRangeMax(unsigned short level) {
2245 return itsAmrObject_mp->Geom(level).ProbHi(1);
2246}
2247
2248
2249double AmrMultiGrid::getZRangeMin(unsigned short level) {
2250 return itsAmrObject_mp->Geom(level).ProbLo(2);
2251}
2252
2253
2254double AmrMultiGrid::getZRangeMax(unsigned short level) {
2255 return itsAmrObject_mp->Geom(level).ProbHi(2);
2256}
2257
2258
2260 os << "* ********************* A M R M u l t i G r i d ********************** " << endl
2261 //FIXME
2262 << "* ******************************************************************** " << endl;
2263 return os;
2264}
PartBunchBase< T, Dim >::ConstIterator end(PartBunchBase< T, Dim > const &bunch)
PartBunchBase< T, Dim >::ConstIterator begin(PartBunchBase< T, Dim > const &bunch)
Inform & endl(Inform &inf)
Definition Inform.cpp:42
#define INFOMSG(msg)
Definition IpplInfo.h:348
void serialize(Param_t params, std::ostringstream &os)
serializes params using Boost archive
Definition MPIHelper.cpp:32
STL namespace.
constexpr double s2ns
Definition Units.h:44
std::string getGitRevision()
Definition Util.cpp:33
The global OPAL structure.
Definition OpalData.h:49
bool isInOPALTMode()
Definition OpalData.cpp:276
bool isInOPALCyclMode()
Definition OpalData.cpp:272
std::string getInputFn()
get opals input filename
Definition OpalData.cpp:670
static OpalData * getInstance()
Definition OpalData.cpp:196
void trilinos2amrex_m(const lo_t &level, const lo_t &comp, AmrField_t &mf, const Teuchos::RCP< vector_t > &mv)
void buildCrseBoundaryMatrix_m(const lo_t &level, const go_t &gidx, const AmrIntVect_t &iv, const basefab_t &mfab, const basefab_t &cfab, const scalar_t *invdx2)
amrex::BaseFab< int > basefab_t
void relax_m(const lo_t &level)
void initGuess_m(bool reset)
void close_m(const lo_t &level, const bool &matrices)
AmrMultiGrid(AmrBoxLib *itsAmrObject_p, const std::string &bsolver, const std::string &prec, const bool &rebalance, const std::string &reuse, const std::string &bcx, const std::string &bcy, const std::string &bcz, const std::string &smoother, const std::size_t &nSweeps, const std::string &interp, const std::string &norm)
double getZRangeMax(unsigned short level=0)
amr::Preconditioner Preconditioner
void restrict_m(const lo_t &level)
Boundary
Supported physical boundaries.
Boundary convertToEnumBoundary_m(const std::string &bc)
std::shared_ptr< bsolver_t > solver_mp
bottom solver
BaseSolver convertToEnumBaseSolver_m(const std::string &bsolver)
void setup_m(const amrex::Vector< AmrField_u > &rho, const amrex::Vector< AmrField_u > &phi, const bool &matrices=true)
double getXRangeMin(unsigned short level=0)
void initPhysicalBoundary_m(const Boundary *bc)
AmrMultiGridLevel_t::coefficients_t coefficients_t
Teuchos::RCP< comm_t > comm_mp
communicator
std::string fname_m
SDDS filename.
BaseSolver
Supported bottom solvers.
AmrMultiGridLevel_t::umap_t umap_t
amrex::FArrayBox farraybox_t
std::string snorm_m
norm for convergence criteria
std::size_t nSweeps_m
number of smoothing iterations
void map2vector_m(umap_t &map, indices_t &indices, coefficients_t &values)
void buildFineBoundaryMatrix_m(const lo_t &level, const go_t &gidx, const AmrIntVect_t &iv, const basefab_t &mfab, const basefab_t &rfab)
void writeSDDSHeader_m(std::ofstream &outfile)
std::vector< std::unique_ptr< AmrMultiGridLevel_t > > mglevel_m
container for levels
void amrex2trilinos_m(const lo_t &level, const lo_t &comp, const AmrField_t &mf, Teuchos::RCP< vector_t > &mv)
scalar_t getLevelResidualNorm(lo_t level)
scalar_t evalNorm_m(const Teuchos::RCP< const vector_t > &x)
void buildSingleLevel_m(const amrex::Vector< AmrField_u > &rho, const amrex::Vector< AmrField_u > &phi, const bool &matrices=true)
int lbase_m
base level (currently only 0 supported)
void initPrec_m(const Preconditioner &prec, const std::string &reuse)
void setVerbose(bool verbose)
void initInterpolater_m(const Interpolater &interp)
amr::comm_t comm_t
std::unique_ptr< AmrInterpolater< AmrMultiGridLevel_t > > interp_mp
interpolater without coarse-fine interface
amr::scalar_t scalar_t
scalar_t iterate_m()
void buildGradientMatrix_m(const lo_t &level, const go_t &gidx, const AmrIntVect_t &iv, const basefab_t &mfab, const scalar_t *invdx)
BelosBottomSolver< AmrMultiGridLevel_t > BelosSolver_t
void initCrseFineInterp_m(const Interpolater &interface)
std::shared_ptr< preconditioner_t > prec_mp
preconditioner for bottom solver
amrex::Box box_t
std::unique_ptr< AmrInterpolater< AmrMultiGridLevel_t > > interface_mp
interpolater for coarse-fine interface
Norm convertToEnumNorm_m(const std::string &norm)
Interpolater convertToEnumInterpolater_m(const std::string &interp)
amr::global_ordinal_t go_t
void averageDown_m(const lo_t &level)
amr::matrix_t matrix_t
Smoother convertToEnumSmoother_m(const std::string &smoother)
double getZRangeMin(unsigned short level=0)
amr::vector_t vector_t
Amesos2BottomSolver< AmrMultiGridLevel_t > Amesos2Solver_t
void buildPotentialVector_m(const lo_t &level, const AmrField_t &phi)
void residual_no_fine_m(const lo_t &level, Teuchos::RCP< vector_t > &result, const Teuchos::RCP< vector_t > &rhs, const Teuchos::RCP< vector_t > &crs_rhs, const Teuchos::RCP< vector_t > &b)
scalar_t eps_m
rhs scale for convergence
AmrMultiGridLevel_t::AmrField_t AmrField_t
std::size_t nIter_m
number of iterations till convergence
int lfine_m
fineste level
double getXRangeMax(unsigned short level=0)
Inform & print(Inform &os) const
AmrMultiGridLevel_t::AmrIntVect_t AmrIntVect_t
bool verbose_m
If true, a SDDS file is written.
std::size_t bIter_m
number of iterations of bottom solver
int nlevel_m
number of levelss
amr::local_ordinal_t lo_t
std::ios_base::openmode flag_m
std::ios::out or std::ios::app
void setNumberOfSweeps(const std::size_t &nSweeps)
void initResidual_m(std::vector< scalar_t > &rhsNorms, std::vector< scalar_t > &resNorms)
AmrSmoother::Smoother Smoother
void buildMultiLevel_m(const amrex::Vector< AmrField_u > &rho, const amrex::Vector< AmrField_u > &phi, const bool &matrices=true)
void open_m(const lo_t &level, const bool &matrices)
void buildDensityVector_m(const lo_t &level, const AmrField_t &rho)
MueLuPreconditioner< AmrMultiGridLevel_t > MueLuPreconditioner_t
AmrMultiGridLevel< matrix_t, vector_t > AmrMultiGridLevel_t
std::size_t getNumIters()
std::size_t maxiter_m
maximum number of iterations allowed
void buildNoFinePoissonMatrix_m(const lo_t &level, const go_t &gidx, const AmrIntVect_t &iv, const basefab_t &mfab, const scalar_t *invdx2)
double getYRangeMin(unsigned short level=0)
void writeSDDSData_m(const scalar_t &error)
int nBcPoints_m
maximum number of stencils points for BC
void setMaxNumberOfIterations(const std::size_t &maxiter)
void initBaseSolver_m(const BaseSolver &solver, const bool &rebalance, const std::string &reuse)
Smoother smootherType_m
type of smoother
Norm
Supported convergence criteria.
void initLevels_m(const amrex::Vector< AmrField_u > &rho, const amrex::Vector< AmrGeometry_t > &geom, bool regrid)
void smooth_m(const lo_t &level, Teuchos::RCP< vector_t > &e, Teuchos::RCP< vector_t > &r)
Ifpack2Preconditioner< AmrMultiGridLevel_t > Ifpack2Preconditioner_t
boundary_t bc_m[AMREX_SPACEDIM]
boundary conditions
void residual_m(const lo_t &level, Teuchos::RCP< vector_t > &r, const Teuchos::RCP< vector_t > &b, const Teuchos::RCP< vector_t > &x)
std::vector< std::shared_ptr< AmrSmoother > > smoother_m
error smoother
void buildInterpolationMatrix_m(const lo_t &level, const go_t &gidx, const AmrIntVect_t &iv, const basefab_t &cfab)
void computeEfield_m(AmrVectorFieldContainer_t &efield)
AmrMultiGridLevel_t::indices_t indices_t
Interpolater
Supported interpolaters for prolongation operation.
void prolongate_m(const lo_t &level)
bool isConverged_m(std::vector< scalar_t > &rhsNorms, std::vector< scalar_t > &resNorms)
Norm norm_m
norm for convergence criteria (l1, l2, linf)
double getYRangeMax(unsigned short level=0)
void setTolerance(const scalar_t &eps)
void solve(AmrScalarFieldContainer_t &rho, AmrScalarFieldContainer_t &phi, AmrVectorFieldContainer_t &efield, unsigned short baseLevel, unsigned short finestLevel, bool prevAsGuess=true)
MueLuBottomSolver< AmrMultiGridLevel_t > MueLuSolver_t
void buildRestrictionMatrix_m(const lo_t &level, const go_t &gidx, const AmrIntVect_t &iv, D_DECL(const go_t &ii, const go_t &jj, const go_t &kk), const basefab_t &rfab)
Preconditioner convertToEnumPreconditioner_m(const std::string &prec)
void buildCompositePoissonMatrix_m(const lo_t &level, const go_t &gidx, const AmrIntVect_t &iv, const basefab_t &mfab, const basefab_t &rfab, const scalar_t *invdx2)
amrex::FabArray< basefab_t > mask_t
static Smoother convertToEnumSmoother(const std::string &smoother)
static std::string convertToMueLuReuseOption(const std::string &reuse)
static std::string convertToMueLuReuseOption(const std::string &reuse)
AmrPoissonSolver(AmrBoxLib *itsAmrObject_p)
The base class for all OPAL exceptions.
std::string date() const
Return date.
std::string time() const
Return time.
static TimerRef getTimer(const char *nm)
static void stopTimer(TimerRef t)
static void startTimer(TimerRef t)