1c4762a1bSJed Brown #include "finitevolume1d.h" 2c4762a1bSJed Brown #include <petscdmda.h> 3c4762a1bSJed Brown #include <petscdraw.h> 4c4762a1bSJed Brown #include <petsc/private/tsimpl.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown #include <petsc/private/kernels/blockinvert.h> /* For the Kernel_*_gets_* stuff for BAIJ */ 784ff8c8bSHong Zhang const char *FVBCTypes[] = {"PERIODIC", "OUTFLOW", "INFLOW", "FVBCType", "FVBC_", 0}; 8c4762a1bSJed Brown 9d71ae5a4SJacob Faibussowitsch static inline PetscReal Sgn(PetscReal a) 10d71ae5a4SJacob Faibussowitsch { 119371c9d4SSatish Balay return (a < 0) ? -1 : 1; 129371c9d4SSatish Balay } 13d71ae5a4SJacob Faibussowitsch static inline PetscReal Abs(PetscReal a) 14d71ae5a4SJacob Faibussowitsch { 159371c9d4SSatish Balay return (a < 0) ? 0 : a; 169371c9d4SSatish Balay } 17d71ae5a4SJacob Faibussowitsch static inline PetscReal Sqr(PetscReal a) 18d71ae5a4SJacob Faibussowitsch { 199371c9d4SSatish Balay return a * a; 209371c9d4SSatish Balay } 2184ff8c8bSHong Zhang 22d71ae5a4SJacob Faibussowitsch PETSC_UNUSED static inline PetscReal MinAbs(PetscReal a, PetscReal b) 23d71ae5a4SJacob Faibussowitsch { 249371c9d4SSatish Balay return (PetscAbs(a) < PetscAbs(b)) ? a : b; 259371c9d4SSatish Balay } 26d71ae5a4SJacob Faibussowitsch static inline PetscReal MinMod2(PetscReal a, PetscReal b) 27d71ae5a4SJacob Faibussowitsch { 289371c9d4SSatish Balay return (a * b < 0) ? 0 : Sgn(a) * PetscMin(PetscAbs(a), PetscAbs(b)); 299371c9d4SSatish Balay } 30d71ae5a4SJacob Faibussowitsch static inline PetscReal MaxMod2(PetscReal a, PetscReal b) 31d71ae5a4SJacob Faibussowitsch { 329371c9d4SSatish Balay return (a * b < 0) ? 0 : Sgn(a) * PetscMax(PetscAbs(a), PetscAbs(b)); 339371c9d4SSatish Balay } 34d71ae5a4SJacob Faibussowitsch static inline PetscReal MinMod3(PetscReal a, PetscReal b, PetscReal c) 35d71ae5a4SJacob Faibussowitsch { 369371c9d4SSatish Balay return (a * b < 0 || a * c < 0) ? 0 : Sgn(a) * PetscMin(PetscAbs(a), PetscMin(PetscAbs(b), PetscAbs(c))); 379371c9d4SSatish Balay } 38c4762a1bSJed Brown 39c4762a1bSJed Brown /* ----------------------- Lots of limiters, these could go in a separate library ------------------------- */ 40d71ae5a4SJacob Faibussowitsch void Limit_Upwind(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt) 41d71ae5a4SJacob Faibussowitsch { 42c4762a1bSJed Brown PetscInt i; 43c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = 0; 44c4762a1bSJed Brown } 45d71ae5a4SJacob Faibussowitsch void Limit_LaxWendroff(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt) 46d71ae5a4SJacob Faibussowitsch { 47c4762a1bSJed Brown PetscInt i; 48c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = jR[i]; 49c4762a1bSJed Brown } 50d71ae5a4SJacob Faibussowitsch void Limit_BeamWarming(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt) 51d71ae5a4SJacob Faibussowitsch { 52c4762a1bSJed Brown PetscInt i; 53c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = jL[i]; 54c4762a1bSJed Brown } 55d71ae5a4SJacob Faibussowitsch void Limit_Fromm(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt) 56d71ae5a4SJacob Faibussowitsch { 57c4762a1bSJed Brown PetscInt i; 58c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = 0.5 * (jL[i] + jR[i]); 59c4762a1bSJed Brown } 60d71ae5a4SJacob Faibussowitsch void Limit_Minmod(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt) 61d71ae5a4SJacob Faibussowitsch { 62c4762a1bSJed Brown PetscInt i; 63c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod2(jL[i], jR[i]); 64c4762a1bSJed Brown } 65d71ae5a4SJacob Faibussowitsch void Limit_Superbee(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt) 66d71ae5a4SJacob Faibussowitsch { 67c4762a1bSJed Brown PetscInt i; 68c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MaxMod2(MinMod2(jL[i], 2 * jR[i]), MinMod2(2 * jL[i], jR[i])); 69c4762a1bSJed Brown } 70d71ae5a4SJacob Faibussowitsch void Limit_MC(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt) 71d71ae5a4SJacob Faibussowitsch { 72c4762a1bSJed Brown PetscInt i; 73c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i], 0.5 * (jL[i] + jR[i]), 2 * jR[i]); 74c4762a1bSJed Brown } 75d71ae5a4SJacob Faibussowitsch void Limit_VanLeer(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt) 76d71ae5a4SJacob Faibussowitsch { /* phi = (t + abs(t)) / (1 + abs(t)) */ 77c4762a1bSJed Brown PetscInt i; 78c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = (jL[i] * Abs(jR[i]) + Abs(jL[i]) * jR[i]) / (Abs(jL[i]) + Abs(jR[i]) + 1e-15); 79c4762a1bSJed Brown } 80c4762a1bSJed Brown void Limit_VanAlbada(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt) /* differentiable */ 81c4762a1bSJed Brown { /* phi = (t + t^2) / (1 + t^2) */ 82c4762a1bSJed Brown PetscInt i; 83c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = (jL[i] * Sqr(jR[i]) + Sqr(jL[i]) * jR[i]) / (Sqr(jL[i]) + Sqr(jR[i]) + 1e-15); 84c4762a1bSJed Brown } 85d71ae5a4SJacob Faibussowitsch void Limit_VanAlbadaTVD(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt) 86d71ae5a4SJacob Faibussowitsch { /* phi = (t + t^2) / (1 + t^2) */ 87c4762a1bSJed Brown PetscInt i; 88c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = (jL[i] * jR[i] < 0) ? 0 : (jL[i] * Sqr(jR[i]) + Sqr(jL[i]) * jR[i]) / (Sqr(jL[i]) + Sqr(jR[i]) + 1e-15); 89c4762a1bSJed Brown } 90c4762a1bSJed Brown void Limit_Koren(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt) /* differentiable */ 91c4762a1bSJed Brown { /* phi = (t + 2*t^2) / (2 - t + 2*t^2) */ 92c4762a1bSJed Brown PetscInt i; 93c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = ((jL[i] * Sqr(jR[i]) + 2 * Sqr(jL[i]) * jR[i]) / (2 * Sqr(jL[i]) - jL[i] * jR[i] + 2 * Sqr(jR[i]) + 1e-15)); 94c4762a1bSJed Brown } 95c4762a1bSJed Brown void Limit_KorenSym(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt) /* differentiable */ 96c4762a1bSJed Brown { /* Symmetric version of above */ 97c4762a1bSJed Brown PetscInt i; 98c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = (1.5 * (jL[i] * Sqr(jR[i]) + Sqr(jL[i]) * jR[i]) / (2 * Sqr(jL[i]) - jL[i] * jR[i] + 2 * Sqr(jR[i]) + 1e-15)); 99c4762a1bSJed Brown } 100d71ae5a4SJacob Faibussowitsch void Limit_Koren3(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt) 101d71ae5a4SJacob Faibussowitsch { /* Eq 11 of Cada-Torrilhon 2009 */ 102c4762a1bSJed Brown PetscInt i; 103c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i], (jL[i] + 2 * jR[i]) / 3, 2 * jR[i]); 104c4762a1bSJed Brown } 105d71ae5a4SJacob Faibussowitsch static PetscReal CadaTorrilhonPhiHatR_Eq13(PetscReal L, PetscReal R) 106d71ae5a4SJacob Faibussowitsch { 107c4762a1bSJed Brown return PetscMax(0, PetscMin((L + 2 * R) / 3, PetscMax(-0.5 * L, PetscMin(2 * L, PetscMin((L + 2 * R) / 3, 1.6 * R))))); 108c4762a1bSJed Brown } 109d71ae5a4SJacob Faibussowitsch void Limit_CadaTorrilhon2(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt) 110d71ae5a4SJacob Faibussowitsch { /* Cada-Torrilhon 2009, Eq 13 */ 111c4762a1bSJed Brown PetscInt i; 112c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = CadaTorrilhonPhiHatR_Eq13(jL[i], jR[i]); 113c4762a1bSJed Brown } 114d71ae5a4SJacob Faibussowitsch void Limit_CadaTorrilhon3R(PetscReal r, LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt) 115d71ae5a4SJacob Faibussowitsch { /* Cada-Torrilhon 2009, Eq 22 */ 116c4762a1bSJed Brown /* They recommend 0.001 < r < 1, but larger values are more accurate in smooth regions */ 117c4762a1bSJed Brown const PetscReal eps = 1e-7, hx = info->hx; 118c4762a1bSJed Brown PetscInt i; 119c4762a1bSJed Brown for (i = 0; i < info->m; i++) { 120c4762a1bSJed Brown const PetscReal eta = (Sqr(jL[i]) + Sqr(jR[i])) / Sqr(r * hx); 121c4762a1bSJed Brown lmt[i] = ((eta < 1 - eps) ? (jL[i] + 2 * jR[i]) / 3 : ((eta > 1 + eps) ? CadaTorrilhonPhiHatR_Eq13(jL[i], jR[i]) : 0.5 * ((1 - (eta - 1) / eps) * (jL[i] + 2 * jR[i]) / 3 + (1 + (eta + 1) / eps) * CadaTorrilhonPhiHatR_Eq13(jL[i], jR[i])))); 122c4762a1bSJed Brown } 123c4762a1bSJed Brown } 124d71ae5a4SJacob Faibussowitsch void Limit_CadaTorrilhon3R0p1(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt) 125d71ae5a4SJacob Faibussowitsch { 126c4762a1bSJed Brown Limit_CadaTorrilhon3R(0.1, info, jL, jR, lmt); 127c4762a1bSJed Brown } 128d71ae5a4SJacob Faibussowitsch void Limit_CadaTorrilhon3R1(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt) 129d71ae5a4SJacob Faibussowitsch { 130c4762a1bSJed Brown Limit_CadaTorrilhon3R(1, info, jL, jR, lmt); 131c4762a1bSJed Brown } 132d71ae5a4SJacob Faibussowitsch void Limit_CadaTorrilhon3R10(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt) 133d71ae5a4SJacob Faibussowitsch { 134c4762a1bSJed Brown Limit_CadaTorrilhon3R(10, info, jL, jR, lmt); 135c4762a1bSJed Brown } 136d71ae5a4SJacob Faibussowitsch void Limit_CadaTorrilhon3R100(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, PetscScalar *lmt) 137d71ae5a4SJacob Faibussowitsch { 138c4762a1bSJed Brown Limit_CadaTorrilhon3R(100, info, jL, jR, lmt); 139c4762a1bSJed Brown } 140c4762a1bSJed Brown 141c4762a1bSJed Brown /* ----------------------- Limiters for split systems ------------------------- */ 142c4762a1bSJed Brown 143d71ae5a4SJacob Faibussowitsch void Limit2_Upwind(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, const PetscInt sf, const PetscInt fs, PetscInt n, PetscScalar *lmt) 144d71ae5a4SJacob Faibussowitsch { 145c4762a1bSJed Brown PetscInt i; 146c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = 0; 147c4762a1bSJed Brown } 148d71ae5a4SJacob Faibussowitsch void Limit2_LaxWendroff(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, const PetscInt sf, const PetscInt fs, PetscInt n, PetscScalar *lmt) 149d71ae5a4SJacob Faibussowitsch { 150c4762a1bSJed Brown PetscInt i; 151c4762a1bSJed Brown if (n < sf - 1) { /* slow components */ 152c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = jR[i] / info->hxs; 153c4762a1bSJed Brown } else if (n == sf - 1) { /* slow component which is next to fast components */ 154c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = jR[i] / (info->hxs / 2.0 + info->hxf / 2.0); 155c4762a1bSJed Brown } else if (n == sf) { /* fast component which is next to slow components */ 156c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = jR[i] / info->hxf; 157c4762a1bSJed Brown } else if (n > sf && n < fs - 1) { /* fast components */ 158c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = jR[i] / info->hxf; 159c4762a1bSJed Brown } else if (n == fs - 1) { /* fast component next to slow components */ 160c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = jR[i] / (info->hxf / 2.0 + info->hxs / 2.0); 161c4762a1bSJed Brown } else if (n == fs) { /* slow component next to fast components */ 162c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = jR[i] / info->hxs; 163c4762a1bSJed Brown } else { /* slow components */ 164c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = jR[i] / info->hxs; 165c4762a1bSJed Brown } 166c4762a1bSJed Brown } 167d71ae5a4SJacob Faibussowitsch void Limit2_BeamWarming(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, const PetscInt sf, const PetscInt fs, PetscInt n, PetscScalar *lmt) 168d71ae5a4SJacob Faibussowitsch { 169c4762a1bSJed Brown PetscInt i; 170c4762a1bSJed Brown if (n < sf - 1) { 171c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = jL[i] / info->hxs; 172c4762a1bSJed Brown } else if (n == sf - 1) { 173c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = jL[i] / info->hxs; 174c4762a1bSJed Brown } else if (n == sf) { 175c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = jL[i] / (info->hxs / 2.0 + info->hxf / 2.0); 176c4762a1bSJed Brown } else if (n > sf && n < fs - 1) { 177c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = jL[i] / info->hxf; 178c4762a1bSJed Brown } else if (n == fs - 1) { 179c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = jL[i] / info->hxf; 180c4762a1bSJed Brown } else if (n == fs) { 181c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = jL[i] / (info->hxf / 2.0 + info->hxs / 2.0); 182c4762a1bSJed Brown } else { 183c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = jL[i] / info->hxs; 184c4762a1bSJed Brown } 185c4762a1bSJed Brown } 186d71ae5a4SJacob Faibussowitsch void Limit2_Fromm(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, const PetscInt sf, const PetscInt fs, PetscInt n, PetscScalar *lmt) 187d71ae5a4SJacob Faibussowitsch { 188c4762a1bSJed Brown PetscInt i; 189c4762a1bSJed Brown if (n < sf - 1) { 190c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = 0.5 * (jL[i] + jR[i]) / info->hxs; 191c4762a1bSJed Brown } else if (n == sf - 1) { 192c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = 0.5 * (jL[i] / info->hxs + jR[i] / (info->hxs / 2.0 + info->hxf / 2.0)); 193c4762a1bSJed Brown } else if (n == sf) { 194c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = 0.5 * (jL[i] / (info->hxs / 2.0 + info->hxf / 2.0) + jR[i] / info->hxf); 195c4762a1bSJed Brown } else if (n > sf && n < fs - 1) { 196c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = 0.5 * (jL[i] + jR[i]) / info->hxf; 197c4762a1bSJed Brown } else if (n == fs - 1) { 198c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = 0.5 * (jL[i] / info->hxf + jR[i] / (info->hxf / 2.0 + info->hxs / 2.0)); 199c4762a1bSJed Brown } else if (n == fs) { 200c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = 0.5 * (jL[i] / (info->hxf / 2.0 + info->hxs / 2.0) + jR[i] / info->hxs); 201c4762a1bSJed Brown } else { 202c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = 0.5 * (jL[i] + jR[i]) / info->hxs; 203c4762a1bSJed Brown } 204c4762a1bSJed Brown } 205d71ae5a4SJacob Faibussowitsch void Limit2_Minmod(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, const PetscInt sf, const PetscInt fs, PetscInt n, PetscScalar *lmt) 206d71ae5a4SJacob Faibussowitsch { 207c4762a1bSJed Brown PetscInt i; 208c4762a1bSJed Brown if (n < sf - 1) { 209c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod2(jL[i], jR[i]) / info->hxs; 210c4762a1bSJed Brown } else if (n == sf - 1) { 211c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod2(jL[i] / info->hxs, jR[i] / (info->hxs / 2.0 + info->hxf / 2.0)); 212c4762a1bSJed Brown } else if (n == sf) { 213c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod2(jL[i] / (info->hxs / 2.0 + info->hxf / 2.0), jR[i] / info->hxf); 214c4762a1bSJed Brown } else if (n > sf && n < fs - 1) { 215c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod2(jL[i], jR[i]) / info->hxf; 216c4762a1bSJed Brown } else if (n == fs - 1) { 217c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod2(jL[i] / info->hxf, jR[i] / (info->hxf / 2.0 + info->hxs / 2.0)); 218c4762a1bSJed Brown } else if (n == fs) { 219c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod2(jL[i] / (info->hxf / 2.0 + info->hxs / 2.0), jR[i] / info->hxs); 220c4762a1bSJed Brown } else { 221c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod2(jL[i], jR[i]) / info->hxs; 222c4762a1bSJed Brown } 223c4762a1bSJed Brown } 224d71ae5a4SJacob Faibussowitsch void Limit2_Superbee(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, const PetscInt sf, const PetscInt fs, PetscInt n, PetscScalar *lmt) 225d71ae5a4SJacob Faibussowitsch { 226c4762a1bSJed Brown PetscInt i; 227c4762a1bSJed Brown if (n < sf - 1) { 228c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MaxMod2(MinMod2(jL[i], 2 * jR[i]), MinMod2(2 * jL[i], jR[i])) / info->hxs; 229c4762a1bSJed Brown } else if (n == sf - 1) { 230c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MaxMod2(MinMod2(jL[i] / info->hxs, 2 * jR[i] / (info->hxs / 2.0 + info->hxf / 2.0)), MinMod2(2 * jL[i] / info->hxs, jR[i] / (info->hxs / 2.0 + info->hxf / 2.0))); 231c4762a1bSJed Brown } else if (n == sf) { 232c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MaxMod2(MinMod2(jL[i] / (info->hxs / 2.0 + info->hxf / 2.0), 2 * jR[i] / info->hxf), MinMod2(2 * jL[i] / (info->hxs / 2.0 + info->hxf / 2.0), jR[i] / info->hxf)); 233c4762a1bSJed Brown } else if (n > sf && n < fs - 1) { 234c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod2(MinMod2(jL[i], 2 * jR[i]), MinMod2(2 * jL[i], jR[i])) / info->hxf; 235c4762a1bSJed Brown } else if (n == fs - 1) { 236c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod2(MinMod2(jL[i] / info->hxf, 2 * jR[i] / (info->hxf / 2.0 + info->hxs / 2.0)), MinMod2(2 * jL[i] / info->hxf, jR[i] / (info->hxf / 2.0 + info->hxs / 2.0))); 237c4762a1bSJed Brown } else if (n == fs) { 238c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod2(MinMod2(jL[i] / (info->hxf / 2.0 + info->hxs / 2.0), 2 * jR[i] / info->hxs), MinMod2(2 * jL[i] / (info->hxf / 2.0 + info->hxs / 2.0), jR[i] / info->hxs)); 239c4762a1bSJed Brown } else { 240c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod2(MinMod2(jL[i], 2 * jR[i]), MinMod2(2 * jL[i], jR[i])) / info->hxs; 241c4762a1bSJed Brown } 242c4762a1bSJed Brown } 243d71ae5a4SJacob Faibussowitsch void Limit2_MC(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, const PetscInt sf, const PetscInt fs, PetscInt n, PetscScalar *lmt) 244d71ae5a4SJacob Faibussowitsch { 245c4762a1bSJed Brown PetscInt i; 246c4762a1bSJed Brown if (n < sf - 1) { 247c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i], 0.5 * (jL[i] + jR[i]), 2 * jR[i]) / info->hxs; 248c4762a1bSJed Brown } else if (n == sf - 1) { 249c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i] / info->hxs, 0.5 * (jL[i] / info->hxs + jR[i] / (info->hxf / 2.0 + info->hxs / 2.0)), 2 * jR[i] / (info->hxf / 2.0 + info->hxs / 2.0)); 250c4762a1bSJed Brown } else if (n == sf) { 251c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i] / (info->hxs / 2.0 + info->hxf / 2.0), 0.5 * (jL[i] / (info->hxs / 2.0 + info->hxf / 2.0) + jR[i] / info->hxf), 2 * jR[i] / info->hxf); 252c4762a1bSJed Brown } else if (n > sf && n < fs - 1) { 253c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i], 0.5 * (jL[i] + jR[i]), 2 * jR[i]) / info->hxf; 254c4762a1bSJed Brown } else if (n == fs - 1) { 255c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i] / info->hxf, 0.5 * (jL[i] / info->hxf + jR[i] / (info->hxf / 2.0 + info->hxs / 2.0)), 2 * jR[i] / (info->hxf / 2.0 + info->hxs / 2.0)); 256c4762a1bSJed Brown } else if (n == fs) { 257c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i] / (info->hxf / 2.0 + info->hxs / 2.0), 0.5 * (jL[i] / (info->hxf / 2.0 + info->hxs / 2.0) + jR[i] / info->hxs), 2 * jR[i] / info->hxs); 258c4762a1bSJed Brown } else { 259c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i], 0.5 * (jL[i] + jR[i]), 2 * jR[i]) / info->hxs; 260c4762a1bSJed Brown } 261c4762a1bSJed Brown } 262d71ae5a4SJacob Faibussowitsch void Limit2_Koren3(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, const PetscInt sf, const PetscInt fs, PetscInt n, PetscScalar *lmt) 263d71ae5a4SJacob Faibussowitsch { /* Eq 11 of Cada-Torrilhon 2009 */ 264c4762a1bSJed Brown PetscInt i; 265c4762a1bSJed Brown if (n < sf - 1) { 266c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i], (jL[i] + 2 * jR[i]) / 3, 2 * jR[i]) / info->hxs; 267c4762a1bSJed Brown } else if (n == sf - 1) { 268c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i] / info->hxs, (jL[i] / info->hxs + 2 * jR[i] / (info->hxf / 2.0 + info->hxs / 2.0)) / 3, 2 * jR[i] / (info->hxf / 2.0 + info->hxs / 2.0)); 269c4762a1bSJed Brown } else if (n == sf) { 270c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i] / (info->hxs / 2.0 + info->hxf / 2.0), (jL[i] / (info->hxs / 2.0 + info->hxf / 2.0) + 2 * jR[i] / info->hxf) / 3, 2 * jR[i] / info->hxf); 271c4762a1bSJed Brown } else if (n > sf && n < fs - 1) { 272c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i], (jL[i] + 2 * jR[i]) / 3, 2 * jR[i]) / info->hxf; 273c4762a1bSJed Brown } else if (n == fs - 1) { 274c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i] / info->hxf, (jL[i] / info->hxf + 2 * jR[i] / (info->hxf / 2.0 + info->hxs / 2.0)) / 3, 2 * jR[i] / (info->hxf / 2.0 + info->hxs / 2.0)); 275c4762a1bSJed Brown } else if (n == fs) { 276c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i] / (info->hxf / 2.0 + info->hxs / 2.0), (jL[i] / (info->hxf / 2.0 + info->hxs / 2.0) + 2 * jR[i] / info->hxs) / 3, 2 * jR[i] / info->hxs); 277c4762a1bSJed Brown } else { 278c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i], (jL[i] + 2 * jR[i]) / 3, 2 * jR[i]) / info->hxs; 279c4762a1bSJed Brown } 280c4762a1bSJed Brown } 281c4762a1bSJed Brown 282c4762a1bSJed Brown /* ---- Three-way splitting ---- */ 283d71ae5a4SJacob Faibussowitsch void Limit3_Upwind(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, const PetscInt sm, const PetscInt mf, const PetscInt fm, const PetscInt ms, PetscInt n, PetscScalar *lmt) 284d71ae5a4SJacob Faibussowitsch { 285c4762a1bSJed Brown PetscInt i; 286c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = 0; 287c4762a1bSJed Brown } 288d71ae5a4SJacob Faibussowitsch void Limit3_LaxWendroff(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, const PetscInt sm, const PetscInt mf, const PetscInt fm, const PetscInt ms, PetscInt n, PetscScalar *lmt) 289d71ae5a4SJacob Faibussowitsch { 290c4762a1bSJed Brown PetscInt i; 291c4762a1bSJed Brown if (n < sm - 1 || n > ms) { /* slow components */ 292c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = jR[i] / info->hxs; 293c4762a1bSJed Brown } else if (n == sm - 1 || n == ms - 1) { /* slow component which is next to medium components */ 294c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = jR[i] / (info->hxs / 2.0 + info->hxm / 2.0); 295c4762a1bSJed Brown } else if (n < mf - 1 || n > fm) { /* medium components */ 296c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = jR[i] / info->hxm; 297c4762a1bSJed Brown } else if (n == mf - 1 || n == fm) { /* medium component next to fast components */ 298c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = jR[i] / (info->hxm / 2.0 + info->hxf / 2.0); 299c4762a1bSJed Brown } else { /* fast components */ 300c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = jR[i] / info->hxf; 301c4762a1bSJed Brown } 302c4762a1bSJed Brown } 303d71ae5a4SJacob Faibussowitsch void Limit3_BeamWarming(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, const PetscInt sm, const PetscInt mf, const PetscInt fm, const PetscInt ms, PetscInt n, PetscScalar *lmt) 304d71ae5a4SJacob Faibussowitsch { 305c4762a1bSJed Brown PetscInt i; 306c4762a1bSJed Brown if (n < sm || n > ms) { 307c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = jL[i] / info->hxs; 308c4762a1bSJed Brown } else if (n == sm || n == ms) { 309c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = jL[i] / (info->hxs / 2.0 + info->hxf / 2.0); 310c4762a1bSJed Brown } else if (n < mf || n > fm) { 311c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = jL[i] / info->hxm; 312c4762a1bSJed Brown } else if (n == mf || n == fm) { 313c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = jL[i] / (info->hxm / 2.0 + info->hxf / 2.0); 314c4762a1bSJed Brown } else { 315c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = jL[i] / info->hxf; 316c4762a1bSJed Brown } 317c4762a1bSJed Brown } 318d71ae5a4SJacob Faibussowitsch void Limit3_Fromm(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, const PetscInt sm, const PetscInt mf, const PetscInt fm, const PetscInt ms, PetscInt n, PetscScalar *lmt) 319d71ae5a4SJacob Faibussowitsch { 320c4762a1bSJed Brown PetscInt i; 321c4762a1bSJed Brown if (n < sm - 1 || n > ms) { 322c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = 0.5 * (jL[i] + jR[i]) / info->hxs; 323c4762a1bSJed Brown } else if (n == sm - 1) { 324c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = 0.5 * (jL[i] / info->hxs + jR[i] / (info->hxs / 2.0 + info->hxf / 2.0)); 325c4762a1bSJed Brown } else if (n == sm) { 326c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = 0.5 * (jL[i] / (info->hxs / 2.0 + info->hxm / 2.0) + jR[i] / info->hxm); 327c4762a1bSJed Brown } else if (n == ms - 1) { 328c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = 0.5 * (jL[i] / info->hxm + jR[i] / (info->hxs / 2.0 + info->hxf / 2.0)); 329c4762a1bSJed Brown } else if (n == ms) { 330c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = 0.5 * (jL[i] / (info->hxm / 2.0 + info->hxs / 2.0) + jR[i] / info->hxs); 331c4762a1bSJed Brown } else if (n < mf - 1 || n > fm) { 332c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = 0.5 * (jL[i] + jR[i]) / info->hxm; 333c4762a1bSJed Brown } else if (n == mf - 1) { 334c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = 0.5 * (jL[i] / info->hxm + jR[i] / (info->hxm / 2.0 + info->hxf / 2.0)); 335c4762a1bSJed Brown } else if (n == mf) { 336c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = 0.5 * (jL[i] / (info->hxm / 2.0 + info->hxf / 2.0) + jR[i] / info->hxf); 337c4762a1bSJed Brown } else if (n == fm - 1) { 338c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = 0.5 * (jL[i] / info->hxf + jR[i] / (info->hxf / 2.0 + info->hxm / 2.0)); 339c4762a1bSJed Brown } else if (n == fm) { 340c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = 0.5 * (jL[i] / (info->hxf / 2.0 + info->hxm / 2.0) + jR[i] / info->hxm); 341c4762a1bSJed Brown } else { 342c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = 0.5 * (jL[i] + jR[i]) / info->hxf; 343c4762a1bSJed Brown } 344c4762a1bSJed Brown } 345d71ae5a4SJacob Faibussowitsch void Limit3_Minmod(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, const PetscInt sm, const PetscInt mf, const PetscInt fm, const PetscInt ms, PetscInt n, PetscScalar *lmt) 346d71ae5a4SJacob Faibussowitsch { 347c4762a1bSJed Brown PetscInt i; 348c4762a1bSJed Brown if (n < sm - 1 || n > ms) { 349c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod2(jL[i], jR[i]) / info->hxs; 350c4762a1bSJed Brown } else if (n == sm - 1) { 351c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod2(jL[i] / info->hxs, jR[i] / (info->hxs / 2.0 + info->hxf / 2.0)); 352c4762a1bSJed Brown } else if (n == sm) { 353c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod2(jL[i] / (info->hxs / 2.0 + info->hxf / 2.0), jR[i] / info->hxf); 354c4762a1bSJed Brown } else if (n == ms - 1) { 355c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod2(jL[i] / info->hxm, jR[i] / (info->hxm / 2.0 + info->hxs / 2.0)); 356c4762a1bSJed Brown } else if (n == ms) { 357c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod2(jL[i] / (info->hxm / 2.0 + info->hxs / 2.0), jR[i] / info->hxs); 358c4762a1bSJed Brown } else if (n < mf - 1 || n > fm) { 359c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod2(jL[i], jR[i]) / info->hxm; 360c4762a1bSJed Brown } else if (n == mf - 1) { 361c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod2(jL[i] / info->hxm, jR[i] / (info->hxm / 2.0 + info->hxf / 2.0)); 362c4762a1bSJed Brown } else if (n == mf) { 363c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod2(jL[i] / (info->hxm / 2.0 + info->hxf / 2.0), jR[i] / info->hxf); 364c4762a1bSJed Brown } else if (n == fm - 1) { 365c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod2(jL[i] / info->hxf, jR[i] / (info->hxf / 2.0 + info->hxm / 2.0)); 366c4762a1bSJed Brown } else if (n == fm) { 367c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod2(jL[i] / (info->hxf / 2.0 + info->hxm / 2.0), jR[i] / info->hxm); 368c4762a1bSJed Brown } else { 369c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = 0.5 * (jL[i] + jR[i]) / info->hxf; 370c4762a1bSJed Brown } 371c4762a1bSJed Brown } 372d71ae5a4SJacob Faibussowitsch void Limit3_Superbee(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, const PetscInt sm, const PetscInt mf, const PetscInt fm, const PetscInt ms, PetscInt n, PetscScalar *lmt) 373d71ae5a4SJacob Faibussowitsch { 374c4762a1bSJed Brown PetscInt i; 375c4762a1bSJed Brown if (n < sm - 1 || n > ms) { 376c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MaxMod2(MinMod2(jL[i], 2 * jR[i]), MinMod2(2 * jL[i], jR[i])) / info->hxs; 377c4762a1bSJed Brown } else if (n == sm - 1) { 378c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MaxMod2(MinMod2(jL[i] / info->hxs, 2 * jR[i] / (info->hxs / 2.0 + info->hxm / 2.0)), MinMod2(2 * jL[i] / info->hxs, jR[i] / (info->hxs / 2.0 + info->hxm / 2.0))); 379c4762a1bSJed Brown } else if (n == sm) { 380c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MaxMod2(MinMod2(jL[i] / (info->hxs / 2.0 + info->hxm / 2.0), 2 * jR[i] / info->hxm), MinMod2(2 * jL[i] / (info->hxs / 2.0 + info->hxm / 2.0), jR[i] / info->hxm)); 381c4762a1bSJed Brown } else if (n == ms - 1) { 382c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod2(MinMod2(jL[i] / info->hxm, 2 * jR[i] / (info->hxm / 2.0 + info->hxs / 2.0)), MinMod2(2 * jL[i] / info->hxm, jR[i] / (info->hxm / 2.0 + info->hxs / 2.0))); 383c4762a1bSJed Brown } else if (n == ms) { 384c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod2(MinMod2(jL[i] / (info->hxm / 2.0 + info->hxs / 2.0), 2 * jR[i] / info->hxs), MinMod2(2 * jL[i] / (info->hxm / 2.0 + info->hxs / 2.0), jR[i] / info->hxs)); 385c4762a1bSJed Brown } else if (n < mf - 1 || n > fm) { 386c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MaxMod2(MinMod2(jL[i], 2 * jR[i]), MinMod2(2 * jL[i], jR[i])) / info->hxm; 387c4762a1bSJed Brown } else if (n == mf - 1) { 388c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MaxMod2(MinMod2(jL[i] / info->hxm, 2 * jR[i] / (info->hxm / 2.0 + info->hxf / 2.0)), MinMod2(2 * jL[i] / info->hxm, jR[i] / (info->hxm / 2.0 + info->hxf / 2.0))); 389c4762a1bSJed Brown } else if (n == mf) { 390c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MaxMod2(MinMod2(jL[i] / (info->hxm / 2.0 + info->hxf / 2.0), 2 * jR[i] / info->hxf), MinMod2(2 * jL[i] / (info->hxm / 2.0 + info->hxf / 2.0), jR[i] / info->hxf)); 391c4762a1bSJed Brown } else if (n == fm - 1) { 392c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod2(MinMod2(jL[i] / info->hxf, 2 * jR[i] / (info->hxf / 2.0 + info->hxm / 2.0)), MinMod2(2 * jL[i] / info->hxf, jR[i] / (info->hxf / 2.0 + info->hxm / 2.0))); 393c4762a1bSJed Brown } else if (n == fm) { 394c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod2(MinMod2(jL[i] / (info->hxf / 2.0 + info->hxm / 2.0), 2 * jR[i] / info->hxm), MinMod2(2 * jL[i] / (info->hxf / 2.0 + info->hxm / 2.0), jR[i] / info->hxm)); 395c4762a1bSJed Brown } else { 396c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MaxMod2(MinMod2(jL[i], 2 * jR[i]), MinMod2(2 * jL[i], jR[i])) / info->hxf; 397c4762a1bSJed Brown } 398c4762a1bSJed Brown } 399d71ae5a4SJacob Faibussowitsch void Limit3_MC(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, const PetscInt sm, const PetscInt mf, const PetscInt fm, const PetscInt ms, PetscInt n, PetscScalar *lmt) 400d71ae5a4SJacob Faibussowitsch { 401c4762a1bSJed Brown PetscInt i; 402c4762a1bSJed Brown if (n < sm - 1 || n > ms) { 403c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i], 0.5 * (jL[i] + jR[i]), 2 * jR[i]) / info->hxs; 404c4762a1bSJed Brown } else if (n == sm - 1) { 405c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i] / info->hxs, 0.5 * (jL[i] / info->hxs + jR[i] / (info->hxs / 2.0 + info->hxm / 2.0)), 2 * jR[i] / (info->hxs / 2.0 + info->hxm / 2.0)); 406c4762a1bSJed Brown } else if (n == sm) { 407c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i] / (info->hxs / 2.0 + info->hxm / 2.0), 0.5 * (jL[i] / (info->hxs / 2.0 + info->hxm / 2.0) + jR[i] / info->hxm), 2 * jR[i] / info->hxm); 408c4762a1bSJed Brown } else if (n == ms - 1) { 409c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i] / info->hxm, 0.5 * (jL[i] / info->hxm + jR[i] / (info->hxm / 2.0 + info->hxs / 2.0)), 2 * jR[i] / (info->hxm / 2.0 + info->hxs / 2.0)); 410c4762a1bSJed Brown } else if (n == ms) { 411c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i] / (info->hxm / 2.0 + info->hxs / 2.0), 0.5 * (jL[i] / (info->hxm / 2.0 + info->hxs / 2.0) + jR[i] / info->hxs), 2 * jR[i] / info->hxs); 412c4762a1bSJed Brown } else if (n < mf - 1 || n > fm) { 413c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i], 0.5 * (jL[i] + jR[i]), 2 * jR[i]) / info->hxm; 414c4762a1bSJed Brown } else if (n == mf - 1) { 415c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i] / info->hxm, 0.5 * (jL[i] / info->hxm + jR[i] / (info->hxm / 2.0 + info->hxf / 2.0)), 2 * jR[i] / (info->hxm / 2.0 + info->hxf / 2.0)); 416c4762a1bSJed Brown } else if (n == mf) { 417c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i] / (info->hxm / 2.0 + info->hxf / 2.0), 0.5 * (jL[i] / (info->hxm / 2.0 + info->hxf / 2.0) + jR[i] / info->hxf), 2 * jR[i] / info->hxf); 418c4762a1bSJed Brown } else if (n == fm - 1) { 419c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i] / info->hxf, 0.5 * (jL[i] / info->hxf + jR[i] / (info->hxf / 2.0 + info->hxm / 2.0)), 2 * jR[i] / (info->hxf / 2.0 + info->hxm / 2.0)); 420c4762a1bSJed Brown } else if (n == fm) { 421c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i] / (info->hxf / 2.0 + info->hxm / 2.0), 0.5 * (jL[i] / (info->hxf / 2.0 + info->hxm / 2.0) + jR[i] / info->hxm), 2 * jR[i] / info->hxm); 422c4762a1bSJed Brown } else { 423c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i], 0.5 * (jL[i] + jR[i]), 2 * jR[i]) / info->hxf; 424c4762a1bSJed Brown } 425c4762a1bSJed Brown } 426d71ae5a4SJacob Faibussowitsch void Limit3_Koren3(LimitInfo info, const PetscScalar *jL, const PetscScalar *jR, const PetscInt sm, const PetscInt mf, const PetscInt fm, const PetscInt ms, PetscInt n, PetscScalar *lmt) 427d71ae5a4SJacob Faibussowitsch { /* Eq 11 of Cada-Torrilhon 2009 */ 428c4762a1bSJed Brown PetscInt i; 429c4762a1bSJed Brown if (n < sm - 1 || n > ms) { 430c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i], (jL[i] + 2 * jR[i]) / 3, 2 * jR[i]) / info->hxs; 431c4762a1bSJed Brown } else if (n == sm - 1) { 432c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i] / info->hxs, (jL[i] / info->hxs + 2 * jR[i] / (info->hxm / 2.0 + info->hxs / 2.0)) / 3, 2 * jR[i] / (info->hxm / 2.0 + info->hxs / 2.0)); 433c4762a1bSJed Brown } else if (n == sm) { 434c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i] / (info->hxs / 2.0 + info->hxm / 2.0), (jL[i] / (info->hxs / 2.0 + info->hxm / 2.0) + 2 * jR[i] / info->hxm) / 3, 2 * jR[i] / info->hxm); 435c4762a1bSJed Brown } else if (n == ms - 1) { 436c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i] / info->hxm, (jL[i] / info->hxm + 2 * jR[i] / (info->hxm / 2.0 + info->hxs / 2.0)) / 3, 2 * jR[i] / (info->hxm / 2.0 + info->hxs / 2.0)); 437c4762a1bSJed Brown } else if (n == ms) { 438c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i] / (info->hxm / 2.0 + info->hxs / 2.0), (jL[i] / (info->hxm / 2.0 + info->hxs / 2.0) + 2 * jR[i] / info->hxs) / 3, 2 * jR[i] / info->hxs); 439c4762a1bSJed Brown } else if (n < mf - 1 || n > fm) { 440c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i], (jL[i] + 2 * jR[i]) / 3, 2 * jR[i]) / info->hxm; 441c4762a1bSJed Brown } else if (n == mf - 1) { 442c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i] / info->hxm, (jL[i] / info->hxm + 2 * jR[i] / (info->hxm / 2.0 + info->hxf / 2.0)) / 3, 2 * jR[i] / (info->hxm / 2.0 + info->hxf / 2.0)); 443c4762a1bSJed Brown } else if (n == mf) { 444c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i] / (info->hxm / 2.0 + info->hxf / 2.0), (jL[i] / (info->hxm / 2.0 + info->hxf / 2.0) + 2 * jR[i] / info->hxf) / 3, 2 * jR[i] / info->hxf); 445c4762a1bSJed Brown } else if (n == fm - 1) { 446c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i] / info->hxf, (jL[i] / info->hxf + 2 * jR[i] / (info->hxf / 2.0 + info->hxm / 2.0)) / 3, 2 * jR[i] / (info->hxf / 2.0 + info->hxm / 2.0)); 447c4762a1bSJed Brown } else if (n == fm) { 448c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i] / (info->hxf / 2.0 + info->hxm / 2.0), (jL[i] / (info->hxf / 2.0 + info->hxm / 2.0) + 2 * jR[i] / info->hxm) / 3, 2 * jR[i] / info->hxm); 449c4762a1bSJed Brown } else { 450c4762a1bSJed Brown for (i = 0; i < info->m; i++) lmt[i] = MinMod3(2 * jL[i], (jL[i] + 2 * jR[i]) / 3, 2 * jR[i]) / info->hxs; 451c4762a1bSJed Brown } 452c4762a1bSJed Brown } 45384ff8c8bSHong Zhang 454d71ae5a4SJacob Faibussowitsch PetscErrorCode RiemannListAdd(PetscFunctionList *flist, const char *name, RiemannFunction rsolve) 455d71ae5a4SJacob Faibussowitsch { 456c4762a1bSJed Brown PetscFunctionBeginUser; 4579566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(flist, name, rsolve)); 4583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 459c4762a1bSJed Brown } 460c4762a1bSJed Brown 461d71ae5a4SJacob Faibussowitsch PetscErrorCode RiemannListFind(PetscFunctionList flist, const char *name, RiemannFunction *rsolve) 462d71ae5a4SJacob Faibussowitsch { 463c4762a1bSJed Brown PetscFunctionBeginUser; 4649566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(flist, name, rsolve)); 4653c633725SBarry Smith PetscCheck(*rsolve, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Riemann solver \"%s\" could not be found", name); 4663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 467c4762a1bSJed Brown } 468c4762a1bSJed Brown 469d71ae5a4SJacob Faibussowitsch PetscErrorCode ReconstructListAdd(PetscFunctionList *flist, const char *name, ReconstructFunction r) 470d71ae5a4SJacob Faibussowitsch { 471c4762a1bSJed Brown PetscFunctionBeginUser; 4729566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(flist, name, r)); 4733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 474c4762a1bSJed Brown } 475c4762a1bSJed Brown 476d71ae5a4SJacob Faibussowitsch PetscErrorCode ReconstructListFind(PetscFunctionList flist, const char *name, ReconstructFunction *r) 477d71ae5a4SJacob Faibussowitsch { 478c4762a1bSJed Brown PetscFunctionBeginUser; 4799566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(flist, name, r)); 4803c633725SBarry Smith PetscCheck(*r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Reconstruction \"%s\" could not be found", name); 4813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 482c4762a1bSJed Brown } 483c4762a1bSJed Brown 484d71ae5a4SJacob Faibussowitsch PetscErrorCode RiemannListAdd_2WaySplit(PetscFunctionList *flist, const char *name, RiemannFunction_2WaySplit rsolve) 485d71ae5a4SJacob Faibussowitsch { 48684ff8c8bSHong Zhang PetscFunctionBeginUser; 4879566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(flist, name, rsolve)); 4883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 48984ff8c8bSHong Zhang } 49084ff8c8bSHong Zhang 491d71ae5a4SJacob Faibussowitsch PetscErrorCode RiemannListFind_2WaySplit(PetscFunctionList flist, const char *name, RiemannFunction_2WaySplit *rsolve) 492d71ae5a4SJacob Faibussowitsch { 49384ff8c8bSHong Zhang PetscFunctionBeginUser; 4949566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(flist, name, rsolve)); 4953c633725SBarry Smith PetscCheck(*rsolve, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Riemann solver \"%s\" could not be found", name); 4963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 49784ff8c8bSHong Zhang } 49884ff8c8bSHong Zhang 499d71ae5a4SJacob Faibussowitsch PetscErrorCode ReconstructListAdd_2WaySplit(PetscFunctionList *flist, const char *name, ReconstructFunction_2WaySplit r) 500d71ae5a4SJacob Faibussowitsch { 50184ff8c8bSHong Zhang PetscFunctionBeginUser; 5029566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(flist, name, r)); 5033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 50484ff8c8bSHong Zhang } 50584ff8c8bSHong Zhang 506d71ae5a4SJacob Faibussowitsch PetscErrorCode ReconstructListFind_2WaySplit(PetscFunctionList flist, const char *name, ReconstructFunction_2WaySplit *r) 507d71ae5a4SJacob Faibussowitsch { 50884ff8c8bSHong Zhang PetscFunctionBeginUser; 5099566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(flist, name, r)); 5103c633725SBarry Smith PetscCheck(*r, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "Reconstruction \"%s\" could not be found", name); 5113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 51284ff8c8bSHong Zhang } 51384ff8c8bSHong Zhang 51484ff8c8bSHong Zhang /* --------------------------------- Physics ------- */ 515d71ae5a4SJacob Faibussowitsch PetscErrorCode PhysicsDestroy_SimpleFree(void *vctx) 516d71ae5a4SJacob Faibussowitsch { 517c4762a1bSJed Brown PetscFunctionBeginUser; 5189566063dSJacob Faibussowitsch PetscCall(PetscFree(vctx)); 5193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 520c4762a1bSJed Brown } 521c4762a1bSJed Brown 52284ff8c8bSHong Zhang /* --------------------------------- Finite Volume Solver --------------- */ 523d71ae5a4SJacob Faibussowitsch PetscErrorCode FVRHSFunction(TS ts, PetscReal time, Vec X, Vec F, void *vctx) 524d71ae5a4SJacob Faibussowitsch { 525c4762a1bSJed Brown FVCtx *ctx = (FVCtx *)vctx; 526c4762a1bSJed Brown PetscInt i, j, k, Mx, dof, xs, xm; 527c4762a1bSJed Brown PetscReal hx, cfl_idt = 0; 528c4762a1bSJed Brown PetscScalar *x, *f, *slope; 529c4762a1bSJed Brown Vec Xloc; 530c4762a1bSJed Brown DM da; 531c4762a1bSJed Brown 532c4762a1bSJed Brown PetscFunctionBeginUser; 5339566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 5349566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &Xloc)); /* Xloc contains ghost points */ 5359566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); /* Mx is the number of center points */ 536c4762a1bSJed Brown hx = (ctx->xmax - ctx->xmin) / Mx; 5379566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); /* X is solution vector which does not contain ghost points */ 5389566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc)); 539*dd8e379bSPierre Jolivet PetscCall(VecZeroEntries(F)); /* F is the right-hand side function corresponds to center points */ 5409566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, Xloc, &x)); 5419566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, F, &f)); 5429566063dSJacob Faibussowitsch PetscCall(DMDAGetArray(da, PETSC_TRUE, &slope)); /* contains ghost points */ 5439566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); 544c4762a1bSJed Brown 545c4762a1bSJed Brown if (ctx->bctype == FVBC_OUTFLOW) { 546c4762a1bSJed Brown for (i = xs - 2; i < 0; i++) { 547c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[j]; 548c4762a1bSJed Brown } 549c4762a1bSJed Brown for (i = Mx; i < xs + xm + 2; i++) { 550c4762a1bSJed Brown for (j = 0; j < dof; j++) x[i * dof + j] = x[(xs + xm - 1) * dof + j]; 551c4762a1bSJed Brown } 552c4762a1bSJed Brown } 55384ff8c8bSHong Zhang 554c4762a1bSJed Brown for (i = xs - 1; i < xs + xm + 1; i++) { 555c4762a1bSJed Brown struct _LimitInfo info; 556c4762a1bSJed Brown PetscScalar *cjmpL, *cjmpR; 557c4762a1bSJed Brown /* Determine the right eigenvectors R, where A = R \Lambda R^{-1} */ 5589566063dSJacob Faibussowitsch PetscCall((*ctx->physics.characteristic)(ctx->physics.user, dof, &x[i * dof], ctx->R, ctx->Rinv, ctx->speeds, ctx->xmin + hx * i)); 559c4762a1bSJed Brown /* Evaluate jumps across interfaces (i-1, i) and (i, i+1), put in characteristic basis */ 5609566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(ctx->cjmpLR, 2 * dof)); 561c4762a1bSJed Brown cjmpL = &ctx->cjmpLR[0]; 562c4762a1bSJed Brown cjmpR = &ctx->cjmpLR[dof]; 563c4762a1bSJed Brown for (j = 0; j < dof; j++) { 564c4762a1bSJed Brown PetscScalar jmpL, jmpR; 565c4762a1bSJed Brown jmpL = x[(i + 0) * dof + j] - x[(i - 1) * dof + j]; 566c4762a1bSJed Brown jmpR = x[(i + 1) * dof + j] - x[(i + 0) * dof + j]; 567c4762a1bSJed Brown for (k = 0; k < dof; k++) { 568c4762a1bSJed Brown cjmpL[k] += ctx->Rinv[k + j * dof] * jmpL; 569c4762a1bSJed Brown cjmpR[k] += ctx->Rinv[k + j * dof] * jmpR; 570c4762a1bSJed Brown } 571c4762a1bSJed Brown } 572c4762a1bSJed Brown /* Apply limiter to the left and right characteristic jumps */ 573c4762a1bSJed Brown info.m = dof; 574c4762a1bSJed Brown info.hx = hx; 575c4762a1bSJed Brown (*ctx->limit)(&info, cjmpL, cjmpR, ctx->cslope); 576c4762a1bSJed Brown for (j = 0; j < dof; j++) ctx->cslope[j] /= hx; /* rescale to a slope */ 577c4762a1bSJed Brown for (j = 0; j < dof; j++) { 578c4762a1bSJed Brown PetscScalar tmp = 0; 579c4762a1bSJed Brown for (k = 0; k < dof; k++) tmp += ctx->R[j + k * dof] * ctx->cslope[k]; 580c4762a1bSJed Brown slope[i * dof + j] = tmp; 581c4762a1bSJed Brown } 582c4762a1bSJed Brown } 583c4762a1bSJed Brown 584c4762a1bSJed Brown for (i = xs; i < xs + xm + 1; i++) { 585c4762a1bSJed Brown PetscReal maxspeed; 586c4762a1bSJed Brown PetscScalar *uL, *uR; 587c4762a1bSJed Brown uL = &ctx->uLR[0]; 588c4762a1bSJed Brown uR = &ctx->uLR[dof]; 589c4762a1bSJed Brown for (j = 0; j < dof; j++) { 590c4762a1bSJed Brown uL[j] = x[(i - 1) * dof + j] + slope[(i - 1) * dof + j] * hx / 2; 591c4762a1bSJed Brown uR[j] = x[(i - 0) * dof + j] - slope[(i - 0) * dof + j] * hx / 2; 592c4762a1bSJed Brown } 5939566063dSJacob Faibussowitsch PetscCall((*ctx->physics.riemann)(ctx->physics.user, dof, uL, uR, ctx->flux, &maxspeed, ctx->xmin + hx * i, ctx->xmin, ctx->xmax)); 594c4762a1bSJed Brown cfl_idt = PetscMax(cfl_idt, PetscAbsScalar(maxspeed / hx)); /* Max allowable value of 1/Delta t */ 595c4762a1bSJed Brown if (i > xs) { 596c4762a1bSJed Brown for (j = 0; j < dof; j++) f[(i - 1) * dof + j] -= ctx->flux[j] / hx; 597c4762a1bSJed Brown } 598c4762a1bSJed Brown if (i < xs + xm) { 599c4762a1bSJed Brown for (j = 0; j < dof; j++) f[i * dof + j] += ctx->flux[j] / hx; 600c4762a1bSJed Brown } 601c4762a1bSJed Brown } 6029566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, Xloc, &x)); 6039566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, F, &f)); 6049566063dSJacob Faibussowitsch PetscCall(DMDARestoreArray(da, PETSC_TRUE, &slope)); 6059566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &Xloc)); 606712fec58SPierre Jolivet PetscCall(MPIU_Allreduce(&cfl_idt, &ctx->cfl_idt, 1, MPIU_REAL, MPIU_MAX, PetscObjectComm((PetscObject)da))); 607c4762a1bSJed Brown if (0) { 608c4762a1bSJed Brown /* We need to a way to inform the TS of a CFL constraint, this is a debugging fragment */ 609c4762a1bSJed Brown PetscReal dt, tnow; 6109566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &dt)); 6119566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &tnow)); 612c4762a1bSJed Brown if (dt > 0.5 / ctx->cfl_idt) { 61384ff8c8bSHong Zhang if (1) { 6149566063dSJacob Faibussowitsch PetscCall(PetscPrintf(ctx->comm, "Stability constraint exceeded at t=%g, dt %g > %g\n", (double)tnow, (double)dt, (double)(0.5 / ctx->cfl_idt))); 61598921bdaSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Stability constraint exceeded, %g > %g", (double)dt, (double)(ctx->cfl / ctx->cfl_idt)); 616c4762a1bSJed Brown } 617c4762a1bSJed Brown } 6183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 619c4762a1bSJed Brown } 620c4762a1bSJed Brown 621d71ae5a4SJacob Faibussowitsch PetscErrorCode FVSample(FVCtx *ctx, DM da, PetscReal time, Vec U) 622d71ae5a4SJacob Faibussowitsch { 623c4762a1bSJed Brown PetscScalar *u, *uj; 624c4762a1bSJed Brown PetscInt i, j, k, dof, xs, xm, Mx; 625c4762a1bSJed Brown 626c4762a1bSJed Brown PetscFunctionBeginUser; 6273c633725SBarry Smith PetscCheck(ctx->physics.sample, PETSC_COMM_SELF, PETSC_ERR_SUP, "Physics has not provided a sampling function"); 6289566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); 6299566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); 6309566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, U, &u)); 6319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dof, &uj)); 632c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 633c4762a1bSJed Brown const PetscReal h = (ctx->xmax - ctx->xmin) / Mx, xi = ctx->xmin + h / 2 + i * h; 634c4762a1bSJed Brown const PetscInt N = 200; 635c4762a1bSJed Brown /* Integrate over cell i using trapezoid rule with N points. */ 636c4762a1bSJed Brown for (k = 0; k < dof; k++) u[i * dof + k] = 0; 637c4762a1bSJed Brown for (j = 0; j < N + 1; j++) { 638c4762a1bSJed Brown PetscScalar xj = xi + h * (j - N / 2) / (PetscReal)N; 6399566063dSJacob Faibussowitsch PetscCall((*ctx->physics.sample)(ctx->physics.user, ctx->initial, ctx->bctype, ctx->xmin, ctx->xmax, time, xj, uj)); 640c4762a1bSJed Brown for (k = 0; k < dof; k++) u[i * dof + k] += ((j == 0 || j == N) ? 0.5 : 1.0) * uj[k] / N; 641c4762a1bSJed Brown } 642c4762a1bSJed Brown } 6439566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, U, &u)); 6449566063dSJacob Faibussowitsch PetscCall(PetscFree(uj)); 6453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 646c4762a1bSJed Brown } 647c4762a1bSJed Brown 648d71ae5a4SJacob Faibussowitsch PetscErrorCode SolutionStatsView(DM da, Vec X, PetscViewer viewer) 649d71ae5a4SJacob Faibussowitsch { 650c4762a1bSJed Brown PetscReal xmin, xmax; 651c4762a1bSJed Brown PetscScalar sum, tvsum, tvgsum; 652c4762a1bSJed Brown const PetscScalar *x; 653c4762a1bSJed Brown PetscInt imin, imax, Mx, i, j, xs, xm, dof; 654c4762a1bSJed Brown Vec Xloc; 655c4762a1bSJed Brown PetscBool iascii; 656c4762a1bSJed Brown 657c4762a1bSJed Brown PetscFunctionBeginUser; 6589566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 659c4762a1bSJed Brown if (iascii) { 660c4762a1bSJed Brown /* PETSc lacks a function to compute total variation norm (difficult in multiple dimensions), we do it here */ 6619566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &Xloc)); 6629566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc)); 6639566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc)); 6649566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, Xloc, (void *)&x)); 6659566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, 0, 0, &xm, 0, 0)); 6669566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, 0, &Mx, 0, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0)); 667c4762a1bSJed Brown tvsum = 0; 668c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 669c4762a1bSJed Brown for (j = 0; j < dof; j++) tvsum += PetscAbsScalar(x[i * dof + j] - x[(i - 1) * dof + j]); 670c4762a1bSJed Brown } 671712fec58SPierre Jolivet PetscCall(MPIU_Allreduce(&tvsum, &tvgsum, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)da))); 6729566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, Xloc, (void *)&x)); 6739566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &Xloc)); 6749566063dSJacob Faibussowitsch PetscCall(VecMin(X, &imin, &xmin)); 6759566063dSJacob Faibussowitsch PetscCall(VecMax(X, &imax, &xmax)); 6769566063dSJacob Faibussowitsch PetscCall(VecSum(X, &sum)); 67763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Solution range [%8.5f,%8.5f] with minimum at %" PetscInt_FMT ", mean %8.5f, ||x||_TV %8.5f\n", (double)xmin, (double)xmax, imin, (double)(sum / Mx), (double)(tvgsum / Mx))); 678c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Viewer type not supported"); 6793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 680c4762a1bSJed Brown } 681