Ice Sheet System Model  4.18
Code documentation
solutionsequence_shakti_nonlinear.cpp
Go to the documentation of this file.
1 
5 #include "./solutionsequences.h"
6 #include "../toolkits/toolkits.h"
7 #include "../classes/classes.h"
8 #include "../shared/shared.h"
9 #include "../modules/modules.h"
10 
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;
26 
27  /*Recover parameters: */
28  femmodel->parameters->FindParam(&max_nonlinear_iterations,StressbalanceMaxiterEnum);
33 
34  int count=0;
35  bool converged=false;
36 
37  /*Start non-linear iteration using input velocity: */
40 
41  for(;;){
42  /*save pointer to old solution*/
43  delete old_uf;old_uf=uf;
44  delete ug;
45 
46  SystemMatricesx(&Kff,&Kfs,&pf,&df,NULL,femmodel);
48  Reduceloadx(pf, Kfs, ys); delete Kfs;
50  Solverx(&uf, Kff, pf, old_uf, df, femmodel->parameters);
53 
54  convergence(&converged,Kff,pf,uf,old_uf,eps_res,eps_rel,eps_abs); delete Kff; delete pf; delete df;
56 
57  /*Increase count: */
58  count++;
59  if(count>=max_nonlinear_iterations){
60  _printf0_(" maximum number of nonlinear iterations (" << max_nonlinear_iterations << ") exceeded\n");
61  converged = true;
62  }
63  if(converged) break;
64  }
65 
66  if(VerboseConvergence()) _printf0_("\n total number of iterations: " << count << "\n");
67 
68  /*clean-up*/
69  delete uf;
70  delete ug;
71  delete old_uf;
72 }
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
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
solutionsequence_shakti_nonlinear
void solutionsequence_shakti_nonlinear(FemModel *femmodel)
Definition: solutionsequence_shakti_nonlinear.cpp:11
FemModel::UpdateConstraintsx
void UpdateConstraintsx(void)
Definition: FemModel.cpp:3027
FemModel::nodes
Nodes * nodes
Definition: FemModel.h:56
solutionsequences.h
StressbalanceAbstolEnum
@ StressbalanceAbstolEnum
Definition: EnumDefinitions.h:404
StressbalanceReltolEnum
@ StressbalanceReltolEnum
Definition: EnumDefinitions.h:410
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
FemModel
Definition: FemModel.h:31
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
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