1c4762a1bSJed Brown static char help[] = "Demonstrates adjoint sensitivity analysis for Reaction-Diffusion Equations.\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /* 4c4762a1bSJed Brown REQUIRES configuration of PETSc with option --download-adolc. 5c4762a1bSJed Brown 6c4762a1bSJed Brown For documentation on ADOL-C, see 7c4762a1bSJed Brown $PETSC_ARCH/externalpackages/ADOL-C-2.6.0/ADOL-C/doc/adolc-manual.pdf 8c4762a1bSJed Brown 9c4762a1bSJed Brown REQUIRES configuration of PETSc with option --download-colpack 10c4762a1bSJed Brown 11c4762a1bSJed Brown For documentation on ColPack, see 12c4762a1bSJed Brown $PETSC_ARCH/externalpackages/git.colpack/README.md 13c4762a1bSJed Brown */ 14c4762a1bSJed Brown /* ------------------------------------------------------------------------ 15c4762a1bSJed Brown See ../advection-diffusion-reaction/ex5 for a description of the problem 16c4762a1bSJed Brown ------------------------------------------------------------------------- */ 17c4762a1bSJed Brown 18c4762a1bSJed Brown /* 19c4762a1bSJed Brown Runtime options: 20c4762a1bSJed Brown 21c4762a1bSJed Brown Solver: 22c4762a1bSJed Brown -forwardonly - Run the forward simulation without adjoint. 23c4762a1bSJed Brown -implicitform - Provide IFunction and IJacobian to TS, if not set, RHSFunction and RHSJacobian will be used. 24c4762a1bSJed Brown -aijpc - Set the preconditioner matrix to be aij (the Jacobian matrix can be of a different type such as ELL). 25c4762a1bSJed Brown 26c4762a1bSJed Brown Jacobian generation: 27c4762a1bSJed Brown -jacobian_by_hand - Use the hand-coded Jacobian of ex5.c, rather than generating it automatically. 28c4762a1bSJed Brown -no_annotation - Do not annotate ADOL-C active variables. (Should be used alongside -jacobian_by_hand.) 29c4762a1bSJed Brown -adolc_sparse - Calculate Jacobian in compressed format, using ADOL-C sparse functionality. 30c4762a1bSJed Brown -adolc_sparse_view - Print sparsity pattern used by -adolc_sparse option. 31c4762a1bSJed Brown */ 32c4762a1bSJed Brown /* 33c4762a1bSJed Brown NOTE: If -adolc_sparse option is used, at least four processors should be used, so that no processor owns boundaries which are 34c4762a1bSJed Brown identified by the periodic boundary conditions. The number of grid points in both x- and y-directions should be multiples 35c4762a1bSJed Brown of 5, in order for the 5-point stencil to be cleanly parallelised. 36c4762a1bSJed Brown */ 37c4762a1bSJed Brown 38c4762a1bSJed Brown #include <petscdmda.h> 39c4762a1bSJed Brown #include <petscts.h> 40c4762a1bSJed Brown #include "adolc-utils/drivers.cxx" 41c4762a1bSJed Brown #include <adolc/adolc.h> 42c4762a1bSJed Brown 43c4762a1bSJed Brown /* (Passive) field for the two variables */ 44c4762a1bSJed Brown typedef struct { 45c4762a1bSJed Brown PetscScalar u, v; 46c4762a1bSJed Brown } Field; 47c4762a1bSJed Brown 48c4762a1bSJed Brown /* Active field for the two variables */ 49c4762a1bSJed Brown typedef struct { 50c4762a1bSJed Brown adouble u, v; 51c4762a1bSJed Brown } AField; 52c4762a1bSJed Brown 53c4762a1bSJed Brown /* Application context */ 54c4762a1bSJed Brown typedef struct { 55c4762a1bSJed Brown PetscReal D1, D2, gamma, kappa; 56c4762a1bSJed Brown AField **u_a, **f_a; 57c4762a1bSJed Brown PetscBool aijpc; 58c4762a1bSJed Brown AdolcCtx *adctx; /* Automatic differentation support */ 59c4762a1bSJed Brown } AppCtx; 60c4762a1bSJed Brown 61c4762a1bSJed Brown extern PetscErrorCode InitialConditions(DM da, Vec U); 62c4762a1bSJed Brown extern PetscErrorCode InitializeLambda(DM da, Vec lambda, PetscReal x, PetscReal y); 63c4762a1bSJed Brown extern PetscErrorCode IFunctionLocalPassive(DMDALocalInfo *info, PetscReal t, Field **u, Field **udot, Field **f, void *ptr); 64c4762a1bSJed Brown extern PetscErrorCode IFunctionActive(TS ts, PetscReal ftime, Vec U, Vec Udot, Vec F, void *ptr); 65c4762a1bSJed Brown extern PetscErrorCode RHSFunctionActive(TS ts, PetscReal ftime, Vec U, Vec F, void *ptr); 66c4762a1bSJed Brown extern PetscErrorCode RHSFunctionPassive(TS ts, PetscReal ftime, Vec U, Vec F, void *ptr); 67c4762a1bSJed Brown extern PetscErrorCode IJacobianByHand(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, void *ctx); 68c4762a1bSJed Brown extern PetscErrorCode IJacobianAdolc(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, void *ctx); 69c4762a1bSJed Brown extern PetscErrorCode RHSJacobianByHand(TS ts, PetscReal t, Vec U, Mat A, Mat B, void *ctx); 70c4762a1bSJed Brown extern PetscErrorCode RHSJacobianAdolc(TS ts, PetscReal t, Vec U, Mat A, Mat B, void *ctx); 71c4762a1bSJed Brown 729371c9d4SSatish Balay int main(int argc, char **argv) { 73410585f6SHong Zhang TS ts; 74410585f6SHong Zhang Vec x, r, xdot; 75c4762a1bSJed Brown DM da; 76c4762a1bSJed Brown AppCtx appctx; 77c4762a1bSJed Brown AdolcCtx *adctx; 78c4762a1bSJed Brown Vec lambda[1]; 79c4762a1bSJed Brown PetscBool forwardonly = PETSC_FALSE, implicitform = PETSC_FALSE, byhand = PETSC_FALSE; 80c4762a1bSJed Brown PetscInt gxm, gym, i, dofs = 2, ctrl[3] = {0, 0, 0}; 81c4762a1bSJed Brown PetscScalar **Seed = NULL, **Rec = NULL, *u_vec; 82c4762a1bSJed Brown unsigned int **JP = NULL; 83c4762a1bSJed Brown ISColoring iscoloring; 84c4762a1bSJed Brown 85327415f7SBarry Smith PetscFunctionBeginUser; 869566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 879566063dSJacob Faibussowitsch PetscCall(PetscNew(&adctx)); 889566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-forwardonly", &forwardonly, NULL)); 899566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-implicitform", &implicitform, NULL)); 909371c9d4SSatish Balay appctx.aijpc = PETSC_FALSE, adctx->no_an = PETSC_FALSE, adctx->sparse = PETSC_FALSE, adctx->sparse_view = PETSC_FALSE; 919371c9d4SSatish Balay adctx->sparse_view_done = PETSC_FALSE; 929566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-aijpc", &appctx.aijpc, NULL)); 939566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-no_annotation", &adctx->no_an, NULL)); 949566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-jacobian_by_hand", &byhand, NULL)); 959566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-adolc_sparse", &adctx->sparse, NULL)); 969566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-adolc_sparse_view", &adctx->sparse_view, NULL)); 97c4762a1bSJed Brown appctx.D1 = 8.0e-5; 98c4762a1bSJed Brown appctx.D2 = 4.0e-5; 99c4762a1bSJed Brown appctx.gamma = .024; 100c4762a1bSJed Brown appctx.kappa = .06; 101c4762a1bSJed Brown appctx.adctx = adctx; 102c4762a1bSJed Brown 103c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 104c4762a1bSJed Brown Create distributed array (DMDA) to manage parallel grid and vectors 105c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1069566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, DM_BOUNDARY_PERIODIC, DMDA_STENCIL_STAR, 65, 65, PETSC_DECIDE, PETSC_DECIDE, 2, 1, NULL, NULL, &da)); 1079566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 1089566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 1099566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 0, "u")); 1109566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, 1, "v")); 111c4762a1bSJed Brown 112c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 113c4762a1bSJed Brown Extract global vectors from DMDA; then duplicate for remaining 114c4762a1bSJed Brown vectors that are the same types 115c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1169566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &x)); 1179566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &r)); 1189566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &xdot)); 119c4762a1bSJed Brown 120c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 121c4762a1bSJed Brown Create timestepping solver context 122c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1239566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 1249566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSCN)); 1259566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, da)); 1269566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 127c4762a1bSJed Brown if (!implicitform) { 1289566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, RHSFunctionPassive, &appctx)); 129c4762a1bSJed Brown } else { 1309566063dSJacob Faibussowitsch PetscCall(DMDATSSetIFunctionLocal(da, INSERT_VALUES, (DMDATSIFunctionLocal)IFunctionLocalPassive, &appctx)); 131c4762a1bSJed Brown } 132c4762a1bSJed Brown 133c4762a1bSJed Brown if (!adctx->no_an) { 1349566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, NULL, NULL, NULL, &gxm, &gym, NULL)); 1359371c9d4SSatish Balay adctx->m = dofs * gxm * gym; 1369371c9d4SSatish Balay adctx->n = dofs * gxm * gym; /* Number of dependent and independent variables */ 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 139c4762a1bSJed Brown Trace function(s) just once 140c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 141c4762a1bSJed Brown if (!implicitform) { 1429566063dSJacob Faibussowitsch PetscCall(RHSFunctionActive(ts, 1.0, x, r, &appctx)); 143c4762a1bSJed Brown } else { 1449566063dSJacob Faibussowitsch PetscCall(IFunctionActive(ts, 1.0, x, xdot, r, &appctx)); 145c4762a1bSJed Brown } 146c4762a1bSJed Brown 147c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 148c4762a1bSJed Brown In the case where ADOL-C generates the Jacobian in compressed format, 149c4762a1bSJed Brown seed and recovery matrices are required. Since the sparsity structure 150c4762a1bSJed Brown of the Jacobian does not change over the course of the time 151c4762a1bSJed Brown integration, we can save computational effort by only generating 152c4762a1bSJed Brown these objects once. 153c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 154c4762a1bSJed Brown if (adctx->sparse) { 155c4762a1bSJed Brown /* 156c4762a1bSJed Brown Generate sparsity pattern 157c4762a1bSJed Brown 158c4762a1bSJed Brown One way to do this is to use the Jacobian pattern driver `jac_pat` 159c4762a1bSJed Brown provided by ColPack. Otherwise, it is the responsibility of the user 160c4762a1bSJed Brown to obtain the sparsity pattern. 161c4762a1bSJed Brown */ 1629566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(adctx->n, &u_vec)); 163c4762a1bSJed Brown JP = (unsigned int **)malloc(adctx->m * sizeof(unsigned int *)); 164c4762a1bSJed Brown jac_pat(1, adctx->m, adctx->n, u_vec, JP, ctrl); 165*48a46eb9SPierre Jolivet if (adctx->sparse_view) PetscCall(PrintSparsity(MPI_COMM_WORLD, adctx->m, JP)); 166c4762a1bSJed Brown 167c4762a1bSJed Brown /* 168c4762a1bSJed Brown Extract a column colouring 169c4762a1bSJed Brown 170c4762a1bSJed Brown For problems using DMDA, colourings can be extracted directly, as 171c4762a1bSJed Brown follows. There are also ColPack drivers available. Otherwise, it is the 172c4762a1bSJed Brown responsibility of the user to obtain a suitable colouring. 173c4762a1bSJed Brown */ 1749566063dSJacob Faibussowitsch PetscCall(DMCreateColoring(da, IS_COLORING_LOCAL, &iscoloring)); 1759566063dSJacob Faibussowitsch PetscCall(ISColoringGetIS(iscoloring, PETSC_USE_POINTER, &adctx->p, NULL)); 176c4762a1bSJed Brown 177c4762a1bSJed Brown /* Generate seed matrix to propagate through the forward mode of AD */ 1789566063dSJacob Faibussowitsch PetscCall(AdolcMalloc2(adctx->n, adctx->p, &Seed)); 1799566063dSJacob Faibussowitsch PetscCall(GenerateSeedMatrix(iscoloring, Seed)); 1809566063dSJacob Faibussowitsch PetscCall(ISColoringDestroy(&iscoloring)); 181c4762a1bSJed Brown 182c4762a1bSJed Brown /* 183c4762a1bSJed Brown Generate recovery matrix, which is used to recover the Jacobian from 184c4762a1bSJed Brown compressed format */ 1859566063dSJacob Faibussowitsch PetscCall(AdolcMalloc2(adctx->m, adctx->p, &Rec)); 1869566063dSJacob Faibussowitsch PetscCall(GetRecoveryMatrix(Seed, JP, adctx->m, adctx->p, Rec)); 187c4762a1bSJed Brown 188c4762a1bSJed Brown /* Store results and free workspace */ 189c4762a1bSJed Brown adctx->Rec = Rec; 1909371c9d4SSatish Balay for (i = 0; i < adctx->m; i++) free(JP[i]); 191c4762a1bSJed Brown free(JP); 1929566063dSJacob Faibussowitsch PetscCall(PetscFree(u_vec)); 193c4762a1bSJed Brown 194c4762a1bSJed Brown } else { 195c4762a1bSJed Brown /* 196c4762a1bSJed Brown In 'full' Jacobian mode, propagate an identity matrix through the 197c4762a1bSJed Brown forward mode of AD. 198c4762a1bSJed Brown */ 199c4762a1bSJed Brown adctx->p = adctx->n; 2009566063dSJacob Faibussowitsch PetscCall(AdolcMalloc2(adctx->n, adctx->p, &Seed)); 2019566063dSJacob Faibussowitsch PetscCall(Identity(adctx->n, Seed)); 202c4762a1bSJed Brown } 203c4762a1bSJed Brown adctx->Seed = Seed; 204c4762a1bSJed Brown } 205c4762a1bSJed Brown 206c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 207c4762a1bSJed Brown Set Jacobian 208c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 209c4762a1bSJed Brown if (!implicitform) { 210c4762a1bSJed Brown if (!byhand) { 2119566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts, NULL, NULL, RHSJacobianAdolc, &appctx)); 212c4762a1bSJed Brown } else { 2139566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts, NULL, NULL, RHSJacobianByHand, &appctx)); 214c4762a1bSJed Brown } 215c4762a1bSJed Brown } else { 216c4762a1bSJed Brown if (appctx.aijpc) { 217c4762a1bSJed Brown Mat A, B; 218c4762a1bSJed Brown 2199566063dSJacob Faibussowitsch PetscCall(DMSetMatType(da, MATSELL)); 2209566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(da, &A)); 2219566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATAIJ, MAT_INITIAL_MATRIX, &B)); 222c4762a1bSJed Brown /* FIXME do we need to change viewer to display matrix in natural ordering as DMCreateMatrix_DA does? */ 223c4762a1bSJed Brown if (!byhand) { 2249566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, A, B, IJacobianAdolc, &appctx)); 225c4762a1bSJed Brown } else { 2269566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, A, B, IJacobianByHand, &appctx)); 227c4762a1bSJed Brown } 2289566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 2299566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 230c4762a1bSJed Brown } else { 231c4762a1bSJed Brown if (!byhand) { 2329566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, NULL, NULL, IJacobianAdolc, &appctx)); 233c4762a1bSJed Brown } else { 2349566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, NULL, NULL, IJacobianByHand, &appctx)); 235c4762a1bSJed Brown } 236c4762a1bSJed Brown } 237c4762a1bSJed Brown } 238c4762a1bSJed Brown 239c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 240c4762a1bSJed Brown Set initial conditions 241c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2429566063dSJacob Faibussowitsch PetscCall(InitialConditions(da, x)); 2439566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, x)); 244c4762a1bSJed Brown 245c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 246c4762a1bSJed Brown Have the TS save its trajectory so that TSAdjointSolve() may be used 247c4762a1bSJed Brown and set solver options 248c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 249c4762a1bSJed Brown if (!forwardonly) { 2509566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(ts)); 2519566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 200.0)); 2529566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 0.5)); 253c4762a1bSJed Brown } else { 2549566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 2000.0)); 2559566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 10)); 256c4762a1bSJed Brown } 2579566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 2589566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 259c4762a1bSJed Brown 260c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 261c4762a1bSJed Brown Solve ODE system 262c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2639566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, x)); 264c4762a1bSJed Brown if (!forwardonly) { 265c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 266c4762a1bSJed Brown Start the Adjoint model 267c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2689566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &lambda[0])); 269c4762a1bSJed Brown /* Reset initial conditions for the adjoint integration */ 2709566063dSJacob Faibussowitsch PetscCall(InitializeLambda(da, lambda[0], 0.5, 0.5)); 2719566063dSJacob Faibussowitsch PetscCall(TSSetCostGradients(ts, 1, lambda, NULL)); 2729566063dSJacob Faibussowitsch PetscCall(TSAdjointSolve(ts)); 2739566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lambda[0])); 274c4762a1bSJed Brown } 275c4762a1bSJed Brown 276c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 277c4762a1bSJed Brown Free work space. 278c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xdot)); 2809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 2819566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 2829566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 283c4762a1bSJed Brown if (!adctx->no_an) { 2849566063dSJacob Faibussowitsch if (adctx->sparse) PetscCall(AdolcFree2(Rec)); 2859566063dSJacob Faibussowitsch PetscCall(AdolcFree2(Seed)); 286c4762a1bSJed Brown } 2879566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 2889566063dSJacob Faibussowitsch PetscCall(PetscFree(adctx)); 2899566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 290b122ec5aSJacob Faibussowitsch return 0; 291c4762a1bSJed Brown } 292c4762a1bSJed Brown 2939371c9d4SSatish Balay PetscErrorCode InitialConditions(DM da, Vec U) { 294c4762a1bSJed Brown PetscInt i, j, xs, ys, xm, ym, Mx, My; 295c4762a1bSJed Brown Field **u; 296c4762a1bSJed Brown PetscReal hx, hy, x, y; 297c4762a1bSJed Brown 298c4762a1bSJed Brown PetscFunctionBegin; 2999566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, &My, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE)); 300c4762a1bSJed Brown 301c4762a1bSJed Brown hx = 2.5 / (PetscReal)Mx; 302c4762a1bSJed Brown hy = 2.5 / (PetscReal)My; 303c4762a1bSJed Brown 304c4762a1bSJed Brown /* 305c4762a1bSJed Brown Get pointers to vector data 306c4762a1bSJed Brown */ 3079566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, U, &u)); 308c4762a1bSJed Brown 309c4762a1bSJed Brown /* 310c4762a1bSJed Brown Get local grid boundaries 311c4762a1bSJed Brown */ 3129566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 313c4762a1bSJed Brown 314c4762a1bSJed Brown /* 315c4762a1bSJed Brown Compute function over the locally owned part of the grid 316c4762a1bSJed Brown */ 317c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 318c4762a1bSJed Brown y = j * hy; 319c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 320c4762a1bSJed Brown x = i * hx; 3219371c9d4SSatish Balay if (PetscApproximateGTE(x, 1.0) && PetscApproximateLTE(x, 1.5) && PetscApproximateGTE(y, 1.0) && PetscApproximateLTE(y, 1.5)) 3229371c9d4SSatish Balay u[j][i].v = PetscPowReal(PetscSinReal(4.0 * PETSC_PI * x), 2.0) * PetscPowReal(PetscSinReal(4.0 * PETSC_PI * y), 2.0) / 4.0; 323c4762a1bSJed Brown else u[j][i].v = 0.0; 324c4762a1bSJed Brown 325c4762a1bSJed Brown u[j][i].u = 1.0 - 2.0 * u[j][i].v; 326c4762a1bSJed Brown } 327c4762a1bSJed Brown } 328c4762a1bSJed Brown 329c4762a1bSJed Brown /* 330c4762a1bSJed Brown Restore vectors 331c4762a1bSJed Brown */ 3329566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, U, &u)); 333c4762a1bSJed Brown PetscFunctionReturn(0); 334c4762a1bSJed Brown } 335c4762a1bSJed Brown 3369371c9d4SSatish Balay PetscErrorCode InitializeLambda(DM da, Vec lambda, PetscReal x, PetscReal y) { 337c4762a1bSJed Brown PetscInt i, j, Mx, My, xs, ys, xm, ym; 338c4762a1bSJed Brown Field **l; 339c4762a1bSJed Brown PetscFunctionBegin; 340c4762a1bSJed Brown 3419566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, &My, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE)); 342c4762a1bSJed Brown /* locate the global i index for x and j index for y */ 343c4762a1bSJed Brown i = (PetscInt)(x * (Mx - 1)); 344c4762a1bSJed Brown j = (PetscInt)(y * (My - 1)); 3459566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 346c4762a1bSJed Brown 347c4762a1bSJed Brown if (xs <= i && i < xs + xm && ys <= j && j < ys + ym) { 348c4762a1bSJed Brown /* the i,j vertex is on this process */ 3499566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, lambda, &l)); 350c4762a1bSJed Brown l[j][i].u = 1.0; 351c4762a1bSJed Brown l[j][i].v = 1.0; 3529566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, lambda, &l)); 353c4762a1bSJed Brown } 354c4762a1bSJed Brown PetscFunctionReturn(0); 355c4762a1bSJed Brown } 356c4762a1bSJed Brown 3579371c9d4SSatish Balay PetscErrorCode IFunctionLocalPassive(DMDALocalInfo *info, PetscReal t, Field **u, Field **udot, Field **f, void *ptr) { 358c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ptr; 359c4762a1bSJed Brown PetscInt i, j, xs, ys, xm, ym; 360c4762a1bSJed Brown PetscReal hx, hy, sx, sy; 361c4762a1bSJed Brown PetscScalar uc, uxx, uyy, vc, vxx, vyy; 362c4762a1bSJed Brown 363c4762a1bSJed Brown PetscFunctionBegin; 3649371c9d4SSatish Balay hx = 2.50 / (PetscReal)(info->mx); 3659371c9d4SSatish Balay sx = 1.0 / (hx * hx); 3669371c9d4SSatish Balay hy = 2.50 / (PetscReal)(info->my); 3679371c9d4SSatish Balay sy = 1.0 / (hy * hy); 368c4762a1bSJed Brown 369c4762a1bSJed Brown /* Get local grid boundaries */ 3709371c9d4SSatish Balay xs = info->xs; 3719371c9d4SSatish Balay xm = info->xm; 3729371c9d4SSatish Balay ys = info->ys; 3739371c9d4SSatish Balay ym = info->ym; 374c4762a1bSJed Brown 375c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 376c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 377c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 378c4762a1bSJed Brown uc = u[j][i].u; 379c4762a1bSJed Brown uxx = (-2.0 * uc + u[j][i - 1].u + u[j][i + 1].u) * sx; 380c4762a1bSJed Brown uyy = (-2.0 * uc + u[j - 1][i].u + u[j + 1][i].u) * sy; 381c4762a1bSJed Brown vc = u[j][i].v; 382c4762a1bSJed Brown vxx = (-2.0 * vc + u[j][i - 1].v + u[j][i + 1].v) * sx; 383c4762a1bSJed Brown vyy = (-2.0 * vc + u[j - 1][i].v + u[j + 1][i].v) * sy; 384c4762a1bSJed Brown f[j][i].u = udot[j][i].u - appctx->D1 * (uxx + uyy) + uc * vc * vc - appctx->gamma * (1.0 - uc); 385c4762a1bSJed Brown f[j][i].v = udot[j][i].v - appctx->D2 * (vxx + vyy) - uc * vc * vc + (appctx->gamma + appctx->kappa) * vc; 386c4762a1bSJed Brown } 387c4762a1bSJed Brown } 3889566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(16.0 * xm * ym)); 389c4762a1bSJed Brown PetscFunctionReturn(0); 390c4762a1bSJed Brown } 391c4762a1bSJed Brown 3929371c9d4SSatish Balay PetscErrorCode IFunctionActive(TS ts, PetscReal ftime, Vec U, Vec Udot, Vec F, void *ptr) { 393c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ptr; 394c4762a1bSJed Brown DM da; 395c4762a1bSJed Brown DMDALocalInfo info; 396c4762a1bSJed Brown Field **u, **f, **udot; 397c4762a1bSJed Brown Vec localU; 398c4762a1bSJed Brown PetscInt i, j, xs, ys, xm, ym, gxs, gys, gxm, gym; 399c4762a1bSJed Brown PetscReal hx, hy, sx, sy; 400c4762a1bSJed Brown adouble uc, uxx, uyy, vc, vxx, vyy; 401c4762a1bSJed Brown AField **f_a, *f_c, **u_a, *u_c; 402c4762a1bSJed Brown PetscScalar dummy; 403c4762a1bSJed Brown 404c4762a1bSJed Brown PetscFunctionBegin; 405c4762a1bSJed Brown 4069566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 4079566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 4089566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localU)); 4099371c9d4SSatish Balay hx = 2.50 / (PetscReal)(info.mx); 4109371c9d4SSatish Balay sx = 1.0 / (hx * hx); 4119371c9d4SSatish Balay hy = 2.50 / (PetscReal)(info.my); 4129371c9d4SSatish Balay sy = 1.0 / (hy * hy); 4139371c9d4SSatish Balay xs = info.xs; 4149371c9d4SSatish Balay xm = info.xm; 4159371c9d4SSatish Balay gxs = info.gxs; 4169371c9d4SSatish Balay gxm = info.gxm; 4179371c9d4SSatish Balay ys = info.ys; 4189371c9d4SSatish Balay ym = info.ym; 4199371c9d4SSatish Balay gys = info.gys; 4209371c9d4SSatish Balay gym = info.gym; 421c4762a1bSJed Brown 422c4762a1bSJed Brown /* 423c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 424c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 425c4762a1bSJed Brown By placing code between these two statements, computations can be 426c4762a1bSJed Brown done while messages are in transition. 427c4762a1bSJed Brown */ 4289566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, localU)); 4299566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, localU)); 430c4762a1bSJed Brown 431c4762a1bSJed Brown /* 432c4762a1bSJed Brown Get pointers to vector data 433c4762a1bSJed Brown */ 4349566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localU, &u)); 4359566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, F, &f)); 4369566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, Udot, &udot)); 437c4762a1bSJed Brown 438c4762a1bSJed Brown /* 439c4762a1bSJed Brown Create contiguous 1-arrays of AFields 440c4762a1bSJed Brown 441c4762a1bSJed Brown NOTE: Memory for ADOL-C active variables (such as adouble and AField) 442c4762a1bSJed Brown cannot be allocated using PetscMalloc, as this does not call the 443c4762a1bSJed Brown relevant class constructor. Instead, we use the C++ keyword `new`. 444c4762a1bSJed Brown */ 445c4762a1bSJed Brown u_c = new AField[info.gxm * info.gym]; 446c4762a1bSJed Brown f_c = new AField[info.gxm * info.gym]; 447c4762a1bSJed Brown 448c4762a1bSJed Brown /* Create corresponding 2-arrays of AFields */ 449c4762a1bSJed Brown u_a = new AField *[info.gym]; 450c4762a1bSJed Brown f_a = new AField *[info.gym]; 451c4762a1bSJed Brown 452c4762a1bSJed Brown /* Align indices between array types to endow 2d array with ghost points */ 4539566063dSJacob Faibussowitsch PetscCall(GiveGhostPoints(da, u_c, &u_a)); 4549566063dSJacob Faibussowitsch PetscCall(GiveGhostPoints(da, f_c, &f_a)); 455c4762a1bSJed Brown 456c4762a1bSJed Brown trace_on(1); /* Start of active section on tape 1 */ 457c4762a1bSJed Brown 458c4762a1bSJed Brown /* 459c4762a1bSJed Brown Mark independence 460c4762a1bSJed Brown 461c4762a1bSJed Brown NOTE: Ghost points are marked as independent, in place of the points they represent on 462c4762a1bSJed Brown other processors / on other boundaries. 463c4762a1bSJed Brown */ 464c4762a1bSJed Brown for (j = gys; j < gys + gym; j++) { 465c4762a1bSJed Brown for (i = gxs; i < gxs + gxm; i++) { 466c4762a1bSJed Brown u_a[j][i].u <<= u[j][i].u; 467c4762a1bSJed Brown u_a[j][i].v <<= u[j][i].v; 468c4762a1bSJed Brown } 469c4762a1bSJed Brown } 470c4762a1bSJed Brown 471c4762a1bSJed Brown /* Compute function over the locally owned part of the grid */ 472c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 473c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 474c4762a1bSJed Brown uc = u_a[j][i].u; 475c4762a1bSJed Brown uxx = (-2.0 * uc + u_a[j][i - 1].u + u_a[j][i + 1].u) * sx; 476c4762a1bSJed Brown uyy = (-2.0 * uc + u_a[j - 1][i].u + u_a[j + 1][i].u) * sy; 477c4762a1bSJed Brown vc = u_a[j][i].v; 478c4762a1bSJed Brown vxx = (-2.0 * vc + u_a[j][i - 1].v + u_a[j][i + 1].v) * sx; 479c4762a1bSJed Brown vyy = (-2.0 * vc + u_a[j - 1][i].v + u_a[j + 1][i].v) * sy; 480c4762a1bSJed Brown f_a[j][i].u = udot[j][i].u - appctx->D1 * (uxx + uyy) + uc * vc * vc - appctx->gamma * (1.0 - uc); 481c4762a1bSJed Brown f_a[j][i].v = udot[j][i].v - appctx->D2 * (vxx + vyy) - uc * vc * vc + (appctx->gamma + appctx->kappa) * vc; 482c4762a1bSJed Brown } 483c4762a1bSJed Brown } 484c4762a1bSJed Brown 485c4762a1bSJed Brown /* 486c4762a1bSJed Brown Mark dependence 487c4762a1bSJed Brown 488c4762a1bSJed Brown NOTE: Marking dependence of dummy variables makes the index notation much simpler when forming 489c4762a1bSJed Brown the Jacobian later. 490c4762a1bSJed Brown */ 491c4762a1bSJed Brown for (j = gys; j < gys + gym; j++) { 492c4762a1bSJed Brown for (i = gxs; i < gxs + gxm; i++) { 493c4762a1bSJed Brown if ((i < xs) || (i >= xs + xm) || (j < ys) || (j >= ys + ym)) { 494c4762a1bSJed Brown f_a[j][i].u >>= dummy; 495c4762a1bSJed Brown f_a[j][i].v >>= dummy; 496c4762a1bSJed Brown } else { 497c4762a1bSJed Brown f_a[j][i].u >>= f[j][i].u; 498c4762a1bSJed Brown f_a[j][i].v >>= f[j][i].v; 499c4762a1bSJed Brown } 500c4762a1bSJed Brown } 501c4762a1bSJed Brown } 502c4762a1bSJed Brown trace_off(); /* End of active section */ 5039566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(16.0 * xm * ym)); 504c4762a1bSJed Brown 505c4762a1bSJed Brown /* Restore vectors */ 5069566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, F, &f)); 5079566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localU, &u)); 5089566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, Udot, &udot)); 509c4762a1bSJed Brown 5109566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localU)); 511410585f6SHong Zhang 512c4762a1bSJed Brown /* Destroy AFields appropriately */ 513c4762a1bSJed Brown f_a += info.gys; 514c4762a1bSJed Brown u_a += info.gys; 515c4762a1bSJed Brown delete[] f_a; 516c4762a1bSJed Brown delete[] u_a; 517c4762a1bSJed Brown delete[] f_c; 518c4762a1bSJed Brown delete[] u_c; 519c4762a1bSJed Brown 520c4762a1bSJed Brown PetscFunctionReturn(0); 521c4762a1bSJed Brown } 522c4762a1bSJed Brown 5239371c9d4SSatish Balay PetscErrorCode RHSFunctionPassive(TS ts, PetscReal ftime, Vec U, Vec F, void *ptr) { 524c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ptr; 525c4762a1bSJed Brown DM da; 526c4762a1bSJed Brown PetscInt i, j, xs, ys, xm, ym, Mx, My; 527c4762a1bSJed Brown PetscReal hx, hy, sx, sy; 528c4762a1bSJed Brown PetscScalar uc, uxx, uyy, vc, vxx, vyy; 529c4762a1bSJed Brown Field **u, **f; 530c4762a1bSJed Brown Vec localU, localF; 531c4762a1bSJed Brown 532c4762a1bSJed Brown PetscFunctionBegin; 5339566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 5349566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, &My, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE)); 5359371c9d4SSatish Balay hx = 2.50 / (PetscReal)(Mx); 5369371c9d4SSatish Balay sx = 1.0 / (hx * hx); 5379371c9d4SSatish Balay hy = 2.50 / (PetscReal)(My); 5389371c9d4SSatish Balay sy = 1.0 / (hy * hy); 5399566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localU)); 5409566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localF)); 541c4762a1bSJed Brown 542c4762a1bSJed Brown /* 543c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 544c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 545c4762a1bSJed Brown By placing code between these two statements, computations can be 546c4762a1bSJed Brown done while messages are in transition. 547c4762a1bSJed Brown */ 5489566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, localU)); 5499566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, localU)); 5509566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); // NOTE (1): See (2) below 5519566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, F, INSERT_VALUES, localF)); 5529566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, F, INSERT_VALUES, localF)); 553c4762a1bSJed Brown 554c4762a1bSJed Brown /* 555c4762a1bSJed Brown Get pointers to vector data 556c4762a1bSJed Brown */ 5579566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localU, &u)); 5589566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, localF, &f)); 559c4762a1bSJed Brown 560c4762a1bSJed Brown /* 561c4762a1bSJed Brown Get local grid boundaries 562c4762a1bSJed Brown */ 5639566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 564c4762a1bSJed Brown 565c4762a1bSJed Brown /* 566c4762a1bSJed Brown Compute function over the locally owned part of the grid 567c4762a1bSJed Brown */ 568c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 569c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 570c4762a1bSJed Brown uc = u[j][i].u; 571c4762a1bSJed Brown uxx = (-2.0 * uc + u[j][i - 1].u + u[j][i + 1].u) * sx; 572c4762a1bSJed Brown uyy = (-2.0 * uc + u[j - 1][i].u + u[j + 1][i].u) * sy; 573c4762a1bSJed Brown vc = u[j][i].v; 574c4762a1bSJed Brown vxx = (-2.0 * vc + u[j][i - 1].v + u[j][i + 1].v) * sx; 575c4762a1bSJed Brown vyy = (-2.0 * vc + u[j - 1][i].v + u[j + 1][i].v) * sy; 576c4762a1bSJed Brown f[j][i].u = appctx->D1 * (uxx + uyy) - uc * vc * vc + appctx->gamma * (1.0 - uc); 577c4762a1bSJed Brown f[j][i].v = appctx->D2 * (vxx + vyy) + uc * vc * vc - (appctx->gamma + appctx->kappa) * vc; 578c4762a1bSJed Brown } 579c4762a1bSJed Brown } 580c4762a1bSJed Brown 581c4762a1bSJed Brown /* 582c4762a1bSJed Brown Gather global vector, using the 2-step process 583c4762a1bSJed Brown DMLocalToGlobalBegin(),DMLocalToGlobalEnd(). 584c4762a1bSJed Brown 585c4762a1bSJed Brown NOTE (2): We need to use ADD_VALUES if boundaries are not of type DM_BOUNDARY_NONE or 586c4762a1bSJed Brown DM_BOUNDARY_GHOSTED, meaning we should also zero F before addition (see (1) above). 587c4762a1bSJed Brown */ 5889566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(da, localF, ADD_VALUES, F)); 5899566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(da, localF, ADD_VALUES, F)); 590c4762a1bSJed Brown 591c4762a1bSJed Brown /* 592c4762a1bSJed Brown Restore vectors 593c4762a1bSJed Brown */ 5949566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, localF, &f)); 5959566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localU, &u)); 5969566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localF)); 5979566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localU)); 5989566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(16.0 * xm * ym)); 599c4762a1bSJed Brown PetscFunctionReturn(0); 600c4762a1bSJed Brown } 601c4762a1bSJed Brown 6029371c9d4SSatish Balay PetscErrorCode RHSFunctionActive(TS ts, PetscReal ftime, Vec U, Vec F, void *ptr) { 603c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ptr; 604c4762a1bSJed Brown DM da; 605c4762a1bSJed Brown PetscInt i, j, xs, ys, xm, ym, gxs, gys, gxm, gym, Mx, My; 606c4762a1bSJed Brown PetscReal hx, hy, sx, sy; 607c4762a1bSJed Brown AField **f_a, *f_c, **u_a, *u_c; 608c4762a1bSJed Brown adouble uc, uxx, uyy, vc, vxx, vyy; 609c4762a1bSJed Brown Field **u, **f; 610c4762a1bSJed Brown Vec localU, localF; 611c4762a1bSJed Brown 612c4762a1bSJed Brown PetscFunctionBegin; 6139566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 6149566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, &My, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE)); 6159371c9d4SSatish Balay hx = 2.50 / (PetscReal)(Mx); 6169371c9d4SSatish Balay sx = 1.0 / (hx * hx); 6179371c9d4SSatish Balay hy = 2.50 / (PetscReal)(My); 6189371c9d4SSatish Balay sy = 1.0 / (hy * hy); 6199566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localU)); 6209566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localF)); 621c4762a1bSJed Brown 622c4762a1bSJed Brown /* 623c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 624c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 625c4762a1bSJed Brown By placing code between these two statements, computations can be 626c4762a1bSJed Brown done while messages are in transition. 627c4762a1bSJed Brown */ 6289566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, localU)); 6299566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, localU)); 6309566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); // NOTE (1): See (2) below 6319566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, F, INSERT_VALUES, localF)); 6329566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, F, INSERT_VALUES, localF)); 633c4762a1bSJed Brown 634c4762a1bSJed Brown /* 635c4762a1bSJed Brown Get pointers to vector data 636c4762a1bSJed Brown */ 6379566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localU, &u)); 6389566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, localF, &f)); 639c4762a1bSJed Brown 640c4762a1bSJed Brown /* 641c4762a1bSJed Brown Get local and ghosted grid boundaries 642c4762a1bSJed Brown */ 6439566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, NULL, &gxm, &gym, NULL)); 6449566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 645c4762a1bSJed Brown 646c4762a1bSJed Brown /* 647c4762a1bSJed Brown Create contiguous 1-arrays of AFields 648c4762a1bSJed Brown 649c4762a1bSJed Brown NOTE: Memory for ADOL-C active variables (such as adouble and AField) 650c4762a1bSJed Brown cannot be allocated using PetscMalloc, as this does not call the 651c4762a1bSJed Brown relevant class constructor. Instead, we use the C++ keyword `new`. 652c4762a1bSJed Brown */ 653c4762a1bSJed Brown u_c = new AField[gxm * gym]; 654c4762a1bSJed Brown f_c = new AField[gxm * gym]; 655c4762a1bSJed Brown 656c4762a1bSJed Brown /* Create corresponding 2-arrays of AFields */ 657c4762a1bSJed Brown u_a = new AField *[gym]; 658c4762a1bSJed Brown f_a = new AField *[gym]; 659c4762a1bSJed Brown 660c4762a1bSJed Brown /* Align indices between array types to endow 2d array with ghost points */ 6619566063dSJacob Faibussowitsch PetscCall(GiveGhostPoints(da, u_c, &u_a)); 6629566063dSJacob Faibussowitsch PetscCall(GiveGhostPoints(da, f_c, &f_a)); 663c4762a1bSJed Brown 664c4762a1bSJed Brown /* 665c4762a1bSJed Brown Compute function over the locally owned part of the grid 666c4762a1bSJed Brown */ 667c4762a1bSJed Brown trace_on(1); // ----------------------------------------------- Start of active section 668c4762a1bSJed Brown 669c4762a1bSJed Brown /* 670c4762a1bSJed Brown Mark independence 671c4762a1bSJed Brown 672c4762a1bSJed Brown NOTE: Ghost points are marked as independent, in place of the points they represent on 673c4762a1bSJed Brown other processors / on other boundaries. 674c4762a1bSJed Brown */ 675c4762a1bSJed Brown for (j = gys; j < gys + gym; j++) { 676c4762a1bSJed Brown for (i = gxs; i < gxs + gxm; i++) { 677c4762a1bSJed Brown u_a[j][i].u <<= u[j][i].u; 678c4762a1bSJed Brown u_a[j][i].v <<= u[j][i].v; 679c4762a1bSJed Brown } 680c4762a1bSJed Brown } 681c4762a1bSJed Brown 682c4762a1bSJed Brown /* 683c4762a1bSJed Brown Compute function over the locally owned part of the grid 684c4762a1bSJed Brown */ 685c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 686c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 687c4762a1bSJed Brown uc = u_a[j][i].u; 688c4762a1bSJed Brown uxx = (-2.0 * uc + u_a[j][i - 1].u + u_a[j][i + 1].u) * sx; 689c4762a1bSJed Brown uyy = (-2.0 * uc + u_a[j - 1][i].u + u_a[j + 1][i].u) * sy; 690c4762a1bSJed Brown vc = u_a[j][i].v; 691c4762a1bSJed Brown vxx = (-2.0 * vc + u_a[j][i - 1].v + u_a[j][i + 1].v) * sx; 692c4762a1bSJed Brown vyy = (-2.0 * vc + u_a[j - 1][i].v + u_a[j + 1][i].v) * sy; 693c4762a1bSJed Brown f_a[j][i].u = appctx->D1 * (uxx + uyy) - uc * vc * vc + appctx->gamma * (1.0 - uc); 694c4762a1bSJed Brown f_a[j][i].v = appctx->D2 * (vxx + vyy) + uc * vc * vc - (appctx->gamma + appctx->kappa) * vc; 695c4762a1bSJed Brown } 696c4762a1bSJed Brown } 697c4762a1bSJed Brown /* 698c4762a1bSJed Brown Mark dependence 699c4762a1bSJed Brown 700c4762a1bSJed Brown NOTE: Ghost points are marked as dependent in order to vastly simplify index notation 701c4762a1bSJed Brown during Jacobian assembly. 702c4762a1bSJed Brown */ 703c4762a1bSJed Brown for (j = gys; j < gys + gym; j++) { 704c4762a1bSJed Brown for (i = gxs; i < gxs + gxm; i++) { 705c4762a1bSJed Brown f_a[j][i].u >>= f[j][i].u; 706c4762a1bSJed Brown f_a[j][i].v >>= f[j][i].v; 707c4762a1bSJed Brown } 708c4762a1bSJed Brown } 709c4762a1bSJed Brown trace_off(); // ----------------------------------------------- End of active section 710c4762a1bSJed Brown 711c4762a1bSJed Brown /* Test zeroth order scalar evaluation in ADOL-C gives the same result */ 712c4762a1bSJed Brown // if (appctx->adctx->zos) { 7139566063dSJacob Faibussowitsch // PetscCall(TestZOS2d(da,f,u,appctx)); // FIXME: Why does this give nonzero? 714c4762a1bSJed Brown // } 715c4762a1bSJed Brown 716c4762a1bSJed Brown /* 717c4762a1bSJed Brown Gather global vector, using the 2-step process 718c4762a1bSJed Brown DMLocalToGlobalBegin(),DMLocalToGlobalEnd(). 719c4762a1bSJed Brown 720c4762a1bSJed Brown NOTE (2): We need to use ADD_VALUES if boundaries are not of type DM_BOUNDARY_NONE or 721c4762a1bSJed Brown DM_BOUNDARY_GHOSTED, meaning we should also zero F before addition (see (1) above). 722c4762a1bSJed Brown */ 7239566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(da, localF, ADD_VALUES, F)); 7249566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(da, localF, ADD_VALUES, F)); 725c4762a1bSJed Brown 726c4762a1bSJed Brown /* 727c4762a1bSJed Brown Restore vectors 728c4762a1bSJed Brown */ 7299566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, localF, &f)); 7309566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localU, &u)); 7319566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localF)); 7329566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localU)); 7339566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(16.0 * xm * ym)); 734c4762a1bSJed Brown 735c4762a1bSJed Brown /* Destroy AFields appropriately */ 736c4762a1bSJed Brown f_a += gys; 737c4762a1bSJed Brown u_a += gys; 738c4762a1bSJed Brown delete[] f_a; 739c4762a1bSJed Brown delete[] u_a; 740c4762a1bSJed Brown delete[] f_c; 741c4762a1bSJed Brown delete[] u_c; 742c4762a1bSJed Brown 743c4762a1bSJed Brown PetscFunctionReturn(0); 744c4762a1bSJed Brown } 745c4762a1bSJed Brown 7469371c9d4SSatish Balay PetscErrorCode IJacobianAdolc(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, void *ctx) { 747c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx; 748c4762a1bSJed Brown DM da; 749a8c08197SHong Zhang const PetscScalar *u_vec; 750c4762a1bSJed Brown Vec localU; 751c4762a1bSJed Brown 752c4762a1bSJed Brown PetscFunctionBegin; 7539566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 7549566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localU)); 755c4762a1bSJed Brown 756c4762a1bSJed Brown /* 757c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 758c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 759c4762a1bSJed Brown By placing code between these two statements, computations can be 760c4762a1bSJed Brown done while messages are in transition. 761c4762a1bSJed Brown */ 7629566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, localU)); 7639566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, localU)); 764c4762a1bSJed Brown 765c4762a1bSJed Brown /* Get pointers to vector data */ 7669566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(localU, &u_vec)); 767c4762a1bSJed Brown 768c4762a1bSJed Brown /* 769c4762a1bSJed Brown Compute Jacobian 770c4762a1bSJed Brown */ 7719566063dSJacob Faibussowitsch PetscCall(PetscAdolcComputeIJacobianLocalIDMass(1, A, u_vec, a, appctx->adctx)); 772c4762a1bSJed Brown 773c4762a1bSJed Brown /* 774c4762a1bSJed Brown Restore vectors 775c4762a1bSJed Brown */ 7769566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(localU, &u_vec)); 7779566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localU)); 778c4762a1bSJed Brown PetscFunctionReturn(0); 779c4762a1bSJed Brown } 780c4762a1bSJed Brown 7819371c9d4SSatish Balay PetscErrorCode IJacobianByHand(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, void *ctx) { 782c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx; /* user-defined application context */ 783c4762a1bSJed Brown DM da; 784c4762a1bSJed Brown PetscInt i, j, Mx, My, xs, ys, xm, ym; 785c4762a1bSJed Brown PetscReal hx, hy, sx, sy; 786c4762a1bSJed Brown PetscScalar uc, vc; 787c4762a1bSJed Brown Field **u; 788c4762a1bSJed Brown Vec localU; 789c4762a1bSJed Brown MatStencil stencil[6], rowstencil; 790c4762a1bSJed Brown PetscScalar entries[6]; 791c4762a1bSJed Brown 792c4762a1bSJed Brown PetscFunctionBegin; 7939566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 7949566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localU)); 7959566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, &My, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE)); 796c4762a1bSJed Brown 7979371c9d4SSatish Balay hx = 2.50 / (PetscReal)Mx; 7989371c9d4SSatish Balay sx = 1.0 / (hx * hx); 7999371c9d4SSatish Balay hy = 2.50 / (PetscReal)My; 8009371c9d4SSatish Balay sy = 1.0 / (hy * hy); 801c4762a1bSJed Brown 802c4762a1bSJed Brown /* 803c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 804c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 805c4762a1bSJed Brown By placing code between these two statements, computations can be 806c4762a1bSJed Brown done while messages are in transition. 807c4762a1bSJed Brown */ 8089566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, localU)); 8099566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, localU)); 810c4762a1bSJed Brown 811c4762a1bSJed Brown /* 812c4762a1bSJed Brown Get pointers to vector data 813c4762a1bSJed Brown */ 8149566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localU, &u)); 815c4762a1bSJed Brown 816c4762a1bSJed Brown /* 817c4762a1bSJed Brown Get local grid boundaries 818c4762a1bSJed Brown */ 8199566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 820c4762a1bSJed Brown 821c4762a1bSJed Brown stencil[0].k = 0; 822c4762a1bSJed Brown stencil[1].k = 0; 823c4762a1bSJed Brown stencil[2].k = 0; 824c4762a1bSJed Brown stencil[3].k = 0; 825c4762a1bSJed Brown stencil[4].k = 0; 826c4762a1bSJed Brown stencil[5].k = 0; 827c4762a1bSJed Brown rowstencil.k = 0; 828c4762a1bSJed Brown /* 829c4762a1bSJed Brown Compute function over the locally owned part of the grid 830c4762a1bSJed Brown */ 831c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 832c4762a1bSJed Brown stencil[0].j = j - 1; 833c4762a1bSJed Brown stencil[1].j = j + 1; 834c4762a1bSJed Brown stencil[2].j = j; 835c4762a1bSJed Brown stencil[3].j = j; 836c4762a1bSJed Brown stencil[4].j = j; 837c4762a1bSJed Brown stencil[5].j = j; 8389371c9d4SSatish Balay rowstencil.k = 0; 8399371c9d4SSatish Balay rowstencil.j = j; 840c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 841c4762a1bSJed Brown uc = u[j][i].u; 842c4762a1bSJed Brown vc = u[j][i].v; 843c4762a1bSJed Brown 844c4762a1bSJed Brown /* uxx = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx; 845c4762a1bSJed Brown uyy = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy; 846c4762a1bSJed Brown 847c4762a1bSJed Brown vxx = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx; 848c4762a1bSJed Brown vyy = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy; 849c4762a1bSJed Brown f[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc);*/ 850c4762a1bSJed Brown 8519371c9d4SSatish Balay stencil[0].i = i; 8529371c9d4SSatish Balay stencil[0].c = 0; 8539371c9d4SSatish Balay entries[0] = -appctx->D1 * sy; 8549371c9d4SSatish Balay stencil[1].i = i; 8559371c9d4SSatish Balay stencil[1].c = 0; 8569371c9d4SSatish Balay entries[1] = -appctx->D1 * sy; 8579371c9d4SSatish Balay stencil[2].i = i - 1; 8589371c9d4SSatish Balay stencil[2].c = 0; 8599371c9d4SSatish Balay entries[2] = -appctx->D1 * sx; 8609371c9d4SSatish Balay stencil[3].i = i + 1; 8619371c9d4SSatish Balay stencil[3].c = 0; 8629371c9d4SSatish Balay entries[3] = -appctx->D1 * sx; 8639371c9d4SSatish Balay stencil[4].i = i; 8649371c9d4SSatish Balay stencil[4].c = 0; 8659371c9d4SSatish Balay entries[4] = 2.0 * appctx->D1 * (sx + sy) + vc * vc + appctx->gamma + a; 8669371c9d4SSatish Balay stencil[5].i = i; 8679371c9d4SSatish Balay stencil[5].c = 1; 8689371c9d4SSatish Balay entries[5] = 2.0 * uc * vc; 8699371c9d4SSatish Balay rowstencil.i = i; 8709371c9d4SSatish Balay rowstencil.c = 0; 871c4762a1bSJed Brown 8729566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(A, 1, &rowstencil, 6, stencil, entries, INSERT_VALUES)); 873*48a46eb9SPierre Jolivet if (appctx->aijpc) PetscCall(MatSetValuesStencil(B, 1, &rowstencil, 6, stencil, entries, INSERT_VALUES)); 8749371c9d4SSatish Balay stencil[0].c = 1; 8759371c9d4SSatish Balay entries[0] = -appctx->D2 * sy; 8769371c9d4SSatish Balay stencil[1].c = 1; 8779371c9d4SSatish Balay entries[1] = -appctx->D2 * sy; 8789371c9d4SSatish Balay stencil[2].c = 1; 8799371c9d4SSatish Balay entries[2] = -appctx->D2 * sx; 8809371c9d4SSatish Balay stencil[3].c = 1; 8819371c9d4SSatish Balay entries[3] = -appctx->D2 * sx; 8829371c9d4SSatish Balay stencil[4].c = 1; 8839371c9d4SSatish Balay entries[4] = 2.0 * appctx->D2 * (sx + sy) - 2.0 * uc * vc + appctx->gamma + appctx->kappa + a; 8849371c9d4SSatish Balay stencil[5].c = 0; 8859371c9d4SSatish Balay entries[5] = -vc * vc; 886c4762a1bSJed Brown 8879566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(A, 1, &rowstencil, 6, stencil, entries, INSERT_VALUES)); 888*48a46eb9SPierre Jolivet if (appctx->aijpc) PetscCall(MatSetValuesStencil(B, 1, &rowstencil, 6, stencil, entries, INSERT_VALUES)); 889c4762a1bSJed Brown /* f[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc; */ 890c4762a1bSJed Brown } 891c4762a1bSJed Brown } 892c4762a1bSJed Brown 893c4762a1bSJed Brown /* 894c4762a1bSJed Brown Restore vectors 895c4762a1bSJed Brown */ 8969566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(19.0 * xm * ym)); 8979566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localU, &u)); 8989566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localU)); 8999566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 9009566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 9019566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE)); 902c4762a1bSJed Brown if (appctx->aijpc) { 9039566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 9049566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 9059566063dSJacob Faibussowitsch PetscCall(MatSetOption(B, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE)); 906c4762a1bSJed Brown } 907c4762a1bSJed Brown PetscFunctionReturn(0); 908c4762a1bSJed Brown } 909c4762a1bSJed Brown 9109371c9d4SSatish Balay PetscErrorCode RHSJacobianByHand(TS ts, PetscReal t, Vec U, Mat A, Mat B, void *ctx) { 911c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx; /* user-defined application context */ 912c4762a1bSJed Brown DM da; 913c4762a1bSJed Brown PetscInt i, j, Mx, My, xs, ys, xm, ym; 914c4762a1bSJed Brown PetscReal hx, hy, sx, sy; 915c4762a1bSJed Brown PetscScalar uc, vc; 916c4762a1bSJed Brown Field **u; 917c4762a1bSJed Brown Vec localU; 918c4762a1bSJed Brown MatStencil stencil[6], rowstencil; 919c4762a1bSJed Brown PetscScalar entries[6]; 920c4762a1bSJed Brown 921c4762a1bSJed Brown PetscFunctionBegin; 9229566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 9239566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localU)); 9249566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, &My, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE)); 925c4762a1bSJed Brown 9269371c9d4SSatish Balay hx = 2.50 / (PetscReal)(Mx); 9279371c9d4SSatish Balay sx = 1.0 / (hx * hx); 9289371c9d4SSatish Balay hy = 2.50 / (PetscReal)(My); 9299371c9d4SSatish Balay sy = 1.0 / (hy * hy); 930c4762a1bSJed Brown 931c4762a1bSJed Brown /* 932c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 933c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 934c4762a1bSJed Brown By placing code between these two statements, computations can be 935c4762a1bSJed Brown done while messages are in transition. 936c4762a1bSJed Brown */ 9379566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, localU)); 9389566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, localU)); 939c4762a1bSJed Brown 940c4762a1bSJed Brown /* 941c4762a1bSJed Brown Get pointers to vector data 942c4762a1bSJed Brown */ 9439566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayRead(da, localU, &u)); 944c4762a1bSJed Brown 945c4762a1bSJed Brown /* 946c4762a1bSJed Brown Get local grid boundaries 947c4762a1bSJed Brown */ 9489566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL)); 949c4762a1bSJed Brown 950c4762a1bSJed Brown stencil[0].k = 0; 951c4762a1bSJed Brown stencil[1].k = 0; 952c4762a1bSJed Brown stencil[2].k = 0; 953c4762a1bSJed Brown stencil[3].k = 0; 954c4762a1bSJed Brown stencil[4].k = 0; 955c4762a1bSJed Brown stencil[5].k = 0; 956c4762a1bSJed Brown rowstencil.k = 0; 957c4762a1bSJed Brown 958c4762a1bSJed Brown /* 959c4762a1bSJed Brown Compute function over the locally owned part of the grid 960c4762a1bSJed Brown */ 961c4762a1bSJed Brown for (j = ys; j < ys + ym; j++) { 962c4762a1bSJed Brown stencil[0].j = j - 1; 963c4762a1bSJed Brown stencil[1].j = j + 1; 964c4762a1bSJed Brown stencil[2].j = j; 965c4762a1bSJed Brown stencil[3].j = j; 966c4762a1bSJed Brown stencil[4].j = j; 967c4762a1bSJed Brown stencil[5].j = j; 9689371c9d4SSatish Balay rowstencil.k = 0; 9699371c9d4SSatish Balay rowstencil.j = j; 970c4762a1bSJed Brown for (i = xs; i < xs + xm; i++) { 971c4762a1bSJed Brown uc = u[j][i].u; 972c4762a1bSJed Brown vc = u[j][i].v; 973c4762a1bSJed Brown 974c4762a1bSJed Brown /* uxx = (-2.0*uc + u[j][i-1].u + u[j][i+1].u)*sx; 975c4762a1bSJed Brown uyy = (-2.0*uc + u[j-1][i].u + u[j+1][i].u)*sy; 976c4762a1bSJed Brown 977c4762a1bSJed Brown vxx = (-2.0*vc + u[j][i-1].v + u[j][i+1].v)*sx; 978c4762a1bSJed Brown vyy = (-2.0*vc + u[j-1][i].v + u[j+1][i].v)*sy; 979c4762a1bSJed Brown f[j][i].u = appctx->D1*(uxx + uyy) - uc*vc*vc + appctx->gamma*(1.0 - uc);*/ 980c4762a1bSJed Brown 9819371c9d4SSatish Balay stencil[0].i = i; 9829371c9d4SSatish Balay stencil[0].c = 0; 9839371c9d4SSatish Balay entries[0] = appctx->D1 * sy; 9849371c9d4SSatish Balay stencil[1].i = i; 9859371c9d4SSatish Balay stencil[1].c = 0; 9869371c9d4SSatish Balay entries[1] = appctx->D1 * sy; 9879371c9d4SSatish Balay stencil[2].i = i - 1; 9889371c9d4SSatish Balay stencil[2].c = 0; 9899371c9d4SSatish Balay entries[2] = appctx->D1 * sx; 9909371c9d4SSatish Balay stencil[3].i = i + 1; 9919371c9d4SSatish Balay stencil[3].c = 0; 9929371c9d4SSatish Balay entries[3] = appctx->D1 * sx; 9939371c9d4SSatish Balay stencil[4].i = i; 9949371c9d4SSatish Balay stencil[4].c = 0; 9959371c9d4SSatish Balay entries[4] = -2.0 * appctx->D1 * (sx + sy) - vc * vc - appctx->gamma; 9969371c9d4SSatish Balay stencil[5].i = i; 9979371c9d4SSatish Balay stencil[5].c = 1; 9989371c9d4SSatish Balay entries[5] = -2.0 * uc * vc; 9999371c9d4SSatish Balay rowstencil.i = i; 10009371c9d4SSatish Balay rowstencil.c = 0; 1001c4762a1bSJed Brown 10029566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(A, 1, &rowstencil, 6, stencil, entries, INSERT_VALUES)); 1003c4762a1bSJed Brown 10049371c9d4SSatish Balay stencil[0].c = 1; 10059371c9d4SSatish Balay entries[0] = appctx->D2 * sy; 10069371c9d4SSatish Balay stencil[1].c = 1; 10079371c9d4SSatish Balay entries[1] = appctx->D2 * sy; 10089371c9d4SSatish Balay stencil[2].c = 1; 10099371c9d4SSatish Balay entries[2] = appctx->D2 * sx; 10109371c9d4SSatish Balay stencil[3].c = 1; 10119371c9d4SSatish Balay entries[3] = appctx->D2 * sx; 10129371c9d4SSatish Balay stencil[4].c = 1; 10139371c9d4SSatish Balay entries[4] = -2.0 * appctx->D2 * (sx + sy) + 2.0 * uc * vc - appctx->gamma - appctx->kappa; 10149371c9d4SSatish Balay stencil[5].c = 0; 10159371c9d4SSatish Balay entries[5] = vc * vc; 1016c4762a1bSJed Brown rowstencil.c = 1; 1017c4762a1bSJed Brown 10189566063dSJacob Faibussowitsch PetscCall(MatSetValuesStencil(A, 1, &rowstencil, 6, stencil, entries, INSERT_VALUES)); 1019c4762a1bSJed Brown /* f[j][i].v = appctx->D2*(vxx + vyy) + uc*vc*vc - (appctx->gamma + appctx->kappa)*vc; */ 1020c4762a1bSJed Brown } 1021c4762a1bSJed Brown } 1022c4762a1bSJed Brown 1023c4762a1bSJed Brown /* 1024c4762a1bSJed Brown Restore vectors 1025c4762a1bSJed Brown */ 10269566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(19.0 * xm * ym)); 10279566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayRead(da, localU, &u)); 10289566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localU)); 10299566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 10309566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 10319566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE)); 1032c4762a1bSJed Brown if (appctx->aijpc) { 10339566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 10349566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 10359566063dSJacob Faibussowitsch PetscCall(MatSetOption(B, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE)); 1036c4762a1bSJed Brown } 1037c4762a1bSJed Brown PetscFunctionReturn(0); 1038c4762a1bSJed Brown } 1039c4762a1bSJed Brown 10409371c9d4SSatish Balay PetscErrorCode RHSJacobianAdolc(TS ts, PetscReal t, Vec U, Mat A, Mat B, void *ctx) { 1041c4762a1bSJed Brown AppCtx *appctx = (AppCtx *)ctx; 1042c4762a1bSJed Brown DM da; 1043c4762a1bSJed Brown PetscScalar *u_vec; 1044c4762a1bSJed Brown Vec localU; 1045c4762a1bSJed Brown 1046c4762a1bSJed Brown PetscFunctionBegin; 10479566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &da)); 10489566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(da, &localU)); 1049c4762a1bSJed Brown 1050c4762a1bSJed Brown /* 1051c4762a1bSJed Brown Scatter ghost points to local vector,using the 2-step process 1052c4762a1bSJed Brown DMGlobalToLocalBegin(),DMGlobalToLocalEnd(). 1053c4762a1bSJed Brown By placing code between these two statements, computations can be 1054c4762a1bSJed Brown done while messages are in transition. 1055c4762a1bSJed Brown */ 10569566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, localU)); 10579566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, localU)); 1058c4762a1bSJed Brown 1059c4762a1bSJed Brown /* Get pointers to vector data */ 10609566063dSJacob Faibussowitsch PetscCall(VecGetArray(localU, &u_vec)); 1061c4762a1bSJed Brown 1062c4762a1bSJed Brown /* 1063c4762a1bSJed Brown Compute Jacobian 1064c4762a1bSJed Brown */ 10659566063dSJacob Faibussowitsch PetscCall(PetscAdolcComputeRHSJacobianLocal(1, A, u_vec, appctx->adctx)); 1066c4762a1bSJed Brown 1067c4762a1bSJed Brown /* 1068c4762a1bSJed Brown Restore vectors 1069c4762a1bSJed Brown */ 10709566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(localU, &u_vec)); 10719566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(da, &localU)); 1072c4762a1bSJed Brown PetscFunctionReturn(0); 1073c4762a1bSJed Brown } 1074c4762a1bSJed Brown 1075c4762a1bSJed Brown /*TEST 1076c4762a1bSJed Brown 1077c4762a1bSJed Brown build: 1078c4762a1bSJed Brown requires: double !complex adolc colpack 1079c4762a1bSJed Brown 1080c4762a1bSJed Brown test: 1081c4762a1bSJed Brown suffix: 1 1082c4762a1bSJed Brown nsize: 1 1083c4762a1bSJed Brown args: -ts_max_steps 1 -da_grid_x 12 -da_grid_y 12 -snes_test_jacobian 1084c4762a1bSJed Brown output_file: output/adr_ex5adj_1.out 1085c4762a1bSJed Brown 1086c4762a1bSJed Brown test: 1087c4762a1bSJed Brown suffix: 2 1088c4762a1bSJed Brown nsize: 1 1089c4762a1bSJed Brown args: -ts_max_steps 1 -da_grid_x 12 -da_grid_y 12 -snes_test_jacobian -implicitform 1090c4762a1bSJed Brown output_file: output/adr_ex5adj_2.out 1091c4762a1bSJed Brown 1092c4762a1bSJed Brown test: 1093c4762a1bSJed Brown suffix: 3 1094c4762a1bSJed Brown nsize: 4 1095c4762a1bSJed Brown args: -ts_max_steps 10 -da_grid_x 12 -da_grid_y 12 -ts_monitor -ts_adjoint_monitor 1096c4762a1bSJed Brown output_file: output/adr_ex5adj_3.out 1097c4762a1bSJed Brown 1098c4762a1bSJed Brown test: 1099c4762a1bSJed Brown suffix: 4 1100c4762a1bSJed Brown nsize: 4 1101c4762a1bSJed Brown args: -ts_max_steps 10 -da_grid_x 12 -da_grid_y 12 -ts_monitor -ts_adjoint_monitor -implicitform 1102c4762a1bSJed Brown output_file: output/adr_ex5adj_4.out 1103c4762a1bSJed Brown 1104c4762a1bSJed Brown testset: 1105c4762a1bSJed Brown suffix: 5 1106c4762a1bSJed Brown nsize: 4 1107c4762a1bSJed Brown args: -ts_max_steps 10 -da_grid_x 15 -da_grid_y 15 -ts_monitor -ts_adjoint_monitor -adolc_sparse 1108c4762a1bSJed Brown output_file: output/adr_ex5adj_5.out 1109c4762a1bSJed Brown 1110c4762a1bSJed Brown testset: 1111c4762a1bSJed Brown suffix: 6 1112c4762a1bSJed Brown nsize: 4 1113c4762a1bSJed Brown args: -ts_max_steps 10 -da_grid_x 15 -da_grid_y 15 -ts_monitor -ts_adjoint_monitor -adolc_sparse -implicitform 1114c4762a1bSJed Brown output_file: output/adr_ex5adj_6.out 1115c4762a1bSJed Brown 1116c4762a1bSJed Brown TEST*/ 1117