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 PetscInt dim=2; 11 TS ts; 12 Mat J; 13 SNES snes; 14 KSP ksp; 15 PC pc; 16 SNESLineSearch linesearch; 17 PetscReal time; 18 19 CHKERRQ(PetscInitialize(&argc, &argv, NULL,help)); 20 CHKERRQ(PetscOptionsGetInt(NULL,NULL, "-dim", &dim, NULL)); 21 /* Create a mesh */ 22 CHKERRQ(DMPlexLandauCreateVelocitySpace(PETSC_COMM_SELF, dim, "", &X, &J, &dm)); 23 CHKERRQ(DMSetUp(dm)); 24 CHKERRQ(VecDuplicate(X,&X_0)); 25 CHKERRQ(VecCopy(X,X_0)); 26 CHKERRQ(DMPlexLandauPrintNorms(X,0)); 27 CHKERRQ(DMSetOutputSequenceNumber(dm, 0, 0.0)); 28 CHKERRQ(DMViewFromOptions(dm,NULL,"-dm_view")); 29 CHKERRQ(VecViewFromOptions(X,NULL,"-vec_view")); 30 /* Create timestepping solver context */ 31 CHKERRQ(TSCreate(PETSC_COMM_SELF,&ts)); 32 CHKERRQ(TSSetDM(ts,dm)); 33 CHKERRQ(TSGetSNES(ts,&snes)); 34 CHKERRQ(SNESGetLineSearch(snes,&linesearch)); 35 CHKERRQ(SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC)); 36 CHKERRQ(TSSetIFunction(ts,NULL,DMPlexLandauIFunction,NULL)); 37 CHKERRQ(TSSetIJacobian(ts,J,J,DMPlexLandauIJacobian,NULL)); 38 CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 39 CHKERRQ(SNESGetKSP(snes,&ksp)); 40 CHKERRQ(KSPGetPC(ksp,&pc)); 41 CHKERRQ(TSSetFromOptions(ts)); 42 CHKERRQ(TSSetSolution(ts,X)); 43 CHKERRQ(TSSolve(ts,X)); 44 CHKERRQ(DMPlexLandauPrintNorms(X,1)); 45 CHKERRQ(TSGetTime(ts, &time)); 46 CHKERRQ(DMSetOutputSequenceNumber(dm, 1, time)); 47 CHKERRQ(VecViewFromOptions(X,NULL,"-vec_view")); 48 CHKERRQ(VecAXPY(X,-1,X_0)); 49 /* clean up */ 50 CHKERRQ(DMPlexLandauDestroyVelocitySpace(&dm)); 51 CHKERRQ(TSDestroy(&ts)); 52 CHKERRQ(VecDestroy(&X)); 53 CHKERRQ(VecDestroy(&X_0)); 54 CHKERRQ(PetscFinalize()); 55 return 0; 56 } 57 58 /*TEST 59 testset: 60 requires: p4est !complex double 61 output_file: output/ex1_0.out 62 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 63 test: 64 suffix: cpu 65 args: -dm_landau_device_type cpu 66 test: 67 suffix: kokkos 68 requires: kokkos_kernels 69 args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos 70 test: 71 suffix: cuda 72 requires: cuda 73 args: -dm_landau_device_type cuda -dm_mat_type aijcusparse -dm_vec_type cuda -mat_cusparse_use_cpu_solve 74 75 TEST*/ 76