File:NewtonHook.f90: Difference between revisions

From openpipeflow.org
Jump to navigation Jump to search
mNo edit summary
mNo edit summary
Line 1: Line 1:
{{latexPreamble}}
{{latexPreamble}}
'''See also Template/Example ''' [[File:Newton_Lorenz.f90]].
 
'''See also Template/Example ''' [[File:Newton_Lorenz.f90]](FORTRAN90), [[File:Newton_Lorenz_m.zip]](MATLAB)




Line 11: Line 12:
'''*** THE DESCRIPTION OF THE METHOD HAS MOVED TO [[Newton-Krylov_method]] ***'''
'''*** THE DESCRIPTION OF THE METHOD HAS MOVED TO [[Newton-Krylov_method]] ***'''


Follow the link just above for MATLAB and FORTRAN90 codes.


== The code ==
== The code ==

Revision as of 08:27, 4 October 2019

$ \renewcommand{\vec}[1]{ {\bf #1} } \newcommand{\bnabla}{ \vec{\nabla} } \newcommand{\Rey}{Re} \def\vechat#1{ \hat{ \vec{#1} } } \def\mat#1{#1} $

See also Template/Example File:Newton Lorenz.f90(FORTRAN90), File:Newton Lorenz m.zip(MATLAB)


The Jacobian-free Newton-Krylov (JFNK) Method

The code above implements the Jacobian-free Newton-Krylov (JFNK) method for solving \begin{equation}\label{eq:ROOTF}\vec{F}(\vec{x})=\vec{0},\end{equation} where $\vec{x}$ and $\vec{F}$ are $n$-vectors, supplemented with a Hookstep--Trust-region approach.

*** THE DESCRIPTION OF THE METHOD HAS MOVED TO Newton-Krylov_method ***


The code

To download, click the links at top of page. A template for calling the routine: File:Newton template.f90. The code also requires File:GMRESm.f90 and LAPACK.

  • Prior to calling the subroutine newtonhook(...), the arrays new_x(:)$\equiv\vec{x}_i$ and new_fx(:)$\equiv\vec{F}(\vec{x}_i)$ need to be allocated the dimension $n$.
  • Also prior to calling newtonhook(...), new_x(:) must be assigned the initial guess $\vec{x}_0$. (As the calculation proceeds, the vector new_x(:) is updated with improved vectors, having smaller corresponding new_fx(:).)
  • In addition to scalar and array variables, newtonhook(...) needs to be passed
    • a function that calculates the dot product of two vectors,
    • a function that calculates $\vec{F}(\vec{x})$ for a given $\vec{\delta x}$,
    • a function that calculates $\left.\frac{\vec{dF}}{\vec{dx}}\right|_{\vec{x}_i}\,\vec{\delta x}$ for a given $\vec{\delta x}$. See possible approximation above, where $\vec{F}(\vec{x}_i)$ is already stored in new_fx(:).
    • a subroutine that replaces a vector $\vec{x}$ with the solution of $M\vec{x}'=\vec{x}$ for $\vec{x}'$, where $M$ is a preconditioner matrix. This may simply be an empty subroutine if no preconditioner is required, i.e. $M=I$.
    • a subroutine that is called at the end of each iteration. This may be used to save the current state after each iteration, if desired.
  • The functions above may require auxiliary data, in addition to the given $\vec{x}$ or $\vec{\delta x}$. A few ways to access the data:
    • place the required data in a module and access via 'use mymodule'.
    • place the data in a common block.
    • write the vector to disk, call a script that runs an external programme, then load the result.
  • Extra constraints may be necessary, e.g. that determine an update to the period of an orbit. The function that evaluates $\left.\frac{\vec{dF}}{\vec{dx}}\right|_{\vec{x}_i}\,\vec{\delta x}$ for a given $\vec{\delta x}$ should append to the result the evaluations of the constraints. Correspondingly, for each constraint an extra 0 should be appended to $\vec{F}(\vec{x}_i)$ (the evaluated constraint should equal zero when converged).

Parallel use

It is NOT necessary to edit this code for parallel (MPI) use:

  • let each thread pass its subsection of the vector $\vec{x}$,
  • make the dot product function mpi_allreduce the result of the dot product.
  • to avoid multiple outputs to the terminal, set info=1 on rank 0 and info=0 for the other ranks.

File history

Click on a date/time to view the file as it appeared at that time.

Date/TimeDimensionsUserComment
current04:50, 26 June 2019 (6 KB)Apwillis (talk | contribs)Minor edits in comments only.
05:32, 13 December 2016 (6 KB)Apwillis (talk | contribs)== The method == Implements the Newton-Raphson method for solving \[\vec{F}(\vec{x})=\vec{0}\] with a Hookstep-Trust-region approach. At each iteration, the linear problem for the step in $\vec{x}$ is solved subject to the condition that the step be ...