Ice Sheet System Model  4.18
Code documentation
Functions
solutionsequence_glads_nonlinear.cpp File Reference
#include "./solutionsequences.h"
#include "../toolkits/toolkits.h"
#include "../classes/classes.h"
#include "../shared/shared.h"
#include "../modules/modules.h"

Go to the source code of this file.

Functions

void solutionsequence_glads_nonlinear (FemModel *femmodel)
 

Function Documentation

◆ solutionsequence_glads_nonlinear()

void solutionsequence_glads_nonlinear ( FemModel femmodel)

Definition at line 11 of file solutionsequence_glads_nonlinear.cpp.

11  {
12 
13  /*intermediary: */
14  Matrix<IssmDouble>* Kff = NULL;
15  Matrix<IssmDouble>* Kfs = NULL;
16  Vector<IssmDouble>* ug = NULL;
17  Vector<IssmDouble>* uf = NULL;
18  Vector<IssmDouble>* old_uf = NULL;
19  Vector<IssmDouble>* pf = NULL;
20  Vector<IssmDouble>* df = NULL;
21  Vector<IssmDouble>* ys = NULL;
22 
23  /*parameters:*/
24  int max_nonlinear_iterations;
25  IssmDouble eps_res,eps_rel,eps_abs;
27 
28  /*Recover parameters (FIXME: from Stress balance for now :( )*/
29  femmodel->parameters->FindParam(&max_nonlinear_iterations,StressbalanceMaxiterEnum);
34 
35  int count_out=0;
36  bool converged_out=false;
37  int count_in=0;
38  bool converged_in=false;
39 
40  /*Start non-linear iteration using input velocity: */
43 
44  while(!converged_out){
45 
46  count_in=0;
47  converged_in=false;
48 
49  while(!converged_in){
50  /*save pointer to old solution*/
51  delete old_uf;old_uf=uf;
52  delete ug;
53 
54  SystemMatricesx(&Kff,&Kfs,&pf,&df,NULL,femmodel);
56  Reduceloadx(pf, Kfs, ys); delete Kfs;
58  Solverx(&uf, Kff, pf, old_uf, df, femmodel->parameters);
61 
62  convergence(&converged_in,Kff,pf,uf,old_uf,eps_res,eps_rel,eps_abs); delete Kff; delete pf; delete df;
64 
65  /*Increase count: */
66  count_in++;
67  if(count_in>=max_nonlinear_iterations && !converged_in){
68  _printf0_(" maximum number of nonlinear iterations of inner loop (" << max_nonlinear_iterations << ") exceeded\n");
69  converged_in = true;
70  }
71  }
72  if(VerboseConvergence()) _printf0_(setw(50) << left << " Inner loop converged in "<<count_in<<" iterations\n");
73 
74  if(VerboseConvergence()) _printf0_(" updating sheet thickness and channels cross section\n");
75  analysis->UpdateSheetThickness(femmodel);
77 
78  /*Converged if inner loop converged in one solution*/
79  if(count_in==1) converged_out = true;
80 
81  /*Increase count: */
82  count_out++;
83  if(count_out>=max_nonlinear_iterations && !converged_out){
84  _printf0_(" maximum number of nonlinear iterations of outer loop (" << max_nonlinear_iterations << ") exceeded\n");
85  converged_out = true;
86  }
87  }
88 
89  if(VerboseConvergence()) _printf0_("\n total number of iterations: " << count_out<<"x"<<count_in<<"\n");
90 
91  /*clean-up*/
92  delete uf;
93  delete ug;
94  delete old_uf;
95  delete analysis;
96 }
Matrix< IssmDouble >
SOLVER
#define SOLVER
Definition: Profiler.h:16
Mergesolutionfromftogx
void Mergesolutionfromftogx(Vector< IssmDouble > **pug, Vector< IssmDouble > *uf, Vector< IssmDouble > *ys, Nodes *nodes, Parameters *parameters, bool flag_ys0)
Definition: Mergesolutionfromftogx.cpp:8
IssmDouble
double IssmDouble
Definition: types.h:37
_printf0_
#define _printf0_(StreamArgs)
Definition: Print.h:29
VerboseConvergence
bool VerboseConvergence(void)
Definition: Verbosity.cpp:26
FemModel::parameters
Parameters * parameters
Definition: FemModel.h:46
HydrologyGlaDSAnalysis::UpdateSheetThickness
void UpdateSheetThickness(FemModel *femmodel)
Definition: HydrologyGlaDSAnalysis.cpp:458
Solverx
void Solverx(Vector< IssmDouble > **puf, Matrix< IssmDouble > *Kff, Vector< IssmDouble > *pf, Vector< IssmDouble > *uf0, Vector< IssmDouble > *df, Parameters *parameters)
Definition: Solverx.cpp:15
StressbalanceMaxiterEnum
@ StressbalanceMaxiterEnum
Definition: EnumDefinitions.h:407
FemModel::UpdateConstraintsx
void UpdateConstraintsx(void)
Definition: FemModel.cpp:3027
FemModel::nodes
Nodes * nodes
Definition: FemModel.h:56
StressbalanceAbstolEnum
@ StressbalanceAbstolEnum
Definition: EnumDefinitions.h:404
StressbalanceReltolEnum
@ StressbalanceReltolEnum
Definition: EnumDefinitions.h:410
HydrologyGlaDSAnalysis::UpdateChannelCrossSection
void UpdateChannelCrossSection(FemModel *femmodel)
Definition: HydrologyGlaDSAnalysis.cpp:607
Profiler::Stop
void Stop(int tagenum, bool dontmpisync=true)
Definition: Profiler.cpp:179
Reducevectorgtofx
void Reducevectorgtofx(Vector< IssmDouble > **puf, Vector< IssmDouble > *ug, Nodes *nodes, Parameters *parameters)
Definition: Reducevectorgtofx.cpp:8
Reduceloadx
void Reduceloadx(Vector< IssmDouble > *pf, Matrix< IssmDouble > *Kfs, Vector< IssmDouble > *y_s, bool flag_ys0)
Definition: Reduceloadx.cpp:14
Profiler::Start
void Start(int tagenum, bool dontmpisync=true)
Definition: Profiler.cpp:139
InputUpdateFromSolutionx
void InputUpdateFromSolutionx(FemModel *femmodel, Vector< IssmDouble > *solution)
Definition: InputUpdateFromSolutionx.cpp:9
StressbalanceRestolEnum
@ StressbalanceRestolEnum
Definition: EnumDefinitions.h:412
GetSolutionFromInputsx
void GetSolutionFromInputsx(Vector< IssmDouble > **psolution, FemModel *femmodel)
Definition: GetSolutionFromInputsx.cpp:9
HydrologyGlaDSAnalysis
Definition: HydrologyGlaDSAnalysis.h:11
Parameters::FindParam
void FindParam(bool *pinteger, int enum_type)
Definition: Parameters.cpp:262
CreateNodalConstraintsx
void CreateNodalConstraintsx(Vector< IssmDouble > **pys, Nodes *nodes)
Definition: CreateNodalConstraintsx.cpp:10
SystemMatricesx
void SystemMatricesx(Matrix< IssmDouble > **pKff, Matrix< IssmDouble > **pKfs, Vector< IssmDouble > **ppf, Vector< IssmDouble > **pdf, IssmDouble *pkmax, FemModel *femmodel, bool isAllocated)
Definition: SystemMatricesx.cpp:10
convergence
void convergence(bool *pconverged, Matrix< IssmDouble > *Kff, Vector< IssmDouble > *pf, Vector< IssmDouble > *uf, Vector< IssmDouble > *old_uf, IssmDouble eps_res, IssmDouble eps_rel, IssmDouble eps_abs)
Definition: convergence.cpp:9
FemModel::profiler
Profiler * profiler
Definition: FemModel.h:42
Vector< IssmDouble >
femmodel
FemModel * femmodel
Definition: esmfbinders.cpp:16