1545aaa6fSHong Zhang #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 2545aaa6fSHong Zhang 3*1d06f6b3SHong Zhang static PetscErrorCode TSRHSSplitGetRHSSplit(TS ts,const char splitname[],TS_RHSSplitLink *isplit) 4*1d06f6b3SHong Zhang { 5*1d06f6b3SHong Zhang PetscBool found = PETSC_FALSE; 6*1d06f6b3SHong Zhang PetscErrorCode ierr; 7*1d06f6b3SHong Zhang 8*1d06f6b3SHong Zhang PetscFunctionBegin; 9*1d06f6b3SHong Zhang *isplit = ts->tsrhssplit; 10*1d06f6b3SHong Zhang /* look up the split */ 11*1d06f6b3SHong Zhang while (*isplit) { 12*1d06f6b3SHong Zhang ierr = PetscStrcmp((*isplit)->splitname,splitname,&found);CHKERRQ(ierr); 13*1d06f6b3SHong Zhang if (found) break; 14*1d06f6b3SHong Zhang *isplit = (*isplit)->next; 15*1d06f6b3SHong Zhang } 16*1d06f6b3SHong Zhang PetscFunctionReturn(0); 17*1d06f6b3SHong Zhang } 18*1d06f6b3SHong Zhang 19545aaa6fSHong Zhang /*@C 20545aaa6fSHong Zhang TSRHSSplitSetIS - Set the index set for the specified split 21545aaa6fSHong Zhang 22545aaa6fSHong Zhang Logically Collective on TS 23545aaa6fSHong Zhang 24545aaa6fSHong Zhang Input Parameters: 25545aaa6fSHong Zhang + ts - the TS context obtained from TSCreate() 26545aaa6fSHong Zhang . splitname - name of this split, if NULL the number of the split is used 27545aaa6fSHong Zhang - is - the index set for part of the solution vector 28545aaa6fSHong Zhang 29545aaa6fSHong Zhang Level: intermediate 30545aaa6fSHong Zhang 31*1d06f6b3SHong Zhang .seealso: TSRHSSplitGetIS() 32545aaa6fSHong Zhang 33545aaa6fSHong Zhang .keywords: TS, TSRHSSplit 34545aaa6fSHong Zhang @*/ 35545aaa6fSHong Zhang PetscErrorCode TSRHSSplitSetIS(TS ts,const char splitname[],IS is) 36545aaa6fSHong Zhang { 37*1d06f6b3SHong Zhang TS_RHSSplitLink newsplit,next = ts->tsrhssplit; 38545aaa6fSHong Zhang char prefix[128]; 39545aaa6fSHong Zhang PetscErrorCode ierr; 40545aaa6fSHong Zhang 41545aaa6fSHong Zhang PetscFunctionBegin; 42545aaa6fSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 43545aaa6fSHong Zhang PetscValidHeaderSpecific(is,IS_CLASSID,3); 44*1d06f6b3SHong Zhang 45545aaa6fSHong Zhang ierr = PetscNew(&newsplit);CHKERRQ(ierr); 46545aaa6fSHong Zhang if (splitname) { 47545aaa6fSHong Zhang ierr = PetscStrallocpy(splitname,&newsplit->splitname);CHKERRQ(ierr); 48545aaa6fSHong Zhang } else { 49545aaa6fSHong Zhang ierr = PetscMalloc1(8,&newsplit->splitname);CHKERRQ(ierr); 50545aaa6fSHong Zhang ierr = PetscSNPrintf(newsplit->splitname,7,"%D",ts->num_rhs_splits);CHKERRQ(ierr); 51545aaa6fSHong Zhang } 52545aaa6fSHong Zhang ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr); 53545aaa6fSHong Zhang newsplit->is = is; 54545aaa6fSHong Zhang ierr = TSCreate(PetscObjectComm((PetscObject)ts),&newsplit->ts);CHKERRQ(ierr); 55545aaa6fSHong Zhang ierr = PetscObjectIncrementTabLevel((PetscObject)newsplit->ts,(PetscObject)ts,1);CHKERRQ(ierr); 56545aaa6fSHong Zhang ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)newsplit->ts);CHKERRQ(ierr); 57545aaa6fSHong Zhang ierr = PetscSNPrintf(prefix,sizeof(prefix),"%srhsplit_%s_",((PetscObject)ts)->prefix ? ((PetscObject)ts)->prefix : "",newsplit->splitname); 58545aaa6fSHong Zhang ierr = TSSetOptionsPrefix(newsplit->ts,prefix);CHKERRQ(ierr); 59*1d06f6b3SHong Zhang if (!next) ts->tsrhssplit = newsplit; 60*1d06f6b3SHong Zhang else { 61*1d06f6b3SHong Zhang while (next->next) next = next->next; 62*1d06f6b3SHong Zhang next->next = newsplit; 63*1d06f6b3SHong Zhang } 64*1d06f6b3SHong Zhang ts->num_rhs_splits++; 65545aaa6fSHong Zhang PetscFunctionReturn(0); 66545aaa6fSHong Zhang } 67545aaa6fSHong Zhang 68545aaa6fSHong Zhang /*@C 69545aaa6fSHong Zhang TSRHSSplitGetIS - Retrieves the elements for a split as an IS 70545aaa6fSHong Zhang 71545aaa6fSHong Zhang Logically Collective on TS 72545aaa6fSHong Zhang 73545aaa6fSHong Zhang Input Parameters: 74545aaa6fSHong Zhang + ts - the TS context obtained from TSCreate() 75545aaa6fSHong Zhang - splitname - name of this split 76545aaa6fSHong Zhang 77545aaa6fSHong Zhang Output Parameters: 78545aaa6fSHong Zhang - is - the index set for part of the solution vector 79545aaa6fSHong Zhang 80545aaa6fSHong Zhang Level: intermediate 81545aaa6fSHong Zhang 82*1d06f6b3SHong Zhang .seealso: TSRHSSplitSetIS() 83545aaa6fSHong Zhang 84545aaa6fSHong Zhang .keywords: TS, TSRHSSplit 85545aaa6fSHong Zhang @*/ 86545aaa6fSHong Zhang PetscErrorCode TSRHSSplitGetIS(TS ts,const char splitname[],IS *is) 87545aaa6fSHong Zhang { 88*1d06f6b3SHong Zhang TS_RHSSplitLink isplit; 89545aaa6fSHong Zhang PetscErrorCode ierr; 90545aaa6fSHong Zhang 91545aaa6fSHong Zhang PetscFunctionBegin; 92545aaa6fSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 93*1d06f6b3SHong Zhang *is = NULL; 94545aaa6fSHong Zhang /* look up the split */ 95*1d06f6b3SHong Zhang ierr = TSRHSSplitGetRHSSplit(ts,splitname,&isplit);CHKERRQ(ierr); 96*1d06f6b3SHong Zhang if (isplit) *is = isplit->is; 97545aaa6fSHong Zhang PetscFunctionReturn(0); 98545aaa6fSHong Zhang } 99545aaa6fSHong Zhang 100545aaa6fSHong Zhang /*@C 101545aaa6fSHong Zhang TSRHSSplitSetRHSFunction - Set the split right-hand-side functions. 102545aaa6fSHong Zhang 103545aaa6fSHong Zhang Logically Collective on TS 104545aaa6fSHong Zhang 105545aaa6fSHong Zhang Input Parameters: 106545aaa6fSHong Zhang + ts - the TS context obtained from TSCreate() 107545aaa6fSHong Zhang . splitname - name of this split 108545aaa6fSHong Zhang . r - vector to hold the residual (or NULL to have it created internally) 109545aaa6fSHong Zhang . rhsfunc - the RHS function evaluation routine 110545aaa6fSHong Zhang - ctx - user-defined context for private data for the split function evaluation routine (may be NULL) 111545aaa6fSHong Zhang 112545aaa6fSHong Zhang Calling sequence of fun: 113545aaa6fSHong Zhang $ rhsfunc(TS ts,PetscReal t,Vec u,Vec f,ctx); 114545aaa6fSHong Zhang 115545aaa6fSHong Zhang + t - time at step/stage being solved 116545aaa6fSHong Zhang . u - state vector 117545aaa6fSHong Zhang . f - function vector 118545aaa6fSHong Zhang - ctx - [optional] user-defined context for matrix evaluation routine (may be NULL) 119545aaa6fSHong Zhang 120545aaa6fSHong Zhang Level: beginner 121545aaa6fSHong Zhang 122545aaa6fSHong Zhang .keywords: TS, timestep, set, ODE, Hamiltonian, Function 123545aaa6fSHong Zhang @*/ 124545aaa6fSHong Zhang PetscErrorCode TSRHSSplitSetRHSFunction(TS ts,const char splitname[],Vec r,TSRHSFunction rhsfunc,void *ctx) 125545aaa6fSHong Zhang { 126*1d06f6b3SHong Zhang TS_RHSSplitLink isplit; 127545aaa6fSHong Zhang Vec subvec,ralloc = NULL; 128545aaa6fSHong Zhang PetscErrorCode ierr; 129545aaa6fSHong Zhang 130545aaa6fSHong Zhang PetscFunctionBegin; 131545aaa6fSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 132545aaa6fSHong Zhang if (r) PetscValidHeaderSpecific(r,VEC_CLASSID,2); 133545aaa6fSHong Zhang 134545aaa6fSHong Zhang /* look up the split */ 135*1d06f6b3SHong Zhang ierr = TSRHSSplitGetRHSSplit(ts,splitname,&isplit);CHKERRQ(ierr); 136*1d06f6b3SHong Zhang if (!isplit) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_USER,"The split %s is not created, check the split name or call TSRHSSplitSetIS() to create one",splitname); 137545aaa6fSHong Zhang 138*1d06f6b3SHong Zhang if (!r && ts->vec_sol) { 139*1d06f6b3SHong Zhang ierr = VecGetSubVector(ts->vec_sol,isplit->is,&subvec);CHKERRQ(ierr); 140545aaa6fSHong Zhang ierr = VecDuplicate(subvec,&ralloc);CHKERRQ(ierr); 141545aaa6fSHong Zhang r = ralloc; 142*1d06f6b3SHong Zhang ierr = VecRestoreSubVector(ts->vec_sol,isplit->is,&subvec);CHKERRQ(ierr); 143545aaa6fSHong Zhang } 144*1d06f6b3SHong Zhang ierr = TSSetRHSFunction(isplit->ts,r,rhsfunc,ctx);CHKERRQ(ierr); 145545aaa6fSHong Zhang ierr = VecDestroy(&ralloc);CHKERRQ(ierr); 146545aaa6fSHong Zhang PetscFunctionReturn(0); 147545aaa6fSHong Zhang } 148545aaa6fSHong Zhang 149545aaa6fSHong Zhang /*@C 150*1d06f6b3SHong Zhang TSRHSSplitGetSubTS - Get the sub-TS by split name. 151*1d06f6b3SHong Zhang 152*1d06f6b3SHong Zhang Logically Collective on TS 153*1d06f6b3SHong Zhang 154*1d06f6b3SHong Zhang Output Parameters: 155*1d06f6b3SHong Zhang + splitname - the number of the split 156*1d06f6b3SHong Zhang - subts - the array of TS contexts 157*1d06f6b3SHong Zhang 158*1d06f6b3SHong Zhang Level: advanced 159*1d06f6b3SHong Zhang 160*1d06f6b3SHong Zhang .seealso: TSGetRHSSplitFunction() 161*1d06f6b3SHong Zhang @*/ 162*1d06f6b3SHong Zhang PetscErrorCode TSRHSSplitGetSubTS(TS ts,const char splitname[],TS *subts) 163*1d06f6b3SHong Zhang { 164*1d06f6b3SHong Zhang TS_RHSSplitLink isplit; 165*1d06f6b3SHong Zhang PetscErrorCode ierr; 166*1d06f6b3SHong Zhang 167*1d06f6b3SHong Zhang PetscFunctionBegin; 168*1d06f6b3SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 169*1d06f6b3SHong Zhang PetscValidPointer(subts,3); 170*1d06f6b3SHong Zhang *subts = NULL; 171*1d06f6b3SHong Zhang /* look up the split */ 172*1d06f6b3SHong Zhang ierr = TSRHSSplitGetRHSSplit(ts,splitname,&isplit);CHKERRQ(ierr); 173*1d06f6b3SHong Zhang if (isplit) *subts = isplit->ts; 174*1d06f6b3SHong Zhang PetscFunctionReturn(0); 175*1d06f6b3SHong Zhang } 176*1d06f6b3SHong Zhang 177*1d06f6b3SHong Zhang /*@C 178*1d06f6b3SHong Zhang TSRHSSplitGetSubTSs - Get an array of all sub-TS contexts. 179545aaa6fSHong Zhang 180545aaa6fSHong Zhang Logically Collective on TS 181545aaa6fSHong Zhang 182545aaa6fSHong Zhang Output Parameters: 183545aaa6fSHong Zhang + n - the number of splits 184545aaa6fSHong Zhang - subksp - the array of TS contexts 185545aaa6fSHong Zhang 186545aaa6fSHong Zhang Note: 187545aaa6fSHong Zhang After TSRHSSplitGetSubTS() the array of TSs is to be freed by the user with PetscFree() 188545aaa6fSHong Zhang (not the TS just the array that contains them). 189545aaa6fSHong Zhang 190545aaa6fSHong Zhang Level: advanced 191545aaa6fSHong Zhang 192545aaa6fSHong Zhang .seealso: TSGetRHSSplitFunction() 193545aaa6fSHong Zhang @*/ 194*1d06f6b3SHong Zhang PetscErrorCode TSRHSSplitGetSubTSs(TS ts,PetscInt *n,TS *subts[]) 195545aaa6fSHong Zhang { 196*1d06f6b3SHong Zhang TS_RHSSplitLink ilink = ts->tsrhssplit; 197545aaa6fSHong Zhang PetscInt i = 0; 198545aaa6fSHong Zhang PetscErrorCode ierr; 199545aaa6fSHong Zhang 200545aaa6fSHong Zhang PetscFunctionBegin; 201545aaa6fSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 202545aaa6fSHong Zhang if (subts) { 203545aaa6fSHong Zhang ierr = PetscMalloc1(ts->num_rhs_splits,subts);CHKERRQ(ierr); 204*1d06f6b3SHong Zhang while (ilink) { 205*1d06f6b3SHong Zhang (*subts)[i++] = ilink->ts; 206*1d06f6b3SHong Zhang ilink = ilink->next; 207545aaa6fSHong Zhang } 208545aaa6fSHong Zhang } 209545aaa6fSHong Zhang if (n) *n = ts->num_rhs_splits; 210545aaa6fSHong Zhang PetscFunctionReturn(0); 211545aaa6fSHong Zhang } 212