11447629fSBarry Smith 21447629fSBarry Smith /* 31447629fSBarry Smith These AO application ordering routines do not require that the input 41447629fSBarry Smith be a permutation, but merely a 1-1 mapping. This implementation still 51447629fSBarry Smith keeps the entire ordering on each processor. 61447629fSBarry Smith */ 71447629fSBarry Smith 81447629fSBarry Smith #include <../src/vec/is/ao/aoimpl.h> /*I "petscao.h" I*/ 91447629fSBarry Smith 101447629fSBarry Smith typedef struct { 111447629fSBarry Smith PetscInt N; 121447629fSBarry Smith PetscInt *app; /* app[i] is the partner for petsc[appPerm[i]] */ 131447629fSBarry Smith PetscInt *appPerm; 141447629fSBarry Smith PetscInt *petsc; /* petsc[j] is the partner for app[petscPerm[j]] */ 151447629fSBarry Smith PetscInt *petscPerm; 161447629fSBarry Smith } AO_Mapping; 171447629fSBarry Smith 181447629fSBarry Smith PetscErrorCode AODestroy_Mapping(AO ao) 191447629fSBarry Smith { 201447629fSBarry Smith AO_Mapping *aomap = (AO_Mapping*) ao->data; 211447629fSBarry Smith 221447629fSBarry Smith PetscFunctionBegin; 239566063dSJacob Faibussowitsch PetscCall(PetscFree4(aomap->app,aomap->appPerm,aomap->petsc,aomap->petscPerm)); 249566063dSJacob Faibussowitsch PetscCall(PetscFree(aomap)); 251447629fSBarry Smith PetscFunctionReturn(0); 261447629fSBarry Smith } 271447629fSBarry Smith 281447629fSBarry Smith PetscErrorCode AOView_Mapping(AO ao, PetscViewer viewer) 291447629fSBarry Smith { 301447629fSBarry Smith AO_Mapping *aomap = (AO_Mapping*) ao->data; 311447629fSBarry Smith PetscMPIInt rank; 321447629fSBarry Smith PetscInt i; 331447629fSBarry Smith PetscBool iascii; 341447629fSBarry Smith 351447629fSBarry Smith PetscFunctionBegin; 369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)ao), &rank)); 371447629fSBarry Smith if (rank) PetscFunctionReturn(0); 389566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii)); 391447629fSBarry Smith if (iascii) { 402abc8c78SJacob Faibussowitsch PetscViewerASCIIPrintf(viewer, "Number of elements in ordering %" PetscInt_FMT "\n", aomap->N); 411447629fSBarry Smith PetscViewerASCIIPrintf(viewer, " App. PETSc\n"); 421447629fSBarry Smith for (i = 0; i < aomap->N; i++) { 432abc8c78SJacob Faibussowitsch PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n", i, aomap->app[i], aomap->petsc[aomap->appPerm[i]]); 441447629fSBarry Smith } 451447629fSBarry Smith } 461447629fSBarry Smith PetscFunctionReturn(0); 471447629fSBarry Smith } 481447629fSBarry Smith 491447629fSBarry Smith PetscErrorCode AOPetscToApplication_Mapping(AO ao, PetscInt n, PetscInt *ia) 501447629fSBarry Smith { 511447629fSBarry Smith AO_Mapping *aomap = (AO_Mapping*) ao->data; 521447629fSBarry Smith PetscInt *app = aomap->app; 531447629fSBarry Smith PetscInt *petsc = aomap->petsc; 541447629fSBarry Smith PetscInt *perm = aomap->petscPerm; 551447629fSBarry Smith PetscInt N = aomap->N; 561447629fSBarry Smith PetscInt low, high, mid=0; 571447629fSBarry Smith PetscInt idex; 581447629fSBarry Smith PetscInt i; 591447629fSBarry Smith 601447629fSBarry Smith /* It would be possible to use a single bisection search, which 611447629fSBarry Smith recursively divided the indices to be converted, and searched 621447629fSBarry Smith partitions which contained an index. This would result in 631447629fSBarry Smith better running times if indices are clustered. 641447629fSBarry Smith */ 651447629fSBarry Smith PetscFunctionBegin; 661447629fSBarry Smith for (i = 0; i < n; i++) { 671447629fSBarry Smith idex = ia[i]; 681447629fSBarry Smith if (idex < 0) continue; 691447629fSBarry Smith /* Use bisection since the array is sorted */ 701447629fSBarry Smith low = 0; 711447629fSBarry Smith high = N - 1; 721447629fSBarry Smith while (low <= high) { 731447629fSBarry Smith mid = (low + high)/2; 741447629fSBarry Smith if (idex == petsc[mid]) break; 751447629fSBarry Smith else if (idex < petsc[mid]) high = mid - 1; 761447629fSBarry Smith else low = mid + 1; 771447629fSBarry Smith } 781447629fSBarry Smith if (low > high) ia[i] = -1; 791447629fSBarry Smith else ia[i] = app[perm[mid]]; 801447629fSBarry Smith } 811447629fSBarry Smith PetscFunctionReturn(0); 821447629fSBarry Smith } 831447629fSBarry Smith 841447629fSBarry Smith PetscErrorCode AOApplicationToPetsc_Mapping(AO ao, PetscInt n, PetscInt *ia) 851447629fSBarry Smith { 861447629fSBarry Smith AO_Mapping *aomap = (AO_Mapping*) ao->data; 871447629fSBarry Smith PetscInt *app = aomap->app; 881447629fSBarry Smith PetscInt *petsc = aomap->petsc; 891447629fSBarry Smith PetscInt *perm = aomap->appPerm; 901447629fSBarry Smith PetscInt N = aomap->N; 911447629fSBarry Smith PetscInt low, high, mid=0; 921447629fSBarry Smith PetscInt idex; 931447629fSBarry Smith PetscInt i; 941447629fSBarry Smith 951447629fSBarry Smith /* It would be possible to use a single bisection search, which 961447629fSBarry Smith recursively divided the indices to be converted, and searched 971447629fSBarry Smith partitions which contained an index. This would result in 981447629fSBarry Smith better running times if indices are clustered. 991447629fSBarry Smith */ 1001447629fSBarry Smith PetscFunctionBegin; 1011447629fSBarry Smith for (i = 0; i < n; i++) { 1021447629fSBarry Smith idex = ia[i]; 1031447629fSBarry Smith if (idex < 0) continue; 1041447629fSBarry Smith /* Use bisection since the array is sorted */ 1051447629fSBarry Smith low = 0; 1061447629fSBarry Smith high = N - 1; 1071447629fSBarry Smith while (low <= high) { 1081447629fSBarry Smith mid = (low + high)/2; 1091447629fSBarry Smith if (idex == app[mid]) break; 1101447629fSBarry Smith else if (idex < app[mid]) high = mid - 1; 1111447629fSBarry Smith else low = mid + 1; 1121447629fSBarry Smith } 1131447629fSBarry Smith if (low > high) ia[i] = -1; 1141447629fSBarry Smith else ia[i] = petsc[perm[mid]]; 1151447629fSBarry Smith } 1161447629fSBarry Smith PetscFunctionReturn(0); 1171447629fSBarry Smith } 1181447629fSBarry Smith 119267267bdSJacob Faibussowitsch static struct _AOOps AOps = { 120267267bdSJacob Faibussowitsch PetscDesignatedInitializer(view,AOView_Mapping), 121267267bdSJacob Faibussowitsch PetscDesignatedInitializer(destroy,AODestroy_Mapping), 122267267bdSJacob Faibussowitsch PetscDesignatedInitializer(petsctoapplication,AOPetscToApplication_Mapping), 123267267bdSJacob Faibussowitsch PetscDesignatedInitializer(applicationtopetsc,AOApplicationToPetsc_Mapping), 124267267bdSJacob Faibussowitsch }; 1251447629fSBarry Smith 1261447629fSBarry Smith /*@C 1271447629fSBarry Smith AOMappingHasApplicationIndex - Searches for the supplied application index. 1281447629fSBarry Smith 1291447629fSBarry Smith Input Parameters: 1301447629fSBarry Smith + ao - The AOMapping 1311447629fSBarry Smith - index - The application index 1321447629fSBarry Smith 1331447629fSBarry Smith Output Parameter: 1341447629fSBarry Smith . hasIndex - Flag is PETSC_TRUE if the index exists 1351447629fSBarry Smith 1361447629fSBarry Smith Level: intermediate 1371447629fSBarry Smith 1381447629fSBarry Smith .seealso: AOMappingHasPetscIndex(), AOCreateMapping() 1391447629fSBarry Smith @*/ 1401447629fSBarry Smith PetscErrorCode AOMappingHasApplicationIndex(AO ao, PetscInt idex, PetscBool *hasIndex) 1411447629fSBarry Smith { 1421447629fSBarry Smith AO_Mapping *aomap; 1431447629fSBarry Smith PetscInt *app; 1441447629fSBarry Smith PetscInt low, high, mid; 1451447629fSBarry Smith 1461447629fSBarry Smith PetscFunctionBegin; 1471447629fSBarry Smith PetscValidHeaderSpecific(ao, AO_CLASSID,1); 148dadcf809SJacob Faibussowitsch PetscValidBoolPointer(hasIndex,3); 1491447629fSBarry Smith aomap = (AO_Mapping*) ao->data; 1501447629fSBarry Smith app = aomap->app; 1511447629fSBarry Smith /* Use bisection since the array is sorted */ 1521447629fSBarry Smith low = 0; 1531447629fSBarry Smith high = aomap->N - 1; 1541447629fSBarry Smith while (low <= high) { 1551447629fSBarry Smith mid = (low + high)/2; 1561447629fSBarry Smith if (idex == app[mid]) break; 1571447629fSBarry Smith else if (idex < app[mid]) high = mid - 1; 1581447629fSBarry Smith else low = mid + 1; 1591447629fSBarry Smith } 1601447629fSBarry Smith if (low > high) *hasIndex = PETSC_FALSE; 1611447629fSBarry Smith else *hasIndex = PETSC_TRUE; 1621447629fSBarry Smith PetscFunctionReturn(0); 1631447629fSBarry Smith } 1641447629fSBarry Smith 1651447629fSBarry Smith /*@ 1661447629fSBarry Smith AOMappingHasPetscIndex - Searches for the supplied petsc index. 1671447629fSBarry Smith 1681447629fSBarry Smith Input Parameters: 1691447629fSBarry Smith + ao - The AOMapping 1701447629fSBarry Smith - index - The petsc index 1711447629fSBarry Smith 1721447629fSBarry Smith Output Parameter: 1731447629fSBarry Smith . hasIndex - Flag is PETSC_TRUE if the index exists 1741447629fSBarry Smith 1751447629fSBarry Smith Level: intermediate 1761447629fSBarry Smith 1771447629fSBarry Smith .seealso: AOMappingHasApplicationIndex(), AOCreateMapping() 1781447629fSBarry Smith @*/ 1791447629fSBarry Smith PetscErrorCode AOMappingHasPetscIndex(AO ao, PetscInt idex, PetscBool *hasIndex) 1801447629fSBarry Smith { 1811447629fSBarry Smith AO_Mapping *aomap; 1821447629fSBarry Smith PetscInt *petsc; 1831447629fSBarry Smith PetscInt low, high, mid; 1841447629fSBarry Smith 1851447629fSBarry Smith PetscFunctionBegin; 1861447629fSBarry Smith PetscValidHeaderSpecific(ao, AO_CLASSID,1); 187dadcf809SJacob Faibussowitsch PetscValidBoolPointer(hasIndex,3); 1881447629fSBarry Smith aomap = (AO_Mapping*) ao->data; 1891447629fSBarry Smith petsc = aomap->petsc; 1901447629fSBarry Smith /* Use bisection since the array is sorted */ 1911447629fSBarry Smith low = 0; 1921447629fSBarry Smith high = aomap->N - 1; 1931447629fSBarry Smith while (low <= high) { 1941447629fSBarry Smith mid = (low + high)/2; 1951447629fSBarry Smith if (idex == petsc[mid]) break; 1961447629fSBarry Smith else if (idex < petsc[mid]) high = mid - 1; 1971447629fSBarry Smith else low = mid + 1; 1981447629fSBarry Smith } 1991447629fSBarry Smith if (low > high) *hasIndex = PETSC_FALSE; 2001447629fSBarry Smith else *hasIndex = PETSC_TRUE; 2011447629fSBarry Smith PetscFunctionReturn(0); 2021447629fSBarry Smith } 2031447629fSBarry Smith 2041447629fSBarry Smith /*@C 2051447629fSBarry Smith AOCreateMapping - Creates a basic application mapping using two integer arrays. 2061447629fSBarry Smith 2071447629fSBarry Smith Input Parameters: 2081447629fSBarry Smith + comm - MPI communicator that is to share AO 2091447629fSBarry Smith . napp - size of integer arrays 2101447629fSBarry Smith . myapp - integer array that defines an ordering 2111447629fSBarry Smith - mypetsc - integer array that defines another ordering (may be NULL to indicate the identity ordering) 2121447629fSBarry Smith 2131447629fSBarry Smith Output Parameter: 2141447629fSBarry Smith . aoout - the new application mapping 2151447629fSBarry Smith 2161447629fSBarry Smith Options Database Key: 217147403d9SBarry Smith . -ao_view - call AOView() at the conclusion of AOCreateMapping() 2181447629fSBarry Smith 2191447629fSBarry Smith Level: beginner 2201447629fSBarry Smith 22195452b02SPatrick Sanan Notes: 22295452b02SPatrick Sanan the arrays myapp and mypetsc need NOT contain the all the integers 0 to napp-1, that is there CAN be "holes" in the indices. 2231447629fSBarry Smith Use AOCreateBasic() or AOCreateBasicIS() if they do not have holes for better performance. 2241447629fSBarry Smith 2251447629fSBarry Smith .seealso: AOCreateBasic(), AOCreateBasic(), AOCreateMappingIS(), AODestroy() 2261447629fSBarry Smith @*/ 2271447629fSBarry Smith PetscErrorCode AOCreateMapping(MPI_Comm comm,PetscInt napp,const PetscInt myapp[],const PetscInt mypetsc[],AO *aoout) 2281447629fSBarry Smith { 2291447629fSBarry Smith AO ao; 2301447629fSBarry Smith AO_Mapping *aomap; 2311447629fSBarry Smith PetscInt *allpetsc, *allapp; 2321447629fSBarry Smith PetscInt *petscPerm, *appPerm; 2331447629fSBarry Smith PetscInt *petsc; 2341447629fSBarry Smith PetscMPIInt size, rank,*lens, *disp,nnapp; 2351447629fSBarry Smith PetscInt N, start; 2361447629fSBarry Smith PetscInt i; 2371447629fSBarry Smith 2381447629fSBarry Smith PetscFunctionBegin; 2391447629fSBarry Smith PetscValidPointer(aoout,5); 2404c8fdceaSLisandro Dalcin *aoout = NULL; 2419566063dSJacob Faibussowitsch PetscCall(AOInitializePackage()); 2421447629fSBarry Smith 2439566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(ao, AO_CLASSID, "AO", "Application Ordering", "AO", comm, AODestroy, AOView)); 2449566063dSJacob Faibussowitsch PetscCall(PetscNewLog(ao,&aomap)); 2459566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(ao->ops, &AOps, sizeof(AOps))); 2461447629fSBarry Smith ao->data = (void*) aomap; 2471447629fSBarry Smith 2481447629fSBarry Smith /* transmit all lengths to all processors */ 2499566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 2509566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 2519566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(size, &lens,size,&disp)); 2521447629fSBarry Smith nnapp = napp; 2539566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&nnapp, 1, MPI_INT, lens, 1, MPI_INT, comm)); 2541447629fSBarry Smith N = 0; 2551447629fSBarry Smith for (i = 0; i < size; i++) { 2561447629fSBarry Smith disp[i] = N; 2571447629fSBarry Smith N += lens[i]; 2581447629fSBarry Smith } 2591447629fSBarry Smith aomap->N = N; 2601447629fSBarry Smith ao->N = N; 2611447629fSBarry Smith ao->n = N; 2621447629fSBarry Smith 2631447629fSBarry Smith /* If mypetsc is 0 then use "natural" numbering */ 2641447629fSBarry Smith if (!mypetsc) { 2651447629fSBarry Smith start = disp[rank]; 2669566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(napp+1, &petsc)); 2671447629fSBarry Smith for (i = 0; i < napp; i++) petsc[i] = start + i; 2681447629fSBarry Smith } else { 2691447629fSBarry Smith petsc = (PetscInt*)mypetsc; 2701447629fSBarry Smith } 2711447629fSBarry Smith 2721447629fSBarry Smith /* get all indices on all processors */ 2739566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(N, &allapp,N,&appPerm,N,&allpetsc,N,&petscPerm)); 2749566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgatherv((void*)myapp, napp, MPIU_INT, allapp, lens, disp, MPIU_INT, comm)); 2759566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgatherv((void*)petsc, napp, MPIU_INT, allpetsc, lens, disp, MPIU_INT, comm)); 2769566063dSJacob Faibussowitsch PetscCall(PetscFree2(lens,disp)); 2771447629fSBarry Smith 2781447629fSBarry Smith /* generate a list of application and PETSc node numbers */ 2799566063dSJacob Faibussowitsch PetscCall(PetscMalloc4(N, &aomap->app,N,&aomap->appPerm,N,&aomap->petsc,N,&aomap->petscPerm)); 2809566063dSJacob Faibussowitsch PetscCall(PetscLogObjectMemory((PetscObject)ao, 4*N * sizeof(PetscInt))); 2811447629fSBarry Smith for (i = 0; i < N; i++) { 2821447629fSBarry Smith appPerm[i] = i; 2831447629fSBarry Smith petscPerm[i] = i; 2841447629fSBarry Smith } 2859566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithPermutation(N, allpetsc, petscPerm)); 2869566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithPermutation(N, allapp, appPerm)); 2871447629fSBarry Smith /* Form sorted arrays of indices */ 2881447629fSBarry Smith for (i = 0; i < N; i++) { 2891447629fSBarry Smith aomap->app[i] = allapp[appPerm[i]]; 2901447629fSBarry Smith aomap->petsc[i] = allpetsc[petscPerm[i]]; 2911447629fSBarry Smith } 2921447629fSBarry Smith /* Invert petscPerm[] into aomap->petscPerm[] */ 2931447629fSBarry Smith for (i = 0; i < N; i++) aomap->petscPerm[petscPerm[i]] = i; 2941447629fSBarry Smith 2951447629fSBarry Smith /* Form map between aomap->app[] and aomap->petsc[] */ 2961447629fSBarry Smith for (i = 0; i < N; i++) aomap->appPerm[i] = aomap->petscPerm[appPerm[i]]; 2971447629fSBarry Smith 2981447629fSBarry Smith /* Invert appPerm[] into allapp[] */ 2991447629fSBarry Smith for (i = 0; i < N; i++) allapp[appPerm[i]] = i; 3001447629fSBarry Smith 3011447629fSBarry Smith /* Form map between aomap->petsc[] and aomap->app[] */ 3021447629fSBarry Smith for (i = 0; i < N; i++) aomap->petscPerm[i] = allapp[petscPerm[i]]; 3031447629fSBarry Smith 30476bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { 3051447629fSBarry Smith /* Check that the permutations are complementary */ 3061447629fSBarry Smith for (i = 0; i < N; i++) { 3072c71b3e2SJacob Faibussowitsch PetscCheckFalse(i != aomap->appPerm[aomap->petscPerm[i]],PETSC_COMM_SELF,PETSC_ERR_PLIB, "Invalid ordering"); 3081447629fSBarry Smith } 30976bd3646SJed Brown } 3101447629fSBarry Smith /* Cleanup */ 3111447629fSBarry Smith if (!mypetsc) { 3129566063dSJacob Faibussowitsch PetscCall(PetscFree(petsc)); 3131447629fSBarry Smith } 3149566063dSJacob Faibussowitsch PetscCall(PetscFree4(allapp,appPerm,allpetsc,petscPerm)); 3151447629fSBarry Smith 3169566063dSJacob Faibussowitsch PetscCall(AOViewFromOptions(ao,NULL,"-ao_view")); 3171447629fSBarry Smith 3181447629fSBarry Smith *aoout = ao; 3191447629fSBarry Smith PetscFunctionReturn(0); 3201447629fSBarry Smith } 3211447629fSBarry Smith 322487a658cSBarry Smith /*@ 3231447629fSBarry Smith AOCreateMappingIS - Creates a basic application ordering using two index sets. 3241447629fSBarry Smith 3251447629fSBarry Smith Input Parameters: 3261447629fSBarry Smith + comm - MPI communicator that is to share AO 3271447629fSBarry Smith . isapp - index set that defines an ordering 3281447629fSBarry Smith - ispetsc - index set that defines another ordering, maybe NULL for identity IS 3291447629fSBarry Smith 3301447629fSBarry Smith Output Parameter: 3311447629fSBarry Smith . aoout - the new application ordering 3321447629fSBarry Smith 3331447629fSBarry Smith Options Database Key: 334147403d9SBarry Smith . -ao_view - call AOView() at the conclusion of AOCreateMappingIS() 3351447629fSBarry Smith 3361447629fSBarry Smith Level: beginner 3371447629fSBarry Smith 33895452b02SPatrick Sanan Notes: 33995452b02SPatrick Sanan the index sets isapp and ispetsc need NOT contain the all the integers 0 to N-1, that is there CAN be "holes" in the indices. 3401447629fSBarry Smith Use AOCreateBasic() or AOCreateBasicIS() if they do not have holes for better performance. 3411447629fSBarry Smith 3421447629fSBarry Smith .seealso: AOCreateBasic(), AOCreateMapping(), AODestroy() 3431447629fSBarry Smith @*/ 3441447629fSBarry Smith PetscErrorCode AOCreateMappingIS(IS isapp, IS ispetsc, AO *aoout) 3451447629fSBarry Smith { 3461447629fSBarry Smith MPI_Comm comm; 3471447629fSBarry Smith const PetscInt *mypetsc, *myapp; 3481447629fSBarry Smith PetscInt napp, npetsc; 3491447629fSBarry Smith 3501447629fSBarry Smith PetscFunctionBegin; 3519566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject) isapp, &comm)); 3529566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isapp, &napp)); 3531447629fSBarry Smith if (ispetsc) { 3549566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(ispetsc, &npetsc)); 355*08401ef6SPierre Jolivet PetscCheck(napp == npetsc,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ, "Local IS lengths must match"); 3569566063dSJacob Faibussowitsch PetscCall(ISGetIndices(ispetsc, &mypetsc)); 3571447629fSBarry Smith } else { 3581447629fSBarry Smith mypetsc = NULL; 3591447629fSBarry Smith } 3609566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isapp, &myapp)); 3611447629fSBarry Smith 3629566063dSJacob Faibussowitsch PetscCall(AOCreateMapping(comm, napp, myapp, mypetsc, aoout)); 3631447629fSBarry Smith 3649566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isapp, &myapp)); 3651447629fSBarry Smith if (ispetsc) { 3669566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(ispetsc, &mypetsc)); 3671447629fSBarry Smith } 3681447629fSBarry Smith PetscFunctionReturn(0); 3691447629fSBarry Smith } 370