1 /* -----------------------------------------------------------------
2  * Programmer(s): Carol Woodward @ LLNL
3  *                Daniel R. Reynolds @ SMU
4  *                David J. Gardner @ LLNL
5  * -----------------------------------------------------------------
6  * SUNDIALS Copyright Start
7  * Copyright (c) 2002-2021, Lawrence Livermore National Security
8  * and Southern Methodist University.
9  * All rights reserved.
10  *
11  * See the top-level LICENSE and NOTICE files for details.
12  *
13  * SPDX-License-Identifier: BSD-3-Clause
14  * SUNDIALS Copyright End
15  * -----------------------------------------------------------------*/
16 
17 #include <stdio.h>
18 #include <stdlib.h>
19 
20 #include "fkinsol.h"
21 #include "kinsol_impl.h"
22 
23 #include <kinsol/kinsol_ls.h>
24 #include <sunmatrix/sunmatrix_sparse.h>
25 
26 /*=============================================================*/
27 
28 /* Prototype of the Fortran routine */
29 
30 #ifdef __cplusplus  /* wrapper to enable C++ usage */
31 extern "C" {
32 #endif
33 
34 extern void FKIN_SPJAC(realtype *Y, realtype *FY, long int *N,
35                        long int *NNZ, realtype *JDATA,
36                        sunindextype *JRVALS, sunindextype *JCPTRS,
37                        realtype *V1, realtype *V2, int *ier);
38 
39 #ifdef __cplusplus
40 }
41 #endif
42 
43 /*=============================================================*/
44 
45 /* Fortran interface to C routine KINSlsSetSparseJacFn; see
46    fkinsol.h for further information */
FKIN_SPARSESETJAC(int * ier)47 void FKIN_SPARSESETJAC(int *ier)
48 {
49 #if defined(SUNDIALS_INT32_T)
50   KINProcessError((KINMem) KIN_kinmem, KIN_ILL_INPUT, "KIN",
51                   "FKINSPARSESETJAC",
52                   "Sparse Fortran users must configure SUNDIALS with 64-bit integers.");
53   *ier = 1;
54 #else
55   *ier = KINSetJacFn(KIN_kinmem, FKINSparseJac);
56 #endif
57 }
58 
59 /*=============================================================*/
60 
61 /* C interface to user-supplied Fortran routine FKINSPJAC; see
62    fkinsol.h for additional information  */
FKINSparseJac(N_Vector y,N_Vector fy,SUNMatrix J,void * user_data,N_Vector vtemp1,N_Vector vtemp2)63 int FKINSparseJac(N_Vector y, N_Vector fy, SUNMatrix J,
64                   void *user_data, N_Vector vtemp1,
65                   N_Vector vtemp2)
66 {
67   int ier;
68   realtype *ydata, *fydata, *v1data, *v2data, *Jdata;
69   long int NP, NNZ;
70   sunindextype *indexvals, *indexptrs;
71 
72   ydata   = N_VGetArrayPointer(y);
73   fydata  = N_VGetArrayPointer(fy);
74   v1data  = N_VGetArrayPointer(vtemp1);
75   v2data  = N_VGetArrayPointer(vtemp2);
76 
77   NP = SUNSparseMatrix_NP(J);
78   NNZ = SUNSparseMatrix_NNZ(J);
79   Jdata = SUNSparseMatrix_Data(J);
80   indexvals = SUNSparseMatrix_IndexValues(J);
81   indexptrs = SUNSparseMatrix_IndexPointers(J);
82 
83   FKIN_SPJAC(ydata, fydata, &NP, &NNZ,
84              Jdata, indexvals, indexptrs,
85              v1data, v2data, &ier);
86   return(ier);
87 }
88 
89