1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Solves 1D heat equation with FEM formulation.\n\ 3c4762a1bSJed Brown Input arguments are\n\ 4c4762a1bSJed Brown -useAlhs: solve Alhs*U' = (Arhs*U + g) \n\ 5c4762a1bSJed Brown otherwise, solve U' = inv(Alhs)*(Arhs*U + g) \n\n"; 6c4762a1bSJed Brown 7c4762a1bSJed Brown /*-------------------------------------------------------------------------- 8c4762a1bSJed Brown Solves 1D heat equation U_t = U_xx with FEM formulation: 9c4762a1bSJed Brown Alhs*U' = rhs (= Arhs*U + g) 10c4762a1bSJed Brown We thank Chris Cox <clcox@clemson.edu> for contributing the original code 11c4762a1bSJed Brown ----------------------------------------------------------------------------*/ 12c4762a1bSJed Brown 13c4762a1bSJed Brown #include <petscksp.h> 14c4762a1bSJed Brown #include <petscts.h> 15c4762a1bSJed Brown 16c4762a1bSJed Brown /* special variable - max size of all arrays */ 17c4762a1bSJed Brown #define num_z 10 18c4762a1bSJed Brown 19c4762a1bSJed Brown /* 20c4762a1bSJed Brown User-defined application context - contains data needed by the 21c4762a1bSJed Brown application-provided call-back routines. 22c4762a1bSJed Brown */ 23c4762a1bSJed Brown typedef struct { 24c4762a1bSJed Brown Mat Amat; /* left hand side matrix */ 25c4762a1bSJed Brown Vec ksp_rhs,ksp_sol; /* working vectors for formulating inv(Alhs)*(Arhs*U+g) */ 26c4762a1bSJed Brown int max_probsz; /* max size of the problem */ 27c4762a1bSJed Brown PetscBool useAlhs; /* flag (1 indicates solving Alhs*U' = Arhs*U+g */ 28c4762a1bSJed Brown int nz; /* total number of grid points */ 29c4762a1bSJed Brown PetscInt m; /* total number of interio grid points */ 30c4762a1bSJed Brown Vec solution; /* global exact ts solution vector */ 31c4762a1bSJed Brown PetscScalar *z; /* array of grid points */ 32c4762a1bSJed Brown PetscBool debug; /* flag (1 indicates activation of debugging printouts) */ 33c4762a1bSJed Brown } AppCtx; 34c4762a1bSJed Brown 35c4762a1bSJed Brown extern PetscScalar exact(PetscScalar,PetscReal); 36c4762a1bSJed Brown extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*); 37c4762a1bSJed Brown extern PetscErrorCode Petsc_KSPSolve(AppCtx*); 38c4762a1bSJed Brown extern PetscScalar bspl(PetscScalar*,PetscScalar,PetscInt,PetscInt,PetscInt[][2],PetscInt); 39c4762a1bSJed Brown extern PetscErrorCode femBg(PetscScalar[][3],PetscScalar*,PetscInt,PetscScalar*,PetscReal); 40c4762a1bSJed Brown extern PetscErrorCode femA(AppCtx*,PetscInt,PetscScalar*); 41c4762a1bSJed Brown extern PetscErrorCode rhs(AppCtx*,PetscScalar*, PetscInt, PetscScalar*,PetscReal); 42c4762a1bSJed Brown extern PetscErrorCode RHSfunction(TS,PetscReal,Vec,Vec,void*); 43c4762a1bSJed Brown 44c4762a1bSJed Brown int main(int argc,char **argv) 45c4762a1bSJed Brown { 46c4762a1bSJed Brown PetscInt i,m,nz,steps,max_steps,k,nphase=1; 47c4762a1bSJed Brown PetscScalar zInitial,zFinal,val,*z; 48c4762a1bSJed Brown PetscReal stepsz[4],T,ftime; 49c4762a1bSJed Brown TS ts; 50c4762a1bSJed Brown SNES snes; 51c4762a1bSJed Brown Mat Jmat; 52c4762a1bSJed Brown AppCtx appctx; /* user-defined application context */ 53c4762a1bSJed Brown Vec init_sol; /* ts solution vector */ 54c4762a1bSJed Brown PetscMPIInt size; 55c4762a1bSJed Brown 56*327415f7SBarry Smith PetscFunctionBeginUser; 579566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 589566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 593c633725SBarry Smith PetscCheck(size == 1,PETSC_COMM_SELF,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only"); 60c4762a1bSJed Brown 61c4762a1bSJed Brown /* initializations */ 62c4762a1bSJed Brown zInitial = 0.0; 63c4762a1bSJed Brown zFinal = 1.0; 64c4762a1bSJed Brown nz = num_z; 65c4762a1bSJed Brown m = nz-2; 66c4762a1bSJed Brown appctx.nz = nz; 67c4762a1bSJed Brown max_steps = (PetscInt)10000; 68c4762a1bSJed Brown 69c4762a1bSJed Brown appctx.m = m; 70c4762a1bSJed Brown appctx.max_probsz = nz; 71c4762a1bSJed Brown appctx.debug = PETSC_FALSE; 72c4762a1bSJed Brown appctx.useAlhs = PETSC_FALSE; 73c4762a1bSJed Brown 74d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"",""); 759566063dSJacob Faibussowitsch PetscCall(PetscOptionsName("-debug",NULL,NULL,&appctx.debug)); 769566063dSJacob Faibussowitsch PetscCall(PetscOptionsName("-useAlhs",NULL,NULL,&appctx.useAlhs)); 779566063dSJacob Faibussowitsch PetscCall(PetscOptionsRangeInt("-nphase",NULL,NULL,nphase,&nphase,NULL,1,3)); 78d0609cedSBarry Smith PetscOptionsEnd(); 79303a5415SBarry Smith T = 0.014/nphase; 80303a5415SBarry Smith 81c4762a1bSJed Brown /* create vector to hold ts solution */ 82c4762a1bSJed Brown /*-----------------------------------*/ 839566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &init_sol)); 849566063dSJacob Faibussowitsch PetscCall(VecSetSizes(init_sol, PETSC_DECIDE, m)); 859566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(init_sol)); 86c4762a1bSJed Brown 87c4762a1bSJed Brown /* create vector to hold true ts soln for comparison */ 889566063dSJacob Faibussowitsch PetscCall(VecDuplicate(init_sol, &appctx.solution)); 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* create LHS matrix Amat */ 91c4762a1bSJed Brown /*------------------------*/ 929566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD, m, m, 3, NULL, &appctx.Amat)); 939566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(appctx.Amat)); 949566063dSJacob Faibussowitsch PetscCall(MatSetUp(appctx.Amat)); 95c4762a1bSJed Brown /* set space grid points - interio points only! */ 969566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nz+1,&z)); 97c4762a1bSJed Brown for (i=0; i<nz; i++) z[i]=(i)*((zFinal-zInitial)/(nz-1)); 98c4762a1bSJed Brown appctx.z = z; 99c4762a1bSJed Brown femA(&appctx,nz,z); 100c4762a1bSJed Brown 101c4762a1bSJed Brown /* create the jacobian matrix */ 102c4762a1bSJed Brown /*----------------------------*/ 1039566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &Jmat)); 1049566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Jmat,PETSC_DECIDE,PETSC_DECIDE,m,m)); 1059566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(Jmat)); 1069566063dSJacob Faibussowitsch PetscCall(MatSetUp(Jmat)); 107c4762a1bSJed Brown 108c4762a1bSJed Brown /* create working vectors for formulating rhs=inv(Alhs)*(Arhs*U + g) */ 1099566063dSJacob Faibussowitsch PetscCall(VecDuplicate(init_sol,&appctx.ksp_rhs)); 1109566063dSJacob Faibussowitsch PetscCall(VecDuplicate(init_sol,&appctx.ksp_sol)); 111c4762a1bSJed Brown 1122d4ee042Sprj- /* set initial guess */ 1132d4ee042Sprj- /*-------------------*/ 114c4762a1bSJed Brown for (i=0; i<nz-2; i++) { 115c4762a1bSJed Brown val = exact(z[i+1], 0.0); 1169566063dSJacob Faibussowitsch PetscCall(VecSetValue(init_sol,i,(PetscScalar)val,INSERT_VALUES)); 117c4762a1bSJed Brown } 1189566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(init_sol)); 1199566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(init_sol)); 120c4762a1bSJed Brown 121c4762a1bSJed Brown /*create a time-stepping context and set the problem type */ 122c4762a1bSJed Brown /*--------------------------------------------------------*/ 1239566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 1249566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts,TS_NONLINEAR)); 125c4762a1bSJed Brown 126c4762a1bSJed Brown /* set time-step method */ 1279566063dSJacob Faibussowitsch PetscCall(TSSetType(ts,TSCN)); 128c4762a1bSJed Brown 129c4762a1bSJed Brown /* Set optional user-defined monitoring routine */ 1309566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts,Monitor,&appctx,NULL)); 131c4762a1bSJed Brown /* set the right hand side of U_t = RHSfunction(U,t) */ 1329566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts,NULL,(PetscErrorCode (*)(TS,PetscScalar,Vec,Vec,void*))RHSfunction,&appctx)); 133c4762a1bSJed Brown 134c4762a1bSJed Brown if (appctx.useAlhs) { 135c4762a1bSJed Brown /* set the left hand side matrix of Amat*U_t = rhs(U,t) */ 136c4762a1bSJed Brown 137c4762a1bSJed Brown /* Note: this approach is incompatible with the finite differenced Jacobian set below because we can't restore the 138c4762a1bSJed Brown * Alhs matrix without making a copy. Either finite difference the entire thing or use analytic Jacobians in both 139c4762a1bSJed Brown * places. 140c4762a1bSJed Brown */ 1419566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts,NULL,TSComputeIFunctionLinear,&appctx)); 1429566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts,appctx.Amat,appctx.Amat,TSComputeIJacobianConstant,&appctx)); 143c4762a1bSJed Brown } 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* use petsc to compute the jacobian by finite differences */ 1469566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts,&snes)); 1479566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes,Jmat,Jmat,SNESComputeJacobianDefault,NULL)); 148c4762a1bSJed Brown 149c4762a1bSJed Brown /* get the command line options if there are any and set them */ 1509566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 151c4762a1bSJed Brown 152e808b789SPatrick Sanan #if defined(PETSC_HAVE_SUNDIALS2) 153c4762a1bSJed Brown { 154c4762a1bSJed Brown TSType type; 155c4762a1bSJed Brown PetscBool sundialstype=PETSC_FALSE; 1569566063dSJacob Faibussowitsch PetscCall(TSGetType(ts,&type)); 1579566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)ts,TSSUNDIALS,&sundialstype)); 1583c633725SBarry Smith PetscCheck(!sundialstype || !appctx.useAlhs,PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use Alhs formulation for TSSUNDIALS type"); 159c4762a1bSJed Brown } 160c4762a1bSJed Brown #endif 161c4762a1bSJed Brown /* Sets the initial solution */ 1629566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts,init_sol)); 163c4762a1bSJed Brown 164c4762a1bSJed Brown stepsz[0] = 1.0/(2.0*(nz-1)*(nz-1)); /* (mesh_size)^2/2.0 */ 165c4762a1bSJed Brown ftime = 0.0; 166c4762a1bSJed Brown for (k=0; k<nphase; k++) { 16763a3b9bcSJacob Faibussowitsch if (nphase > 1) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Phase %" PetscInt_FMT " initial time %g, stepsz %g, duration: %g\n",k,(double)ftime,(double)stepsz[k],(double)((k+1)*T))); 1689566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts,ftime)); 1699566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,stepsz[k])); 1709566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts,max_steps)); 1719566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts,(k+1)*T)); 1729566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 173c4762a1bSJed Brown 174c4762a1bSJed Brown /* loop over time steps */ 175c4762a1bSJed Brown /*----------------------*/ 1769566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,init_sol)); 1779566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts,&ftime)); 1789566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts,&steps)); 179c4762a1bSJed Brown stepsz[k+1] = stepsz[k]*1.5; /* change step size for the next phase */ 180c4762a1bSJed Brown } 181c4762a1bSJed Brown 182c4762a1bSJed Brown /* free space */ 1839566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 1849566063dSJacob Faibussowitsch PetscCall(MatDestroy(&appctx.Amat)); 1859566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Jmat)); 1869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.ksp_rhs)); 1879566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.ksp_sol)); 1889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&init_sol)); 1899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&appctx.solution)); 1909566063dSJacob Faibussowitsch PetscCall(PetscFree(z)); 191c4762a1bSJed Brown 1929566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 193b122ec5aSJacob Faibussowitsch return 0; 194c4762a1bSJed Brown } 195c4762a1bSJed Brown 196c4762a1bSJed Brown /*------------------------------------------------------------------------ 197c4762a1bSJed Brown Set exact solution 198c4762a1bSJed Brown u(z,t) = sin(6*PI*z)*exp(-36.*PI*PI*t) + 3.*sin(2*PI*z)*exp(-4.*PI*PI*t) 199c4762a1bSJed Brown --------------------------------------------------------------------------*/ 200c4762a1bSJed Brown PetscScalar exact(PetscScalar z,PetscReal t) 201c4762a1bSJed Brown { 202c4762a1bSJed Brown PetscScalar val, ex1, ex2; 203c4762a1bSJed Brown 204c4762a1bSJed Brown ex1 = PetscExpReal(-36.*PETSC_PI*PETSC_PI*t); 205c4762a1bSJed Brown ex2 = PetscExpReal(-4.*PETSC_PI*PETSC_PI*t); 206c4762a1bSJed Brown val = PetscSinScalar(6*PETSC_PI*z)*ex1 + 3.*PetscSinScalar(2*PETSC_PI*z)*ex2; 207c4762a1bSJed Brown return val; 208c4762a1bSJed Brown } 209c4762a1bSJed Brown 210c4762a1bSJed Brown /* 211c4762a1bSJed Brown Monitor - User-provided routine to monitor the solution computed at 212c4762a1bSJed Brown each timestep. This example plots the solution and computes the 213c4762a1bSJed Brown error in two different norms. 214c4762a1bSJed Brown 215c4762a1bSJed Brown Input Parameters: 216c4762a1bSJed Brown ts - the timestep context 217c4762a1bSJed Brown step - the count of the current step (with 0 meaning the 218c4762a1bSJed Brown initial condition) 219c4762a1bSJed Brown time - the current time 220c4762a1bSJed Brown u - the solution at this timestep 221c4762a1bSJed Brown ctx - the user-provided context for this monitoring routine. 222c4762a1bSJed Brown In this case we use the application context which contains 223c4762a1bSJed Brown information about the problem size, workspace and the exact 224c4762a1bSJed Brown solution. 225c4762a1bSJed Brown */ 226c4762a1bSJed Brown PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec u,void *ctx) 227c4762a1bSJed Brown { 228c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; 229c4762a1bSJed Brown PetscInt i,m=appctx->m; 230c4762a1bSJed Brown PetscReal norm_2,norm_max,h=1.0/(m+1); 231c4762a1bSJed Brown PetscScalar *u_exact; 232c4762a1bSJed Brown 233c4762a1bSJed Brown /* Compute the exact solution */ 2349566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(appctx->solution,&u_exact)); 235c4762a1bSJed Brown for (i=0; i<m; i++) u_exact[i] = exact(appctx->z[i+1],time); 2369566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(appctx->solution,&u_exact)); 237c4762a1bSJed Brown 238c4762a1bSJed Brown /* Print debugging information if desired */ 239c4762a1bSJed Brown if (appctx->debug) { 2409566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Computed solution vector at time %g\n",(double)time)); 2419566063dSJacob Faibussowitsch PetscCall(VecView(u,PETSC_VIEWER_STDOUT_SELF)); 2429566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Exact solution vector\n")); 2439566063dSJacob Faibussowitsch PetscCall(VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF)); 244c4762a1bSJed Brown } 245c4762a1bSJed Brown 246c4762a1bSJed Brown /* Compute the 2-norm and max-norm of the error */ 2479566063dSJacob Faibussowitsch PetscCall(VecAXPY(appctx->solution,-1.0,u)); 2489566063dSJacob Faibussowitsch PetscCall(VecNorm(appctx->solution,NORM_2,&norm_2)); 249c4762a1bSJed Brown 250c4762a1bSJed Brown norm_2 = PetscSqrtReal(h)*norm_2; 2519566063dSJacob Faibussowitsch PetscCall(VecNorm(appctx->solution,NORM_MAX,&norm_max)); 25263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Timestep %" PetscInt_FMT ": time = %g, 2-norm error = %6.4f, max norm error = %6.4f\n",step,(double)time,(double)norm_2,(double)norm_max)); 253c4762a1bSJed Brown 254c4762a1bSJed Brown /* 255c4762a1bSJed Brown Print debugging information if desired 256c4762a1bSJed Brown */ 257c4762a1bSJed Brown if (appctx->debug) { 2589566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error vector\n")); 2599566063dSJacob Faibussowitsch PetscCall(VecView(appctx->solution,PETSC_VIEWER_STDOUT_SELF)); 260c4762a1bSJed Brown } 261c4762a1bSJed Brown return 0; 262c4762a1bSJed Brown } 263c4762a1bSJed Brown 264c4762a1bSJed Brown /*%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% 2650e3d61c9SBarry Smith Function to solve a linear system using KSP 266c4762a1bSJed Brown %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%*/ 267c4762a1bSJed Brown 268c4762a1bSJed Brown PetscErrorCode Petsc_KSPSolve(AppCtx *obj) 269c4762a1bSJed Brown { 270c4762a1bSJed Brown KSP ksp; 271c4762a1bSJed Brown PC pc; 272c4762a1bSJed Brown 273c4762a1bSJed Brown /*create the ksp context and set the operators,that is, associate the system matrix with it*/ 2749566063dSJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_WORLD,&ksp)); 2759566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(ksp,obj->Amat,obj->Amat)); 276c4762a1bSJed Brown 277c4762a1bSJed Brown /*get the preconditioner context, set its type and the tolerances*/ 2789566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp,&pc)); 2799566063dSJacob Faibussowitsch PetscCall(PCSetType(pc,PCLU)); 2809566063dSJacob Faibussowitsch PetscCall(KSPSetTolerances(ksp,1.e-7,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT)); 281c4762a1bSJed Brown 282c4762a1bSJed Brown /*get the command line options if there are any and set them*/ 2839566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(ksp)); 284c4762a1bSJed Brown 285c4762a1bSJed Brown /*get the linear system (ksp) solve*/ 2869566063dSJacob Faibussowitsch PetscCall(KSPSolve(ksp,obj->ksp_rhs,obj->ksp_sol)); 287c4762a1bSJed Brown 2889566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&ksp)); 289c4762a1bSJed Brown return 0; 290c4762a1bSJed Brown } 291c4762a1bSJed Brown 292c4762a1bSJed Brown /*********************************************************************** 2930e3d61c9SBarry Smith Function to return value of basis function or derivative of basis function. 294c4762a1bSJed Brown *********************************************************************** 2950e3d61c9SBarry Smith 2960e3d61c9SBarry Smith Arguments: 2970e3d61c9SBarry Smith x = array of xpoints or nodal values 2980e3d61c9SBarry Smith xx = point at which the basis function is to be 2990e3d61c9SBarry Smith evaluated. 3000e3d61c9SBarry Smith il = interval containing xx. 3010e3d61c9SBarry Smith iq = indicates which of the two basis functions in 3020e3d61c9SBarry Smith interval intrvl should be used 3030e3d61c9SBarry Smith nll = array containing the endpoints of each interval. 3040e3d61c9SBarry Smith id = If id ~= 2, the value of the basis function 3050e3d61c9SBarry Smith is calculated; if id = 2, the value of the 3060e3d61c9SBarry Smith derivative of the basis function is returned. 307c4762a1bSJed Brown ***********************************************************************/ 308c4762a1bSJed Brown 309c4762a1bSJed Brown PetscScalar bspl(PetscScalar *x, PetscScalar xx,PetscInt il,PetscInt iq,PetscInt nll[][2],PetscInt id) 310c4762a1bSJed Brown { 311c4762a1bSJed Brown PetscScalar x1,x2,bfcn; 312c4762a1bSJed Brown PetscInt i1,i2,iq1,iq2; 313c4762a1bSJed Brown 3140e3d61c9SBarry Smith /* Determine which basis function in interval intrvl is to be used in */ 315c4762a1bSJed Brown iq1 = iq; 316c4762a1bSJed Brown if (iq1==0) iq2 = 1; 317c4762a1bSJed Brown else iq2 = 0; 318c4762a1bSJed Brown 3190e3d61c9SBarry Smith /* Determine endpoint of the interval intrvl */ 320c4762a1bSJed Brown i1=nll[il][iq1]; 321c4762a1bSJed Brown i2=nll[il][iq2]; 322c4762a1bSJed Brown 3230e3d61c9SBarry Smith /* Determine nodal values at the endpoints of the interval intrvl */ 324c4762a1bSJed Brown x1=x[i1]; 325c4762a1bSJed Brown x2=x[i2]; 326303a5415SBarry Smith 3270e3d61c9SBarry Smith /* Evaluate basis function */ 328c4762a1bSJed Brown if (id == 2) bfcn=(1.0)/(x1-x2); 329c4762a1bSJed Brown else bfcn=(xx-x2)/(x1-x2); 330c4762a1bSJed Brown return bfcn; 331c4762a1bSJed Brown } 332c4762a1bSJed Brown 333c4762a1bSJed Brown /*--------------------------------------------------------- 334c4762a1bSJed Brown Function called by rhs function to get B and g 335c4762a1bSJed Brown ---------------------------------------------------------*/ 336c4762a1bSJed Brown PetscErrorCode femBg(PetscScalar btri[][3],PetscScalar *f,PetscInt nz,PetscScalar *z, PetscReal t) 337c4762a1bSJed Brown { 338c4762a1bSJed Brown PetscInt i,j,jj,il,ip,ipp,ipq,iq,iquad,iqq; 339c4762a1bSJed Brown PetscInt nli[num_z][2],indx[num_z]; 340c4762a1bSJed Brown PetscScalar dd,dl,zip,zipq,zz,b_z,bb_z,bij; 341c4762a1bSJed Brown PetscScalar zquad[num_z][3],dlen[num_z],qdwt[3]; 342c4762a1bSJed Brown 343c4762a1bSJed Brown /* initializing everything - btri and f are initialized in rhs.c */ 344c4762a1bSJed Brown for (i=0; i < nz; i++) { 345c4762a1bSJed Brown nli[i][0] = 0; 346c4762a1bSJed Brown nli[i][1] = 0; 347c4762a1bSJed Brown indx[i] = 0; 348c4762a1bSJed Brown zquad[i][0] = 0.0; 349c4762a1bSJed Brown zquad[i][1] = 0.0; 350c4762a1bSJed Brown zquad[i][2] = 0.0; 351c4762a1bSJed Brown dlen[i] = 0.0; 352c4762a1bSJed Brown } /*end for (i)*/ 353c4762a1bSJed Brown 354c4762a1bSJed Brown /* quadrature weights */ 355c4762a1bSJed Brown qdwt[0] = 1.0/6.0; 356c4762a1bSJed Brown qdwt[1] = 4.0/6.0; 357c4762a1bSJed Brown qdwt[2] = 1.0/6.0; 358c4762a1bSJed Brown 359c4762a1bSJed Brown /* 1st and last nodes have Dirichlet boundary condition - 360c4762a1bSJed Brown set indices there to -1 */ 361c4762a1bSJed Brown 362c4762a1bSJed Brown for (i=0; i < nz-1; i++) indx[i] = i-1; 363c4762a1bSJed Brown indx[nz-1] = -1; 364c4762a1bSJed Brown 365c4762a1bSJed Brown ipq = 0; 366c4762a1bSJed Brown for (il=0; il < nz-1; il++) { 367c4762a1bSJed Brown ip = ipq; 368c4762a1bSJed Brown ipq = ip+1; 369c4762a1bSJed Brown zip = z[ip]; 370c4762a1bSJed Brown zipq = z[ipq]; 371c4762a1bSJed Brown dl = zipq-zip; 372c4762a1bSJed Brown zquad[il][0] = zip; 373c4762a1bSJed Brown zquad[il][1] = (0.5)*(zip+zipq); 374c4762a1bSJed Brown zquad[il][2] = zipq; 375c4762a1bSJed Brown dlen[il] = PetscAbsScalar(dl); 376c4762a1bSJed Brown nli[il][0] = ip; 377c4762a1bSJed Brown nli[il][1] = ipq; 378c4762a1bSJed Brown } 379c4762a1bSJed Brown 380c4762a1bSJed Brown for (il=0; il < nz-1; il++) { 381c4762a1bSJed Brown for (iquad=0; iquad < 3; iquad++) { 382c4762a1bSJed Brown dd = (dlen[il])*(qdwt[iquad]); 383c4762a1bSJed Brown zz = zquad[il][iquad]; 384c4762a1bSJed Brown 385c4762a1bSJed Brown for (iq=0; iq < 2; iq++) { 386c4762a1bSJed Brown ip = nli[il][iq]; 387c4762a1bSJed Brown b_z = bspl(z,zz,il,iq,nli,2); 388c4762a1bSJed Brown i = indx[ip]; 389c4762a1bSJed Brown 390c4762a1bSJed Brown if (i > -1) { 391c4762a1bSJed Brown for (iqq=0; iqq < 2; iqq++) { 392c4762a1bSJed Brown ipp = nli[il][iqq]; 393c4762a1bSJed Brown bb_z = bspl(z,zz,il,iqq,nli,2); 394c4762a1bSJed Brown j = indx[ipp]; 395c4762a1bSJed Brown bij = -b_z*bb_z; 396c4762a1bSJed Brown 397c4762a1bSJed Brown if (j > -1) { 398c4762a1bSJed Brown jj = 1+j-i; 399c4762a1bSJed Brown btri[i][jj] += bij*dd; 400c4762a1bSJed Brown } else { 401c4762a1bSJed Brown f[i] += bij*dd*exact(z[ipp], t); 402c4762a1bSJed Brown /* f[i] += 0.0; */ 403c4762a1bSJed Brown /* if (il==0 && j==-1) { */ 404c4762a1bSJed Brown /* f[i] += bij*dd*exact(zz,t); */ 405c4762a1bSJed Brown /* }*/ /*end if*/ 406c4762a1bSJed Brown } /*end else*/ 407c4762a1bSJed Brown } /*end for (iqq)*/ 408c4762a1bSJed Brown } /*end if (i>0)*/ 409c4762a1bSJed Brown } /*end for (iq)*/ 410c4762a1bSJed Brown } /*end for (iquad)*/ 411c4762a1bSJed Brown } /*end for (il)*/ 412c4762a1bSJed Brown return 0; 413c4762a1bSJed Brown } 414c4762a1bSJed Brown 415c4762a1bSJed Brown PetscErrorCode femA(AppCtx *obj,PetscInt nz,PetscScalar *z) 416c4762a1bSJed Brown { 417c4762a1bSJed Brown PetscInt i,j,il,ip,ipp,ipq,iq,iquad,iqq; 418c4762a1bSJed Brown PetscInt nli[num_z][2],indx[num_z]; 419c4762a1bSJed Brown PetscScalar dd,dl,zip,zipq,zz,bb,bbb,aij; 420c4762a1bSJed Brown PetscScalar rquad[num_z][3],dlen[num_z],qdwt[3],add_term; 421c4762a1bSJed Brown 422c4762a1bSJed Brown /* initializing everything */ 423c4762a1bSJed Brown for (i=0; i < nz; i++) { 424c4762a1bSJed Brown nli[i][0] = 0; 425c4762a1bSJed Brown nli[i][1] = 0; 426c4762a1bSJed Brown indx[i] = 0; 427c4762a1bSJed Brown rquad[i][0] = 0.0; 428c4762a1bSJed Brown rquad[i][1] = 0.0; 429c4762a1bSJed Brown rquad[i][2] = 0.0; 430c4762a1bSJed Brown dlen[i] = 0.0; 431c4762a1bSJed Brown } /*end for (i)*/ 432c4762a1bSJed Brown 433c4762a1bSJed Brown /* quadrature weights */ 434c4762a1bSJed Brown qdwt[0] = 1.0/6.0; 435c4762a1bSJed Brown qdwt[1] = 4.0/6.0; 436c4762a1bSJed Brown qdwt[2] = 1.0/6.0; 437c4762a1bSJed Brown 438c4762a1bSJed Brown /* 1st and last nodes have Dirichlet boundary condition - 439c4762a1bSJed Brown set indices there to -1 */ 440c4762a1bSJed Brown 441c4762a1bSJed Brown for (i=0; i < nz-1; i++) indx[i]=i-1; 442c4762a1bSJed Brown indx[nz-1]=-1; 443c4762a1bSJed Brown 444c4762a1bSJed Brown ipq = 0; 445c4762a1bSJed Brown 446c4762a1bSJed Brown for (il=0; il < nz-1; il++) { 447c4762a1bSJed Brown ip = ipq; 448c4762a1bSJed Brown ipq = ip+1; 449c4762a1bSJed Brown zip = z[ip]; 450c4762a1bSJed Brown zipq = z[ipq]; 451c4762a1bSJed Brown dl = zipq-zip; 452c4762a1bSJed Brown rquad[il][0] = zip; 453c4762a1bSJed Brown rquad[il][1] = (0.5)*(zip+zipq); 454c4762a1bSJed Brown rquad[il][2] = zipq; 455c4762a1bSJed Brown dlen[il] = PetscAbsScalar(dl); 456c4762a1bSJed Brown nli[il][0] = ip; 457c4762a1bSJed Brown nli[il][1] = ipq; 458c4762a1bSJed Brown } /*end for (il)*/ 459c4762a1bSJed Brown 460c4762a1bSJed Brown for (il=0; il < nz-1; il++) { 461c4762a1bSJed Brown for (iquad=0; iquad < 3; iquad++) { 462c4762a1bSJed Brown dd = (dlen[il])*(qdwt[iquad]); 463c4762a1bSJed Brown zz = rquad[il][iquad]; 464c4762a1bSJed Brown 465c4762a1bSJed Brown for (iq=0; iq < 2; iq++) { 466c4762a1bSJed Brown ip = nli[il][iq]; 467c4762a1bSJed Brown bb = bspl(z,zz,il,iq,nli,1); 468c4762a1bSJed Brown i = indx[ip]; 469c4762a1bSJed Brown if (i > -1) { 470c4762a1bSJed Brown for (iqq=0; iqq < 2; iqq++) { 471c4762a1bSJed Brown ipp = nli[il][iqq]; 472c4762a1bSJed Brown bbb = bspl(z,zz,il,iqq,nli,1); 473c4762a1bSJed Brown j = indx[ipp]; 474c4762a1bSJed Brown aij = bb*bbb; 475c4762a1bSJed Brown if (j > -1) { 476c4762a1bSJed Brown add_term = aij*dd; 4779566063dSJacob Faibussowitsch PetscCall(MatSetValue(obj->Amat,i,j,add_term,ADD_VALUES)); 478c4762a1bSJed Brown }/*endif*/ 479c4762a1bSJed Brown } /*end for (iqq)*/ 480c4762a1bSJed Brown } /*end if (i>0)*/ 481c4762a1bSJed Brown } /*end for (iq)*/ 482c4762a1bSJed Brown } /*end for (iquad)*/ 483c4762a1bSJed Brown } /*end for (il)*/ 4849566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(obj->Amat,MAT_FINAL_ASSEMBLY)); 4859566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(obj->Amat,MAT_FINAL_ASSEMBLY)); 486c4762a1bSJed Brown return 0; 487c4762a1bSJed Brown } 488c4762a1bSJed Brown 489c4762a1bSJed Brown /*--------------------------------------------------------- 490c4762a1bSJed Brown Function to fill the rhs vector with 491c4762a1bSJed Brown By + g values **** 492c4762a1bSJed Brown ---------------------------------------------------------*/ 493c4762a1bSJed Brown PetscErrorCode rhs(AppCtx *obj,PetscScalar *y, PetscInt nz, PetscScalar *z, PetscReal t) 494c4762a1bSJed Brown { 495c4762a1bSJed Brown PetscInt i,j,js,je,jj; 496c4762a1bSJed Brown PetscScalar val,g[num_z],btri[num_z][3],add_term; 497c4762a1bSJed Brown 498c4762a1bSJed Brown for (i=0; i < nz-2; i++) { 499c4762a1bSJed Brown for (j=0; j <= 2; j++) btri[i][j]=0.0; 500c4762a1bSJed Brown g[i] = 0.0; 501c4762a1bSJed Brown } 502c4762a1bSJed Brown 503c4762a1bSJed Brown /* call femBg to set the tri-diagonal b matrix and vector g */ 504c4762a1bSJed Brown femBg(btri,g,nz,z,t); 505c4762a1bSJed Brown 506c4762a1bSJed Brown /* setting the entries of the right hand side vector */ 507c4762a1bSJed Brown for (i=0; i < nz-2; i++) { 508c4762a1bSJed Brown val = 0.0; 509c4762a1bSJed Brown js = 0; 510c4762a1bSJed Brown if (i == 0) js = 1; 511c4762a1bSJed Brown je = 2; 512c4762a1bSJed Brown if (i == nz-2) je = 1; 513c4762a1bSJed Brown 514c4762a1bSJed Brown for (jj=js; jj <= je; jj++) { 515c4762a1bSJed Brown j = i+jj-1; 516c4762a1bSJed Brown val += (btri[i][jj])*(y[j]); 517c4762a1bSJed Brown } 518c4762a1bSJed Brown add_term = val + g[i]; 5199566063dSJacob Faibussowitsch PetscCall(VecSetValue(obj->ksp_rhs,(PetscInt)i,(PetscScalar)add_term,INSERT_VALUES)); 520c4762a1bSJed Brown } 5219566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(obj->ksp_rhs)); 5229566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(obj->ksp_rhs)); 523c4762a1bSJed Brown return 0; 524c4762a1bSJed Brown } 525c4762a1bSJed Brown 526c4762a1bSJed Brown /*%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% 527c4762a1bSJed Brown %% Function to form the right hand side of the time-stepping problem. %% 528c4762a1bSJed Brown %% -------------------------------------------------------------------------------------------%% 529c4762a1bSJed Brown if (useAlhs): 530c4762a1bSJed Brown globalout = By+g 531c4762a1bSJed Brown else if (!useAlhs): 532c4762a1bSJed Brown globalout = f(y,t)=Ainv(By+g), 533c4762a1bSJed Brown in which the ksp solver to transform the problem A*ydot=By+g 534c4762a1bSJed Brown to the problem ydot=f(y,t)=inv(A)*(By+g) 535c4762a1bSJed Brown %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%*/ 536c4762a1bSJed Brown 537c4762a1bSJed Brown PetscErrorCode RHSfunction(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx) 538c4762a1bSJed Brown { 539c4762a1bSJed Brown AppCtx *obj = (AppCtx*)ctx; 540c4762a1bSJed Brown PetscScalar soln[num_z]; 541c4762a1bSJed Brown const PetscScalar *soln_ptr; 542c4762a1bSJed Brown PetscInt i,nz=obj->nz; 543c4762a1bSJed Brown PetscReal time; 544c4762a1bSJed Brown 545c4762a1bSJed Brown /* get the previous solution to compute updated system */ 5469566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(globalin,&soln_ptr)); 547c4762a1bSJed Brown for (i=0; i < num_z-2; i++) soln[i] = soln_ptr[i]; 5489566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(globalin,&soln_ptr)); 549c4762a1bSJed Brown soln[num_z-1] = 0.0; 550c4762a1bSJed Brown soln[num_z-2] = 0.0; 551c4762a1bSJed Brown 552c4762a1bSJed Brown /* clear out the matrix and rhs for ksp to keep things straight */ 5539566063dSJacob Faibussowitsch PetscCall(VecSet(obj->ksp_rhs,(PetscScalar)0.0)); 554c4762a1bSJed Brown 555c4762a1bSJed Brown time = t; 556c4762a1bSJed Brown /* get the updated system */ 557c4762a1bSJed Brown rhs(obj,soln,nz,obj->z,time); /* setup of the By+g rhs */ 558c4762a1bSJed Brown 559c4762a1bSJed Brown /* do a ksp solve to get the rhs for the ts problem */ 560c4762a1bSJed Brown if (obj->useAlhs) { 561c4762a1bSJed Brown /* ksp_sol = ksp_rhs */ 5629566063dSJacob Faibussowitsch PetscCall(VecCopy(obj->ksp_rhs,globalout)); 563c4762a1bSJed Brown } else { 564c4762a1bSJed Brown /* ksp_sol = inv(Amat)*ksp_rhs */ 5659566063dSJacob Faibussowitsch PetscCall(Petsc_KSPSolve(obj)); 5669566063dSJacob Faibussowitsch PetscCall(VecCopy(obj->ksp_sol,globalout)); 567c4762a1bSJed Brown } 568c4762a1bSJed Brown return 0; 569c4762a1bSJed Brown } 570c4762a1bSJed Brown 571c4762a1bSJed Brown /*TEST 572c4762a1bSJed Brown 573c4762a1bSJed Brown build: 574c4762a1bSJed Brown requires: !complex 575c4762a1bSJed Brown 576c4762a1bSJed Brown test: 577c4762a1bSJed Brown suffix: euler 578c4762a1bSJed Brown output_file: output/ex3.out 579c4762a1bSJed Brown 580c4762a1bSJed Brown test: 581c4762a1bSJed Brown suffix: 2 582c4762a1bSJed Brown args: -useAlhs 583c4762a1bSJed Brown output_file: output/ex3.out 584c4762a1bSJed Brown TODO: Broken because SNESComputeJacobianDefault is incompatible with TSComputeIJacobianConstant 585c4762a1bSJed Brown 586c4762a1bSJed Brown TEST*/ 587