17d0a6c19SBarry Smith 2e5c89e4eSSatish Balay /* 3e5c89e4eSSatish Balay This file contains routines for sorting doubles. Values are sorted in place. 4e5c89e4eSSatish Balay These are provided because the general sort routines incur a great deal 5a5b23f4aSJose E. Roman of overhead in calling the comparison routines. 6e5c89e4eSSatish Balay 7e5c89e4eSSatish Balay */ 8c6db04a5SJed Brown #include <petscsys.h> /*I "petscsys.h" I*/ 939f41f7fSStefano Zampini #include <petsc/private/petscimpl.h> 10e5c89e4eSSatish Balay 119371c9d4SSatish Balay #define SWAP(a, b, t) \ 129371c9d4SSatish Balay { \ 139371c9d4SSatish Balay t = a; \ 149371c9d4SSatish Balay a = b; \ 159371c9d4SSatish Balay b = t; \ 169371c9d4SSatish Balay } 17e5c89e4eSSatish Balay 186a7c706bSVaclav Hapla /*@ 19811af0c4SBarry Smith PetscSortedReal - Determines whether the array of `PetscReal` is sorted. 206a7c706bSVaclav Hapla 216a7c706bSVaclav Hapla Not Collective 226a7c706bSVaclav Hapla 236a7c706bSVaclav Hapla Input Parameters: 246a7c706bSVaclav Hapla + n - number of values 256a7c706bSVaclav Hapla - X - array of integers 266a7c706bSVaclav Hapla 27*2fe279fdSBarry Smith Output Parameter: 286a7c706bSVaclav Hapla . sorted - flag whether the array is sorted 296a7c706bSVaclav Hapla 306a7c706bSVaclav Hapla Level: intermediate 316a7c706bSVaclav Hapla 32db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortedInt()`, `PetscSortedMPIInt()` 336a7c706bSVaclav Hapla @*/ 34d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortedReal(PetscInt n, const PetscReal X[], PetscBool *sorted) 35d71ae5a4SJacob Faibussowitsch { 366a7c706bSVaclav Hapla PetscFunctionBegin; 376a7c706bSVaclav Hapla PetscSorted(n, X, *sorted); 383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 396a7c706bSVaclav Hapla } 406a7c706bSVaclav Hapla 41e5c89e4eSSatish Balay /* A simple version of quicksort; taken from Kernighan and Ritchie, page 87 */ 42d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSortReal_Private(PetscReal *v, PetscInt right) 43d71ae5a4SJacob Faibussowitsch { 44e5c89e4eSSatish Balay PetscInt i, last; 45e5c89e4eSSatish Balay PetscReal vl, tmp; 46e5c89e4eSSatish Balay 47e5c89e4eSSatish Balay PetscFunctionBegin; 48e5c89e4eSSatish Balay if (right <= 1) { 49e5c89e4eSSatish Balay if (right == 1) { 50e5c89e4eSSatish Balay if (v[0] > v[1]) SWAP(v[0], v[1], tmp); 51e5c89e4eSSatish Balay } 523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 53e5c89e4eSSatish Balay } 54e5c89e4eSSatish Balay SWAP(v[0], v[right / 2], tmp); 55e5c89e4eSSatish Balay vl = v[0]; 56e5c89e4eSSatish Balay last = 0; 57e5c89e4eSSatish Balay for (i = 1; i <= right; i++) { 589371c9d4SSatish Balay if (v[i] < vl) { 599371c9d4SSatish Balay last++; 609371c9d4SSatish Balay SWAP(v[last], v[i], tmp); 619371c9d4SSatish Balay } 62e5c89e4eSSatish Balay } 63e5c89e4eSSatish Balay SWAP(v[0], v[last], tmp); 643ba16761SJacob Faibussowitsch PetscCall(PetscSortReal_Private(v, last - 1)); 653ba16761SJacob Faibussowitsch PetscCall(PetscSortReal_Private(v + last + 1, right - (last + 1))); 663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 67e5c89e4eSSatish Balay } 68e5c89e4eSSatish Balay 69e5c89e4eSSatish Balay /*@ 70811af0c4SBarry Smith PetscSortReal - Sorts an array of `PetscReal` in place in increasing order. 71e5c89e4eSSatish Balay 72e5c89e4eSSatish Balay Not Collective 73e5c89e4eSSatish Balay 74e5c89e4eSSatish Balay Input Parameters: 75e5c89e4eSSatish Balay + n - number of values 76e5c89e4eSSatish Balay - v - array of doubles 77e5c89e4eSSatish Balay 78667f096bSBarry Smith Level: intermediate 79667f096bSBarry Smith 80811af0c4SBarry Smith Note: 81811af0c4SBarry Smith This function serves as an alternative to `PetscRealSortSemiOrdered()`, and may perform faster especially if the array 82a5b23f4aSJose E. Roman is completely random. There are exceptions to this and so it is __highly__ recommended that the user benchmark their 83676f2a66SJacob Faibussowitsch code to see which routine is fastest. 84676f2a66SJacob Faibussowitsch 85db781477SPatrick Sanan .seealso: `PetscRealSortSemiOrdered()`, `PetscSortInt()`, `PetscSortRealWithPermutation()`, `PetscSortRealWithArrayInt()` 86e5c89e4eSSatish Balay @*/ 87d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortReal(PetscInt n, PetscReal v[]) 88d71ae5a4SJacob Faibussowitsch { 89e5c89e4eSSatish Balay PetscInt j, k; 90e5c89e4eSSatish Balay PetscReal tmp, vk; 91e5c89e4eSSatish Balay 92e5c89e4eSSatish Balay PetscFunctionBegin; 93dadcf809SJacob Faibussowitsch PetscValidRealPointer(v, 2); 94e5c89e4eSSatish Balay if (n < 8) { 95e5c89e4eSSatish Balay for (k = 0; k < n; k++) { 96e5c89e4eSSatish Balay vk = v[k]; 97e5c89e4eSSatish Balay for (j = k + 1; j < n; j++) { 98e5c89e4eSSatish Balay if (vk > v[j]) { 99e5c89e4eSSatish Balay SWAP(v[k], v[j], tmp); 100e5c89e4eSSatish Balay vk = v[k]; 101e5c89e4eSSatish Balay } 102e5c89e4eSSatish Balay } 103e5c89e4eSSatish Balay } 1043ba16761SJacob Faibussowitsch } else PetscCall(PetscSortReal_Private(v, n - 1)); 1053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 106e5c89e4eSSatish Balay } 107e5c89e4eSSatish Balay 1089371c9d4SSatish Balay #define SWAP2ri(a, b, c, d, rt, it) \ 1099371c9d4SSatish Balay { \ 1109371c9d4SSatish Balay rt = a; \ 1119371c9d4SSatish Balay a = b; \ 1129371c9d4SSatish Balay b = rt; \ 1139371c9d4SSatish Balay it = c; \ 1149371c9d4SSatish Balay c = d; \ 1159371c9d4SSatish Balay d = it; \ 1169371c9d4SSatish Balay } 11739f41f7fSStefano Zampini 11839f41f7fSStefano Zampini /* modified from PetscSortIntWithArray_Private */ 119d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSortRealWithArrayInt_Private(PetscReal *v, PetscInt *V, PetscInt right) 120d71ae5a4SJacob Faibussowitsch { 12139f41f7fSStefano Zampini PetscInt i, last, itmp; 12239f41f7fSStefano Zampini PetscReal rvl, rtmp; 12339f41f7fSStefano Zampini 12439f41f7fSStefano Zampini PetscFunctionBegin; 12539f41f7fSStefano Zampini if (right <= 1) { 12639f41f7fSStefano Zampini if (right == 1) { 12739f41f7fSStefano Zampini if (v[0] > v[1]) SWAP2ri(v[0], v[1], V[0], V[1], rtmp, itmp); 12839f41f7fSStefano Zampini } 1293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13039f41f7fSStefano Zampini } 13139f41f7fSStefano Zampini SWAP2ri(v[0], v[right / 2], V[0], V[right / 2], rtmp, itmp); 13239f41f7fSStefano Zampini rvl = v[0]; 13339f41f7fSStefano Zampini last = 0; 13439f41f7fSStefano Zampini for (i = 1; i <= right; i++) { 1359371c9d4SSatish Balay if (v[i] < rvl) { 1369371c9d4SSatish Balay last++; 1379371c9d4SSatish Balay SWAP2ri(v[last], v[i], V[last], V[i], rtmp, itmp); 1389371c9d4SSatish Balay } 13939f41f7fSStefano Zampini } 14039f41f7fSStefano Zampini SWAP2ri(v[0], v[last], V[0], V[last], rtmp, itmp); 1419566063dSJacob Faibussowitsch PetscCall(PetscSortRealWithArrayInt_Private(v, V, last - 1)); 1429566063dSJacob Faibussowitsch PetscCall(PetscSortRealWithArrayInt_Private(v + last + 1, V + last + 1, right - (last + 1))); 1433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14439f41f7fSStefano Zampini } 14539f41f7fSStefano Zampini /*@ 146811af0c4SBarry Smith PetscSortRealWithArrayInt - Sorts an array of `PetscReal` in place in increasing order; 147811af0c4SBarry Smith changes a second `PetscInt` array to match the sorted first array. 14839f41f7fSStefano Zampini 14939f41f7fSStefano Zampini Not Collective 15039f41f7fSStefano Zampini 15139f41f7fSStefano Zampini Input Parameters: 15239f41f7fSStefano Zampini + n - number of values 15339f41f7fSStefano Zampini . i - array of integers 15439f41f7fSStefano Zampini - I - second array of integers 15539f41f7fSStefano Zampini 15639f41f7fSStefano Zampini Level: intermediate 15739f41f7fSStefano Zampini 158db781477SPatrick Sanan .seealso: `PetscSortReal()` 15939f41f7fSStefano Zampini @*/ 160d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortRealWithArrayInt(PetscInt n, PetscReal r[], PetscInt Ii[]) 161d71ae5a4SJacob Faibussowitsch { 16239f41f7fSStefano Zampini PetscInt j, k, itmp; 16339f41f7fSStefano Zampini PetscReal rk, rtmp; 16439f41f7fSStefano Zampini 16539f41f7fSStefano Zampini PetscFunctionBegin; 166dadcf809SJacob Faibussowitsch PetscValidRealPointer(r, 2); 167dadcf809SJacob Faibussowitsch PetscValidIntPointer(Ii, 3); 16839f41f7fSStefano Zampini if (n < 8) { 16939f41f7fSStefano Zampini for (k = 0; k < n; k++) { 17039f41f7fSStefano Zampini rk = r[k]; 17139f41f7fSStefano Zampini for (j = k + 1; j < n; j++) { 17239f41f7fSStefano Zampini if (rk > r[j]) { 17339f41f7fSStefano Zampini SWAP2ri(r[k], r[j], Ii[k], Ii[j], rtmp, itmp); 17439f41f7fSStefano Zampini rk = r[k]; 17539f41f7fSStefano Zampini } 17639f41f7fSStefano Zampini } 17739f41f7fSStefano Zampini } 17839f41f7fSStefano Zampini } else { 1799566063dSJacob Faibussowitsch PetscCall(PetscSortRealWithArrayInt_Private(r, Ii, n - 1)); 18039f41f7fSStefano Zampini } 1813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 18239f41f7fSStefano Zampini } 18339f41f7fSStefano Zampini 18439f41f7fSStefano Zampini /*@ 185811af0c4SBarry Smith PetscFindReal - Finds a PetscReal` in a sorted array of `PetscReal`s 18639f41f7fSStefano Zampini 18739f41f7fSStefano Zampini Not Collective 18839f41f7fSStefano Zampini 18939f41f7fSStefano Zampini Input Parameters: 19039f41f7fSStefano Zampini + key - the value to locate 19139f41f7fSStefano Zampini . n - number of values in the array 19239f41f7fSStefano Zampini . ii - array of values 19339f41f7fSStefano Zampini - eps - tolerance used to compare 19439f41f7fSStefano Zampini 19539f41f7fSStefano Zampini Output Parameter: 19639f41f7fSStefano Zampini . loc - the location if found, otherwise -(slot+1) where slot is the place the value would go 19739f41f7fSStefano Zampini 19839f41f7fSStefano Zampini Level: intermediate 19939f41f7fSStefano Zampini 200db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortRealWithArrayInt()` 20139f41f7fSStefano Zampini @*/ 202d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscFindReal(PetscReal key, PetscInt n, const PetscReal t[], PetscReal eps, PetscInt *loc) 203d71ae5a4SJacob Faibussowitsch { 20439f41f7fSStefano Zampini PetscInt lo = 0, hi = n; 20539f41f7fSStefano Zampini 20639f41f7fSStefano Zampini PetscFunctionBegin; 207dadcf809SJacob Faibussowitsch PetscValidIntPointer(loc, 5); 2089371c9d4SSatish Balay if (!n) { 2099371c9d4SSatish Balay *loc = -1; 2103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2119371c9d4SSatish Balay } 212dadcf809SJacob Faibussowitsch PetscValidRealPointer(t, 3); 2136a7c706bSVaclav Hapla PetscCheckSorted(n, t); 21439f41f7fSStefano Zampini while (hi - lo > 1) { 21539f41f7fSStefano Zampini PetscInt mid = lo + (hi - lo) / 2; 21639f41f7fSStefano Zampini if (key < t[mid]) hi = mid; 21739f41f7fSStefano Zampini else lo = mid; 21839f41f7fSStefano Zampini } 21939f41f7fSStefano Zampini *loc = (PetscAbsReal(key - t[lo]) < eps) ? lo : -(lo + (key > t[lo]) + 1); 2203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 22139f41f7fSStefano Zampini } 222745b41b2SMatthew G. Knepley 223745b41b2SMatthew G. Knepley /*@ 224811af0c4SBarry Smith PetscSortRemoveDupsReal - Sorts an array of `PetscReal` in place in increasing order and removes all duplicate entries 225745b41b2SMatthew G. Knepley 226745b41b2SMatthew G. Knepley Not Collective 227745b41b2SMatthew G. Knepley 228745b41b2SMatthew G. Knepley Input Parameters: 229745b41b2SMatthew G. Knepley + n - number of values 230745b41b2SMatthew G. Knepley - v - array of doubles 231745b41b2SMatthew G. Knepley 232745b41b2SMatthew G. Knepley Output Parameter: 233745b41b2SMatthew G. Knepley . n - number of non-redundant values 234745b41b2SMatthew G. Knepley 235745b41b2SMatthew G. Knepley Level: intermediate 236745b41b2SMatthew G. Knepley 237db781477SPatrick Sanan .seealso: `PetscSortReal()`, `PetscSortRemoveDupsInt()` 238745b41b2SMatthew G. Knepley @*/ 239d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortRemoveDupsReal(PetscInt *n, PetscReal v[]) 240d71ae5a4SJacob Faibussowitsch { 241745b41b2SMatthew G. Knepley PetscInt i, s = 0, N = *n, b = 0; 242745b41b2SMatthew G. Knepley 243745b41b2SMatthew G. Knepley PetscFunctionBegin; 2449566063dSJacob Faibussowitsch PetscCall(PetscSortReal(N, v)); 245745b41b2SMatthew G. Knepley for (i = 0; i < N - 1; i++) { 246745b41b2SMatthew G. Knepley if (v[b + s + 1] != v[b]) { 2479371c9d4SSatish Balay v[b + 1] = v[b + s + 1]; 2489371c9d4SSatish Balay b++; 249745b41b2SMatthew G. Knepley } else s++; 250745b41b2SMatthew G. Knepley } 251745b41b2SMatthew G. Knepley *n = N - s; 2523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 253745b41b2SMatthew G. Knepley } 254745b41b2SMatthew G. Knepley 255d97c5584SHong Zhang /*@ 256811af0c4SBarry Smith PetscSortSplit - Quick-sort split of an array of `PetscScalar`s in place. 257d97c5584SHong Zhang 258d97c5584SHong Zhang Not Collective 259d97c5584SHong Zhang 260d97c5584SHong Zhang Input Parameters: 261d97c5584SHong Zhang + ncut - splitig index 2626b867d5aSJose E. Roman - n - number of values to sort 263d97c5584SHong Zhang 2646b867d5aSJose E. Roman Input/Output Parameters: 2656b867d5aSJose E. Roman + a - array of values, on output the values are permuted such that its elements satisfy: 266d97c5584SHong Zhang abs(a[i]) >= abs(a[ncut-1]) for i < ncut and 267d97c5584SHong Zhang abs(a[i]) <= abs(a[ncut-1]) for i >= ncut 2686b867d5aSJose E. Roman - idx - index for array a, on output permuted accordingly 269d97c5584SHong Zhang 270d97c5584SHong Zhang Level: intermediate 271d97c5584SHong Zhang 272db781477SPatrick Sanan .seealso: `PetscSortInt()`, `PetscSortRealWithPermutation()` 273d97c5584SHong Zhang @*/ 274d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortSplit(PetscInt ncut, PetscInt n, PetscScalar a[], PetscInt idx[]) 275d71ae5a4SJacob Faibussowitsch { 276d97c5584SHong Zhang PetscInt i, mid, last, itmp, j, first; 277d97c5584SHong Zhang PetscScalar d, tmp; 278d97c5584SHong Zhang PetscReal abskey; 279d97c5584SHong Zhang 280d97c5584SHong Zhang PetscFunctionBegin; 281d97c5584SHong Zhang first = 0; 282d97c5584SHong Zhang last = n - 1; 2833ba16761SJacob Faibussowitsch if (ncut < first || ncut > last) PetscFunctionReturn(PETSC_SUCCESS); 284d97c5584SHong Zhang 285d97c5584SHong Zhang while (1) { 286d97c5584SHong Zhang mid = first; 2872a274a78SSatish Balay d = a[mid]; 2882a274a78SSatish Balay abskey = PetscAbsScalar(d); 289d97c5584SHong Zhang i = last; 290d97c5584SHong Zhang for (j = first + 1; j <= i; ++j) { 2912a274a78SSatish Balay d = a[j]; 2922a274a78SSatish Balay if (PetscAbsScalar(d) >= abskey) { 293d97c5584SHong Zhang ++mid; 294d97c5584SHong Zhang /* interchange */ 2959371c9d4SSatish Balay tmp = a[mid]; 2969371c9d4SSatish Balay itmp = idx[mid]; 2979371c9d4SSatish Balay a[mid] = a[j]; 2989371c9d4SSatish Balay idx[mid] = idx[j]; 2999371c9d4SSatish Balay a[j] = tmp; 3009371c9d4SSatish Balay idx[j] = itmp; 301d97c5584SHong Zhang } 302d97c5584SHong Zhang } 303d97c5584SHong Zhang 304d97c5584SHong Zhang /* interchange */ 3059371c9d4SSatish Balay tmp = a[mid]; 3069371c9d4SSatish Balay itmp = idx[mid]; 3079371c9d4SSatish Balay a[mid] = a[first]; 3089371c9d4SSatish Balay idx[mid] = idx[first]; 3099371c9d4SSatish Balay a[first] = tmp; 3109371c9d4SSatish Balay idx[first] = itmp; 311d97c5584SHong Zhang 312d97c5584SHong Zhang /* test for while loop */ 313a297a907SKarl Rupp if (mid == ncut) break; 314a297a907SKarl Rupp else if (mid > ncut) last = mid - 1; 315a297a907SKarl Rupp else first = mid + 1; 316d97c5584SHong Zhang } 3173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 318d97c5584SHong Zhang } 319021822bcSHong Zhang 320021822bcSHong Zhang /*@ 321811af0c4SBarry Smith PetscSortSplitReal - Quick-sort split of an array of `PetscReal`s in place. 322021822bcSHong Zhang 323021822bcSHong Zhang Not Collective 324021822bcSHong Zhang 325021822bcSHong Zhang Input Parameters: 326021822bcSHong Zhang + ncut - splitig index 3276b867d5aSJose E. Roman - n - number of values to sort 328021822bcSHong Zhang 3296b867d5aSJose E. Roman Input/Output Parameters: 3306b867d5aSJose E. Roman + a - array of values, on output the values are permuted such that its elements satisfy: 331021822bcSHong Zhang abs(a[i]) >= abs(a[ncut-1]) for i < ncut and 332021822bcSHong Zhang abs(a[i]) <= abs(a[ncut-1]) for i >= ncut 3336b867d5aSJose E. Roman - idx - index for array a, on output permuted accordingly 334021822bcSHong Zhang 335021822bcSHong Zhang Level: intermediate 336021822bcSHong Zhang 337db781477SPatrick Sanan .seealso: `PetscSortInt()`, `PetscSortRealWithPermutation()` 338021822bcSHong Zhang @*/ 339d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSortSplitReal(PetscInt ncut, PetscInt n, PetscReal a[], PetscInt idx[]) 340d71ae5a4SJacob Faibussowitsch { 341021822bcSHong Zhang PetscInt i, mid, last, itmp, j, first; 342021822bcSHong Zhang PetscReal d, tmp; 343021822bcSHong Zhang PetscReal abskey; 344021822bcSHong Zhang 345021822bcSHong Zhang PetscFunctionBegin; 346021822bcSHong Zhang first = 0; 347021822bcSHong Zhang last = n - 1; 3483ba16761SJacob Faibussowitsch if (ncut < first || ncut > last) PetscFunctionReturn(PETSC_SUCCESS); 349021822bcSHong Zhang 350021822bcSHong Zhang while (1) { 351021822bcSHong Zhang mid = first; 3522a274a78SSatish Balay d = a[mid]; 3532a274a78SSatish Balay abskey = PetscAbsReal(d); 354021822bcSHong Zhang i = last; 355021822bcSHong Zhang for (j = first + 1; j <= i; ++j) { 3562a274a78SSatish Balay d = a[j]; 3572a274a78SSatish Balay if (PetscAbsReal(d) >= abskey) { 358021822bcSHong Zhang ++mid; 359021822bcSHong Zhang /* interchange */ 3609371c9d4SSatish Balay tmp = a[mid]; 3619371c9d4SSatish Balay itmp = idx[mid]; 3629371c9d4SSatish Balay a[mid] = a[j]; 3639371c9d4SSatish Balay idx[mid] = idx[j]; 3649371c9d4SSatish Balay a[j] = tmp; 3659371c9d4SSatish Balay idx[j] = itmp; 366021822bcSHong Zhang } 367021822bcSHong Zhang } 368021822bcSHong Zhang 369021822bcSHong Zhang /* interchange */ 3709371c9d4SSatish Balay tmp = a[mid]; 3719371c9d4SSatish Balay itmp = idx[mid]; 3729371c9d4SSatish Balay a[mid] = a[first]; 3739371c9d4SSatish Balay idx[mid] = idx[first]; 3749371c9d4SSatish Balay a[first] = tmp; 3759371c9d4SSatish Balay idx[first] = itmp; 376021822bcSHong Zhang 377021822bcSHong Zhang /* test for while loop */ 378a297a907SKarl Rupp if (mid == ncut) break; 379a297a907SKarl Rupp else if (mid > ncut) last = mid - 1; 380a297a907SKarl Rupp else first = mid + 1; 381021822bcSHong Zhang } 3823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 383021822bcSHong Zhang } 384