Cabana 0.8.0-dev
 
Loading...
Searching...
No Matches
Cabana_Grid_HypreStructuredSolver.hpp
Go to the documentation of this file.
1/****************************************************************************
2 * Copyright (c) 2018-2023 by the Cabana authors *
3 * All rights reserved. *
4 * *
5 * This file is part of the Cabana library. Cabana is distributed under a *
6 * BSD 3-clause license. For the licensing terms see the LICENSE file in *
7 * the top-level directory. *
8 * *
9 * SPDX-License-Identifier: BSD-3-Clause *
10 ****************************************************************************/
11
16#ifndef CABANA_GRID_HYPRESTRUCTUREDSOLVER_HPP
17#define CABANA_GRID_HYPRESTRUCTUREDSOLVER_HPP
18
19#include <Cabana_Grid_Array.hpp>
21#include <Cabana_Grid_Hypre.hpp>
24#include <Cabana_Grid_Types.hpp>
25
26#include <HYPRE_config.h>
27#include <HYPRE_struct_ls.h>
28#include <HYPRE_struct_mv.h>
29
30#include <Kokkos_Core.hpp>
31#include <Kokkos_Profiling_ScopedRegion.hpp>
32
33#include <array>
34#include <memory>
35#include <numeric>
36#include <sstream>
37#include <string>
38#include <type_traits>
39#include <vector>
40
41namespace Cabana
42{
43namespace Grid
44{
45//---------------------------------------------------------------------------//
47template <class Scalar, class EntityType, class MemorySpace>
49{
50 public:
52 using entity_type = EntityType;
54 using memory_space = MemorySpace;
56 using value_type = Scalar;
59 "HYPRE not compatible with solver memory space" );
60
67 template <class ArrayLayout_t>
68 HypreStructuredSolver( const ArrayLayout_t& layout,
69 const bool is_preconditioner = false )
70 : _comm( layout.localGrid()->globalGrid().comm() )
71 , _is_preconditioner( is_preconditioner )
72 {
74 "Must use an array layout" );
75 static_assert(
76 std::is_same<typename ArrayLayout_t::entity_type,
77 entity_type>::value,
78 "Array layout entity type mush match solver entity type" );
79
80 // Spatial dimension.
81 const std::size_t num_space_dim = ArrayLayout_t::num_space_dim;
82
83 // Only create data structures if this is not a preconditioner.
84 if ( !_is_preconditioner )
85 {
86 // Create the grid.
87 auto error = HYPRE_StructGridCreate( _comm, num_space_dim, &_grid );
88 checkHypreError( error );
89
90 // Get the global index space spanned by the local grid on this
91 // rank. Note that the upper bound is not a bound but rather the
92 // last index as this is what Hypre wants. Note that we reordered
93 // this to KJI from IJK to be consistent with HYPRE ordering. By
94 // setting up the grid like this, HYPRE will then want layout-right
95 // data indexed as (i,j,k) or (i,j,k,l) which will allow us to
96 // directly use Kokkos::deep_copy to move data between arrays and
97 // HYPRE data structures.
98 auto global_space = layout.indexSpace( Own(), Global() );
99 _lower.resize( num_space_dim );
100 _upper.resize( num_space_dim );
101 for ( std::size_t d = 0; d < num_space_dim; ++d )
102 {
103 _lower[d] = static_cast<HYPRE_Int>(
104 global_space.min( num_space_dim - d - 1 ) );
105 _upper[d] = static_cast<HYPRE_Int>(
106 global_space.max( num_space_dim - d - 1 ) - 1 );
107 }
108 error = HYPRE_StructGridSetExtents( _grid, _lower.data(),
109 _upper.data() );
110 checkHypreError( error );
111
112 // Get periodicity. Note we invert the order of this to KJI as well.
113 const auto& global_grid = layout.localGrid()->globalGrid();
114 HYPRE_Int periodic[num_space_dim];
115 for ( std::size_t d = 0; d < num_space_dim; ++d )
116 periodic[num_space_dim - 1 - d] =
117 global_grid.isPeriodic( d )
118 ? layout.localGrid()->globalGrid().globalNumEntity(
119 EntityType(), d )
120 : 0;
121 error = HYPRE_StructGridSetPeriodic( _grid, periodic );
122 checkHypreError( error );
123
124 // Assemble the grid.
125 error = HYPRE_StructGridAssemble( _grid );
126 checkHypreError( error );
127
128 // Allocate LHS and RHS vectors and initialize to zero. Note that we
129 // are fixing the views under these vectors to layout-right.
130 std::array<long, num_space_dim> reorder_size;
131 for ( std::size_t d = 0; d < num_space_dim; ++d )
132 {
133 reorder_size[d] = global_space.extent( d );
134 }
135 IndexSpace<num_space_dim> reorder_space( reorder_size );
136 auto vector_values =
138 "vector_values", reorder_space );
139 Kokkos::deep_copy( vector_values, 0.0 );
140
141 error = HYPRE_StructVectorCreate( _comm, _grid, &_b );
142 checkHypreError( error );
143 error = HYPRE_StructVectorInitialize( _b );
144 checkHypreError( error );
145 error = HYPRE_StructVectorSetBoxValues(
146 _b, _lower.data(), _upper.data(), vector_values.data() );
147 checkHypreError( error );
148 error = HYPRE_StructVectorAssemble( _b );
149 checkHypreError( error );
150
151 error = HYPRE_StructVectorCreate( _comm, _grid, &_x );
152 checkHypreError( error );
153 error = HYPRE_StructVectorInitialize( _x );
154 checkHypreError( error );
155 error = HYPRE_StructVectorSetBoxValues(
156 _x, _lower.data(), _upper.data(), vector_values.data() );
157 checkHypreError( error );
158 error = HYPRE_StructVectorAssemble( _x );
159 checkHypreError( error );
160 }
161 }
162
163 // Destructor.
164 virtual ~HypreStructuredSolver()
165 {
166 // We only make data if this is not a preconditioner.
167 if ( !_is_preconditioner )
168 {
169 HYPRE_StructVectorDestroy( _x );
170 HYPRE_StructVectorDestroy( _b );
171 HYPRE_StructMatrixDestroy( _A );
172 HYPRE_StructStencilDestroy( _stencil );
173 HYPRE_StructGridDestroy( _grid );
174 }
175 }
176
178 bool isPreconditioner() const { return _is_preconditioner; }
179
188 template <std::size_t NumSpaceDim>
189 void
190 setMatrixStencil( const std::vector<std::array<int, NumSpaceDim>>& stencil,
191 const bool is_symmetric = false )
192 {
193 // This function is only valid for non-preconditioners.
194 if ( _is_preconditioner )
195 throw std::logic_error(
196 "Cannot call setMatrixStencil() on preconditioners" );
197
198 // Create the stencil.
199 _stencil_size = stencil.size();
200 auto error =
201 HYPRE_StructStencilCreate( NumSpaceDim, _stencil_size, &_stencil );
202 checkHypreError( error );
203 std::array<HYPRE_Int, NumSpaceDim> offset;
204 for ( unsigned n = 0; n < stencil.size(); ++n )
205 {
206 for ( std::size_t d = 0; d < NumSpaceDim; ++d )
207 offset[d] = stencil[n][d];
208 error = HYPRE_StructStencilSetElement( _stencil, n, offset.data() );
209 checkHypreError( error );
210 }
211
212 // Create the matrix object. Must be done after the stencil is setup
213 error = HYPRE_StructMatrixCreate( _comm, _grid, _stencil, &_A );
214 checkHypreError( error );
215 error = HYPRE_StructMatrixSetSymmetric( _A, is_symmetric );
216 checkHypreError( error );
217 error = HYPRE_StructMatrixInitialize( _A );
218 checkHypreError( error );
219 }
220
229 template <class Array_t>
230 void setMatrixValues( const Array_t& values )
231 {
232 static_assert( is_array<Array_t>::value, "Must use an array" );
233 static_assert(
234 std::is_same<typename Array_t::entity_type, entity_type>::value,
235 "Array entity type mush match solver entity type" );
236 static_assert(
237 std::is_same<typename Array_t::memory_space, MemorySpace>::value,
238 "Array memory space and solver memory space are different." );
239
240 static_assert(
241 std::is_same<typename Array_t::value_type, value_type>::value,
242 "Array value type and solver value type are different." );
243
244 // This function is only valid for non-preconditioners.
245 if ( _is_preconditioner )
246 throw std::logic_error(
247 "Cannot call setMatrixValues() on preconditioners" );
248
249 if ( values.layout()->dofsPerEntity() !=
250 static_cast<int>( _stencil_size ) )
251 throw std::runtime_error(
252 "Number of matrix values does not match stencil size" );
253
254 // Spatial dimension.
255 const std::size_t num_space_dim = Array_t::num_space_dim;
256
257 // Copy the matrix entries into HYPRE. The HYPRE layout is fixed as
258 // layout-right.
259 auto owned_space = values.layout()->indexSpace( Own(), Local() );
260 std::array<long, num_space_dim + 1> reorder_size;
261 for ( std::size_t d = 0; d < num_space_dim; ++d )
262 {
263 reorder_size[d] = owned_space.extent( d );
264 }
265 reorder_size.back() = _stencil_size;
266 IndexSpace<num_space_dim + 1> reorder_space( reorder_size );
267 auto a_values =
269 "a_values", reorder_space );
270 auto values_subv = createSubview( values.view(), owned_space );
271 Kokkos::deep_copy( a_values, values_subv );
272
273 // Insert values into the HYPRE matrix.
274 std::vector<HYPRE_Int> indices( _stencil_size );
275 std::iota( indices.begin(), indices.end(), 0 );
276 auto error = HYPRE_StructMatrixSetBoxValues(
277 _A, _lower.data(), _upper.data(), indices.size(), indices.data(),
278 a_values.data() );
279 checkHypreError( error );
280 error = HYPRE_StructMatrixAssemble( _A );
281 checkHypreError( error );
282 }
283
288 void printMatrix( const char* prefix )
289 {
290 HYPRE_StructMatrixPrint( prefix, _A, 0 );
291 }
292
297 void printLHS( const char* prefix )
298 {
299 HYPRE_StructVectorPrint( prefix, _x, 0 );
300 }
301
306 void printRHS( const char* prefix )
307 {
308 HYPRE_StructVectorPrint( prefix, _b, 0 );
309 }
310
312 void setTolerance( const double tol ) { this->setToleranceImpl( tol ); }
313
315 void setMaxIter( const int max_iter ) { this->setMaxIterImpl( max_iter ); }
316
318 void setPrintLevel( const int print_level )
319 {
320 this->setPrintLevelImpl( print_level );
321 }
322
324 void
326 Scalar, EntityType, MemorySpace>>& preconditioner )
327 {
328 // This function is only valid for non-preconditioners.
329 if ( _is_preconditioner )
330 throw std::logic_error(
331 "Cannot call setPreconditioner() on a preconditioner" );
332
333 // Only a preconditioner can be used as a preconditioner.
334 if ( !preconditioner->isPreconditioner() )
335 throw std::logic_error( "Not a preconditioner" );
336
337 _preconditioner = preconditioner;
338 this->setPreconditionerImpl( *_preconditioner );
339 }
340
342 void setup()
343 {
344 // This function is only valid for non-preconditioners.
345 if ( _is_preconditioner )
346 throw std::logic_error( "Cannot call setup() on preconditioners" );
347
348 // FIXME: appears to be a memory issue in the call to this function
349 this->setupImpl();
350 }
351
357 template <class Array_t>
358 void solve( const Array_t& b, Array_t& x )
359 {
360 Kokkos::Profiling::ScopedRegion region(
361 "Cabana::Grid::HypreStructuredSolver::solve" );
362
363 static_assert( is_array<Array_t>::value, "Must use an array" );
364 static_assert(
365 std::is_same<typename Array_t::entity_type, entity_type>::value,
366 "Array entity type mush match solver entity type" );
367 static_assert(
368 std::is_same<typename Array_t::memory_space, MemorySpace>::value,
369 "Array memory space and solver memory space are different." );
370
371 static_assert(
372 std::is_same<typename Array_t::value_type, value_type>::value,
373 "Array value type and solver value type are different." );
374
375 // This function is only valid for non-preconditioners.
376 if ( _is_preconditioner )
377 throw std::logic_error( "Cannot call solve() on preconditioners" );
378
379 if ( b.layout()->dofsPerEntity() != 1 ||
380 x.layout()->dofsPerEntity() != 1 )
381 throw std::runtime_error(
382 "Structured solver only for scalar fields" );
383
384 // Spatial dimension.
385 const std::size_t num_space_dim = Array_t::num_space_dim;
386
387 // Copy the RHS into HYPRE. The HYPRE layout is fixed as layout-right.
388 auto owned_space = b.layout()->indexSpace( Own(), Local() );
389 std::array<long, num_space_dim + 1> reorder_size;
390 for ( std::size_t d = 0; d < num_space_dim; ++d )
391 {
392 reorder_size[d] = owned_space.extent( d );
393 }
394 reorder_size.back() = 1;
395 IndexSpace<num_space_dim + 1> reorder_space( reorder_size );
396 auto vector_values =
398 "vector_values", reorder_space );
399 auto b_subv = createSubview( b.view(), owned_space );
400 Kokkos::deep_copy( vector_values, b_subv );
401
402 // Insert b values into the HYPRE vector.
403 auto error = HYPRE_StructVectorSetBoxValues(
404 _b, _lower.data(), _upper.data(), vector_values.data() );
405 checkHypreError( error );
406 error = HYPRE_StructVectorAssemble( _b );
407 checkHypreError( error );
408
409 // Solve the problem
410 this->solveImpl();
411
412 // Extract the solution from the LHS
413 error = HYPRE_StructVectorGetBoxValues(
414 _x, _lower.data(), _upper.data(), vector_values.data() );
415 checkHypreError( error );
416
417 // Copy the HYPRE solution to the LHS.
418 auto x_subv = createSubview( x.view(), owned_space );
419 Kokkos::deep_copy( x_subv, vector_values );
420 }
421
423 int getNumIter() { return this->getNumIterImpl(); }
424
427 {
429 }
430
432 virtual HYPRE_StructSolver getHypreSolver() const = 0;
434 virtual HYPRE_PtrToStructSolverFcn getHypreSetupFunction() const = 0;
436 virtual HYPRE_PtrToStructSolverFcn getHypreSolveFunction() const = 0;
437
438 protected:
440 virtual void setToleranceImpl( const double tol ) = 0;
441
443 virtual void setMaxIterImpl( const int max_iter ) = 0;
444
446 virtual void setPrintLevelImpl( const int print_level ) = 0;
447
449 virtual void setupImpl() = 0;
450
452 virtual void solveImpl() = 0;
453
455 virtual int getNumIterImpl() = 0;
456
459
463 preconditioner ) = 0;
464
466 void checkHypreError( const int error ) const
467 {
468 if ( error > 0 )
469 {
470 char error_msg[256];
471 HYPRE_DescribeError( error, error_msg );
472 std::stringstream out;
473 out << "HYPRE structured solver error: ";
474 out << error << " " << error_msg;
475 HYPRE_ClearError( error );
476 throw std::runtime_error( out.str() );
477 }
478 }
479
480 protected:
482 HYPRE_StructMatrix _A;
484 HYPRE_StructVector _b;
486 HYPRE_StructVector _x;
487
488 private:
489 MPI_Comm _comm;
490 bool _is_preconditioner;
491 HYPRE_StructGrid _grid;
492 std::vector<HYPRE_Int> _lower;
493 std::vector<HYPRE_Int> _upper;
494 HYPRE_StructStencil _stencil;
495 unsigned _stencil_size;
496 std::shared_ptr<HypreStructuredSolver<Scalar, EntityType, MemorySpace>>
497 _preconditioner;
498};
499
500//---------------------------------------------------------------------------//
502template <class Scalar, class EntityType, class MemorySpace>
504 : public HypreStructuredSolver<Scalar, EntityType, MemorySpace>
505{
506 public:
510 template <class ArrayLayout_t>
511 HypreStructPCG( const ArrayLayout_t& layout,
512 const bool is_preconditioner = false )
513 : base_type( layout, is_preconditioner )
514 {
515 if ( is_preconditioner )
516 throw std::logic_error(
517 "HYPRE PCG cannot be used as a preconditioner" );
518
519 auto error = HYPRE_StructPCGCreate(
520 layout.localGrid()->globalGrid().comm(), &_solver );
521 this->checkHypreError( error );
522
523 HYPRE_StructPCGSetTwoNorm( _solver, 1 );
524 }
525
526 ~HypreStructPCG() { HYPRE_StructPCGDestroy( _solver ); }
527
528 // PCG SETTINGS
529
531 void setAbsoluteTol( const double tol )
532 {
533 auto error = HYPRE_StructPCGSetAbsoluteTol( _solver, tol );
534 this->checkHypreError( error );
535 }
536
539 void setRelChange( const int rel_change )
540 {
541 auto error = HYPRE_StructPCGSetRelChange( _solver, rel_change );
542 this->checkHypreError( error );
543 }
544
546 void setLogging( const int logging )
547 {
548 auto error = HYPRE_StructPCGSetLogging( _solver, logging );
549 this->checkHypreError( error );
550 }
551
552 HYPRE_StructSolver getHypreSolver() const override { return _solver; }
553 HYPRE_PtrToStructSolverFcn getHypreSetupFunction() const override
554 {
555 return HYPRE_StructPCGSetup;
556 }
557 HYPRE_PtrToStructSolverFcn getHypreSolveFunction() const override
558 {
559 return HYPRE_StructPCGSolve;
560 }
561
562 protected:
563 void setToleranceImpl( const double tol ) override
564 {
565 auto error = HYPRE_StructPCGSetTol( _solver, tol );
566 this->checkHypreError( error );
567 }
568
569 void setMaxIterImpl( const int max_iter ) override
570 {
571 auto error = HYPRE_StructPCGSetMaxIter( _solver, max_iter );
572 this->checkHypreError( error );
573 }
574
575 void setPrintLevelImpl( const int print_level ) override
576 {
577 auto error = HYPRE_StructPCGSetPrintLevel( _solver, print_level );
578 this->checkHypreError( error );
579 }
580
581 void setupImpl() override
582 {
583 auto error = HYPRE_StructPCGSetup( _solver, _A, _b, _x );
584 this->checkHypreError( error );
585 }
586
587 void solveImpl() override
588 {
589 auto error = HYPRE_StructPCGSolve( _solver, _A, _b, _x );
590 this->checkHypreError( error );
591 }
592
593 int getNumIterImpl() override
594 {
595 HYPRE_Int num_iter;
596 auto error = HYPRE_StructPCGGetNumIterations( _solver, &num_iter );
597 this->checkHypreError( error );
598 return num_iter;
599 }
600
602 {
603 HYPRE_Real norm;
604 auto error =
605 HYPRE_StructPCGGetFinalRelativeResidualNorm( _solver, &norm );
606 this->checkHypreError( error );
607 return norm;
608 }
609
612 preconditioner ) override
613 {
614 auto error = HYPRE_StructPCGSetPrecond(
615 _solver, preconditioner.getHypreSolveFunction(),
616 preconditioner.getHypreSetupFunction(),
617 preconditioner.getHypreSolver() );
618 this->checkHypreError( error );
619 }
620
621 private:
622 HYPRE_StructSolver _solver;
623 using base_type::_A;
624 using base_type::_b;
625 using base_type::_x;
626};
627
628//---------------------------------------------------------------------------//
630template <class Scalar, class EntityType, class MemorySpace>
632 : public HypreStructuredSolver<Scalar, EntityType, MemorySpace>
633{
634 public:
638 template <class ArrayLayout_t>
639 HypreStructGMRES( const ArrayLayout_t& layout,
640 const bool is_preconditioner = false )
641 : base_type( layout, is_preconditioner )
642 {
643 if ( is_preconditioner )
644 throw std::logic_error(
645 "HYPRE GMRES cannot be used as a preconditioner" );
646
647 auto error = HYPRE_StructGMRESCreate(
648 layout.localGrid()->globalGrid().comm(), &_solver );
649 this->checkHypreError( error );
650 }
651
652 ~HypreStructGMRES() { HYPRE_StructGMRESDestroy( _solver ); }
653
654 // GMRES SETTINGS
655
657 void setAbsoluteTol( const double tol )
658 {
659 auto error = HYPRE_StructGMRESSetAbsoluteTol( _solver, tol );
660 this->checkHypreError( error );
661 }
662
664 void setKDim( const int k_dim )
665 {
666 auto error = HYPRE_StructGMRESSetKDim( _solver, k_dim );
667 this->checkHypreError( error );
668 }
669
671 void setLogging( const int logging )
672 {
673 auto error = HYPRE_StructGMRESSetLogging( _solver, logging );
674 this->checkHypreError( error );
675 }
676
677 HYPRE_StructSolver getHypreSolver() const override { return _solver; }
678 HYPRE_PtrToStructSolverFcn getHypreSetupFunction() const override
679 {
680 return HYPRE_StructGMRESSetup;
681 }
682 HYPRE_PtrToStructSolverFcn getHypreSolveFunction() const override
683 {
684 return HYPRE_StructGMRESSolve;
685 }
686
687 protected:
688 void setToleranceImpl( const double tol ) override
689 {
690 auto error = HYPRE_StructGMRESSetTol( _solver, tol );
691 this->checkHypreError( error );
692 }
693
694 void setMaxIterImpl( const int max_iter ) override
695 {
696 auto error = HYPRE_StructGMRESSetMaxIter( _solver, max_iter );
697 this->checkHypreError( error );
698 }
699
700 void setPrintLevelImpl( const int print_level ) override
701 {
702 auto error = HYPRE_StructGMRESSetPrintLevel( _solver, print_level );
703 this->checkHypreError( error );
704 }
705
706 void setupImpl() override
707 {
708 auto error = HYPRE_StructGMRESSetup( _solver, _A, _b, _x );
709 this->checkHypreError( error );
710 }
711
712 void solveImpl() override
713 {
714 auto error = HYPRE_StructGMRESSolve( _solver, _A, _b, _x );
715 this->checkHypreError( error );
716 }
717
718 int getNumIterImpl() override
719 {
720 HYPRE_Int num_iter;
721 auto error = HYPRE_StructGMRESGetNumIterations( _solver, &num_iter );
722 this->checkHypreError( error );
723 return num_iter;
724 }
725
727 {
728 HYPRE_Real norm;
729 auto error =
730 HYPRE_StructGMRESGetFinalRelativeResidualNorm( _solver, &norm );
731 this->checkHypreError( error );
732 return norm;
733 }
734
737 preconditioner ) override
738 {
739 auto error = HYPRE_StructGMRESSetPrecond(
740 _solver, preconditioner.getHypreSolveFunction(),
741 preconditioner.getHypreSetupFunction(),
742 preconditioner.getHypreSolver() );
743 this->checkHypreError( error );
744 }
745
746 private:
747 HYPRE_StructSolver _solver;
748 using base_type::_A;
749 using base_type::_b;
750 using base_type::_x;
751};
752
753//---------------------------------------------------------------------------//
755template <class Scalar, class EntityType, class MemorySpace>
757 : public HypreStructuredSolver<Scalar, EntityType, MemorySpace>
758{
759 public:
763 template <class ArrayLayout_t>
764 HypreStructBiCGSTAB( const ArrayLayout_t& layout,
765 const bool is_preconditioner = false )
766 : base_type( layout, is_preconditioner )
767 {
768 if ( is_preconditioner )
769 throw std::logic_error(
770 "HYPRE BiCGSTAB cannot be used as a preconditioner" );
771
772 auto error = HYPRE_StructBiCGSTABCreate(
773 layout.localGrid()->globalGrid().comm(), &_solver );
774 this->checkHypreError( error );
775 }
776
777 ~HypreStructBiCGSTAB() { HYPRE_StructBiCGSTABDestroy( _solver ); }
778
779 // BiCGSTAB SETTINGS
780
782 void setAbsoluteTol( const double tol )
783 {
784 auto error = HYPRE_StructBiCGSTABSetAbsoluteTol( _solver, tol );
785 this->checkHypreError( error );
786 }
787
789 void setLogging( const int logging )
790 {
791 auto error = HYPRE_StructBiCGSTABSetLogging( _solver, logging );
792 this->checkHypreError( error );
793 }
794
795 HYPRE_StructSolver getHypreSolver() const override { return _solver; }
796 HYPRE_PtrToStructSolverFcn getHypreSetupFunction() const override
797 {
798 return HYPRE_StructBiCGSTABSetup;
799 }
800 HYPRE_PtrToStructSolverFcn getHypreSolveFunction() const override
801 {
802 return HYPRE_StructBiCGSTABSolve;
803 }
804
805 protected:
806 void setToleranceImpl( const double tol ) override
807 {
808 auto error = HYPRE_StructBiCGSTABSetTol( _solver, tol );
809 this->checkHypreError( error );
810 }
811
812 void setMaxIterImpl( const int max_iter ) override
813 {
814 auto error = HYPRE_StructBiCGSTABSetMaxIter( _solver, max_iter );
815 this->checkHypreError( error );
816 }
817
818 void setPrintLevelImpl( const int print_level ) override
819 {
820 auto error = HYPRE_StructBiCGSTABSetPrintLevel( _solver, print_level );
821 this->checkHypreError( error );
822 }
823
824 void setupImpl() override
825 {
826 auto error = HYPRE_StructBiCGSTABSetup( _solver, _A, _b, _x );
827 this->checkHypreError( error );
828 }
829
830 void solveImpl() override
831 {
832 auto error = HYPRE_StructBiCGSTABSolve( _solver, _A, _b, _x );
833 this->checkHypreError( error );
834 }
835
836 int getNumIterImpl() override
837 {
838 HYPRE_Int num_iter;
839 auto error = HYPRE_StructBiCGSTABGetNumIterations( _solver, &num_iter );
840 this->checkHypreError( error );
841 return num_iter;
842 }
843
845 {
846 HYPRE_Real norm;
847 auto error =
848 HYPRE_StructBiCGSTABGetFinalRelativeResidualNorm( _solver, &norm );
849 this->checkHypreError( error );
850 return norm;
851 }
852
855 preconditioner ) override
856 {
857 auto error = HYPRE_StructBiCGSTABSetPrecond(
858 _solver, preconditioner.getHypreSolveFunction(),
859 preconditioner.getHypreSetupFunction(),
860 preconditioner.getHypreSolver() );
861 this->checkHypreError( error );
862 }
863
864 private:
865 HYPRE_StructSolver _solver;
866 using base_type::_A;
867 using base_type::_b;
868 using base_type::_x;
869};
870
871//---------------------------------------------------------------------------//
873template <class Scalar, class EntityType, class MemorySpace>
875 : public HypreStructuredSolver<Scalar, EntityType, MemorySpace>
876{
877 public:
881 template <class ArrayLayout_t>
882 HypreStructPFMG( const ArrayLayout_t& layout,
883 const bool is_preconditioner = false )
884 : base_type( layout, is_preconditioner )
885 {
886 auto error = HYPRE_StructPFMGCreate(
887 layout.localGrid()->globalGrid().comm(), &_solver );
888 this->checkHypreError( error );
889
890 if ( is_preconditioner )
891 {
892 error = HYPRE_StructPFMGSetZeroGuess( _solver );
893 this->checkHypreError( error );
894 }
895 }
896
897 ~HypreStructPFMG() { HYPRE_StructPFMGDestroy( _solver ); }
898
899 // PFMG SETTINGS
900
902 void setMaxLevels( const int max_levels )
903 {
904 auto error = HYPRE_StructPFMGSetMaxLevels( _solver, max_levels );
905 this->checkHypreError( error );
906 }
907
910 void setRelChange( const int rel_change )
911 {
912 auto error = HYPRE_StructPFMGSetRelChange( _solver, rel_change );
913 this->checkHypreError( error );
914 }
915
925 void setRelaxType( const int relax_type )
926 {
927 auto error = HYPRE_StructPFMGSetRelaxType( _solver, relax_type );
928 this->checkHypreError( error );
929 }
930
932 void setJacobiWeight( const double weight )
933 {
934 auto error = HYPRE_StructPFMGSetJacobiWeight( _solver, weight );
935 this->checkHypreError( error );
936 }
937
948 void setRAPType( const int rap_type )
949 {
950 auto error = HYPRE_StructPFMGSetRAPType( _solver, rap_type );
951 this->checkHypreError( error );
952 }
953
955 void setNumPreRelax( const int num_pre_relax )
956 {
957 auto error = HYPRE_StructPFMGSetNumPreRelax( _solver, num_pre_relax );
958 this->checkHypreError( error );
959 }
960
962 void setNumPostRelax( const int num_post_relax )
963 {
964 auto error = HYPRE_StructPFMGSetNumPostRelax( _solver, num_post_relax );
965 this->checkHypreError( error );
966 }
967
971 void setSkipRelax( const int skip_relax )
972 {
973 auto error = HYPRE_StructPFMGSetSkipRelax( _solver, skip_relax );
974 this->checkHypreError( error );
975 }
976
978 void setLogging( const int logging )
979 {
980 auto error = HYPRE_StructPFMGSetLogging( _solver, logging );
981 this->checkHypreError( error );
982 }
983
984 HYPRE_StructSolver getHypreSolver() const override { return _solver; }
985 HYPRE_PtrToStructSolverFcn getHypreSetupFunction() const override
986 {
987 return HYPRE_StructPFMGSetup;
988 }
989 HYPRE_PtrToStructSolverFcn getHypreSolveFunction() const override
990 {
991 return HYPRE_StructPFMGSolve;
992 }
993
994 protected:
995 void setToleranceImpl( const double tol ) override
996 {
997 auto error = HYPRE_StructPFMGSetTol( _solver, tol );
998 this->checkHypreError( error );
999 }
1000
1001 void setMaxIterImpl( const int max_iter ) override
1002 {
1003 auto error = HYPRE_StructPFMGSetMaxIter( _solver, max_iter );
1004 this->checkHypreError( error );
1005 }
1006
1007 void setPrintLevelImpl( const int print_level ) override
1008 {
1009 auto error = HYPRE_StructPFMGSetPrintLevel( _solver, print_level );
1010 this->checkHypreError( error );
1011 }
1012
1013 void setupImpl() override
1014 {
1015 auto error = HYPRE_StructPFMGSetup( _solver, _A, _b, _x );
1016 this->checkHypreError( error );
1017 }
1018
1019 void solveImpl() override
1020 {
1021 auto error = HYPRE_StructPFMGSolve( _solver, _A, _b, _x );
1022 this->checkHypreError( error );
1023 }
1024
1025 int getNumIterImpl() override
1026 {
1027 HYPRE_Int num_iter;
1028 auto error = HYPRE_StructPFMGGetNumIterations( _solver, &num_iter );
1029 this->checkHypreError( error );
1030 return num_iter;
1031 }
1032
1034 {
1035 HYPRE_Real norm;
1036 auto error =
1037 HYPRE_StructPFMGGetFinalRelativeResidualNorm( _solver, &norm );
1038 this->checkHypreError( error );
1039 return norm;
1040 }
1041
1044 {
1045 throw std::logic_error(
1046 "HYPRE PFMG solver does not support preconditioning." );
1047 }
1048
1049 private:
1050 HYPRE_StructSolver _solver;
1051 using base_type::_A;
1052 using base_type::_b;
1053 using base_type::_x;
1054};
1055
1056//---------------------------------------------------------------------------//
1058template <class Scalar, class EntityType, class MemorySpace>
1060 : public HypreStructuredSolver<Scalar, EntityType, MemorySpace>
1061{
1062 public:
1066 template <class ArrayLayout_t>
1067 HypreStructSMG( const ArrayLayout_t& layout,
1068 const bool is_preconditioner = false )
1069 : base_type( layout, is_preconditioner )
1070 {
1071 auto error = HYPRE_StructSMGCreate(
1072 layout.localGrid()->globalGrid().comm(), &_solver );
1073 this->checkHypreError( error );
1074
1075 if ( is_preconditioner )
1076 {
1077 error = HYPRE_StructSMGSetZeroGuess( _solver );
1078 this->checkHypreError( error );
1079 }
1080 }
1081
1082 ~HypreStructSMG() { HYPRE_StructSMGDestroy( _solver ); }
1083
1084 // SMG Settings
1085
1088 void setRelChange( const int rel_change )
1089 {
1090 auto error = HYPRE_StructSMGSetRelChange( _solver, rel_change );
1091 this->checkHypreError( error );
1092 }
1093
1095 void setNumPreRelax( const int num_pre_relax )
1096 {
1097 auto error = HYPRE_StructSMGSetNumPreRelax( _solver, num_pre_relax );
1098 this->checkHypreError( error );
1099 }
1100
1102 void setNumPostRelax( const int num_post_relax )
1103 {
1104 auto error = HYPRE_StructSMGSetNumPostRelax( _solver, num_post_relax );
1105 this->checkHypreError( error );
1106 }
1107
1109 void setLogging( const int logging )
1110 {
1111 auto error = HYPRE_StructSMGSetLogging( _solver, logging );
1112 this->checkHypreError( error );
1113 }
1114
1115 HYPRE_StructSolver getHypreSolver() const override { return _solver; }
1116 HYPRE_PtrToStructSolverFcn getHypreSetupFunction() const override
1117 {
1118 return HYPRE_StructSMGSetup;
1119 }
1120 HYPRE_PtrToStructSolverFcn getHypreSolveFunction() const override
1121 {
1122 return HYPRE_StructSMGSolve;
1123 }
1124
1125 protected:
1126 void setToleranceImpl( const double tol ) override
1127 {
1128 auto error = HYPRE_StructSMGSetTol( _solver, tol );
1129 this->checkHypreError( error );
1130 }
1131
1132 void setMaxIterImpl( const int max_iter ) override
1133 {
1134 auto error = HYPRE_StructSMGSetMaxIter( _solver, max_iter );
1135 this->checkHypreError( error );
1136 }
1137
1138 void setPrintLevelImpl( const int print_level ) override
1139 {
1140 auto error = HYPRE_StructSMGSetPrintLevel( _solver, print_level );
1141 this->checkHypreError( error );
1142 }
1143
1144 void setupImpl() override
1145 {
1146 auto error = HYPRE_StructSMGSetup( _solver, _A, _b, _x );
1147 this->checkHypreError( error );
1148 }
1149
1150 void solveImpl() override
1151 {
1152 auto error = HYPRE_StructSMGSolve( _solver, _A, _b, _x );
1153 this->checkHypreError( error );
1154 }
1155
1156 int getNumIterImpl() override
1157 {
1158 HYPRE_Int num_iter;
1159 auto error = HYPRE_StructSMGGetNumIterations( _solver, &num_iter );
1160 this->checkHypreError( error );
1161 return num_iter;
1162 }
1163
1165 {
1166 HYPRE_Real norm;
1167 auto error =
1168 HYPRE_StructSMGGetFinalRelativeResidualNorm( _solver, &norm );
1169 this->checkHypreError( error );
1170 return norm;
1171 }
1172
1175 {
1176 throw std::logic_error(
1177 "HYPRE SMG solver does not support preconditioning." );
1178 }
1179
1180 private:
1181 HYPRE_StructSolver _solver;
1182 using base_type::_A;
1183 using base_type::_b;
1184 using base_type::_x;
1185};
1186
1187//---------------------------------------------------------------------------//
1189template <class Scalar, class EntityType, class MemorySpace>
1191 : public HypreStructuredSolver<Scalar, EntityType, MemorySpace>
1192{
1193 public:
1197 template <class ArrayLayout_t>
1198 HypreStructJacobi( const ArrayLayout_t& layout,
1199 const bool is_preconditioner = false )
1200 : base_type( layout, is_preconditioner )
1201 {
1202 auto error = HYPRE_StructJacobiCreate(
1203 layout.localGrid()->globalGrid().comm(), &_solver );
1204 this->checkHypreError( error );
1205
1206 if ( is_preconditioner )
1207 {
1208 error = HYPRE_StructJacobiSetZeroGuess( _solver );
1209 this->checkHypreError( error );
1210 }
1211 }
1212
1213 ~HypreStructJacobi() { HYPRE_StructJacobiDestroy( _solver ); }
1214
1215 HYPRE_StructSolver getHypreSolver() const override { return _solver; }
1216 HYPRE_PtrToStructSolverFcn getHypreSetupFunction() const override
1217 {
1218 return HYPRE_StructJacobiSetup;
1219 }
1220 HYPRE_PtrToStructSolverFcn getHypreSolveFunction() const override
1221 {
1222 return HYPRE_StructJacobiSolve;
1223 }
1224
1225 protected:
1226 void setToleranceImpl( const double tol ) override
1227 {
1228 auto error = HYPRE_StructJacobiSetTol( _solver, tol );
1229 this->checkHypreError( error );
1230 }
1231
1232 void setMaxIterImpl( const int max_iter ) override
1233 {
1234 auto error = HYPRE_StructJacobiSetMaxIter( _solver, max_iter );
1235 this->checkHypreError( error );
1236 }
1237
1238 void setPrintLevelImpl( const int ) override
1239 {
1240 // The Jacobi solver does not support a print level.
1241 }
1242
1243 void setupImpl() override
1244 {
1245 auto error = HYPRE_StructJacobiSetup( _solver, _A, _b, _x );
1246 this->checkHypreError( error );
1247 }
1248
1249 void solveImpl() override
1250 {
1251 auto error = HYPRE_StructJacobiSolve( _solver, _A, _b, _x );
1252 this->checkHypreError( error );
1253 }
1254
1255 int getNumIterImpl() override
1256 {
1257 HYPRE_Int num_iter;
1258 auto error = HYPRE_StructJacobiGetNumIterations( _solver, &num_iter );
1259 this->checkHypreError( error );
1260 return num_iter;
1261 }
1262
1264 {
1265 HYPRE_Real norm;
1266 auto error =
1267 HYPRE_StructJacobiGetFinalRelativeResidualNorm( _solver, &norm );
1268 this->checkHypreError( error );
1269 return norm;
1270 }
1271
1274 {
1275 throw std::logic_error(
1276 "HYPRE Jacobi solver does not support preconditioning." );
1277 }
1278
1279 private:
1280 HYPRE_StructSolver _solver;
1281 using base_type::_A;
1282 using base_type::_b;
1283 using base_type::_x;
1284};
1285
1286//---------------------------------------------------------------------------//
1288template <class Scalar, class EntityType, class MemorySpace>
1290 : public HypreStructuredSolver<Scalar, EntityType, MemorySpace>
1291{
1292 public:
1296 template <class ArrayLayout_t>
1297 HypreStructDiagonal( const ArrayLayout_t& layout,
1298 const bool is_preconditioner = false )
1299 : base_type( layout, is_preconditioner )
1300 {
1301 if ( !is_preconditioner )
1302 throw std::logic_error(
1303 "Diagonal preconditioner cannot be used as a solver" );
1304 }
1305
1306 HYPRE_StructSolver getHypreSolver() const override { return nullptr; }
1307 HYPRE_PtrToStructSolverFcn getHypreSetupFunction() const override
1308 {
1309 return HYPRE_StructDiagScaleSetup;
1310 }
1311 HYPRE_PtrToStructSolverFcn getHypreSolveFunction() const override
1312 {
1313 return HYPRE_StructDiagScale;
1314 }
1315
1316 protected:
1317 void setToleranceImpl( const double ) override
1318 {
1319 throw std::logic_error(
1320 "Diagonal preconditioner cannot be used as a solver" );
1321 }
1322
1323 void setMaxIterImpl( const int ) override
1324 {
1325 throw std::logic_error(
1326 "Diagonal preconditioner cannot be used as a solver" );
1327 }
1328
1329 void setPrintLevelImpl( const int ) override
1330 {
1331 throw std::logic_error(
1332 "Diagonal preconditioner cannot be used as a solver" );
1333 }
1334
1335 void setupImpl() override
1336 {
1337 throw std::logic_error(
1338 "Diagonal preconditioner cannot be used as a solver" );
1339 }
1340
1341 void solveImpl() override
1342 {
1343 throw std::logic_error(
1344 "Diagonal preconditioner cannot be used as a solver" );
1345 }
1346
1347 int getNumIterImpl() override
1348 {
1349 throw std::logic_error(
1350 "Diagonal preconditioner cannot be used as a solver" );
1351 }
1352
1354 {
1355 throw std::logic_error(
1356 "Diagonal preconditioner cannot be used as a solver" );
1357 }
1358
1361 {
1362 throw std::logic_error(
1363 "Diagonal preconditioner does not support preconditioning." );
1364 }
1365};
1366
1367//---------------------------------------------------------------------------//
1368// Builders
1369//---------------------------------------------------------------------------//
1372template <class Scalar, class MemorySpace, class ArrayLayout_t>
1373std::shared_ptr<
1374 HypreStructPCG<Scalar, typename ArrayLayout_t::entity_type, MemorySpace>>
1375createHypreStructPCG( const ArrayLayout_t& layout,
1376 const bool is_preconditioner = false )
1377{
1379 "Must use an array layout" );
1380 return std::make_shared<HypreStructPCG<
1381 Scalar, typename ArrayLayout_t::entity_type, MemorySpace>>(
1382 layout, is_preconditioner );
1383}
1384
1387template <class Scalar, class MemorySpace, class ArrayLayout_t>
1388std::shared_ptr<
1389 HypreStructGMRES<Scalar, typename ArrayLayout_t::entity_type, MemorySpace>>
1390createHypreStructGMRES( const ArrayLayout_t& layout,
1391 const bool is_preconditioner = false )
1392{
1394 "Must use an array layout" );
1395 return std::make_shared<HypreStructGMRES<
1396 Scalar, typename ArrayLayout_t::entity_type, MemorySpace>>(
1397 layout, is_preconditioner );
1398}
1399
1402template <class Scalar, class MemorySpace, class ArrayLayout_t>
1403std::shared_ptr<HypreStructBiCGSTAB<Scalar, typename ArrayLayout_t::entity_type,
1404 MemorySpace>>
1405createHypreStructBiCGSTAB( const ArrayLayout_t& layout,
1406 const bool is_preconditioner = false )
1407{
1409 "Must use an array layout" );
1410 return std::make_shared<HypreStructBiCGSTAB<
1411 Scalar, typename ArrayLayout_t::entity_type, MemorySpace>>(
1412 layout, is_preconditioner );
1413}
1414
1417template <class Scalar, class MemorySpace, class ArrayLayout_t>
1418std::shared_ptr<
1419 HypreStructPFMG<Scalar, typename ArrayLayout_t::entity_type, MemorySpace>>
1420createHypreStructPFMG( const ArrayLayout_t& layout,
1421 const bool is_preconditioner = false )
1422{
1424 "Must use an array layout" );
1425 return std::make_shared<HypreStructPFMG<
1426 Scalar, typename ArrayLayout_t::entity_type, MemorySpace>>(
1427 layout, is_preconditioner );
1428}
1429
1432template <class Scalar, class MemorySpace, class ArrayLayout_t>
1433std::shared_ptr<
1434 HypreStructSMG<Scalar, typename ArrayLayout_t::entity_type, MemorySpace>>
1435createHypreStructSMG( const ArrayLayout_t& layout,
1436 const bool is_preconditioner = false )
1437{
1439 "Must use an array layout" );
1440 return std::make_shared<HypreStructSMG<
1441 Scalar, typename ArrayLayout_t::entity_type, MemorySpace>>(
1442 layout, is_preconditioner );
1443}
1444
1447template <class Scalar, class MemorySpace, class ArrayLayout_t>
1448std::shared_ptr<
1449 HypreStructJacobi<Scalar, typename ArrayLayout_t::entity_type, MemorySpace>>
1450createHypreStructJacobi( const ArrayLayout_t& layout,
1451 const bool is_preconditioner = false )
1452{
1454 "Must use an array layout" );
1455 return std::make_shared<HypreStructJacobi<
1456 Scalar, typename ArrayLayout_t::entity_type, MemorySpace>>(
1457 layout, is_preconditioner );
1458}
1459
1462template <class Scalar, class MemorySpace, class ArrayLayout_t>
1463std::shared_ptr<HypreStructDiagonal<Scalar, typename ArrayLayout_t::entity_type,
1464 MemorySpace>>
1465createHypreStructDiagonal( const ArrayLayout_t& layout,
1466 const bool is_preconditioner = false )
1467{
1469 "Must use an array layout" );
1470 return std::make_shared<HypreStructDiagonal<
1471 Scalar, typename ArrayLayout_t::entity_type, MemorySpace>>(
1472 layout, is_preconditioner );
1473}
1474
1475//---------------------------------------------------------------------------//
1476// Factory
1477//---------------------------------------------------------------------------//
1486template <class Scalar, class MemorySpace, class ArrayLayout_t>
1487std::shared_ptr<HypreStructuredSolver<
1488 Scalar, typename ArrayLayout_t::entity_type, MemorySpace>>
1489createHypreStructuredSolver( const std::string& solver_type,
1490 const ArrayLayout_t& layout,
1491 const bool is_preconditioner = false )
1492{
1494 "Must use an array layout" );
1495
1496 if ( "PCG" == solver_type )
1498 is_preconditioner );
1499 else if ( "GMRES" == solver_type )
1501 is_preconditioner );
1502 else if ( "BiCGSTAB" == solver_type )
1504 layout, is_preconditioner );
1505 else if ( "PFMG" == solver_type )
1507 is_preconditioner );
1508 else if ( "SMG" == solver_type )
1510 is_preconditioner );
1511 else if ( "Jacobi" == solver_type )
1513 layout, is_preconditioner );
1514 else if ( "Diagonal" == solver_type )
1516 layout, is_preconditioner );
1517 else
1518 throw std::runtime_error( "Invalid solver type" );
1519}
1520
1521//---------------------------------------------------------------------------//
1522
1523} // namespace Grid
1524} // namespace Cabana
1525
1526#endif // end CABANA_GRID_HYPRESTRUCTUREDSOLVER_HPP
Grid field arrays.
std::shared_ptr< HypreStructGMRES< Scalar, typename ArrayLayout_t::entity_type, MemorySpace > > createHypreStructGMRES(const ArrayLayout_t &layout, const bool is_preconditioner=false)
Definition Cabana_Grid_HypreStructuredSolver.hpp:1390
std::shared_ptr< HypreStructJacobi< Scalar, typename ArrayLayout_t::entity_type, MemorySpace > > createHypreStructJacobi(const ArrayLayout_t &layout, const bool is_preconditioner=false)
Definition Cabana_Grid_HypreStructuredSolver.hpp:1450
std::shared_ptr< HypreStructBiCGSTAB< Scalar, typename ArrayLayout_t::entity_type, MemorySpace > > createHypreStructBiCGSTAB(const ArrayLayout_t &layout, const bool is_preconditioner=false)
Definition Cabana_Grid_HypreStructuredSolver.hpp:1405
std::shared_ptr< HypreStructDiagonal< Scalar, typename ArrayLayout_t::entity_type, MemorySpace > > createHypreStructDiagonal(const ArrayLayout_t &layout, const bool is_preconditioner=false)
Definition Cabana_Grid_HypreStructuredSolver.hpp:1465
std::shared_ptr< HypreStructPCG< Scalar, typename ArrayLayout_t::entity_type, MemorySpace > > createHypreStructPCG(const ArrayLayout_t &layout, const bool is_preconditioner=false)
Definition Cabana_Grid_HypreStructuredSolver.hpp:1375
std::shared_ptr< HypreStructPFMG< Scalar, typename ArrayLayout_t::entity_type, MemorySpace > > createHypreStructPFMG(const ArrayLayout_t &layout, const bool is_preconditioner=false)
Definition Cabana_Grid_HypreStructuredSolver.hpp:1420
std::shared_ptr< HypreStructSMG< Scalar, typename ArrayLayout_t::entity_type, MemorySpace > > createHypreStructSMG(const ArrayLayout_t &layout, const bool is_preconditioner=false)
Definition Cabana_Grid_HypreStructuredSolver.hpp:1435
std::shared_ptr< HypreStructuredSolver< Scalar, typename ArrayLayout_t::entity_type, MemorySpace > > createHypreStructuredSolver(const std::string &solver_type, const ArrayLayout_t &layout, const bool is_preconditioner=false)
Create a HYPRE structured solver.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1489
HYPRE memory space handling.
Logical grid indexing.
KOKKOS_INLINE_FUNCTION auto createSubview(const ViewType &view, const IndexSpace< 1 > &index_space) -> decltype(Kokkos::subview(view, index_space.range(0)))
Given a view create a subview over the given index space.
Definition Cabana_Grid_IndexSpace.hpp:369
Kokkos::View< Scalar *, Params... > createView(const std::string &label, const IndexSpace< 1 > &index_space)
Given an index space create a view over the extent of that index space.
Definition Cabana_Grid_IndexSpace.hpp:235
Grid type tags.
BiCGSTAB solver.
Definition Cabana_Grid_HypreStructuredSolver.hpp:758
double getFinalRelativeResidualNormImpl() override
Get the relative residual norm achieved on the last solve.
Definition Cabana_Grid_HypreStructuredSolver.hpp:844
void solveImpl() override
Solver implementation.
Definition Cabana_Grid_HypreStructuredSolver.hpp:830
void setLogging(const int logging)
Set the amount of logging to do.
Definition Cabana_Grid_HypreStructuredSolver.hpp:789
void setAbsoluteTol(const double tol)
Set the absolute tolerance.
Definition Cabana_Grid_HypreStructuredSolver.hpp:782
HYPRE_StructSolver getHypreSolver() const override
Get the preconditioner.
Definition Cabana_Grid_HypreStructuredSolver.hpp:795
int getNumIterImpl() override
Get the number of iterations taken on the last solve.
Definition Cabana_Grid_HypreStructuredSolver.hpp:836
HYPRE_PtrToStructSolverFcn getHypreSetupFunction() const override
Get the preconditioner setup function.
Definition Cabana_Grid_HypreStructuredSolver.hpp:796
HypreStructuredSolver< Scalar, EntityType, MemorySpace > base_type
Base HYPRE structured solver type.
Definition Cabana_Grid_HypreStructuredSolver.hpp:761
void setMaxIterImpl(const int max_iter) override
Set maximum iteration implementation.
Definition Cabana_Grid_HypreStructuredSolver.hpp:812
void setPreconditionerImpl(const HypreStructuredSolver< Scalar, EntityType, MemorySpace > &preconditioner) override
Set a preconditioner.
Definition Cabana_Grid_HypreStructuredSolver.hpp:853
void setPrintLevelImpl(const int print_level) override
Set the output level.
Definition Cabana_Grid_HypreStructuredSolver.hpp:818
HYPRE_PtrToStructSolverFcn getHypreSolveFunction() const override
Get the preconditioner solve function.
Definition Cabana_Grid_HypreStructuredSolver.hpp:800
void setupImpl() override
Setup implementation.
Definition Cabana_Grid_HypreStructuredSolver.hpp:824
void setToleranceImpl(const double tol) override
Set convergence tolerance implementation.
Definition Cabana_Grid_HypreStructuredSolver.hpp:806
HypreStructBiCGSTAB(const ArrayLayout_t &layout, const bool is_preconditioner=false)
Constructor.
Definition Cabana_Grid_HypreStructuredSolver.hpp:764
Diagonal preconditioner.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1291
void solveImpl() override
Solver implementation.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1341
double getFinalRelativeResidualNormImpl() override
Get the relative residual norm achieved on the last solve.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1353
int getNumIterImpl() override
Get the number of iterations taken on the last solve.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1347
void setupImpl() override
Setup implementation.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1335
void setPrintLevelImpl(const int) override
Set the output level.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1329
HYPRE_PtrToStructSolverFcn getHypreSetupFunction() const override
Get the preconditioner setup function.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1307
HYPRE_StructSolver getHypreSolver() const override
Get the preconditioner.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1306
HypreStructuredSolver< Scalar, EntityType, MemorySpace > base_type
Base HYPRE structured solver type.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1294
void setToleranceImpl(const double) override
Set convergence tolerance implementation.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1317
HypreStructDiagonal(const ArrayLayout_t &layout, const bool is_preconditioner=false)
Constructor.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1297
void setPreconditionerImpl(const HypreStructuredSolver< Scalar, EntityType, MemorySpace > &) override
Set a preconditioner.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1359
void setMaxIterImpl(const int) override
Set maximum iteration implementation.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1323
HYPRE_PtrToStructSolverFcn getHypreSolveFunction() const override
Get the preconditioner solve function.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1311
GMRES solver.
Definition Cabana_Grid_HypreStructuredSolver.hpp:633
void setToleranceImpl(const double tol) override
Set convergence tolerance implementation.
Definition Cabana_Grid_HypreStructuredSolver.hpp:688
void solveImpl() override
Solver implementation.
Definition Cabana_Grid_HypreStructuredSolver.hpp:712
HypreStructuredSolver< Scalar, EntityType, MemorySpace > base_type
Base HYPRE structured solver type.
Definition Cabana_Grid_HypreStructuredSolver.hpp:636
HYPRE_PtrToStructSolverFcn getHypreSetupFunction() const override
Get the preconditioner setup function.
Definition Cabana_Grid_HypreStructuredSolver.hpp:678
HypreStructGMRES(const ArrayLayout_t &layout, const bool is_preconditioner=false)
Constructor.
Definition Cabana_Grid_HypreStructuredSolver.hpp:639
int getNumIterImpl() override
Get the number of iterations taken on the last solve.
Definition Cabana_Grid_HypreStructuredSolver.hpp:718
void setPreconditionerImpl(const HypreStructuredSolver< Scalar, EntityType, MemorySpace > &preconditioner) override
Set a preconditioner.
Definition Cabana_Grid_HypreStructuredSolver.hpp:735
void setMaxIterImpl(const int max_iter) override
Set maximum iteration implementation.
Definition Cabana_Grid_HypreStructuredSolver.hpp:694
HYPRE_PtrToStructSolverFcn getHypreSolveFunction() const override
Get the preconditioner solve function.
Definition Cabana_Grid_HypreStructuredSolver.hpp:682
HYPRE_StructSolver getHypreSolver() const override
Get the preconditioner.
Definition Cabana_Grid_HypreStructuredSolver.hpp:677
void setKDim(const int k_dim)
Set the max size of the Krylov space.
Definition Cabana_Grid_HypreStructuredSolver.hpp:664
double getFinalRelativeResidualNormImpl() override
Get the relative residual norm achieved on the last solve.
Definition Cabana_Grid_HypreStructuredSolver.hpp:726
void setLogging(const int logging)
Set the amount of logging to do.
Definition Cabana_Grid_HypreStructuredSolver.hpp:671
void setAbsoluteTol(const double tol)
Set the absolute tolerance.
Definition Cabana_Grid_HypreStructuredSolver.hpp:657
void setupImpl() override
Setup implementation.
Definition Cabana_Grid_HypreStructuredSolver.hpp:706
void setPrintLevelImpl(const int print_level) override
Set the output level.
Definition Cabana_Grid_HypreStructuredSolver.hpp:700
Jacobi solver.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1192
void setPreconditionerImpl(const HypreStructuredSolver< Scalar, EntityType, MemorySpace > &) override
Set a preconditioner.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1272
void setToleranceImpl(const double tol) override
Set convergence tolerance implementation.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1226
double getFinalRelativeResidualNormImpl() override
Get the relative residual norm achieved on the last solve.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1263
HypreStructuredSolver< Scalar, EntityType, MemorySpace > base_type
Base HYPRE structured solver type.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1195
HYPRE_PtrToStructSolverFcn getHypreSolveFunction() const override
Get the preconditioner solve function.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1220
void setupImpl() override
Setup implementation.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1243
void setMaxIterImpl(const int max_iter) override
Set maximum iteration implementation.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1232
int getNumIterImpl() override
Get the number of iterations taken on the last solve.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1255
void setPrintLevelImpl(const int) override
Set the output level.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1238
HYPRE_PtrToStructSolverFcn getHypreSetupFunction() const override
Get the preconditioner setup function.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1216
void solveImpl() override
Solver implementation.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1249
HYPRE_StructSolver getHypreSolver() const override
Get the preconditioner.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1215
HypreStructJacobi(const ArrayLayout_t &layout, const bool is_preconditioner=false)
Constructor.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1198
PCG solver.
Definition Cabana_Grid_HypreStructuredSolver.hpp:505
void setRelChange(const int rel_change)
Definition Cabana_Grid_HypreStructuredSolver.hpp:539
void setToleranceImpl(const double tol) override
Set convergence tolerance implementation.
Definition Cabana_Grid_HypreStructuredSolver.hpp:563
void setupImpl() override
Setup implementation.
Definition Cabana_Grid_HypreStructuredSolver.hpp:581
HYPRE_StructSolver getHypreSolver() const override
Get the preconditioner.
Definition Cabana_Grid_HypreStructuredSolver.hpp:552
HypreStructuredSolver< Scalar, EntityType, MemorySpace > base_type
Base HYPRE structured solver type.
Definition Cabana_Grid_HypreStructuredSolver.hpp:508
HypreStructPCG(const ArrayLayout_t &layout, const bool is_preconditioner=false)
Constructor.
Definition Cabana_Grid_HypreStructuredSolver.hpp:511
HYPRE_PtrToStructSolverFcn getHypreSetupFunction() const override
Get the preconditioner setup function.
Definition Cabana_Grid_HypreStructuredSolver.hpp:553
void solveImpl() override
Solver implementation.
Definition Cabana_Grid_HypreStructuredSolver.hpp:587
int getNumIterImpl() override
Get the number of iterations taken on the last solve.
Definition Cabana_Grid_HypreStructuredSolver.hpp:593
void setAbsoluteTol(const double tol)
Set the absolute tolerance.
Definition Cabana_Grid_HypreStructuredSolver.hpp:531
double getFinalRelativeResidualNormImpl() override
Get the relative residual norm achieved on the last solve.
Definition Cabana_Grid_HypreStructuredSolver.hpp:601
void setMaxIterImpl(const int max_iter) override
Set maximum iteration implementation.
Definition Cabana_Grid_HypreStructuredSolver.hpp:569
HYPRE_PtrToStructSolverFcn getHypreSolveFunction() const override
Get the preconditioner solve function.
Definition Cabana_Grid_HypreStructuredSolver.hpp:557
void setPrintLevelImpl(const int print_level) override
Set the output level.
Definition Cabana_Grid_HypreStructuredSolver.hpp:575
void setLogging(const int logging)
Set the amount of logging to do.
Definition Cabana_Grid_HypreStructuredSolver.hpp:546
void setPreconditionerImpl(const HypreStructuredSolver< Scalar, EntityType, MemorySpace > &preconditioner) override
Set a preconditioner.
Definition Cabana_Grid_HypreStructuredSolver.hpp:610
PFMG solver.
Definition Cabana_Grid_HypreStructuredSolver.hpp:876
HypreStructuredSolver< Scalar, EntityType, MemorySpace > base_type
Base HYPRE structured solver type.
Definition Cabana_Grid_HypreStructuredSolver.hpp:879
HYPRE_PtrToStructSolverFcn getHypreSolveFunction() const override
Get the preconditioner solve function.
Definition Cabana_Grid_HypreStructuredSolver.hpp:989
void setRelaxType(const int relax_type)
Set relaxation type.
Definition Cabana_Grid_HypreStructuredSolver.hpp:925
void setPreconditionerImpl(const HypreStructuredSolver< Scalar, EntityType, MemorySpace > &) override
Set a preconditioner.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1042
void setSkipRelax(const int skip_relax)
Definition Cabana_Grid_HypreStructuredSolver.hpp:971
void setPrintLevelImpl(const int print_level) override
Set the output level.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1007
double getFinalRelativeResidualNormImpl() override
Get the relative residual norm achieved on the last solve.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1033
void setLogging(const int logging)
Set the amount of logging to do.
Definition Cabana_Grid_HypreStructuredSolver.hpp:978
void setNumPreRelax(const int num_pre_relax)
Set number of relaxation sweeps before coarse-grid correction.
Definition Cabana_Grid_HypreStructuredSolver.hpp:955
void setJacobiWeight(const double weight)
Set the Jacobi weight.
Definition Cabana_Grid_HypreStructuredSolver.hpp:932
HYPRE_StructSolver getHypreSolver() const override
Get the preconditioner.
Definition Cabana_Grid_HypreStructuredSolver.hpp:984
void solveImpl() override
Solver implementation.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1019
HypreStructPFMG(const ArrayLayout_t &layout, const bool is_preconditioner=false)
Constructor.
Definition Cabana_Grid_HypreStructuredSolver.hpp:882
HYPRE_PtrToStructSolverFcn getHypreSetupFunction() const override
Get the preconditioner setup function.
Definition Cabana_Grid_HypreStructuredSolver.hpp:985
void setupImpl() override
Setup implementation.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1013
int getNumIterImpl() override
Get the number of iterations taken on the last solve.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1025
void setMaxIterImpl(const int max_iter) override
Set maximum iteration implementation.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1001
void setToleranceImpl(const double tol) override
Set convergence tolerance implementation.
Definition Cabana_Grid_HypreStructuredSolver.hpp:995
void setRelChange(const int rel_change)
Definition Cabana_Grid_HypreStructuredSolver.hpp:910
void setNumPostRelax(const int num_post_relax)
Set number of relaxation sweeps before coarse-grid correction.
Definition Cabana_Grid_HypreStructuredSolver.hpp:962
void setMaxLevels(const int max_levels)
Set the maximum number of multigrid levels.
Definition Cabana_Grid_HypreStructuredSolver.hpp:902
void setRAPType(const int rap_type)
Set type of coarse-grid operator to use.
Definition Cabana_Grid_HypreStructuredSolver.hpp:948
SMG solver.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1061
HypreStructSMG(const ArrayLayout_t &layout, const bool is_preconditioner=false)
Constructor.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1067
void solveImpl() override
Solver implementation.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1150
void setPreconditionerImpl(const HypreStructuredSolver< Scalar, EntityType, MemorySpace > &) override
Set a preconditioner.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1173
void setToleranceImpl(const double tol) override
Set convergence tolerance implementation.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1126
void setRelChange(const int rel_change)
Definition Cabana_Grid_HypreStructuredSolver.hpp:1088
HypreStructuredSolver< Scalar, EntityType, MemorySpace > base_type
Base HYPRE structured solver type.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1064
HYPRE_PtrToStructSolverFcn getHypreSetupFunction() const override
Get the preconditioner setup function.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1116
HYPRE_PtrToStructSolverFcn getHypreSolveFunction() const override
Get the preconditioner solve function.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1120
void setLogging(const int logging)
Set the amount of logging to do.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1109
void setMaxIterImpl(const int max_iter) override
Set maximum iteration implementation.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1132
int getNumIterImpl() override
Get the number of iterations taken on the last solve.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1156
void setPrintLevelImpl(const int print_level) override
Set the output level.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1138
void setupImpl() override
Setup implementation.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1144
double getFinalRelativeResidualNormImpl() override
Get the relative residual norm achieved on the last solve.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1164
void setNumPostRelax(const int num_post_relax)
Set number of relaxation sweeps before coarse-grid correction.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1102
void setNumPreRelax(const int num_pre_relax)
Set number of relaxation sweeps before coarse-grid correction.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1095
HYPRE_StructSolver getHypreSolver() const override
Get the preconditioner.
Definition Cabana_Grid_HypreStructuredSolver.hpp:1115
Hypre structured solver interface for scalar fields.
Definition Cabana_Grid_HypreStructuredSolver.hpp:49
bool isPreconditioner() const
Return if this solver is a preconditioner.
Definition Cabana_Grid_HypreStructuredSolver.hpp:178
Scalar value_type
Scalar value type.
Definition Cabana_Grid_HypreStructuredSolver.hpp:56
virtual HYPRE_PtrToStructSolverFcn getHypreSetupFunction() const =0
Get the preconditioner setup function.
void printMatrix(const char *prefix)
Print the hypre matrix to output file.
Definition Cabana_Grid_HypreStructuredSolver.hpp:288
HYPRE_StructMatrix _A
Matrix for the problem Ax = b.
Definition Cabana_Grid_HypreStructuredSolver.hpp:482
HYPRE_StructVector _b
Forcing term for the problem Ax = b.
Definition Cabana_Grid_HypreStructuredSolver.hpp:484
double getFinalRelativeResidualNorm()
Get the relative residual norm achieved on the last solve.
Definition Cabana_Grid_HypreStructuredSolver.hpp:426
void setup()
Setup the problem.
Definition Cabana_Grid_HypreStructuredSolver.hpp:342
virtual void setMaxIterImpl(const int max_iter)=0
Set maximum iteration implementation.
HypreStructuredSolver(const ArrayLayout_t &layout, const bool is_preconditioner=false)
Hypre memory space compatibility check.
Definition Cabana_Grid_HypreStructuredSolver.hpp:68
virtual void setPreconditionerImpl(const HypreStructuredSolver< Scalar, EntityType, MemorySpace > &preconditioner)=0
Set a preconditioner.
virtual void setToleranceImpl(const double tol)=0
Set convergence tolerance implementation.
void setMatrixStencil(const std::vector< std::array< int, NumSpaceDim > > &stencil, const bool is_symmetric=false)
Set the operator stencil.
Definition Cabana_Grid_HypreStructuredSolver.hpp:190
EntityType entity_type
Entity type.
Definition Cabana_Grid_HypreStructuredSolver.hpp:52
virtual double getFinalRelativeResidualNormImpl()=0
Get the relative residual norm achieved on the last solve.
void setPrintLevel(const int print_level)
Set the output level.
Definition Cabana_Grid_HypreStructuredSolver.hpp:318
virtual void solveImpl()=0
Solver implementation.
virtual int getNumIterImpl()=0
Get the number of iterations taken on the last solve.
MemorySpace memory_space
Kokkos memory space..
Definition Cabana_Grid_HypreStructuredSolver.hpp:54
virtual HYPRE_StructSolver getHypreSolver() const =0
Get the preconditioner.
void setMatrixValues(const Array_t &values)
Set the matrix values.
Definition Cabana_Grid_HypreStructuredSolver.hpp:230
void printLHS(const char *prefix)
Print the hypre LHS to output file.
Definition Cabana_Grid_HypreStructuredSolver.hpp:297
void printRHS(const char *prefix)
Print the hypre RHS to output file.
Definition Cabana_Grid_HypreStructuredSolver.hpp:306
void checkHypreError(const int error) const
Check a hypre error.
Definition Cabana_Grid_HypreStructuredSolver.hpp:466
int getNumIter()
Get the number of iterations taken on the last solve.
Definition Cabana_Grid_HypreStructuredSolver.hpp:423
virtual HYPRE_PtrToStructSolverFcn getHypreSolveFunction() const =0
Get the preconditioner solve function.
virtual void setupImpl()=0
Setup implementation.
HYPRE_StructVector _x
Solution to the problem Ax = b.
Definition Cabana_Grid_HypreStructuredSolver.hpp:486
void setMaxIter(const int max_iter)
Set maximum iteration implementation.
Definition Cabana_Grid_HypreStructuredSolver.hpp:315
void solve(const Array_t &b, Array_t &x)
Solve the problem Ax = b for x.
Definition Cabana_Grid_HypreStructuredSolver.hpp:358
void setTolerance(const double tol)
Set convergence tolerance implementation.
Definition Cabana_Grid_HypreStructuredSolver.hpp:312
void setPreconditioner(const std::shared_ptr< HypreStructuredSolver< Scalar, EntityType, MemorySpace > > &preconditioner)
Set a preconditioner.
Definition Cabana_Grid_HypreStructuredSolver.hpp:325
virtual void setPrintLevelImpl(const int print_level)=0
Set the output level.
Structured index space.
Definition Cabana_Grid_IndexSpace.hpp:37
Core: particle data structures and algorithms.
Definition Cabana_AoSoA.hpp:36
Global index tag.
Definition Cabana_Grid_Types.hpp:215
Hypre device compatibility check.
Definition Cabana_Grid_Hypre.hpp:39
Local index tag.
Definition Cabana_Grid_Types.hpp:208
Owned decomposition tag.
Definition Cabana_Grid_Types.hpp:190
Array static type checker.
Definition Cabana_Grid_Array.hpp:160
Definition Cabana_Grid_Array.hpp:322