5.2.5.1. Teko Solver for OMD with Porous Flows Example
This example shows how to set up a Teko solver for an OMD example problem with porous flows.
5.2.5.1.1. Problem Files
The files required for this example can be downloaded here, or found at $TESTS/aria_rtest/example/training/omd_fic_enclosure (where TESTS=/projects/sierra/tests in a CEE environment, or your Sierra test distribution otherwise).
This example problem simulates a two dimensional foam-in-a-can experimental setup with enclosure radiation. There are three materials in this example: a fluid-only region (block 6), foam - a porous medium containing both fluid and solid phases (block 5), and a solid material for the can (blocks 1-4). For this particular case, a single Aria region encompassing the entire domain is used. A single, monolithic equation system is used for this simulation. Our main concern in this example, therefore, is demonstrating how to construct a Teko block preconditioner using the advice from the Road-map for Constructing Teko Solvers for the monolithic, multiphysics system.
5.2.5.1.2. Setting up a Teko Solver
Here we follow steps 1-4 as outlined in the Road-map for Constructing Teko Solvers.
We start by gather the DOFs in the monolithic linear system, which in alphabetical order are:
ENTHALPY GAS_PHASE_MASS_FRACTION_CO2 GAS_PHASE_MASS_FRACTION_Organic GAS_PHASE_PRESSURE GAS_PHASE_TEMPERATURE MASS_FRACTION_CO2 MASS_FRACTION_Organic PRESSURE SOLID_PHASE_DENSITY SOLID_PHASE_MASS_FRACTION_FoamA SOLID_PHASE_MASS_FRACTION_FoamB SOLID_PHASE_TEMPERATURE TEMPERATURE VELOCITY
With the DOFs identified above, we start by assigning each DOF to its own block. We employ a Block Gauss-Seidel approach with a direct solver for each diagonal sub-block, as recommended in step 2. The sample input deck block is shown below:
BEGIN TPETRA EQUATION SOLVER GMRES_TEKO BEGIN GMRES SOLVER CONVERGENCE TOLERANCE = 1e-8 BEGIN EXPERT TEKO PRECONDITIONER INVERSE METHOD = OMD_SOLVER XML FILE = teko.xml DEFINE MATRIX SUBBLOCK NUMBER=1 DOF=ENTHALPY DEFINE MATRIX SUBBLOCK NUMBER=2 DOF=GAS_PHASE_MASS_FRACTION_CO2 DEFINE MATRIX SUBBLOCK NUMBER=3 DOF=GAS_PHASE_MASS_FRACTION_Organic DEFINE MATRIX SUBBLOCK NUMBER=4 DOF=GAS_PHASE_PRESSURE DEFINE MATRIX SUBBLOCK NUMBER=5 DOF=GAS_PHASE_TEMPERATURE DEFINE MATRIX SUBBLOCK NUMBER=6 DOF=MASS_FRACTION_CO2 DEFINE MATRIX SUBBLOCK NUMBER=7 DOF=MASS_FRACTION_Organic DEFINE MATRIX SUBBLOCK NUMBER=8 DOF=PRESSURE DEFINE MATRIX SUBBLOCK NUMBER=9 DOF=SOLID_PHASE_DENSITY DEFINE MATRIX SUBBLOCK NUMBER=10 DOF=SOLID_PHASE_MASS_FRACTION_FoamA DEFINE MATRIX SUBBLOCK NUMBER=11 DOF=SOLID_PHASE_MASS_FRACTION_FoamB DEFINE MATRIX SUBBLOCK NUMBER=12 DOF=SOLID_PHASE_TEMPERATURE DEFINE MATRIX SUBBLOCK NUMBER=13 DOF=TEMPERATURE DEFINE MATRIX SUBBLOCK NUMBER=14 DOF=VELOCITY END END END TPETRA EQUATION SOLVERThe minimal input for
teko.xmlis:<ParameterList> <ParameterList name="OMD_SOLVER"> <Parameter name="Type" type="string" value="Block Gauss-Seidel"/> <Parameter name="Use Upper Triangle" type="bool" value="true"/> <Parameter name="Inverse Type" type="string" value="Amesos2"/> </ParameterList> </ParameterList>
Now we need to identify the DOF-to-block mapping. We start with the one-to-one mapping from step 2 and run a single time-step and observe the solve:
Equation System all->main: * Step : Transient, Strategy: LINE_SEARCH, Time: 1.00e-04, Step: 1.00e-04 * Matrix: Solver: "gmres_teko", Unknowns: 12270, Nonzeros: 688594 * Mesh : Processor 0 of 1: 2396 of 2396 elems, 2514 of 2514 nodes * Computing View Factors for enclosure encl_Enc1 N O N L I N E A R L I N E A R T I M E ------------------------------- -------------------- ----------------------- Step SubItns Resid Delta Itns Status Resid Asm/Slv/SubItn Time ---- ------- -------- -------- ---- ------ -------- ----------------------- 1 1 2.16e+02 1.63e+01 170 ok 9.72e-09 7.2e-02/1.7e+00/2.5e-02 2 1 7.13e-03 6.00e+00 143 ok 8.06e-09 7.4e-02/1.6e+00/2.4e-02 3 1 2.10e-05 2.71e-04 206 ok 8.19e-09 7.2e-02/1.8e+00/2.4e-02 4 0 2.51e-06 9.47e-07 222 ok 9.71e-09 7.2e-02/1.9e+00/3.0e-04 Termination reason: 9.46792e-07 < nonlinear_correction_tolerance(1e-06)Considering that we are employing a direct solver for each of the sub-blocks, these iterations counts are too high. Keeping usage of our direct solver, we can start grouping different DOFs to the same sub-block. Through some trial-and-error and physical intuition, we realize that:
The enthalpy and pressure equations are strongly coupled
The solid phase and gas phase temperatures, along with the gas phase pressure, are strongly coupled
All other blocks, in comparison, are only ‘weakly’ coupled
This is incorporated into the following input deck:
BEGIN TPETRA EQUATION SOLVER GMRES_TEKO BEGIN GMRES SOLVER CONVERGENCE TOLERANCE = 1e-8 # in this case, for fast nonlinear convergence BEGIN EXPERT TEKO PRECONDITIONER INVERSE METHOD = OMD_SOLVER XML FILE = teko.xml DEFINE MATRIX SUBBLOCK NUMBER=1 DOF=PRESSURE DEFINE MATRIX SUBBLOCK NUMBER=1 DOF=ENTHALPY DEFINE MATRIX SUBBLOCK NUMBER=2 DOF=SOLID_PHASE_TEMPERATURE DEFINE MATRIX SUBBLOCK NUMBER=2 DOF=GAS_PHASE_TEMPERATURE DEFINE MATRIX SUBBLOCK NUMBER=2 DOF=GAS_PHASE_PRESSURE DEFINE MATRIX SUBBLOCK NUMBER=3 DOF=VELOCITY DEFINE MATRIX SUBBLOCK NUMBER=4 DOF=TEMPERATURE DEFINE MATRIX SUBBLOCK NUMBER=5 DOF=GAS_PHASE_MASS_FRACTION_CO2 DEFINE MATRIX SUBBLOCK NUMBER=6 DOF=GAS_PHASE_MASS_FRACTION_Organic DEFINE MATRIX SUBBLOCK NUMBER=7 DOF=MASS_FRACTION_CO2 DEFINE MATRIX SUBBLOCK NUMBER=8 DOF=MASS_FRACTION_Organic DEFINE MATRIX SUBBLOCK NUMBER=9 DOF=SOLID_PHASE_DENSITY DEFINE MATRIX SUBBLOCK NUMBER=10 DOF=SOLID_PHASE_MASS_FRACTION_FoamA DEFINE MATRIX SUBBLOCK NUMBER=11 DOF=SOLID_PHASE_MASS_FRACTION_FoamB END END END TPETRA EQUATION SOLVERUsing the input deck from above, we see that the iteration counts are significantly improved:
Equation System all->main: * Step : Transient, Strategy: LINE_SEARCH, Time: 1.00e-04, Step: 1.00e-04 * Matrix: Solver: "gmres_teko", Unknowns: 12270, Nonzeros: 688594 * Mesh : Processor 0 of 1: 2396 of 2396 elems, 2514 of 2514 nodes * Computing View Factors for enclosure encl_Enc1 N O N L I N E A R L I N E A R T I M E ------------------------------- -------------------- ----------------------- Step SubItns Resid Delta Itns Status Resid Asm/Slv/SubItn Time ---- ------- -------- -------- ---- ------ -------- ----------------------- 1 1 2.16e+02 1.48e+01 1 ok 4.37e-09 6.9e-02/1.4e+00/2.3e-02 2 1 7.16e-03 3.38e-04 2 ok 2.59e-10 6.8e-02/7.8e-01/2.3e-02 3 0 5.38e-05 4.78e-08 2 ok 1.83e-10 6.8e-02/8.0e-01/2.8e-04 Termination reason: 4.77609e-08 < nonlinear_correction_tolerance(1e-06)
Now that we can reach convergence in a few iterations with Block Gauss-Seidel, we must now try to remove the direct solver usage for each sub-block. Converting each sub-block one-at-a-time, we observe that the solid phase density and all mass fraction terms can be readily solved via Jacobi. That is, the following XML yields a solver with the same iteration counts as the direct solver for those sub-blocks
<ParameterList> <ParameterList name="OMD_SOLVER"> <Parameter name="Type" type="string" value="Block Gauss-Seidel"/> <Parameter name="Use Upper Triangle" type="bool" value="true"/> <Parameter name="Inverse Type" type="string" value="Amesos2"/> <!-- Solid phase mass fractions, solid phase density --> <Parameter name="Inverse Type 5" type="string" value="Jacobi"/> <Parameter name="Inverse Type 6" type="string" value="Jacobi"/> <Parameter name="Inverse Type 7" type="string" value="Jacobi"/> <Parameter name="Inverse Type 8" type="string" value="Jacobi"/> <Parameter name="Inverse Type 9" type="string" value="Jacobi"/> <Parameter name="Inverse Type 10" type="string" value="Jacobi"/> <Parameter name="Inverse Type 11" type="string" value="Jacobi"/> </ParameterList> <ParameterList name="Jacobi"> <Parameter name="Type" type="string" value="Ifpack2"/> <Parameter name="Prec Type" type="string" value="RELAXATION"/> <ParameterList name="Ifpack2 Settings"> <Parameter name="relaxation: type" type="string" value="Jacobi"/> </ParameterList> </ParameterList> </ParameterList>With the following results from the log file:
Equation System all->main: * Step : Transient, Strategy: LINE_SEARCH, Time: 1.00e-04, Step: 1.00e-04 * Matrix: Solver: "gmres_teko", Unknowns: 12270, Nonzeros: 688594 * Mesh : Processor 0 of 1: 2396 of 2396 elems, 2514 of 2514 nodes * Computing View Factors for enclosure encl_Enc1 N O N L I N E A R L I N E A R T I M E ------------------------------- -------------------- ----------------------- Step SubItns Resid Delta Itns Status Resid Asm/Slv/SubItn Time ---- ------- -------- -------- ---- ------ -------- ----------------------- 1 1 2.16e+02 1.48e+01 1 ok 4.37e-09 6.8e-02/7.5e-01/2.3e-02 2 1 7.16e-03 3.38e-04 2 ok 2.59e-10 6.6e-02/7.4e-01/2.2e-02 3 0 5.38e-05 4.78e-08 2 ok 1.83e-10 6.6e-02/7.4e-01/2.6e-04 Termination reason: 4.77609e-08 < nonlinear_correction_tolerance(1e-06)When this occurs, you can generally combine all of the blocks using Jacobi into a single, large block solves via Jacobi. This typically leads to the same iteration count and a faster solver time-to-solution. In addition, we assume that DD-ILU can be easily applied to the velocity and temperature equations. The same cannot be said, however, for the combined pressure-enthalpy block and the combined solid/gas phase temperature-gas phase pressure block. Using a diagnostic preconditioner as recommended in the Road-map for Constructing Teko Solvers, we observe that those two sub-blocks require either a better preconditioner or a sub-iteration scheme. Here, we choose the latter option. Since we plan on using an iterative solver inside our preconditioner, we must switch to a flexible GMRES variant provided through the
EXPERT BELOS SOLVERblock. Combining all of these elements together, we construct our Teko preconditioner. The solver from the input file is shown below:BEGIN TPETRA EQUATION SOLVER GMRES_TEKO BEGIN EXPERT BELOS SOLVER XML FILE = belos.xml BEGIN EXPERT TEKO PRECONDITIONER INVERSE METHOD = OMD_SOLVER XML FILE = teko.xml DEFINE MATRIX SUBBLOCK NUMBER=1 DOF=PRESSURE DEFINE MATRIX SUBBLOCK NUMBER=1 DOF=ENTHALPY DEFINE MATRIX SUBBLOCK NUMBER=2 DOF=SOLID_PHASE_TEMPERATURE DEFINE MATRIX SUBBLOCK NUMBER=2 DOF=GAS_PHASE_TEMPERATURE DEFINE MATRIX SUBBLOCK NUMBER=2 DOF=GAS_PHASE_PRESSURE DEFINE MATRIX SUBBLOCK NUMBER=3 DOF=VELOCITY DEFINE MATRIX SUBBLOCK NUMBER=4 DOF=TEMPERATURE DEFINE MATRIX SUBBLOCK NUMBER=5 DOF=GAS_PHASE_MASS_FRACTION_CO2 DEFINE MATRIX SUBBLOCK NUMBER=5 DOF=GAS_PHASE_MASS_FRACTION_Organic DEFINE MATRIX SUBBLOCK NUMBER=5 DOF=MASS_FRACTION_CO2 DEFINE MATRIX SUBBLOCK NUMBER=5 DOF=MASS_FRACTION_Organic DEFINE MATRIX SUBBLOCK NUMBER=5 DOF=SOLID_PHASE_DENSITY DEFINE MATRIX SUBBLOCK NUMBER=5 DOF=SOLID_PHASE_MASS_FRACTION_FoamA DEFINE MATRIX SUBBLOCK NUMBER=5 DOF=SOLID_PHASE_MASS_FRACTION_FoamB END END END TPETRA EQUATION SOLVERThe contents of
teko.xmlare as follows:<ParameterList> <ParameterList name="OMD_SOLVER"> <Parameter name="Type" type="string" value="Block Gauss-Seidel"/> <Parameter name="Use Upper Triangle" type="bool" value="true"/> <!-- Combined energy/pressure block --> <Parameter name="Inverse Type 1" type="string" value="gmres"/> <Parameter name="Preconditioner Type 1" type="string" value="ddilu"/> <!-- Combined gas/solid temperature, gas pressure --> <Parameter name="Inverse Type 2" type="string" value="gmres"/> <Parameter name="Preconditioner Type 2" type="string" value="ddilu"/> <!-- Velocity --> <Parameter name="Inverse Type 3" type="string" value="ddilu"/> <!-- Temperature --> <Parameter name="Inverse Type 4" type="string" value="ddilu"/> <!-- Solid density, mass fractions --> <Parameter name="Inverse Type 5" type="string" value="Jacobi"/> </ParameterList> <ParameterList name="Jacobi"> <Parameter name="Type" type="string" value="Ifpack2"/> <Parameter name="Prec Type" type="string" value="RELAXATION"/> <ParameterList name="Ifpack2 Settings"> <Parameter name="relaxation: type" type="string" value="Jacobi"/> </ParameterList> </ParameterList> <ParameterList name="ddilu"> <Parameter name="Type" type="string" value="Ifpack2"/> <Parameter name="Prec Type" type="string" value="SCHWARZ"/> <ParameterList name="Ifpack2 Settings"> <Parameter name="schwarz: use reordering" type="bool" value="true"/> <Parameter name="subdomain solver name" type="string" value="RILUK"/> </ParameterList> </ParameterList> <ParameterList name="gmres"> <Parameter name="Type" type="string" value="Belos"/> <Parameter name="Solver Type" type="string" value="Pseudo Block GMRES"/> <ParameterList name="Solver Types"> <ParameterList name="Pseudo Block GMRES"> <Parameter name="Maximum Iterations" type="int" value="30"/> <Parameter name="Num Blocks" type="int" value="30"/> <Parameter name="Convergence Tolerance" type="double" value="1e-2"/> <Parameter name="Orthogonalization" type="string" value="ICGS"/> </ParameterList> </ParameterList> </ParameterList> </ParameterList>Finally, the contents of
belos.xmlrequired to use flexible GMRES are shown below:<ParameterList name="OMD_linear_solver"> <Parameter name="solution method" type="string" value="Flexible GMRES"/> <ParameterList name="Belos-Tpetra"> <Parameter name="Flexible Gmres" type="bool" value="true"/> <Parameter name="Maximum Iterations" type="int" value="100"/> <Parameter name="Num Blocks" type="int" value="100"/> <Parameter name="Implicit Residual Scaling" type="string" value="Norm of Initial Residual"/> <Parameter name="Explicit Residual Scaling" type="string" value="Norm of Initial Residual"/> <Parameter name="Convergence Tolerance" type="double" value="1e-8"/> <Parameter name="Orthogonalization" type="string" value="ICGS"/> </ParameterList> </ParameterList>Running the first step and inspecting the logfile, we observe that we have a fast, robust block preconditioner:
Equation System all->main: * Step : Transient, Strategy: LINE_SEARCH, Time: 1.00e-04, Step: 1.00e-04 * Matrix: Solver: "gmres_teko", Unknowns: 12270, Nonzeros: 688594 * Mesh : Processor 0 of 1: 2396 of 2396 elems, 2514 of 2514 nodes * Computing View Factors for enclosure encl_Enc1 N O N L I N E A R L I N E A R T I M E ------------------------------- -------------------- ----------------------- Step SubItns Resid Delta Itns Status Resid Asm/Slv/SubItn Time ---- ------- -------- -------- ---- ------ -------- ----------------------- 1 1 2.16e+02 1.48e+01 4 ok 6.08e-10 6.9e-02/4.2e-01/2.3e-02 2 1 7.13e-03 3.36e-04 4 ok 1.91e-09 6.7e-02/4.0e-01/2.2e-02 3 0 2.44e-06 9.08e-09 4 ok 1.39e-09 6.6e-02/4.0e-01/2.6e-04 Termination reason: 9.0814e-09 < nonlinear_correction_tolerance(1e-06)
Now that we have completed steps 1-4 as outlined in the Road-map for Constructing Teko Solvers, we have constructed a good Teko preconditioner for our OMD multiphysics problem.