1 static char help[] = "Landau collision operator driver\n\n"; 2 3 #include <petscts.h> 4 #include <petsclandau.h> 5 6 int main(int argc, char **argv) 7 { 8 DM dm; 9 Vec X,X_0; 10 PetscErrorCode ierr; 11 PetscInt dim=2; 12 TS ts; 13 Mat J; 14 SNES snes; 15 KSP ksp; 16 PC pc; 17 SNESLineSearch linesearch; 18 PetscReal time; 19 20 ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 21 ierr = PetscOptionsGetInt(NULL,NULL, "-dim", &dim, NULL);CHKERRQ(ierr); 22 /* Create a mesh */ 23 ierr = DMPlexLandauCreateVelocitySpace(PETSC_COMM_SELF, dim, "", &X, &J, &dm);CHKERRQ(ierr); 24 ierr = DMSetUp(dm);CHKERRQ(ierr); 25 ierr = VecDuplicate(X,&X_0);CHKERRQ(ierr); 26 ierr = VecCopy(X,X_0);CHKERRQ(ierr); 27 ierr = DMPlexLandauPrintNorms(X,0);CHKERRQ(ierr); 28 ierr = DMSetOutputSequenceNumber(dm, 0, 0.0);CHKERRQ(ierr); 29 ierr = DMViewFromOptions(dm,NULL,"-dm_view");CHKERRQ(ierr); 30 ierr = VecViewFromOptions(X,NULL,"-vec_view");CHKERRQ(ierr); 31 /* Create timestepping solver context */ 32 ierr = TSCreate(PETSC_COMM_SELF,&ts);CHKERRQ(ierr); 33 ierr = TSSetDM(ts,dm);CHKERRQ(ierr); 34 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 35 ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr); 36 ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr); 37 ierr = TSSetIFunction(ts,NULL,DMPlexLandauIFunction,NULL);CHKERRQ(ierr); 38 ierr = TSSetIJacobian(ts,J,J,DMPlexLandauIJacobian,NULL);CHKERRQ(ierr); 39 ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr); 40 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 41 ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 42 ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 43 ierr = TSSetSolution(ts,X);CHKERRQ(ierr); 44 ierr = TSSolve(ts,X);CHKERRQ(ierr); 45 ierr = DMPlexLandauPrintNorms(X,1);CHKERRQ(ierr); 46 ierr = TSGetTime(ts, &time);CHKERRQ(ierr); 47 ierr = DMSetOutputSequenceNumber(dm, 1, time);CHKERRQ(ierr); 48 ierr = VecViewFromOptions(X,NULL,"-vec_view");CHKERRQ(ierr); 49 ierr = VecAXPY(X,-1,X_0);CHKERRQ(ierr); 50 /* clean up */ 51 ierr = DMPlexLandauDestroyVelocitySpace(&dm);CHKERRQ(ierr); 52 ierr = TSDestroy(&ts);CHKERRQ(ierr); 53 ierr = VecDestroy(&X);CHKERRQ(ierr); 54 ierr = VecDestroy(&X_0);CHKERRQ(ierr); 55 ierr = PetscFinalize(); 56 return ierr; 57 } 58 59 /*TEST 60 testset: 61 requires: p4est !complex double 62 output_file: output/ex1_0.out 63 args: -dm_landau_num_species_grid 1,2 -petscspace_degree 3 -petscspace_poly_tensor 1 -dm_landau_type p4est -dm_landau_ion_masses 2,4 -dm_landau_ion_charges 1,18 -dm_landau_thermal_temps 5,5,.5 -dm_landau_n 1.00018,1,1e-5 -dm_landau_n_0 1e20 -ts_monitor -snes_rtol 1.e-14 -snes_stol 1.e-14 -snes_monitor -snes_converged_reason -ts_type arkimex -ts_arkimex_type 1bee -ts_max_snes_failures -1 -ts_rtol 1e-1 -ts_dt 1.e-1 -ts_max_time 1 -ts_adapt_clip .5,1.25 -ts_adapt_scale_solve_failed 0.75 -ts_adapt_time_step_increase_delay 5 -ts_max_steps 1 -pc_type lu -ksp_type preonly -dm_landau_amr_levels_max 2,1 64 test: 65 suffix: cpu 66 args: -dm_landau_device_type cpu 67 test: 68 suffix: kokkos 69 requires: kokkos_kernels 70 args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos 71 test: 72 suffix: cuda 73 requires: cuda 74 args: -dm_landau_device_type cuda -dm_mat_type aijcusparse -dm_vec_type cuda -mat_cusparse_use_cpu_solve 75 76 TEST*/ 77