1 /* 2 Objects which encapsulate discretizations+continuum residuals 3 */ 4 #if !defined(__PETSCDS_H) 5 #define __PETSCDS_H 6 #include <petscfe.h> 7 #include <petscfv.h> 8 #include <petscdstypes.h> 9 10 PETSC_EXTERN PetscErrorCode PetscDSInitializePackage(void); 11 12 PETSC_EXTERN PetscClassId PETSCDS_CLASSID; 13 14 /*J 15 PetscDSType - String with the name of a PETSc discrete system 16 17 Level: beginner 18 19 .seealso: PetscDSSetType(), PetscDS 20 J*/ 21 typedef const char *PetscDSType; 22 #define PETSCDSBASIC "basic" 23 24 typedef void (*PetscPointFunc)(PetscInt, PetscInt, PetscInt, 25 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 26 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 27 PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]); 28 typedef void (*PetscPointJac)(PetscInt, PetscInt, PetscInt, 29 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 30 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 31 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]); 32 typedef void (*PetscBdPointFunc)(PetscInt, PetscInt, PetscInt, 33 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 34 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 35 PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]); 36 typedef void (*PetscBdPointJac)(PetscInt, PetscInt, PetscInt, 37 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 38 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 39 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]); 40 typedef void (*PetscRiemannFunc)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *); 41 typedef PetscErrorCode (*PetscSimplePointFunc)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *); 42 43 PETSC_EXTERN PetscFunctionList PetscDSList; 44 PETSC_EXTERN PetscErrorCode PetscDSCreate(MPI_Comm, PetscDS *); 45 PETSC_EXTERN PetscErrorCode PetscDSDestroy(PetscDS *); 46 PETSC_EXTERN PetscErrorCode PetscDSSetType(PetscDS, PetscDSType); 47 PETSC_EXTERN PetscErrorCode PetscDSGetType(PetscDS, PetscDSType *); 48 PETSC_EXTERN PetscErrorCode PetscDSSetUp(PetscDS); 49 PETSC_EXTERN PetscErrorCode PetscDSSetFromOptions(PetscDS); 50 PETSC_STATIC_INLINE PetscErrorCode PetscDSViewFromOptions(PetscDS A,PetscObject B,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,B,name);} 51 52 PETSC_EXTERN PetscErrorCode PetscDSView(PetscDS,PetscViewer); 53 PETSC_EXTERN PetscErrorCode PetscDSRegister(const char [], PetscErrorCode (*)(PetscDS)); 54 PETSC_EXTERN PetscErrorCode PetscDSRegisterDestroy(void); 55 56 PETSC_EXTERN PetscErrorCode PetscDSGetSpatialDimension(PetscDS, PetscInt *); 57 PETSC_EXTERN PetscErrorCode PetscDSGetNumFields(PetscDS, PetscInt *); 58 PETSC_EXTERN PetscErrorCode PetscDSGetTotalDimension(PetscDS, PetscInt *); 59 PETSC_EXTERN PetscErrorCode PetscDSGetTotalComponents(PetscDS, PetscInt *); 60 PETSC_EXTERN PetscErrorCode PetscDSGetFieldIndex(PetscDS, PetscObject, PetscInt *); 61 PETSC_EXTERN PetscErrorCode PetscDSGetFieldSize(PetscDS, PetscInt, PetscInt *); 62 PETSC_EXTERN PetscErrorCode PetscDSGetFieldOffset(PetscDS, PetscInt, PetscInt *); 63 PETSC_EXTERN PetscErrorCode PetscDSGetDimensions(PetscDS, PetscInt *[]); 64 PETSC_EXTERN PetscErrorCode PetscDSGetComponents(PetscDS, PetscInt *[]); 65 PETSC_EXTERN PetscErrorCode PetscDSGetComponentOffset(PetscDS, PetscInt, PetscInt *); 66 PETSC_EXTERN PetscErrorCode PetscDSGetComponentOffsets(PetscDS, PetscInt *[]); 67 PETSC_EXTERN PetscErrorCode PetscDSGetComponentDerivativeOffsets(PetscDS, PetscInt *[]); 68 69 PETSC_EXTERN PetscErrorCode PetscDSGetDiscretization(PetscDS, PetscInt, PetscObject *); 70 PETSC_EXTERN PetscErrorCode PetscDSSetDiscretization(PetscDS, PetscInt, PetscObject); 71 PETSC_EXTERN PetscErrorCode PetscDSAddDiscretization(PetscDS, PetscObject); 72 PETSC_EXTERN PetscErrorCode PetscDSGetImplicit(PetscDS, PetscInt, PetscBool*); 73 PETSC_EXTERN PetscErrorCode PetscDSSetImplicit(PetscDS, PetscInt, PetscBool); 74 PETSC_EXTERN PetscErrorCode PetscDSGetAdjacency(PetscDS, PetscInt, PetscBool*, PetscBool*); 75 PETSC_EXTERN PetscErrorCode PetscDSSetAdjacency(PetscDS, PetscInt, PetscBool, PetscBool); 76 PETSC_EXTERN PetscErrorCode PetscDSGetConstants(PetscDS, PetscInt *, const PetscScalar *[]); 77 PETSC_EXTERN PetscErrorCode PetscDSSetConstants(PetscDS, PetscInt, PetscScalar[]); 78 PETSC_EXTERN PetscErrorCode PetscDSGetObjective(PetscDS, PetscInt, 79 void (**)(PetscInt, PetscInt, PetscInt, 80 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 81 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 82 PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])); 83 PETSC_EXTERN PetscErrorCode PetscDSSetObjective(PetscDS, PetscInt, 84 void (*)(PetscInt, PetscInt, PetscInt, 85 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 86 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 87 PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])); 88 PETSC_EXTERN PetscErrorCode PetscDSGetResidual(PetscDS, PetscInt, 89 void (**)(PetscInt, PetscInt, PetscInt, 90 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 91 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 92 PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 93 void (**)(PetscInt, PetscInt, PetscInt, 94 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 95 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 96 PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])); 97 PETSC_EXTERN PetscErrorCode PetscDSSetResidual(PetscDS, PetscInt, 98 void (*)(PetscInt, PetscInt, PetscInt, 99 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 100 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 101 PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 102 void (*)(PetscInt, PetscInt, PetscInt, 103 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 104 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 105 PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])); 106 PETSC_EXTERN PetscErrorCode PetscDSHasJacobian(PetscDS, PetscBool *); 107 PETSC_EXTERN PetscErrorCode PetscDSGetJacobian(PetscDS, PetscInt, PetscInt, 108 void (**)(PetscInt, PetscInt, PetscInt, 109 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 110 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 111 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 112 void (**)(PetscInt, PetscInt, PetscInt, 113 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 114 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 115 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 116 void (**)(PetscInt, PetscInt, PetscInt, 117 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 118 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 119 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 120 void (**)(PetscInt, PetscInt, PetscInt, 121 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 122 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 123 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])); 124 PETSC_EXTERN PetscErrorCode PetscDSSetJacobian(PetscDS, PetscInt, PetscInt, 125 void (*)(PetscInt, PetscInt, PetscInt, 126 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 127 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 128 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 129 void (*)(PetscInt, PetscInt, PetscInt, 130 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 131 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 132 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 133 void (*)(PetscInt, PetscInt, PetscInt, 134 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 135 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 136 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 137 void (*)(PetscInt, PetscInt, PetscInt, 138 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 139 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 140 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])); 141 PETSC_EXTERN PetscErrorCode PetscDSHasJacobianPreconditioner(PetscDS, PetscBool *); 142 PETSC_EXTERN PetscErrorCode PetscDSGetJacobianPreconditioner(PetscDS, PetscInt, PetscInt, 143 void (**)(PetscInt, PetscInt, PetscInt, 144 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 145 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 146 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 147 void (**)(PetscInt, PetscInt, PetscInt, 148 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 149 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 150 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 151 void (**)(PetscInt, PetscInt, PetscInt, 152 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 153 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 154 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 155 void (**)(PetscInt, PetscInt, PetscInt, 156 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 157 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 158 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])); 159 PETSC_EXTERN PetscErrorCode PetscDSSetJacobianPreconditioner(PetscDS, PetscInt, PetscInt, 160 void (*)(PetscInt, PetscInt, PetscInt, 161 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 162 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 163 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 164 void (*)(PetscInt, PetscInt, PetscInt, 165 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 166 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 167 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 168 void (*)(PetscInt, PetscInt, PetscInt, 169 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 170 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 171 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 172 void (*)(PetscInt, PetscInt, PetscInt, 173 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 174 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 175 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])); 176 PETSC_EXTERN PetscErrorCode PetscDSHasDynamicJacobian(PetscDS, PetscBool *); 177 PETSC_EXTERN PetscErrorCode PetscDSGetDynamicJacobian(PetscDS, PetscInt, PetscInt, 178 void (**)(PetscInt, PetscInt, PetscInt, 179 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 180 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 181 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 182 void (**)(PetscInt, PetscInt, PetscInt, 183 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 184 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 185 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 186 void (**)(PetscInt, PetscInt, PetscInt, 187 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 188 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 189 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 190 void (**)(PetscInt, PetscInt, PetscInt, 191 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 192 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 193 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])); 194 PETSC_EXTERN PetscErrorCode PetscDSSetDynamicJacobian(PetscDS, PetscInt, PetscInt, 195 void (*)(PetscInt, PetscInt, PetscInt, 196 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 197 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 198 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 199 void (*)(PetscInt, PetscInt, PetscInt, 200 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 201 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 202 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 203 void (*)(PetscInt, PetscInt, PetscInt, 204 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 205 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 206 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 207 void (*)(PetscInt, PetscInt, PetscInt, 208 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 209 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 210 PetscReal, PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])); 211 PETSC_EXTERN PetscErrorCode PetscDSGetRiemannSolver(PetscDS, PetscInt, 212 void (**)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *)); 213 PETSC_EXTERN PetscErrorCode PetscDSSetRiemannSolver(PetscDS, PetscInt, 214 void (*)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *)); 215 PETSC_EXTERN PetscErrorCode PetscDSGetUpdate(PetscDS, PetscInt, 216 void (**)(PetscInt, PetscInt, PetscInt, 217 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 218 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 219 PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])); 220 PETSC_EXTERN PetscErrorCode PetscDSSetUpdate(PetscDS, PetscInt, 221 void (*)(PetscInt, PetscInt, PetscInt, 222 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 223 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 224 PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])); 225 PETSC_EXTERN PetscErrorCode PetscDSGetContext(PetscDS, PetscInt, void **); 226 PETSC_EXTERN PetscErrorCode PetscDSSetContext(PetscDS, PetscInt, void *); 227 PETSC_EXTERN PetscErrorCode PetscDSGetBdResidual(PetscDS, PetscInt, 228 void (**)(PetscInt, PetscInt, PetscInt, 229 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 230 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 231 PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 232 void (**)(PetscInt, PetscInt, PetscInt, 233 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 234 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 235 PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])); 236 PETSC_EXTERN PetscErrorCode PetscDSSetBdResidual(PetscDS, PetscInt, 237 void (*)(PetscInt, PetscInt, PetscInt, 238 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 239 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 240 PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 241 void (*)(PetscInt, PetscInt, PetscInt, 242 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 243 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 244 PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])); 245 PETSC_EXTERN PetscErrorCode PetscDSGetBdJacobian(PetscDS, PetscInt, PetscInt, 246 void (**)(PetscInt, PetscInt, PetscInt, 247 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 248 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 249 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 250 void (**)(PetscInt, PetscInt, PetscInt, 251 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 252 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 253 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 254 void (**)(PetscInt, PetscInt, PetscInt, 255 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 256 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 257 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 258 void (**)(PetscInt, PetscInt, PetscInt, 259 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 260 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 261 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])); 262 PETSC_EXTERN PetscErrorCode PetscDSSetBdJacobian(PetscDS, PetscInt, PetscInt, 263 void (*)(PetscInt, PetscInt, PetscInt, 264 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 265 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 266 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 267 void (*)(PetscInt, PetscInt, PetscInt, 268 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 269 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 270 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 271 void (*)(PetscInt, PetscInt, PetscInt, 272 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 273 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 274 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), 275 void (*)(PetscInt, PetscInt, PetscInt, 276 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 277 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 278 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[])); 279 PETSC_EXTERN PetscErrorCode PetscDSGetExactSolution(PetscDS, PetscInt, PetscErrorCode (**)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *)); 280 PETSC_EXTERN PetscErrorCode PetscDSSetExactSolution(PetscDS, PetscInt, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *)); 281 PETSC_EXTERN PetscErrorCode PetscDSGetTabulation(PetscDS, PetscReal ***, PetscReal ***); 282 PETSC_EXTERN PetscErrorCode PetscDSGetFaceTabulation(PetscDS, PetscReal ***, PetscReal ***); 283 PETSC_EXTERN PetscErrorCode PetscDSGetEvaluationArrays(PetscDS, PetscScalar **, PetscScalar **, PetscScalar **); 284 PETSC_EXTERN PetscErrorCode PetscDSGetWeakFormArrays(PetscDS, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **); 285 PETSC_EXTERN PetscErrorCode PetscDSGetRefCoordArrays(PetscDS, PetscReal **, PetscScalar **); 286 PETSC_EXTERN PetscErrorCode PetscDSCopyEquations(PetscDS, PetscDS); 287 PETSC_EXTERN PetscErrorCode PetscDSAddBoundary(PetscDS, DMBoundaryConditionType, const char[], const char[], PetscInt, PetscInt, const PetscInt *, void (*)(void), PetscInt, const PetscInt *, void *); 288 PETSC_EXTERN PetscErrorCode PetscDSGetNumBoundary(PetscDS, PetscInt *); 289 PETSC_EXTERN PetscErrorCode PetscDSGetBoundary(PetscDS, PetscInt, DMBoundaryConditionType *, const char **, const char **, PetscInt *, PetscInt *, const PetscInt **, void (**)(void), PetscInt *, const PetscInt **, void **); 290 PETSC_EXTERN PetscErrorCode PetscDSCopyBoundary(PetscDS, PetscDS); 291 292 #endif 293