1d2567f34SHong Zhang /* 2d2567f34SHong Zhang Code for timestepping with implicit Runge-Kutta method 3d2567f34SHong Zhang 4d2567f34SHong Zhang Notes: 5d2567f34SHong Zhang The general system is written as 6d2567f34SHong Zhang 7d2567f34SHong Zhang F(t,U,Udot) = 0 8d2567f34SHong Zhang 9d2567f34SHong Zhang */ 10d2567f34SHong Zhang #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 11d2567f34SHong Zhang #include <petscdm.h> 12d2567f34SHong Zhang #include <petscdt.h> 13d2567f34SHong Zhang 14d2567f34SHong Zhang static TSIRKType TSIRKDefault = TSIRKGAUSS; 15d2567f34SHong Zhang static PetscBool TSIRKRegisterAllCalled; 16d2567f34SHong Zhang static PetscBool TSIRKPackageInitialized; 17d2567f34SHong Zhang static PetscFunctionList TSIRKList; 18d2567f34SHong Zhang 19d2567f34SHong Zhang struct _IRKTableau{ 20d2567f34SHong Zhang PetscReal *A,*b,*c; 21d2567f34SHong Zhang PetscScalar *A_inv,*A_inv_rowsum,*I_s; 22d2567f34SHong Zhang PetscReal *binterp; /* Dense output formula */ 23d2567f34SHong Zhang }; 24d2567f34SHong Zhang 25d2567f34SHong Zhang typedef struct _IRKTableau *IRKTableau; 26d2567f34SHong Zhang 27d2567f34SHong Zhang typedef struct { 28d2567f34SHong Zhang char *method_name; 29d2567f34SHong Zhang PetscInt order; /* Classical approximation order of the method */ 30d2567f34SHong Zhang PetscInt nstages; /* Number of stages */ 31d2567f34SHong Zhang PetscBool stiffly_accurate; 32d2567f34SHong Zhang PetscInt pinterp; /* Interpolation order */ 33d2567f34SHong Zhang IRKTableau tableau; 34d2567f34SHong Zhang Vec U0; /* Backup vector */ 35d2567f34SHong Zhang Vec Z; /* Combined stage vector */ 36d2567f34SHong Zhang Vec *Y; /* States computed during the step */ 37d2567f34SHong Zhang Vec Ydot; /* Work vector holding time derivatives during residual evaluation */ 38d2567f34SHong Zhang Vec U; /* U is used to compute Ydot = shift(Y-U) */ 39d2567f34SHong Zhang Vec *YdotI; /* Work vectors to hold the residual evaluation */ 40d2567f34SHong Zhang Mat TJ; /* KAIJ matrix for the Jacobian of the combined system */ 41d2567f34SHong Zhang PetscScalar *work; /* Scalar work */ 42d2567f34SHong Zhang TSStepStatus status; 43d2567f34SHong Zhang PetscBool rebuild_completion; 44d2567f34SHong Zhang PetscReal ccfl; 45d2567f34SHong Zhang } TS_IRK; 46d2567f34SHong Zhang 47d2567f34SHong Zhang /*@C 48d2567f34SHong Zhang TSIRKTableauCreate - create the tableau for TSIRK and provide the entries 49d2567f34SHong Zhang 50d2567f34SHong Zhang Not Collective 51d2567f34SHong Zhang 52d2567f34SHong Zhang Input Parameters: 53d2567f34SHong Zhang + ts - timestepping context 54d2567f34SHong Zhang . nstages - number of stages, this is the dimension of the matrices below 55d2567f34SHong Zhang . A - stage coefficients (dimension nstages*nstages, row-major) 56d2567f34SHong Zhang . b - step completion table (dimension nstages) 57d2567f34SHong Zhang . c - abscissa (dimension nstages) 5802024617SSatish Balay . binterp - coefficients of the interpolation formula (dimension nstages) 59d2567f34SHong Zhang . A_inv - inverse of A (dimension nstages*nstages, row-major) 60d2567f34SHong Zhang . A_inv_rowsum - row sum of the inverse of A (dimension nstages) 6102024617SSatish Balay - I_s - identity matrix (dimension nstages*nstages) 62d2567f34SHong Zhang 63d2567f34SHong Zhang Level: advanced 64d2567f34SHong Zhang 65d2567f34SHong Zhang .seealso: TSIRK, TSIRKRegister() 66d2567f34SHong Zhang @*/ 67d2567f34SHong Zhang PetscErrorCode TSIRKTableauCreate(TS ts,PetscInt nstages,const PetscReal *A,const PetscReal *b,const PetscReal *c,const PetscReal *binterp,const PetscScalar *A_inv,const PetscScalar *A_inv_rowsum,const PetscScalar *I_s) 68d2567f34SHong Zhang { 69d2567f34SHong Zhang TS_IRK *irk = (TS_IRK*)ts->data; 70d2567f34SHong Zhang IRKTableau tab = irk->tableau; 71d2567f34SHong Zhang PetscErrorCode ierr; 72d2567f34SHong Zhang 73d2567f34SHong Zhang PetscFunctionBegin; 74d2567f34SHong Zhang irk->order = nstages; 75d2567f34SHong Zhang ierr = PetscMalloc3(PetscSqr(nstages),&tab->A,PetscSqr(nstages),&tab->A_inv,PetscSqr(nstages),&tab->I_s);CHKERRQ(ierr); 76d2567f34SHong Zhang ierr = PetscMalloc4(nstages,&tab->b,nstages,&tab->c,nstages,&tab->binterp,nstages,&tab->A_inv_rowsum);CHKERRQ(ierr); 77d2567f34SHong Zhang ierr = PetscArraycpy(tab->A,A,PetscSqr(nstages));CHKERRQ(ierr); 78d2567f34SHong Zhang ierr = PetscArraycpy(tab->b,b,nstages);CHKERRQ(ierr); 79d2567f34SHong Zhang ierr = PetscArraycpy(tab->c,c,nstages);CHKERRQ(ierr); 80d2567f34SHong Zhang /* optional coefficient arrays */ 81d2567f34SHong Zhang if (binterp) { 82d2567f34SHong Zhang ierr = PetscArraycpy(tab->binterp,binterp,nstages);CHKERRQ(ierr); 83d2567f34SHong Zhang } 84d2567f34SHong Zhang if (A_inv) { 85d2567f34SHong Zhang ierr = PetscArraycpy(tab->A_inv,A_inv,PetscSqr(nstages));CHKERRQ(ierr); 86d2567f34SHong Zhang } 87d2567f34SHong Zhang if (A_inv_rowsum) { 88d2567f34SHong Zhang ierr = PetscArraycpy(tab->A_inv_rowsum,A_inv_rowsum,nstages);CHKERRQ(ierr); 89d2567f34SHong Zhang } 90d2567f34SHong Zhang if (I_s) { 91d2567f34SHong Zhang ierr = PetscArraycpy(tab->I_s,I_s,PetscSqr(nstages));CHKERRQ(ierr); 92d2567f34SHong Zhang } 93d2567f34SHong Zhang PetscFunctionReturn(0); 94d2567f34SHong Zhang } 95d2567f34SHong Zhang 96d2567f34SHong Zhang /* Arrays should be freed with PetscFree3(A,b,c) */ 97d2567f34SHong Zhang static PetscErrorCode TSIRKCreate_Gauss(TS ts) 98d2567f34SHong Zhang { 99d2567f34SHong Zhang PetscInt nstages; 100d2567f34SHong Zhang PetscReal *gauss_A_real,*gauss_b,*b,*gauss_c; 101d2567f34SHong Zhang PetscScalar *gauss_A,*gauss_A_inv,*gauss_A_inv_rowsum,*I_s; 102d2567f34SHong Zhang PetscScalar *G0,*G1; 103d2567f34SHong Zhang PetscInt i,j; 104d2567f34SHong Zhang Mat G0mat,G1mat,Amat; 105d2567f34SHong Zhang PetscErrorCode ierr; 106d2567f34SHong Zhang 107d2567f34SHong Zhang PetscFunctionBegin; 108d2567f34SHong Zhang ierr = TSIRKGetNumStages(ts,&nstages);CHKERRQ(ierr); 109d2567f34SHong Zhang ierr = PetscMalloc3(PetscSqr(nstages),&gauss_A_real,nstages,&gauss_b,nstages,&gauss_c);CHKERRQ(ierr); 110d2567f34SHong Zhang ierr = PetscMalloc4(PetscSqr(nstages),&gauss_A,PetscSqr(nstages),&gauss_A_inv,nstages,&gauss_A_inv_rowsum,PetscSqr(nstages),&I_s);CHKERRQ(ierr); 111d2567f34SHong Zhang ierr = PetscMalloc3(nstages,&b,PetscSqr(nstages),&G0,PetscSqr(nstages),&G1);CHKERRQ(ierr); 112d2567f34SHong Zhang ierr = PetscDTGaussQuadrature(nstages,0.,1.,gauss_c,b);CHKERRQ(ierr); 113d2567f34SHong Zhang for (i=0; i<nstages; i++) gauss_b[i] = b[i]; /* copy to possibly-complex array */ 114d2567f34SHong Zhang 115d2567f34SHong Zhang /* A^T = G0^{-1} G1 */ 116d2567f34SHong Zhang for (i=0; i<nstages; i++) { 117d2567f34SHong Zhang for (j=0; j<nstages; j++) { 118d2567f34SHong Zhang G0[i*nstages+j] = PetscPowRealInt(gauss_c[i],j); 119d2567f34SHong Zhang G1[i*nstages+j] = PetscPowRealInt(gauss_c[i],j+1)/(j+1); 120d2567f34SHong Zhang } 121d2567f34SHong Zhang } 122d2567f34SHong Zhang /* The arrays above are row-aligned, but we create dense matrices as the transpose */ 123d2567f34SHong Zhang ierr = MatCreateSeqDense(PETSC_COMM_SELF,nstages,nstages,G0,&G0mat);CHKERRQ(ierr); 124d2567f34SHong Zhang ierr = MatCreateSeqDense(PETSC_COMM_SELF,nstages,nstages,G1,&G1mat);CHKERRQ(ierr); 125d2567f34SHong Zhang ierr = MatCreateSeqDense(PETSC_COMM_SELF,nstages,nstages,gauss_A,&Amat);CHKERRQ(ierr); 126d2567f34SHong Zhang ierr = MatLUFactor(G0mat,NULL,NULL,NULL);CHKERRQ(ierr); 127d2567f34SHong Zhang ierr = MatMatSolve(G0mat,G1mat,Amat);CHKERRQ(ierr); 128d2567f34SHong Zhang ierr = MatTranspose(Amat,MAT_INPLACE_MATRIX,&Amat);CHKERRQ(ierr); 129d2567f34SHong Zhang for (i=0; i<nstages; i++) 130d2567f34SHong Zhang for (j=0; j<nstages; j++) 131d2567f34SHong Zhang gauss_A_real[i*nstages+j] = PetscRealPart(gauss_A[i*nstages+j]); 132d2567f34SHong Zhang 133d2567f34SHong Zhang ierr = MatDestroy(&G0mat);CHKERRQ(ierr); 134d2567f34SHong Zhang ierr = MatDestroy(&G1mat);CHKERRQ(ierr); 135d2567f34SHong Zhang ierr = MatDestroy(&Amat);CHKERRQ(ierr); 136d2567f34SHong Zhang ierr = PetscFree3(b,G0,G1);CHKERRQ(ierr); 137d2567f34SHong Zhang 138d2567f34SHong Zhang {/* Invert A */ 139d2567f34SHong Zhang /* PETSc does not provide a routine to calculate the inverse of a general matrix. 140d2567f34SHong Zhang * To get the inverse of A, we form a sequential BAIJ matrix from it, consisting of a single block with block size 141d2567f34SHong Zhang * equal to the dimension of A, and then use MatInvertBlockDiagonal(). */ 142d2567f34SHong Zhang Mat A_baij; 143d2567f34SHong Zhang PetscInt idxm[1]={0},idxn[1]={0}; 144d2567f34SHong Zhang const PetscScalar *A_inv; 145d2567f34SHong Zhang 146d2567f34SHong Zhang ierr = MatCreateSeqBAIJ(PETSC_COMM_SELF,nstages,nstages,nstages,1,NULL,&A_baij);CHKERRQ(ierr); 147d2567f34SHong Zhang ierr = MatSetOption(A_baij,MAT_ROW_ORIENTED,PETSC_FALSE);CHKERRQ(ierr); 148d2567f34SHong Zhang ierr = MatSetValuesBlocked(A_baij,1,idxm,1,idxn,gauss_A,INSERT_VALUES);CHKERRQ(ierr); 149d2567f34SHong Zhang ierr = MatAssemblyBegin(A_baij,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 150d2567f34SHong Zhang ierr = MatAssemblyEnd(A_baij,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 151d2567f34SHong Zhang ierr = MatInvertBlockDiagonal(A_baij,&A_inv);CHKERRQ(ierr); 152d2567f34SHong Zhang ierr = PetscMemcpy(gauss_A_inv,A_inv,nstages*nstages*sizeof(PetscScalar));CHKERRQ(ierr); 153d2567f34SHong Zhang ierr = MatDestroy(&A_baij);CHKERRQ(ierr); 154d2567f34SHong Zhang } 155d2567f34SHong Zhang 156d2567f34SHong Zhang /* Compute row sums A_inv_rowsum and identity I_s */ 157d2567f34SHong Zhang for (i=0; i<nstages; i++) { 158d2567f34SHong Zhang gauss_A_inv_rowsum[i] = 0; 159d2567f34SHong Zhang for (j=0; j<nstages; j++) { 160d2567f34SHong Zhang gauss_A_inv_rowsum[i] += gauss_A_inv[i+nstages*j]; 161d2567f34SHong Zhang I_s[i+nstages*j] = 1.*(i == j); 162d2567f34SHong Zhang } 163d2567f34SHong Zhang } 164d2567f34SHong Zhang ierr = TSIRKTableauCreate(ts,nstages,gauss_A_real,gauss_b,gauss_c,NULL,gauss_A_inv,gauss_A_inv_rowsum,I_s);CHKERRQ(ierr); 165d2567f34SHong Zhang ierr = PetscFree3(gauss_A_real,gauss_b,gauss_c);CHKERRQ(ierr); 166d2567f34SHong Zhang ierr = PetscFree4(gauss_A,gauss_A_inv,gauss_A_inv_rowsum,I_s);CHKERRQ(ierr); 167d2567f34SHong Zhang PetscFunctionReturn(0); 168d2567f34SHong Zhang } 169d2567f34SHong Zhang 170d2567f34SHong Zhang /*@C 171d2567f34SHong Zhang TSIRKRegister - adds a TSIRK implementation 172d2567f34SHong Zhang 173d2567f34SHong Zhang Not Collective 174d2567f34SHong Zhang 175d2567f34SHong Zhang Input Parameters: 176d2567f34SHong Zhang + sname - name of user-defined IRK scheme 177d2567f34SHong Zhang - function - function to create method context 178d2567f34SHong Zhang 179d2567f34SHong Zhang Notes: 180d2567f34SHong Zhang TSIRKRegister() may be called multiple times to add several user-defined families. 181d2567f34SHong Zhang 182d2567f34SHong Zhang Sample usage: 183d2567f34SHong Zhang .vb 184d2567f34SHong Zhang TSIRKRegister("my_scheme",MySchemeCreate); 185d2567f34SHong Zhang .ve 186d2567f34SHong Zhang 187d2567f34SHong Zhang Then, your scheme can be chosen with the procedural interface via 188d2567f34SHong Zhang $ TSIRKSetType(ts,"my_scheme") 189d2567f34SHong Zhang or at runtime via the option 190d2567f34SHong Zhang $ -ts_irk_type my_scheme 191d2567f34SHong Zhang 192d2567f34SHong Zhang Level: advanced 193d2567f34SHong Zhang 194d2567f34SHong Zhang .seealso: TSIRKRegisterAll() 195d2567f34SHong Zhang @*/ 196d2567f34SHong Zhang PetscErrorCode TSIRKRegister(const char sname[],PetscErrorCode (*function)(TS)) 197d2567f34SHong Zhang { 198d2567f34SHong Zhang PetscErrorCode ierr; 199d2567f34SHong Zhang 200d2567f34SHong Zhang PetscFunctionBegin; 201d2567f34SHong Zhang ierr = TSIRKInitializePackage();CHKERRQ(ierr); 202d2567f34SHong Zhang ierr = PetscFunctionListAdd(&TSIRKList,sname,function);CHKERRQ(ierr); 203d2567f34SHong Zhang PetscFunctionReturn(0); 204d2567f34SHong Zhang } 205d2567f34SHong Zhang 206d2567f34SHong Zhang /*@C 207d2567f34SHong Zhang TSIRKRegisterAll - Registers all of the implicit Runge-Kutta methods in TSIRK 208d2567f34SHong Zhang 209d2567f34SHong Zhang Not Collective, but should be called by all processes which will need the schemes to be registered 210d2567f34SHong Zhang 211d2567f34SHong Zhang Level: advanced 212d2567f34SHong Zhang 213d2567f34SHong Zhang .seealso: TSIRKRegisterDestroy() 214d2567f34SHong Zhang @*/ 215d2567f34SHong Zhang PetscErrorCode TSIRKRegisterAll(void) 216d2567f34SHong Zhang { 217d2567f34SHong Zhang PetscErrorCode ierr; 218d2567f34SHong Zhang 219d2567f34SHong Zhang PetscFunctionBegin; 220d2567f34SHong Zhang if (TSIRKRegisterAllCalled) PetscFunctionReturn(0); 221d2567f34SHong Zhang TSIRKRegisterAllCalled = PETSC_TRUE; 222d2567f34SHong Zhang 223d2567f34SHong Zhang ierr = TSIRKRegister(TSIRKGAUSS,TSIRKCreate_Gauss);CHKERRQ(ierr); 224d2567f34SHong Zhang PetscFunctionReturn(0); 225d2567f34SHong Zhang } 226d2567f34SHong Zhang 227d2567f34SHong Zhang /*@C 228d2567f34SHong Zhang TSIRKRegisterDestroy - Frees the list of schemes that were registered by TSIRKRegister(). 229d2567f34SHong Zhang 230d2567f34SHong Zhang Not Collective 231d2567f34SHong Zhang 232d2567f34SHong Zhang Level: advanced 233d2567f34SHong Zhang 234d2567f34SHong Zhang .seealso: TSIRKRegister(), TSIRKRegisterAll() 235d2567f34SHong Zhang @*/ 236d2567f34SHong Zhang PetscErrorCode TSIRKRegisterDestroy(void) 237d2567f34SHong Zhang { 238d2567f34SHong Zhang PetscFunctionBegin; 239d2567f34SHong Zhang TSIRKRegisterAllCalled = PETSC_FALSE; 240d2567f34SHong Zhang PetscFunctionReturn(0); 241d2567f34SHong Zhang } 242d2567f34SHong Zhang 243d2567f34SHong Zhang /*@C 244d2567f34SHong Zhang TSIRKInitializePackage - This function initializes everything in the TSIRK package. It is called 245d2567f34SHong Zhang from TSInitializePackage(). 246d2567f34SHong Zhang 247d2567f34SHong Zhang Level: developer 248d2567f34SHong Zhang 249d2567f34SHong Zhang .seealso: PetscInitialize() 250d2567f34SHong Zhang @*/ 251d2567f34SHong Zhang PetscErrorCode TSIRKInitializePackage(void) 252d2567f34SHong Zhang { 253d2567f34SHong Zhang PetscErrorCode ierr; 254d2567f34SHong Zhang 255d2567f34SHong Zhang PetscFunctionBegin; 256d2567f34SHong Zhang if (TSIRKPackageInitialized) PetscFunctionReturn(0); 257d2567f34SHong Zhang TSIRKPackageInitialized = PETSC_TRUE; 258d2567f34SHong Zhang ierr = TSIRKRegisterAll();CHKERRQ(ierr); 259d2567f34SHong Zhang ierr = PetscRegisterFinalize(TSIRKFinalizePackage);CHKERRQ(ierr); 260d2567f34SHong Zhang PetscFunctionReturn(0); 261d2567f34SHong Zhang } 262d2567f34SHong Zhang 263d2567f34SHong Zhang /*@C 264d2567f34SHong Zhang TSIRKFinalizePackage - This function destroys everything in the TSIRK package. It is 265d2567f34SHong Zhang called from PetscFinalize(). 266d2567f34SHong Zhang 267d2567f34SHong Zhang Level: developer 268d2567f34SHong Zhang 269d2567f34SHong Zhang .seealso: PetscFinalize() 270d2567f34SHong Zhang @*/ 271d2567f34SHong Zhang PetscErrorCode TSIRKFinalizePackage(void) 272d2567f34SHong Zhang { 273d2567f34SHong Zhang PetscErrorCode ierr; 274d2567f34SHong Zhang 275d2567f34SHong Zhang PetscFunctionBegin; 276d2567f34SHong Zhang ierr = PetscFunctionListDestroy(&TSIRKList);CHKERRQ(ierr); 277d2567f34SHong Zhang TSIRKPackageInitialized = PETSC_FALSE; 278d2567f34SHong Zhang PetscFunctionReturn(0); 279d2567f34SHong Zhang } 280d2567f34SHong Zhang 281d2567f34SHong Zhang /* 282d2567f34SHong Zhang This function can be called before or after ts->vec_sol has been updated. 283d2567f34SHong Zhang */ 284d2567f34SHong Zhang static PetscErrorCode TSEvaluateStep_IRK(TS ts,PetscInt order,Vec U,PetscBool *done) 285d2567f34SHong Zhang { 286d2567f34SHong Zhang TS_IRK *irk = (TS_IRK*)ts->data; 287d2567f34SHong Zhang IRKTableau tab = irk->tableau; 288d2567f34SHong Zhang Vec *YdotI = irk->YdotI; 289d2567f34SHong Zhang PetscScalar *w = irk->work; 290d2567f34SHong Zhang PetscReal h; 291d2567f34SHong Zhang PetscInt j; 292d2567f34SHong Zhang PetscErrorCode ierr; 293d2567f34SHong Zhang 294d2567f34SHong Zhang PetscFunctionBegin; 295d2567f34SHong Zhang switch (irk->status) { 296d2567f34SHong Zhang case TS_STEP_INCOMPLETE: 297d2567f34SHong Zhang case TS_STEP_PENDING: 298d2567f34SHong Zhang h = ts->time_step; break; 299d2567f34SHong Zhang case TS_STEP_COMPLETE: 300d2567f34SHong Zhang h = ts->ptime - ts->ptime_prev; break; 301d2567f34SHong Zhang default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 302d2567f34SHong Zhang } 303d2567f34SHong Zhang 304d2567f34SHong Zhang ierr = VecCopy(ts->vec_sol,U);CHKERRQ(ierr); 305d2567f34SHong Zhang for (j=0; j<irk->nstages; j++) w[j] = h*tab->b[j]; 306d2567f34SHong Zhang ierr = VecMAXPY(U,irk->nstages,w,YdotI);CHKERRQ(ierr); 307d2567f34SHong Zhang PetscFunctionReturn(0); 308d2567f34SHong Zhang } 309d2567f34SHong Zhang 310d2567f34SHong Zhang static PetscErrorCode TSRollBack_IRK(TS ts) 311d2567f34SHong Zhang { 312d2567f34SHong Zhang TS_IRK *irk = (TS_IRK*)ts->data; 313d2567f34SHong Zhang PetscErrorCode ierr; 314d2567f34SHong Zhang 315d2567f34SHong Zhang PetscFunctionBegin; 316d2567f34SHong Zhang ierr = VecCopy(irk->U0,ts->vec_sol);CHKERRQ(ierr); 317d2567f34SHong Zhang PetscFunctionReturn(0); 318d2567f34SHong Zhang } 319d2567f34SHong Zhang 320d2567f34SHong Zhang static PetscErrorCode TSStep_IRK(TS ts) 321d2567f34SHong Zhang { 322d2567f34SHong Zhang TS_IRK *irk = (TS_IRK*)ts->data; 323d2567f34SHong Zhang IRKTableau tab = irk->tableau; 324d2567f34SHong Zhang PetscScalar *A_inv = tab->A_inv,*A_inv_rowsum = tab->A_inv_rowsum; 325d2567f34SHong Zhang const PetscInt nstages = irk->nstages; 326d2567f34SHong Zhang SNES snes; 327d2567f34SHong Zhang PetscInt i,j,its,lits,bs; 328d2567f34SHong Zhang TSAdapt adapt; 329d2567f34SHong Zhang PetscInt rejections = 0; 330d2567f34SHong Zhang PetscBool accept = PETSC_TRUE; 331d2567f34SHong Zhang PetscReal next_time_step = ts->time_step; 332d2567f34SHong Zhang PetscErrorCode ierr; 333d2567f34SHong Zhang 334d2567f34SHong Zhang PetscFunctionBegin; 335d2567f34SHong Zhang if (!ts->steprollback) { 336d2567f34SHong Zhang ierr = VecCopy(ts->vec_sol,irk->U0);CHKERRQ(ierr); 337d2567f34SHong Zhang } 338d2567f34SHong Zhang ierr = VecGetBlockSize(ts->vec_sol,&bs);CHKERRQ(ierr); 339d2567f34SHong Zhang for (i=0; i<nstages; i++) { 340d2567f34SHong Zhang ierr = VecStrideScatter(ts->vec_sol,i*bs,irk->Z,INSERT_VALUES);CHKERRQ(ierr); 341d2567f34SHong Zhang } 342d2567f34SHong Zhang 343d2567f34SHong Zhang irk->status = TS_STEP_INCOMPLETE; 344d2567f34SHong Zhang while (!ts->reason && irk->status != TS_STEP_COMPLETE) { 345d2567f34SHong Zhang ierr = VecCopy(ts->vec_sol,irk->U);CHKERRQ(ierr); 346d2567f34SHong Zhang ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 347d2567f34SHong Zhang ierr = SNESSolve(snes,NULL,irk->Z);CHKERRQ(ierr); 348d2567f34SHong Zhang ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr); 349d2567f34SHong Zhang ierr = SNESGetLinearSolveIterations(snes,&lits);CHKERRQ(ierr); 350d2567f34SHong Zhang ts->snes_its += its; ts->ksp_its += lits; 351d2567f34SHong Zhang ierr = VecStrideGatherAll(irk->Z,irk->Y,INSERT_VALUES);CHKERRQ(ierr); 352d2567f34SHong Zhang for (i=0; i<nstages; i++) { 353d2567f34SHong Zhang ierr = VecZeroEntries(irk->YdotI[i]);CHKERRQ(ierr); 354d2567f34SHong Zhang for (j=0; j<nstages; j++) { 355d2567f34SHong Zhang ierr = VecAXPY(irk->YdotI[i],A_inv[i+j*nstages]/ts->time_step,irk->Y[j]);CHKERRQ(ierr); 356d2567f34SHong Zhang } 357d2567f34SHong Zhang ierr = VecAXPY(irk->YdotI[i],-A_inv_rowsum[i]/ts->time_step,irk->U);CHKERRQ(ierr); 358d2567f34SHong Zhang } 359d2567f34SHong Zhang irk->status = TS_STEP_INCOMPLETE; 360d2567f34SHong Zhang ierr = TSEvaluateStep_IRK(ts,irk->order,ts->vec_sol,NULL);CHKERRQ(ierr); 361d2567f34SHong Zhang irk->status = TS_STEP_PENDING; 362d2567f34SHong Zhang ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 363d2567f34SHong Zhang ierr = TSAdaptChoose(adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 364d2567f34SHong Zhang irk->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 365d2567f34SHong Zhang if (!accept) { 366d2567f34SHong Zhang ierr = TSRollBack_IRK(ts);CHKERRQ(ierr); 367d2567f34SHong Zhang ts->time_step = next_time_step; 368d2567f34SHong Zhang goto reject_step; 369d2567f34SHong Zhang } 370d2567f34SHong Zhang 371d2567f34SHong Zhang ts->ptime += ts->time_step; 372d2567f34SHong Zhang ts->time_step = next_time_step; 373d2567f34SHong Zhang break; 374d2567f34SHong Zhang reject_step: 375d2567f34SHong Zhang ts->reject++; accept = PETSC_FALSE; 376d2567f34SHong Zhang if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 377d2567f34SHong Zhang ts->reason = TS_DIVERGED_STEP_REJECTED; 378d2567f34SHong Zhang ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 379d2567f34SHong Zhang } 380d2567f34SHong Zhang } 381d2567f34SHong Zhang PetscFunctionReturn(0); 382d2567f34SHong Zhang } 383d2567f34SHong Zhang 384d2567f34SHong Zhang static PetscErrorCode TSInterpolate_IRK(TS ts,PetscReal itime,Vec U) 385d2567f34SHong Zhang { 386d2567f34SHong Zhang TS_IRK *irk = (TS_IRK*)ts->data; 387d2567f34SHong Zhang PetscInt nstages = irk->nstages,pinterp = irk->pinterp,i,j; 388d2567f34SHong Zhang PetscReal h; 389d2567f34SHong Zhang PetscReal tt,t; 390d2567f34SHong Zhang PetscScalar *bt; 391d2567f34SHong Zhang const PetscReal *B = irk->tableau->binterp; 392d2567f34SHong Zhang PetscErrorCode ierr; 393d2567f34SHong Zhang 394d2567f34SHong Zhang PetscFunctionBegin; 395*98921bdaSJacob Faibussowitsch if (!B) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSIRK %s does not have an interpolation formula",irk->method_name); 396d2567f34SHong Zhang switch (irk->status) { 397d2567f34SHong Zhang case TS_STEP_INCOMPLETE: 398d2567f34SHong Zhang case TS_STEP_PENDING: 399d2567f34SHong Zhang h = ts->time_step; 400d2567f34SHong Zhang t = (itime - ts->ptime)/h; 401d2567f34SHong Zhang break; 402d2567f34SHong Zhang case TS_STEP_COMPLETE: 403d2567f34SHong Zhang h = ts->ptime - ts->ptime_prev; 404d2567f34SHong Zhang t = (itime - ts->ptime)/h + 1; /* In the interval [0,1] */ 405d2567f34SHong Zhang break; 406d2567f34SHong Zhang default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Invalid TSStepStatus"); 407d2567f34SHong Zhang } 408d2567f34SHong Zhang ierr = PetscMalloc1(nstages,&bt);CHKERRQ(ierr); 409d2567f34SHong Zhang for (i=0; i<nstages; i++) bt[i] = 0; 410d2567f34SHong Zhang for (j=0,tt=t; j<pinterp; j++,tt*=t) { 411d2567f34SHong Zhang for (i=0; i<nstages; i++) { 412d2567f34SHong Zhang bt[i] += h * B[i*pinterp+j] * tt; 413d2567f34SHong Zhang } 414d2567f34SHong Zhang } 415d2567f34SHong Zhang ierr = VecMAXPY(U,nstages,bt,irk->YdotI);CHKERRQ(ierr); 416d2567f34SHong Zhang PetscFunctionReturn(0); 417d2567f34SHong Zhang } 418d2567f34SHong Zhang 419d2567f34SHong Zhang static PetscErrorCode TSIRKTableauReset(TS ts) 420d2567f34SHong Zhang { 421d2567f34SHong Zhang TS_IRK *irk = (TS_IRK*)ts->data; 422d2567f34SHong Zhang IRKTableau tab = irk->tableau; 423d2567f34SHong Zhang PetscErrorCode ierr; 424d2567f34SHong Zhang 425d2567f34SHong Zhang PetscFunctionBegin; 426d2567f34SHong Zhang if (!tab) PetscFunctionReturn(0); 427d2567f34SHong Zhang ierr = PetscFree3(tab->A,tab->A_inv,tab->I_s);CHKERRQ(ierr); 428d2567f34SHong Zhang ierr = PetscFree4(tab->b,tab->c,tab->binterp,tab->A_inv_rowsum);CHKERRQ(ierr); 429d2567f34SHong Zhang PetscFunctionReturn(0); 430d2567f34SHong Zhang } 431d2567f34SHong Zhang 432d2567f34SHong Zhang static PetscErrorCode TSReset_IRK(TS ts) 433d2567f34SHong Zhang { 434d2567f34SHong Zhang TS_IRK *irk = (TS_IRK*)ts->data; 435d2567f34SHong Zhang PetscErrorCode ierr; 436d2567f34SHong Zhang 437d2567f34SHong Zhang PetscFunctionBegin; 438d2567f34SHong Zhang ierr = TSIRKTableauReset(ts);CHKERRQ(ierr); 439d2567f34SHong Zhang if (irk->tableau) { 440d2567f34SHong Zhang ierr = PetscFree(irk->tableau);CHKERRQ(ierr); 441d2567f34SHong Zhang } 442d2567f34SHong Zhang if (irk->method_name) { 443d2567f34SHong Zhang ierr = PetscFree(irk->method_name);CHKERRQ(ierr); 444d2567f34SHong Zhang } 445d2567f34SHong Zhang if (irk->work) { 446d2567f34SHong Zhang ierr = PetscFree(irk->work);CHKERRQ(ierr); 447d2567f34SHong Zhang } 448d2567f34SHong Zhang ierr = VecDestroyVecs(irk->nstages,&irk->Y);CHKERRQ(ierr); 449d2567f34SHong Zhang ierr = VecDestroyVecs(irk->nstages,&irk->YdotI);CHKERRQ(ierr); 450d2567f34SHong Zhang ierr = VecDestroy(&irk->Ydot);CHKERRQ(ierr); 451d2567f34SHong Zhang ierr = VecDestroy(&irk->Z);CHKERRQ(ierr); 452d2567f34SHong Zhang ierr = VecDestroy(&irk->U);CHKERRQ(ierr); 453d2567f34SHong Zhang ierr = VecDestroy(&irk->U0);CHKERRQ(ierr); 454d2567f34SHong Zhang ierr = MatDestroy(&irk->TJ);CHKERRQ(ierr); 455d2567f34SHong Zhang PetscFunctionReturn(0); 456d2567f34SHong Zhang } 457d2567f34SHong Zhang 458d2567f34SHong Zhang static PetscErrorCode TSIRKGetVecs(TS ts,DM dm,Vec *U) 459d2567f34SHong Zhang { 460d2567f34SHong Zhang TS_IRK *irk = (TS_IRK*)ts->data; 461d2567f34SHong Zhang PetscErrorCode ierr; 462d2567f34SHong Zhang 463d2567f34SHong Zhang PetscFunctionBegin; 464d2567f34SHong Zhang if (U) { 465d2567f34SHong Zhang if (dm && dm != ts->dm) { 466d2567f34SHong Zhang ierr = DMGetNamedGlobalVector(dm,"TSIRK_U",U);CHKERRQ(ierr); 467d2567f34SHong Zhang } else *U = irk->U; 468d2567f34SHong Zhang } 469d2567f34SHong Zhang PetscFunctionReturn(0); 470d2567f34SHong Zhang } 471d2567f34SHong Zhang 472d2567f34SHong Zhang static PetscErrorCode TSIRKRestoreVecs(TS ts,DM dm,Vec *U) 473d2567f34SHong Zhang { 474d2567f34SHong Zhang PetscErrorCode ierr; 475d2567f34SHong Zhang 476d2567f34SHong Zhang PetscFunctionBegin; 477d2567f34SHong Zhang if (U) { 478d2567f34SHong Zhang if (dm && dm != ts->dm) { 479d2567f34SHong Zhang ierr = DMRestoreNamedGlobalVector(dm,"TSIRK_U",U);CHKERRQ(ierr); 480d2567f34SHong Zhang } 481d2567f34SHong Zhang } 482d2567f34SHong Zhang PetscFunctionReturn(0); 483d2567f34SHong Zhang } 484d2567f34SHong Zhang 485d2567f34SHong Zhang /* 486d2567f34SHong Zhang This defines the nonlinear equations that is to be solved with SNES 487d2567f34SHong Zhang G[e\otimes t + C*dt, Z, Zdot] = 0 488d2567f34SHong Zhang Zdot = (In \otimes S)*Z - (In \otimes Se) U 489d2567f34SHong Zhang where S = 1/(dt*A) 490d2567f34SHong Zhang */ 491d2567f34SHong Zhang static PetscErrorCode SNESTSFormFunction_IRK(SNES snes,Vec ZC,Vec FC,TS ts) 492d2567f34SHong Zhang { 493d2567f34SHong Zhang TS_IRK *irk = (TS_IRK*)ts->data; 494d2567f34SHong Zhang IRKTableau tab = irk->tableau; 495d2567f34SHong Zhang const PetscInt nstages = irk->nstages; 496d2567f34SHong Zhang const PetscReal *c = tab->c; 497d2567f34SHong Zhang const PetscScalar *A_inv = tab->A_inv,*A_inv_rowsum = tab->A_inv_rowsum; 498d2567f34SHong Zhang DM dm,dmsave; 499d2567f34SHong Zhang Vec U,*YdotI = irk->YdotI,Ydot = irk->Ydot,*Y = irk->Y; 500d2567f34SHong Zhang PetscReal h = ts->time_step; 501d2567f34SHong Zhang PetscInt i,j; 502d2567f34SHong Zhang PetscErrorCode ierr; 503d2567f34SHong Zhang 504d2567f34SHong Zhang PetscFunctionBegin; 505d2567f34SHong Zhang ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 506d2567f34SHong Zhang ierr = TSIRKGetVecs(ts,dm,&U);CHKERRQ(ierr); 507d2567f34SHong Zhang ierr = VecStrideGatherAll(ZC,Y,INSERT_VALUES);CHKERRQ(ierr); 508d2567f34SHong Zhang dmsave = ts->dm; 509d2567f34SHong Zhang ts->dm = dm; 510d2567f34SHong Zhang for (i=0; i<nstages; i++) { 511d2567f34SHong Zhang ierr = VecZeroEntries(Ydot);CHKERRQ(ierr); 512d2567f34SHong Zhang for (j=0; j<nstages; j++) { 513d2567f34SHong Zhang ierr = VecAXPY(Ydot,A_inv[j*nstages+i]/h,Y[j]);CHKERRQ(ierr); 514d2567f34SHong Zhang } 515d2567f34SHong Zhang ierr = VecAXPY(Ydot,-A_inv_rowsum[i]/h,U);CHKERRQ(ierr); /* Ydot = (S \otimes In)*Z - (Se \otimes In) U */ 516d2567f34SHong Zhang ierr = TSComputeIFunction(ts,ts->ptime+ts->time_step*c[i],Y[i],Ydot,YdotI[i],PETSC_FALSE);CHKERRQ(ierr); 517d2567f34SHong Zhang } 518d2567f34SHong Zhang ierr = VecStrideScatterAll(YdotI,FC,INSERT_VALUES);CHKERRQ(ierr); 519d2567f34SHong Zhang ts->dm = dmsave; 520d2567f34SHong Zhang ierr = TSIRKRestoreVecs(ts,dm,&U);CHKERRQ(ierr); 521d2567f34SHong Zhang PetscFunctionReturn(0); 522d2567f34SHong Zhang } 523d2567f34SHong Zhang 524d2567f34SHong Zhang /* 525d2567f34SHong Zhang For explicit ODE, the Jacobian is 526d2567f34SHong Zhang JC = I_n \otimes S - J \otimes I_s 527d2567f34SHong Zhang For DAE, the Jacobian is 528d2567f34SHong Zhang JC = M_n \otimes S - J \otimes I_s 529d2567f34SHong Zhang */ 530d2567f34SHong Zhang static PetscErrorCode SNESTSFormJacobian_IRK(SNES snes,Vec ZC,Mat JC,Mat JCpre,TS ts) 531d2567f34SHong Zhang { 532d2567f34SHong Zhang TS_IRK *irk = (TS_IRK*)ts->data; 533d2567f34SHong Zhang IRKTableau tab = irk->tableau; 534d2567f34SHong Zhang const PetscInt nstages = irk->nstages; 535d2567f34SHong Zhang const PetscReal *c = tab->c; 536d2567f34SHong Zhang DM dm,dmsave; 537d2567f34SHong Zhang Vec *Y = irk->Y,Ydot = irk->Ydot; 538d2567f34SHong Zhang Mat J; 539d2567f34SHong Zhang PetscScalar *S; 540d2567f34SHong Zhang PetscInt i,j,bs; 541d2567f34SHong Zhang PetscErrorCode ierr; 542d2567f34SHong Zhang 543d2567f34SHong Zhang PetscFunctionBegin; 544d2567f34SHong Zhang ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 545d2567f34SHong Zhang /* irk->Ydot has already been computed in SNESTSFormFunction_IRK (SNES guarantees this) */ 546d2567f34SHong Zhang dmsave = ts->dm; 547d2567f34SHong Zhang ts->dm = dm; 548d2567f34SHong Zhang ierr = VecGetBlockSize(Y[nstages-1],&bs);CHKERRQ(ierr); 549d2567f34SHong Zhang if (ts->equation_type <= TS_EQ_ODE_EXPLICIT) { /* Support explicit formulas only */ 550d2567f34SHong Zhang ierr = VecStrideGather(ZC,(nstages-1)*bs,Y[nstages-1],INSERT_VALUES);CHKERRQ(ierr); 551d2567f34SHong Zhang ierr = MatKAIJGetAIJ(JC,&J);CHKERRQ(ierr); 552d2567f34SHong Zhang ierr = TSComputeIJacobian(ts,ts->ptime+ts->time_step*c[nstages-1],Y[nstages-1],Ydot,0,J,J,PETSC_FALSE);CHKERRQ(ierr); 553d2567f34SHong Zhang ierr = MatKAIJGetS(JC,NULL,NULL,&S);CHKERRQ(ierr); 554d2567f34SHong Zhang for (i=0; i<nstages; i++) 555d2567f34SHong Zhang for (j=0; j<nstages; j++) 556d2567f34SHong Zhang S[i+nstages*j] = tab->A_inv[i+nstages*j]/ts->time_step; 557d2567f34SHong Zhang ierr = MatKAIJRestoreS(JC,&S);CHKERRQ(ierr); 558*98921bdaSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"TSIRK %s does not support implicit formula",irk->method_name); /* TODO: need the mass matrix for DAE */ 559d2567f34SHong Zhang ts->dm = dmsave; 560d2567f34SHong Zhang PetscFunctionReturn(0); 561d2567f34SHong Zhang } 562d2567f34SHong Zhang 563d2567f34SHong Zhang static PetscErrorCode DMCoarsenHook_TSIRK(DM fine,DM coarse,void *ctx) 564d2567f34SHong Zhang { 565d2567f34SHong Zhang PetscFunctionBegin; 566d2567f34SHong Zhang PetscFunctionReturn(0); 567d2567f34SHong Zhang } 568d2567f34SHong Zhang 569d2567f34SHong Zhang static PetscErrorCode DMRestrictHook_TSIRK(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 570d2567f34SHong Zhang { 571d2567f34SHong Zhang TS ts = (TS)ctx; 572d2567f34SHong Zhang Vec U,U_c; 573d2567f34SHong Zhang PetscErrorCode ierr; 574d2567f34SHong Zhang 575d2567f34SHong Zhang PetscFunctionBegin; 576d2567f34SHong Zhang ierr = TSIRKGetVecs(ts,fine,&U);CHKERRQ(ierr); 577d2567f34SHong Zhang ierr = TSIRKGetVecs(ts,coarse,&U_c);CHKERRQ(ierr); 578d2567f34SHong Zhang ierr = MatRestrict(restrct,U,U_c);CHKERRQ(ierr); 579d2567f34SHong Zhang ierr = VecPointwiseMult(U_c,rscale,U_c);CHKERRQ(ierr); 580d2567f34SHong Zhang ierr = TSIRKRestoreVecs(ts,fine,&U);CHKERRQ(ierr); 581d2567f34SHong Zhang ierr = TSIRKRestoreVecs(ts,coarse,&U_c);CHKERRQ(ierr); 582d2567f34SHong Zhang PetscFunctionReturn(0); 583d2567f34SHong Zhang } 584d2567f34SHong Zhang 585d2567f34SHong Zhang static PetscErrorCode DMSubDomainHook_TSIRK(DM dm,DM subdm,void *ctx) 586d2567f34SHong Zhang { 587d2567f34SHong Zhang PetscFunctionBegin; 588d2567f34SHong Zhang PetscFunctionReturn(0); 589d2567f34SHong Zhang } 590d2567f34SHong Zhang 591d2567f34SHong Zhang static PetscErrorCode DMSubDomainRestrictHook_TSIRK(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 592d2567f34SHong Zhang { 593d2567f34SHong Zhang TS ts = (TS)ctx; 594d2567f34SHong Zhang Vec U,U_c; 595d2567f34SHong Zhang PetscErrorCode ierr; 596d2567f34SHong Zhang 597d2567f34SHong Zhang PetscFunctionBegin; 598d2567f34SHong Zhang ierr = TSIRKGetVecs(ts,dm,&U);CHKERRQ(ierr); 599d2567f34SHong Zhang ierr = TSIRKGetVecs(ts,subdm,&U_c);CHKERRQ(ierr); 600d2567f34SHong Zhang 601d2567f34SHong Zhang ierr = VecScatterBegin(gscat,U,U_c,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 602d2567f34SHong Zhang ierr = VecScatterEnd(gscat,U,U_c,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 603d2567f34SHong Zhang 604d2567f34SHong Zhang ierr = TSIRKRestoreVecs(ts,dm,&U);CHKERRQ(ierr); 605d2567f34SHong Zhang ierr = TSIRKRestoreVecs(ts,subdm,&U_c);CHKERRQ(ierr); 606d2567f34SHong Zhang PetscFunctionReturn(0); 607d2567f34SHong Zhang } 608d2567f34SHong Zhang 609d2567f34SHong Zhang static PetscErrorCode TSSetUp_IRK(TS ts) 610d2567f34SHong Zhang { 611d2567f34SHong Zhang TS_IRK *irk = (TS_IRK*)ts->data; 612d2567f34SHong Zhang IRKTableau tab = irk->tableau; 613d2567f34SHong Zhang DM dm; 614d2567f34SHong Zhang Mat J; 615d2567f34SHong Zhang Vec R; 616d2567f34SHong Zhang const PetscInt nstages = irk->nstages; 617d2567f34SHong Zhang PetscInt vsize,bs; 618d2567f34SHong Zhang PetscErrorCode ierr; 619d2567f34SHong Zhang 620d2567f34SHong Zhang PetscFunctionBegin; 621d2567f34SHong Zhang if (!irk->work) { 622d2567f34SHong Zhang ierr = PetscMalloc1(irk->nstages,&irk->work);CHKERRQ(ierr); 623d2567f34SHong Zhang } 624d2567f34SHong Zhang if (!irk->Y) { 625d2567f34SHong Zhang ierr = VecDuplicateVecs(ts->vec_sol,irk->nstages,&irk->Y);CHKERRQ(ierr); 626d2567f34SHong Zhang } 627d2567f34SHong Zhang if (!irk->YdotI) { 628d2567f34SHong Zhang ierr = VecDuplicateVecs(ts->vec_sol,irk->nstages,&irk->YdotI);CHKERRQ(ierr); 629d2567f34SHong Zhang } 630d2567f34SHong Zhang if (!irk->Ydot) { 631d2567f34SHong Zhang ierr = VecDuplicate(ts->vec_sol,&irk->Ydot);CHKERRQ(ierr); 632d2567f34SHong Zhang } 633d2567f34SHong Zhang if (!irk->U) { 634d2567f34SHong Zhang ierr = VecDuplicate(ts->vec_sol,&irk->U);CHKERRQ(ierr); 635d2567f34SHong Zhang } 636d2567f34SHong Zhang if (!irk->U0) { 637d2567f34SHong Zhang ierr = VecDuplicate(ts->vec_sol,&irk->U0);CHKERRQ(ierr); 638d2567f34SHong Zhang } 639d2567f34SHong Zhang if (!irk->Z) { 640d2567f34SHong Zhang ierr = VecCreate(PetscObjectComm((PetscObject)ts->vec_sol),&irk->Z);CHKERRQ(ierr); 641d2567f34SHong Zhang ierr = VecGetSize(ts->vec_sol,&vsize);CHKERRQ(ierr); 642d2567f34SHong Zhang ierr = VecSetSizes(irk->Z,PETSC_DECIDE,vsize*irk->nstages);CHKERRQ(ierr); 643d2567f34SHong Zhang ierr = VecGetBlockSize(ts->vec_sol,&bs);CHKERRQ(ierr); 644d2567f34SHong Zhang ierr = VecSetBlockSize(irk->Z,irk->nstages*bs);CHKERRQ(ierr); 645d2567f34SHong Zhang ierr = VecSetFromOptions(irk->Z);CHKERRQ(ierr); 646d2567f34SHong Zhang } 647d2567f34SHong Zhang ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 648d2567f34SHong Zhang ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_TSIRK,DMRestrictHook_TSIRK,ts);CHKERRQ(ierr); 649d2567f34SHong Zhang ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_TSIRK,DMSubDomainRestrictHook_TSIRK,ts);CHKERRQ(ierr); 650d2567f34SHong Zhang 651d2567f34SHong Zhang ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr); 652d2567f34SHong Zhang ierr = VecDuplicate(irk->Z,&R);CHKERRQ(ierr); 653d2567f34SHong Zhang ierr = SNESSetFunction(ts->snes,R,SNESTSFormFunction,ts);CHKERRQ(ierr); 6543397a780SHong Zhang ierr = TSGetIJacobian(ts,&J,NULL,NULL,NULL);CHKERRQ(ierr); 655d2567f34SHong Zhang if (!irk->TJ) { 656d2567f34SHong Zhang /* Create the KAIJ matrix for solving the stages */ 657d2567f34SHong Zhang ierr = MatCreateKAIJ(J,nstages,nstages,tab->A_inv,tab->I_s,&irk->TJ);CHKERRQ(ierr); 658d2567f34SHong Zhang } 659d2567f34SHong Zhang ierr = SNESSetJacobian(ts->snes,irk->TJ,irk->TJ,SNESTSFormJacobian,ts);CHKERRQ(ierr); 660d2567f34SHong Zhang ierr = VecDestroy(&R);CHKERRQ(ierr); 661d2567f34SHong Zhang PetscFunctionReturn(0); 662d2567f34SHong Zhang } 663d2567f34SHong Zhang 664d2567f34SHong Zhang static PetscErrorCode TSSetFromOptions_IRK(PetscOptionItems *PetscOptionsObject,TS ts) 665d2567f34SHong Zhang { 666d2567f34SHong Zhang TS_IRK *irk = (TS_IRK*)ts->data; 667d2567f34SHong Zhang char tname[256] = TSIRKGAUSS; 668d2567f34SHong Zhang PetscErrorCode ierr; 669d2567f34SHong Zhang 670d2567f34SHong Zhang PetscFunctionBegin; 671d2567f34SHong Zhang ierr = PetscOptionsHead(PetscOptionsObject,"IRK ODE solver options");CHKERRQ(ierr); 672d2567f34SHong Zhang { 673d2567f34SHong Zhang PetscBool flg1,flg2; 674d2567f34SHong Zhang ierr = PetscOptionsInt("-ts_irk_nstages","Stages of the IRK method","TSIRKSetNumStages",irk->nstages,&irk->nstages,&flg1);CHKERRQ(ierr); 675d2567f34SHong Zhang ierr = PetscOptionsFList("-ts_irk_type","Type of IRK method","TSIRKSetType",TSIRKList,irk->method_name[0] ? irk->method_name : tname,tname,sizeof(tname),&flg2);CHKERRQ(ierr); 676d2567f34SHong Zhang if (flg1 ||flg2 || !irk->method_name[0]) { /* Create the method tableau after nstages or method is set */ 677d2567f34SHong Zhang ierr = TSIRKSetType(ts,tname);CHKERRQ(ierr); 678d2567f34SHong Zhang } 679d2567f34SHong Zhang } 680d2567f34SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 681d2567f34SHong Zhang PetscFunctionReturn(0); 682d2567f34SHong Zhang } 683d2567f34SHong Zhang 684d2567f34SHong Zhang static PetscErrorCode TSView_IRK(TS ts,PetscViewer viewer) 685d2567f34SHong Zhang { 686d2567f34SHong Zhang TS_IRK *irk = (TS_IRK*)ts->data; 687d2567f34SHong Zhang PetscBool iascii; 688d2567f34SHong Zhang PetscErrorCode ierr; 689d2567f34SHong Zhang 690d2567f34SHong Zhang PetscFunctionBegin; 691d2567f34SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 692d2567f34SHong Zhang if (iascii) { 693d2567f34SHong Zhang IRKTableau tab = irk->tableau; 694d2567f34SHong Zhang TSIRKType irktype; 695d2567f34SHong Zhang char buf[512]; 696d2567f34SHong Zhang 697d2567f34SHong Zhang ierr = TSIRKGetType(ts,&irktype);CHKERRQ(ierr); 698d2567f34SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," IRK type %s\n",irktype);CHKERRQ(ierr); 699d2567f34SHong Zhang ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",irk->nstages,tab->c);CHKERRQ(ierr); 700d2567f34SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Abscissa c = %s\n",buf);CHKERRQ(ierr); 701d2567f34SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"Stiffly accurate: %s\n",irk->stiffly_accurate ? "yes" : "no");CHKERRQ(ierr); 702d2567f34SHong Zhang ierr = PetscFormatRealArray(buf,sizeof(buf),"% 8.6f",PetscSqr(irk->nstages),tab->A);CHKERRQ(ierr); 703d2567f34SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," A coefficients A = %s\n",buf);CHKERRQ(ierr); 704d2567f34SHong Zhang } 705d2567f34SHong Zhang PetscFunctionReturn(0); 706d2567f34SHong Zhang } 707d2567f34SHong Zhang 708d2567f34SHong Zhang static PetscErrorCode TSLoad_IRK(TS ts,PetscViewer viewer) 709d2567f34SHong Zhang { 710d2567f34SHong Zhang PetscErrorCode ierr; 711d2567f34SHong Zhang SNES snes; 712d2567f34SHong Zhang TSAdapt adapt; 713d2567f34SHong Zhang 714d2567f34SHong Zhang PetscFunctionBegin; 715d2567f34SHong Zhang ierr = TSGetAdapt(ts,&adapt);CHKERRQ(ierr); 716d2567f34SHong Zhang ierr = TSAdaptLoad(adapt,viewer);CHKERRQ(ierr); 717d2567f34SHong Zhang ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 718d2567f34SHong Zhang ierr = SNESLoad(snes,viewer);CHKERRQ(ierr); 719d2567f34SHong Zhang /* function and Jacobian context for SNES when used with TS is always ts object */ 720d2567f34SHong Zhang ierr = SNESSetFunction(snes,NULL,NULL,ts);CHKERRQ(ierr); 721d2567f34SHong Zhang ierr = SNESSetJacobian(snes,NULL,NULL,NULL,ts);CHKERRQ(ierr); 722d2567f34SHong Zhang PetscFunctionReturn(0); 723d2567f34SHong Zhang } 724d2567f34SHong Zhang 725d2567f34SHong Zhang /*@C 726d2567f34SHong Zhang TSIRKSetType - Set the type of IRK scheme 727d2567f34SHong Zhang 728d2567f34SHong Zhang Logically collective 729d2567f34SHong Zhang 73002024617SSatish Balay Input Parameters: 731d2567f34SHong Zhang + ts - timestepping context 732d2567f34SHong Zhang - irktype - type of IRK scheme 733d2567f34SHong Zhang 734d2567f34SHong Zhang Options Database: 735d2567f34SHong Zhang . -ts_irk_type <gauss> 736d2567f34SHong Zhang 737d2567f34SHong Zhang Level: intermediate 738d2567f34SHong Zhang 739d2567f34SHong Zhang .seealso: TSIRKGetType(), TSIRK, TSIRKType, TSIRKGAUSS 740d2567f34SHong Zhang @*/ 741d2567f34SHong Zhang PetscErrorCode TSIRKSetType(TS ts,TSIRKType irktype) 742d2567f34SHong Zhang { 743d2567f34SHong Zhang PetscErrorCode ierr; 744d2567f34SHong Zhang 745d2567f34SHong Zhang PetscFunctionBegin; 746d2567f34SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 747d2567f34SHong Zhang PetscValidCharPointer(irktype,2); 748d2567f34SHong Zhang ierr = PetscTryMethod(ts,"TSIRKSetType_C",(TS,TSIRKType),(ts,irktype));CHKERRQ(ierr); 749d2567f34SHong Zhang PetscFunctionReturn(0); 750d2567f34SHong Zhang } 751d2567f34SHong Zhang 752d2567f34SHong Zhang /*@C 753d2567f34SHong Zhang TSIRKGetType - Get the type of IRK IMEX scheme 754d2567f34SHong Zhang 755d2567f34SHong Zhang Logically collective 756d2567f34SHong Zhang 757d2567f34SHong Zhang Input Parameter: 758d2567f34SHong Zhang . ts - timestepping context 759d2567f34SHong Zhang 760d2567f34SHong Zhang Output Parameter: 761d2567f34SHong Zhang . irktype - type of IRK-IMEX scheme 762d2567f34SHong Zhang 763d2567f34SHong Zhang Level: intermediate 764d2567f34SHong Zhang 765d2567f34SHong Zhang .seealso: TSIRKGetType() 766d2567f34SHong Zhang @*/ 767d2567f34SHong Zhang PetscErrorCode TSIRKGetType(TS ts,TSIRKType *irktype) 768d2567f34SHong Zhang { 769d2567f34SHong Zhang PetscErrorCode ierr; 770d2567f34SHong Zhang 771d2567f34SHong Zhang PetscFunctionBegin; 772d2567f34SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 773d2567f34SHong Zhang ierr = PetscUseMethod(ts,"TSIRKGetType_C",(TS,TSIRKType*),(ts,irktype));CHKERRQ(ierr); 774d2567f34SHong Zhang PetscFunctionReturn(0); 775d2567f34SHong Zhang } 776d2567f34SHong Zhang 777d2567f34SHong Zhang /*@C 778d2567f34SHong Zhang TSIRKSetNumStages - Set the number of stages of IRK scheme 779d2567f34SHong Zhang 780d2567f34SHong Zhang Logically collective 781d2567f34SHong Zhang 78202024617SSatish Balay Input Parameters: 783d2567f34SHong Zhang + ts - timestepping context 784d2567f34SHong Zhang - nstages - number of stages of IRK scheme 785d2567f34SHong Zhang 786d2567f34SHong Zhang Options Database: 787d2567f34SHong Zhang . -ts_irk_nstages <int>: Number of stages 788d2567f34SHong Zhang 789d2567f34SHong Zhang Level: intermediate 790d2567f34SHong Zhang 791d2567f34SHong Zhang .seealso: TSIRKGetNumStages(), TSIRK 792d2567f34SHong Zhang @*/ 793d2567f34SHong Zhang PetscErrorCode TSIRKSetNumStages(TS ts,PetscInt nstages) 794d2567f34SHong Zhang { 795d2567f34SHong Zhang PetscErrorCode ierr; 796d2567f34SHong Zhang 797d2567f34SHong Zhang PetscFunctionBegin; 798d2567f34SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 799d2567f34SHong Zhang ierr = PetscTryMethod(ts,"TSIRKSetNumStages_C",(TS,PetscInt),(ts,nstages));CHKERRQ(ierr); 800d2567f34SHong Zhang PetscFunctionReturn(0); 801d2567f34SHong Zhang } 802d2567f34SHong Zhang 803d2567f34SHong Zhang /*@C 804d2567f34SHong Zhang TSIRKGetNumStages - Get the number of stages of IRK scheme 805d2567f34SHong Zhang 806d2567f34SHong Zhang Logically collective 807d2567f34SHong Zhang 80802024617SSatish Balay Input Parameters: 809d2567f34SHong Zhang + ts - timestepping context 810d2567f34SHong Zhang - nstages - number of stages of IRK scheme 811d2567f34SHong Zhang 812d2567f34SHong Zhang Level: intermediate 813d2567f34SHong Zhang 814d2567f34SHong Zhang .seealso: TSIRKSetNumStages(), TSIRK 815d2567f34SHong Zhang @*/ 816d2567f34SHong Zhang PetscErrorCode TSIRKGetNumStages(TS ts,PetscInt *nstages) 817d2567f34SHong Zhang { 818d2567f34SHong Zhang PetscErrorCode ierr; 819d2567f34SHong Zhang 820d2567f34SHong Zhang PetscFunctionBegin; 821d2567f34SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 822d2567f34SHong Zhang PetscValidIntPointer(nstages,2); 823d2567f34SHong Zhang ierr = PetscTryMethod(ts,"TSIRKGetNumStages_C",(TS,PetscInt*),(ts,nstages));CHKERRQ(ierr); 824d2567f34SHong Zhang PetscFunctionReturn(0); 825d2567f34SHong Zhang } 826d2567f34SHong Zhang 827d2567f34SHong Zhang static PetscErrorCode TSIRKGetType_IRK(TS ts,TSIRKType *irktype) 828d2567f34SHong Zhang { 829d2567f34SHong Zhang TS_IRK *irk = (TS_IRK*)ts->data; 830d2567f34SHong Zhang 831d2567f34SHong Zhang PetscFunctionBegin; 832d2567f34SHong Zhang *irktype = irk->method_name; 833d2567f34SHong Zhang PetscFunctionReturn(0); 834d2567f34SHong Zhang } 835d2567f34SHong Zhang 836d2567f34SHong Zhang static PetscErrorCode TSIRKSetType_IRK(TS ts,TSIRKType irktype) 837d2567f34SHong Zhang { 838d2567f34SHong Zhang TS_IRK *irk = (TS_IRK*)ts->data; 839d2567f34SHong Zhang PetscErrorCode ierr,(*irkcreate)(TS); 840d2567f34SHong Zhang 841d2567f34SHong Zhang PetscFunctionBegin; 842d2567f34SHong Zhang if (irk->method_name) { 843d2567f34SHong Zhang ierr = PetscFree(irk->method_name);CHKERRQ(ierr); 844d2567f34SHong Zhang ierr = TSIRKTableauReset(ts);CHKERRQ(ierr); 845d2567f34SHong Zhang } 846d2567f34SHong Zhang ierr = PetscFunctionListFind(TSIRKList,irktype,&irkcreate);CHKERRQ(ierr); 847*98921bdaSJacob Faibussowitsch if (!irkcreate) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSIRK type \"%s\" given",irktype); 848d2567f34SHong Zhang ierr = (*irkcreate)(ts);CHKERRQ(ierr); 849d2567f34SHong Zhang ierr = PetscStrallocpy(irktype,&irk->method_name);CHKERRQ(ierr); 850d2567f34SHong Zhang PetscFunctionReturn(0); 851d2567f34SHong Zhang } 852d2567f34SHong Zhang 853d2567f34SHong Zhang static PetscErrorCode TSIRKSetNumStages_IRK(TS ts,PetscInt nstages) 854d2567f34SHong Zhang { 855d2567f34SHong Zhang TS_IRK *irk = (TS_IRK*)ts->data; 856d2567f34SHong Zhang 857d2567f34SHong Zhang PetscFunctionBegin; 858*98921bdaSJacob Faibussowitsch if (nstages<=0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"input argument, %d, out of range",nstages); 859d2567f34SHong Zhang irk->nstages = nstages; 860d2567f34SHong Zhang PetscFunctionReturn(0); 861d2567f34SHong Zhang } 862d2567f34SHong Zhang 863d2567f34SHong Zhang static PetscErrorCode TSIRKGetNumStages_IRK(TS ts,PetscInt *nstages) 864d2567f34SHong Zhang { 865d2567f34SHong Zhang TS_IRK *irk = (TS_IRK*)ts->data; 866d2567f34SHong Zhang 867d2567f34SHong Zhang PetscFunctionBegin; 868d2567f34SHong Zhang PetscValidIntPointer(nstages,2); 869d2567f34SHong Zhang *nstages = irk->nstages; 870d2567f34SHong Zhang PetscFunctionReturn(0); 871d2567f34SHong Zhang } 872d2567f34SHong Zhang 873d2567f34SHong Zhang static PetscErrorCode TSDestroy_IRK(TS ts) 874d2567f34SHong Zhang { 875d2567f34SHong Zhang PetscErrorCode ierr; 876d2567f34SHong Zhang 877d2567f34SHong Zhang PetscFunctionBegin; 878d2567f34SHong Zhang ierr = TSReset_IRK(ts);CHKERRQ(ierr); 879d2567f34SHong Zhang if (ts->dm) { 880d2567f34SHong Zhang ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSIRK,DMRestrictHook_TSIRK,ts);CHKERRQ(ierr); 881d2567f34SHong Zhang ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSIRK,DMSubDomainRestrictHook_TSIRK,ts);CHKERRQ(ierr); 882d2567f34SHong Zhang } 883d2567f34SHong Zhang ierr = PetscFree(ts->data);CHKERRQ(ierr); 884d2567f34SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSIRKSetType_C",NULL);CHKERRQ(ierr); 885d2567f34SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSIRKGetType_C",NULL);CHKERRQ(ierr); 886d2567f34SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSIRKSetNumStages_C",NULL);CHKERRQ(ierr); 887d2567f34SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSIRKGetNumStages_C",NULL);CHKERRQ(ierr); 888d2567f34SHong Zhang PetscFunctionReturn(0); 889d2567f34SHong Zhang } 890d2567f34SHong Zhang 891d2567f34SHong Zhang /*MC 892d2567f34SHong Zhang TSIRK - ODE and DAE solver using Implicit Runge-Kutta schemes 893d2567f34SHong Zhang 894d2567f34SHong Zhang Notes: 895d2567f34SHong Zhang 896d2567f34SHong Zhang TSIRK uses the sparse Kronecker product matrix implementation of MATKAIJ to achieve good arithmetic intensity. 897d2567f34SHong Zhang 898517f343bSHong Zhang Gauss-Legrendre methods are currently supported. These are A-stable symplectic methods with an arbitrary number of stages. The order of accuracy is 2s when using s stages. The default method uses three stages and thus has an order of six. The number of stages (thus order) can be set with -ts_irk_nstages or TSIRKSetNumStages(). 899d2567f34SHong Zhang 900d2567f34SHong Zhang Level: beginner 901d2567f34SHong Zhang 902517f343bSHong Zhang .seealso: TSCreate(), TS, TSSetType(), TSIRKSetType(), TSIRKGetType(), TSIRKGAUSS, TSIRKRegister(), TSIRKSetNumStages() 903d2567f34SHong Zhang 904d2567f34SHong Zhang M*/ 905d2567f34SHong Zhang PETSC_EXTERN PetscErrorCode TSCreate_IRK(TS ts) 906d2567f34SHong Zhang { 907d2567f34SHong Zhang TS_IRK *irk; 908d2567f34SHong Zhang PetscErrorCode ierr; 909d2567f34SHong Zhang 910d2567f34SHong Zhang PetscFunctionBegin; 911d2567f34SHong Zhang ierr = TSIRKInitializePackage();CHKERRQ(ierr); 912d2567f34SHong Zhang 913d2567f34SHong Zhang ts->ops->reset = TSReset_IRK; 914d2567f34SHong Zhang ts->ops->destroy = TSDestroy_IRK; 915d2567f34SHong Zhang ts->ops->view = TSView_IRK; 916d2567f34SHong Zhang ts->ops->load = TSLoad_IRK; 917d2567f34SHong Zhang ts->ops->setup = TSSetUp_IRK; 918d2567f34SHong Zhang ts->ops->step = TSStep_IRK; 919d2567f34SHong Zhang ts->ops->interpolate = TSInterpolate_IRK; 920d2567f34SHong Zhang ts->ops->evaluatestep = TSEvaluateStep_IRK; 921d2567f34SHong Zhang ts->ops->rollback = TSRollBack_IRK; 922d2567f34SHong Zhang ts->ops->setfromoptions = TSSetFromOptions_IRK; 923d2567f34SHong Zhang ts->ops->snesfunction = SNESTSFormFunction_IRK; 924d2567f34SHong Zhang ts->ops->snesjacobian = SNESTSFormJacobian_IRK; 925d2567f34SHong Zhang 926d2567f34SHong Zhang ts->usessnes = PETSC_TRUE; 927d2567f34SHong Zhang 928d2567f34SHong Zhang ierr = PetscNewLog(ts,&irk);CHKERRQ(ierr); 929d2567f34SHong Zhang ts->data = (void*)irk; 930d2567f34SHong Zhang 931d2567f34SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSIRKSetType_C",TSIRKSetType_IRK);CHKERRQ(ierr); 932d2567f34SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSIRKGetType_C",TSIRKGetType_IRK);CHKERRQ(ierr); 933d2567f34SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSIRKSetNumStages_C",TSIRKSetNumStages_IRK);CHKERRQ(ierr); 934d2567f34SHong Zhang ierr = PetscObjectComposeFunction((PetscObject)ts,"TSIRKGetNumStages_C",TSIRKGetNumStages_IRK);CHKERRQ(ierr); 935d2567f34SHong Zhang /* 3-stage IRK_Gauss is the default */ 936d2567f34SHong Zhang ierr = PetscNew(&irk->tableau);CHKERRQ(ierr); 937d2567f34SHong Zhang irk->nstages = 3; 938d2567f34SHong Zhang ierr = TSIRKSetType(ts,TSIRKDefault);CHKERRQ(ierr); 939d2567f34SHong Zhang PetscFunctionReturn(0); 940d2567f34SHong Zhang } 941