-
Notifications
You must be signed in to change notification settings - Fork 3
/
newton.f90
638 lines (571 loc) · 21.6 KB
/
newton.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
MODULE newton
USE par_solve_mumps
USE qv_sp
USE qc_sp_M
USE qs_L_sp
USE dynamic_structures
USE sparse_matrix_profiles
USE sparse_matrix_operations
USE global_variables
USE miscellaneous_subroutines
USE prep_mesh_p1p2_sp
USE start_sparse_kit
USE Dirichlet_Neumann
USE axisym_boundary_values
USE case_dependent
!
USE loca_types
USE loca_pd
USE loca_bord
CONTAINS
!=======
FUNCTION nonlinear_solver_conwrap (x_vec, step_num, lambda, delta_s, con_ptr) &
RESULT(num_newt_its)
!
! Put the call to your nonlinear solver here.
! Input:
! x_vec solution vector
! con_ptr pointer to continuation structure, cast to (void *)
! must be passed to nonlinear solver and then passed
! to bordering algorithms.
! step_num Continuation step number
!
! Output:
! x_vec solution vector
!
! Return Value:
! num_newt_its Number of Newton iterations needed for
! convergence, used to pick next step size.
! Negative value means nonlinear solver didn't converge.
!
IMPLICIT NONE
! input variables
REAL(KIND=8), DIMENSION(0:Nx-1) :: x_vec ! FixMe check (or change) the indices of the vector
INTEGER :: step_num
REAL(KIND=8) :: lambda
REAL(KIND=8) :: delta_s
TYPE(con_struct), OPTIONAL :: con_ptr
! output variables
INTEGER :: num_newt_its
! common variables used
! Jacobian
! p_in
! velCmpnnts, np, np_L, Nx
! x0, u0, p0
! mm, jj, jj_L, js_D, zero_bvs_D
! DESINGULARIZE
! Re
! local variables
INTEGER :: n
REAL(KIND=8) :: residual
REAL(KIND=8), DIMENSION(:,:), ALLOCATABLE :: vv
REAL(KIND=8), DIMENSION(:), ALLOCATABLE :: dx, ww, rhs, dx0
INTEGER :: continuation_converged=0
REAL(KIND=8), DIMENSION(0:Nx-1) :: x_hook ! FixMe check (or change) the indices of the vector
REAL(KIND=8), DIMENSION(0:Nx-1) :: delta_x_hook ! FixMe check (or change) the indices of the vector
REAL(KIND=8) :: Reltol = 1d-3
REAL(KIND=8) :: Abstol = 1d-8
! (Reltol=1.0e-3, Abstol=1.0e-8 are good defaults.)
! executable statements
WRITE(*,*)
WRITE(*,*) '+++++++++++++++++++++++++++++++++++++'
WRITE(*,*) '--> CALL to nonlinear_solver_conwrap'
WRITE(*,*)
ALLOCATE ( vv (velCmpnnts, np) )
ALLOCATE ( ww (np_L) )
ALLOCATE ( dx (Nx) )
ALLOCATE ( rhs (Nx) )
ALLOCATE ( dx0 (Nx) )
dx0 = 1.29
x0 = x_vec
xx = x_vec
!************************************************************************************
!*** Newton't method in bi-incremental form is reported at the end of this module ***
!************************************************************************************
!=============================================================
! start of NEWTON'S ITERATIONS IN NON-INCREMENTAL FORM
!
WRITE(*,*) ' Start of Newton''s iterations'
WRITE(*,*) ' in NON-INCREMENTAL form'
!
DO n = 1, p_in%nwtn_maxite
WRITE(*,*)
WRITE(*,*) ' n = ', n
IF (n == p_in%nwtn_maxite) THEN
WRITE(*,*) ' ************************************************'
WRITE(*,*) ' *Maximum number of Newton''s iterations reached:'
WRITE(*,*) ' *ite_max = ', p_in%nwtn_maxite
WRITE(*,*) ' ************************************************'
WRITE(*,*)
ENDIF
CALL extract (x0, u0)
! CALL extract (x0, u0, p0)
! CALL vtk_plot_P2 (rr, jj, jj_L, u0, p0, trim(p_in%plot_directory) // 'iteSol.vtk')
!------------------------------------------------------------------
! call case dependent subroutine
!
CALL case_newton_iteprocess(n, continuation_converged)
!------------------------------------------------------------------
!-------------GENERATION OF THE RIGHT-HAND SIDE--------------------
!
! NON-INCREMENTAL FORM
! rhs <--- (u0 \dot \nabla)u0 + f
! 0
vv = 0
! CALL extract (x0, u0)
CALL qv_0y01_sp (mm, jj, u0, vv)
! CALL qc_ty0_sp_s (ms_2, jjs, iis, c_2, vv) ! cumulative
! CALL qc_ny0_sp_s (ms_3, jjs, iis, -q_3, vv) ! cumulative
u0(1,:) = volumeForcing(1,1)
u0(2,:) = volumeForcing(2,1)
u0(3,:) = volumeForcing(3,1)
CALL qv_0y0_sp (mm, jj, u0, 1d0, vv)
u0(1,:) = volumeForcing(1,2)
u0(2,:) = volumeForcing(2,2)
u0(3,:) = volumeForcing(3,2)
CALL qv_0y0_dR_sp(mm, jj, u0, 1d0, vv)
ww = 0
CALL collect (vv, ww, dx) ! here dx is the RHS
!------------------------------------------------------------------
!-------------ENFORCING DIRICHLET BOUNDARY CONDITIONS ON THE RHS---
!
! NON-INCREMENTAL FORM
! non-homogeneous boundary conditions
!
CALL Dirichlet_c (np, js_Axis, js_D, bvs_D, dx)
IF (DESINGULARIZE) dx(Nx) = 0d0
!------------------------------------------------------------------
!-------------COMPUTE RESIDUAL-------------------------------------
! rhs <-- (du0 \dot \nabla)du0
! div{u0}
vv = 0
CALL extract (dx0, u0)
CALL qv_0y01_sp (mm, jj, u0, vv)
ww = 0
CALL collect (vv, ww, rhs)
CALL Dirichlet_c (np, js_Axis, js_D, zero_bvs_D, rhs)
IF (DESINGULARIZE) rhs(Nx) = 0d0
residual = MAXVAL(ABS(rhs))
WRITE(*,*) ' |res|_L-infty = ', residual
!===================================================================
!===================================================================
IF (residual < p_in%nwtn_tol .AND. continuation_converged == 1) THEN
x_vec = x0
WRITE (*,*)
WRITE (*,*) ' Residual on the converged solution:'
WRITE (*,*) ' |res|_L-infty = ', residual
WRITE (*,*) ' End of Newton''s iterations.'
WRITE (*,*)
EXIT
ENDIF
!===================================================================
!===================================================================
!------------------------------------------------------------------
!-------------GENERATION OF THE JACOBIAN MATRIX--------------------
!-------------OF THE COUPLED EQUATION SYSTEM-----------------------
! Jacobian <--- [(u0.V)_ + (_.V)u0)] + K_ + V_ (ibp)
!
#if DEBUG > 1
WRITE(*,*) '*check*'
WRITE(*,*) ' Re = ', Re
#endif
CALL extract (x0, u0)
CALL ComputeJacobianMatrix (np, mm, jj, jj_L, js_Axis, js_D, DESINGULARIZE, Jacobian, Re, u0)
CALL par_mumps_master (NUMER_FACTOR, 1, Jacobian, 0)
!------------------------------------------------------------------
!-------------DIRECT SOLUTION OF THE COUPLED SYSTEM----------------
CALL par_mumps_master (DIRECT_SOLUTION, 1, Jacobian, 0, dx)
! as Newton's method is implemented in NON-INCREMENTAL FORM,
! we actually compute x_n even though we call it dx,
! dx has then to be computed as dx = x_n - x_n-1
!
dx = dx - x0
!------------------------------------------------------------------
!-------------LOCA'S STUFF-----------------------------------------
IF ( PRESENT(con_ptr) ) THEN
WRITE(*,*)
WRITE(*,*) ' --> CALL to continuation_hook'
#if DEBUG > 2
WRITE(*,*) ' |x0|_L-infty = ', MAXVAL(ABS(x0))
WRITE(*,*) ' |dx|_L-infty = ', MAXVAL(ABS(dx))
#endif
WRITE(*,*)
!-------------------------------------
! WARNING: continuation_hook expects the solution of
! J(-x) = + R and not
! J( x) = - R
!-------------------------------------
x_hook = x0
delta_x_hook = - dx
continuation_converged = continuation_hook(x_hook, delta_x_hook, &
con_ptr, Reltol, Abstol);
dx = - delta_x_hook
WRITE(*,*) ' done.'
WRITE(*,*)
ELSE
continuation_converged = 1
ENDIF
!------------------------------------------------------------------
!-------------UPDATE SOLUTION VECTOR-------------------------------
x_vec = x0 + dx
x0 = x_vec
dx0 = dx
ENDDO
!
! end of NEWTON'S ITERATIONS IN NON-INCREMENTAL FORM
!=============================================================
!------------------------------------------------------------------
!-------------UPDATE EVERYTHING------------------------------------
x0 = x_vec
xx = x_vec
pd%x => xx
CALL extract (x0, u0, p0)
CALL extract (xx, uu, pp)
IF ( n <= p_in%nwtn_maxite ) THEN
num_newt_its = n
ELSE
num_newt_its = -1
ENDIF
!------------------------------------------------------------------
! call case dependent subroutine
!
CALL case_newton_postprocess()
DEALLOCATE( vv, ww, dx, dx0, rhs )
END FUNCTION nonlinear_solver_conwrap
!==============================================================================
END MODULE newton
!------------------------------------------------------------------------------
! qui sotto e` riportato qualcosa che non deve essre perso
!! FUNCTION nonlinear_solver_conwrap (x_vec, con_ptr, step_num, lambda, delta_s) &
!! RESULT(num_newt_its) BIND(C, NAME='nonlinear_solver_conwrap')
!! !
!! ! Put the call to your nonlinear solver here.
!! ! Input:
!! ! x_vec solution vector
!! ! con_ptr pointer to continuation structure, cast to (void *)
!! ! must be passed to nonlinear solver and then passed
!! ! to bordering algorithms.
!! ! step_num Continuation step number
!! !
!! ! Output:
!! ! x_vec solution vector
!! !
!! ! Return Value:
!! ! num_newt_its Number of Newton iterations needed for
!! ! convergence, used to pick next step size.
!! ! Negative value means nonlinear solver didn't converge.
!! !
!!
!! USE ISO_C_BINDING
!!
!! IMPLICIT NONE
!!
!! INTERFACE
!! FUNCTION continuation_hook(x_hook, delta_x_hook, con_ptr, Reltol, Abstol) &
!! RESULT(continuation_converged) BIND(C, NAME='continuation_hook')
!! USE ISO_C_BINDING
!! REAL(KIND=C_DOUBLE) :: x_hook
!! REAL(KIND=C_DOUBLE) :: delta_x_hook
!! TYPE(C_PTR), VALUE :: con_ptr
!! REAL(KIND=C_DOUBLE), VALUE :: Reltol
!! REAL(KIND=C_DOUBLE), VALUE :: Abstol
!! INTEGER(KIND=C_INT) :: continuation_converged
!! END FUNCTION continuation_hook
!! END INTERFACE
!!
!! ! input variables
!! REAL(KIND=C_DOUBLE), DIMENSION(Nx) :: x_vec
!! TYPE(C_PTR), VALUE :: con_ptr
!! INTEGER(KIND=C_INT), VALUE :: step_num
!! REAL(KIND=C_DOUBLE), VALUE :: lambda
!! REAL(KIND=C_DOUBLE), VALUE :: delta_s
!! ! output variables
!! INTEGER(KIND=C_INT) :: num_newt_its
!!
!! ! common variables used
!! ! Jacobian
!! ! p_in
!! ! velCmpnnts, np, np_L, Nx
!! ! x0, u0, p0
!! ! mm, jj, jj_L, js_D, zero_bvs_D
!! ! DESINGULARIZE
!! ! Re
!!
!! ! local variables
!! INTEGER :: n
!! REAL(KIND=8) :: residual
!! REAL(KIND=8), DIMENSION(:,:), ALLOCATABLE :: vv, du0
!! REAL(KIND=8), DIMENSION(:), ALLOCATABLE :: dx, ww, dx0
!!
!!
!! INTEGER(KIND=C_INT) :: continuation_converged=0
!! REAL(KIND=C_DOUBLE), DIMENSION(Nx) :: x_hook
!! REAL(KIND=C_DOUBLE), DIMENSION(Nx) :: delta_x_hook
!! REAL(KIND=C_DOUBLE) :: Reltol = 1d-3
!! REAL(KIND=C_DOUBLE) :: Abstol = 1d-8
!! ! (Reltol=1.0e-3, Abstol=1.0e-8 are good defaults.)
!!
!!
!! ! executable statements
!!
!! WRITE(*,*)
!! WRITE(*,*) '+++++++++++++++++++++++++++++++++++++'
!! WRITE(*,*) '--> CALL to nonlinear_solver_conwrap'
!! WRITE(*,*)
!!
!!
!! ALLOCATE ( vv (velCmpnnts, np) )
!! ALLOCATE ( ww (np_L) )
!! ALLOCATE ( dx (Nx) )
!! ALLOCATE ( dx0(Nx) )
!! ALLOCATE ( du0(velCmpnnts, np) )
!!
!!
!! x0 = x_vec
!! xx = x_vec
!!
!!
!! !============================================
!! ! start of FIRST STEP IN NON-INCREMENTAL FORM
!! !
!! WRITE(*,*) ' First step in NON-INCREMENTAL form'
!! !
!! !------------------------------------------------------------------
!! !-------------GENERATION OF THE RIGHT-HAND SIDE--------------------
!! !
!! ! NON-INCREMENTAL FORM
!! ! rhs <--- (u0 \dot \nabla)u0
!! ! 0
!! vv = 0
!! ! CALL extract (x0, u0, p0)
!! CALL extract (x0, u0)
!! CALL qv_0y01_sp (mm, jj, u0, vv)
!! ! CALL qc_ty0_sp_s (ms_2, jjs, iis, c_2, vv) ! cumulative
!! ! CALL qc_ny0_sp_s (ms_3, jjs, iis, -q_3, vv) ! cumulative
!!
!! u0(1,:) = volumeForcing(1,1)
!! u0(2,:) = volumeForcing(2,1)
!! u0(3,:) = volumeForcing(3,1)
!! CALL qv_0y0_sp (mm, jj, u0, 1d0, vv)
!!
!! u0(1,:) = volumeForcing(1,2)
!! u0(2,:) = volumeForcing(2,2)
!! u0(3,:) = volumeForcing(3,2)
!! CALL qv_0y0_dR_sp(mm, jj, u0, 1d0, vv)
!!
!! ww = 0
!! CALL collect (vv, ww, dx) ! here dx is the RHS
!! !------------------------------------------------------------------
!! !-------------ENFORCING DIRICHLET BOUNDARY CONDITIONS ON THE RHS---
!! !
!! ! NON-INCREMENTAL FORM
!! ! non-homogeneous boundary conditions
!! !
!! CALL Dirichlet_c (np, js_Axis, js_D, bvs_D, dx)
!! IF (DESINGULARIZE) dx(Nx) = 0d0
!! !------------------------------------------------------------------
!! !-------------GENERATION OF THE JACOBIAN MATRIX--------------------
!! !-------------OF THE COUPLED EQUATION SYSTEM-----------------------
!! ! Jacobian <--- [(u0.V)_ + (_.V)u0)] + K_ + V_ (ibp)
!! !
!! !write(*,*) '*check*'
!! !write(*,*) ' Re = ', Re
!! CALL extract (x0, u0)
!! CALL ComputeJacobianMatrix (np, mm, jj, jj_L, js_Axis, js_D, DESINGULARIZE, Jacobian, Re, u0)
!! CALL par_mumps_master (NUMER_FACTOR, 1, Jacobian, 0)
!! !------------------------------------------------------------------
!! !-------------DIRECT SOLUTION OF THE COUPLED SYSTEM----------------
!! CALL par_mumps_master (DIRECT_SOLUTION, 1, Jacobian, 0, dx)
!! ! as the first step is done in NON-INCREMENTAL FORM,
!! ! we actually compute x_1 even though we call it dx,
!! ! dx has then to be computed as dx = x_1 - x_0
!! !
!! dx = dx - x0
!! !------------------------------------------------------------------
!! !-------------UPDATE SOLUTION VECTOR-------------------------------
!! x_vec = x0 + dx
!! x0 = x_vec
!! dx0 = dx
!! !
!! ! end of FIRST STEP IN NON-INCREMENTAL FORM
!! !============================================
!!
!!
!! WRITE(*,*)
!!
!!
!! !=============================================================
!! ! start of NEWTON'S ITERATIONS IN BI-INCREMENTAL FORM
!! !
!! WRITE(*,*) ' Start of Newton''s iterations'
!! WRITE(*,*) ' in BI-INCREMENTAL form'
!! !
!! DO n = 1, p_in%nwtn_maxite
!!
!! WRITE(*,*)
!! WRITE(*,*) ' n = ', n
!!
!! IF (n == p_in%nwtn_maxite) THEN
!!
!! WRITE(*,*) ' ************************************************'
!! WRITE(*,*) ' *Maximum number of Newton''s iterations reached:'
!! WRITE(*,*) ' *ite_max = ', p_in%nwtn_maxite
!! WRITE(*,*) ' ************************************************'
!! WRITE(*,*)
!!
!! ENDIF
!!
!! CALL extract (x0, u0)
!!
!! ! CALL extract (x0, u0, p0)
!! ! CALL vtk_plot_P2 (rr, jj, jj_L, u0, p0, trim(p_in%plot_directory) // 'iteSol.vtk')
!!
!! !------------------------------------------------------------------
!! ! call case dependent subroutine
!! !
!! CALL case_newton_iteprocess(n, continuation_converged)
!!
!! !------------------------------------------------------------------
!! !-------------GENERATION OF THE RIGHT-HAND SIDE--------------------
!! !
!! ! BI-INCREMENTAL FORM
!! ! rhs <--- - (du0 \dot \nabla)du0
!! ! 0
!! !
!! vv = 0
!! CALL extract (dx0, du0)
!! CALL qv_0y01_sp (mm, jj, du0, vv)
!!
!! ww = 0
!!
!! CALL collect (-vv, ww, dx) ! here dx is the RHS
!!
!! !------------------------------------------------------------------
!! !-------------ENFORCING DIRICHLET BOUNDARY CONDITIONS ON THE RHS---
!! !
!! ! BI-INCREMENTAL FORM
!! ! differential type boundary conditions
!! ! (homogeneous if not calling LOCA)
!! !
!! CALL extract_Dirichlet_c (np, js_Axis, js_D, x0, old_bvs_D)
!! CALL Dirichlet_c_DIFF (np, js_Axis, js_D, bvs_D, old_bvs_D, dx)
!!
!! IF (DESINGULARIZE) dx(Nx) = 0d0
!!
!! !write(*,*) '*check*'
!! !do ii = 1, SIZE(du0, 1)
!! ! write(*,*) 'MAXdelta_bvs_D = ', MAXVAL(bvs_D(ii)%DRL-old_bvs_D(ii)%DRL)
!! ! write(*,*) 'MINdelta_bvs_D = ', MINVAL(bvs_D(ii)%DRL-old_bvs_D(ii)%DRL)
!! !enddo
!!
!! !------------------------------------------------------------------
!! !-------------COMPUTE RESIDUAL-------------------------------------
!!
!! residual = MAXVAL(ABS(dx))
!!
!! WRITE(*,*) ' |res|_L-infty = ', residual
!!
!! !===================================================================
!! !===================================================================
!! IF (residual < p_in%nwtn_tol .AND. continuation_converged == 1) THEN
!! x_vec = x0
!! WRITE (*,*)
!! WRITE (*,*) ' Residual on the converged solution:'
!! WRITE (*,*) ' |res|_L-infty = ', residual
!! WRITE (*,*) ' End of Newton''s iterations.'
!! WRITE (*,*)
!! EXIT
!! ENDIF
!! !===================================================================
!! !===================================================================
!!
!! !------------------------------------------------------------------
!! !-------------GENERATION OF THE JACOBIAN MATRIX--------------------
!! !-------------OF THE COUPLED EQUATION SYSTEM-----------------------
!! ! Jacobian <--- [(u0.V)_ + (_.V)u0)] + K_ + V_ (ibp)
!! !
!!
!! !write(*,*) '*check*'
!! !write(*,*) ' Re = ', Re
!! CALL extract (x0, u0)
!! CALL ComputeJacobianMatrix (np, mm, jj, jj_L, js_Axis, js_D, DESINGULARIZE, Jacobian, Re, u0)
!! CALL par_mumps_master (NUMER_FACTOR, 1, Jacobian, 0)
!!
!! !------------------------------------------------------------------
!! !-------------DIRECT SOLUTION OF THE COUPLED SYSTEM----------------
!!
!! CALL par_mumps_master (DIRECT_SOLUTION, 1, Jacobian, 0, dx)
!!
!! !------------------------------------------------------------------
!! !-------------LOCA'S STUFF-----------------------------------------
!!
!! IF ( C_ASSOCIATED(con_ptr) ) THEN
!!
!! WRITE(*,*)
!! WRITE(*,*) ' --> CALL to continuation_hook'
!! !WRITE(*,*) ' |x0|_L-infty = ', MAXVAL(ABS(x0))
!! !WRITE(*,*) ' |dx|_L-infty = ', MAXVAL(ABS(dx))
!! WRITE(*,*)
!! !-------------------------------------
!! ! WARNING: continuation_hook expects the solution of
!! ! J(-x) = + R and not
!! ! J( x) = - R
!! !-------------------------------------
!!
!! x_hook = x0
!! delta_x_hook = - dx
!!
!! continuation_converged = continuation_hook(x_hook(1), delta_x_hook(1), &
!! con_ptr, Reltol, Abstol);
!! dx = - delta_x_hook
!!
!! WRITE(*,*) ' done.'
!! WRITE(*,*)
!!
!! ELSE
!! continuation_converged = 1
!! ENDIF
!!
!! !------------------------------------------------------------------
!! !-------------UPDATE SOLUTION VECTOR-------------------------------
!!
!! x_vec = x0 + dx
!!
!! x0 = x_vec
!! dx0 = dx
!!
!! ENDDO
!! !
!! ! end of NEWTON'S ITERATIONS IN BI-INCREMENTAL FORM
!! !=============================================================
!!
!!
!! !------------------------------------------------------------------
!! !-------------UPDATE EVERYTHING------------------------------------
!!
!! x0 = x_vec
!! xx = x_vec
!! pd%x = C_LOC(xx)
!! CALL extract (x0, u0, p0)
!! CALL extract (xx, uu, pp)
!!
!!
!! IF ( n <= p_in%nwtn_maxite ) THEN
!! num_newt_its = n
!! ELSE
!! num_newt_its = -1
!! ENDIF
!!
!! !------------------------------------------------------------------
!! ! call case dependent subroutine
!! !
!! CALL case_newton_postprocess()
!!
!!
!!
!! DEALLOCATE( vv, ww, dx, dx0, du0 )
!!
!!
!! END FUNCTION nonlinear_solver_conwrap
!!
!! !------------------------------------------------------------------------------