1*2d3f70b5SBarry Smith #ifndef lint 2*2d3f70b5SBarry Smith static char vcid[] = "$Id: posindep.c,v 1.5 1996/09/14 12:37:20 bsmith Exp bsmith $"; 3*2d3f70b5SBarry Smith #endif 4*2d3f70b5SBarry Smith /* 5*2d3f70b5SBarry Smith Code for Time Stepping with implicit backwards Euler. 6*2d3f70b5SBarry Smith */ 7*2d3f70b5SBarry Smith #include <math.h> 8*2d3f70b5SBarry Smith #include "src/ts/tsimpl.h" /*I "ts.h" I*/ 9*2d3f70b5SBarry Smith #include "pinclude/pviewer.h" 10*2d3f70b5SBarry Smith 11*2d3f70b5SBarry Smith 12*2d3f70b5SBarry Smith typedef struct { 13*2d3f70b5SBarry Smith Vec update; /* work vector where new solution is formed */ 14*2d3f70b5SBarry Smith Vec func; /* work vector where F(t[i],u[i]) is stored */ 15*2d3f70b5SBarry Smith Vec rhs; /* work vector for RHS; vec_sol/dt */ 16*2d3f70b5SBarry Smith 17*2d3f70b5SBarry Smith /* information used for Pseudo-timestepping */ 18*2d3f70b5SBarry Smith 19*2d3f70b5SBarry Smith int (*dt)(TS,double*,void*); 20*2d3f70b5SBarry Smith void *dtctx; 21*2d3f70b5SBarry Smith double initial_fnorm,fnorm; /* original and current norm of F(u) */ 22*2d3f70b5SBarry Smith 23*2d3f70b5SBarry Smith } TS_PosIndep; 24*2d3f70b5SBarry Smith 25*2d3f70b5SBarry Smith /* ------------------------------------------------------------------------------*/ 26*2d3f70b5SBarry Smith /*@C 27*2d3f70b5SBarry Smith TSPseudoDefaultPosIndTimeStep - Default code to compute 28*2d3f70b5SBarry Smith pseudo timestepping. Use with TSPseudoSetPosIndTimeStep(). 29*2d3f70b5SBarry Smith 30*2d3f70b5SBarry Smith Input Parameters: 31*2d3f70b5SBarry Smith . ts - the time step context 32*2d3f70b5SBarry Smith . dtctx - 33*2d3f70b5SBarry Smith 34*2d3f70b5SBarry Smith Output Parameter: 35*2d3f70b5SBarry Smith . newdt - the time step to use for the next step 36*2d3f70b5SBarry Smith 37*2d3f70b5SBarry Smith @*/ 38*2d3f70b5SBarry Smith int TSPseudoDefaultPosIndTimeStep(TS ts,double* newdt,void* dtctx) 39*2d3f70b5SBarry Smith { 40*2d3f70b5SBarry Smith TS_PosIndep *posindep = (TS_PosIndep*) ts->data; 41*2d3f70b5SBarry Smith int ierr; 42*2d3f70b5SBarry Smith 43*2d3f70b5SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,ts->vec_sol,posindep->func);CHKERRQ(ierr); 44*2d3f70b5SBarry Smith ierr = VecNorm(posindep->func,NORM_2,&posindep->fnorm); CHKERRQ(ierr); 45*2d3f70b5SBarry Smith if (posindep->initial_fnorm == 0.0) { 46*2d3f70b5SBarry Smith /* first time through so compute initial function norm */ 47*2d3f70b5SBarry Smith posindep->initial_fnorm = posindep->fnorm; 48*2d3f70b5SBarry Smith } 49*2d3f70b5SBarry Smith if (posindep->fnorm == 0.0) *newdt = 1.e12*1.1*ts->time_step; 50*2d3f70b5SBarry Smith else *newdt = 1.1*ts->time_step*posindep->initial_fnorm/posindep->fnorm; 51*2d3f70b5SBarry Smith return 0; 52*2d3f70b5SBarry Smith } 53*2d3f70b5SBarry Smith 54*2d3f70b5SBarry Smith /*@ 55*2d3f70b5SBarry Smith TSPseudoSetPosIndTimeStep - Sets the user routine to be 56*2d3f70b5SBarry Smith called at each pseudo-time-step to update the time-step. 57*2d3f70b5SBarry Smith 58*2d3f70b5SBarry Smith Input Parameters: 59*2d3f70b5SBarry Smith . ts - time step context 60*2d3f70b5SBarry Smith . dt - function to compute timestep 61*2d3f70b5SBarry Smith . ctx - [optional] context required by function 62*2d3f70b5SBarry Smith 63*2d3f70b5SBarry Smith @*/ 64*2d3f70b5SBarry Smith int TSPseudoSetPosIndTimeStep(TS ts,int (*dt)(TS,double*,void*),void* ctx) 65*2d3f70b5SBarry Smith { 66*2d3f70b5SBarry Smith TS_PosIndep *posindep; 67*2d3f70b5SBarry Smith 68*2d3f70b5SBarry Smith PetscValidHeaderSpecific(ts,TS_COOKIE); 69*2d3f70b5SBarry Smith if (ts->type != TS_PSEUDO_POSIND) return 0; 70*2d3f70b5SBarry Smith posindep = (TS_PosIndep*) ts->data; 71*2d3f70b5SBarry Smith posindep->dt = dt; 72*2d3f70b5SBarry Smith posindep->dtctx = ctx; 73*2d3f70b5SBarry Smith return 0; 74*2d3f70b5SBarry Smith } 75*2d3f70b5SBarry Smith 76*2d3f70b5SBarry Smith /* 77*2d3f70b5SBarry Smith 78*2d3f70b5SBarry Smith */ 79*2d3f70b5SBarry Smith static int TSStep_PosIndep(TS ts,int *steps,double *time) 80*2d3f70b5SBarry Smith { 81*2d3f70b5SBarry Smith Vec sol = ts->vec_sol; 82*2d3f70b5SBarry Smith int ierr,i,max_steps = ts->max_steps,its; 83*2d3f70b5SBarry Smith TS_PosIndep *posindep = (TS_PosIndep*) ts->data; 84*2d3f70b5SBarry Smith 85*2d3f70b5SBarry Smith *steps = -ts->steps; 86*2d3f70b5SBarry Smith 87*2d3f70b5SBarry Smith ierr = VecCopy(sol,posindep->update); CHKERRQ(ierr); 88*2d3f70b5SBarry Smith for ( i=0; i<max_steps && ts->ptime < ts->max_time; i++ ) { 89*2d3f70b5SBarry Smith ierr = (*posindep->dt)(ts,&ts->time_step, posindep->dtctx); CHKERRQ(ierr); 90*2d3f70b5SBarry Smith ts->ptime += ts->time_step; 91*2d3f70b5SBarry Smith ierr = SNESSolve(ts->snes,posindep->update,&its); CHKERRQ(ierr); 92*2d3f70b5SBarry Smith ierr = VecCopy(posindep->update,sol); CHKERRQ(ierr); 93*2d3f70b5SBarry Smith ts->steps++; 94*2d3f70b5SBarry Smith ierr = TSMonitor(ts,ts->steps,ts->ptime,sol);CHKERRQ(ierr); 95*2d3f70b5SBarry Smith } 96*2d3f70b5SBarry Smith 97*2d3f70b5SBarry Smith *steps += ts->steps; 98*2d3f70b5SBarry Smith *time = ts->ptime; 99*2d3f70b5SBarry Smith return 0; 100*2d3f70b5SBarry Smith } 101*2d3f70b5SBarry Smith 102*2d3f70b5SBarry Smith /*------------------------------------------------------------*/ 103*2d3f70b5SBarry Smith static int TSDestroy_PosIndep(PetscObject obj ) 104*2d3f70b5SBarry Smith { 105*2d3f70b5SBarry Smith TS ts = (TS) obj; 106*2d3f70b5SBarry Smith TS_PosIndep *posindep = (TS_PosIndep*) ts->data; 107*2d3f70b5SBarry Smith int ierr; 108*2d3f70b5SBarry Smith 109*2d3f70b5SBarry Smith ierr = VecDestroy(posindep->update); CHKERRQ(ierr); 110*2d3f70b5SBarry Smith if (posindep->func) {ierr = VecDestroy(posindep->func);CHKERRQ(ierr);} 111*2d3f70b5SBarry Smith if (posindep->rhs) {ierr = VecDestroy(posindep->rhs);CHKERRQ(ierr);} 112*2d3f70b5SBarry Smith if (ts->Ashell) {ierr = MatDestroy(ts->A); CHKERRQ(ierr);} 113*2d3f70b5SBarry Smith PetscFree(posindep); 114*2d3f70b5SBarry Smith return 0; 115*2d3f70b5SBarry Smith } 116*2d3f70b5SBarry Smith 117*2d3f70b5SBarry Smith 118*2d3f70b5SBarry Smith /*------------------------------------------------------------*/ 119*2d3f70b5SBarry Smith /* 120*2d3f70b5SBarry Smith This matrix shell multiply where user provided Shell matrix 121*2d3f70b5SBarry Smith */ 122*2d3f70b5SBarry Smith 123*2d3f70b5SBarry Smith int TSPosIndepMatMult(Mat mat,Vec x,Vec y) 124*2d3f70b5SBarry Smith { 125*2d3f70b5SBarry Smith TS ts; 126*2d3f70b5SBarry Smith Scalar mdt,mone = -1.0; 127*2d3f70b5SBarry Smith int ierr; 128*2d3f70b5SBarry Smith 129*2d3f70b5SBarry Smith MatShellGetContext(mat,(void **)&ts); 130*2d3f70b5SBarry Smith mdt = 1.0/ts->time_step; 131*2d3f70b5SBarry Smith 132*2d3f70b5SBarry Smith /* apply user provided function */ 133*2d3f70b5SBarry Smith ierr = MatMult(ts->Ashell,x,y); CHKERRQ(ierr); 134*2d3f70b5SBarry Smith /* shift and scale by 1/dt - F */ 135*2d3f70b5SBarry Smith ierr = VecAXPBY(&mdt,&mone,x,y); CHKERRQ(ierr); 136*2d3f70b5SBarry Smith return 0; 137*2d3f70b5SBarry Smith } 138*2d3f70b5SBarry Smith 139*2d3f70b5SBarry Smith /* 140*2d3f70b5SBarry Smith This defines the nonlinear equation that is to be solved with SNES 141*2d3f70b5SBarry Smith 142*2d3f70b5SBarry Smith (U^{n+1} - U^{n})/dt - F(U^{n+1}) 143*2d3f70b5SBarry Smith */ 144*2d3f70b5SBarry Smith int TSPosIndepFunction(SNES snes,Vec x,Vec y,void *ctx) 145*2d3f70b5SBarry Smith { 146*2d3f70b5SBarry Smith TS ts = (TS) ctx; 147*2d3f70b5SBarry Smith Scalar mdt = 1.0/ts->time_step,*unp1,*un,*Funp1; 148*2d3f70b5SBarry Smith int ierr,i,n; 149*2d3f70b5SBarry Smith 150*2d3f70b5SBarry Smith /* apply user provided function */ 151*2d3f70b5SBarry Smith ierr = TSComputeRHSFunction(ts,ts->ptime,x,y); CHKERRQ(ierr); 152*2d3f70b5SBarry Smith /* (u^{n+1) - u^{n})/dt - F(u^{n+1}) */ 153*2d3f70b5SBarry Smith ierr = VecGetArray(ts->vec_sol,&un); CHKERRQ(ierr); 154*2d3f70b5SBarry Smith ierr = VecGetArray(x,&unp1); CHKERRQ(ierr); 155*2d3f70b5SBarry Smith ierr = VecGetArray(y,&Funp1); CHKERRQ(ierr); 156*2d3f70b5SBarry Smith ierr = VecGetLocalSize(x,&n); CHKERRQ(ierr); 157*2d3f70b5SBarry Smith for ( i=0; i<n; i++ ) { 158*2d3f70b5SBarry Smith Funp1[i] = mdt*(unp1[i] - un[i]) - Funp1[i]; 159*2d3f70b5SBarry Smith } 160*2d3f70b5SBarry Smith ierr = VecRestoreArray(ts->vec_sol,&un); 161*2d3f70b5SBarry Smith ierr = VecRestoreArray(x,&unp1); 162*2d3f70b5SBarry Smith ierr = VecRestoreArray(y,&Funp1); 163*2d3f70b5SBarry Smith 164*2d3f70b5SBarry Smith return 0; 165*2d3f70b5SBarry Smith } 166*2d3f70b5SBarry Smith 167*2d3f70b5SBarry Smith /* 168*2d3f70b5SBarry Smith This constructs the Jacobian needed for SNES 169*2d3f70b5SBarry Smith 170*2d3f70b5SBarry Smith J = I/dt - J_{F} where J_{F} is the given Jacobian of F. 171*2d3f70b5SBarry Smith */ 172*2d3f70b5SBarry Smith int TSPosIndepJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx) 173*2d3f70b5SBarry Smith { 174*2d3f70b5SBarry Smith TS ts = (TS) ctx; 175*2d3f70b5SBarry Smith int ierr; 176*2d3f70b5SBarry Smith Scalar mone = -1.0, mdt = 1.0/ts->time_step; 177*2d3f70b5SBarry Smith MatType mtype; 178*2d3f70b5SBarry Smith 179*2d3f70b5SBarry Smith /* construct users Jacobian */ 180*2d3f70b5SBarry Smith if (ts->rhsjacobian) { 181*2d3f70b5SBarry Smith ierr = (*ts->rhsjacobian)(ts,ts->ptime,x,AA,BB,str,ts->jacP);CHKERRQ(ierr); 182*2d3f70b5SBarry Smith } 183*2d3f70b5SBarry Smith 184*2d3f70b5SBarry Smith /* shift and scale Jacobian, if not a shell matrix */ 185*2d3f70b5SBarry Smith ierr = MatGetType(*AA,&mtype,PETSC_NULL); 186*2d3f70b5SBarry Smith if (mtype != MATSHELL) { 187*2d3f70b5SBarry Smith ierr = MatScale(&mone,*AA); CHKERRQ(ierr); 188*2d3f70b5SBarry Smith ierr = MatShift(&mdt,*AA); CHKERRQ(ierr); 189*2d3f70b5SBarry Smith } 190*2d3f70b5SBarry Smith ierr = MatGetType(*BB,&mtype,PETSC_NULL); 191*2d3f70b5SBarry Smith if (*BB != *AA && *str != SAME_PRECONDITIONER && mtype != MATSHELL) { 192*2d3f70b5SBarry Smith ierr = MatScale(&mone,*BB); CHKERRQ(ierr); 193*2d3f70b5SBarry Smith ierr = MatShift(&mdt,*BB); CHKERRQ(ierr); 194*2d3f70b5SBarry Smith } 195*2d3f70b5SBarry Smith 196*2d3f70b5SBarry Smith return 0; 197*2d3f70b5SBarry Smith } 198*2d3f70b5SBarry Smith 199*2d3f70b5SBarry Smith 200*2d3f70b5SBarry Smith static int TSSetUp_PosIndep(TS ts) 201*2d3f70b5SBarry Smith { 202*2d3f70b5SBarry Smith TS_PosIndep *posindep = (TS_PosIndep*) ts->data; 203*2d3f70b5SBarry Smith int ierr, M, m; 204*2d3f70b5SBarry Smith 205*2d3f70b5SBarry Smith ierr = VecDuplicate(ts->vec_sol,&posindep->update); CHKERRQ(ierr); 206*2d3f70b5SBarry Smith ierr = VecDuplicate(ts->vec_sol,&posindep->func); CHKERRQ(ierr); 207*2d3f70b5SBarry Smith ierr = SNESSetFunction(ts->snes,posindep->func,TSPosIndepFunction,ts);CHKERRQ(ierr); 208*2d3f70b5SBarry Smith if (ts->Ashell) { /* construct new shell matrix */ 209*2d3f70b5SBarry Smith ierr = VecGetSize(ts->vec_sol,&M); CHKERRQ(ierr); 210*2d3f70b5SBarry Smith ierr = VecGetLocalSize(ts->vec_sol,&m); CHKERRQ(ierr); 211*2d3f70b5SBarry Smith ierr = MatCreateShell(ts->comm,m,M,M,M,ts,&ts->A); CHKERRQ(ierr); 212*2d3f70b5SBarry Smith ierr = MatShellSetOperation(ts->A,MAT_MULT,(void*)TSPosIndepMatMult);CHKERRQ(ierr); 213*2d3f70b5SBarry Smith } 214*2d3f70b5SBarry Smith ierr = SNESSetJacobian(ts->snes,ts->A,ts->B,TSPosIndepJacobian,ts);CHKERRQ(ierr); 215*2d3f70b5SBarry Smith return 0; 216*2d3f70b5SBarry Smith } 217*2d3f70b5SBarry Smith /*------------------------------------------------------------*/ 218*2d3f70b5SBarry Smith 219*2d3f70b5SBarry Smith int TSPseudoDefaultMonitor(TS ts, int step, double time,Vec v, void *ctx) 220*2d3f70b5SBarry Smith { 221*2d3f70b5SBarry Smith TS_PosIndep *posindep = (TS_PosIndep*) ts->data; 222*2d3f70b5SBarry Smith 223*2d3f70b5SBarry Smith PetscPrintf(ts->comm,"TS %d dt %g time %g fnorm %g\n",step,ts->time_step,time, 224*2d3f70b5SBarry Smith posindep->fnorm); 225*2d3f70b5SBarry Smith return 0; 226*2d3f70b5SBarry Smith } 227*2d3f70b5SBarry Smith 228*2d3f70b5SBarry Smith 229*2d3f70b5SBarry Smith static int TSSetFromOptions_PosIndep(TS ts) 230*2d3f70b5SBarry Smith { 231*2d3f70b5SBarry Smith int ierr,flg; 232*2d3f70b5SBarry Smith 233*2d3f70b5SBarry Smith ierr = SNESSetFromOptions(ts->snes); CHKERRQ(ierr); 234*2d3f70b5SBarry Smith 235*2d3f70b5SBarry Smith ierr = OptionsHasName(ts->prefix,"-ts_monitor",&flg); CHKERRQ(ierr); 236*2d3f70b5SBarry Smith if (flg) { 237*2d3f70b5SBarry Smith TSSetMonitor(ts,TSPseudoDefaultMonitor,0); 238*2d3f70b5SBarry Smith } 239*2d3f70b5SBarry Smith return 0; 240*2d3f70b5SBarry Smith } 241*2d3f70b5SBarry Smith 242*2d3f70b5SBarry Smith static int TSPrintHelp_PosIndep(TS ts) 243*2d3f70b5SBarry Smith { 244*2d3f70b5SBarry Smith 245*2d3f70b5SBarry Smith return 0; 246*2d3f70b5SBarry Smith } 247*2d3f70b5SBarry Smith 248*2d3f70b5SBarry Smith static int TSView_PosIndep(PetscObject obj,Viewer viewer) 249*2d3f70b5SBarry Smith { 250*2d3f70b5SBarry Smith return 0; 251*2d3f70b5SBarry Smith } 252*2d3f70b5SBarry Smith 253*2d3f70b5SBarry Smith /* ------------------------------------------------------------ */ 254*2d3f70b5SBarry Smith int TSCreate_PosIndep(TS ts ) 255*2d3f70b5SBarry Smith { 256*2d3f70b5SBarry Smith TS_PosIndep *posindep; 257*2d3f70b5SBarry Smith int ierr; 258*2d3f70b5SBarry Smith MatType mtype; 259*2d3f70b5SBarry Smith 260*2d3f70b5SBarry Smith ts->type = TS_PSEUDO_POSIND; 261*2d3f70b5SBarry Smith ts->destroy = TSDestroy_PosIndep; 262*2d3f70b5SBarry Smith ts->printhelp = TSPrintHelp_PosIndep; 263*2d3f70b5SBarry Smith ts->view = TSView_PosIndep; 264*2d3f70b5SBarry Smith 265*2d3f70b5SBarry Smith if (ts->problem_type == TS_LINEAR) { 266*2d3f70b5SBarry Smith SETERRQ(1,"TSCreate_PosIndep:Only for nonlinear problems"); 267*2d3f70b5SBarry Smith } 268*2d3f70b5SBarry Smith if (!ts->A) { 269*2d3f70b5SBarry Smith SETERRQ(1,"TSCreate_PosIndep:Must set Jacobian"); 270*2d3f70b5SBarry Smith } 271*2d3f70b5SBarry Smith ierr = MatGetType(ts->A,&mtype,PETSC_NULL); 272*2d3f70b5SBarry Smith if (mtype == MATSHELL) { 273*2d3f70b5SBarry Smith ts->Ashell = ts->A; 274*2d3f70b5SBarry Smith } 275*2d3f70b5SBarry Smith ts->setup = TSSetUp_PosIndep; 276*2d3f70b5SBarry Smith ts->step = TSStep_PosIndep; 277*2d3f70b5SBarry Smith ts->setfromoptions = TSSetFromOptions_PosIndep; 278*2d3f70b5SBarry Smith ierr = SNESCreate(ts->comm,SNES_NONLINEAR_EQUATIONS,&ts->snes);CHKERRQ(ierr); 279*2d3f70b5SBarry Smith 280*2d3f70b5SBarry Smith posindep = PetscNew(TS_PosIndep); CHKPTRQ(posindep); 281*2d3f70b5SBarry Smith PetscMemzero(posindep,sizeof(TS_PosIndep)); 282*2d3f70b5SBarry Smith ts->data = (void *) posindep; 283*2d3f70b5SBarry Smith 284*2d3f70b5SBarry Smith return 0; 285*2d3f70b5SBarry Smith } 286*2d3f70b5SBarry Smith 287*2d3f70b5SBarry Smith 288*2d3f70b5SBarry Smith 289*2d3f70b5SBarry Smith 290*2d3f70b5SBarry Smith 291