source: issm/trunk-jpl/src/c/solvers/solver_newton.cpp@ 14416

Last change on this file since 14416 was 14416, checked in by Mathieu Morlighem, 12 years ago

FIX: commited wrong file

File size: 4.3 KB
Line 
1/*!\file: solver_nonlinear.cpp
2 * \brief: core of a non-linear solution, using fixed-point method
3 */
4
5#include "../toolkits/toolkits.h"
6#include "../classes/objects/objects.h"
7#include "../io/io.h"
8#include "../EnumDefinitions/EnumDefinitions.h"
9#include "../modules/modules.h"
10#include "../solutions/solutions.h"
11#include "./solvers.h"
12
13void solver_newton(FemModel* femmodel){
14
15 /*intermediary: */
16 bool converged;
17 int count,newton;
18 IssmDouble kmax;
19 Matrix<IssmDouble>* Kff = NULL;
20 Matrix<IssmDouble>* Kfs = NULL;
21 Matrix<IssmDouble>* Jff = NULL;
22 Vector<IssmDouble>* ug = NULL;
23 Vector<IssmDouble>* old_ug = NULL;
24 Vector<IssmDouble>* uf = NULL;
25 Vector<IssmDouble>* old_uf = NULL;
26 Vector<IssmDouble>* duf = NULL;
27 Vector<IssmDouble>* pf = NULL;
28 Vector<IssmDouble>* pJf = NULL;
29 Vector<IssmDouble>* df = NULL;
30 Vector<IssmDouble>* ys = NULL;
31
32 /*parameters:*/
33 int max_nonlinear_iterations;
34 int configuration_type;
35
36 /*Recover parameters: */
37 femmodel->parameters->FindParam(&max_nonlinear_iterations,DiagnosticMaxiterEnum);
38 femmodel->parameters->FindParam(&configuration_type,ConfigurationTypeEnum);
39 femmodel->parameters->FindParam(&newton,DiagnosticIsnewtonEnum);
40 femmodel->UpdateConstraintsx();
41
42 count=1;
43 converged=false;
44
45 /*Start non-linear iteration using input velocity: */
46 GetSolutionFromInputsx(&ug,femmodel->elements,femmodel->nodes,femmodel->vertices,femmodel->loads,femmodel->materials,femmodel->parameters);
47 Reducevectorgtofx(&uf,ug,femmodel->nodes,femmodel->parameters);
48
49 //Update once again the solution to make sure that vx and vxold are similar (for next step in transient or steadystate)
50 InputUpdateFromConstantx(femmodel->elements,femmodel->nodes,femmodel->vertices,femmodel->loads,femmodel->materials,femmodel->parameters,converged,ConvergedEnum);
51 InputUpdateFromSolutionx(femmodel->elements,femmodel->nodes,femmodel->vertices,femmodel->loads,femmodel->materials,femmodel->parameters,ug);
52
53 for(;;){
54
55 xdelete(&old_ug);old_ug=ug;
56 xdelete(&old_uf);old_uf=uf;
57
58 /*Solver forward model*/
59 if(count==1 || newton==2){
60 femmodel->SystemMatricesx(&Kff, &Kfs, &pf, &df, NULL);
61 CreateNodalConstraintsx(&ys,femmodel->nodes,configuration_type);
62 Reduceloadx(pf,Kfs,ys);xdelete(&Kfs);
63 Solverx(&uf,Kff,pf,old_uf,df,femmodel->parameters);xdelete(&df);
64 Mergesolutionfromftogx(&ug,uf,ys,femmodel->nodes,femmodel->parameters);xdelete(&ys);
65 InputUpdateFromSolutionx(femmodel->elements,femmodel->nodes,femmodel->vertices,femmodel->loads,femmodel->materials,femmodel->parameters,ug);xdelete(&ug);
66 xdelete(&old_ug);old_ug=ug;
67 xdelete(&old_uf);old_uf=uf;
68 }
69 uf=old_uf->Duplicate(); old_uf->Copy(uf);
70
71 /*Prepare next iteration using Newton's method*/
72 femmodel->SystemMatricesx(&Kff, &Kfs, &pf, &df, NULL);xdelete(&df);
73 CreateNodalConstraintsx(&ys,femmodel->nodes,configuration_type);
74 Reduceloadx(pf,Kfs,ys);xdelete(&Kfs);
75
76 pJf=pf->Duplicate();
77 Kff->MatMult(uf,pJf);// xdelete(&Kff);
78 pJf->Scale(-1.0); pJf->AXPY(pf,+1.0); //xdelete(&pf);
79
80 femmodel->CreateJacobianMatrixx(&Jff,kmax);
81 Solverx(&duf,Jff,pJf,NULL,NULL,femmodel->parameters); xdelete(&Jff); xdelete(&pJf);
82 uf->AXPY(duf, 1.0); xdelete(&duf);
83 Mergesolutionfromftogx(&ug,uf,ys,femmodel->nodes,femmodel->parameters);xdelete(&ys);
84 InputUpdateFromSolutionx(femmodel->elements,femmodel->nodes,femmodel->vertices,femmodel->loads,femmodel->materials,femmodel->parameters,ug);
85
86 /*Check convergence*/
87 convergence(&converged,Kff,pf,uf,old_uf,femmodel->parameters);
88 xdelete(&Kff); xdelete(&pf);
89 if(converged==true){
90 bool max_iteration_state=false;
91 int tempStep=1;
92 IssmDouble tempTime=1.0;
93 femmodel->results->AddObject(new GenericExternalResult<bool>(femmodel->results->Size()+1, MaxIterationConvergenceFlagEnum, max_iteration_state, tempStep, tempTime));
94 break;
95 }
96 if(count>=max_nonlinear_iterations){
97 _pprintLine_(" maximum number of Newton iterations (" << max_nonlinear_iterations << ") exceeded");
98 bool max_iteration_state=true;
99 int tempStep=1;
100 IssmDouble tempTime=1.0;
101 femmodel->results->AddObject(new GenericExternalResult<bool>(femmodel->results->Size()+1, MaxIterationConvergenceFlagEnum, max_iteration_state, tempStep, tempTime));
102 break;
103 }
104
105 count++;
106 }
107
108 if(VerboseConvergence()) _pprintLine_("\n total number of iterations: " << count-1);
109
110 /*clean-up*/
111 xdelete(&uf);
112 xdelete(&ug);
113 xdelete(&old_ug);
114 xdelete(&old_uf);
115}
Note: See TracBrowser for help on using the repository browser.