GEOS
CoupledSolver.hpp
Go to the documentation of this file.
1 /*
2  * ------------------------------------------------------------------------------------------------------------
3  * SPDX-License-Identifier: LGPL-2.1-only
4  *
5  * Copyright (c) 2016-2024 Lawrence Livermore National Security LLC
6  * Copyright (c) 2018-2024 Total, S.A
7  * Copyright (c) 2018-2024 The Board of Trustees of the Leland Stanford Junior University
8  * Copyright (c) 2023-2024 Chevron
9  * Copyright (c) 2019- GEOS/GEOSX Contributors
10  * All rights reserved
11  *
12  * See top level LICENSE, COPYRIGHT, CONTRIBUTORS, NOTICE, and ACKNOWLEDGEMENTS files for details.
13  * ------------------------------------------------------------------------------------------------------------
14  */
15 
21 #ifndef GEOS_PHYSICSSOLVERS_MULTIPHYSICS_COUPLEDSOLVER_HPP_
22 #define GEOS_PHYSICSSOLVERS_MULTIPHYSICS_COUPLEDSOLVER_HPP_
23 
26 
27 #include <tuple>
28 
29 namespace geos
30 {
31 
32 template< typename ... SOLVERS >
34 {
35 
36 public:
37 
43  CoupledSolver( const string & name,
44  Group * const parent )
45  : PhysicsSolverBase( name, parent )
46  {
47  forEachArgInTuple( m_solvers, [&]( auto solver, auto idx )
48  {
49  using SolverType = TYPEOFPTR( solver );
50  string const key = SolverType::coupledSolverAttributePrefix() + "SolverName";
51  registerWrapper( key, &m_names[idx()] ).
52  setRTTypeName( rtTypes::CustomTypes::groupNameRef ).
54  setDescription( "Name of the " + SolverType::coupledSolverAttributePrefix() + " solver used by the coupled solver" );
55  } );
56 
57  this->getWrapper< string >( PhysicsSolverBase::viewKeyStruct::discretizationString() ).
58  setInputFlag( dataRepository::InputFlags::FALSE );
59 
60  addLogLevel< logInfo::Coupling >();
61  }
62 
64  CoupledSolver( CoupledSolver const & ) = delete;
65 
67  CoupledSolver( CoupledSolver && ) = default;
68 
70  CoupledSolver & operator=( CoupledSolver const & ) = delete;
71 
74 
75 
79  void
81  {
82  forEachArgInTuple( m_solvers, [&]( auto & solver, auto idx )
83  {
84  using SolverPtr = TYPEOFREF( solver );
85  using SolverType = TYPEOFPTR( SolverPtr {} );
86  auto const & solverName = m_names[idx()];
87  auto const & solverType = LvArray::system::demangleType< SolverType >();
88  solver = this->getParent().template getGroupPointer< SolverType >( solverName );
89  GEOS_THROW_IF( solver == nullptr,
90  GEOS_FMT( "{}: Could not find solver '{}' of type {}",
92  solverName, solverType ),
93  InputError );
94  GEOS_LOG_LEVEL_INFO_RANK_0( logInfo::Coupling, GEOS_FMT( "{}: found {} solver named {}", getName(), solver->getCatalogName(), solverName ) );
95  } );
96  }
97 
98 
104  virtual void
106  DofManager & dofManager ) const
107  { GEOS_UNUSED_VAR( domain, dofManager ); }
108 
118  virtual void
120  real64 const dt,
121  DomainPartition const & domain,
122  DofManager const & dofManager,
123  CRSMatrixView< real64, globalIndex const > const & localMatrix,
124  arrayView1d< real64 > const & localRhs )
125  { GEOS_UNUSED_VAR( time_n, dt, domain, dofManager, localMatrix, localRhs ); }
126 
134  void
135  setupDofs( DomainPartition const & domain,
136  DofManager & dofManager ) const override
137  {
138  forEachArgInTuple( m_solvers, [&]( auto & solver, auto )
139  {
140  solver->setupDofs( domain, dofManager );
141  } );
142 
143  setupCoupling( domain, dofManager );
144  }
145 
146  virtual void
147  implicitStepSetup( real64 const & time_n,
148  real64 const & dt,
149  DomainPartition & domain ) override
150  {
151  forEachArgInTuple( m_solvers, [&]( auto & solver, auto )
152  {
153  solver->implicitStepSetup( time_n, dt, domain );
154  } );
155  }
156 
157  virtual void
158  implicitStepComplete( real64 const & time_n,
159  real64 const & dt,
160  DomainPartition & domain ) override
161  {
162  forEachArgInTuple( m_solvers, [&]( auto & solver, auto )
163  {
164  solver->implicitStepComplete( time_n, dt, domain );
165  } );
166  }
167 
168  // general version of assembleSystem function, keep in mind many solvers will override it
169  virtual void
170  assembleSystem( real64 const time_n,
171  real64 const dt,
172  DomainPartition & domain,
173  DofManager const & dofManager,
174  CRSMatrixView< real64, globalIndex const > const & localMatrix,
175  arrayView1d< real64 > const & localRhs ) override
176  {
178 
179  // 1. Assemble matrix blocks of each individual solver
180  forEachArgInTuple( m_solvers, [&]( auto & solver, auto )
181  {
182  solver->assembleSystem( time_n, dt, domain, dofManager, localMatrix, localRhs );
183  } );
184 
185  // 2. Assemble coupling blocks
186  assembleCouplingTerms( time_n, dt, domain, dofManager, localMatrix, localRhs );
187  }
188 
189  virtual void
190  applySystemSolution( DofManager const & dofManager,
191  arrayView1d< real64 const > const & localSolution,
192  real64 const scalingFactor,
193  real64 const dt,
194  DomainPartition & domain ) override
195  {
196  forEachArgInTuple( m_solvers, [&]( auto & solver, auto )
197  {
198  solver->applySystemSolution( dofManager, localSolution, scalingFactor, dt, domain );
199  } );
200  }
201 
202  virtual void
203  updateState( DomainPartition & domain ) override
204  {
205  forEachArgInTuple( m_solvers, [&]( auto & solver, auto )
206  {
207  solver->updateState( domain );
208  } );
209  }
210 
211  virtual void
213  {
214  forEachArgInTuple( m_solvers, [&]( auto & solver, auto )
215  {
216  solver->resetStateToBeginningOfStep( domain );
217  } );
218  }
219 
222  real64
223  solverStep( real64 const & time_n,
224  real64 const & dt,
225  int const cycleNumber,
226  DomainPartition & domain ) override final
227  {
229 
231  {
232  return fullyCoupledSolverStep( time_n, dt, cycleNumber, domain );
233  }
235  {
236  return sequentiallyCoupledSolverStep( time_n, dt, cycleNumber, domain );
237  }
238  else
239  {
240  GEOS_ERROR( getDataContext() << ": Invalid coupling type option." );
241  return 0;
242  }
243 
244  }
245 
246 
247  virtual real64
248  calculateResidualNorm( real64 const & time_n,
249  real64 const & dt,
250  DomainPartition const & domain,
251  DofManager const & dofManager,
252  arrayView1d< real64 const > const & localRhs ) override
253  {
254  real64 norm = 0.0;
255  forEachArgInTuple( m_solvers, [&]( auto & solver, auto )
256  {
257  real64 const singlePhysicsNorm = solver->calculateResidualNorm( time_n, dt, domain, dofManager, localRhs );
258  norm += singlePhysicsNorm * singlePhysicsNorm;
259  } );
260 
261  return sqrt( norm );
262  }
263 
264  virtual void
266  real64 const dt,
267  DomainPartition & domain,
268  DofManager const & dofManager,
269  CRSMatrixView< real64, globalIndex const > const & localMatrix,
270  arrayView1d< real64 > const & localRhs ) override
271  {
272  forEachArgInTuple( m_solvers, [&]( auto & solver, auto )
273  {
274  solver->applyBoundaryConditions( time_n, dt, domain, dofManager, localMatrix, localRhs );
275  } );
276  }
277 
278  virtual bool
280  DofManager const & dofManager,
281  arrayView1d< real64 const > const & localSolution,
282  real64 const scalingFactor ) override
283  {
284  bool validSolution = true;
285  forEachArgInTuple( m_solvers, [&]( auto & solver, auto )
286  {
287  bool const validSinglePhysicsSolution = solver->checkSystemSolution( domain, dofManager, localSolution, scalingFactor );
288  if( !validSinglePhysicsSolution )
289  {
290  GEOS_LOG_RANK_0( GEOS_FMT( " {}/{}: Solution check failed. Newton loop terminated.", getName(), solver->getName()) );
291  }
292  validSolution = validSolution && validSinglePhysicsSolution;
293  } );
294  return validSolution;
295  }
296 
297  virtual real64
299  DofManager const & dofManager,
300  arrayView1d< real64 const > const & localSolution ) override
301  {
302  real64 scalingFactor = 1e9;
303  forEachArgInTuple( m_solvers, [&]( auto & solver, auto )
304  {
305  real64 const singlePhysicsScalingFactor = solver->scalingForSystemSolution( domain, dofManager, localSolution );
306  scalingFactor = LvArray::math::min( scalingFactor, singlePhysicsScalingFactor );
307  } );
308  return scalingFactor;
309  }
310 
311  virtual real64
313  DomainPartition & domain ) override
314  {
315  real64 nextDt = LvArray::NumericLimits< real64 >::max;
316  forEachArgInTuple( m_solvers, [&]( auto & solver, auto )
317  {
318  real64 const singlePhysicsNextDt =
319  solver->setNextDtBasedOnStateChange( currentDt, domain );
320  nextDt = LvArray::math::min( singlePhysicsNextDt, nextDt );
321  } );
322  return nextDt;
323  }
324 
325  virtual real64 setNextDtBasedOnNewtonIter( real64 const & currentDt ) override
326  {
328  forEachArgInTuple( m_solvers, [&]( auto & solver, auto )
329  {
330  real64 const singlePhysicsNextDt =
331  solver->setNextDtBasedOnNewtonIter( currentDt );
332  nextDt = LvArray::math::min( singlePhysicsNextDt, nextDt );
333  } );
334  return nextDt;
335  }
336 
337  virtual void cleanup( real64 const time_n,
338  integer const cycleNumber,
339  integer const eventCounter,
340  real64 const eventProgress,
341  DomainPartition & domain ) override
342  {
343  forEachArgInTuple( m_solvers, [&]( auto & solver, auto )
344  {
345  solver->cleanup( time_n, cycleNumber, eventCounter, eventProgress, domain );
346  } );
347  PhysicsSolverBase::cleanup( time_n, cycleNumber, eventCounter, eventProgress, domain );
348  }
349 
352  virtual bool checkSequentialSolutionIncrements( DomainPartition & domain ) const override
353  {
354  bool isConverged = true;
355  forEachArgInTuple( m_solvers, [&]( auto & solver, auto )
356  {
357  isConverged &= solver->checkSequentialSolutionIncrements( domain );
358  } );
359  return isConverged;
360  }
361 
362  virtual bool updateConfiguration( DomainPartition & domain ) override
363  {
364  bool result = true;
365  forEachArgInTuple( m_solvers, [&]( auto & solver, auto )
366  {
367  result &= solver->updateConfiguration( domain );
368  } );
369  return result;
370  }
371 
372  virtual void outputConfigurationStatistics( DomainPartition const & domain ) const override
373  {
374  forEachArgInTuple( m_solvers, [&]( auto & solver, auto )
375  {
376  solver->outputConfigurationStatistics( domain );
377  } );
378  }
379 
380  virtual void resetConfigurationToBeginningOfStep( DomainPartition & domain ) override
381  {
382  forEachArgInTuple( m_solvers, [&]( auto & solver, auto )
383  {
384  solver->resetConfigurationToBeginningOfStep( domain );
385  } );
386  }
387 
388  virtual bool resetConfigurationToDefault( DomainPartition & domain ) const override
389  {
390  bool result = true;
391  forEachArgInTuple( m_solvers, [&]( auto & solver, auto )
392  {
393  result &= solver->resetConfigurationToDefault( domain );
394  } );
395  return result;
396  }
397 
398 protected:
399 
409  virtual real64 fullyCoupledSolverStep( real64 const & time_n,
410  real64 const & dt,
411  int const cycleNumber,
412  DomainPartition & domain )
413  {
414  return PhysicsSolverBase::solverStep( time_n, dt, cycleNumber, domain );
415  }
416 
427  virtual real64 sequentiallyCoupledSolverStep( real64 const & time_n,
428  real64 const & dt,
429  int const cycleNumber,
430  DomainPartition & domain )
431  {
433 
434  // Only build the sparsity pattern if the mesh has changed
435  Timestamp const meshModificationTimestamp = getMeshModificationTimestamp( domain );
436  forEachArgInTuple( m_solvers, [&]( auto & solver, auto )
437  {
438  if( meshModificationTimestamp > solver->getSystemSetupTimestamp() )
439  {
440  solver->setupSystem( domain,
441  solver->getDofManager(),
442  solver->getLocalMatrix(),
443  solver->getSystemRhs(),
444  solver->getSystemSolution() );
445  solver->setSystemSetupTimestamp( meshModificationTimestamp );
446  }
447  } );
448 
449  implicitStepSetup( time_n, dt, domain );
450 
452  integer const maxNumberDtCuts = solverParams.m_maxTimeStepCuts;
453  real64 const dtCutFactor = solverParams.m_timeStepCutFactor;
454  integer & dtAttempt = solverParams.m_numTimeStepAttempts;
455 
456  bool isConverged = false;
457  // dt may be cut during the course of this step, so we are keeping a local
458  // value to track the achieved dt for this step.
459  real64 stepDt = dt;
460 
461  // outer loop attempts to apply full timestep, and managed the cutting of the timestep if
462  // required.
463  for( dtAttempt = 0; dtAttempt < maxNumberDtCuts; ++dtAttempt )
464  {
465  // TODO configuration loop
466 
467  // Reset the states of all solvers if any of them had to restart
468  forEachArgInTuple( m_solvers, [&]( auto & solver, auto )
469  {
470  solver->resetStateToBeginningOfStep( domain );
471  solver->getSolverStatistics().initializeTimeStepStatistics(); // initialize counters for subsolvers
472  } );
473  resetStateToBeginningOfStep( domain );
474 
475  integer & iter = solverParams.m_numNewtonIterations;
476 
478  for( iter = 0; iter < solverParams.m_maxIterNewton; iter++ )
479  {
480  // Increment the solver statistics for reporting purposes
481  // Pass a "0" as argument (0 linear iteration) to skip the output of linear iteration stats at the end
483 
484  startSequentialIteration( iter, domain );
485 
486  // Solve the subproblems nonlinearly
487  forEachArgInTuple( m_solvers, [&]( auto & solver, auto idx )
488  {
489  GEOS_LOG_LEVEL_INFO_RANK_0( logInfo::NonlinearSolver, GEOS_FMT( " Iteration {:2}: {}", iter + 1, solver->getName() ) );
490  real64 solverDt = solver->nonlinearImplicitStep( time_n,
491  stepDt,
492  cycleNumber,
493  domain );
494 
495  // save fields (e.g. pressure and temperature) after inner solve
496  solver->saveSequentialIterationState( domain );
497 
498  mapSolutionBetweenSolvers( domain, idx() );
499 
500  if( solverDt < stepDt ) // subsolver had to cut the time step
501  {
502  iter = 0; // restart outer loop
503  stepDt = solverDt; // sync time step
504  }
505  } );
506 
507  // Check convergence of the outer loop
508  isConverged = checkSequentialConvergence( iter,
509  time_n,
510  stepDt,
511  domain );
512 
513  if( isConverged )
514  {
515  // we still want to count current iteration
516  ++iter;
517  // exit outer loop
518  break;
519  }
520  else
521  {
522  finishSequentialIteration( iter, domain );
523  }
524  }
525 
526  if( isConverged )
527  {
528  // Save time step statistics for the subsolvers
529  forEachArgInTuple( m_solvers, [&]( auto & solver, auto )
530  {
531  solver->getSolverStatistics().saveTimeStepStatistics();
532  } );
533  // get out of the time loop
534  break;
535  }
536  else
537  {
538  // cut timestep, go back to beginning of step and restart the Newton loop
539  stepDt *= dtCutFactor;
540  GEOS_LOG_LEVEL_INFO_RANK_0( logInfo::TimeStep, GEOS_FMT( "New dt = {}", stepDt ) );
541 
542  // notify the solver statistics counter that this is a time step cut
544  forEachArgInTuple( m_solvers, [&]( auto & solver, auto )
545  {
546  solver->getSolverStatistics().logTimeStepCut();
547  } );
548  }
549  }
550 
551  if( !isConverged )
552  {
553  GEOS_LOG_RANK_0( "Convergence not achieved." );
554 
556  {
557  GEOS_LOG_RANK_0( "The accepted solution may be inaccurate." );
558  }
559  else
560  {
561  GEOS_ERROR( "Nonconverged solutions not allowed. Terminating..." );
562  }
563  }
564 
565  implicitStepComplete( time_n, stepDt, domain );
566 
567  return stepDt;
568  }
569 
577  integer const solverType )
578  {
579  GEOS_UNUSED_VAR( domain, solverType );
580  }
581 
582  virtual bool checkSequentialConvergence( int const & iter,
583  real64 const & time_n,
584  real64 const & dt,
585  DomainPartition & domain )
586  {
588  bool isConverged = true;
589 
590  if( params.m_subcyclingOption == 0 )
591  {
592  GEOS_LOG_LEVEL_INFO_RANK_0( logInfo::Convergence, "***** Single Pass solver, no subcycling *****" );
593  }
594  else
595  {
596  GEOS_LOG_LEVEL_INFO_RANK_0( logInfo::Convergence, GEOS_FMT( " Iteration {:2}: outer-loop convergence check", iter + 1 ) );
597 
599  {
600  real64 residualNorm = 0;
601 
602  // loop over all the single-physics solvers
603  forEachArgInTuple( m_solvers, [&]( auto & solver, auto )
604  {
605 
606  solver->getLocalMatrix().toViewConstSizes().zero();
607  solver->getSystemRhs().zero();
608  arrayView1d< real64 > const localRhs = solver->getSystemRhs().open();
609 
610  // for each solver, we have to recompute the residual (and Jacobian, although not necessary)
611  solver->assembleSystem( time_n,
612  dt,
613  domain,
614  solver->getDofManager(),
615  solver->getLocalMatrix().toViewConstSizes(),
616  localRhs );
617  solver->applyBoundaryConditions( time_n,
618  dt,
619  domain,
620  solver->getDofManager(),
621  solver->getLocalMatrix().toViewConstSizes(),
622  localRhs );
623  solver->getSystemRhs().close();
624 
625  // once this is done, we recompute the single-physics residual
626  real64 const singlePhysicsNorm =
627  solver->calculateResidualNorm( time_n,
628  dt,
629  domain,
630  solver->getDofManager(),
631  solver->getSystemRhs().values() );
632  residualNorm += singlePhysicsNorm * singlePhysicsNorm;
633  } );
634 
635  // finally, we perform the convergence check on the multiphysics residual
636  residualNorm = sqrt( residualNorm );
637  GEOS_LOG_LEVEL_INFO_RANK_0( logInfo::Convergence, GEOS_FMT( " ( R ) = ( {:4.2e} )", residualNorm ) );
638  isConverged = ( residualNorm < params.m_newtonTol );
639 
640  }
642  {
643  // TODO also make recursive?
644  forEachArgInTuple( m_solvers, [&]( auto & solver, auto )
645  {
646  NonlinearSolverParameters const & singlePhysicsParams = solver->getNonlinearSolverParameters();
647  if( singlePhysicsParams.m_numNewtonIterations > singlePhysicsParams.m_minIterNewton )
648  {
649  isConverged = false;
650  }
651  } );
652  }
654  {
655  isConverged = checkSequentialSolutionIncrements( domain );
656  }
657  else
658  {
659  GEOS_ERROR( getDataContext() << ": Invalid sequential convergence criterion." );
660  }
661 
662  if( isConverged )
663  {
664  GEOS_LOG_LEVEL_INFO_RANK_0( logInfo::Convergence, GEOS_FMT( "***** The iterative coupling has converged in {} iteration(s) *****", iter + 1 ) );
665  }
666  }
667  return isConverged;
668  }
669 
670  virtual void
672  {
673  setSubSolvers();
674 
677  GEOS_THROW_IF( isSequential && usesLineSearch,
678  GEOS_FMT( "{}: line search is not supported by the coupled solver when {} is set to `{}`. Please set {} to `{}` to remove this error",
679  getNonlinearSolverParameters().getWrapperDataContext( NonlinearSolverParameters::viewKeysStruct::couplingTypeString() ),
680  NonlinearSolverParameters::viewKeysStruct::couplingTypeString(),
682  NonlinearSolverParameters::viewKeysStruct::lineSearchActionString(),
684  InputError );
685 
686  if( !isSequential )
687  {
689  }
690 
692  validateNonlinearAcceleration();
693  }
694 
695  virtual void validateNonlinearAcceleration()
696  {
697  GEOS_THROW ( GEOS_FMT( "{}: Nonlinear acceleration {} is not supported by {} solver '{}'",
698  getWrapperDataContext( NonlinearSolverParameters::viewKeysStruct::nonlinearAccelerationTypeString() ),
700  getCatalogName(), getName()),
701  InputError );
702  }
703 
704  virtual void
706  {
707  forEachArgInTuple( m_solvers, [&]( auto & solver, auto )
708  {
709  solver->getNonlinearSolverParameters() = getNonlinearSolverParameters();
710  } );
711  }
712 
713  virtual void startSequentialIteration( integer const & iter,
714  DomainPartition & domain )
715  {
716  GEOS_UNUSED_VAR( iter, domain );
717  }
718 
719  virtual void finishSequentialIteration( integer const & iter,
720  DomainPartition & domain )
721  {
722  GEOS_UNUSED_VAR( iter, domain );
723  }
724 
725 protected:
726 
728  std::tuple< SOLVERS *... > m_solvers;
729 
731  std::array< string, sizeof...( SOLVERS ) > m_names;
732 };
733 
734 } /* namespace geos */
735 
736 #endif /* GEOS_PHYSICSSOLVERS_MULTIPHYSICS_COUPLEDSOLVER_HPP_ */
#define GEOS_UNUSED_VAR(...)
Mark an unused variable and silence compiler warnings.
Definition: GeosxMacros.hpp:84
#define GEOS_THROW(msg, TYPE)
Throw an exception.
Definition: Logger.hpp:164
#define GEOS_ERROR(msg)
Raise a hard error and terminate the program.
Definition: Logger.hpp:157
#define GEOS_LOG_RANK_0(msg)
Log a message on screen on rank 0.
Definition: Logger.hpp:101
#define GEOS_THROW_IF(EXP, msg, TYPE)
Conditionally throw an exception.
Definition: Logger.hpp:151
#define GEOS_MARK_FUNCTION
Mark function with both Caliper and NVTX if enabled.
virtual void assembleCouplingTerms(real64 const time_n, real64 const dt, DomainPartition const &domain, DofManager const &dofManager, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs)
Utility function to compute coupling terms.
CoupledSolver & operator=(CoupledSolver const &)=delete
deleted assignment operator
virtual real64 fullyCoupledSolverStep(real64 const &time_n, real64 const &dt, int const cycleNumber, DomainPartition &domain)
Fully coupled solution approach solution step.
virtual bool resetConfigurationToDefault(DomainPartition &domain) const override
resets the configuration to the default value.
CoupledSolver(const string &name, Group *const parent)
main constructor for CoupledSolver Objects
virtual void outputConfigurationStatistics(DomainPartition const &domain) const override
CoupledSolver(CoupledSolver &&)=default
default move constructor
virtual void resetConfigurationToBeginningOfStep(DomainPartition &domain) override
resets the configuration to the beginning of the time-step.
CoupledSolver(CoupledSolver const &)=delete
deleted copy constructor
virtual void synchronizeNonlinearSolverParameters() override
syncronize the nonlinear solver parameters.
virtual void setupCoupling(DomainPartition const &domain, DofManager &dofManager) const
Utility function to set the coupling between degrees of freedom.
virtual real64 sequentiallyCoupledSolverStep(real64 const &time_n, real64 const &dt, int const cycleNumber, DomainPartition &domain)
Sequentially coupled solver step. It solves a nonlinear system of equations using a sequential approa...
virtual bool checkSequentialSolutionIncrements(DomainPartition &domain) const override
Check if the solution increments are ok to use.
virtual void postInputInitialization() override
virtual bool updateConfiguration(DomainPartition &domain) override
updates the configuration (if needed) based on the state after a converged Newton loop.
std::array< string, sizeof...(SOLVERS) > m_names
Names of the single-physics solvers.
std::tuple< SOLVERS *... > m_solvers
Pointers of the single-physics solvers.
CoupledSolver & operator=(CoupledSolver &&)=delete
deleted move operator
void setSubSolvers()
Utility function to set the subsolvers pointers using the names provided by the user.
virtual void mapSolutionBetweenSolvers(DomainPartition &domain, integer const solverType)
Maps the solution obtained from one solver to the fields used by the other solver(s)
The DoFManager is responsible for allocating global dofs, constructing sparsity patterns,...
Definition: DofManager.hpp:44
Partition of the decomposed physical domain. It also manages the connexion information to its neighbo...
SequentialConvergenceCriterion sequentialConvergenceCriterion() const
Getter for the sequential convergence criterion.
integer m_allowNonConverged
Flag to allow for a non-converged nonlinear solution and continue with the problem.
real64 m_newtonTol
The tolerance for the nonlinear convergence check.
NonlinearAccelerationType m_nonlinearAccelerationType
Type of nonlinear acceleration for sequential solver.
integer m_maxIterNewton
The maximum number of nonlinear iterations that are allowed.
real64 m_timeStepCutFactor
Factor by which the time step will be cut if a timestep cut is required.
integer m_numNewtonIterations
The number of nonlinear iterations that have been exectued.
integer m_numTimeStepAttempts
Number of times that the time-step had to be cut.
integer m_maxTimeStepCuts
Max number of time step cuts.
CouplingType couplingType() const
Getter for the coupling type.
LineSearchAction m_lineSearchAction
Flag to apply a line search.
@ ResidualNorm
convergence achieved when the residual drops below a given norm
@ NumberOfNonlinearIterations
convergence achieved when the subproblems convergence is achieved in less than minNewtonIteration
@ SolutionIncrements
convergence achieved when the solution increments are small enough
integer m_subcyclingOption
Flag to specify whether subcycling is allowed or not in sequential schemes.
Base class for all physics solvers.
virtual string getCatalogName() const =0
SolverStatistics m_solverStatistics
Solver statistics.
virtual void cleanup(real64 const time_n, integer const cycleNumber, integer const eventCounter, real64 const eventProgress, DomainPartition &domain) override
Called as the code exits the main run loop.
Timestamp getMeshModificationTimestamp(DomainPartition &domain) const
getter for the timestamp of the mesh modification on the mesh levels
NonlinearSolverParameters & getNonlinearSolverParameters()
accessor for the nonlinear solver parameters.
NonlinearSolverParameters m_nonlinearSolverParameters
Nonlinear solver parameters.
void logTimeStepCut()
Tell the solverStatistics that there is a time step cut.
void logNonlinearIteration(integer const numLinearIterations)
Tell the solverStatistics that we are doing a nonlinear iteration.
Wrapper< TBASE > & registerWrapper(string const &name, wrapperMap::KeyIndex::index_type *const rkey=nullptr)
Create and register a Wrapper around a new object.
DataContext const & getDataContext() const
Definition: Group.hpp:1343
string const & getName() const
Get group name.
Definition: Group.hpp:1329
Group & getParent()
Access the group's parent.
Definition: Group.hpp:1362
DataContext const & getWrapperDataContext(KEY key) const
Definition: Group.hpp:1354
#define GEOS_LOG_LEVEL_INFO_RANK_0(logInfoStruct, msg)
Output messages (only on rank 0) based on current Group's log level.
virtual void implicitStepSetup(real64 const &time_n, real64 const &dt, DomainPartition &domain) override
function to perform setup for implicit timestep
virtual void cleanup(real64 const time_n, integer const cycleNumber, integer const eventCounter, real64 const eventProgress, DomainPartition &domain) override
Called as the code exits the main run loop.
virtual real64 scalingForSystemSolution(DomainPartition &domain, DofManager const &dofManager, arrayView1d< real64 const > const &localSolution) override
Function to determine if the solution vector should be scaled back in order to maintain a known const...
virtual void updateState(DomainPartition &domain) override
Recompute all dependent quantities from primary variables (including constitutive models)
virtual void implicitStepComplete(real64 const &time_n, real64 const &dt, DomainPartition &domain) override
perform cleanup for implicit timestep
virtual void applyBoundaryConditions(real64 const time_n, real64 const dt, DomainPartition &domain, DofManager const &dofManager, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs) override
apply boundary condition to system
virtual real64 calculateResidualNorm(real64 const &time_n, real64 const &dt, DomainPartition const &domain, DofManager const &dofManager, arrayView1d< real64 const > const &localRhs) override
calculate the norm of the global system residual
real64 solverStep(real64 const &time_n, real64 const &dt, int const cycleNumber, DomainPartition &domain) override final
virtual real64 setNextDtBasedOnStateChange(real64 const &currentDt, DomainPartition &domain) override
function to set the next dt based on state change
virtual real64 solverStep(real64 const &time_n, real64 const &dt, integer const cycleNumber, DomainPartition &domain)
entry function to perform a solver step
void setupDofs(DomainPartition const &domain, DofManager &dofManager) const override
Populate degree-of-freedom manager with fields relevant to this solver.
virtual void applySystemSolution(DofManager const &dofManager, arrayView1d< real64 const > const &localSolution, real64 const scalingFactor, real64 const dt, DomainPartition &domain) override
Function to apply the solution vector to the state.
virtual real64 setNextDtBasedOnNewtonIter(real64 const &currentDt) override
function to set the next time step size based on Newton convergence
virtual void assembleSystem(real64 const time_n, real64 const dt, DomainPartition &domain, DofManager const &dofManager, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs) override
function to assemble the linear system matrix and rhs
virtual real64 setNextDtBasedOnNewtonIter(real64 const &currentDt)
function to set the next time step size based on Newton convergence
virtual bool checkSystemSolution(DomainPartition &domain, DofManager const &dofManager, arrayView1d< real64 const > const &localSolution, real64 const scalingFactor) override
Function to check system solution for physical consistency and constraint violation.
virtual void resetStateToBeginningOfStep(DomainPartition &domain) override
reset state of physics back to the beginning of the step.
@ FALSE
Not read from input.
@ REQUIRED
Required in input.
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
Definition: DataTypes.hpp:180
unsigned long long int Timestamp
Timestamp type (used to perform actions such a sparsity pattern computation after mesh modifications)
Definition: DataTypes.hpp:127
LvArray::CRSMatrixView< T, COL_INDEX, localIndex const, LvArray::ChaiBuffer > CRSMatrixView
Alias for CRS Matrix View.
Definition: DataTypes.hpp:310
std::string string
String type.
Definition: DataTypes.hpp:91
double real64
64-bit floating point type.
Definition: DataTypes.hpp:99
std::int32_t integer
Signed integer type.
Definition: DataTypes.hpp:82
Provides enum <-> string conversion facilities.
Exception class used to report errors in user input.
Definition: Logger.hpp:502
static constexpr char const * discretizationString()