1545aaa6fSHong Zhang #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 2424fc49dSJoe Pusztay #include <petscdm.h> 31d06f6b3SHong Zhang static PetscErrorCode TSRHSSplitGetRHSSplit(TS ts,const char splitname[],TS_RHSSplitLink *isplit) 41d06f6b3SHong Zhang { 51d06f6b3SHong Zhang PetscBool found = PETSC_FALSE; 61d06f6b3SHong Zhang 71d06f6b3SHong Zhang PetscFunctionBegin; 81d06f6b3SHong Zhang *isplit = ts->tsrhssplit; 91d06f6b3SHong Zhang /* look up the split */ 101d06f6b3SHong Zhang while (*isplit) { 119566063dSJacob Faibussowitsch PetscCall(PetscStrcmp((*isplit)->splitname,splitname,&found)); 121d06f6b3SHong Zhang if (found) break; 131d06f6b3SHong Zhang *isplit = (*isplit)->next; 141d06f6b3SHong Zhang } 151d06f6b3SHong Zhang PetscFunctionReturn(0); 161d06f6b3SHong Zhang } 171d06f6b3SHong Zhang 18545aaa6fSHong Zhang /*@C 19545aaa6fSHong Zhang TSRHSSplitSetIS - Set the index set for the specified split 20545aaa6fSHong Zhang 21545aaa6fSHong Zhang Logically Collective on TS 22545aaa6fSHong Zhang 23545aaa6fSHong Zhang Input Parameters: 24545aaa6fSHong Zhang + ts - the TS context obtained from TSCreate() 25545aaa6fSHong Zhang . splitname - name of this split, if NULL the number of the split is used 26545aaa6fSHong Zhang - is - the index set for part of the solution vector 27545aaa6fSHong Zhang 28545aaa6fSHong Zhang Level: intermediate 29545aaa6fSHong Zhang 301d06f6b3SHong Zhang .seealso: TSRHSSplitGetIS() 31545aaa6fSHong Zhang 32545aaa6fSHong Zhang @*/ 33545aaa6fSHong Zhang PetscErrorCode TSRHSSplitSetIS(TS ts,const char splitname[],IS is) 34545aaa6fSHong Zhang { 351d06f6b3SHong Zhang TS_RHSSplitLink newsplit,next = ts->tsrhssplit; 36545aaa6fSHong Zhang char prefix[128]; 37545aaa6fSHong Zhang 38545aaa6fSHong Zhang PetscFunctionBegin; 39545aaa6fSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 40545aaa6fSHong Zhang PetscValidHeaderSpecific(is,IS_CLASSID,3); 411d06f6b3SHong Zhang 429566063dSJacob Faibussowitsch PetscCall(PetscNew(&newsplit)); 43545aaa6fSHong Zhang if (splitname) { 449566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(splitname,&newsplit->splitname)); 45545aaa6fSHong Zhang } else { 469566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(8,&newsplit->splitname)); 47*63a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(newsplit->splitname,7,"%" PetscInt_FMT,ts->num_rhs_splits)); 48545aaa6fSHong Zhang } 499566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)is)); 50545aaa6fSHong Zhang newsplit->is = is; 519566063dSJacob Faibussowitsch PetscCall(TSCreate(PetscObjectComm((PetscObject)ts),&newsplit->ts)); 52109dc152SHong Zhang 539566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)newsplit->ts,(PetscObject)ts,1)); 549566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)ts,(PetscObject)newsplit->ts)); 559566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix,sizeof(prefix),"%srhsplit_%s_",((PetscObject)ts)->prefix ? ((PetscObject)ts)->prefix : "",newsplit->splitname)); 569566063dSJacob Faibussowitsch PetscCall(TSSetOptionsPrefix(newsplit->ts,prefix)); 571d06f6b3SHong Zhang if (!next) ts->tsrhssplit = newsplit; 581d06f6b3SHong Zhang else { 591d06f6b3SHong Zhang while (next->next) next = next->next; 601d06f6b3SHong Zhang next->next = newsplit; 611d06f6b3SHong Zhang } 621d06f6b3SHong Zhang ts->num_rhs_splits++; 63545aaa6fSHong Zhang PetscFunctionReturn(0); 64545aaa6fSHong Zhang } 65545aaa6fSHong Zhang 66545aaa6fSHong Zhang /*@C 67545aaa6fSHong Zhang TSRHSSplitGetIS - Retrieves the elements for a split as an IS 68545aaa6fSHong Zhang 69545aaa6fSHong Zhang Logically Collective on TS 70545aaa6fSHong Zhang 71545aaa6fSHong Zhang Input Parameters: 72545aaa6fSHong Zhang + ts - the TS context obtained from TSCreate() 73545aaa6fSHong Zhang - splitname - name of this split 74545aaa6fSHong Zhang 75545aaa6fSHong Zhang Output Parameters: 76545aaa6fSHong Zhang - is - the index set for part of the solution vector 77545aaa6fSHong Zhang 78545aaa6fSHong Zhang Level: intermediate 79545aaa6fSHong Zhang 801d06f6b3SHong Zhang .seealso: TSRHSSplitSetIS() 81545aaa6fSHong Zhang 82545aaa6fSHong Zhang @*/ 83545aaa6fSHong Zhang PetscErrorCode TSRHSSplitGetIS(TS ts,const char splitname[],IS *is) 84545aaa6fSHong Zhang { 851d06f6b3SHong Zhang TS_RHSSplitLink isplit; 86545aaa6fSHong Zhang 87545aaa6fSHong Zhang PetscFunctionBegin; 88545aaa6fSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 891d06f6b3SHong Zhang *is = NULL; 90545aaa6fSHong Zhang /* look up the split */ 919566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetRHSSplit(ts,splitname,&isplit)); 921d06f6b3SHong Zhang if (isplit) *is = isplit->is; 93545aaa6fSHong Zhang PetscFunctionReturn(0); 94545aaa6fSHong Zhang } 95545aaa6fSHong Zhang 96545aaa6fSHong Zhang /*@C 97545aaa6fSHong Zhang TSRHSSplitSetRHSFunction - Set the split right-hand-side functions. 98545aaa6fSHong Zhang 99545aaa6fSHong Zhang Logically Collective on TS 100545aaa6fSHong Zhang 101545aaa6fSHong Zhang Input Parameters: 102545aaa6fSHong Zhang + ts - the TS context obtained from TSCreate() 103545aaa6fSHong Zhang . splitname - name of this split 104545aaa6fSHong Zhang . r - vector to hold the residual (or NULL to have it created internally) 105545aaa6fSHong Zhang . rhsfunc - the RHS function evaluation routine 106545aaa6fSHong Zhang - ctx - user-defined context for private data for the split function evaluation routine (may be NULL) 107545aaa6fSHong Zhang 108545aaa6fSHong Zhang Calling sequence of fun: 109545aaa6fSHong Zhang $ rhsfunc(TS ts,PetscReal t,Vec u,Vec f,ctx); 110545aaa6fSHong Zhang 111545aaa6fSHong Zhang + t - time at step/stage being solved 112545aaa6fSHong Zhang . u - state vector 113545aaa6fSHong Zhang . f - function vector 114545aaa6fSHong Zhang - ctx - [optional] user-defined context for matrix evaluation routine (may be NULL) 115545aaa6fSHong Zhang 116545aaa6fSHong Zhang Level: beginner 117545aaa6fSHong Zhang 118545aaa6fSHong Zhang @*/ 119545aaa6fSHong Zhang PetscErrorCode TSRHSSplitSetRHSFunction(TS ts,const char splitname[],Vec r,TSRHSFunction rhsfunc,void *ctx) 120545aaa6fSHong Zhang { 1211d06f6b3SHong Zhang TS_RHSSplitLink isplit; 122424fc49dSJoe Pusztay DM dmc; 123545aaa6fSHong Zhang Vec subvec,ralloc = NULL; 124545aaa6fSHong Zhang 125545aaa6fSHong Zhang PetscFunctionBegin; 126545aaa6fSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 127064a246eSJacob Faibussowitsch if (r) PetscValidHeaderSpecific(r,VEC_CLASSID,3); 128545aaa6fSHong Zhang 129545aaa6fSHong Zhang /* look up the split */ 1309566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetRHSSplit(ts,splitname,&isplit)); 1313c633725SBarry Smith PetscCheck(isplit,PETSC_COMM_SELF,PETSC_ERR_USER,"The split %s is not created, check the split name or call TSRHSSplitSetIS() to create one",splitname); 132545aaa6fSHong Zhang 1331d06f6b3SHong Zhang if (!r && ts->vec_sol) { 1349566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(ts->vec_sol,isplit->is,&subvec)); 1359566063dSJacob Faibussowitsch PetscCall(VecDuplicate(subvec,&ralloc)); 136545aaa6fSHong Zhang r = ralloc; 1379566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(ts->vec_sol,isplit->is,&subvec)); 138545aaa6fSHong Zhang } 139424fc49dSJoe Pusztay 140424fc49dSJoe Pusztay if (ts->dm) { 141424fc49dSJoe Pusztay PetscInt dim; 142424fc49dSJoe Pusztay 1439566063dSJacob Faibussowitsch PetscCall(DMGetDimension(ts->dm, &dim)); 144424fc49dSJoe Pusztay if (dim != -1) { 1459566063dSJacob Faibussowitsch PetscCall(DMClone(ts->dm, &dmc)); 1469566063dSJacob Faibussowitsch PetscCall(TSSetDM(isplit->ts, dmc)); 1479566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dmc)); 148424fc49dSJoe Pusztay } 149424fc49dSJoe Pusztay } 150424fc49dSJoe Pusztay 1519566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(isplit->ts,r,rhsfunc,ctx)); 1529566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ralloc)); 153545aaa6fSHong Zhang PetscFunctionReturn(0); 154545aaa6fSHong Zhang } 155545aaa6fSHong Zhang 156545aaa6fSHong Zhang /*@C 1571d06f6b3SHong Zhang TSRHSSplitGetSubTS - Get the sub-TS by split name. 1581d06f6b3SHong Zhang 1591d06f6b3SHong Zhang Logically Collective on TS 1601d06f6b3SHong Zhang 1616b867d5aSJose E. Roman Input Parameter: 1626b867d5aSJose E. Roman . ts - the TS context obtained from TSCreate() 1636b867d5aSJose E. Roman 1641d06f6b3SHong Zhang Output Parameters: 1651d06f6b3SHong Zhang + splitname - the number of the split 1661d06f6b3SHong Zhang - subts - the array of TS contexts 1671d06f6b3SHong Zhang 1681d06f6b3SHong Zhang Level: advanced 1691d06f6b3SHong Zhang 1701d06f6b3SHong Zhang .seealso: TSGetRHSSplitFunction() 1711d06f6b3SHong Zhang @*/ 1721d06f6b3SHong Zhang PetscErrorCode TSRHSSplitGetSubTS(TS ts,const char splitname[],TS *subts) 1731d06f6b3SHong Zhang { 1741d06f6b3SHong Zhang TS_RHSSplitLink isplit; 1751d06f6b3SHong Zhang 1761d06f6b3SHong Zhang PetscFunctionBegin; 1771d06f6b3SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1781d06f6b3SHong Zhang PetscValidPointer(subts,3); 1791d06f6b3SHong Zhang *subts = NULL; 1801d06f6b3SHong Zhang /* look up the split */ 1819566063dSJacob Faibussowitsch PetscCall(TSRHSSplitGetRHSSplit(ts,splitname,&isplit)); 1821d06f6b3SHong Zhang if (isplit) *subts = isplit->ts; 1831d06f6b3SHong Zhang PetscFunctionReturn(0); 1841d06f6b3SHong Zhang } 1851d06f6b3SHong Zhang 1861d06f6b3SHong Zhang /*@C 1871d06f6b3SHong Zhang TSRHSSplitGetSubTSs - Get an array of all sub-TS contexts. 188545aaa6fSHong Zhang 189545aaa6fSHong Zhang Logically Collective on TS 190545aaa6fSHong Zhang 1916b867d5aSJose E. Roman Input Parameter: 1926b867d5aSJose E. Roman . ts - the TS context obtained from TSCreate() 1936b867d5aSJose E. Roman 194545aaa6fSHong Zhang Output Parameters: 195545aaa6fSHong Zhang + n - the number of splits 196545aaa6fSHong Zhang - subksp - the array of TS contexts 197545aaa6fSHong Zhang 198545aaa6fSHong Zhang Note: 199545aaa6fSHong Zhang After TSRHSSplitGetSubTS() the array of TSs is to be freed by the user with PetscFree() 200545aaa6fSHong Zhang (not the TS just the array that contains them). 201545aaa6fSHong Zhang 202545aaa6fSHong Zhang Level: advanced 203545aaa6fSHong Zhang 204545aaa6fSHong Zhang .seealso: TSGetRHSSplitFunction() 205545aaa6fSHong Zhang @*/ 2061d06f6b3SHong Zhang PetscErrorCode TSRHSSplitGetSubTSs(TS ts,PetscInt *n,TS *subts[]) 207545aaa6fSHong Zhang { 2081d06f6b3SHong Zhang TS_RHSSplitLink ilink = ts->tsrhssplit; 209545aaa6fSHong Zhang PetscInt i = 0; 210545aaa6fSHong Zhang 211545aaa6fSHong Zhang PetscFunctionBegin; 212545aaa6fSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 213545aaa6fSHong Zhang if (subts) { 2149566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ts->num_rhs_splits,subts)); 2151d06f6b3SHong Zhang while (ilink) { 2161d06f6b3SHong Zhang (*subts)[i++] = ilink->ts; 2171d06f6b3SHong Zhang ilink = ilink->next; 218545aaa6fSHong Zhang } 219545aaa6fSHong Zhang } 220545aaa6fSHong Zhang if (n) *n = ts->num_rhs_splits; 221545aaa6fSHong Zhang PetscFunctionReturn(0); 222545aaa6fSHong Zhang } 223