/* nag_roots_sys_func_easy (c05qbc) Example Program.
*
* Copyright 2024 Numerical Algorithms Group.
*
* Mark 30.0, 2024.
*/
#include <math.h>
#include <nag.h>
#include <stdio.h>
#ifdef __cplusplus
extern "C" {
#endif
static void NAG_CALL fcn(Integer n, const double x[], double fvec[],
Nag_Comm *comm, Integer *iflag);
#ifdef __cplusplus
}
#endif
int main(void) {
static double ruser[1] = {-1.0};
Integer exit_status = 0, i, n = 9;
double *fvec = 0, *x = 0, xtol;
/* Nag Types */
NagError fail;
Nag_Comm comm;
INIT_FAIL(fail);
printf("nag_roots_sys_func_easy (c05qbc) Example Program Results\n");
/* For communication with user-supplied functions: */
comm.user = ruser;
if (n > 0) {
if (!(fvec = NAG_ALLOC(n, double)) || !(x = NAG_ALLOC(n, double))) {
printf("Allocation failure\n");
exit_status = -1;
goto END;
}
} else {
printf("Invalid n.\n");
exit_status = 1;
goto END;
}
/* The following starting values provide a rough solution. */
for (i = 0; i < n; i++)
x[i] = -1.0;
/* nag_machine_precision (x02ajc).
* The machine precision
*/
xtol = sqrt(nag_machine_precision);
/* nag_roots_sys_func_easy (c05qbc).
* Solution of a system of nonlinear equations (function
* values only)
*/
nag_roots_sys_func_easy(fcn, n, x, fvec, xtol, &comm, &fail);
if (fail.code != NE_NOERROR) {
printf("Error from nag_roots_sys_func_easy (c05qbc).\n%s\n", fail.message);
exit_status = 1;
if (fail.code != NE_TOO_MANY_FEVALS && fail.code != NE_TOO_SMALL &&
fail.code != NE_NO_IMPROVEMENT)
goto END;
}
printf(fail.code == NE_NOERROR ? "Final approximate" : "Approximate");
printf(" solution\n\n");
for (i = 0; i < n; i++)
printf("%12.4f%s", x[i], (i % 3 == 2 || i == n - 1) ? "\n" : " ");
if (fail.code != NE_NOERROR)
exit_status = 2;
END:
NAG_FREE(fvec);
NAG_FREE(x);
return exit_status;
}
static void NAG_CALL fcn(Integer n, const double x[], double fvec[],
Nag_Comm *comm, Integer *iflag) {
Integer k;
if (comm->user[0] == -1.0) {
printf("(User-supplied callback fcn, first invocation.)\n");
comm->user[0] = 0.0;
}
for (k = 0; k < n; ++k) {
fvec[k] = (3.0 - x[k] * 2.0) * x[k] + 1.0;
if (k > 0)
fvec[k] -= x[k - 1];
if (k < n - 1)
fvec[k] -= x[k + 1] * 2.0;
}
/* Set iflag negative to terminate execution for any reason. */
*iflag = 0;
}