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[], 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[], 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[], 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[], PetscScalar[]); 40 typedef void (*PetscRiemannFunc)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[], void *); 41 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 PetscDSGetTotalBdDimension(PetscDS, PetscInt *); 60 PETSC_EXTERN PetscErrorCode PetscDSGetTotalComponents(PetscDS, PetscInt *); 61 PETSC_EXTERN PetscErrorCode PetscDSGetFieldIndex(PetscDS, PetscObject, PetscInt *); 62 PETSC_EXTERN PetscErrorCode PetscDSGetFieldSize(PetscDS, PetscInt, PetscInt *); 63 PETSC_EXTERN PetscErrorCode PetscDSGetFieldOffset(PetscDS, PetscInt, PetscInt *); 64 PETSC_EXTERN PetscErrorCode PetscDSGetBdFieldOffset(PetscDS, PetscInt, PetscInt *); 65 PETSC_EXTERN PetscErrorCode PetscDSGetComponentOffset(PetscDS, PetscInt, PetscInt *); 66 PETSC_EXTERN PetscErrorCode PetscDSGetComponentOffsets(PetscDS, PetscInt *[]); 67 PETSC_EXTERN PetscErrorCode PetscDSGetComponentBdOffsets(PetscDS, PetscInt *[]); 68 PETSC_EXTERN PetscErrorCode PetscDSGetComponentDerivativeOffsets(PetscDS, PetscInt *[]); 69 PETSC_EXTERN PetscErrorCode PetscDSGetComponentBdDerivativeOffsets(PetscDS, PetscInt *[]); 70 71 PETSC_EXTERN PetscErrorCode PetscDSGetDiscretization(PetscDS, PetscInt, PetscObject *); 72 PETSC_EXTERN PetscErrorCode PetscDSSetDiscretization(PetscDS, PetscInt, PetscObject); 73 PETSC_EXTERN PetscErrorCode PetscDSAddDiscretization(PetscDS, PetscObject); 74 PETSC_EXTERN PetscErrorCode PetscDSGetBdDiscretization(PetscDS, PetscInt, PetscObject *); 75 PETSC_EXTERN PetscErrorCode PetscDSSetBdDiscretization(PetscDS, PetscInt, PetscObject); 76 PETSC_EXTERN PetscErrorCode PetscDSAddBdDiscretization(PetscDS, PetscObject); 77 PETSC_EXTERN PetscErrorCode PetscDSGetImplicit(PetscDS, PetscInt, PetscBool*); 78 PETSC_EXTERN PetscErrorCode PetscDSSetImplicit(PetscDS, PetscInt, PetscBool); 79 PETSC_EXTERN PetscErrorCode PetscDSGetAdjacency(PetscDS, PetscInt, PetscBool*, PetscBool*); 80 PETSC_EXTERN PetscErrorCode PetscDSSetAdjacency(PetscDS, PetscInt, PetscBool, PetscBool); 81 PETSC_EXTERN PetscErrorCode PetscDSGetObjective(PetscDS, PetscInt, 82 void (**)(PetscInt, PetscInt, PetscInt, 83 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 84 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 85 PetscReal, const PetscReal[], PetscScalar[])); 86 PETSC_EXTERN PetscErrorCode PetscDSSetObjective(PetscDS, PetscInt, 87 void (*)(PetscInt, PetscInt, PetscInt, 88 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 89 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 90 PetscReal, const PetscReal[], PetscScalar[])); 91 PETSC_EXTERN PetscErrorCode PetscDSGetResidual(PetscDS, PetscInt, 92 void (**)(PetscInt, PetscInt, PetscInt, 93 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 94 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 95 PetscReal, const PetscReal[], PetscScalar[]), 96 void (**)(PetscInt, PetscInt, PetscInt, 97 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 98 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 99 PetscReal, const PetscReal[], PetscScalar[])); 100 PETSC_EXTERN PetscErrorCode PetscDSSetResidual(PetscDS, PetscInt, 101 void (*)(PetscInt, PetscInt, PetscInt, 102 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 103 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 104 PetscReal, const PetscReal[], PetscScalar[]), 105 void (*)(PetscInt, PetscInt, PetscInt, 106 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 107 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 108 PetscReal, const PetscReal[], PetscScalar[])); 109 PETSC_EXTERN PetscErrorCode PetscDSHasJacobian(PetscDS, PetscBool *); 110 PETSC_EXTERN PetscErrorCode PetscDSGetJacobian(PetscDS, PetscInt, PetscInt, 111 void (**)(PetscInt, PetscInt, PetscInt, 112 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 113 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 114 PetscReal, PetscReal, const PetscReal[], PetscScalar[]), 115 void (**)(PetscInt, PetscInt, PetscInt, 116 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 117 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 118 PetscReal, PetscReal, const PetscReal[], PetscScalar[]), 119 void (**)(PetscInt, PetscInt, PetscInt, 120 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 121 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 122 PetscReal, PetscReal, const PetscReal[], PetscScalar[]), 123 void (**)(PetscInt, PetscInt, PetscInt, 124 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 125 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 126 PetscReal, PetscReal, const PetscReal[], PetscScalar[])); 127 PETSC_EXTERN PetscErrorCode PetscDSSetJacobian(PetscDS, PetscInt, PetscInt, 128 void (*)(PetscInt, PetscInt, PetscInt, 129 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 130 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 131 PetscReal, PetscReal, const PetscReal[], PetscScalar[]), 132 void (*)(PetscInt, PetscInt, PetscInt, 133 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 134 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 135 PetscReal, PetscReal, const PetscReal[], PetscScalar[]), 136 void (*)(PetscInt, PetscInt, PetscInt, 137 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 138 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 139 PetscReal, PetscReal, const PetscReal[], PetscScalar[]), 140 void (*)(PetscInt, PetscInt, PetscInt, 141 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 142 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 143 PetscReal, PetscReal, const PetscReal[], PetscScalar[])); 144 PETSC_EXTERN PetscErrorCode PetscDSHasJacobianPreconditioner(PetscDS, PetscBool *); 145 PETSC_EXTERN PetscErrorCode PetscDSGetJacobianPreconditioner(PetscDS, PetscInt, PetscInt, 146 void (**)(PetscInt, PetscInt, PetscInt, 147 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 148 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 149 PetscReal, PetscReal, const PetscReal[], PetscScalar[]), 150 void (**)(PetscInt, PetscInt, PetscInt, 151 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 152 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 153 PetscReal, PetscReal, const PetscReal[], PetscScalar[]), 154 void (**)(PetscInt, PetscInt, PetscInt, 155 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 156 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 157 PetscReal, PetscReal, const PetscReal[], PetscScalar[]), 158 void (**)(PetscInt, PetscInt, PetscInt, 159 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 160 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 161 PetscReal, PetscReal, const PetscReal[], PetscScalar[])); 162 PETSC_EXTERN PetscErrorCode PetscDSSetJacobianPreconditioner(PetscDS, PetscInt, PetscInt, 163 void (*)(PetscInt, PetscInt, PetscInt, 164 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 165 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 166 PetscReal, PetscReal, const PetscReal[], PetscScalar[]), 167 void (*)(PetscInt, PetscInt, PetscInt, 168 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 169 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 170 PetscReal, PetscReal, const PetscReal[], PetscScalar[]), 171 void (*)(PetscInt, PetscInt, PetscInt, 172 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 173 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 174 PetscReal, PetscReal, const PetscReal[], PetscScalar[]), 175 void (*)(PetscInt, PetscInt, PetscInt, 176 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 177 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 178 PetscReal, PetscReal, const PetscReal[], PetscScalar[])); 179 PETSC_EXTERN PetscErrorCode PetscDSHasDynamicJacobian(PetscDS, PetscBool *); 180 PETSC_EXTERN PetscErrorCode PetscDSGetDynamicJacobian(PetscDS, PetscInt, PetscInt, 181 void (**)(PetscInt, PetscInt, PetscInt, 182 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 183 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 184 PetscReal, PetscReal, const PetscReal[], PetscScalar[]), 185 void (**)(PetscInt, PetscInt, PetscInt, 186 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 187 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 188 PetscReal, PetscReal, const PetscReal[], PetscScalar[]), 189 void (**)(PetscInt, PetscInt, PetscInt, 190 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 191 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 192 PetscReal, PetscReal, const PetscReal[], PetscScalar[]), 193 void (**)(PetscInt, PetscInt, PetscInt, 194 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 195 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 196 PetscReal, PetscReal, const PetscReal[], PetscScalar[])); 197 PETSC_EXTERN PetscErrorCode PetscDSSetDynamicJacobian(PetscDS, PetscInt, PetscInt, 198 void (*)(PetscInt, PetscInt, PetscInt, 199 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 200 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 201 PetscReal, PetscReal, const PetscReal[], PetscScalar[]), 202 void (*)(PetscInt, PetscInt, PetscInt, 203 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 204 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 205 PetscReal, PetscReal, const PetscReal[], PetscScalar[]), 206 void (*)(PetscInt, PetscInt, PetscInt, 207 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 208 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 209 PetscReal, PetscReal, const PetscReal[], PetscScalar[]), 210 void (*)(PetscInt, PetscInt, PetscInt, 211 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 212 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 213 PetscReal, PetscReal, const PetscReal[], PetscScalar[])); 214 PETSC_EXTERN PetscErrorCode PetscDSGetRiemannSolver(PetscDS, PetscInt, 215 void (**)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[], void *)); 216 PETSC_EXTERN PetscErrorCode PetscDSSetRiemannSolver(PetscDS, PetscInt, 217 void (*)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[], void *)); 218 PETSC_EXTERN PetscErrorCode PetscDSGetContext(PetscDS, PetscInt, void **); 219 PETSC_EXTERN PetscErrorCode PetscDSSetContext(PetscDS, PetscInt, void *); 220 PETSC_EXTERN PetscErrorCode PetscDSGetBdResidual(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[], const PetscReal[], PetscScalar[]), 225 void (**)(PetscInt, PetscInt, PetscInt, 226 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 227 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 228 PetscReal, const PetscReal[], const PetscReal[], PetscScalar[])); 229 PETSC_EXTERN PetscErrorCode PetscDSSetBdResidual(PetscDS, PetscInt, 230 void (*)(PetscInt, PetscInt, PetscInt, 231 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 232 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 233 PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]), 234 void (*)(PetscInt, PetscInt, PetscInt, 235 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 236 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 237 PetscReal, const PetscReal[], const PetscReal[], PetscScalar[])); 238 PETSC_EXTERN PetscErrorCode PetscDSGetBdJacobian(PetscDS, PetscInt, PetscInt, 239 void (**)(PetscInt, PetscInt, PetscInt, 240 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 241 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 242 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]), 243 void (**)(PetscInt, PetscInt, PetscInt, 244 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 245 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 246 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]), 247 void (**)(PetscInt, PetscInt, PetscInt, 248 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 249 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 250 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]), 251 void (**)(PetscInt, PetscInt, PetscInt, 252 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 253 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 254 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[])); 255 PETSC_EXTERN PetscErrorCode PetscDSSetBdJacobian(PetscDS, PetscInt, PetscInt, 256 void (*)(PetscInt, PetscInt, PetscInt, 257 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 258 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 259 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]), 260 void (*)(PetscInt, PetscInt, PetscInt, 261 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 262 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 263 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]), 264 void (*)(PetscInt, PetscInt, PetscInt, 265 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 266 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 267 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[]), 268 void (*)(PetscInt, PetscInt, PetscInt, 269 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 270 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 271 PetscReal, PetscReal, const PetscReal[], const PetscReal[], PetscScalar[])); 272 PETSC_EXTERN PetscErrorCode PetscDSGetTabulation(PetscDS, PetscReal ***, PetscReal ***); 273 PETSC_EXTERN PetscErrorCode PetscDSGetBdTabulation(PetscDS, PetscReal ***, PetscReal ***); 274 PETSC_EXTERN PetscErrorCode PetscDSGetEvaluationArrays(PetscDS, PetscScalar **, PetscScalar **, PetscScalar **); 275 PETSC_EXTERN PetscErrorCode PetscDSGetWeakFormArrays(PetscDS, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **, PetscScalar **); 276 PETSC_EXTERN PetscErrorCode PetscDSGetRefCoordArrays(PetscDS, PetscReal **, PetscScalar **); 277 PETSC_EXTERN PetscErrorCode PetscDSCopyEquations(PetscDS, PetscDS); 278 279 #endif 280