1 /* 2 Code for timestepping with implicit Theta method 3 */ 4 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 5 #include <petscsnes.h> 6 #include <petscdm.h> 7 #include <petscmat.h> 8 9 typedef struct { 10 /* context for time stepping */ 11 PetscReal stage_time; 12 Vec X0,X,Xdot; /* Storage for stages and time derivative */ 13 Vec affine; /* Affine vector needed for residual at beginning of step in endpoint formulation */ 14 PetscReal Theta; 15 PetscReal ptime; 16 PetscReal time_step; 17 PetscInt order; 18 PetscBool endpoint; 19 PetscBool extrapolate; 20 TSStepStatus status; 21 Vec VecCostIntegral0; /* Backup for roll-backs due to events */ 22 23 /* context for sensitivity analysis */ 24 PetscInt num_tlm; /* Total number of tangent linear equations */ 25 Vec *VecsDeltaLam; /* Increment of the adjoint sensitivity w.r.t IC at stage */ 26 Vec *VecsDeltaMu; /* Increment of the adjoint sensitivity w.r.t P at stage */ 27 Vec *VecsSensiTemp; /* Vector to be multiplied with Jacobian transpose */ 28 Mat MatDeltaFwdSensip; /* Increment of the forward sensitivity at stage */ 29 Vec VecDeltaFwdSensipCol; /* Working vector for holding one column of the sensitivity matrix */ 30 Mat MatFwdSensip0; /* backup for roll-backs due to events */ 31 Vec VecIntegralSensipTemp; /* Working vector for forward integral sensitivity */ 32 Vec *VecsIntegralSensip0; /* backup for roll-backs due to events */ 33 Vec *VecsDeltaLam2; /* Increment of the 2nd-order adjoint sensitivity w.r.t IC at stage */ 34 Vec *VecsDeltaMu2; /* Increment of the 2nd-order adjoint sensitivity w.r.t P at stage */ 35 Vec *VecsSensi2Temp; /* Working vectors that holds the residual for the second-order adjoint */ 36 Vec *VecsAffine; /* Working vectors to store residuals */ 37 /* context for error estimation */ 38 Vec vec_sol_prev; 39 Vec vec_lte_work; 40 } TS_Theta; 41 42 static PetscErrorCode TSThetaGetX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 43 { 44 TS_Theta *th = (TS_Theta*)ts->data; 45 PetscErrorCode ierr; 46 47 PetscFunctionBegin; 48 if (X0) { 49 if (dm && dm != ts->dm) { 50 ierr = DMGetNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 51 } else *X0 = ts->vec_sol; 52 } 53 if (Xdot) { 54 if (dm && dm != ts->dm) { 55 ierr = DMGetNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 56 } else *Xdot = th->Xdot; 57 } 58 PetscFunctionReturn(0); 59 } 60 61 static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts,DM dm,Vec *X0,Vec *Xdot) 62 { 63 PetscErrorCode ierr; 64 65 PetscFunctionBegin; 66 if (X0) { 67 if (dm && dm != ts->dm) { 68 ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_X0",X0);CHKERRQ(ierr); 69 } 70 } 71 if (Xdot) { 72 if (dm && dm != ts->dm) { 73 ierr = DMRestoreNamedGlobalVector(dm,"TSTheta_Xdot",Xdot);CHKERRQ(ierr); 74 } 75 } 76 PetscFunctionReturn(0); 77 } 78 79 static PetscErrorCode DMCoarsenHook_TSTheta(DM fine,DM coarse,void *ctx) 80 { 81 PetscFunctionBegin; 82 PetscFunctionReturn(0); 83 } 84 85 static PetscErrorCode DMRestrictHook_TSTheta(DM fine,Mat restrct,Vec rscale,Mat inject,DM coarse,void *ctx) 86 { 87 TS ts = (TS)ctx; 88 PetscErrorCode ierr; 89 Vec X0,Xdot,X0_c,Xdot_c; 90 91 PetscFunctionBegin; 92 ierr = TSThetaGetX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr); 93 ierr = TSThetaGetX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr); 94 ierr = MatRestrict(restrct,X0,X0_c);CHKERRQ(ierr); 95 ierr = MatRestrict(restrct,Xdot,Xdot_c);CHKERRQ(ierr); 96 ierr = VecPointwiseMult(X0_c,rscale,X0_c);CHKERRQ(ierr); 97 ierr = VecPointwiseMult(Xdot_c,rscale,Xdot_c);CHKERRQ(ierr); 98 ierr = TSThetaRestoreX0AndXdot(ts,fine,&X0,&Xdot);CHKERRQ(ierr); 99 ierr = TSThetaRestoreX0AndXdot(ts,coarse,&X0_c,&Xdot_c);CHKERRQ(ierr); 100 PetscFunctionReturn(0); 101 } 102 103 static PetscErrorCode DMSubDomainHook_TSTheta(DM dm,DM subdm,void *ctx) 104 { 105 PetscFunctionBegin; 106 PetscFunctionReturn(0); 107 } 108 109 static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 110 { 111 TS ts = (TS)ctx; 112 PetscErrorCode ierr; 113 Vec X0,Xdot,X0_sub,Xdot_sub; 114 115 PetscFunctionBegin; 116 ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 117 ierr = TSThetaGetX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr); 118 119 ierr = VecScatterBegin(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 120 ierr = VecScatterEnd(gscat,X0,X0_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 121 122 ierr = VecScatterBegin(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 123 ierr = VecScatterEnd(gscat,Xdot,Xdot_sub,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 124 125 ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 126 ierr = TSThetaRestoreX0AndXdot(ts,subdm,&X0_sub,&Xdot_sub);CHKERRQ(ierr); 127 PetscFunctionReturn(0); 128 } 129 130 static PetscErrorCode TSThetaEvaluateCostIntegral(TS ts) 131 { 132 TS_Theta *th = (TS_Theta*)ts->data; 133 PetscErrorCode ierr; 134 135 PetscFunctionBegin; 136 if (th->endpoint) { 137 /* Evolve ts->vec_costintegral to compute integrals */ 138 if (th->Theta!=1.0) { 139 ierr = TSComputeCostIntegrand(ts,th->ptime,th->X0,ts->vec_costintegrand);CHKERRQ(ierr); 140 ierr = VecAXPY(ts->vec_costintegral,th->time_step*(1.0-th->Theta),ts->vec_costintegrand);CHKERRQ(ierr); 141 } 142 ierr = TSComputeCostIntegrand(ts,ts->ptime,ts->vec_sol,ts->vec_costintegrand);CHKERRQ(ierr); 143 ierr = VecAXPY(ts->vec_costintegral,th->time_step*th->Theta,ts->vec_costintegrand);CHKERRQ(ierr); 144 } else { 145 ierr = TSComputeCostIntegrand(ts,th->stage_time,th->X,ts->vec_costintegrand);CHKERRQ(ierr); 146 ierr = VecAXPY(ts->vec_costintegral,th->time_step,ts->vec_costintegrand);CHKERRQ(ierr); 147 } 148 PetscFunctionReturn(0); 149 } 150 151 static PetscErrorCode TSForwardCostIntegral_Theta(TS ts) 152 { 153 TS_Theta *th = (TS_Theta*)ts->data; 154 PetscErrorCode ierr; 155 156 PetscFunctionBegin; 157 /* backup cost integral */ 158 ierr = VecCopy(ts->vec_costintegral,th->VecCostIntegral0);CHKERRQ(ierr); 159 ierr = TSThetaEvaluateCostIntegral(ts);CHKERRQ(ierr); 160 PetscFunctionReturn(0); 161 } 162 163 static PetscErrorCode TSAdjointCostIntegral_Theta(TS ts) 164 { 165 PetscErrorCode ierr; 166 167 PetscFunctionBegin; 168 ierr = TSThetaEvaluateCostIntegral(ts);CHKERRQ(ierr); 169 PetscFunctionReturn(0); 170 } 171 172 static PetscErrorCode TSTheta_SNESSolve(TS ts,Vec b,Vec x) 173 { 174 PetscInt nits,lits; 175 PetscErrorCode ierr; 176 177 PetscFunctionBegin; 178 ierr = SNESSolve(ts->snes,b,x);CHKERRQ(ierr); 179 ierr = SNESGetIterationNumber(ts->snes,&nits);CHKERRQ(ierr); 180 ierr = SNESGetLinearSolveIterations(ts->snes,&lits);CHKERRQ(ierr); 181 ts->snes_its += nits; ts->ksp_its += lits; 182 PetscFunctionReturn(0); 183 } 184 185 static PetscErrorCode TSStep_Theta(TS ts) 186 { 187 TS_Theta *th = (TS_Theta*)ts->data; 188 PetscInt rejections = 0; 189 PetscBool stageok,accept = PETSC_TRUE; 190 PetscReal next_time_step = ts->time_step; 191 PetscErrorCode ierr; 192 193 PetscFunctionBegin; 194 if (!ts->steprollback) { 195 if (th->vec_sol_prev) { ierr = VecCopy(th->X0,th->vec_sol_prev);CHKERRQ(ierr); } 196 ierr = VecCopy(ts->vec_sol,th->X0);CHKERRQ(ierr); 197 } 198 199 th->status = TS_STEP_INCOMPLETE; 200 while (!ts->reason && th->status != TS_STEP_COMPLETE) { 201 202 PetscReal shift = 1/(th->Theta*ts->time_step); 203 th->stage_time = ts->ptime + (th->endpoint ? (PetscReal)1 : th->Theta)*ts->time_step; 204 205 ierr = VecCopy(th->X0,th->X);CHKERRQ(ierr); 206 if (th->extrapolate && !ts->steprestart) { 207 ierr = VecAXPY(th->X,1/shift,th->Xdot);CHKERRQ(ierr); 208 } 209 if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */ 210 if (!th->affine) {ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr);} 211 ierr = VecZeroEntries(th->Xdot);CHKERRQ(ierr); 212 ierr = TSComputeIFunction(ts,ts->ptime,th->X0,th->Xdot,th->affine,PETSC_FALSE);CHKERRQ(ierr); 213 ierr = VecScale(th->affine,(th->Theta-1)/th->Theta);CHKERRQ(ierr); 214 } else if (th->affine) { /* Just in case th->endpoint is changed between calls to TSStep_Theta() */ 215 ierr = VecZeroEntries(th->affine);CHKERRQ(ierr); 216 } 217 ierr = TSPreStage(ts,th->stage_time);CHKERRQ(ierr); 218 ierr = TSTheta_SNESSolve(ts,th->affine,th->X);CHKERRQ(ierr); 219 ierr = TSPostStage(ts,th->stage_time,0,&th->X);CHKERRQ(ierr); 220 ierr = TSAdaptCheckStage(ts->adapt,ts,th->stage_time,th->X,&stageok);CHKERRQ(ierr); 221 if (!stageok) goto reject_step; 222 223 th->status = TS_STEP_PENDING; 224 if (th->endpoint) { 225 ierr = VecCopy(th->X,ts->vec_sol);CHKERRQ(ierr); 226 } else { 227 ierr = VecAXPBYPCZ(th->Xdot,-shift,shift,0,th->X0,th->X);CHKERRQ(ierr); 228 ierr = VecAXPY(ts->vec_sol,ts->time_step,th->Xdot);CHKERRQ(ierr); 229 } 230 ierr = TSAdaptChoose(ts->adapt,ts,ts->time_step,NULL,&next_time_step,&accept);CHKERRQ(ierr); 231 th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE; 232 if (!accept) { 233 ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 234 ts->time_step = next_time_step; 235 goto reject_step; 236 } 237 238 if (ts->forward_solve || ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */ 239 th->ptime = ts->ptime; 240 th->time_step = ts->time_step; 241 } 242 243 ts->ptime += ts->time_step; 244 ts->time_step = next_time_step; 245 break; 246 247 reject_step: 248 ts->reject++; accept = PETSC_FALSE; 249 if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) { 250 ts->reason = TS_DIVERGED_STEP_REJECTED; 251 ierr = PetscInfo2(ts,"Step=%D, step rejections %D greater than current TS allowed, stopping solve\n",ts->steps,rejections);CHKERRQ(ierr); 252 } 253 } 254 PetscFunctionReturn(0); 255 } 256 257 static PetscErrorCode TSAdjointStep_Theta(TS ts) 258 { 259 TS_Theta *th = (TS_Theta*)ts->data; 260 Vec *VecsDeltaLam = th->VecsDeltaLam,*VecsDeltaMu = th->VecsDeltaMu,*VecsSensiTemp = th->VecsSensiTemp; 261 Vec *VecsDeltaLam2 = th->VecsDeltaLam2,*VecsDeltaMu2 = th->VecsDeltaMu2,*VecsSensi2Temp = th->VecsSensi2Temp; 262 PetscInt nadj; 263 Mat J,Jp; 264 KSP ksp; 265 PetscReal shift; 266 PetscScalar *xarr; 267 PetscErrorCode ierr; 268 269 PetscFunctionBegin; 270 th->status = TS_STEP_INCOMPLETE; 271 ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 272 ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr); 273 274 /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */ 275 th->stage_time = th->endpoint ? ts->ptime : (ts->ptime+(1.-th->Theta)*ts->time_step); /* time_step is negative*/ 276 th->ptime = ts->ptime + ts->time_step; 277 th->time_step = -ts->time_step; 278 279 /* Build RHS for first-order adjoint */ 280 if (ts->vec_costintegral) { /* Cost function has an integral term */ 281 if (th->endpoint) { 282 ierr = TSComputeDRDUFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdu);CHKERRQ(ierr); 283 } else { 284 ierr = TSComputeDRDUFunction(ts,th->stage_time,th->X,ts->vecs_drdu);CHKERRQ(ierr); 285 } 286 } 287 for (nadj=0; nadj<ts->numcost; nadj++) { 288 ierr = VecCopy(ts->vecs_sensi[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 289 ierr = VecScale(VecsSensiTemp[nadj],1./(th->Theta*th->time_step));CHKERRQ(ierr); 290 if (ts->vec_costintegral) { 291 ierr = VecAXPY(VecsSensiTemp[nadj],1.,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 292 } 293 } 294 295 /* Build LHS for first-order adjoint */ 296 shift = 1./(th->Theta*th->time_step); 297 if (th->endpoint) { 298 ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 299 } else { 300 ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 301 } 302 ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr); 303 304 /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */ 305 for (nadj=0; nadj<ts->numcost; nadj++) { 306 ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam[nadj]);CHKERRQ(ierr); 307 } 308 309 if (ts->vecs_sensi2) { /* U_{n+1} */ 310 if (!th->endpoint) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"Operation not implemented in TS_Theta"); 311 /* Get w1 at t_{n+1} from TLM matrix */ 312 ierr = MatDenseGetColumn(ts->mat_sensip,0,&xarr);CHKERRQ(ierr); 313 ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 314 /* lambda_s^T F_UU w_1 */ 315 ierr = TSComputeIHessianProductFunction1(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr); 316 ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 317 ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 318 if (ts->vecs_fup) { 319 /* lambda_s^T F_UP w_2 */ 320 ierr = TSComputeIHessianProductFunction2(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr); 321 } 322 for (nadj=0; nadj<ts->numcost; nadj++) { /* compute the residual */ 323 ierr = VecCopy(ts->vecs_sensi2[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr); 324 ierr = VecScale(VecsSensi2Temp[nadj],1./shift);CHKERRQ(ierr); 325 ierr = VecAXPY(VecsSensi2Temp[nadj],1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 326 ierr = VecAXPY(VecsSensi2Temp[nadj],1.,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 327 if (ts->vecs_fup) { 328 ierr = VecAXPY(VecsSensi2Temp[nadj],1.,ts->vecs_fup[nadj]);CHKERRQ(ierr); 329 } 330 if (ts->vec_costintegral) { 331 ierr = VecAXPY(VecsSensi2Temp[nadj],1.,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 332 } 333 } 334 /* Solve stage equation LHS X = RHS for second-order adjoint */ 335 for (nadj=0; nadj<ts->numcost; nadj++) { 336 ierr = KSPSolveTranspose(ksp,VecsSensiTemp[nadj],VecsDeltaLam2[nadj]);CHKERRQ(ierr); 337 } 338 } 339 340 /* Update sensitivities, and evaluate integrals if there is any */ 341 if(th->endpoint) { /* two-stage Theta methods */ 342 if (th->Theta!=1.) { /* general case */ 343 shift = 1./((th->Theta-1.)*th->time_step); 344 ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 345 if (ts->vec_costintegral) { /* R_U at t_n */ 346 ierr = TSComputeDRDUFunction(ts,th->ptime,th->X0,ts->vecs_drdu);CHKERRQ(ierr); 347 } 348 for (nadj=0; nadj<ts->numcost; nadj++) { 349 ierr = MatMultTranspose(J,VecsDeltaLam[nadj],ts->vecs_sensi[nadj]);CHKERRQ(ierr); 350 ierr = VecScale(ts->vecs_sensi[nadj],1./shift);CHKERRQ(ierr); 351 if (ts->vec_costintegral) { 352 ierr = VecAXPY(ts->vecs_sensi[nadj],-1./shift,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 353 } 354 } 355 if (ts->vecs_sensi2) { /* second-order */ 356 /* Get w1 at t_n from TLM matrix */ 357 ierr = MatDenseGetColumn(th->MatFwdSensip0,0,&xarr);CHKERRQ(ierr); 358 ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 359 /* lambda_s^T F_UU w_1 */ 360 ierr = TSComputeIHessianProductFunction1(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fuu);CHKERRQ(ierr); 361 ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 362 ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 363 if (ts->vecs_fup) { 364 /* lambda_s^T F_UU w_2 */ 365 ierr = TSComputeIHessianProductFunction2(ts,th->ptime,th->X0,VecsDeltaLam,ts->vec_dir,ts->vecs_fup);CHKERRQ(ierr); 366 } 367 for (nadj=0; nadj<ts->numcost; nadj++) { 368 /* M^T Lambda_s + h(1-theta) F_U^T Lambda_s + h(1-theta) R_U */ 369 ierr = MatMultTranspose(J,VecsDeltaLam2[nadj],ts->vecs_sensi2[nadj]);CHKERRQ(ierr); 370 ierr = VecScale(ts->vecs_sensi2[nadj],1./shift);CHKERRQ(ierr); 371 ierr = VecAXPY(ts->vecs_sensi2[nadj],-1./shift,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 372 ierr = VecAXPY(ts->vecs_sensi2[nadj],-1./shift,ts->vecs_fuu[nadj]);CHKERRQ(ierr); 373 if (ts->vecs_fup) { 374 ierr = VecAXPY(ts->vecs_sensi2[nadj],-1./shift,ts->vecs_fup[nadj]);CHKERRQ(ierr); 375 } 376 if (ts->vec_costintegral) { 377 ierr = VecAXPY(ts->vecs_sensi2[nadj],-1./shift,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 378 } 379 } 380 } 381 } else { /* backward Euler */ 382 shift = 0.0; 383 ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_u */ 384 for (nadj=0; nadj<ts->numcost; nadj++) { 385 ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 386 ierr = VecAXPY(ts->vecs_sensi[nadj],-th->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr); 387 if (ts->vec_costintegral) { /* wrong? */ 388 ierr = VecAXPY(ts->vecs_sensi[nadj],th->time_step,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 389 } 390 } 391 if (ts->vecs_sensi2) { 392 for (nadj=0; nadj<ts->numcost; nadj++) { 393 ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensi2Temp[nadj]);CHKERRQ(ierr); 394 ierr = VecAXPY(ts->vecs_sensi2[nadj],-th->time_step,VecsSensi2Temp[nadj]);CHKERRQ(ierr); 395 } 396 } 397 } 398 399 if (ts->vecs_sensip) { /* sensitivities wrt parameters */ 400 /* U_{n+1} */ 401 ierr = TSComputeRHSJacobianP(ts,th->stage_time,ts->vec_sol,ts->Jacp);CHKERRQ(ierr); 402 if (ts->vec_costintegral) { 403 ierr = TSComputeDRDPFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdp);CHKERRQ(ierr); 404 } 405 for (nadj=0; nadj<ts->numcost; nadj++) { 406 ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 407 ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step*th->Theta,VecsDeltaMu[nadj]);CHKERRQ(ierr); 408 } 409 if (ts->vecs_sensip2) { /* second-order */ 410 /* lambda_s^T F_PU w_1 */ 411 ierr = TSComputeIHessianProductFunction3(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr); 412 /* lambda_s^T F_PP w_2 */ 413 ierr = TSComputeIHessianProductFunction4(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr); 414 for (nadj=0; nadj<ts->numcost; nadj++) { 415 ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr); 416 ierr = VecAXPY(ts->vecs_sensip2[nadj],th->time_step*th->Theta,VecsDeltaMu2[nadj]);CHKERRQ(ierr); 417 if (ts->vecs_fpu) { 418 ierr = VecAXPY(ts->vecs_sensi2[nadj],th->time_step*th->Theta,ts->vecs_fpu[nadj]);CHKERRQ(ierr); 419 } 420 if (ts->vecs_fpp) { 421 ierr = VecAXPY(ts->vecs_sensi2[nadj],th->time_step*th->Theta,ts->vecs_fpp[nadj]);CHKERRQ(ierr); 422 } 423 if (ts->vec_costintegral) { 424 ierr = VecAXPY(ts->vecs_sensip2[nadj],th->time_step*th->Theta,ts->vecs_drdp[nadj]);CHKERRQ(ierr); 425 } 426 } 427 } 428 429 /* U_s */ 430 if (th->Theta!=1.) { 431 ierr = TSComputeRHSJacobianP(ts,th->ptime,th->X0,ts->Jacp);CHKERRQ(ierr); 432 if (ts->vec_costintegral) { 433 ierr = TSComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr); 434 } 435 for (nadj=0; nadj<ts->numcost; nadj++) { 436 ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 437 ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step*(1.-th->Theta),VecsDeltaMu[nadj]);CHKERRQ(ierr); 438 if (ts->vecs_sensip2) { /* second-order */ 439 /* lambda_s^T F_PU w_1 */ 440 ierr = TSComputeIHessianProductFunction3(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_sensip_col,ts->vecs_fpu);CHKERRQ(ierr); 441 /* lambda_s^T F_PP w_2 */ 442 ierr = TSComputeIHessianProductFunction4(ts,th->stage_time,ts->vec_sol,VecsDeltaLam,ts->vec_dir,ts->vecs_fpp);CHKERRQ(ierr); 443 for (nadj=0; nadj<ts->numcost; nadj++) { 444 ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam2[nadj],VecsDeltaMu2[nadj]);CHKERRQ(ierr); 445 ierr = VecAXPY(ts->vecs_sensip2[nadj],th->time_step*(1.-th->Theta),VecsDeltaMu2[nadj]);CHKERRQ(ierr); 446 if (ts->vecs_fpu) { 447 ierr = VecAXPY(ts->vecs_sensi2[nadj],th->time_step*(1.-th->Theta),ts->vecs_fpu[nadj]);CHKERRQ(ierr); 448 } 449 if (ts->vecs_fpp) { 450 ierr = VecAXPY(ts->vecs_sensi2[nadj],th->time_step*(1.-th->Theta),ts->vecs_fpp[nadj]);CHKERRQ(ierr); 451 } 452 if (ts->vec_costintegral) { 453 ierr = VecAXPY(ts->vecs_sensip2[nadj],th->time_step*(1.-th->Theta),ts->vecs_drdp[nadj]);CHKERRQ(ierr); 454 } 455 } 456 } 457 } 458 } 459 } 460 } else { /* one-stage case */ 461 shift = 0.0; 462 ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); /* get -f_y */ 463 if (ts->vec_costintegral) { 464 ierr = TSComputeDRDUFunction(ts,th->stage_time,th->X,ts->vecs_drdu);CHKERRQ(ierr); 465 } 466 for (nadj=0; nadj<ts->numcost; nadj++) { 467 ierr = MatMultTranspose(J,VecsDeltaLam[nadj],VecsSensiTemp[nadj]);CHKERRQ(ierr); 468 ierr = VecAXPY(ts->vecs_sensi[nadj],-th->time_step,VecsSensiTemp[nadj]);CHKERRQ(ierr); 469 if (ts->vec_costintegral) { 470 ierr = VecAXPY(ts->vecs_sensi[nadj],th->time_step,ts->vecs_drdu[nadj]);CHKERRQ(ierr); 471 } 472 } 473 if (ts->vecs_sensip) { 474 ierr = TSComputeRHSJacobianP(ts,th->stage_time,th->X,ts->Jacp);CHKERRQ(ierr); 475 for (nadj=0; nadj<ts->numcost; nadj++) { 476 ierr = MatMultTranspose(ts->Jacp,VecsDeltaLam[nadj],VecsDeltaMu[nadj]);CHKERRQ(ierr); 477 ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,VecsDeltaMu[nadj]);CHKERRQ(ierr); 478 } 479 if (ts->vec_costintegral) { 480 ierr = TSComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr); 481 for (nadj=0; nadj<ts->numcost; nadj++) { 482 ierr = VecAXPY(ts->vecs_sensip[nadj],th->time_step,ts->vecs_drdp[nadj]);CHKERRQ(ierr); 483 } 484 } 485 } 486 } 487 488 th->status = TS_STEP_COMPLETE; 489 PetscFunctionReturn(0); 490 } 491 492 static PetscErrorCode TSInterpolate_Theta(TS ts,PetscReal t,Vec X) 493 { 494 TS_Theta *th = (TS_Theta*)ts->data; 495 PetscReal dt = t - ts->ptime; 496 PetscErrorCode ierr; 497 498 PetscFunctionBegin; 499 ierr = VecCopy(ts->vec_sol,th->X);CHKERRQ(ierr); 500 if (th->endpoint) dt *= th->Theta; 501 ierr = VecWAXPY(X,dt,th->Xdot,th->X);CHKERRQ(ierr); 502 PetscFunctionReturn(0); 503 } 504 505 static PetscErrorCode TSEvaluateWLTE_Theta(TS ts,NormType wnormtype,PetscInt *order,PetscReal *wlte) 506 { 507 TS_Theta *th = (TS_Theta*)ts->data; 508 Vec X = ts->vec_sol; /* X = solution */ 509 Vec Y = th->vec_lte_work; /* Y = X + LTE */ 510 PetscReal wltea,wlter; 511 PetscErrorCode ierr; 512 513 PetscFunctionBegin; 514 if (!th->vec_sol_prev) {*wlte = -1; PetscFunctionReturn(0);} 515 /* Cannot compute LTE in first step or in restart after event */ 516 if (ts->steprestart) {*wlte = -1; PetscFunctionReturn(0);} 517 /* Compute LTE using backward differences with non-constant time step */ 518 { 519 PetscReal h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev; 520 PetscReal a = 1 + h_prev/h; 521 PetscScalar scal[3]; Vec vecs[3]; 522 scal[0] = +1/a; scal[1] = -1/(a-1); scal[2] = +1/(a*(a-1)); 523 vecs[0] = X; vecs[1] = th->X0; vecs[2] = th->vec_sol_prev; 524 ierr = VecCopy(X,Y);CHKERRQ(ierr); 525 ierr = VecMAXPY(Y,3,scal,vecs);CHKERRQ(ierr); 526 ierr = TSErrorWeightedNorm(ts,X,Y,wnormtype,wlte,&wltea,&wlter);CHKERRQ(ierr); 527 } 528 if (order) *order = 2; 529 PetscFunctionReturn(0); 530 } 531 532 static PetscErrorCode TSRollBack_Theta(TS ts) 533 { 534 TS_Theta *th = (TS_Theta*)ts->data; 535 PetscInt ncost; 536 PetscErrorCode ierr; 537 538 PetscFunctionBegin; 539 ierr = VecCopy(th->X0,ts->vec_sol);CHKERRQ(ierr); 540 if (ts->vec_costintegral && ts->costintegralfwd) { 541 ierr = VecCopy(th->VecCostIntegral0,ts->vec_costintegral);CHKERRQ(ierr); 542 } 543 th->status = TS_STEP_INCOMPLETE; 544 if (ts->mat_sensip) { 545 ierr = MatCopy(th->MatFwdSensip0,ts->mat_sensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 546 } 547 if (ts->vecs_integral_sensip) { 548 for (ncost=0;ncost<ts->numcost;ncost++) { 549 ierr = VecCopy(th->VecsIntegralSensip0[ncost],ts->vecs_integral_sensip[ncost]);CHKERRQ(ierr); 550 } 551 } 552 PetscFunctionReturn(0); 553 } 554 555 static PetscErrorCode TSForwardStep_Theta(TS ts) 556 { 557 TS_Theta *th = (TS_Theta*)ts->data; 558 Mat MatDeltaFwdSensip = th->MatDeltaFwdSensip; 559 Vec VecDeltaFwdSensipCol = th->VecDeltaFwdSensipCol; 560 PetscInt ncost,ntlm; 561 KSP ksp; 562 Mat J,Jp; 563 PetscReal shift; 564 PetscScalar *barr,*xarr; 565 PetscErrorCode ierr; 566 567 PetscFunctionBegin; 568 ierr = MatCopy(ts->mat_sensip,th->MatFwdSensip0,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 569 570 for (ncost=0; ncost<ts->numcost; ncost++) { 571 if (ts->vecs_integral_sensip) { 572 ierr = VecCopy(ts->vecs_integral_sensip[ncost],th->VecsIntegralSensip0[ncost]);CHKERRQ(ierr); 573 } 574 } 575 576 ierr = SNESGetKSP(ts->snes,&ksp);CHKERRQ(ierr); 577 ierr = TSGetIJacobian(ts,&J,&Jp,NULL,NULL);CHKERRQ(ierr); 578 579 /* Build RHS */ 580 if (th->endpoint) { /* 2-stage method*/ 581 shift = 1./((th->Theta-1.)*th->time_step); 582 ierr = TSComputeIJacobian(ts,th->ptime,th->X0,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 583 ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr); 584 ierr = MatScale(MatDeltaFwdSensip,(th->Theta-1.)/th->Theta);CHKERRQ(ierr); 585 586 /* Add the f_p forcing terms */ 587 if (ts->Jacp) { 588 ierr = TSComputeRHSJacobianP(ts,th->ptime,th->X0,ts->Jacp);CHKERRQ(ierr); 589 ierr = MatAXPY(MatDeltaFwdSensip,(1.-th->Theta)/th->Theta,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 590 ierr = TSComputeRHSJacobianP(ts,th->stage_time,ts->vec_sol,ts->Jacp);CHKERRQ(ierr); 591 ierr = MatAXPY(MatDeltaFwdSensip,1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 592 } 593 } else { /* 1-stage method */ 594 shift = 0.0; 595 ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 596 ierr = MatMatMult(J,ts->mat_sensip,MAT_REUSE_MATRIX,PETSC_DEFAULT,&MatDeltaFwdSensip);CHKERRQ(ierr); 597 ierr = MatScale(MatDeltaFwdSensip,-1.);CHKERRQ(ierr); 598 599 /* Add the f_p forcing terms */ 600 if (ts->Jacp) { 601 ierr = TSComputeRHSJacobianP(ts,th->stage_time,th->X,ts->Jacp);CHKERRQ(ierr); 602 ierr = MatAXPY(MatDeltaFwdSensip,1.,ts->Jacp,SUBSET_NONZERO_PATTERN);CHKERRQ(ierr); 603 } 604 } 605 606 /* Build LHS */ 607 shift = 1/(th->Theta*th->time_step); 608 if (th->endpoint) { 609 ierr = TSComputeIJacobian(ts,th->stage_time,ts->vec_sol,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 610 } else { 611 ierr = TSComputeIJacobian(ts,th->stage_time,th->X,th->Xdot,shift,J,Jp,PETSC_FALSE);CHKERRQ(ierr); 612 } 613 ierr = KSPSetOperators(ksp,J,Jp);CHKERRQ(ierr); 614 615 /* 616 Evaluate the first stage of integral gradients with the 2-stage method: 617 drdu|t_n*S(t_n) + drdp|t_n 618 This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1}) 619 */ 620 if (th->endpoint) { /* 2-stage method only */ 621 if (ts->vecs_integral_sensip) { 622 ierr = TSComputeDRDUFunction(ts,th->ptime,th->X0,ts->vecs_drdu);CHKERRQ(ierr); 623 if (ts->vecs_drdp) { 624 ierr = TSComputeDRDPFunction(ts,th->ptime,th->X0,ts->vecs_drdp);CHKERRQ(ierr); 625 } 626 for (ncost=0; ncost<ts->numcost; ncost++) { 627 ierr = MatMultTranspose(ts->mat_sensip,ts->vecs_drdu[ncost],th->VecIntegralSensipTemp);CHKERRQ(ierr); 628 if (ts->vecs_drdp) { 629 ierr = VecAXPY(th->VecIntegralSensipTemp,1,ts->vecs_drdp[ncost]);CHKERRQ(ierr); 630 } 631 ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step*(1.-th->Theta),th->VecIntegralSensipTemp);CHKERRQ(ierr); 632 } 633 } 634 } 635 636 /* Solve the tangent linear equation for forward sensitivities to parameters */ 637 for (ntlm=0; ntlm<th->num_tlm; ntlm++) { 638 KSPConvergedReason kspreason; 639 ierr = MatDenseGetColumn(MatDeltaFwdSensip,ntlm,&barr);CHKERRQ(ierr); 640 ierr = VecPlaceArray(VecDeltaFwdSensipCol,barr);CHKERRQ(ierr); 641 if (th->endpoint) { 642 ierr = MatDenseGetColumn(ts->mat_sensip,ntlm,&xarr);CHKERRQ(ierr); 643 ierr = VecPlaceArray(ts->vec_sensip_col,xarr);CHKERRQ(ierr); 644 ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,ts->vec_sensip_col);CHKERRQ(ierr); 645 ierr = VecResetArray(ts->vec_sensip_col);CHKERRQ(ierr); 646 ierr = MatDenseRestoreColumn(ts->mat_sensip,&xarr);CHKERRQ(ierr); 647 } else { 648 ierr = KSPSolve(ksp,VecDeltaFwdSensipCol,VecDeltaFwdSensipCol);CHKERRQ(ierr); 649 } 650 ierr = KSPGetConvergedReason(ksp,&kspreason);CHKERRQ(ierr); 651 if (kspreason < 0) { 652 ts->reason = TSFORWARD_DIVERGED_LINEAR_SOLVE; 653 ierr = PetscInfo2(ts,"Step=%D, %Dth tangent linear solve, linear solve fails, stopping tangent linear solve\n",ts->steps,ntlm);CHKERRQ(ierr); 654 } 655 ierr = VecResetArray(VecDeltaFwdSensipCol);CHKERRQ(ierr); 656 ierr = MatDenseRestoreColumn(MatDeltaFwdSensip,&barr);CHKERRQ(ierr); 657 } 658 659 660 /* 661 Evaluate the second stage of integral gradients with the 2-stage method: 662 drdu|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1} 663 */ 664 if (ts->vecs_integral_sensip) { 665 if (!th->endpoint) { 666 ierr = MatAXPY(ts->mat_sensip,1,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 667 ierr = TSComputeDRDUFunction(ts,th->stage_time,th->X,ts->vecs_drdu);CHKERRQ(ierr); 668 if (ts->vecs_drdp) { 669 ierr = TSComputeDRDPFunction(ts,th->stage_time,th->X,ts->vecs_drdp);CHKERRQ(ierr); 670 } 671 for (ncost=0; ncost<ts->numcost; ncost++) { 672 ierr = MatMultTranspose(ts->mat_sensip,ts->vecs_drdu[ncost],th->VecIntegralSensipTemp);CHKERRQ(ierr); 673 if (ts->vecs_drdp) { 674 ierr = VecAXPY(th->VecIntegralSensipTemp,1,ts->vecs_drdp[ncost]);CHKERRQ(ierr); 675 } 676 ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step,th->VecIntegralSensipTemp);CHKERRQ(ierr); 677 } 678 ierr = MatAXPY(ts->mat_sensip,(1.-th->Theta)/th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 679 } else { 680 ierr = TSComputeDRDUFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdu);CHKERRQ(ierr); 681 if (ts->vecs_drdp) { 682 ierr = TSComputeDRDPFunction(ts,th->stage_time,ts->vec_sol,ts->vecs_drdp);CHKERRQ(ierr); 683 } 684 for (ncost=0; ncost<ts->numcost; ncost++) { 685 ierr = MatMultTranspose(ts->mat_sensip,ts->vecs_drdu[ncost],th->VecIntegralSensipTemp);CHKERRQ(ierr); 686 if (ts->vecs_drdp) { 687 ierr = VecAXPY(th->VecIntegralSensipTemp,1,ts->vecs_drdp[ncost]);CHKERRQ(ierr); 688 } 689 ierr = VecAXPY(ts->vecs_integral_sensip[ncost],th->time_step*th->Theta,th->VecIntegralSensipTemp);CHKERRQ(ierr); 690 } 691 } 692 } else { 693 if (!th->endpoint) { 694 ierr = MatAXPY(ts->mat_sensip,1./th->Theta,MatDeltaFwdSensip,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 695 } 696 } 697 PetscFunctionReturn(0); 698 } 699 700 static PetscErrorCode TSForwardGetStages_Theta(TS ts,PetscInt *ns,Mat **stagesensip) 701 { 702 TS_Theta *th = (TS_Theta*)ts->data; 703 704 PetscFunctionBegin; 705 if (ns) *ns = 1; 706 if (stagesensip) *stagesensip = th->endpoint ? &(th->MatFwdSensip0) : &(th->MatDeltaFwdSensip); 707 PetscFunctionReturn(0); 708 } 709 710 /*------------------------------------------------------------*/ 711 static PetscErrorCode TSReset_Theta(TS ts) 712 { 713 TS_Theta *th = (TS_Theta*)ts->data; 714 PetscErrorCode ierr; 715 716 PetscFunctionBegin; 717 ierr = VecDestroy(&th->X);CHKERRQ(ierr); 718 ierr = VecDestroy(&th->Xdot);CHKERRQ(ierr); 719 ierr = VecDestroy(&th->X0);CHKERRQ(ierr); 720 ierr = VecDestroy(&th->affine);CHKERRQ(ierr); 721 722 ierr = VecDestroy(&th->vec_sol_prev);CHKERRQ(ierr); 723 ierr = VecDestroy(&th->vec_lte_work);CHKERRQ(ierr); 724 725 ierr = VecDestroy(&th->VecCostIntegral0);CHKERRQ(ierr); 726 if (ts->forward_solve) { 727 if (ts->vecs_integral_sensip) { 728 ierr = VecDestroy(&th->VecIntegralSensipTemp);CHKERRQ(ierr); 729 ierr = VecDestroyVecs(ts->numcost,&th->VecsIntegralSensip0);CHKERRQ(ierr); 730 } 731 ierr = VecDestroy(&th->VecDeltaFwdSensipCol);CHKERRQ(ierr); 732 ierr = MatDestroy(&th->MatDeltaFwdSensip);CHKERRQ(ierr); 733 ierr = MatDestroy(&th->MatFwdSensip0);CHKERRQ(ierr); 734 } 735 ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 736 ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 737 ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr); 738 ierr = VecDestroyVecs(ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr); 739 ierr = VecDestroyVecs(ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 740 ierr = VecDestroyVecs(ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr); 741 742 PetscFunctionReturn(0); 743 } 744 745 static PetscErrorCode TSDestroy_Theta(TS ts) 746 { 747 PetscErrorCode ierr; 748 749 PetscFunctionBegin; 750 ierr = TSReset_Theta(ts);CHKERRQ(ierr); 751 if (ts->dm) { 752 ierr = DMCoarsenHookRemove(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 753 ierr = DMSubDomainHookRemove(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 754 } 755 ierr = PetscFree(ts->data);CHKERRQ(ierr); 756 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",NULL);CHKERRQ(ierr); 757 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",NULL);CHKERRQ(ierr); 758 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",NULL);CHKERRQ(ierr); 759 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",NULL);CHKERRQ(ierr); 760 PetscFunctionReturn(0); 761 } 762 763 /* 764 This defines the nonlinear equation that is to be solved with SNES 765 G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0 766 */ 767 static PetscErrorCode SNESTSFormFunction_Theta(SNES snes,Vec x,Vec y,TS ts) 768 { 769 TS_Theta *th = (TS_Theta*)ts->data; 770 PetscErrorCode ierr; 771 Vec X0,Xdot; 772 DM dm,dmsave; 773 PetscReal shift = 1/(th->Theta*ts->time_step); 774 775 PetscFunctionBegin; 776 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 777 /* When using the endpoint variant, this is actually 1/Theta * Xdot */ 778 ierr = TSThetaGetX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 779 ierr = VecAXPBYPCZ(Xdot,-shift,shift,0,X0,x);CHKERRQ(ierr); 780 781 /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */ 782 dmsave = ts->dm; 783 ts->dm = dm; 784 ierr = TSComputeIFunction(ts,th->stage_time,x,Xdot,y,PETSC_FALSE);CHKERRQ(ierr); 785 ts->dm = dmsave; 786 ierr = TSThetaRestoreX0AndXdot(ts,dm,&X0,&Xdot);CHKERRQ(ierr); 787 PetscFunctionReturn(0); 788 } 789 790 static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes,Vec x,Mat A,Mat B,TS ts) 791 { 792 TS_Theta *th = (TS_Theta*)ts->data; 793 PetscErrorCode ierr; 794 Vec Xdot; 795 DM dm,dmsave; 796 PetscReal shift = 1/(th->Theta*ts->time_step); 797 798 PetscFunctionBegin; 799 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 800 /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */ 801 ierr = TSThetaGetX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 802 803 dmsave = ts->dm; 804 ts->dm = dm; 805 ierr = TSComputeIJacobian(ts,th->stage_time,x,Xdot,shift,A,B,PETSC_FALSE);CHKERRQ(ierr); 806 ts->dm = dmsave; 807 ierr = TSThetaRestoreX0AndXdot(ts,dm,NULL,&Xdot);CHKERRQ(ierr); 808 PetscFunctionReturn(0); 809 } 810 811 static PetscErrorCode TSForwardSetUp_Theta(TS ts) 812 { 813 TS_Theta *th = (TS_Theta*)ts->data; 814 PetscErrorCode ierr; 815 816 PetscFunctionBegin; 817 /* combine sensitivities to parameters and sensitivities to initial values into one array */ 818 th->num_tlm = ts->num_parameters; 819 ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatDeltaFwdSensip);CHKERRQ(ierr); 820 if (ts->vecs_integral_sensip) { 821 ierr = VecDuplicate(ts->vecs_integral_sensip[0],&th->VecIntegralSensipTemp);CHKERRQ(ierr); 822 } 823 /* backup sensitivity results for roll-backs */ 824 ierr = MatDuplicate(ts->mat_sensip,MAT_DO_NOT_COPY_VALUES,&th->MatFwdSensip0);CHKERRQ(ierr); 825 826 if (ts->vecs_integral_sensip) { 827 ierr = VecDuplicateVecs(ts->vecs_integral_sensip[0],ts->numcost,&th->VecsIntegralSensip0);CHKERRQ(ierr); 828 } 829 ierr = VecDuplicate(ts->vec_sol,&th->VecDeltaFwdSensipCol);CHKERRQ(ierr); 830 ierr = VecDuplicate(ts->vec_sol,&ts->vec_sensip_col);CHKERRQ(ierr); 831 PetscFunctionReturn(0); 832 } 833 834 static PetscErrorCode TSSetUp_Theta(TS ts) 835 { 836 TS_Theta *th = (TS_Theta*)ts->data; 837 PetscBool match; 838 PetscErrorCode ierr; 839 840 PetscFunctionBegin; 841 if (!th->VecCostIntegral0 && ts->vec_costintegral && ts->costintegralfwd) { /* back up cost integral */ 842 ierr = VecDuplicate(ts->vec_costintegral,&th->VecCostIntegral0);CHKERRQ(ierr); 843 } 844 if (!th->X) { 845 ierr = VecDuplicate(ts->vec_sol,&th->X);CHKERRQ(ierr); 846 } 847 if (!th->Xdot) { 848 ierr = VecDuplicate(ts->vec_sol,&th->Xdot);CHKERRQ(ierr); 849 } 850 if (!th->X0) { 851 ierr = VecDuplicate(ts->vec_sol,&th->X0);CHKERRQ(ierr); 852 } 853 if (th->endpoint) { 854 ierr = VecDuplicate(ts->vec_sol,&th->affine);CHKERRQ(ierr); 855 } 856 857 th->order = (th->Theta == 0.5) ? 2 : 1; 858 859 ierr = TSGetDM(ts,&ts->dm);CHKERRQ(ierr); 860 ierr = DMCoarsenHookAdd(ts->dm,DMCoarsenHook_TSTheta,DMRestrictHook_TSTheta,ts);CHKERRQ(ierr); 861 ierr = DMSubDomainHookAdd(ts->dm,DMSubDomainHook_TSTheta,DMSubDomainRestrictHook_TSTheta,ts);CHKERRQ(ierr); 862 863 ierr = TSGetAdapt(ts,&ts->adapt);CHKERRQ(ierr); 864 ierr = TSAdaptCandidatesClear(ts->adapt);CHKERRQ(ierr); 865 ierr = PetscObjectTypeCompare((PetscObject)ts->adapt,TSADAPTNONE,&match);CHKERRQ(ierr); 866 if (!match) { 867 ierr = VecDuplicate(ts->vec_sol,&th->vec_sol_prev);CHKERRQ(ierr); 868 ierr = VecDuplicate(ts->vec_sol,&th->vec_lte_work);CHKERRQ(ierr); 869 } 870 ierr = TSGetSNES(ts,&ts->snes);CHKERRQ(ierr); 871 PetscFunctionReturn(0); 872 } 873 874 /*------------------------------------------------------------*/ 875 876 static PetscErrorCode TSAdjointSetUp_Theta(TS ts) 877 { 878 TS_Theta *th = (TS_Theta*)ts->data; 879 PetscErrorCode ierr; 880 881 PetscFunctionBegin; 882 ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam);CHKERRQ(ierr); 883 ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsSensiTemp);CHKERRQ(ierr); 884 if (ts->vecs_sensip) { 885 ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&th->VecsDeltaMu);CHKERRQ(ierr); 886 } 887 if (ts->vecs_sensi2) { 888 ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&th->VecsDeltaLam2);CHKERRQ(ierr); 889 ierr = VecDuplicateVecs(ts->vecs_sensi2[0],ts->numcost,&th->VecsSensi2Temp);CHKERRQ(ierr); 890 } 891 if (ts->vecs_sensip2) { 892 ierr = VecDuplicateVecs(ts->vecs_sensi2[0],ts->numcost,&th->VecsDeltaMu2);CHKERRQ(ierr); 893 } 894 PetscFunctionReturn(0); 895 } 896 897 static PetscErrorCode TSSetFromOptions_Theta(PetscOptionItems *PetscOptionsObject,TS ts) 898 { 899 TS_Theta *th = (TS_Theta*)ts->data; 900 PetscErrorCode ierr; 901 902 PetscFunctionBegin; 903 ierr = PetscOptionsHead(PetscOptionsObject,"Theta ODE solver options");CHKERRQ(ierr); 904 { 905 ierr = PetscOptionsReal("-ts_theta_theta","Location of stage (0<Theta<=1)","TSThetaSetTheta",th->Theta,&th->Theta,NULL);CHKERRQ(ierr); 906 ierr = PetscOptionsBool("-ts_theta_endpoint","Use the endpoint instead of midpoint form of the Theta method","TSThetaSetEndpoint",th->endpoint,&th->endpoint,NULL);CHKERRQ(ierr); 907 ierr = PetscOptionsBool("-ts_theta_initial_guess_extrapolate","Extrapolate stage initial guess from previous solution (sometimes unstable)","TSThetaSetExtrapolate",th->extrapolate,&th->extrapolate,NULL);CHKERRQ(ierr); 908 } 909 ierr = PetscOptionsTail();CHKERRQ(ierr); 910 PetscFunctionReturn(0); 911 } 912 913 static PetscErrorCode TSView_Theta(TS ts,PetscViewer viewer) 914 { 915 TS_Theta *th = (TS_Theta*)ts->data; 916 PetscBool iascii; 917 PetscErrorCode ierr; 918 919 PetscFunctionBegin; 920 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 921 if (iascii) { 922 ierr = PetscViewerASCIIPrintf(viewer," Theta=%g\n",(double)th->Theta);CHKERRQ(ierr); 923 ierr = PetscViewerASCIIPrintf(viewer," Extrapolation=%s\n",th->extrapolate ? "yes" : "no");CHKERRQ(ierr); 924 } 925 PetscFunctionReturn(0); 926 } 927 928 static PetscErrorCode TSThetaGetTheta_Theta(TS ts,PetscReal *theta) 929 { 930 TS_Theta *th = (TS_Theta*)ts->data; 931 932 PetscFunctionBegin; 933 *theta = th->Theta; 934 PetscFunctionReturn(0); 935 } 936 937 static PetscErrorCode TSThetaSetTheta_Theta(TS ts,PetscReal theta) 938 { 939 TS_Theta *th = (TS_Theta*)ts->data; 940 941 PetscFunctionBegin; 942 if (theta <= 0 || 1 < theta) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Theta %g not in range (0,1]",(double)theta); 943 th->Theta = theta; 944 th->order = (th->Theta == 0.5) ? 2 : 1; 945 PetscFunctionReturn(0); 946 } 947 948 static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts,PetscBool *endpoint) 949 { 950 TS_Theta *th = (TS_Theta*)ts->data; 951 952 PetscFunctionBegin; 953 *endpoint = th->endpoint; 954 PetscFunctionReturn(0); 955 } 956 957 static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts,PetscBool flg) 958 { 959 TS_Theta *th = (TS_Theta*)ts->data; 960 961 PetscFunctionBegin; 962 th->endpoint = flg; 963 PetscFunctionReturn(0); 964 } 965 966 #if defined(PETSC_HAVE_COMPLEX) 967 static PetscErrorCode TSComputeLinearStability_Theta(TS ts,PetscReal xr,PetscReal xi,PetscReal *yr,PetscReal *yi) 968 { 969 PetscComplex z = xr + xi*PETSC_i,f; 970 TS_Theta *th = (TS_Theta*)ts->data; 971 const PetscReal one = 1.0; 972 973 PetscFunctionBegin; 974 f = (one + (one - th->Theta)*z)/(one - th->Theta*z); 975 *yr = PetscRealPartComplex(f); 976 *yi = PetscImaginaryPartComplex(f); 977 PetscFunctionReturn(0); 978 } 979 #endif 980 981 static PetscErrorCode TSGetStages_Theta(TS ts,PetscInt *ns,Vec **Y) 982 { 983 TS_Theta *th = (TS_Theta*)ts->data; 984 985 PetscFunctionBegin; 986 if (ns) *ns = 1; 987 if (Y) *Y = th->endpoint ? &(th->X0) : &(th->X); 988 PetscFunctionReturn(0); 989 } 990 991 /* ------------------------------------------------------------ */ 992 /*MC 993 TSTHETA - DAE solver using the implicit Theta method 994 995 Level: beginner 996 997 Options Database: 998 + -ts_theta_theta <Theta> - Location of stage (0<Theta<=1) 999 . -ts_theta_endpoint <flag> - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method 1000 - -ts_theta_initial_guess_extrapolate <flg> - Extrapolate stage initial guess from previous solution (sometimes unstable) 1001 1002 Notes: 1003 $ -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER) 1004 $ -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule 1005 $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN) 1006 1007 This method can be applied to DAE. 1008 1009 This method is cast as a 1-stage implicit Runge-Kutta method. 1010 1011 .vb 1012 Theta | Theta 1013 ------------- 1014 | 1 1015 .ve 1016 1017 For the default Theta=0.5, this is also known as the implicit midpoint rule. 1018 1019 When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit: 1020 1021 .vb 1022 0 | 0 0 1023 1 | 1-Theta Theta 1024 ------------------- 1025 | 1-Theta Theta 1026 .ve 1027 1028 For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN). 1029 1030 To apply a diagonally implicit RK method to DAE, the stage formula 1031 1032 $ Y_i = X + h sum_j a_ij Y'_j 1033 1034 is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i) 1035 1036 .seealso: TSCreate(), TS, TSSetType(), TSCN, TSBEULER, TSThetaSetTheta(), TSThetaSetEndpoint() 1037 1038 M*/ 1039 PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts) 1040 { 1041 TS_Theta *th; 1042 PetscErrorCode ierr; 1043 1044 PetscFunctionBegin; 1045 ts->ops->reset = TSReset_Theta; 1046 ts->ops->destroy = TSDestroy_Theta; 1047 ts->ops->view = TSView_Theta; 1048 ts->ops->setup = TSSetUp_Theta; 1049 ts->ops->adjointsetup = TSAdjointSetUp_Theta; 1050 ts->ops->step = TSStep_Theta; 1051 ts->ops->interpolate = TSInterpolate_Theta; 1052 ts->ops->evaluatewlte = TSEvaluateWLTE_Theta; 1053 ts->ops->rollback = TSRollBack_Theta; 1054 ts->ops->setfromoptions = TSSetFromOptions_Theta; 1055 ts->ops->snesfunction = SNESTSFormFunction_Theta; 1056 ts->ops->snesjacobian = SNESTSFormJacobian_Theta; 1057 #if defined(PETSC_HAVE_COMPLEX) 1058 ts->ops->linearstability = TSComputeLinearStability_Theta; 1059 #endif 1060 ts->ops->getstages = TSGetStages_Theta; 1061 ts->ops->adjointstep = TSAdjointStep_Theta; 1062 ts->ops->adjointintegral = TSAdjointCostIntegral_Theta; 1063 ts->ops->forwardintegral = TSForwardCostIntegral_Theta; 1064 ts->default_adapt_type = TSADAPTNONE; 1065 1066 ts->ops->forwardsetup = TSForwardSetUp_Theta; 1067 ts->ops->forwardstep = TSForwardStep_Theta; 1068 ts->ops->forwardgetstages = TSForwardGetStages_Theta; 1069 1070 ts->usessnes = PETSC_TRUE; 1071 1072 ierr = PetscNewLog(ts,&th);CHKERRQ(ierr); 1073 ts->data = (void*)th; 1074 1075 th->VecsDeltaLam = NULL; 1076 th->VecsDeltaMu = NULL; 1077 th->VecsSensiTemp = NULL; 1078 th->VecsSensi2Temp = NULL; 1079 1080 th->extrapolate = PETSC_FALSE; 1081 th->Theta = 0.5; 1082 th->order = 2; 1083 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetTheta_C",TSThetaGetTheta_Theta);CHKERRQ(ierr); 1084 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetTheta_C",TSThetaSetTheta_Theta);CHKERRQ(ierr); 1085 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaGetEndpoint_C",TSThetaGetEndpoint_Theta);CHKERRQ(ierr); 1086 ierr = PetscObjectComposeFunction((PetscObject)ts,"TSThetaSetEndpoint_C",TSThetaSetEndpoint_Theta);CHKERRQ(ierr); 1087 PetscFunctionReturn(0); 1088 } 1089 1090 /*@ 1091 TSThetaGetTheta - Get the abscissa of the stage in (0,1]. 1092 1093 Not Collective 1094 1095 Input Parameter: 1096 . ts - timestepping context 1097 1098 Output Parameter: 1099 . theta - stage abscissa 1100 1101 Note: 1102 Use of this function is normally only required to hack TSTHETA to use a modified integration scheme. 1103 1104 Level: Advanced 1105 1106 .seealso: TSThetaSetTheta() 1107 @*/ 1108 PetscErrorCode TSThetaGetTheta(TS ts,PetscReal *theta) 1109 { 1110 PetscErrorCode ierr; 1111 1112 PetscFunctionBegin; 1113 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1114 PetscValidPointer(theta,2); 1115 ierr = PetscUseMethod(ts,"TSThetaGetTheta_C",(TS,PetscReal*),(ts,theta));CHKERRQ(ierr); 1116 PetscFunctionReturn(0); 1117 } 1118 1119 /*@ 1120 TSThetaSetTheta - Set the abscissa of the stage in (0,1]. 1121 1122 Not Collective 1123 1124 Input Parameter: 1125 + ts - timestepping context 1126 - theta - stage abscissa 1127 1128 Options Database: 1129 . -ts_theta_theta <theta> 1130 1131 Level: Intermediate 1132 1133 .seealso: TSThetaGetTheta() 1134 @*/ 1135 PetscErrorCode TSThetaSetTheta(TS ts,PetscReal theta) 1136 { 1137 PetscErrorCode ierr; 1138 1139 PetscFunctionBegin; 1140 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1141 ierr = PetscTryMethod(ts,"TSThetaSetTheta_C",(TS,PetscReal),(ts,theta));CHKERRQ(ierr); 1142 PetscFunctionReturn(0); 1143 } 1144 1145 /*@ 1146 TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 1147 1148 Not Collective 1149 1150 Input Parameter: 1151 . ts - timestepping context 1152 1153 Output Parameter: 1154 . endpoint - PETSC_TRUE when using the endpoint variant 1155 1156 Level: Advanced 1157 1158 .seealso: TSThetaSetEndpoint(), TSTHETA, TSCN 1159 @*/ 1160 PetscErrorCode TSThetaGetEndpoint(TS ts,PetscBool *endpoint) 1161 { 1162 PetscErrorCode ierr; 1163 1164 PetscFunctionBegin; 1165 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1166 PetscValidPointer(endpoint,2); 1167 ierr = PetscUseMethod(ts,"TSThetaGetEndpoint_C",(TS,PetscBool*),(ts,endpoint));CHKERRQ(ierr); 1168 PetscFunctionReturn(0); 1169 } 1170 1171 /*@ 1172 TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule). 1173 1174 Not Collective 1175 1176 Input Parameter: 1177 + ts - timestepping context 1178 - flg - PETSC_TRUE to use the endpoint variant 1179 1180 Options Database: 1181 . -ts_theta_endpoint <flg> 1182 1183 Level: Intermediate 1184 1185 .seealso: TSTHETA, TSCN 1186 @*/ 1187 PetscErrorCode TSThetaSetEndpoint(TS ts,PetscBool flg) 1188 { 1189 PetscErrorCode ierr; 1190 1191 PetscFunctionBegin; 1192 PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1193 ierr = PetscTryMethod(ts,"TSThetaSetEndpoint_C",(TS,PetscBool),(ts,flg));CHKERRQ(ierr); 1194 PetscFunctionReturn(0); 1195 } 1196 1197 /* 1198 * TSBEULER and TSCN are straightforward specializations of TSTHETA. 1199 * The creation functions for these specializations are below. 1200 */ 1201 1202 static PetscErrorCode TSSetUp_BEuler(TS ts) 1203 { 1204 TS_Theta *th = (TS_Theta*)ts->data; 1205 PetscErrorCode ierr; 1206 1207 PetscFunctionBegin; 1208 if (th->Theta != 1.0) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_OPT_OVERWRITE,"Can not change the default value (1) of theta when using backward Euler\n"); 1209 if (th->endpoint) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_OPT_OVERWRITE,"Can not change to the endpoint form of the Theta methods when using backward Euler\n"); 1210 ierr = TSSetUp_Theta(ts);CHKERRQ(ierr); 1211 PetscFunctionReturn(0); 1212 } 1213 1214 static PetscErrorCode TSView_BEuler(TS ts,PetscViewer viewer) 1215 { 1216 PetscFunctionBegin; 1217 PetscFunctionReturn(0); 1218 } 1219 1220 /*MC 1221 TSBEULER - ODE solver using the implicit backward Euler method 1222 1223 Level: beginner 1224 1225 Notes: 1226 TSBEULER is equivalent to TSTHETA with Theta=1.0 1227 1228 $ -ts_type theta -ts_theta_theta 1.0 1229 1230 .seealso: TSCreate(), TS, TSSetType(), TSEULER, TSCN, TSTHETA 1231 1232 M*/ 1233 PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts) 1234 { 1235 PetscErrorCode ierr; 1236 1237 PetscFunctionBegin; 1238 ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 1239 ierr = TSThetaSetTheta(ts,1.0);CHKERRQ(ierr); 1240 ierr = TSThetaSetEndpoint(ts,PETSC_FALSE);CHKERRQ(ierr); 1241 ts->ops->setup = TSSetUp_BEuler; 1242 ts->ops->view = TSView_BEuler; 1243 PetscFunctionReturn(0); 1244 } 1245 1246 static PetscErrorCode TSSetUp_CN(TS ts) 1247 { 1248 TS_Theta *th = (TS_Theta*)ts->data; 1249 PetscErrorCode ierr; 1250 1251 PetscFunctionBegin; 1252 if (th->Theta != 0.5) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_OPT_OVERWRITE,"Can not change the default value (0.5) of theta when using Crank-Nicolson\n"); 1253 if (!th->endpoint) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_OPT_OVERWRITE,"Can not change to the midpoint form of the Theta methods when using Crank-Nicolson\n"); 1254 ierr = TSSetUp_Theta(ts);CHKERRQ(ierr); 1255 PetscFunctionReturn(0); 1256 } 1257 1258 static PetscErrorCode TSView_CN(TS ts,PetscViewer viewer) 1259 { 1260 PetscFunctionBegin; 1261 PetscFunctionReturn(0); 1262 } 1263 1264 /*MC 1265 TSCN - ODE solver using the implicit Crank-Nicolson method. 1266 1267 Level: beginner 1268 1269 Notes: 1270 TSCN is equivalent to TSTHETA with Theta=0.5 and the "endpoint" option set. I.e. 1271 1272 $ -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint 1273 1274 .seealso: TSCreate(), TS, TSSetType(), TSBEULER, TSTHETA 1275 1276 M*/ 1277 PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts) 1278 { 1279 PetscErrorCode ierr; 1280 1281 PetscFunctionBegin; 1282 ierr = TSCreate_Theta(ts);CHKERRQ(ierr); 1283 ierr = TSThetaSetTheta(ts,0.5);CHKERRQ(ierr); 1284 ierr = TSThetaSetEndpoint(ts,PETSC_TRUE);CHKERRQ(ierr); 1285 ts->ops->setup = TSSetUp_CN; 1286 ts->ops->view = TSView_CN; 1287 PetscFunctionReturn(0); 1288 } 1289