| 1 | function solutionsequence_nonlinear(femmodel::FemModel,analysis::Analysis,maxiter::Int64,restol::Float64,reltol::Float64,abstol::Float64) # {{{
|
|---|
| 2 |
|
|---|
| 3 | #Initialize number of iterations
|
|---|
| 4 | count = 0
|
|---|
| 5 | converged = false
|
|---|
| 6 |
|
|---|
| 7 | #Get existing solution
|
|---|
| 8 | ug = GetSolutionFromInputsx(analysis,femmodel)
|
|---|
| 9 | uf = Reducevectorgtofx(ug,femmodel.nodes)
|
|---|
| 10 |
|
|---|
| 11 | #Update once again the solution to make sure that vx and vxold are similar (for next step in transient or steadystate)
|
|---|
| 12 | InputUpdateFromSolutionx(analysis,ug,femmodel)
|
|---|
| 13 |
|
|---|
| 14 | #Loop until we reach convergence
|
|---|
| 15 | while(~converged)
|
|---|
| 16 |
|
|---|
| 17 | #Get new matrices
|
|---|
| 18 | Kff, Kfs, pf = SystemMatricesx(femmodel,analysis)
|
|---|
| 19 |
|
|---|
| 20 | #Enforce constraints
|
|---|
| 21 | ys = CreateNodalConstraintsx(femmodel.nodes)
|
|---|
| 22 | Reduceloadx!(pf, Kfs, ys)
|
|---|
| 23 |
|
|---|
| 24 | #Solve!
|
|---|
| 25 | old_uf = uf
|
|---|
| 26 | uf = Solverx(Kff, pf, old_uf)
|
|---|
| 27 |
|
|---|
| 28 | #Merge uf with ys
|
|---|
| 29 | Mergesolutionfromftogx(ug, uf, ys, femmodel.nodes)
|
|---|
| 30 |
|
|---|
| 31 | #Check for convergence
|
|---|
| 32 | converged = convergence(Kff,pf,uf,old_uf,restol,reltol,abstol)
|
|---|
| 33 | InputUpdateFromSolutionx(analysis,ug,femmodel)
|
|---|
| 34 |
|
|---|
| 35 | #Increase count
|
|---|
| 36 | count += 1
|
|---|
| 37 | if(count>=maxiter)
|
|---|
| 38 | println(" maximum number of nonlinear iterations (",maxiter,") exceeded")
|
|---|
| 39 | converged = true
|
|---|
| 40 | end
|
|---|
| 41 | end
|
|---|
| 42 |
|
|---|
| 43 | print("\n total number of iterations: ", count, "\n")
|
|---|
| 44 |
|
|---|
| 45 | end# }}}
|
|---|
| 46 | function convergence(Kff::IssmMatrix, pf::IssmVector, uf::IssmVector, old_uf::IssmVector, restol::Float64, reltol::Float64, abstol::Float64)#{{{
|
|---|
| 47 |
|
|---|
| 48 | print(" checking convergence\n");
|
|---|
| 49 |
|
|---|
| 50 | #If solution vector is empty, return true
|
|---|
| 51 | if(IsEmpty(uf))
|
|---|
| 52 | return true
|
|---|
| 53 | end
|
|---|
| 54 |
|
|---|
| 55 | #Convergence criterion #1: force equilibrium (Mandatory)
|
|---|
| 56 | #compute K[n]U[n-1] - F
|
|---|
| 57 | KUold = Duplicate(uf); MatMult!(Kff,old_uf,KUold)
|
|---|
| 58 | KUoldF = Duplicate(KUold); VecCopy!(KUold, KUoldF); AXPY!(KUoldF, -1.0, pf)
|
|---|
| 59 | nKUoldF = Norm(KUoldF,2)
|
|---|
| 60 | nF = Norm(pf,2)
|
|---|
| 61 | res = nKUoldF/nF
|
|---|
| 62 | if ~isfinite(res)
|
|---|
| 63 | println("norm nf = ", nF, " and norm kuold = ",nKUoldF)
|
|---|
| 64 | error("mechanical equilibrium convergence criterion is not finite!")
|
|---|
| 65 | end
|
|---|
| 66 | if(res<restol)
|
|---|
| 67 | print(" mechanical equilibrium convergence criterion ", res*100, " < ", restol*100, " %\n")
|
|---|
| 68 | converged=true
|
|---|
| 69 | else
|
|---|
| 70 | print(" mechanical equilibrium convergence criterion ", res*100, " > ", restol*100, " %\n")
|
|---|
| 71 | converged=false;
|
|---|
| 72 | end
|
|---|
| 73 |
|
|---|
| 74 | #Convergence criterion #2: norm(du)/norm(u)
|
|---|
| 75 | if ~isnan(reltol)
|
|---|
| 76 | duf = Duplicate(old_uf); VecCopy!(old_uf,duf); AXPY!(duf, -1.0, uf)
|
|---|
| 77 | ndu = Norm(duf, 2); nu = Norm(old_uf, 2)
|
|---|
| 78 | if ~isfinite(ndu) | ~isfinite(nu)
|
|---|
| 79 | error("convergence criterion is not finite!")
|
|---|
| 80 | end
|
|---|
| 81 | if((ndu/nu)<reltol)
|
|---|
| 82 | print(" Convergence criterion: norm(du)/norm(u) ", ndu/nu*100, " < ", reltol*100, " %\n")
|
|---|
| 83 | else
|
|---|
| 84 | print(" Convergence criterion: norm(du)/norm(u) ", ndu/nu*100, " > ", reltol*100, " %\n")
|
|---|
| 85 | converged=false;
|
|---|
| 86 | end
|
|---|
| 87 | end
|
|---|
| 88 |
|
|---|
| 89 | #Convergence criterion #3: max(du)
|
|---|
| 90 | if ~isnan(abstol)
|
|---|
| 91 | duf = Duplicate(old_uf); VecCopy!(old_uf,duf); AXPY!(duf, -1.0, uf)
|
|---|
| 92 | nduinf= Norm(duf, 3)
|
|---|
| 93 | if ~isfinite(nduinf)
|
|---|
| 94 | error("convergence criterion is not finite!")
|
|---|
| 95 | end
|
|---|
| 96 | if(nduinf<abstol)
|
|---|
| 97 | print(" Convergence criterion: max(du) ", nduinf, " < ", abstol, "\n")
|
|---|
| 98 | else
|
|---|
| 99 | print(" Convergence criterion: max(du) ", nduinf, " > ", abstol, "\n")
|
|---|
| 100 | converged=false;
|
|---|
| 101 | end
|
|---|
| 102 | end
|
|---|
| 103 |
|
|---|
| 104 | return converged
|
|---|
| 105 |
|
|---|
| 106 | end#}}}
|
|---|