Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027 #ifndef TSF_NVECTOR_HPP
00028 #define TSF_NVECTOR_HPP
00029
00030 #include "SundanceDefs.hpp"
00031 #include "TSFVectorDecl.hpp"
00032 #include "TSFVectorSpaceDecl.hpp"
00033
00034 #ifdef HAVE_SUNDIALS
00035 #include "nvector.h"
00036 #include "sundialstypes.h"
00037
00038
00039
00040
00041
00042 class _TrilinosVectorContent
00043 {
00044 public:
00045 _TrilinosVectorContent() : data() {;}
00046 TSFExtended::Vector<realtype> data;
00047 };
00048
00049 typedef struct _TrilinosVectorContent* TrilinosVectorContent;
00050
00051
00052 extern "C"
00053 {
00054
00055 #define NV_CONTENT_Trilinos(v) ( (TrilinosVectorContent)(v->content) );
00056
00057
00058
00059
00060 inline Vector<realtype>& toTrilinos(N_Vector x)
00061 {
00062 TrilinosVectorContent xPtr = reinterpret_cast<TrilinosVectorContent>(x->content);
00063
00064 return xPtr->data;
00065 }
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075 N_Vector N_VNew_Trilinos(const VectorSpace<realtype>& space);
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085 void N_VPrint_Trilinos(N_Vector v);
00086
00087
00088
00089
00090
00091
00092
00093 N_Vector N_VCloneEmpty_Trilinos(N_Vector w);
00094 N_Vector N_VClone_Trilinos(N_Vector w);
00095 void N_VDestroy_Trilinos(N_Vector v);
00096 void N_VSpace_Trilinos(N_Vector v, long int *lrw, long int *liw);
00097 realtype *N_VGetArrayPointer_Trilinos(N_Vector v);
00098 void N_VSetArrayPointer_Trilinos(realtype *v_data, N_Vector v);
00099 void N_VLinearSum_Trilinos(realtype a, N_Vector x, realtype b,
00100 N_Vector y, N_Vector z);
00101 void N_VConst_Trilinos(realtype c, N_Vector z);
00102 void N_VProd_Trilinos(N_Vector x, N_Vector y, N_Vector z);
00103 void N_VDiv_Trilinos(N_Vector x, N_Vector y, N_Vector z);
00104 void N_VScale_Trilinos(realtype c, N_Vector x, N_Vector z);
00105 void N_VAbs_Trilinos(N_Vector x, N_Vector z);
00106 void N_VInv_Trilinos(N_Vector x, N_Vector z);
00107 void N_VAddConst_Trilinos(N_Vector x, realtype b, N_Vector z);
00108 realtype N_VDotProd_Trilinos(N_Vector x, N_Vector y);
00109 realtype N_VMaxNorm_Trilinos(N_Vector x);
00110 realtype N_VWrmsNorm_Trilinos(N_Vector x, N_Vector w);
00111 realtype N_VWrmsNormMask_Trilinos(N_Vector x, N_Vector w, N_Vector id);
00112 realtype N_VMin_Trilinos(N_Vector x);
00113 realtype N_VWL2Norm_Trilinos(N_Vector x, N_Vector w);
00114 realtype N_VL1Norm_Trilinos(N_Vector x);
00115 void N_VCompare_Trilinos(realtype c, N_Vector x, N_Vector z);
00116 booleantype N_VInvTest_Trilinos(N_Vector x, N_Vector z);
00117 booleantype N_VConstrMask_Trilinos(N_Vector c, N_Vector x, N_Vector m);
00118 realtype N_VMinQuotient_Trilinos(N_Vector num, N_Vector denom);
00119 }
00120
00121 #endif // HAVE_SUNDIALS
00122
00123 #endif