2025-03-18 11:18:46 -04:00
/* ***********************************************************************
//
// Copyright (C) 2025 -- The 4D-STAR Collaboration
// File Author: Emily Boudreaux
2025-04-21 08:56:45 -04:00
// Last Modified: April 21, 2025
2025-03-18 11:18:46 -04:00
//
// 4DSSE is free software; you can use it and/or modify
// it under the terms and restrictions the GNU General Library Public
// License version 3 (GPLv3) as published by the Free Software Foundation.
//
// 4DSSE is distributed in the hope that it will be useful,
// but WITHOUT ANY WARRANTY; without even the implied warranty of
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
// See the GNU Library General Public License for more details.
//
// You should have received a copy of the GNU Library General Public License
// along with this software; if not, write to the Free Software
// Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
//
// *********************************************************************** */
2025-04-21 08:54:59 -04:00
# include "polySolver.h"
2025-02-19 14:35:15 -05:00
2025-03-19 13:49:21 -04:00
# include <memory>
2025-03-19 10:09:37 -04:00
# include <stdexcept>
2025-04-09 15:17:55 -04:00
# include <string>
2025-04-02 14:57:37 -04:00
# include <utility>
2025-03-18 10:15:51 -04:00
2025-04-21 08:54:59 -04:00
# include "mfem.hpp"
2025-04-09 15:17:55 -04:00
# include "4DSTARTypes.h"
# include "config.h"
2025-04-02 14:57:37 -04:00
# include "integrators.h"
2025-04-21 08:54:59 -04:00
# include "mfem.hpp"
2025-04-09 15:17:55 -04:00
# include "operator.h"
2025-02-19 14:35:15 -05:00
# include "polyCoeff.h"
2025-04-02 14:57:37 -04:00
# include "probe.h"
# include "resourceManager.h"
# include "resourceManagerTypes.h"
2025-03-03 09:54:13 -05:00
# include "quill/LogMacros.h"
2025-02-19 14:35:15 -05:00
2025-02-20 16:05:02 -05:00
2025-03-03 09:54:13 -05:00
namespace laneEmden {
2025-04-21 09:09:09 -04:00
double a ( const int k , const double n ) { // NOLINT(*-no-recursion)
2025-03-03 09:54:13 -05:00
if ( k = = 0 ) { return 1 ; }
if ( k = = 1 ) { return 0 ; }
else { return - ( c ( k - 2 , n ) / ( std : : pow ( k , 2 ) + k ) ) ; }
}
2025-04-21 09:09:09 -04:00
double c ( const int m , const double n ) { // NOLINT(*-no-recursion)
2025-03-03 09:54:13 -05:00
if ( m = = 0 ) { return std : : pow ( a ( 0 , n ) , n ) ; }
else {
double termOne = 1.0 / ( m * a ( 0 , n ) ) ;
double acc = 0 ;
for ( int k = 1 ; k < = m ; k + + ) {
acc + = ( k * n - m + k ) * a ( k , n ) * c ( m - k , n ) ;
}
return termOne * acc ;
}
}
2025-04-21 09:09:09 -04:00
double thetaSeriesExpansion ( const double xi , const double n , const int order ) {
2025-03-03 09:54:13 -05:00
double acc = 0 ;
for ( int k = 0 ; k < order ; k + + ) {
acc + = a ( k , n ) * std : : pow ( xi , k ) ;
}
return acc ;
}
}
2025-02-19 14:35:15 -05:00
2025-04-21 08:35:29 -04:00
2025-04-21 09:22:21 -04:00
PolySolver : : PolySolver ( const double n , const double order ) {
2025-03-05 12:55:53 -05:00
2025-03-19 10:09:37 -04:00
// --- Check the polytropic index ---
if ( n > 4.99 | | n < 0.0 ) {
2025-04-02 14:57:37 -04:00
LOG_ERROR ( m_logger , " The polytropic index n must be less than 5.0 and greater than 0.0. Currently it is {} " , m_polytropicIndex ) ;
throw std : : runtime_error ( " The polytropic index n must be less than 5.0 and greater than 0.0. Currently it is " + std : : to_string ( m_polytropicIndex ) ) ;
2025-03-19 10:09:37 -04:00
}
2025-04-02 14:57:37 -04:00
m_polytropicIndex = n ;
m_feOrder = order ;
2025-03-05 12:55:53 -05:00
2025-04-21 09:22:21 -04:00
// Load and rescale the mesh
const ResourceManager & rm = ResourceManager : : getInstance ( ) ;
const Resource & genericResource = rm . getResource ( " mesh:polySphere " ) ;
const auto & meshIO = std : : get < std : : unique_ptr < MeshIO > > ( genericResource ) ;
2025-04-02 14:57:37 -04:00
meshIO - > LinearRescale ( polycoeff : : x1 ( m_polytropicIndex ) ) ;
m_mesh = std : : make_unique < mfem : : Mesh > ( meshIO - > GetMesh ( ) ) ;
2025-03-03 09:54:13 -05:00
2025-04-02 14:57:37 -04:00
// Use feOrder - 1 for the RT space to satisfy Brezzi-Babuska condition
2025-04-21 09:22:21 -04:00
// for the H1 and RT [H(div)] spaces
2025-04-02 14:57:37 -04:00
m_fecH1 = std : : make_unique < mfem : : H1_FECollection > ( m_feOrder , m_mesh - > SpaceDimension ( ) ) ;
m_fecRT = std : : make_unique < mfem : : RT_FECollection > ( m_feOrder - 1 , m_mesh - > SpaceDimension ( ) ) ;
2025-02-19 14:35:15 -05:00
2025-04-02 14:57:37 -04:00
m_feTheta = std : : make_unique < mfem : : FiniteElementSpace > ( m_mesh . get ( ) , m_fecH1 . get ( ) ) ;
m_fePhi = std : : make_unique < mfem : : FiniteElementSpace > ( m_mesh . get ( ) , m_fecRT . get ( ) ) ;
2025-02-20 15:28:00 -05:00
2025-04-02 14:57:37 -04:00
m_theta = std : : make_unique < mfem : : GridFunction > ( m_feTheta . get ( ) ) ;
m_phi = std : : make_unique < mfem : : GridFunction > ( m_fePhi . get ( ) ) ;
2025-02-19 14:35:15 -05:00
2025-04-02 14:57:37 -04:00
assembleBlockSystem ( ) ;
2025-02-19 14:35:15 -05:00
2025-04-02 14:57:37 -04:00
}
2025-04-09 15:17:55 -04:00
PolySolver : : ~ PolySolver ( ) = default ;
2025-04-02 14:57:37 -04:00
void PolySolver : : assembleBlockSystem ( ) {
// Start by defining the block structure of the system
// Block 0: Theta (scalar space, uses m_feTheta)
// Block 1: Phi (vector space, uses m_fePhi)
mfem : : Array < mfem : : FiniteElementSpace * > feSpaces ;
feSpaces . Append ( m_feTheta . get ( ) ) ;
feSpaces . Append ( m_fePhi . get ( ) ) ;
// Create the block offsets. These define the start of each block in the combined vector.
// Block offsets will be [0, thetaDofs, thetaDofs + phiDofs]
mfem : : Array < int > blockOffsets ;
blockOffsets . SetSize ( 3 ) ;
blockOffsets [ 0 ] = 0 ;
blockOffsets [ 1 ] = feSpaces [ 0 ] - > GetVSize ( ) ;
blockOffsets [ 2 ] = feSpaces [ 1 ] - > GetVSize ( ) ;
2025-04-21 09:22:21 -04:00
blockOffsets . PartialSum ( ) ; // Cumulative sum to get the offsets
2025-04-02 14:57:37 -04:00
// Add integrators to block form one by one
2025-04-09 15:17:55 -04:00
// We add integrators corresponding to each term in the weak form
2025-04-02 14:57:37 -04:00
// The block form of the residual matrix
// ⎡ 0 -M ⎤ ⎡ θ ⎤ + ⎡f(θ)⎤ = ⎡ 0 ⎤ = R(X)
// ⎣ -Q D ⎦ ⎣ Φ ⎦ ⎣ 0 ⎦ ⎣ 0 ⎦
// This then simplifies to
// ⎡f(θ) - MΘ ⎤ = ⎡ 0 ⎤ = R
// ⎣ -Qɸ Dθ ⎦ ⎣ 0 ⎦
// Here M, Q, and D are
// M = ∫∇ψᶿ·Nᵠ dV (MixedVectorWeakDivergenceIntegrator)
// D = ∫ψᵠ·Nᵠ dV (VectorFEMassIntegrator)
// Q = ∫ψᵠ·∇Nᶿ dV (MixedVectorGradientIntegrator)
// f(θ) = ∫ψᶿ·θⁿ dV (NonlinearPowerIntegrator)
// Here ψᶿ and ψᵠ are the test functions for the theta and phi spaces, respectively
// Nᵠ and Nᶿ are the basis functions for the theta and phi spaces, respectively
// A full derivation of the weak form can be found in the 4DSSE documentation
// --- Assemble the MixedBilinear and Bilinear forms (M, D, and Q) ---
2025-04-03 11:14:50 -04:00
auto Mform = std : : make_unique < mfem : : MixedBilinearForm > ( m_fePhi . get ( ) , m_feTheta . get ( ) ) ;
auto Qform = std : : make_unique < mfem : : MixedBilinearForm > ( m_feTheta . get ( ) , m_fePhi . get ( ) ) ;
2025-04-02 14:57:37 -04:00
auto Dform = std : : make_unique < mfem : : BilinearForm > ( m_fePhi . get ( ) ) ;
// TODO: Check the sign on all of the integrators
2025-04-03 11:14:50 -04:00
Mform - > AddDomainIntegrator ( new mfem : : MixedVectorWeakDivergenceIntegrator ( ) ) ;
Qform - > AddDomainIntegrator ( new mfem : : MixedVectorGradientIntegrator ( ) ) ;
Dform - > AddDomainIntegrator ( new mfem : : VectorFEMassIntegrator ( ) ) ;
2025-04-02 14:57:37 -04:00
Mform - > Assemble ( ) ;
Mform - > Finalize ( ) ;
Qform - > Assemble ( ) ;
Qform - > Finalize ( ) ;
Dform - > Assemble ( ) ;
Dform - > Finalize ( ) ;
// --- Assemble the NonlinearForm (f) ---
2025-04-21 09:22:21 -04:00
// Note that the nonlinear form is built here but its essential true dofs (boundary conditions) are
// not set until later (when solve is called). They are applied through the operator rather than directly.
2025-04-02 14:57:37 -04:00
auto fform = std : : make_unique < mfem : : NonlinearForm > ( m_feTheta . get ( ) ) ;
2025-04-21 08:54:59 -04:00
fform - > AddDomainIntegrator ( new polyMFEMUtils : : NonlinearPowerIntegrator ( m_polytropicIndex ) ) ;
2025-04-02 14:57:37 -04:00
// -- Build the BlockOperator --
m_polytropOperator = std : : make_unique < PolytropeOperator > (
std : : move ( Mform ) ,
std : : move ( Qform ) ,
std : : move ( Dform ) ,
std : : move ( fform ) ,
blockOffsets
) ;
2025-03-19 13:49:21 -04:00
2025-02-19 14:35:15 -05:00
}
2025-04-09 15:17:55 -04:00
void PolySolver : : solve ( ) const {
2025-02-19 14:35:15 -05:00
// --- Set the initial guess for the solution ---
2025-04-02 14:57:37 -04:00
setInitialGuess ( ) ;
2025-03-19 10:09:37 -04:00
2025-04-21 09:05:34 -04:00
setOperatorEssentialTrueDofs ( ) ;
2025-04-02 14:57:37 -04:00
2025-04-03 11:14:50 -04:00
// It's safer to get the offsets directly from the operator after finalization
const mfem : : Array < int > & block_offsets = m_polytropOperator - > GetBlockOffsets ( ) ; // Assuming a getter exists or accessing member if public/friend
mfem : : BlockVector state_vector ( block_offsets ) ;
2025-04-21 09:05:34 -04:00
state_vector . GetBlock ( 0 ) = static_cast < mfem : : Vector > ( * m_theta ) ; // NOLINT(*-slicing)
state_vector . GetBlock ( 1 ) = static_cast < mfem : : Vector > ( * m_phi ) ; // NOLINT(*-slicing)
2025-04-02 14:57:37 -04:00
2025-04-03 11:14:50 -04:00
mfem : : Vector zero_rhs ( block_offsets . Last ( ) ) ;
zero_rhs = 0.0 ;
2025-04-02 14:57:37 -04:00
2025-04-21 09:22:21 -04:00
const solverBundle sb = setupNewtonSolver ( ) ;
2025-03-03 09:54:13 -05:00
2025-04-09 15:17:55 -04:00
// EMB 2025: Calling Mult gets the gradient of the operator for the NewtonSolver
// This then is set as the operator for the solver for NewtonSolver
// The solver (assuming it is an iterative solver) then sets the
// operator for its preconditioner to this.
// What this means is that there is no need to manually deal
// with updating the preconditioner at every newton step as the
// changes to the jacobian are automatically propagated through the
// solving chain. This is at least true with MFEM 4.8-rc0
2025-04-21 08:35:29 -04:00
sb . newton . Mult ( zero_rhs , state_vector ) ;
2025-04-02 14:57:37 -04:00
// --- Save and view the solution ---
2025-04-03 11:14:50 -04:00
saveAndViewSolution ( state_vector ) ;
2025-04-02 14:57:37 -04:00
}
2025-04-09 15:17:55 -04:00
SSE : : MFEMArrayPairSet PolySolver : : getEssentialTrueDof ( ) const {
2025-04-02 14:57:37 -04:00
mfem : : Array < int > theta_ess_tdof_list ;
mfem : : Array < int > phi_ess_tdof_list ;
2025-03-03 09:54:13 -05:00
2025-04-09 15:17:55 -04:00
mfem : : Array < int > thetaCenterDofs , phiCenterDofs ;
mfem : : Array < double > thetaCenterVals , phiCenterVals ;
std : : tie ( thetaCenterDofs , phi_ess_tdof_list ) = findCenterElement ( ) ;
thetaCenterVals . SetSize ( thetaCenterDofs . Size ( ) ) ;
phiCenterVals . SetSize ( phi_ess_tdof_list . Size ( ) ) ;
2025-04-02 14:57:37 -04:00
2025-04-09 15:17:55 -04:00
thetaCenterVals = 1.0 ;
phiCenterVals = 0.0 ;
2025-04-02 14:57:37 -04:00
mfem : : Array < int > ess_brd ( m_mesh - > bdr_attributes . Max ( ) ) ;
ess_brd = 1 ;
2025-04-09 15:17:55 -04:00
mfem : : Array < double > thetaSurfaceVals ;
2025-04-02 14:57:37 -04:00
m_feTheta - > GetEssentialTrueDofs ( ess_brd , theta_ess_tdof_list ) ;
2025-04-09 15:17:55 -04:00
thetaSurfaceVals . SetSize ( theta_ess_tdof_list . Size ( ) ) ;
thetaSurfaceVals = 0.0 ;
2025-04-02 14:57:37 -04:00
// combine the essential dofs with the center dofs
2025-04-09 15:17:55 -04:00
theta_ess_tdof_list . Append ( thetaCenterDofs ) ;
thetaSurfaceVals . Append ( thetaCenterVals ) ;
SSE : : MFEMArrayPair thetaPair = std : : make_pair ( theta_ess_tdof_list , thetaSurfaceVals ) ;
SSE : : MFEMArrayPair phiPair = std : : make_pair ( phi_ess_tdof_list , phiCenterVals ) ;
SSE : : MFEMArrayPairSet pairSet = std : : make_pair ( thetaPair , phiPair ) ;
2025-04-02 14:57:37 -04:00
2025-04-09 15:17:55 -04:00
return pairSet ;
2025-04-02 14:57:37 -04:00
}
2025-04-09 15:17:55 -04:00
std : : pair < mfem : : Array < int > , mfem : : Array < int > > PolySolver : : findCenterElement ( ) const {
mfem : : Array < int > thetaCenterDofs ;
mfem : : Array < int > phiCenterDofs ;
2025-04-02 14:57:37 -04:00
mfem : : DenseMatrix centerPoint ( m_mesh - > SpaceDimension ( ) , 1 ) ;
centerPoint ( 0 , 0 ) = 0.0 ;
centerPoint ( 1 , 0 ) = 0.0 ;
centerPoint ( 2 , 0 ) = 0.0 ;
mfem : : Array < int > elementIDs ;
mfem : : Array < mfem : : IntegrationPoint > ips ;
m_mesh - > FindPoints ( centerPoint , elementIDs , ips ) ;
mfem : : Array < int > tempDofs ;
for ( int i = 0 ; i < elementIDs . Size ( ) ; i + + ) {
m_feTheta - > GetElementDofs ( elementIDs [ i ] , tempDofs ) ;
2025-04-09 15:17:55 -04:00
thetaCenterDofs . Append ( tempDofs ) ;
m_fePhi - > GetElementDofs ( elementIDs [ i ] , tempDofs ) ;
phiCenterDofs . Append ( tempDofs ) ;
2025-04-02 14:57:37 -04:00
}
2025-04-09 15:17:55 -04:00
return std : : make_pair ( thetaCenterDofs , phiCenterDofs ) ;
2025-04-02 14:57:37 -04:00
}
2025-04-09 15:17:55 -04:00
void PolySolver : : setInitialGuess ( ) const {
2025-04-02 14:57:37 -04:00
// --- Set the initial guess for the solution ---
mfem : : FunctionCoefficient thetaInitGuess (
[ this ] ( const mfem : : Vector & x ) {
2025-04-21 09:22:21 -04:00
const double r = x . Norml2 ( ) ;
const double radius = Probe : : getMeshRadius ( * m_mesh ) ;
const double u = 1 / radius ;
2025-04-02 14:57:37 -04:00
2025-04-21 09:22:21 -04:00
return - std : : pow ( ( u * r ) , 2 ) + 1.0 ; // The series expansion is a better guess; however, this is cheaper and ensures that the value at the surface is very close to zero in a way that the series expansion does not
2025-04-02 14:57:37 -04:00
}
) ;
m_theta - > ProjectCoefficient ( thetaInitGuess ) ;
2025-04-18 11:18:55 -04:00
mfem : : GradientGridFunctionCoefficient phiInitGuess ( m_theta . get ( ) ) ;
2025-04-02 14:57:37 -04:00
m_phi - > ProjectCoefficient ( phiInitGuess ) ;
if ( m_config . get < bool > ( " Poly:Solver:ViewInitialGuess " , false ) ) {
2025-04-09 15:17:55 -04:00
Probe : : glVisView ( * m_theta , * m_mesh , " θ init " ) ;
Probe : : glVisView ( * m_phi , * m_mesh , " ɸ init " ) ;
2025-04-02 14:57:37 -04:00
}
}
2025-04-09 15:17:55 -04:00
void PolySolver : : saveAndViewSolution ( const mfem : : BlockVector & state_vector ) const {
2025-04-03 11:14:50 -04:00
mfem : : BlockVector x_block ( const_cast < mfem : : BlockVector & > ( state_vector ) , m_polytropOperator - > GetBlockOffsets ( ) ) ;
mfem : : Vector & x_theta = x_block . GetBlock ( 0 ) ;
mfem : : Vector & x_phi = x_block . GetBlock ( 1 ) ;
2025-04-09 15:17:55 -04:00
if ( m_config . get < bool > ( " Poly:Output:View " , false ) ) {
2025-04-03 11:14:50 -04:00
Probe : : glVisView ( x_theta , * m_feTheta , " θ Solution " ) ;
Probe : : glVisView ( x_phi , * m_fePhi , " ɸ Solution " ) ;
2025-04-02 14:57:37 -04:00
}
2025-03-03 09:54:13 -05:00
2025-03-05 12:55:53 -05:00
// --- Extract the Solution ---
2025-04-21 08:54:59 -04:00
if ( m_config . get < bool > ( " Poly:Output:1D:Save " , true ) ) {
2025-04-09 15:17:55 -04:00
auto solutionPath = m_config . get < std : : string > ( " Poly:Output:1D:Path " , " polytropeSolution_1D.csv " ) ;
auto derivSolPath = " d " + solutionPath ;
2025-04-03 11:14:50 -04:00
2025-04-09 15:17:55 -04:00
auto rayCoLatitude = m_config . get < double > ( " Poly:Output:1D:RayCoLatitude " , 0.0 ) ;
auto rayLongitude = m_config . get < double > ( " Poly:Output:1D:RayLongitude " , 0.0 ) ;
auto raySamples = m_config . get < int > ( " Poly:Output:1D:RaySamples " , 100 ) ;
2025-03-03 09:54:13 -05:00
2025-03-05 12:55:53 -05:00
std : : vector rayDirection = { rayCoLatitude , rayLongitude } ;
2025-03-03 09:54:13 -05:00
2025-04-03 11:14:50 -04:00
Probe : : getRaySolution ( x_theta , * m_feTheta , rayDirection , raySamples , solutionPath ) ;
// Probe::getRaySolution(x_phi, *m_fePhi, rayDirection, raySamples, derivSolPath);
2025-03-03 09:54:13 -05:00
}
2025-04-03 11:14:50 -04:00
}
2025-04-21 09:05:34 -04:00
void PolySolver : : setOperatorEssentialTrueDofs ( ) const {
2025-04-09 15:17:55 -04:00
SSE : : MFEMArrayPairSet ess_tdof_pair_set = getEssentialTrueDof ( ) ;
m_polytropOperator - > SetEssentialTrueDofs ( ess_tdof_pair_set ) ;
2025-04-03 11:14:50 -04:00
// -- Finalize the operator --
m_polytropOperator - > finalize ( ) ;
if ( ! m_polytropOperator - > isFinalized ( ) ) {
LOG_ERROR ( m_logger , " PolytropeOperator is not finalized. Cannot solve. " ) ;
throw std : : runtime_error ( " PolytropeOperator is not finalized. Cannot solve. " ) ;
}
}
2025-04-09 15:17:55 -04:00
void PolySolver : : LoadSolverUserParams ( double & newtonRelTol , double & newtonAbsTol , int & newtonMaxIter , int & newtonPrintLevel , double & gmresRelTol , double & gmresAbsTol , int & gmresMaxIter , int & gmresPrintLevel ) const {
newtonRelTol = m_config . get < double > ( " Poly:Solver:Newton:RelTol " , 1e-7 ) ;
newtonAbsTol = m_config . get < double > ( " Poly:Solver:Newton:AbsTol " , 1e-7 ) ;
newtonMaxIter = m_config . get < int > ( " Poly:Solver:Newton:MaxIter " , 200 ) ;
newtonPrintLevel = m_config . get < int > ( " Poly:Solver:Newton:PrintLevel " , 1 ) ;
2025-04-03 11:14:50 -04:00
2025-04-09 15:17:55 -04:00
gmresRelTol = m_config . get < double > ( " Poly:Solver:GMRES:RelTol " , 1e-10 ) ;
gmresAbsTol = m_config . get < double > ( " Poly:Solver:GMRES:AbsTol " , 1e-12 ) ;
gmresMaxIter = m_config . get < int > ( " Poly:Solver:GMRES:MaxIter " , 2000 ) ;
gmresPrintLevel = m_config . get < int > ( " Poly:Solver:GMRES:PrintLevel " , 0 ) ;
2025-04-03 11:14:50 -04:00
LOG_DEBUG ( m_logger , " Newton Solver (relTol: {:0.2E}, absTol: {:0.2E}, maxIter: {}, printLevel: {}) " , newtonRelTol , newtonAbsTol , newtonMaxIter , newtonPrintLevel ) ;
LOG_DEBUG ( m_logger , " GMRES Solver (relTol: {:0.2E}, absTol: {:0.2E}, maxIter: {}, printLevel: {}) " , gmresRelTol , gmresAbsTol , gmresMaxIter , gmresPrintLevel ) ;
2025-04-09 15:17:55 -04:00
}
2025-04-21 08:35:29 -04:00
solverBundle PolySolver : : setupNewtonSolver ( ) const {
2025-04-09 15:17:55 -04:00
// --- Load configuration parameters ---
double newtonRelTol , newtonAbsTol , gmresRelTol , gmresAbsTol ;
int newtonMaxIter , newtonPrintLevel , gmresMaxIter , gmresPrintLevel ;
LoadSolverUserParams ( newtonRelTol , newtonAbsTol , newtonMaxIter , newtonPrintLevel , gmresRelTol , gmresAbsTol ,
gmresMaxIter , gmresPrintLevel ) ;
2025-04-03 11:14:50 -04:00
2025-04-21 08:54:59 -04:00
solverBundle solver ; // Use this solver bundle to ensure lifetime safety
2025-04-21 08:35:29 -04:00
solver . solver . SetRelTol ( gmresRelTol ) ;
solver . solver . SetAbsTol ( gmresAbsTol ) ;
solver . solver . SetMaxIter ( gmresMaxIter ) ;
solver . solver . SetPrintLevel ( gmresPrintLevel ) ;
solver . solver . SetPreconditioner ( m_polytropOperator - > GetPreconditioner ( ) ) ;
2025-04-03 11:14:50 -04:00
// --- Set up the Newton solver ---
2025-04-21 08:35:29 -04:00
solver . newton . SetRelTol ( newtonRelTol ) ;
solver . newton . SetAbsTol ( newtonAbsTol ) ;
solver . newton . SetMaxIter ( newtonMaxIter ) ;
solver . newton . SetPrintLevel ( newtonPrintLevel ) ;
solver . newton . SetOperator ( * m_polytropOperator ) ;
2025-04-09 15:17:55 -04:00
// --- Created the linear solver which is used to invert the jacobian ---
2025-04-21 08:35:29 -04:00
solver . newton . SetSolver ( solver . solver ) ;
2025-04-03 11:14:50 -04:00
2025-04-21 08:35:29 -04:00
return solver ;
2025-04-09 15:17:55 -04:00
}